MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
maia::dg::mortar Namespace Reference

Functions

void calcMortarProjectionMatrixP (const MInt sourcePolyDeg, const MFloat *const sourceNodes, const MFloat *const sourceBaryWeights, const MInt targetPolyDeg, const MFloat *const targetNodes, MFloat *const matrix)
 
void calcMortarProjectionMatrixHForward (const MInt polyDeg, const MFloat *const nodes, const MFloat *const wBary, const MInt position, MFloat *const matrix)
 
void calcMortarProjectionMatrixHReverse (const MInt polyDeg, const MFloat *const nodes, const MFloat *const wBary, const MInt position, MFloat *const matrix)
 

Variables

const MBool forward = true
 
const MBool reverse = false
 
const MInt lower = 0
 
const MInt upper = 1
 

Function Documentation

◆ calcMortarProjectionMatrixHForward()

void maia::dg::mortar::calcMortarProjectionMatrixHForward ( const MInt  polyDeg,
const MFloat *const  nodes,
const MFloat *const  wBary,
const MInt  position,
MFloat *const  matrix 
)
inline

Calculate h-refinement forward projection matrix from a coarse element to the mortar element (2D).

Author
Sven Berger, Michael Schlottke
Date
2015-02-27
Parameters
[in]polyDegPolynomial degree of element/mortar.
[in]nodesNode locations of Lagrange polynomial in [-1,1].
[in]wBaryBarycentric weights of Lagrange polynomial.
[in]positionDetermines lower or upper element.
[out]matrixProjection matrix of size (polyDeg + 1) x (polyDeg + 1).

Projection: lower: upper: | _ | | | | /| | | | / | \ | | _| | | | | |

See also pp. 25-26 in Sven Berger: Implementation and validation of an adaptive hp-refinement method for the discontinuous Galerkin spectral element method. Master thesis, RWTH Aachen University, 2014.

Definition at line 106 of file dgcartesianmortar.h.

110 {
111 const MInt noNodes = polyDeg + 1;
112
113 // Forward matrix is the Lagrange polynomials evaluated at other nodes
114 const MFloat shift = (position == dg::mortar::lower) ? -0.5 : 0.5;
115 for(MInt i = 0; i < noNodes; i++) {
117 0.5 * nodes[i] + shift, polyDeg, nodes, wBary, &matrix[i * noNodes]);
118 }
119}
int32_t MInt
Definition: maiatypes.h:62
double MFloat
Definition: maiatypes.h:52
void calcLagrangeInterpolatingPolynomials(const MFloat x, const MInt polyDeg, const MFloat *nodes, const MFloat *wBary, MFloat *polynomials)
Calculates the values of the Lagrangian polynomials l_j for a given point x in [-1,...

◆ calcMortarProjectionMatrixHReverse()

void maia::dg::mortar::calcMortarProjectionMatrixHReverse ( const MInt  polyDeg,
const MFloat *const  nodes,
const MFloat *const  wBary,
const MInt  position,
MFloat *const  matrix 
)
inline

Calculate h-refinement reverse projection matrix from a mortar element to the coarse element (2D).

Author
Sven Berger, Michael Schlottke
Date
2015-02-27
Parameters
[in]polyDegPolynomial degree of element/mortar.
[in]nodesNode locations of Lagrange polynomial in [-1,1].
[in]wBaryBarycentric weights of Lagrange polynomial.
[in]positionDetermines lower or upper element.
[out]matrixProjection matrix of size (polyDeg + 1) x (polyDeg + 1).

Projection: lower: upper: | | | | | / | | _ | |/_ | |\ | | \ | | | | |

See also pp. 25-26 in Sven Berger: Implementation and validation of an adaptive hp-refinement method for the discontinuous Galerkin spectral element method. Master thesis, RWTH Aachen University, 2014. and pp. 11-14 in Patrick Antony: Development of a coupled discontinuous Galerkin -finite volume scheme, Bachelor thesis, RWTH Aachen University, 2018.

Contrary to the method described as in the first reference by Berger, the reverse projection operator performs the following steps (see second reference by Antony): 1) Transform input from any quadrature node types to Gauss nodes 2) Apply mortar projection on Gauss nodes (identical to the original formulation in Berger) 3) Transform result from Gauss nodes to original quadrature node types

Therefore, the projection operator contructed here is a matrix multiplication of three matrices. If the integration nodes are Gauss nodes as well, this is equivalent to the method in the first reference.

Definition at line 160 of file dgcartesianmortar.h.

164 {
165 const MInt noNodes = polyDeg + 1;
166
167 // Calculate nodes and quadrature weights for projection (always Gauss)
168 MFloatVector projectionNodes;
169 MFloatVector projectionWInt;
170 projectionNodes.resize(noNodes);
171 projectionWInt.resize(noNodes);
172 maia::dg::interpolation::calcLegendreGaussNodesAndWeights(polyDeg, &projectionNodes[0], &projectionWInt[0]);
173
174 // Calculate barycentric weights for projection nodes
175 MFloatVector projectionWBary;
176 projectionWBary.resize(noNodes);
177 maia::dg::interpolation::calcBarycentricWeights(polyDeg, &projectionNodes[0], &projectionWBary[0]);
178
179 // Calculate projection matrix from mortar element to coarse element in projection nodes
180 const MFloat shift = (position == dg::mortar::lower) ? -0.5 : 0.5;
181 MFloatTensor mortarToCoarseElement(noNodes, noNodes);
182 for(MInt i = 0; i < noNodes; i++) {
184 polyDeg,
185 &projectionNodes[0],
186 &projectionWBary[0],
187 &mortarToCoarseElement(i, 0));
188 }
189
190 // Transpose and multiply
191 for(MInt i = 0; i < noNodes; i++) {
192 for(MInt j = 0; j < i + 1; j++) {
193 const MFloat temp = mortarToCoarseElement(i, j);
194 mortarToCoarseElement(i, j) = 0.5 * mortarToCoarseElement(j, i) * projectionWInt[j] / projectionWInt[i];
195 mortarToCoarseElement(j, i) = 0.5 * temp * projectionWInt[i] / projectionWInt[j];
196 }
197 }
198
199
200 // Transformation matrix from grid to projection nodes (this is the identity when using Gauss
201 // nodes anyways)
202 MFloatTensor p2dg(noNodes, noNodes);
203 calcMortarProjectionMatrixP(polyDeg, &projectionNodes[0], &projectionWBary[0], polyDeg, nodes, &p2dg[0]);
204
205 // Transformation matrix from projection to grid nodes (this is the identity when using Gauss
206 // nodes anyways)
207 MFloatTensor dg2p(noNodes, noNodes);
208 calcMortarProjectionMatrixP(polyDeg, nodes, wBary, polyDeg, &projectionNodes[0], &dg2p[0]);
209
210 // The total projection matrix m is a product of the transformation matrices to(dg2p) and
211 // from(p2dg) quadrature nodes and the projection matrix in quadrature nodes, in this order: m =
212 // p2dg * mortarToCoarseElement * dg2p (where "*" denotes a matrix-matrix product)
213 MFloatTensor m(matrix, noNodes, noNodes);
214 m.set(0.0);
215 for(MInt q = 0; q < noNodes; q++) {
216 for(MInt j = 0; j < noNodes; j++) {
217 for(MInt i = 0; i < noNodes; i++) {
218 for(MInt n = 0; n < noNodes; n++) {
219 m(q, j) += p2dg(q, i) * mortarToCoarseElement(i, n) * dg2p(n, j);
220 }
221 }
222 }
223 }
224}
size_type resize(size_type n0, size_type n1=1, size_type n2=1, size_type n3=1, size_type n4=1)
Deletes the old data structure and creates a new one with the requested dimensions.
Definition: tensor.h:230
void calcBarycentricWeights(MInt Nmax, const MFloat *nodes, MFloat *weights)
Calculates the barycentric weights for Lagrange interpolation at thei specified nodes.
void calcLegendreGaussNodesAndWeights(MInt Nmax, MFloat *nodes, MFloat *wInt)
Calculate the Gauss integration nodes and weight for the Legendre polynomials on the interval [-1,...

◆ calcMortarProjectionMatrixP()

void maia::dg::mortar::calcMortarProjectionMatrixP ( const MInt  sourcePolyDeg,
const MFloat *const  sourceNodes,
const MFloat *const  sourceBaryWeights,
const MInt  targetPolyDeg,
const MFloat *const  targetNodes,
MFloat *const  matrix 
)
inline

Calculate p-refinement mortar projection matrix from source to target polynomial degree.

Author
Sven Berger, Michael Schlottke
Date
2015-02-27
Parameters
[in]sourcePolyDegSource polynomial degree.
[in]sourceNodesNode locations in [-1,1] of source polynomial degree.
[in]sourceBaryWeightsBarycentric weights of source polynomial degree.
[in]targetPolyDegTarget polynomial degree.
[in]targetNodesNode locations in [-1,1] of target polynomial degree.
[out]matrixProjection matrix of size (targetPolyDeg + 1) x (sourcePolyDeg + 1).

This method can be used to calculate the projection matrix between any two polynomial degrees. For each forward projection, the corresponding reverse projection can be found by reversing the source and target polynomial degrees.

When projecting from polynomial degree N_s to N_t, then sourcePolyDeg = N_s and targetPolyDeg = N_t. The reverse projection can then be obtained by setting sourcePolyDeg = N_t and targetPolyDeg = N_s.

The content of the projection matrix P_ij is from polynomial degree N_s to N_t is as follows:

P_ij = l_j(x_i)

where l_j are the N_s + 1 Lagrange polynomials of degree N_s while x_i are the N_t + 1 Gauss nodes of degree N_t.

See also p. 22, Eq. 3.33 and Eq. 3.34 in Sven Berger: Implementation and validation of an adaptive hp-refinement method for the discontinuous Galerkin spectral element method. Master thesis, RWTH Aachen University, 2014.

Definition at line 63 of file dgcartesianmortar.h.

68 {
69 const MInt sourceNoNodes = sourcePolyDeg + 1;
70 const MInt targetNoNodes = targetPolyDeg + 1;
71
72 // Matrix is Lagrange polynomials of source polynomial degree evaluated at
73 // nodes of target polynomial degree
74 for(MInt i = 0; i < targetNoNodes; i++) {
75 interpolation::calcLagrangeInterpolatingPolynomials(
76 targetNodes[i], sourcePolyDeg, sourceNodes, sourceBaryWeights, &matrix[i * sourceNoNodes]);
77 }
78}

Variable Documentation

◆ forward

const MBool maia::dg::mortar::forward = true

Definition at line 21 of file dgcartesianmortar.h.

◆ lower

const MInt maia::dg::mortar::lower = 0

Definition at line 23 of file dgcartesianmortar.h.

◆ reverse

const MBool maia::dg::mortar::reverse = false

Definition at line 22 of file dgcartesianmortar.h.

◆ upper

const MInt maia::dg::mortar::upper = 1

Definition at line 24 of file dgcartesianmortar.h.