Peano
Loading...
Searching...
No Matches
DGUtils.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 "tarch/la/la.h"
9#include "tarch/la/Vector.h"
11#include "tarch/Assertions.h"
12
16#include "peano4/utils/Loop.h"
17
18#include "Functors.h"
19
21
22
23
24namespace exahype2 {
25 namespace dg {
46 template<typename T>
51 int polynomialOrder,
52 const T* __restrict__ quadraturePoints
53 );
54
69 const tarch::la::Vector<3,double>& cellSize,
70 const tarch::la::Vector<3,int>& index,
71 const double* __restrict__ quadratureWeights
72 );
73
74
82 int getNodesPerCell(int nodesPerAxis);
83
84
85 /* returns a vector with the powers of nodesPerAxis up to the Dimensions - 1
86 ** e.g. for dims = 2, nodesPerAxis = 2 will return (1, 2)
87 ** for dims = 3, nodesPerAxis = 4 will return (1, 4, 16)
88 **
89 ** This is used to iterate over the nodes of a cell in a given dimension,
90 ** since the cells are all in the same array one after another advancing by
91 ** n*strides[i] will move you to the cell which is n steps in direction i
92 ** from the current cell.
93 ** i.e. node+1 would be the "right" neighbour, node+nodesPerAxis would be the
94 ** "upper" neighbour, and so on.
95 */
97
98
99 /*for a given linear node number, returns the index of that node
100 **in cartesian coordinate system.
101 **e.g. for dims = 2, nodesPerAxis = 4, there will be 16 total nodes
102 ** so node 0 is [0,0], 1 is [1,0] and so on until 15 which is [3,3]
103 ** for dims = 3, nodesPerAxis = 4, there will be 16*4 = 64 total nodes
104 ** ranging from 0 which is [0,0,0] to 63 which is [3,3,3]
105 */
107 int node,
109 );
110
111 /*
112 * for a given node index within the cell, returns the index of the neighbouring
113 * node on the surface in a given direction with a given orientation.
114 *
115 * Ordering on the surface is the same as for the previous functions, but only
116 * considering one "side" on each face
117 *
118 * So with the same notation as the previous function, if given a 2*2 cell
119 * in 2 dimensions, there are a total of 4 faces (left,right,down,up), each
120 * with 2 nodes.
121 * If given index which corresponds to [0,1],
122 * direction 0, orientation 0 would return node 1 on the left face
123 * direction 1, orientation 1 would return node 0 on the upper face
124 * 6 7
125 * 1 [ o o ] 3
126 * 0 [ o o ] 2
127 * 4 5
128 */
130 const tarch::la::Vector<Dimensions,int>& indexCell,
131 const int direction,
132 const int orientation,
133 const int nodesPerAxis
134 );
135
136 /*
137 * Computes the gradient of Q in a given direction for a given node
138 *
139 * Mostly used for the non-conservative product.
140 *
141 * Since the solution is represented by the sum of basis functions, its
142 * derivative at a given node is the sum of the derivatives of the basis
143 * functions evaluated at that node.
144 * The derivative operator contains the value of the derivatives evaluated
145 * at the quadrature nodes, therefore it only needs to be multiplied by
146 * their coefficients.
147 *
148 *
149 * u(x) = \sum_i \phi_i(x) * u_i => du(x)/dx = d\phi_i(x_i)/dx * u_i
150 *
151 * @param invDx: the inverse of the size of the cell in any one dimension, this
152 * assumes that the cells are cubic
153 *
154 * @param derivativeOperator: array containing the value of the derivatives of
155 * the n base functions of the dg polynomial evaluated at the n quadrature points
156 * (in 1 dimension)
157 *
158 * @param scalarIndex: index of the node at which the gradient should be evaluated
159 *
160 * @param gradQ: where the gradient values should be placed, has size unknowns*Dimensions
161 * since there is a gradient in each direction and the gradient is only nonzero for the
162 * unknowns
163 *
164 * @param strideQ: total number of variables over which the gradient should be computed
165 */
166 void computeGradient(
167 const double* __restrict__ const QCell,
168 const double* __restrict__ const derivativeOperator,
169 const double invDx,
170 const int nodesPerAxis,
171 const int strideQ,
172 const int scalarIndex,
173 double* __restrict__ gradQ
174 );
175
176 /*
177 * Substracts contents of one cell array from another.
178 */
179 void subtractCell(
180 double* __restrict__ QOut,
181 const double* __restrict__ Qsubstract,
182 const int order,
183 const int unknowns,
184 const int auxiliaryVariables
185 );
186
187
188 /*
189 * Plot patch
190 *
191 * usually used for debugging
192 */
193
194 std::string plotCell(
195 const double* __restrict__ Q,
196 const int order,
197 const int unknowns,
198 const int auxiliaryVariables
199 );
200
201 /*
202 * Plot Face
203 *
204 * Usually used for debugging. We plot all unknowns, i.e. both the unknowns
205 * form the left and the right adjacent side. The order is linear, i.e. you
206 * have to take the normal into account when you want to find out how the
207 * unknowns are spatially arranged.
208 *
209 * See the discussion around ExaHyPE's lexicographic ordering.
210 */
211 std::string plotFace(
212 const double* __restrict__ Q,
213 const int order,
214 const int unknowns,
215 const int auxiliaryVariables,
216 int normal,
217 int numberOfQuantitiesProjectedOntoFace
218 );
219
257 template<typename QStoreType>
258 QStoreType evaluatePolynomial(
260 int order,
261 const double* __restrict__ QuadratureNodes1d,
262 int unknownsPerDoF,
263 const QStoreType* __restrict__ Q,
265 int unknown
266 );
267
272 int unknownsPlusAuxiliaryVariables,
273 int order,
274 int numberOfProjectedQuantities,
275 int normal,
276 int isRightFaceHalf,
277 const double* __restrict__ srcQ,
278 double* __restrict__ destQ
279 );
280 }
281}
282
283
284#include "DGUtils.cpph"
Definition dg.py:1
std::string plotFace(const double *__restrict__ Q, const int order, const int unknowns, const int auxiliaryVariables, int normal, int numberOfQuantitiesProjectedOntoFace)
Definition DGUtils.cpp:176
tarch::la::Vector< Dimensions, int > getStrides(int nodesPerAxis)
Definition DGUtils.cpp:47
tarch::la::Vector< Dimensions, double > getQuadraturePoint(const tarch::la::Vector< Dimensions, double > &cellCentre, const tarch::la::Vector< Dimensions, double > &cellSize, const tarch::la::Vector< Dimensions, int > &index, int polynomialOrder, const T *__restrict__ quadraturePoints)
Construct location of a quadrature point.
Definition DGUtils.cpph:33
void computeGradient(const double *__restrict__ const QCell, const double *__restrict__ const derivativeOperator, const double invDx, const int nodesPerAxis, const int strideQ, const int scalarIndex, double *__restrict__ gradQ)
Definition DGUtils.cpp:110
double getQuadratureWeight(const tarch::la::Vector< 3, double > &cellSize, const tarch::la::Vector< 3, int > &index, const double *__restrict__ quadratureWeights)
Compute integral over shape function over cell defined by index.
Definition DGUtils.cpp:25
tarch::la::Vector< Dimensions, int > getIndex(int node, tarch::la::Vector< Dimensions, int > strides)
Definition DGUtils.cpp:55
int getNodesPerCell(int nodesPerAxis)
The number of nodes in a cell is basically the input to the power of d.
Definition DGUtils.cpp:39
void subtractCell(double *__restrict__ QOut, const double *__restrict__ Qsubstract, const int order, const int unknowns, const int auxiliaryVariables)
Definition DGUtils.cpp:141
QStoreType evaluatePolynomial(const peano4::datamanagement::CellMarker &marker, int order, const double *__restrict__ QuadratureNodes1d, int unknownsPerDoF, const QStoreType *__restrict__ Q, const tarch::la::Vector< Dimensions, double > &x, int unknown)
Evaluate the DG polynomial.
Definition DGUtils.cpph:3
void copyOneSideOfFaceProjection(int unknownsPlusAuxiliaryVariables, int order, int numberOfProjectedQuantities, int normal, int isRightFaceHalf, const double *__restrict__ srcQ, double *__restrict__ destQ)
Delegate to PatchUtils of the Finite Volume scheme.
Definition DGUtils.cpp:4
std::string plotCell(const double *__restrict__ Q, const int order, const int unknowns, const int auxiliaryVariables)
Definition DGUtils.cpp:160
int cellIndexToHullIndex(const tarch::la::Vector< Dimensions, int > &indexCell, const int direction, const int orientation, const int nodesPerAxis)
Definition DGUtils.cpp:71
For the generic kernels that I use here most of the time.
Definition CellAccess.h:13
Simple vector class.
Definition Vector.h:150