mnsolver.cpp 12.5 KB
Newer Older
1
#include "solvers/mnsolver.h"
2
3
#include "common/config.h"
#include "common/io.h"
4
#include "common/mesh.h"
5
6
7
#include "entropies/entropybase.h"
#include "fluxes/numericalflux.h"
#include "optimizers/optimizerbase.h"
8
#include "problems/problembase.h"
9
10
#include "quadratures/quadraturebase.h"
#include "solvers/sphericalharmonics.h"
11
#include "toolboxes/errormessages.h"
12
#include "toolboxes/textprocessingtoolbox.h"
13

14
// externals
15
#include "spdlog/spdlog.h"
16
#include <mpi.h>
steffen.schotthoefer's avatar
steffen.schotthoefer committed
17

18
#include <fstream>
19
20
//#include <chrono>

21
MNSolver::MNSolver( Config* settings ) : Solver( settings ) {
steffen.schotthoefer's avatar
steffen.schotthoefer committed
22

23
    // Is this good (fast) code using a constructor list?
24
25
    _LMaxDegree    = settings->GetMaxMomentDegree();
    _nTotalEntries = GlobalIndex( _LMaxDegree, int( _LMaxDegree ) ) + 1;
26
27
28
29
30
31
32

    // build quadrature object and store quadrature points and weights
    _quadPoints       = _quadrature->GetPoints();
    _weights          = _quadrature->GetWeights();
    _nq               = _quadrature->GetNq();
    _quadPointsSphere = _quadrature->GetPointsSphere();
    _settings->SetNQuadPoints( _nq );
33

steffen.schotthoefer's avatar
steffen.schotthoefer committed
34
35
36
37
38
39
    // transform sigmaT and sigmaS in sigmaA.
    _sigmaA = VectorVector( _nEnergies, Vector( _nCells, 0 ) );    // Get rid of this extra vektor!

    for( unsigned n = 0; n < _nEnergies; n++ ) {
        for( unsigned j = 0; j < _nCells; j++ ) {
            _sigmaA[n][j] = 0;    //_sigmaT[n][j] - _sigmaS[n][j];
40
            _sigmaS[n][j] = 0;
steffen.schotthoefer's avatar
steffen.schotthoefer committed
41
42
43
44
45
46
47
        }
    }

    // Initialize Scatter Matrix
    _scatterMatDiag    = Vector( _nTotalEntries, 1.0 );
    _scatterMatDiag[0] = 0.0;    // First entry is zero by construction.

48
49
50
51
52
    // Initialize Entropy
    _entropy = EntropyBase::Create( _settings );

    // Initialize Optimizer
    _optimizer = OptimizerBase::Create( _settings );
steffen.schotthoefer's avatar
steffen.schotthoefer committed
53

54
55
    // Initialize lagrange Multiplier
    _alpha = VectorVector( _nCells, Vector( _nTotalEntries, 0.0 ) );
56
57

    // Initialize and Pre-Compute Moments at quadrature points
58
    _basis = new SphericalHarmonics( _LMaxDegree );
59

60
    _moments = VectorVector( _nq, Vector( _nTotalEntries, 0.0 ) );
61
    ComputeMoments();
62
63

    // Solver output
64
    PrepareOutputFields();
steffen.schotthoefer's avatar
steffen.schotthoefer committed
65
66
}

67
68
69
MNSolver::~MNSolver() {
    delete _entropy;
    delete _optimizer;
70
    delete _basis;
71
}
72

steffen.schotthoefer's avatar
steffen.schotthoefer committed
73
74
75
76
77
78
int MNSolver::GlobalIndex( int l, int k ) const {
    int numIndicesPrevLevel  = l * l;    // number of previous indices untill level l-1
    int prevIndicesThisLevel = k + l;    // number of previous indices in current level
    return numIndicesPrevLevel + prevIndicesThisLevel;
}

79
void MNSolver::ComputeMoments() {
steffen.schotthoefer's avatar
steffen.schotthoefer committed
80
    double my, phi;
81
82

    for( unsigned idx_quad = 0; idx_quad < _nq; idx_quad++ ) {
83
84
        my  = _quadPointsSphere[idx_quad][0];
        phi = _quadPointsSphere[idx_quad][1];
85

86
        _moments[idx_quad] = _basis->ComputeSphericalBasis( my, phi );
87
    }
88
89
90
91
92
}

Vector MNSolver::ConstructFlux( unsigned idx_cell ) {

    // ---- Integration of Moment of flux ----
93
    double entropyL, entropyR, entropyFlux;
94

95
    Vector flux( _nTotalEntries, 0.0 );
96

97
    for( unsigned idx_quad = 0; idx_quad < _nq; idx_quad++ ) {
98
99
100

        entropyFlux = 0.0;    // Reset temorary flux

101
102
        entropyL = _entropy->EntropyPrimeDual( blaze::dot( _alpha[idx_cell], _moments[idx_quad] ) );

103
104
105
106
107
        for( unsigned idx_neigh = 0; idx_neigh < _neighbors[idx_cell].size(); idx_neigh++ ) {
            // Store fluxes in psiNew, to save memory
            if( _boundaryCells[idx_cell] == BOUNDARY_TYPE::NEUMANN && _neighbors[idx_cell][idx_neigh] == _nCells )
                entropyR = entropyL;
            else {
108
                entropyR = _entropy->EntropyPrimeDual( blaze::dot( _alpha[_neighbors[idx_cell][idx_neigh]], _moments[idx_quad] ) );
109
            }
110
            entropyFlux += _g->Flux( _quadPoints[idx_quad], entropyL, entropyR, _normals[idx_cell][idx_neigh] );
111
        }
112
        flux += _moments[idx_quad] * ( _weights[idx_quad] * entropyFlux );
113
114

        // ------- Relizablity Reconstruction Step ----
115
    }
116
    return flux;
117
118
}

119
120
void MNSolver::ComputeRealizableSolution( unsigned idx_cell ) {
    double entropyReconstruction = 0.0;
121
    _sol[idx_cell]               = 0;
122
    for( unsigned idx_quad = 0; idx_quad < _nq; idx_quad++ ) {
123
        // Make entropyReconstruction a member vector, s.t. it does not have to be re-evaluated in ConstructFlux
124
        entropyReconstruction = _entropy->EntropyPrimeDual( blaze::dot( _alpha[idx_cell], _moments[idx_quad] ) );
125
        _sol[idx_cell] += _moments[idx_quad] * ( _weights[idx_quad] * entropyReconstruction );
126
127
128
    }
}

steffen.schotthoefer's avatar
steffen.schotthoefer committed
129
130
131
132
133
134
135
136
void MNSolver::Solve() {

    int rank;
    MPI_Comm_rank( MPI_COMM_WORLD, &rank );

    auto log = spdlog::get( "event" );

    // angular flux at next time step (maybe store angular flux at all time steps, since time becomes energy?)
137
    VectorVector psiNew = _sol;
steffen.schotthoefer's avatar
steffen.schotthoefer committed
138

139
    double mass = 0;
140
141
142

    Save( -1 );

143
144
145
    if( rank == 0 ) log->info( "{:10}  {:10}", "t", "mass" );

    // if( rank == 0 ) log->info( "{:03.8f}   {:01.5e} {:01.5e}", -1.0, dFlux, mass1 ); Should be deleted
steffen.schotthoefer's avatar
steffen.schotthoefer committed
146
147

    // Loop over energies (pseudo-time of continuous slowing down approach)
148
    for( unsigned idx_energy = 0; idx_energy < _nEnergies; idx_energy++ ) {
steffen.schotthoefer's avatar
steffen.schotthoefer committed
149

150
        // ------- Reconstruction Step ----------------
151

152
        _optimizer->SolveMultiCell( _alpha, _sol, _moments );
153

154
155
        // Loop over the grid cells
        for( unsigned idx_cell = 0; idx_cell < _nCells; idx_cell++ ) {
156

157
158
            // ------- Relizablity Reconstruction Step ----

159
            ComputeRealizableSolution( idx_cell );
160
161

            // ------- Flux Computation Step --------------
162

163
164
            // Dirichlet Boundaries are finished now
            if( _boundaryCells[idx_cell] == BOUNDARY_TYPE::DIRICHLET ) continue;
steffen.schotthoefer's avatar
steffen.schotthoefer committed
165

166
            psiNew[idx_cell] = ConstructFlux( idx_cell );
167

168
            // ------ Finite Volume Update Step ------
169

170
171
172
            // NEED TO VECTORIZE
            for( unsigned idx_system = 0; idx_system < _nTotalEntries; idx_system++ ) {

173
                psiNew[idx_cell][idx_system] = _sol[idx_cell][idx_system] -
174
                                               ( _dE / _areas[idx_cell] ) * psiNew[idx_cell][idx_system] /* cell averaged flux */
175
                                               - _dE * _sol[idx_cell][idx_system] *
176
177
                                                     ( _sigmaA[idx_energy][idx_cell]                                    /* absorbtion influence */
                                                       + _sigmaS[idx_energy][idx_cell] * _scatterMatDiag[idx_system] ); /* scattering influence */
178
            }
steffen.schotthoefer's avatar
steffen.schotthoefer committed
179
        }
180

181
182
        // Update Solution
        _sol = psiNew;
183

184
        // --- VTK and CSV Output ---
185
        mass = WriteOutputFields( idx_energy );
186
        Save( idx_energy );
187
        WriteNNTrainingData( idx_energy );
188

189
        // --- Screen Output ---
190
        if( rank == 0 ) log->info( "{:03.8f}   {:01.5e}", _energies[idx_energy], mass );
steffen.schotthoefer's avatar
steffen.schotthoefer committed
191
192
    }
}
193

194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
void MNSolver::PrepareOutputFields() {
    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 );

211
212
                _outputFields[idx_group][0].resize( _nCells );
                _outputFieldNames[idx_group][0] = "radiation flux density";
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
                break;

            case MOMENTS:
                // As many entries as there are moments in the system
                _outputFields[idx_group].resize( _nTotalEntries );
                _outputFieldNames[idx_group].resize( _nTotalEntries );

                for( int idx_l = 0; idx_l <= (int)_LMaxDegree; idx_l++ ) {
                    for( int idx_k = -idx_l; idx_k <= idx_l; idx_k++ ) {
                        _outputFields[idx_group][GlobalIndex( idx_l, idx_k )].resize( _nCells );

                        _outputFieldNames[idx_group][GlobalIndex( idx_l, idx_k )] =
                            std::string( "u_" + std::to_string( idx_l ) + "^" + std::to_string( idx_k ) );
                    }
                }
                break;
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252

            case DUAL_MOMENTS:
                // As many entries as there are moments in the system
                _outputFields[idx_group].resize( _nTotalEntries );
                _outputFieldNames[idx_group].resize( _nTotalEntries );

                for( int idx_l = 0; idx_l <= (int)_LMaxDegree; idx_l++ ) {
                    for( int idx_k = -idx_l; idx_k <= idx_l; idx_k++ ) {
                        _outputFields[idx_group][GlobalIndex( idx_l, idx_k )].resize( _nCells );

                        _outputFieldNames[idx_group][GlobalIndex( idx_l, idx_k )] =
                            std::string( "alpha_" + std::to_string( idx_l ) + "^" + std::to_string( idx_k ) );
                    }
                }
                break;

            case ANALYTIC:
                // one entry per cell
                _outputFields[idx_group].resize( 1 );
                _outputFieldNames[idx_group].resize( 1 );
                _outputFields[idx_group][0].resize( _nCells );
                _outputFieldNames[idx_group][0] = std::string( "analytic radiation flux density" );
                break;

253
254
255
256
257
            default: ErrorMessages::Error( "Volume Output Group not defined for MN Solver!", CURRENT_FUNCTION ); break;
        }
    }
}

258
double MNSolver::WriteOutputFields( unsigned idx_pseudoTime ) {
259
260
    double mass      = 0.0;
    unsigned nGroups = (unsigned)_settings->GetNVolumeOutput();
261

262
263
264
    // Compute total "mass" of the system ==> to check conservation properties
    for( unsigned idx_cell = 0; idx_cell < _nCells; ++idx_cell ) {
        mass += _sol[idx_cell][0] * _areas[idx_cell];    // Should probably go to postprocessing
265
266
    }

267
268
269
270
271
272
273
274
275
276
277
278
279
280
    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] = _sol[idx_cell][0];
                }
                break;
            case MOMENTS:
                for( unsigned idx_sys = 0; idx_sys < _nTotalEntries; idx_sys++ ) {
                    for( unsigned idx_cell = 0; idx_cell < _nCells; ++idx_cell ) {
                        _outputFields[idx_group][idx_sys][idx_cell] = _sol[idx_cell][idx_sys];
                    }
                }
                break;
281
282
283
284
285
286
287
288
289
290
291
292
            case DUAL_MOMENTS:
                for( unsigned idx_sys = 0; idx_sys < _nTotalEntries; idx_sys++ ) {
                    for( unsigned idx_cell = 0; idx_cell < _nCells; ++idx_cell ) {
                        _outputFields[idx_group][idx_sys][idx_cell] = _alpha[idx_cell][idx_sys];
                    }
                }
                break;
            case ANALYTIC:
                // Compute total "mass" of the system ==> to check conservation properties
                for( unsigned idx_cell = 0; idx_cell < _nCells; ++idx_cell ) {

                    double time  = idx_pseudoTime * _dE;
293
                    double sigma = 1;
294
295
296
297
298
299

                    _outputFields[idx_group][0][idx_cell] = _problem->GetAnalyticalSolution(
                        _mesh->GetCellMidPoints()[idx_cell][0], _mesh->GetCellMidPoints()[idx_cell][1], time, sigma );
                }
                break;

300
301
302
            default: ErrorMessages::Error( "Volume Output Group not defined for MN Solver!", CURRENT_FUNCTION ); break;
        }
    }
303
304
    return mass;
}
305

306
void MNSolver::Save() const { ExportVTK( _settings->GetOutputFile(), _outputFields, _outputFieldNames, _mesh ); }
307

308
void MNSolver::Save( int currEnergy ) const {
309
    ExportVTK( _settings->GetOutputFile() + "_" + std::to_string( currEnergy ), _outputFields, _outputFieldNames, _mesh );
310
}
311
312
313
314
315
316
317

void MNSolver::WriteNNTrainingData( unsigned idx_pseudoTime ) {
    std::string filename = "trainNN.csv";
    std::ofstream myfile;
    myfile.open( filename, std::ofstream::app );

    for( unsigned idx_cell = 0; idx_cell < _nCells; idx_cell++ ) {
318
        myfile << 0 << ", " << _nTotalEntries << "," << idx_pseudoTime;
319
        for( unsigned idx_sys = 0; idx_sys < _nTotalEntries; idx_sys++ ) {
320
            myfile << "," << _sol[idx_cell][idx_sys];
321
        }
322
        myfile << " \n" << 1 << ", " << _nTotalEntries << "," << idx_pseudoTime;
323
        for( unsigned idx_sys = 0; idx_sys < _nTotalEntries; idx_sys++ ) {
324
            myfile << "," << _alpha[idx_cell][idx_sys];
325
        }
326
        myfile << "\n";
327
328
329
    }
    myfile.close();
}