mnsolver.cpp 8 KB
Newer Older
1
#include "solvers/mnsolver.h"
2
3
4
5
6
7
8
#include "entropies/entropybase.h"
#include "fluxes/numericalflux.h"
#include "io.h"
#include "optimizers/optimizerbase.h"
#include "quadratures/quadraturebase.h"
#include "settings/config.h"
#include "solvers/sphericalharmonics.h"
9
#include "toolboxes/textprocessingtoolbox.h"
10

11
#include <mpi.h>
steffen.schotthoefer's avatar
steffen.schotthoefer committed
12

13
14
//#include <chrono>

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

17
18
19
    // Is this good (fast) code using a constructor list?
    _nMaxMomentsOrder = settings->GetMaxMomentDegree();
    _nTotalEntries    = GlobalIndex( _nMaxMomentsOrder, int( _nMaxMomentsOrder ) ) + 1;
20
21
22
23
24
25
26

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

steffen.schotthoefer's avatar
steffen.schotthoefer committed
28
29
30
31
32
33
34
35
36
37
38
39
40
41
    // 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];
            _sigmaS[n][j] = 1;
        }
    }

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

42
43
44
45
46
    // Initialize Entropy
    _entropy = EntropyBase::Create( _settings );

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

48
49
    // Initialize lagrange Multiplier
    _alpha = VectorVector( _nCells, Vector( _nTotalEntries, 0.0 ) );
50
51

    // Initialize and Pre-Compute Moments at quadrature points
52
53
    _basis = new SphericalHarmonics( _nMaxMomentsOrder );

54
    _moments = VectorVector( _nq, Vector( _nTotalEntries, 0.0 ) );
55
    ComputeMoments();
steffen.schotthoefer's avatar
steffen.schotthoefer committed
56
57
}

58
59
60
MNSolver::~MNSolver() {
    delete _entropy;
    delete _optimizer;
61
    delete _basis;
62
}
63

steffen.schotthoefer's avatar
steffen.schotthoefer committed
64
65
66
67
68
69
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;
}

70
void MNSolver::ComputeMoments() {
steffen.schotthoefer's avatar
steffen.schotthoefer committed
71
    double my, phi;
72
73

    for( unsigned idx_quad = 0; idx_quad < _nq; idx_quad++ ) {
74
75
        my  = _quadPointsSphere[idx_quad][0];
        phi = _quadPointsSphere[idx_quad][1];
76

77
        _moments[idx_quad] = _basis->ComputeSphericalBasis( my, phi );
78
    }
79
80
81
82
83
}

Vector MNSolver::ConstructFlux( unsigned idx_cell ) {

    // ---- Integration of Moment of flux ----
84
    double entropyL, entropyR, entropyFlux;
85

86
    Vector flux( _nTotalEntries, 0.0 );
87

88
    for( unsigned idx_quad = 0; idx_quad < _nq; idx_quad++ ) {
89
90
91

        entropyFlux = 0.0;    // Reset temorary flux

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

94
95
96
97
98
        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 {
99
                entropyR = _entropy->EntropyPrimeDual( blaze::dot( _alpha[_neighbors[idx_cell][idx_neigh]], _moments[idx_quad] ) );
100
            }
101
            entropyFlux += _g->Flux( _quadPoints[idx_quad], entropyL, entropyR, _normals[idx_cell][idx_neigh] );
102
        }
103
        flux += _moments[idx_quad] * ( _weights[idx_quad] * entropyFlux );
104
105

        // ------- Relizablity Reconstruction Step ----
106
    }
107
    return flux;
108
109
}

110
111
void MNSolver::ComputeRealizableSolution( unsigned idx_cell ) {
    double entropyReconstruction = 0.0;
112
    _sol[idx_cell]               = 0;
113
    for( unsigned idx_quad = 0; idx_quad < _nq; idx_quad++ ) {
114
        // Make entropyReconstruction a member vector, s.t. it does not have to be re-evaluated in ConstructFlux
115
        entropyReconstruction = _entropy->EntropyPrimeDual( blaze::dot( _alpha[idx_cell], _moments[idx_quad] ) );
116
        _sol[idx_cell] += _moments[idx_quad] * ( _weights[idx_quad] * entropyReconstruction );
117
118
119
    }
}

steffen.schotthoefer's avatar
steffen.schotthoefer committed
120
121
122
123
124
125
126
127
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?)
128
    VectorVector psiNew = _sol;
steffen.schotthoefer's avatar
steffen.schotthoefer committed
129
130
131
    double dFlux        = 1e10;
    Vector fluxNew( _nCells, 0.0 );
    Vector fluxOld( _nCells, 0.0 );
steffen.schotthoefer's avatar
steffen.schotthoefer committed
132
133

    double mass1 = 0;
steffen.schotthoefer's avatar
steffen.schotthoefer committed
134
    for( unsigned i = 0; i < _nCells; ++i ) {
135
136
        _solverOutput[i] = _sol[i][0];
        mass1 += _sol[i][0] * _areas[i];
steffen.schotthoefer's avatar
steffen.schotthoefer committed
137
138
    }

139
140
141
142
143
    dFlux   = blaze::l2Norm( fluxNew - fluxOld );
    fluxOld = fluxNew;

    Save( -1 );

steffen.schotthoefer's avatar
steffen.schotthoefer committed
144
145
146
    if( rank == 0 ) log->info( "{:10}   {:10}", "t", "dFlux" );
    if( rank == 0 ) log->info( "{:03.8f}   {:01.5e} {:01.5e}", -1.0, dFlux, mass1 );

147
148
149
150
    // Time measurement
    // auto start = chrono::steady_clock::now();
    // auto end   = chrono::steady_clock::now();

steffen.schotthoefer's avatar
steffen.schotthoefer committed
151
152
    // Loop over energies (pseudo-time of continuous slowing down approach)

153
    for( unsigned idx_energy = 0; idx_energy < _nEnergies; idx_energy++ ) {
steffen.schotthoefer's avatar
steffen.schotthoefer committed
154

155
156
        // Loop over the grid cells
        for( unsigned idx_cell = 0; idx_cell < _nCells; idx_cell++ ) {
steffen.schotthoefer's avatar
steffen.schotthoefer committed
157

158
            // ------- Reconstruction Step ----------------
159

160
            _optimizer->Solve( _alpha[idx_cell], _sol[idx_cell], _moments, idx_cell );
161

162
163
            // ------- Relizablity Reconstruction Step ----

164
            // ComputeRealizableSolution( idx_cell );
165
166

            // ------- Flux Computation Step --------------
167

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

171
            psiNew[idx_cell] = ConstructFlux( idx_cell );
172

173
            // ------ Finite Volume Update Step ------
174

175
176
177
            // NEED TO VECTORIZE
            for( unsigned idx_system = 0; idx_system < _nTotalEntries; idx_system++ ) {

178
                psiNew[idx_cell][idx_system] = _sol[idx_cell][idx_system] -
179
                                               ( _dE / _areas[idx_cell] ) * psiNew[idx_cell][idx_system] /* cell averaged flux */
180
                                               - _dE * _sol[idx_cell][idx_system] *
181
182
                                                     ( _sigmaA[idx_energy][idx_cell]                                    /* absorbtion influence */
                                                       + _sigmaS[idx_energy][idx_cell] * _scatterMatDiag[idx_system] ); /* scattering influence */
183
            }
steffen.schotthoefer's avatar
steffen.schotthoefer committed
184
        }
185
186

        // Update Solution
187
        _sol = psiNew;
188
189

        // pseudo time iteration output
190
191
        double mass = 0.0;
        for( unsigned idx_cell = 0; idx_cell < _nCells; ++idx_cell ) {
192
193
194
            fluxNew[idx_cell]       = _sol[idx_cell][0];    // zeroth moment is raditation densitiy we are interested in
            _solverOutput[idx_cell] = _sol[idx_cell][0];
            mass += _sol[idx_cell][0] * _areas[idx_cell];
steffen.schotthoefer's avatar
steffen.schotthoefer committed
195
        }
196

steffen.schotthoefer's avatar
steffen.schotthoefer committed
197
198
        dFlux   = blaze::l2Norm( fluxNew - fluxOld );
        fluxOld = fluxNew;
199
        if( rank == 0 ) log->info( "{:03.8f}   {:01.5e} {:01.5e}", _energies[idx_energy], dFlux, mass );
200
        Save( idx_energy );
steffen.schotthoefer's avatar
steffen.schotthoefer committed
201
202
    }
}
203
204
205
206
207
208
209

void MNSolver::Save() const {
    std::vector<std::string> fieldNames{ "flux" };
    std::vector<double> flux;
    flux.resize( _nCells );

    for( unsigned i = 0; i < _nCells; ++i ) {
210
        flux[i] = _sol[i][0];
211
212
213
214
215
216
217
218
219
220
221
222
    }
    std::vector<std::vector<double>> scalarField( 1, flux );
    std::vector<std::vector<std::vector<double>>> results{ scalarField };
    ExportVTK( _settings->GetOutputFile(), results, fieldNames, _mesh );
}

void MNSolver::Save( int currEnergy ) const {
    std::vector<std::string> fieldNames{ "flux" };
    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, fieldNames, _mesh );
}