Peano
Loading...
Searching...
No Matches
Rusanov.h
Go to the documentation of this file.
1// This file is part of the ExaHyPE2 project. For conditions of distribution and
2// use, please see the copyright notice at www.peano-framework.org
3#pragma once
4
5
6#include <vector>
7
8#include "../DGUtils.h"
9#include "../Functors.h"
10
14#include "../CellIntegral.h"
15#include "../Riemann.h"
16
17#include "tarch/la/Vector.h"
19
22
23#if defined(GPUOffloadingOMP)
25#endif
26
27
28namespace exahype2 {
29 namespace dg {
30 namespace rusanov {
103 ::exahype2::dg::MaxEigenvalue maxEigenvalue,
104 const tarch::la::Vector<Dimensions,double>& faceCentre,
106 double t,
107 double dt,
108 int order,
109 int unknowns,
110 int auxiliaryVariables,
111 int faceNumber,
112 const double* __restrict__ quadraturePoints,
113 bool useFlux,
114 bool useNcp,
115 const double* __restrict__ projectedValues,
116 double* __restrict__ solution
117 );
118
119
127 ::exahype2::dg::MaxEigenvalue maxEigenvalue,
128 const tarch::la::Vector<Dimensions,double>& faceCentre,
130 double t,
131 double dt,
132 int order,
133 int unknowns,
134 int auxiliaryVariables,
135 int faceNumber,
136 const double* __restrict__ quadraturePoints,
137 bool useFlux,
138 bool useNcp,
139 const double* __restrict__ projectedValues,
140 double* __restrict__ solution
141 );
142
143
151 const int order,
152 const int unknowns,
153 const int auxiliaryVariables,
154 Flux flux,
155 NonConservativeProduct nonConservativeProduct,
156 Source source,
157 PointSources pointSources,
158 const double* __restrict__ QuadratureNodes1d,
159 const double* __restrict__ MassMatrixDiagonal1d,
160 const double* __restrict__ StiffnessMatrix1d,
161 const double* __restrict__ DerivativeOperator1d,
162 bool evaluateFlux,
163 bool evaluateNonconservativeProduct,
164 bool evaluateSource,
165 bool evaluatePointSources
166 ) {
168 cellData,
169 order, unknowns, auxiliaryVariables,
170 flux, nonConservativeProduct, source, pointSources,
171 QuadratureNodes1d,
172 MassMatrixDiagonal1d,
173 StiffnessMatrix1d,
174 DerivativeOperator1d,
175 evaluateFlux, evaluateNonconservativeProduct, evaluateSource, evaluatePointSources
176 );
177 }
178
179
187 const int order,
188 const int unknowns,
189 const int auxiliaryVariables,
190 const double* __restrict__ MassMatrixDiagonal1d
191 ) {
193 cellData,
194 order, unknowns, auxiliaryVariables, MassMatrixDiagonal1d
195 );
196 }
197
198
205 const double* const __restrict__ faceQLeft,
206 const double* const __restrict__ faceQRight,
207 const double* const __restrict__ faceQBottom,
208 const double* const __restrict__ faceQUp,
209 int order,
210 int unknowns,
211 const int auxiliaryVariables,
212 const tarch::la::Vector<2,double>& cellSize,
213 const double* const __restrict__ BasisFunctionValuesLeft,
214 const double* __restrict__ MassMatrixDiagonal1d,
215 double* __restrict__ cellQ
216 ) {
218 faceQLeft, faceQRight, faceQBottom, faceQUp,
219 order, unknowns, auxiliaryVariables,
220 cellSize,
221 BasisFunctionValuesLeft,
222 MassMatrixDiagonal1d,
223 cellQ
224 );
225 }
226
227
234 const double* const __restrict__ faceQLeft,
235 const double* const __restrict__ faceQRight,
236 const double* const __restrict__ faceQBottom,
237 const double* const __restrict__ faceQUp,
238 const double* const __restrict__ faceQFront,
239 const double* const __restrict__ faceQBack,
240 int order,
241 int unknowns,
242 const int auxiliaryVariables,
243 const tarch::la::Vector<3,double>& cellSize,
244 const double* const __restrict__ BasisFunctionValuesLeft,
245 const double* __restrict__ MassMatrixDiagonal1d,
246 double* __restrict__ cellQ
247 ) {
249 faceQLeft, faceQRight, faceQBottom, faceQUp, faceQFront, faceQBack,
250 order, unknowns, auxiliaryVariables,
251 cellSize,
252 BasisFunctionValuesLeft,
253 MassMatrixDiagonal1d,
254 cellQ
255 );
256 }
257
258
259 template <
260 typename Solver,
261 int order,
262 int unknowns,
263 int auxiliaryVariables
264 >
267 bool evaluateFlux,
268 bool evaluateNonconservativeProduct,
269 bool evaluateSource,
270 bool evaluatePointSources
271 ) {
272 ::exahype2::dg::cellIntegral_patchwise_in_situ_GaussLegendre<Solver,order,unknowns,auxiliaryVariables>(
273 cellData,
274 evaluateFlux,
275 evaluateNonconservativeProduct,
276 evaluateSource,
277 evaluatePointSources
278 );
279 }
280 }
281 }
282}
283
Definition dg.py:1
void solveRiemannProblem_pointwise_in_situ(::exahype2::dg::Flux flux, ::exahype2::dg::NonConservativeProduct ncp, ::exahype2::dg::MaxEigenvalue maxEigenvalue, const tarch::la::Vector< Dimensions, double > &faceCentre, const tarch::la::Vector< Dimensions, double > &cellSize, double t, double dt, int order, int unknowns, int auxiliaryVariables, int faceNumber, const double *__restrict__ quadraturePoints, bool useFlux, bool useNcp, const double *__restrict__ projectedValues, double *__restrict__ solution)
Solve Riemann problem on face.
Definition Rusanov.cpp:6
static void integrateOverRiemannSolutionsAndAddToVolume_GaussLegendre(const double *const __restrict__ faceQLeft, const double *const __restrict__ faceQRight, const double *const __restrict__ faceQBottom, const double *const __restrict__ faceQUp, int order, int unknowns, const int auxiliaryVariables, const tarch::la::Vector< 2, double > &cellSize, const double *const __restrict__ BasisFunctionValuesLeft, const double *__restrict__ MassMatrixDiagonal1d, double *__restrict__ cellQ)
Delegate to generic implementation.
Definition Rusanov.h:204
void cellIntegral_patchwise_in_situ_GaussLegendre(::exahype2::CellData< double, double > &cellData, bool evaluateFlux, bool evaluateNonconservativeProduct, bool evaluateSource, bool evaluatePointSources)
Definition Rusanov.h:265
static void cellIntegral_patchwise_in_situ_GaussLegendre_functors(::exahype2::CellData< double, double > &cellData, const int order, const int unknowns, const int auxiliaryVariables, Flux flux, NonConservativeProduct nonConservativeProduct, Source source, PointSources pointSources, const double *__restrict__ QuadratureNodes1d, const double *__restrict__ MassMatrixDiagonal1d, const double *__restrict__ StiffnessMatrix1d, const double *__restrict__ DerivativeOperator1d, bool evaluateFlux, bool evaluateNonconservativeProduct, bool evaluateSource, bool evaluatePointSources)
Delegate to generic implementation in parent directory.
Definition Rusanov.h:149
static void multiplyWithInvertedMassMatrix_GaussLegendre(::exahype2::CellData< double, double > &cellData, const int order, const int unknowns, const int auxiliaryVariables, const double *__restrict__ MassMatrixDiagonal1d)
Delegate to generic implementation.
Definition Rusanov.h:185
void solveRiemannProblem_pointwise_in_situ_with_gradient_projection(::exahype2::dg::Flux flux, ::exahype2::dg::NonConservativeProduct ncp, ::exahype2::dg::MaxEigenvalue maxEigenvalue, const tarch::la::Vector< Dimensions, double > &faceCentre, const tarch::la::Vector< Dimensions, double > &cellSize, double t, double dt, int order, int unknowns, int auxiliaryVariables, int faceNumber, const double *__restrict__ quadraturePoints, bool useFlux, bool useNcp, const double *__restrict__ projectedValues, double *__restrict__ solution)
Definition Rusanov.cpp:127
std::function< double(const double *__restrict__ Q, const tarch::la::Vector< Dimensions, double > &x, double t, double dt, int normal) MaxEigenvalue)
Definition Functors.h:95
std::function< void(const double *__restrict__ Q, const double *__restrict__ dQdx, const tarch::la::Vector< Dimensions, double > &x, double t, double dt, int normal, double *__restrict__ F) NonConservativeProduct)
Definition Functors.h:54
void multiplyWithInvertedMassMatrix_GaussLegendre(::exahype2::CellData< double, double > &cellData, const int order, const int unknowns, const int auxiliaryVariables, const double *__restrict__ MassMatrixDiagonal1d)
Final step of DG algorithm.
std::function< std::vector< PointSource >(const double *__restrict__ Q, const tarch::la::Vector< Dimensions, double > &cellCentre, const tarch::la::Vector< Dimensions, double > &h, double t, double dt) PointSources)
This is the only routine within the DG framework which accepts the dimensions of the underlying cell ...
Definition Functors.h:87
void cellIntegral_patchwise_in_situ_GaussLegendre_functors(::exahype2::CellData< double, double > &cellData, const int order, const int unknowns, const int auxiliaryVariables, Flux flux, NonConservativeProduct nonconservativeProduct, Source source, PointSources pointSources, const double *__restrict__ QuadratureNodes1d, const double *__restrict__ MassMatrixDiagonal1d, const double *__restrict__ StiffnessMatrix1d, const double *__restrict__ DerivativeOperator1d, bool evaluateFlux, bool evaluateNonconservativeProduct, bool evaluateSource, bool evaluatePointSources)
std::function< void(const double *__restrict__ Q, const tarch::la::Vector< Dimensions, double > &x, double t, double dt, double *__restrict__ S) Source)
Source functor.
Definition Functors.h:75
std::function< void(const double *__restrict__ Q, const tarch::la::Vector< Dimensions, double > &x, double t, double dt, int normal, double *__restrict__ F) Flux)
Flux functor.
Definition Functors.h:44
void integrateOverRiemannSolutionsAndAddToVolume_GaussLegendre(const double *const __restrict__ faceQLeft, const double *const __restrict__ faceQRight, const double *const __restrict__ faceQBottom, const double *const __restrict__ faceQUp, int order, int unknowns, const int auxiliaryVariables, const tarch::la::Vector< 2, double > &cellSize, const double *const __restrict__ BasisFunctionValuesLeft, const double *__restrict__ MassMatrixDiagonal1d, double *__restrict__ cellQ)
Given a numerical flux at the various faces, this computes and adds the Riemann integral of this flux...
Definition Riemann.cpp:332
For the generic kernels that I use here most of the time.
Definition CellAccess.h:13
Representation of a number of cells which contains all information that's required to process the sto...
Definition CellData.h:77
Simple vector class.
Definition Vector.h:150