MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
DgSysEqnLinearScalarAdv< nDim > Class Template Reference

#include <dgcartesiansyseqnlinearscalaradv.h>

Inheritance diagram for DgSysEqnLinearScalarAdv< nDim >:
[legend]
Collaboration diagram for DgSysEqnLinearScalarAdv< nDim >:
[legend]

Public Member Functions

 DgSysEqnLinearScalarAdv (MInt solverId)
 Constructor calls parent constructor & loads all necessary properties for this equation. More...
 
void calcInitialCondition (const MFloat t, const MFloat *x, MFloat *nodeVars, MFloat *u) const
 Calculate the initial condition for a certain point in space. More...
 
void calcFlux (const MFloat *nodeVars, const MFloat *u, const MInt noNodes1D, MFloat *flux) const
 Calculates the physical fluxes in all dimensions for all integrations points within a cell. More...
 
void calcSource (const MFloat *nodeVars, const MFloat *u, const MInt noNodes1D, const MFloat t, const MFloat *x, MFloat *src) const
 Calculates the source terms for all integration points within a cell. More...
 
void calcSpongeSource (const MFloat *nodeVars, const MFloat *u, const MInt noNodes1D, const MFloat *eta, MFloat *src) const
 Calculates the sponge source terms for all integration points within a cell. More...
 
MFloat getTimeStep (const MFloat *nodeVars, const MFloat *u, const MInt noNodes1D, const MFloat invJacobian, const MInt sbpMode) const
 Calculate the time step for an explicit time stepping scheme for a given element. More...
 
void calcRiemann (const MFloat *nodeVarsL, const MFloat *nodeVarsR, const MFloat *stateL, const MFloat *stateR, const MInt noNodes1D, const MInt dirId, MFloat *flux) const
 Calculates the numerical flux at a surface given two states (left and right). More...
 
void primToCons (const MFloat *prim, MFloat *cons) const
 Calculates a set of conservative variables from a set of primitive variables. More...
 
void consToPrim (const MFloat *cons, MFloat *prim) const
 Calculates a set of primitive variables from a set of conservative variables. More...
 
MFloat cflFactor (const MInt polyDeg) const
 
void getDefaultNodeVars (MFloat *const NotUsed(nodeVars)) const
 
MBool extendNodeVar (const MInt NotUsed(varId)) const
 
- Public Member Functions inherited from DgSysEqn< nDim, DgSysEqnLinearScalarAdv< nDim > >
MFloat cfl () const
 
MFloat cflScaled (const MInt polyDeg) const
 

Private Attributes

MFloat m_advectionVelocity [nDim] {}
 
MString m_dgIntegrationMethod
 
MInt m_dgTimeIntegrationScheme
 
std::array< std::array< MFloat, s_maxPolyDeg+1 >, 2 > m_cflFactor
 

Static Private Attributes

static const MString s_sysEqnName = "DG_SYSEQN_LINEARSCALARADV"
 
static const MInt s_noVariables = 1
 
static const MString s_consVarNames [s_noVariables] = {"scalar"}
 
static const MString s_primVarNames [s_noVariables] = {"scalar"}
 
static const MInt s_noNodeVars = 0
 
static const MBool s_hasTimeDependentNodeVars = false
 
static constexpr const MStrings_nodeVarNames = nullptr
 
static const MInt s_maxPolyDeg = 31
 

Friends

class DgSysEqn< nDim, DgSysEqnLinearScalarAdv >
 

Additional Inherited Members

- Static Public Member Functions inherited from DgSysEqn< nDim, DgSysEqnLinearScalarAdv< nDim > >
static constexpr MInt noVars ()
 
static constexpr MInt noNodeVars ()
 
static constexpr MBool hasTimeDependentNodeVars ()
 
static const MStringsysEqnName ()
 
static const MStringconsVarNames (MInt i)
 
static const MStringprimVarNames (MInt i)
 
static const MStringnodeVarNames (MInt i)
 
- Public Attributes inherited from DgSysEqn< nDim, DgSysEqnLinearScalarAdv< nDim > >
MInt m_initialCondition
 
MFloat m_initialNumberWaves
 
MInt m_sourceTerm
 
MInt m_riemannSolver
 
- Protected Member Functions inherited from DgSysEqn< nDim, DgSysEqnLinearScalarAdv< nDim > >
 DgSysEqn (MInt solverId)
 
- Protected Attributes inherited from DgSysEqn< nDim, DgSysEqnLinearScalarAdv< nDim > >
const MInt m_solverId
 

Detailed Description

template<MInt nDim>
class DgSysEqnLinearScalarAdv< nDim >

Definition at line 19 of file dgcartesiansyseqnlinearscalaradv.h.

Constructor & Destructor Documentation

◆ DgSysEqnLinearScalarAdv()

template<MInt nDim>
DgSysEqnLinearScalarAdv< nDim >::DgSysEqnLinearScalarAdv ( MInt  solverId)
inlineexplicit
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-18
Parameters
[in]solverIdCurrent solver id.

Definition at line 91 of file dgcartesiansyseqnlinearscalaradv.h.

Member Function Documentation

◆ calcFlux()

template<MInt nDim>
void DgSysEqnLinearScalarAdv< nDim >::calcFlux ( const MFloat nodeVars,
const MFloat u,
const MInt  noNodes1D,
MFloat flux 
) const
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-18
Parameters
[in]uPointer to solution variables (size: noNodes1D^nDim*noVars)).
[in]noNodes1DNumber of nodes 1D in the cell.
[out]fluxCalculated flux (size: noNodes1D^nDim*noVars*nDim)).

Definition at line 311 of file dgcartesiansyseqnlinearscalaradv.h.

312 {
313 const MInt noNodes1D3 = (nDim == 3) ? noNodes1D : 1;
314 const MFloatTensor U(const_cast<MFloat*>(u), noNodes1D, noNodes1D, noNodes1D3);
315 MFloatTensor f(flux, noNodes1D, noNodes1D, noNodes1D3, nDim);
316
317 for(MInt i = 0; i < noNodes1D; i++) {
318 for(MInt j = 0; j < noNodes1D; j++) {
319 for(MInt k = 0; k < noNodes1D3; k++) {
320 for(MInt l = 0; l < nDim; l++) {
321 f(i, j, k, l) = U(i, j, k) * m_advectionVelocity[l];
322 }
323 }
324 }
325 }
326}
int32_t MInt
Definition: maiatypes.h:62
double MFloat
Definition: maiatypes.h:52

◆ calcInitialCondition()

template<MInt nDim>
void DgSysEqnLinearScalarAdv< nDim >::calcInitialCondition ( const MFloat  t,
const MFloat x,
MFloat nodeVars,
MFloat u 
) const
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2012-11-28
Parameters
[in]tTime at which the initial condition should be applied.
[in]xA pointer to the node coordinates.
[out]uA pointer to the solution storage.

Definition at line 237 of file dgcartesiansyseqnlinearscalaradv.h.

◆ calcRiemann()

template<MInt nDim>
void DgSysEqnLinearScalarAdv< nDim >::calcRiemann ( const MFloat nodeVarsL,
const MFloat nodeVarsR,
const MFloat stateL,
const MFloat stateR,
const MInt  noNodes1D,
const MInt  dirId,
MFloat flux 
) const
inline
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-18
Parameters
[in]stateLSurface variables on the left side of the surface (i.e. in -x/-y/-z-direction).
[in]stateRSurface variables on the right side of the surface (i.e. in +x/+y/+z-direction).
[in]noNodes1DNumber of nodes 1D.
[in]dirIdDirection (0 = x-direction, 1 = y-direction, 2 = z-direction).
[out]fluxCalculated Riemann flux (size: noNodes1D^(nDim-1)*noVars).

Definition at line 461 of file dgcartesiansyseqnlinearscalaradv.h.

467 {
468 const MInt noNodes1D3 = (nDim == 3) ? noNodes1D : 1;
469 const MFloatTensor uL(const_cast<MFloat*>(stateL), noNodes1D, noNodes1D3);
470 const MFloatTensor uR(const_cast<MFloat*>(stateR), noNodes1D, noNodes1D3);
471
472 // Calculate maximum eigenvalue
473 const MFloat maxLambda = m_advectionVelocity[dirId];
474
475 // Classical upwind flux...
476 MFloatTensor riemann(flux, noNodes1D, noNodes1D3);
477 for(MInt i = 0; i < noNodes1D; i++) {
478 for(MInt j = 0; j < noNodes1D3; j++) {
479 riemann(i, j) = F1B2 * ((maxLambda + fabs(maxLambda)) * uL(i, j) + (maxLambda - fabs(maxLambda)) * uR(i, j));
480 }
481 }
482}

◆ calcSource()

template<MInt nDim>
void DgSysEqnLinearScalarAdv< nDim >::calcSource ( const MFloat nodeVars,
const MFloat u,
const MInt  noNodes1D,
const MFloat  t,
const MFloat x,
MFloat src 
) const
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-18
Parameters
[in]uPointer to solution variables (size: noNodes1D^nDim*noVars)).
[in]noNodes1DNumber of nodes 1D in the cell.
[in]tCurrent time.
[in]xCoordinates of the integration points of the current cell.
[out]srcPointer to where source terms are stored (size: noNodes1D^nDim*noVars)).

Definition at line 373 of file dgcartesiansyseqnlinearscalaradv.h.

◆ calcSpongeSource()

template<MInt nDim>
void DgSysEqnLinearScalarAdv< nDim >::calcSpongeSource ( const MFloat nodeVars,
const MFloat u,
const MInt  noNodes1D,
const MFloat eta,
MFloat src 
) const
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-18
Parameters
[in]nodeVarsPointer to node variables (size: noNodes1D^nDim*noNodeVars)).
[in]uPointer to solution variables (size: noNodes1D^nDim*noVars)
[in]noNodes1DNumber of nodes 1D in the cell.
[in]tCurrent time.
[in]xCoordinates of the integration points of the current cell.
[out]srcPointer to where source terms are stored (size: noNodes1D^nDim*noVars)).

Definition at line 402 of file dgcartesiansyseqnlinearscalaradv.h.

406 {
407 mTerm(1, AT_, "Sponge: Implementation for linear scalar adv. is missing");
408}
void mTerm(const MInt errorCode, const MString &location, const MString &message)
Definition: functions.cpp:29

◆ cflFactor()

template<MInt nDim>
MFloat DgSysEqnLinearScalarAdv< nDim >::cflFactor ( const MInt  polyDeg) const

Yield an empirically-derived maximum stable CFL number. The user should be able to specify a CFL number of "1.0" in virtually all cases.

Author
Rodrigo Miguez (rodrigo) rodri.nosp@m.go.m.nosp@m.iguez.nosp@m.@rwt.nosp@m.h-aac.nosp@m.hen..nosp@m.de
Date
2018-02-28
Parameters
[in]polyDegPolynomial degree in the cell.
[out]cflfactor obtained from m_cflFactor.

Definition at line 341 of file dgcartesiansyseqnlinearscalaradv.h.

341 {
342 if(polyDeg > s_maxPolyDeg) {
343 TERMM(1, "Polynomial degree exceeds maximum supported (" + std::to_string(s_maxPolyDeg) + ")!");
344 }
345 // Note: The factor "0.95" is used to avoid numerical instabilities when using the full
346 // cfl factor value.
347 const MFloat factor = 0.95 * m_cflFactor[nDim - 2][polyDeg];
348
349 return factor;
350}
std::array< std::array< MFloat, s_maxPolyDeg+1 >, 2 > m_cflFactor

◆ consToPrim()

template<MInt nDim>
void DgSysEqnLinearScalarAdv< nDim >::consToPrim ( const MFloat cons,
MFloat prim 
) const
inline
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-18
Parameters
[in]consConservative variable storage (size: noVars).
[out]primPrimitive variable storage (size: noVars).

Definition at line 496 of file dgcartesiansyseqnlinearscalaradv.h.

496 {
497 // Very easy...
498 std::copy(cons, cons + 1, prim);
499}

◆ extendNodeVar()

template<MInt nDim>
MBool DgSysEqnLinearScalarAdv< nDim >::extendNodeVar ( const MInt   NotUsedvarId) const
inline

Definition at line 43 of file dgcartesiansyseqnlinearscalaradv.h.

43{ return true; };

◆ getDefaultNodeVars()

template<MInt nDim>
void DgSysEqnLinearScalarAdv< nDim >::getDefaultNodeVars ( MFloat *const   NotUsednodeVars) const
inline

Definition at line 42 of file dgcartesiansyseqnlinearscalaradv.h.

42{};

◆ getTimeStep()

template<MInt nDim>
MFloat DgSysEqnLinearScalarAdv< nDim >::getTimeStep ( const MFloat nodeVars,
const MFloat u,
const MInt  noNodes1D,
const MFloat  invJacobian,
const MInt  sbpMode 
) const
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-18
Parameters
[in]uPointer to element variables.
[in]noNodes1DCurrent polynomial degree.
[in]invJacobianInverse Jacobian of element.
Returns
Calculated time step for the element.

Definition at line 423 of file dgcartesiansyseqnlinearscalaradv.h.

425 {
426 // Calculate maximum propagation speed (lambda)
427 MFloat maxLambda = 0;
428 for(MInt n = 0; n < nDim; n++) {
429 maxLambda += fabs(m_advectionVelocity[n]);
430 }
431
432 // Calculate time step for current element
433 MFloat dt;
434 if(sbpMode) {
435 dt = this->cfl() * F2 / (invJacobian * maxLambda * (noNodes1D - 1));
436 } else {
437 const MInt polyDeg = noNodes1D - 1;
438 dt = this->cflScaled(polyDeg) * cflFactor(polyDeg) * F2 / (invJacobian * maxLambda);
439 }
440 return dt;
441}
MFloat cflFactor(const MInt polyDeg) const

◆ primToCons()

template<MInt nDim>
void DgSysEqnLinearScalarAdv< nDim >::primToCons ( const MFloat prim,
MFloat cons 
) const
inline
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-18
Parameters
[in]primPrimitive variable storage (size: noVars).
[out]consConservative variable storage (size: noVars).

Definition at line 513 of file dgcartesiansyseqnlinearscalaradv.h.

513 {
514 // Very easy...
515 std::copy(prim, prim + 1, cons);
516}

Friends And Related Function Documentation

◆ DgSysEqn< nDim, DgSysEqnLinearScalarAdv >

template<MInt nDim>
friend class DgSysEqn< nDim, DgSysEqnLinearScalarAdv >
friend

Definition at line 513 of file dgcartesiansyseqnlinearscalaradv.h.

Member Data Documentation

◆ m_advectionVelocity

template<MInt nDim>
MFloat DgSysEqnLinearScalarAdv< nDim >::m_advectionVelocity[nDim] {}
private

Definition at line 57 of file dgcartesiansyseqnlinearscalaradv.h.

◆ m_cflFactor

template<MInt nDim>
std::array<std::array<MFloat, s_maxPolyDeg + 1>, 2> DgSysEqnLinearScalarAdv< nDim >::m_cflFactor
private

Definition at line 62 of file dgcartesiansyseqnlinearscalaradv.h.

◆ m_dgIntegrationMethod

template<MInt nDim>
MString DgSysEqnLinearScalarAdv< nDim >::m_dgIntegrationMethod
private

Definition at line 59 of file dgcartesiansyseqnlinearscalaradv.h.

◆ m_dgTimeIntegrationScheme

template<MInt nDim>
MInt DgSysEqnLinearScalarAdv< nDim >::m_dgTimeIntegrationScheme
private

Definition at line 60 of file dgcartesiansyseqnlinearscalaradv.h.

◆ s_consVarNames

template<MInt nDim>
const MString DgSysEqnLinearScalarAdv< nDim >::s_consVarNames = {"scalar"}
staticprivate

Definition at line 50 of file dgcartesiansyseqnlinearscalaradv.h.

◆ s_hasTimeDependentNodeVars

template<MInt nDim>
const MBool DgSysEqnLinearScalarAdv< nDim >::s_hasTimeDependentNodeVars = false
staticprivate

Definition at line 54 of file dgcartesiansyseqnlinearscalaradv.h.

◆ s_maxPolyDeg

template<MInt nDim>
const MInt DgSysEqnLinearScalarAdv< nDim >::s_maxPolyDeg = 31
staticprivate

Definition at line 61 of file dgcartesiansyseqnlinearscalaradv.h.

◆ s_nodeVarNames

template<MInt nDim>
constexpr const MString* DgSysEqnLinearScalarAdv< nDim >::s_nodeVarNames = nullptr
staticconstexprprivate

Definition at line 55 of file dgcartesiansyseqnlinearscalaradv.h.

◆ s_noNodeVars

template<MInt nDim>
const MInt DgSysEqnLinearScalarAdv< nDim >::s_noNodeVars = 0
staticprivate

Definition at line 53 of file dgcartesiansyseqnlinearscalaradv.h.

◆ s_noVariables

template<MInt nDim>
const MInt DgSysEqnLinearScalarAdv< nDim >::s_noVariables = 1
staticprivate

Definition at line 49 of file dgcartesiansyseqnlinearscalaradv.h.

◆ s_primVarNames

template<MInt nDim>
const MString DgSysEqnLinearScalarAdv< nDim >::s_primVarNames = {"scalar"}
staticprivate

Definition at line 51 of file dgcartesiansyseqnlinearscalaradv.h.

◆ s_sysEqnName

template<MInt nDim>
const MString DgSysEqnLinearScalarAdv< nDim >::s_sysEqnName = "DG_SYSEQN_LINEARSCALARADV"
staticprivate

Definition at line 47 of file dgcartesiansyseqnlinearscalaradv.h.


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