csdsolvertrafofpsh2d.cpp 21 KB
Newer Older
1
#include "solvers/csdsolvertrafofpsh2d.h"
2
#include "problems/icru.h"
3
4
5
6
7
8
9

CSDSolverTrafoFPSH2D::CSDSolverTrafoFPSH2D( Config* settings ) : SNSolver( settings ) {
    _dose = std::vector<double>( _settings->GetNCells(), 0.0 );

    // Set angle and energies
    _energies  = Vector( _nEnergies, 0.0 );    // equidistant
    _energyMin = 1e-4 * 0.511;
10
    _energyMax = _settings->GetMaxEnergyCSD();    // 0.1;
11
12
13
14
15

    // write equidistant energy grid (false) or refined grid (true)
    GenerateEnergyGrid( false );

    // create quadrature
16
17
    unsigned order = _quadrature->GetOrder();
    // unsigned nq       = _settings->GetNQuadPoints();
18
19
20
21
22
23
24
25
26
27
28
    _quadPoints       = _quadrature->GetPoints();
    _weights          = _quadrature->GetWeights();
    _quadPointsSphere = _quadrature->GetPointsSphere();

    // transform structured quadrature
    _mu  = Vector( 2 * order );
    _phi = Vector( 2 * order );
    _wp  = Vector( 2 * order );
    _wa  = Vector( 2 * order );

    // create quadrature 1D to compute mu grid
29
    QuadratureBase* quad1D = QuadratureBase::Create( QUAD_GaussLegendre1D, 2 * order );
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
    Vector w               = quad1D->GetWeights();
    VectorVector muVec     = quad1D->GetPoints();
    for( unsigned k = 0; k < 2 * order; ++k ) {
        _mu[k] = muVec[k][0];
        _wp[k] = w[k];
    }

    for( unsigned i = 0; i < 2 * order; ++i ) {
        _phi[i] = ( i + 0.5 ) * M_PI / order;
        _wa[i]  = M_PI / order;
    }

    _FPMethod = 2;

    double DMinus = 0.0;
    double DPlus  = 0.0;

    Vector gamma( 2 * order, 0.0 );
    double dPlus;
    double c, K;

    double dMinus = 0.0;
    DPlus         = DMinus - 2 * _mu[0] * w[0];
    for( unsigned j = 0; j < 2 * order - 1; ++j ) {
        DMinus   = DPlus;
        DPlus    = DMinus - 2 * _mu[j] * w[j];
        dPlus    = ( sqrt( 1 - _mu[j + 1] * _mu[j + 1] ) - sqrt( 1 - _mu[j] * _mu[j] ) ) / ( _mu[j + 1] - _mu[j] );
        c        = ( DPlus * dPlus - DMinus * dMinus ) / _wp[j];
        K        = 2 * ( 1 - _mu[j] * _mu[j] ) + c * sqrt( 1 - _mu[j] * _mu[j] );
        gamma[j] = M_PI * M_PI * K / ( 2 * order * ( 1 - std::cos( M_PI / order ) ) );
        dMinus   = dPlus;
    }

    // Heney-Greenstein parameter
    double g = 0.8;

    // determine momente of Heney-Greenstein
    _xi1 = Vector( _nEnergies, 1.0 - g );    // paper Olbrant, Frank (11)
    _xi2 = Vector( _nEnergies, 4.0 / 3.0 - 2.0 * g + 2.0 / 3.0 * g * g );
    _xi  = Matrix( 4, _nEnergies );
    for( unsigned n = 0; n < _nEnergies; ++n ) {
        _xi( 1, n ) = 1.0 - g;
        _xi( 2, n ) = 4.0 / 3.0 - 2.0 * g + 2.0 / 3.0 * g * g;
    }

    // initialize stopping power vector
    _s = Vector( _nEnergies, 1.0 );

    _RT = true;

    // read in medical data if radiation therapy option selected
    if( _RT ) {
82
        ICRU database( _mu, _energies, _settings );
83
84
85
86
87
88
        database.GetTransportCoefficients( _xi );
        database.GetStoppingPower( _s );
    }

    _quadPointsSphere = _quadrature->GetPointsSphere();

89
90
91
92
    unsigned orderSph = order;
    SphericalHarmonics sph( orderSph );
    unsigned nSph = orderSph * orderSph + 2 * orderSph + 1;
    Matrix Y      = blaze::zero<double>( nSph, _nq );
93
94
95
96
97

    for( unsigned q = 0; q < _nq; q++ ) {
        blaze::column( Y, q ) = sph.ComputeSphericalBasis( _quadPointsSphere[q][0], _quadPointsSphere[q][1] );
    }

98
    //_O = Y;
99
    _O = blaze::trans( Y );
100
    _M = _O;    // transposed to simplify weight multiplication. We will transpose _M later again
101

102
    for( unsigned j = 0; j < nSph; ++j ) {
103
104
105
106
107
        blaze::column( _M, j ) *= _weights;
    }

    _M.transpose();

108
    _S               = Matrix( nSph, nSph, 0.0 );
109
    unsigned counter = 0;
110
    for( int l = 0; l <= (int)orderSph; ++l ) {
111
112
        for( int m = -l; m <= l; ++m ) {
            _S( counter, counter ) = double( -l * ( l + 1 ) );
113
114
115
116
            counter++;
        }
    }

117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
    //_density = std::vector<double>( _nCells, 1.0 );

    _L = _O * _S * _M;

    // check against FD
    // setup Laplace Beltrami matrix L in slab geometry
    Matrix LFD( _nq, _nq, 0.0 );

    _FPMethod = 2;

    DMinus = 0.0;
    DPlus  = 0.0;

    dMinus = 0.0;
    DPlus  = DMinus - 2 * _mu[0] * w[0];
    for( unsigned j = 0; j < order - 1; ++j ) {
        DMinus   = DPlus;
        DPlus    = DMinus - 2 * _mu[j] * w[j];
        dPlus    = ( sqrt( 1 - _mu[j + 1] * _mu[j + 1] ) - sqrt( 1 - _mu[j] * _mu[j] ) ) / ( _mu[j + 1] - _mu[j] );
        c        = ( DPlus * dPlus - DMinus * dMinus ) / _wp[j];
        K        = 2 * ( 1 - _mu[j] * _mu[j] ) + c * sqrt( 1 - _mu[j] * _mu[j] );
        gamma[j] = M_PI * M_PI * K / ( 2 * order * ( 1 - std::cos( M_PI / order ) ) );
        dMinus   = dPlus;
140
141
    }

142
    DPlus = 0.0;
143

144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
    // implementation of 2D spherical Laplacian according book "Advances in Discrete Ordinates Methodology", equation (1.136)
    for( unsigned j = 0; j < order; ++j ) {
        DMinus = DPlus;
        DPlus  = DMinus - 2 * _mu[j] * _wp[j];
        for( unsigned i = 0; i < 2 * order; ++i ) {
            if( j > 0 ) {
                LFD( j * 2 * order + i, ( j - 1 ) * 2 * order + i ) = DMinus / ( _mu[j] - _mu[j - 1] ) / _wp[j];
                LFD( j * 2 * order + i, j * 2 * order + i )         = -DMinus / ( _mu[j] - _mu[j - 1] ) / _wp[j];
            }
            if( i > 0 ) {
                LFD( j * 2 * order + i, j * 2 * order + i - 1 ) = 1.0 / ( 1 - _mu[j] * _mu[j] ) * gamma[j] / ( _phi[i] - _phi[i - 1] ) / _wa[i];
                LFD( j * 2 * order + i, j * 2 * order + i ) += -1.0 / ( 1 - _mu[j] * _mu[j] ) * gamma[j] / ( _phi[i] - _phi[i - 1] ) / _wa[i];
            }
            if( j < order - 1 ) {
                LFD( j * 2 * order + i, ( j + 1 ) * 2 * order + i ) = DPlus / ( _mu[j + 1] - _mu[j] ) / _wp[j];
                LFD( j * 2 * order + i, j * 2 * order + i ) += -DPlus / ( _mu[j + 1] - _mu[j] ) / _wp[j];
            }
            if( i < 2 * order - 1 ) {
                LFD( j * 2 * order + i, j * 2 * order + i + 1 ) = 1.0 / ( 1 - _mu[j] * _mu[j] ) * gamma[j] / ( _phi[i + 1] - _phi[i] ) / _wa[i];
                LFD( j * 2 * order + i, j * 2 * order + i ) += -1.0 / ( 1 - _mu[j] * _mu[j] ) * gamma[j] / ( _phi[i + 1] - _phi[i] ) / _wa[i];
            }
        }
    }

    // std::cout << _L << std::endl << std::endl;
    // std::cout << LFD << std::endl;
    // std::cout << _S << std::endl;
}

// void CSDSolverTrafoFPSH2D::Solve() {
//
//    // Prepare Solver output
//    PrepareVolumeOutput();
//
//    SolverPreprocessing();
//
//    for( unsigned n = 0; n < _nEnergies - 1; ++n ) {
//        IterPreprocessing( n );
//
//        FluxUpdate();
//
//        FVMUpdate( n );
//
//        IterPostprocessing();
//
//        WriteVolumeOutput( n );
//        WriteScalarOutput( n );
//        PrintScreenOutput( n );
//        PrintHistoryOutput( n );
//        PrintVolumeOutput( n );
//
//        //#pragma omp parallel for
//
//        //       for( unsigned j = 0; j < _nCells; ++j ) {
//        //           fluxNew[j] = dot( _sol[j], _weights );
//        //           if( n > 0 ) {
//        //               _dose[j] += 0.5 * _dE * ( fluxNew[j] * _s[_nEnergies - n - 1] + fluxOld[j] * _s[_nEnergies - n] ) /
//        //                           _density[j];    // update dose with trapezoidal rule
//        //           }
//        //           else {
//        //               _dose[j] += _dE * fluxNew[j] * _s[_nEnergies - n - 1] / _density[j];
//        //           }
//        //           _solverOutput[j] = fluxNew[j];
//        //       }
//        //
//        //       // Save( n );
//        //       std::cout << _energiesOrig[_nEnergies - n - 1] << std::endl;
//        //       dFlux   = blaze::l2Norm( fluxNew - fluxOld );
//        //       fluxOld = fluxNew;
//        //       if( rank == 0 )
//        //           log->info( "{:03.8f}  {:03.8f}  {:01.5e}  {:01.5e}", _energiesOrig[_nEnergies - n - 1], _energies[n], _dE / _densityMin, dFlux
//        );
//        //       if( std::isinf( dFlux ) || std::isnan( dFlux ) ) break;
//        //   }
//        //   for( unsigned j = 0; j < _nCells; ++j ) _solverOutput[j] = _density[j];
//        //   Save( 1 );
//        //   Save();
//    }
//}

// void CSDSolverTrafoFPSH2D::Save() const {
//    std::vector<std::vector<std::vector<double>>> outputFields;
//    std::vector<std::vector<std::string>> outputFieldNames;
//
//    // one group
//    outputFields.resize( 1 );
//    outputFieldNames.resize( 1 );
//
//    // 2 fields for this group
//    outputFields[0].resize( 1 );
//    outputFieldNames[0].resize( 1 );
//
//    // Fill names
//    outputFieldNames[0][0] = "dose";
//    // outputFieldNames[0][1] = "normalized dose";
//
//    // Fill fields
//    outputFields[0][0].resize( _nCells );
//    // outputFields[0][1].resize( _nCells );
//
//    // std::vector<double> normalizedDose = _dose;
//    // double maxDose                     = *std::max_element( _dose.begin(), _dose.end() );
//    // for( unsigned i = 0; i < _dose.size(); ++i ) normalizedDose[i] /= maxDose;
//
//    for( unsigned idx_cell = 0; idx_cell < _nCells; idx_cell++ ) {
//        outputFields[0][0][idx_cell] = _dose[idx_cell];
//        // outputFields[0][1][idx_cell] = normalizedDose[idx_cell];
//    }
//
//    // debug
//    // for( unsigned idx_cell = 0; idx_cell < _nCells; idx_cell++ ) {
//    //    std::cout << _dose[idx_cell] << "\n";
//    //}
//    std::cout << _settings->GetOutputFile();
//    ExportVTK( _settings->GetOutputFile() + "_TEST", outputFields, outputFieldNames, _mesh );
//}
//
// void CSDSolverTrafoFPSH2D::Save( int currEnergy ) const {
//    std::vector<std::string> fieldNames{ "flux" };
//    std::vector<std::vector<std::string>> fieldNamesWrapper{ fieldNames };
//
//    std::vector<std::vector<double>> scalarField( 1, _solverOutput );
//    std::vector<std::vector<std::vector<double>>> results{ scalarField };
//    ExportVTK( _settings->GetOutputFile() + "_" + std::to_string( currEnergy ), results, fieldNamesWrapper, _mesh );
//}

void CSDSolverTrafoFPSH2D::GenerateEnergyGrid( bool refinement ) {
    _dE = ComputeTimeStep( _settings->GetCFL() );
    if( !refinement ) {
        _nEnergies = unsigned( ( _energyMax - _energyMin ) / _dE );
        _energies.resize( _nEnergies );
        for( unsigned n = 0; n < _nEnergies; ++n ) {
            _energies[n] = _energyMin + ( _energyMax - _energyMin ) / ( _nEnergies - 1 ) * n;
        }
    }
    else {
        // hard-coded positions for energy-grid refinement
        double energySwitch    = 0.58;
        double energySwitchMin = 0.03;

        // number of energies per intervals [E_min,energySwitchMin], [energySwitchMin,energySwitch], [energySwitch,E_Max]
        unsigned nEnergies1 = unsigned( ( _energyMax - energySwitch ) / _dE );
        unsigned nEnergies2 = unsigned( ( energySwitch - energySwitchMin ) / ( _dE / 2 ) );
        unsigned nEnergies3 = unsigned( ( energySwitchMin - _energyMin ) / ( _dE / 3 ) );
        _nEnergies          = nEnergies1 + nEnergies2 + nEnergies3 - 2;
        _energies.resize( _nEnergies );

        // write equidistant energy grid in each interval
        for( unsigned n = 0; n < nEnergies3; ++n ) {
            _energies[n] = _energyMin + ( energySwitchMin - _energyMin ) / ( nEnergies3 - 1 ) * n;
        }
        for( unsigned n = 1; n < nEnergies2; ++n ) {
            _energies[n + nEnergies3 - 1] = energySwitchMin + ( energySwitch - energySwitchMin ) / ( nEnergies2 - 1 ) * n;
        }
        for( unsigned n = 1; n < nEnergies1; ++n ) {
            _energies[n + nEnergies3 + nEnergies2 - 2] = energySwitch + ( _energyMax - energySwitch ) / ( nEnergies1 - 1 ) * n;
        }
    }
}

// IO
void CSDSolverTrafoFPSH2D::PrepareVolumeOutput() {
    unsigned nGroups = (unsigned)_settings->GetNVolumeOutput();

    _outputFieldNames.resize( nGroups );
    _outputFields.resize( nGroups );

    // Prepare all OutputGroups ==> Specified in option VOLUME_OUTPUT
    for( unsigned idx_group = 0; idx_group < nGroups; idx_group++ ) {
        // Prepare all Output Fields per group

        // Different procedure, depending on the Group...
        switch( _settings->GetVolumeOutput()[idx_group] ) {
            case MINIMAL:
                // Currently only one entry ==> rad flux
                _outputFields[idx_group].resize( 1 );
                _outputFieldNames[idx_group].resize( 1 );

                _outputFields[idx_group][0].resize( _nCells );
                _outputFieldNames[idx_group][0] = "radiation flux density";
                break;

            case MEDICAL:
                _outputFields[idx_group].resize( 2 );
                _outputFieldNames[idx_group].resize( 2 );

                // Dose
                _outputFields[idx_group][0].resize( _nCells );
                _outputFieldNames[idx_group][0] = "dose";
                // Normalized Dose
                _outputFields[idx_group][1].resize( _nCells );
                _outputFieldNames[idx_group][1] = "normalized dose";
                break;

            default: ErrorMessages::Error( "Volume Output Group not defined for CSD_SN_FP_TRAFO Solver!", CURRENT_FUNCTION ); break;
        }
    }
}

void CSDSolverTrafoFPSH2D::WriteVolumeOutput( unsigned idx_pseudoTime ) {
    unsigned nGroups = (unsigned)_settings->GetNVolumeOutput();
    double maxDose;
    if( ( _settings->GetVolumeOutputFrequency() != 0 && idx_pseudoTime % (unsigned)_settings->GetVolumeOutputFrequency() == 0 ) ||
        ( idx_pseudoTime == _maxIter - 1 ) /* need sol at last iteration */ ) {

        for( unsigned idx_group = 0; idx_group < nGroups; idx_group++ ) {
            switch( _settings->GetVolumeOutput()[idx_group] ) {
                case MINIMAL:
                    for( unsigned idx_cell = 0; idx_cell < _nCells; ++idx_cell ) {
                        _outputFields[idx_group][0][idx_cell] = _fluxNew[idx_cell];
                    }
                    break;

                case MEDICAL:
                    // Compute Dose
                    for( unsigned idx_cell = 0; idx_cell < _nCells; ++idx_cell ) {
jonas.kusch's avatar
jonas.kusch committed
360
                        _outputFields[idx_group][0][idx_cell] += _dose[idx_cell];
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
                    }
                    // Compute normalized dose
                    _outputFields[idx_group][1] = _outputFields[idx_group][0];

                    maxDose = *std::max_element( _outputFields[idx_group][0].begin(), _outputFields[idx_group][0].end() );

                    for( unsigned idx_cell = 0; idx_cell < _nCells; ++idx_cell ) {
                        _outputFields[idx_group][1][idx_cell] /= maxDose;
                    }
                    break;

                default: ErrorMessages::Error( "Volume Output Group not defined for CSD_SN_FP_TRAFO Solver!", CURRENT_FUNCTION ); break;
            }
        }
    }
}

// Solver
379
void CSDSolverTrafoFPSH2D::FVMUpdate( unsigned /*idx_energy*/ ) {
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
// loop over all spatial cells
#pragma omp parallel for
    for( unsigned j = 0; j < _nCells; ++j ) {
        if( _boundaryCells[j] == BOUNDARY_TYPE::DIRICHLET ) continue;
        // loop over all ordinates
        for( unsigned i = 0; i < _nq; ++i ) {
            // time update angular flux with numerical flux and total scattering cross section
            _solNew[j][i] = _sol[j][i] - _dE * _solNew[j][i];
        }
    }
}

void CSDSolverTrafoFPSH2D::FluxUpdate() {
#pragma omp parallel for
    for( unsigned j = 0; j < _nCells; ++j ) {
        if( _boundaryCells[j] == BOUNDARY_TYPE::DIRICHLET ) continue;
        // loop over all ordinates
        for( unsigned i = 0; i < _nq; ++i ) {
            _solNew[j][i] = 0.0;
            // loop over all neighbor cells (edges) of cell j and compute numerical fluxes
            for( unsigned idx_neighbor = 0; idx_neighbor < _neighbors[j].size(); ++idx_neighbor ) {
                // store flux contribution on psiNew_sigmaS to save memory
                if( _boundaryCells[j] == BOUNDARY_TYPE::NEUMANN && _neighbors[j][idx_neighbor] == _nCells )
                    continue;    // adiabatic wall, add nothing
                else
                    _solNew[j][i] += _g->Flux( _quadPoints[i],
                                               _sol[j][i] / _density[j],
                                               _sol[_neighbors[j][idx_neighbor]][i] / _density[_neighbors[j][idx_neighbor]],
                                               _normals[j][idx_neighbor] ) /
                                     _areas[j];
            }
        }
    }
}

void CSDSolverTrafoFPSH2D::IterPreprocessing( unsigned idx_pseudotime ) {
    unsigned n = idx_pseudotime;

    _dE = fabs( _energies[n + 1] - _energies[n] );    // is the sign correct here?

    double xi1 = _xi( 1, _nEnergies - n - 1 );
    double xi2 = _xi( 2, _nEnergies - n - 1 );
    double xi3 = _xi( 3, _nEnergies - n - 1 );

    // setup coefficients in FP step
    if( _FPMethod == 1 ) {
        _alpha  = 0.0;
        _alpha2 = xi1 / 2.0;
        _beta   = 0.0;
    }
    else if( _FPMethod == 2 ) {
        _alpha  = xi1 / 2.0 + xi2 / 8.0;
        _alpha2 = 0.0;
        _beta   = xi2 / 8.0 / xi1;
    }
    else if( _FPMethod == 3 ) {
        _alpha  = xi2 * ( 27.0 * xi2 * xi2 + 5.0 * xi3 * xi3 - 24.0 * xi2 * xi3 ) / ( 8.0 * xi3 * ( 3.0 * xi2 - 2.0 * xi3 ) );
        _beta   = xi3 / ( 6.0 * ( 3.0 * xi2 - 2.0 * xi3 ) );
        _alpha2 = xi1 / 2.0 - 9.0 / 8.0 * xi2 * xi2 / xi3 + 3.0 / 8.0 * xi2;
    }

    _IL = _identity - _beta * _L;

    // write BC for water phantom
    if( _RT && false ) {
        for( unsigned k = 0; k < _nq; ++k ) {
            if( _quadPoints[k][0] > 0 ) {
                _sol[0][k] = 1e5 * exp( -200.0 * pow( 1.0 - _quadPoints[k][0], 2 ) ) *
                             exp( -50.0 * pow( _energyMax - _energiesOrig[_nEnergies - n - 1], 2 ) ) * _density[0] * _s[_nEnergies - n - 1];
            }
        }
    }

// add FP scattering term implicitly
#pragma omp parallel for
    for( unsigned j = 0; j < _nCells; ++j ) {
        if( _boundaryCells[j] == BOUNDARY_TYPE::DIRICHLET ) continue;
        //_sol[j] = blaze::solve( identity - _dE * _alpha2 * _L, psiNew[j] );
        _sol[j] = _IL * blaze::solve( _IL - _dE * _alpha * _L, _sol[j] );
    }
}

jonas.kusch's avatar
jonas.kusch committed
462
463
void CSDSolverTrafoFPSH2D::IterPostprocessing( unsigned idx_pseudotime ) {
    unsigned n = idx_pseudotime;
464
465
466
467
468
469
470
    for( unsigned j = 0; j < _nCells; ++j ) {
        if( _boundaryCells[j] == BOUNDARY_TYPE::DIRICHLET ) continue;
        _sol[j] = _solNew[j];
    }

    // --- Compute Flux for solution and Screen Output ---
    ComputeRadFlux();
jonas.kusch's avatar
jonas.kusch committed
471
472
473
474
475
476
477
478
479
480
481
482

    for( unsigned j = 0; j < _nCells; ++j ) {
        _fluxNew[j] = dot( _sol[j], _weights );
        if( n > 0 ) {
            _dose[j] += 0.5 * _dE * ( _fluxNew[j] * _s[_nEnergies - n - 1] + _flux[j] * _s[_nEnergies - n] ) /
                        _density[j];    // update dose with trapezoidal rule
        }
        else {
            _dose[j] += _dE * _fluxNew[j] * _s[_nEnergies - n - 1] / _density[j];
        }
        _solverOutput[j] = _fluxNew[j];
    }
483
484
}

485
void CSDSolverTrafoFPSH2D::SolverPreprocessing() {
486
487
    auto log = spdlog::get( "event" );

488
489
490
491
492
    _densityMin = 0.1;
    for( unsigned j = 0; j < _nCells; ++j ) {
        if( _density[j] < _densityMin ) _density[j] = _densityMin;
    }

493
    // save original energy field for boundary conditions
494
    _energiesOrig = _energies;
495
496

    // setup incoming BC on left
497
498
499
500
    //_sol = VectorVector( _density.size(), Vector( _settings->GetNQuadPoints(), 0.0 ) );    // hard coded IC, needs to be changed
    // for( unsigned k = 0; k < _nq; ++k ) {
    //    if( _quadPoints[k][0] > 0 && !_RT ) _sol[0][k] = 1e5 * exp( -10.0 * pow( 1.0 - _quadPoints[k][0], 2 ) );
    //}
501
    // hard coded boundary type for 1D testcases (otherwise cells will be NEUMANN)
502
503
    //_boundaryCells[0]           = BOUNDARY_TYPE::DIRICHLET;
    //_boundaryCells[_nCells - 1] = BOUNDARY_TYPE::DIRICHLET;
504
505

    // setup identity matrix for FP scattering
506
507
    _identity = Matrix( _nq, _nq, 0.0 );
    for( unsigned k = 0; k < _nq; ++k ) _identity( k, k ) = 1.0;
508
509
510
511
512

    int rank;
    MPI_Comm_rank( MPI_COMM_WORLD, &rank );
    if( rank == 0 ) log->info( "{:10}   {:10}", "E", "dFlux" );

513
514
// do substitution from psi to psiTildeHat (cf. Dissertation Kerstion Kuepper, Eq. 1.23)
#pragma omp parallel for
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
    for( unsigned j = 0; j < _nCells; ++j ) {
        for( unsigned k = 0; k < _nq; ++k ) {
            _sol[j][k] = _sol[j][k] * _density[j] * _s[_nEnergies - 1];    // note that _s[_nEnergies - 1] is stopping power at highest energy
        }
    }

    // store transformed energies ETilde instead of E in _energies vector (cf. Dissertation Kerstion Kuepper, Eq. 1.25)
    double tmp   = 0.0;
    _energies[0] = 0.0;
    for( unsigned n = 1; n < _nEnergies; ++n ) {
        tmp          = tmp + _dE * 0.5 * ( 1.0 / _s[n] + 1.0 / _s[n - 1] );
        _energies[n] = tmp;
    }

    // store transformed energies ETildeTilde instead of ETilde in _energies vector (cf. Dissertation Kerstion Kuepper, Eq. 1.25)
    for( unsigned n = 0; n < _nEnergies; ++n ) {
        _energies[n] = _energies[_nEnergies - 1] - _energies[n];
    }

    // determine minimal density for CFL computation
535
    _densityMin = _density[0];
536
    for( unsigned j = 1; j < _nCells; ++j ) {
537
        if( _densityMin > _density[j] ) _densityMin = _density[j];
538
539
540
541
    }

    // cross sections do not need to be transformed to ETilde energy grid since e.g. TildeSigmaT(ETilde) = SigmaT(E(ETilde))

542
    // std::cout << "nEnergies = " << _nEnergies << std::endl;
543
544
    // loop over energies (pseudo-time)
}