MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
DgInterpolation Class Reference

Class stores precalculated values for interpolation & integration on the reference interval [-1,1]. More...

#include <dgcartesianinterpolation.h>

Collaboration diagram for DgInterpolation:
[legend]

Public Member Functions

 DgInterpolation ()
 Default constructor only sets default values for member variables. More...
 
 DgInterpolation (const MInt polyDeg, const DgPolynomialType polyType, const MInt noNodes, const DgIntegrationMethod intMethod, const MBool sbpMode, const MString sbpOperator)
 Constructor passes arguments to init(). More...
 
 ~DgInterpolation ()
 Destructor clears all member variables. More...
 
void init (const MInt polyDeg, const DgPolynomialType polyType, const MInt noNodes, const DgIntegrationMethod intMethod, const MBool sbpMode, const MString sbpOperator)
 Sets the member variables and calls the appropriate functions to calculate the nodes and weights etc. More...
 
void initInterpolation (const MInt polyDeg, const DgPolynomialType polyType, const MInt noNodes, const DgIntegrationMethod intMethod, const MBool sbpMode, const MString sbpOperator)
 Sets the member variables neccessary for the interpolation between sets of nodes. More...
 

Public Attributes

MBool m_sbpMode
 
MString m_sbpOperator
 
MInt m_polyDeg
 
MInt m_noNodes
 
MInt m_polyType
 
MInt m_intMethod
 
MFloatVector m_nodes
 
MFloatVector m_wInt
 
MFloatVector m_wBary
 
MFloatMatrix m_Dhat
 
MFloatVector m_LFace [2]
 
MFloatVector m_LhatFace [2]
 

Detailed Description

Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2012-11-05

The following two sources may be cited:

Kopriva09: David A. Kopriva, Implementing Spectral Methods for Partial Differential Equations, Springer, 2009
HesthavenWarburton08: J. S. Hesthaven and T. Warburton, Nodal Discontinuous Galerkin Methods, Springer, 2008

Definition at line 35 of file dgcartesianinterpolation.h.

Constructor & Destructor Documentation

◆ DgInterpolation() [1/2]

DgInterpolation::DgInterpolation ( )
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2012-11-07

Definition at line 24 of file dgcartesianinterpolation.cpp.

◆ DgInterpolation() [2/2]

DgInterpolation::DgInterpolation ( const MInt  polyDeg,
const DgPolynomialType  polyType,
const MInt  noNodes1D,
const DgIntegrationMethod  intMethod,
const MBool  sbpMode,
const MString  sbpOperator 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2012-11-07
Parameters
[in]polyDegMaximum polynomial degree.
[in]polyTypePolynomial type.
[in]intMethodIntegration method.

Definition at line 41 of file dgcartesianinterpolation.cpp.

47 : m_sbpMode(false),
48 m_sbpOperator(""),
49 m_polyDeg(-1),
50 m_noNodes(-1),
51 m_polyType(-1),
52 m_intMethod(-1),
53 m_nodes(),
54 m_wInt(),
55 m_wBary() {
56 TRACE();
57
58 init(polyDeg, polyType, noNodes1D, intMethod, sbpMode, sbpOperator);
59}
void init(const MInt polyDeg, const DgPolynomialType polyType, const MInt noNodes, const DgIntegrationMethod intMethod, const MBool sbpMode, const MString sbpOperator)
Sets the member variables and calls the appropriate functions to calculate the nodes and weights etc.

◆ ~DgInterpolation()

DgInterpolation::~DgInterpolation ( )
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2012-11-07

Definition at line 68 of file dgcartesianinterpolation.cpp.

68 {
69 TRACE();
70 m_sbpMode = false;
71 m_sbpOperator = "";
72 m_polyDeg = -1;
73 m_polyType = -1;
74 m_intMethod = -1;
75 m_nodes.clear();
76 m_wInt.clear();
77 m_wBary.clear();
78 m_Dhat.clear();
79 m_LFace[0].clear();
80 m_LFace[1].clear();
81 m_LhatFace[0].clear();
82 m_LhatFace[1].clear();
83}
MFloatVector m_LhatFace[2]
void clear()
Deletes the data (if internal storage was used) and resets dimensions to zero.
Definition: tensor.h:283

Member Function Documentation

◆ init()

void DgInterpolation::init ( const MInt  polyDeg,
const DgPolynomialType  polyType,
const MInt  noNodes,
const DgIntegrationMethod  intMethod,
const MBool  sbpMode,
const MString  sbpOperator 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2012-11-07
Parameters
[in]polyDegMaximum polynomial degree.
[in]polyTypePolynomial type.
[in]intMethodIntegration method.
[in]noNodesNumber of nodes
[in]sbpModeSBP mode
[in]sbpOperatorUsed SBP operator

Definition at line 100 of file dgcartesianinterpolation.cpp.

101 {
102 TRACE();
103
104 // Set member variables
105 m_polyDeg = polyDeg;
106 m_noNodes = noNodes;
107 m_polyType = polyType;
108 m_intMethod = intMethod;
109 m_sbpMode = sbpMode;
110 m_sbpOperator = sbpOperator;
111
112 if(!m_sbpMode) {
113 ASSERT(polyDeg + 1 == noNodes, "polyDeg+1 not equal to noNodes in DG mode! (" << polyDeg << " " << noNodes << ")");
114 }
115
116 m_nodes.resize(noNodes);
117 m_wInt.resize(noNodes);
118 m_wBary.resize(noNodes);
119 m_Dhat.resize(noNodes, noNodes);
120
121 if(m_sbpMode) {
122 // DDRP case: Operators are loaded as full matrix and not as set of coefficients
123 // from which these are constructed from ...
124 if(m_sbpOperator == "go4/DDRP307_N32") {
125 readDDRP(&m_nodes[0], &m_wInt[0], &m_Dhat[0]);
126 } else {
127 // Get Coefficients
128 MFloatVector sbpA, sbpP, sbpQ;
129 // Check if available in sbpoperator.h otherwise try to read from file
130 MBool exists = getSBPOperator(noNodes, sbpOperator, sbpA, sbpP, sbpQ);
131 if(!exists) {
132 readSbpOperator(noNodes, sbpOperator, sbpA, sbpP, sbpQ);
133 }
134
135 // TODO labels:DG,toremove Remove, dg_integration method not needed for SBP
137 mTerm(1, "LG NOT AVAILABLE FOR SBP MODE");
139 calcSBPNodes(noNodes, &m_nodes[0]);
140 calcSBPWeights(noNodes, sbpP, &m_wInt[0]);
141 }
142
143 // Calculate the polynomial derivative matrix D, and normalize it with the
144 // integration weights to get Dhat.
145 calcDhatSBP(noNodes, sbpA, sbpQ, &m_wInt[0], &m_Dhat[0]);
146 }
147 } else {
148 // Go through all possible combinations of the integration method and the
149 // polynomial type, and call the corresponding function to calculate the
150 // integration nodes and weights.
155 }
156
157 // Calculate and set the barycentric weights for Lagrange interpolation
159
160 // Calculate the polynomial derivative matrix D, and normalize it with the
161 // integration weights to get Dhat.
162 calcDhat(m_polyDeg, &m_nodes[0], &m_wInt[0], &m_wBary[0], &m_Dhat[0]);
163 }
164
165 // Calculate the values of the Lagrange polynomials at the boundaries of the
166 // interval [-1,1]
167 m_LFace[0].resize(noNodes);
168 m_LFace[1].resize(noNodes);
169
170 if(m_sbpMode) {
171 calcLinearInterpolationBase(-F1, noNodes, &m_nodes[0], &m_LFace[0][0]);
172 calcLinearInterpolationBase(+F1, noNodes, &m_nodes[0], &m_LFace[1][0]);
173 } else {
176 }
177
178 // Calculate the values of the Lagrange polynomials at the boundaries of the
179 // interval [-1,1], and normalize them with the integration weights.
180 m_LhatFace[0].resize(noNodes);
181 m_LhatFace[1].resize(noNodes);
182
183 if(m_sbpMode) {
184 for(MInt i = 0; i < noNodes; i++) {
185 m_LhatFace[0][i] = m_LFace[0][i] / m_wInt[i];
186 m_LhatFace[1][i] = m_LFace[1][i] / m_wInt[i];
187 }
188 } else {
189 calcLhat(-F1, m_polyDeg, &m_nodes[0], &m_wInt[0], &m_wBary[0], &m_LhatFace[0][0]);
190 calcLhat(+F1, m_polyDeg, &m_nodes[0], &m_wInt[0], &m_wBary[0], &m_LhatFace[1][0]);
191 }
192}
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
@ DG_POLY_LEGENDRE
Definition: enums.h:310
@ DG_INTEGRATE_GAUSS
Definition: enums.h:313
@ DG_INTEGRATE_GAUSS_LOBATTO
Definition: enums.h:313
void mTerm(const MInt errorCode, const MString &location, const MString &message)
Definition: functions.cpp:29
int32_t MInt
Definition: maiatypes.h:62
bool MBool
Definition: maiatypes.h:58
void readSbpOperator(MInt noNodes1D, MString opName, MFloatVector &sbpA, MFloatVector &sbpP, MFloatVector &sbpQ)
Reads all coefficient files to contruct SBP operator.
void calcLegendreGaussLobattoNodesAndWeights(MInt Nmax, MFloat *nodes, MFloat *wInt)
Calculate the Gauss-Lobatto integration nodes and weight for the Legendre polynomials on the interval...
void calcLinearInterpolationBase(const MFloat x, const MInt noNodes, const MFloat *nodes, MFloat *const polynomials)
Calculates linear interpolation base for point (1D).
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,...
void calcSBPNodes(const MInt noNodes, MFloat *nodes)
Generates an equidistant node distribution.
void readDDRP(MFloat *nodesVector, MFloat *wIntVector, MFloat *dhatMatrix)
Reads DDRP operator and nodal distribution from file.
void calcBarycentricWeights(MInt Nmax, const MFloat *nodes, MFloat *weights)
Calculates the barycentric weights for Lagrange interpolation at thei specified nodes.
void calcLhat(const MFloat x, const MInt polyDeg, const MFloat *nodes, const MFloat *wInt, const MFloat *wBary, MFloat *Lhat)
Calculates the Lagrange polynomials evaluated at point x in [-1,1] and pre-divides them by the integr...
void calcDhatSBP(const MInt noNodes, MFloatVector sbpA, MFloatVector sbpQ, const MFloat *wInt, MFloat *dhatMatrix)
Calculates SBP operator from coefficients.
void calcSBPWeights(const MInt noNodes, MFloatVector sbpP, MFloat *wInt)
Calulates the diagonal weight matrix (H) entries.
void calcLegendreGaussNodesAndWeights(MInt Nmax, MFloat *nodes, MFloat *wInt)
Calculate the Gauss integration nodes and weight for the Legendre polynomials on the interval [-1,...
void calcDhat(MInt Nmax, const MFloat *nodes, const MFloat *wInt, const MFloat *wBary, MFloat *dhatMatrix)
Calculates the polynomial derivative matrix and normalizes it using the integration weights.
MBool getSBPOperator(const MInt noNodes1D, const MString opName, MFloatVector &sbpA, MFloatVector &sbpP, MFloatVector &sbpQ)
Gets SBP Coefficients from corresponding header if existent.

◆ initInterpolation()

void DgInterpolation::initInterpolation ( const MInt  polyDeg,
const DgPolynomialType  polyType,
const MInt  noNodes,
const DgIntegrationMethod  intMethod,
const MBool  sbpMode,
const MString  sbpOperator 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2012-11-07
Parameters
[in]polyDegMaximum polynomial degree.
[in]polyTypePolynomial type.
[in]intMethodIntegration method.
[in]noNodesNumber of nodes
[in]sbpModeSBP mode
[in]sbpOperatorUsed SBP operator

Definition at line 207 of file dgcartesianinterpolation.cpp.

209 {
210 TRACE();
211
212 if(sbpMode) {
213 // Set member variables
214 m_polyDeg = polyDeg;
215 m_noNodes = noNodes;
216 // m_polyType = polyType;
217 // m_intMethod = intMethod;
218 // m_sbpMode = sbpMode;
219 // m_sbpOperator = sbpOperator;
220
221 m_nodes.resize(noNodes);
222 // m_wInt.resize(noNodes);
223 // m_wBary.resize(noNodes);
224 // m_Dhat.resize(noNodes, noNodes);
225
226 calcSBPNodes(noNodes, &m_nodes[0]);
227 } else {
228 init(polyDeg, polyType, noNodes, intMethod, sbpMode, sbpOperator);
229 }
230}

Member Data Documentation

◆ m_Dhat

MFloatMatrix DgInterpolation::m_Dhat

Definition at line 70 of file dgcartesianinterpolation.h.

◆ m_intMethod

MInt DgInterpolation::m_intMethod

Definition at line 62 of file dgcartesianinterpolation.h.

◆ m_LFace

MFloatVector DgInterpolation::m_LFace[2]

Definition at line 72 of file dgcartesianinterpolation.h.

◆ m_LhatFace

MFloatVector DgInterpolation::m_LhatFace[2]

Definition at line 75 of file dgcartesianinterpolation.h.

◆ m_nodes

MFloatVector DgInterpolation::m_nodes

Definition at line 64 of file dgcartesianinterpolation.h.

◆ m_noNodes

MInt DgInterpolation::m_noNodes

Definition at line 58 of file dgcartesianinterpolation.h.

◆ m_polyDeg

MInt DgInterpolation::m_polyDeg

Definition at line 56 of file dgcartesianinterpolation.h.

◆ m_polyType

MInt DgInterpolation::m_polyType

Definition at line 60 of file dgcartesianinterpolation.h.

◆ m_sbpMode

MBool DgInterpolation::m_sbpMode

Definition at line 52 of file dgcartesianinterpolation.h.

◆ m_sbpOperator

MString DgInterpolation::m_sbpOperator

Definition at line 54 of file dgcartesianinterpolation.h.

◆ m_wBary

MFloatVector DgInterpolation::m_wBary

Definition at line 68 of file dgcartesianinterpolation.h.

◆ m_wInt

MFloatVector DgInterpolation::m_wInt

Definition at line 66 of file dgcartesianinterpolation.h.


The documentation for this class was generated from the following files: