Peano
Loading...
Searching...
No Matches
CellIntegral.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#include <algorithm>
8
9#include "DGUtils.h"
10#include "Functors.h"
11
12#include "tarch/la/Vector.h"
14
17
18#include "exahype2/CellData.h"
19
20
21namespace exahype2 {
22 namespace dg {
23 /*
24 * Volumetric kernel for Gauss-Legendre shape functions using functors
25 *
26 * Computes the @f$ \hat Q @f$ from the generic Discontinuous Galerkin
27 * explanation. That is the volumetric part of the weak formulation. It
28 * does not include the previous time step's solution, and it is not (yet)
29 * multiplied with the time step size or the inverse of the mass matrix.
30 *
31 * The whole implementation works cell-wisely, i.e. runs through cell by
32 * cell. Within each cell, it works dof-wisely, i.e. it runs through each
33 * individual dof and accumulates all the PDE contriutions to this dof.
34 * This second loop is for loop over node, and the d-dimensional index
35 * of node if given by index.
36 *
37 *
38 * ## Solver flavour
39 *
40 * We rely on Gauss-Legendre quadrature nodes on tensor product spaces, i.e.
41 * we exploit the fact that the mass matrix has diagonal form. Consequently,
42 * we can construct all involved cell-local operators on-the-fly by
43 * multiplying their 1d counterparts.
44 *
45 * We run over the patches on by one (patchwise). We do not exploit that we
46 * have to do the same operation for all dofs of all patches.
47 *
48 * We evaluate all operators in-situ, i.e. get their output and apply it
49 * immediately. See the discussion below that each dof's flux evaluation
50 * needs all the fluxes in that row of dofs. Consequently, we evaluate
51 * each flux p times if we have p quadrature points per axis. This is
52 * redundant work. In return, in-situ implies that we don't have to store
53 * data temporarily.
54 *
55 *
56 * ## Flux evaluation
57 *
58 * This version of the volumetric kernel uses the Gauss Legendre property,
59 * i.e. the orthogonality of the polynomials, plus the tensor product style
60 * of the shape functions aggressively.
61 *
62 * For the flux, we first loop over the dimensions (dim) as we have
63 * fluxes in all d directions. We are interested how the dof in the
64 * position index changes due to all the fluxes around. So we want to
65 * find out what happens if we test against the test function @f$ \phi @f$
66 * in the index dof.
67 *
68 * For the flux, tensor products and orthogonality means that the flux
69 * in x-direction in any point really only interferes with the unknowns
70 * along the x-direction along this same line: We test @f$ F(Q) \partial_x \phi @f$
71 * and assume that @f$ \phi = \phi _{1d}(x)\phi _{1d}(y)\phi _{1d}(z) @f$.
72 * Furthermore, lets represent @f$ F(Q) @f$ with the same polynomials,
73 * too. In the product of the flux with the test function (which now
74 * consists of 2d ingredients), all the polynomials besides the @f$ \partial _x \phi(x) @f$
75 * are orthogonal - the fact that the polynomials have been chosen as
76 * orthogonal ones means that their derivatives are obviously not automatically
77 * orthogonal.
78 *
79 * So all the shape functions spanning F where the y and z coordinates do
80 * not match cancel out in the example above if the test functions don't
81 * have the same coordinates. We are left with all the shapes along the
82 * x-axis within the cubic cell. Analogous arguments hold for y and z.
83 *
84 *
85 * ## Algebraic source
86 *
87 * The treatment of the source term is close to trivial: We assume that
88 * we can represent the source @f$ S(Q) @f$ in the shape space. Again,
89 * we exploit the orthogonality. If we want to test what the impact of
90 * the source is if we test against the index-th test function, the
91 * orthogonality implies that we only have take S(Q) in index multiplied
92 * with its shape function. We get a mass matrix entry.
93 *
94 *
95 * ## Non-conservative product
96 *
97 * The non-conservative product once again relies on a very crude numerical
98 * approximation: We assume that @f$ B_x(Q) \partial _x Q @f$ can be
99 * represented within our orthogonal shape space. Consequently, it is
100 * tried exactly in the same way as the algebraic source.
101 *
102 * This approach seems to be super inaccurate, but we have to keep in mind
103 * that B can be highly non-linear. So we don't know anything about the
104 * analytical shape of @f$ B_x(Q) \partial _x Q @f$ anyway: It could be
105 * any type of polynomial. So it is convenient to assume that it is once
106 * again from the same space as @f$ \phi @f$.
107 *
108 *
109 * ## Signs in front of the individual terms
110 *
111 * We rely on the formulation as discussed on the page Discontinuous
112 * Galerkin. See dg.h or the doxygen html pages for details. One important
113 * fact to keep in mind is that this routine computes $\partial _t Q$.
114 * However, all PDE terms besides the source arise on the left-hand side of
115 * the PDE, i.e. on the same side as the time derivative. We therefore have
116 * to bring these terms to the right which gives us an additional minus
117 * sign.
118 *
119 * Here's the rationale behind the volumetric terms:
120 *
121 * - The source term is already on the right-hand side in our PDE
122 * formulation. Therefore, we can simply add it.
123 * - The flux is on the left-hand side of the PDE formulation where we find
124 * a term @f$ div F(Q) @f$. We bring the derivative
125 * over to the test function and then bring the whole term to the
126 * right-hand side Both steps yield a minus sign, which means that the
127 * minus signs cancel out.
128 * - The ncp enters the volumetric term as it is. No partial derivative is
129 * exploited. It hence enters the solution with its original sign.
130 * However, we have to bring it to the right-hand side for the time
131 * derivatives. Therefore, we need a minus sign here.
132 *
133 *
134 * ## Arguments
135 *
136 * @param cellData The routine computes the DG cell contributions over
137 * multiple cells in one go. cellData is the wrapper around these cells,
138 * i.e. holds the pointers to input and ouput data, time step sizes, ...
139 * However, all cells carry the same polynomial degree and solve the
140 * same PDE.
141 *
142 * @param order Order of underlying dg polynomials. This specifies the number
143 * of nodes in each dimension, which is order+1 (you need two points to specify
144 * a linear polynomial for example). The order has to be the same for all cells
145 * handed over via cellData.
146 *
147 * @param unknowns Mumber of unknowns whose values must be updated. If you
148 * solve a PDE with n eqations over a function Q with n+m entries, then n
149 * is the number of unknonws that evolve according to the hyperbolic system
150 * of PDEs, and m are auxiliary (material) parameters. They do not change due
151 * to the PDE (though any user function might want to alter them).
152 *
153 * @param auxiliaryVariables. Auxiliary variables. See description of unknowns.
154 *
155 * @param quadratureNodes Location of quadrature nodes along one dimension in a
156 * reference element (length). We usually expect the nodes to be
157 * Gauss-Legendre or Gauss-Lobatto quadrature points, but we don't really
158 * bake these considerations into the present routine. However, we do
159 * hardcode the information that we work with a Cartesian tensor product
160 * layout of the dofs, i.e. the (x,y,z) position can be computed per
161 * entry with y and z being independent of x. See getQuadraturePoint()
162 * for some documentation.
163 *
164 */
167 const int order,
168 const int unknowns,
169 const int auxiliaryVariables,
170 Flux flux,
171 NonConservativeProduct nonconservativeProduct,
172 Source source,
173 PointSources pointSources,
174 const double* __restrict__ QuadratureNodes1d,
175 const double* __restrict__ MassMatrixDiagonal1d,
176 const double* __restrict__ StiffnessMatrix1d,
177 const double* __restrict__ DerivativeOperator1d,
178 bool evaluateFlux,
179 bool evaluateNonconservativeProduct,
180 bool evaluateSource,
181 bool evaluatePointSources
182 );
183
184
185 template <
186 typename Solver,
187 int order,
188 int unknowns,
189 int auxiliaryVariables
190 >
193 bool evaluateFlux,
194 bool evaluateNonconservativeProduct,
195 bool evaluateSource,
196 bool evaluatePointSources
197 );
198
199
228 const int order,
229 const int unknowns,
230 const int auxiliaryVariables,
231 const double* __restrict__ MassMatrixDiagonal1d
232 );
233
234
241 const int order,
242 const int unknowns,
243 const int auxiliaryVariables,
244 MaxEigenvalue maxEigenvalue,
245 const double* __restrict__ QuadratureNodes1d
246 );
247
248 namespace internal {
259 void copySolution(
260 const double* __restrict__ cellQin,
261 const int order,
262 const int unknowns,
263 const int auxiliaryVariables,
265 double* __restrict__ cellQout
266 );
267 }
268 }
269}
270
271
272#include "CellIntegral.cpph"
273
274
Definition dg.py:1
void copySolution(const double *__restrict__ cellQin, const int order, const int unknowns, const int auxiliaryVariables, const tarch::la::Vector< Dimensions, int > &node, double *__restrict__ cellQout)
Copy solution over for one node.
void cellIntegral_patchwise_in_situ_GaussLegendre(::exahype2::CellData< double, double > &cellData, bool evaluateFlux, bool evaluateNonconservativeProduct, bool evaluateSource, bool evaluatePointSources)
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 reduceMaxEigenvalue_patchwise_functors(::exahype2::CellData< double, double > &cellData, const int order, const int unknowns, const int auxiliaryVariables, MaxEigenvalue maxEigenvalue, const double *__restrict__ QuadratureNodes1d)
Compute the maximum eigenvalues over a sequence of cells and store the result in the respective CellD...
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
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