MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
DgBoundaryCondition< nDim, SysEqn_ > Class Template Referenceabstract

#include <dgcartesianboundarycondition.h>

Inheritance diagram for DgBoundaryCondition< nDim, SysEqn_ >:
[legend]
Collaboration diagram for DgBoundaryCondition< nDim, SysEqn_ >:
[legend]

Public Types

using SysEqn = SysEqn_
 
using SolverType = DgCartesianSolver< nDim, SysEqn >
 
using ElementCollector = maia::dg::collector::ElementCollector< nDim, SysEqn >
 
using HElementCollector = maia::dg::collector::HElementCollector< nDim, SysEqn >
 
using SurfaceCollector = maia::dg::collector::SurfaceCollector< nDim, SysEqn >
 

Public Member Functions

virtual ~DgBoundaryCondition ()
 Destructor must be virtual. More...
 
void init (const MInt begin_, const MInt end_)
 Init method to initialize boundary condition for range of surfaces. More...
 
virtual void apply (const MFloat time)=0
 Apply method to apply boundary condition. More...
 
virtual MString name () const =0
 Returns name of boundary condition. More...
 
 DgBoundaryCondition (SolverType &solver_, MInt bcId)
 Constructor saves arguments to member variables. More...
 
MInt id () const
 Return boundary condition if of this boundary condition. More...
 
MInt begin () const
 Return index of first surface. More...
 
MInt end () const
 Return index of one-past-last surface. More...
 
MInt count () const
 Return number of boundary surfaces. More...
 
virtual MInt noRestartVars () const
 
virtual MInt getLocalNoNodes () const
 Return local number of nodes. More...
 
virtual MString restartVarName (const MInt NotUsed(id)) const
 Return name of restart variable. More...
 
virtual void setRestartVariable (const MInt NotUsed(id), const MFloat *const NotUsed(data))
 Copy restart variable data from pointer to boundary condition class. More...
 
virtual void getRestartVariable (const MInt NotUsed(id), MFloat *const NotUsed(data)) const
 Copy restart variable data from boundary condition class to pointer. More...
 
virtual MInt noBcElements () const
 
virtual MBool hasBcElement (const MInt NotUsed(elementId)) const
 
virtual MInt noCellDataDlb () const
 
virtual MInt cellDataTypeDlb (const MInt NotUsed(dataId)) const
 
virtual MInt cellDataSizeDlb (const MInt NotUsed(dataId), const MInt NotUsed(cellId)) const
 
virtual void getCellDataDlb (const MInt NotUsed(dataId), MFloat *const NotUsed(data)) const
 
virtual void getCellDataDlb (const MInt NotUsed(dataId), MInt *const NotUsed(data)) const
 
virtual void setCellDataDlb (const MInt NotUsed(dataId), const MFloat *const NotUsed(data))
 
virtual void setCellDataDlb (const MInt NotUsed(dataId), const MInt *const NotUsed(data))
 

Protected Member Functions

SolverTypesolver ()
 Return reference to solver. More...
 
SysEqnsysEqn ()
 Return reference to SysEqn object. More...
 
MFloatflux (const MInt i)
 Return pointer to surface flux. More...
 
ElementCollectorelements ()
 Return reference to elements. More...
 
MInt getElementByCellId (const MInt cellId) const
 Return element id corresponding to given cell id. More...
 
HElementCollectorhelements ()
 Return reference to h-elements. More...
 
SurfaceCollectorsurfaces ()
 Return reference to surfaces. More...
 
MBool isMpiSurface (const MInt id_) const
 Return true if surface is a MPI surface. More...
 
MBool needHElementForCell (const MInt cellId)
 Return if h-element is needed for given cell. More...
 
MInt getHElementId (const MInt elementId)
 Return h-element id for an element. More...
 
const DgInterpolationinterpolation (const MInt polyDeg, const MInt noNodes1D) const
 Return interpolation. More...
 
MInt integrationMethod () const
 Return integration method. More...
 
MInt timeIntegrationScheme () const
 Return time integration scheme. More...
 
MInt maxPolyDeg () const
 Return maximum polynomial degree. More...
 
MInt maxNoNodes1D () const
 Return maximum number of nodes. More...
 
MFloat dt () const
 Return current time step size. More...
 
MBool isRestart () const
 Return if a restart is performed. More...
 
void subTimeStepRk (const MFloat dt_, const MInt stage, const MInt totalSize, const MFloat *const rhs, MFloat *const variables, MFloat *const timeIntStorage)
 Access to time integration method. More...
 
template<class F >
void calcVolumeIntegral (const MInt noElements, ElementCollector &elem, F &fluxFct)
 
void resetBuffer (const MInt totalSize, MFloat *const buffer)
 
void applyJacobian (const MInt noElements, ElementCollector &elem)
 
template<class F >
void calcSourceTerms (const MFloat t, const MInt noElements, ElementCollector &elem, F &sourceFct)
 
void calcSurfaceIntegral (const MInt begin_, const MInt end_, ElementCollector &elem, SurfaceCollector &surf, HElementCollector &helem, const MInt noHElements)
 
template<class F >
void calcRegularSurfaceFlux (const MInt begin_, const MInt end_, SurfaceCollector &surf, F &riemannFct)
 

Private Member Functions

virtual void init ()
 

Private Attributes

SolverTypem_solver
 Store a reference to the solver. More...
 
const MInt m_bcId
 The boundary condition id of this boundary condition. More...
 
MInt m_begin
 Id of first surface of this boundary condition. More...
 
MInt m_end
 Id of one-past-last surface of this boundary condition. More...
 

Detailed Description

template<MInt nDim, class SysEqn_>
class DgBoundaryCondition< nDim, SysEqn_ >

Definition at line 21 of file dgcartesianboundarycondition.h.

Member Typedef Documentation

◆ ElementCollector

template<MInt nDim, class SysEqn_ >
using DgBoundaryCondition< nDim, SysEqn_ >::ElementCollector = maia::dg::collector::ElementCollector<nDim, SysEqn>

Definition at line 26 of file dgcartesianboundarycondition.h.

◆ HElementCollector

template<MInt nDim, class SysEqn_ >
using DgBoundaryCondition< nDim, SysEqn_ >::HElementCollector = maia::dg::collector::HElementCollector<nDim, SysEqn>

Definition at line 27 of file dgcartesianboundarycondition.h.

◆ SolverType

template<MInt nDim, class SysEqn_ >
using DgBoundaryCondition< nDim, SysEqn_ >::SolverType = DgCartesianSolver<nDim, SysEqn>

Definition at line 25 of file dgcartesianboundarycondition.h.

◆ SurfaceCollector

template<MInt nDim, class SysEqn_ >
using DgBoundaryCondition< nDim, SysEqn_ >::SurfaceCollector = maia::dg::collector::SurfaceCollector<nDim, SysEqn>

Definition at line 28 of file dgcartesianboundarycondition.h.

◆ SysEqn

template<MInt nDim, class SysEqn_ >
using DgBoundaryCondition< nDim, SysEqn_ >::SysEqn = SysEqn_

Definition at line 24 of file dgcartesianboundarycondition.h.

Constructor & Destructor Documentation

◆ ~DgBoundaryCondition()

template<MInt nDim, class SysEqn_ >
virtual DgBoundaryCondition< nDim, SysEqn_ >::~DgBoundaryCondition ( )
inlinevirtual

Definition at line 33 of file dgcartesianboundarycondition.h.

33{}

◆ DgBoundaryCondition()

template<MInt nDim, class SysEqn_ >
DgBoundaryCondition< nDim, SysEqn_ >::DgBoundaryCondition ( SolverType solver_,
MInt  bcId 
)
inline

Definition at line 49 of file dgcartesianboundarycondition.h.

49: m_solver(solver_), m_bcId(bcId) {}
SolverType & m_solver
Store a reference to the solver.
const MInt m_bcId
The boundary condition id of this boundary condition.

Member Function Documentation

◆ apply()

◆ applyJacobian()

template<MInt nDim, class SysEqn_ >
void DgBoundaryCondition< nDim, SysEqn_ >::applyJacobian ( const MInt  noElements,
ElementCollector elem 
)
inlineprotected

Definition at line 174 of file dgcartesianboundarycondition.h.

174{ m_solver.applyJacobian(noElements, elem); }
void applyJacobian()
Adds the negative of the inverse Jacobian to the time derivative.

◆ begin()

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::begin ( ) const
inline

Definition at line 55 of file dgcartesianboundarycondition.h.

55{ return m_begin; }
MInt m_begin
Id of first surface of this boundary condition.

◆ calcRegularSurfaceFlux()

template<MInt nDim, class SysEqn_ >
template<class F >
void DgBoundaryCondition< nDim, SysEqn_ >::calcRegularSurfaceFlux ( const MInt  begin_,
const MInt  end_,
SurfaceCollector surf,
F &  riemannFct 
)
inlineprotected

Definition at line 187 of file dgcartesianboundarycondition.h.

187 {
188 m_solver.calcRegularSurfaceFlux(begin_, end_, surf, riemannFct);
189 }
void calcRegularSurfaceFlux(const MInt begin, const MInt end, SurfaceCollector &surf, F &riemannFct)
Calculate the numerical flux for a regular (i.e. inner or MPI) surface.

◆ calcSourceTerms()

template<MInt nDim, class SysEqn_ >
template<class F >
void DgBoundaryCondition< nDim, SysEqn_ >::calcSourceTerms ( const MFloat  t,
const MInt  noElements,
ElementCollector elem,
F &  sourceFct 
)
inlineprotected

Definition at line 177 of file dgcartesianboundarycondition.h.

177 {
178 m_solver.calcSourceTerms(t, noElements, elem, sourceFct);
179 }
void calcSourceTerms(MFloat t)
Calculates the source terms for each node and adds them to the time derivative of the conservative va...

◆ calcSurfaceIntegral()

template<MInt nDim, class SysEqn_ >
void DgBoundaryCondition< nDim, SysEqn_ >::calcSurfaceIntegral ( const MInt  begin_,
const MInt  end_,
ElementCollector elem,
SurfaceCollector surf,
HElementCollector helem,
const MInt  noHElements 
)
inlineprotected

Definition at line 181 of file dgcartesianboundarycondition.h.

182 {
183 m_solver.calcSurfaceIntegral(begin_, end_, elem, surf, helem, noHElements);
184 }
void calcSurfaceIntegral()
Calculate the surface integral for all faces of element and update dU/dt.

◆ calcVolumeIntegral()

template<MInt nDim, class SysEqn_ >
template<class F >
void DgBoundaryCondition< nDim, SysEqn_ >::calcVolumeIntegral ( const MInt  noElements,
ElementCollector elem,
F &  fluxFct 
)
inlineprotected

Definition at line 168 of file dgcartesianboundarycondition.h.

168 {
169 m_solver.calcVolumeIntegral(noElements, elem, fluxFct);
170 }
void calcVolumeIntegral()
Calculate the volume integral for all elements and update m_rightHandSide.

◆ cellDataSizeDlb()

template<MInt nDim, class SysEqn_ >
virtual MInt DgBoundaryCondition< nDim, SysEqn_ >::cellDataSizeDlb ( const MInt   NotUseddataId,
const MInt   NotUsedcellId 
) const
inlinevirtual

Definition at line 90 of file dgcartesianboundarycondition.h.

90{ return 0; }

◆ cellDataTypeDlb()

template<MInt nDim, class SysEqn_ >
virtual MInt DgBoundaryCondition< nDim, SysEqn_ >::cellDataTypeDlb ( const MInt   NotUseddataId) const
inlinevirtual

Reimplemented in DgBcAcousticPerturbRBC< nDim >.

Definition at line 89 of file dgcartesianboundarycondition.h.

89{ return -1; }

◆ count()

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::count ( ) const
inline

Definition at line 61 of file dgcartesianboundarycondition.h.

61{ return end() - begin(); }
MInt begin() const
Return index of first surface.
MInt end() const
Return index of one-past-last surface.

◆ dt()

template<MInt nDim, class SysEqn_ >
MFloat DgBoundaryCondition< nDim, SysEqn_ >::dt ( ) const
inlineprotected

Definition at line 154 of file dgcartesianboundarycondition.h.

154{ return m_solver.m_dt; }

◆ elements()

template<MInt nDim, class SysEqn_ >
ElementCollector & DgBoundaryCondition< nDim, SysEqn_ >::elements ( )
inlineprotected

Definition at line 116 of file dgcartesianboundarycondition.h.

116{ return m_solver.m_elements; }
ElementCollector m_elements

◆ end()

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::end ( ) const
inline

Definition at line 58 of file dgcartesianboundarycondition.h.

58{ return m_end; }
MInt m_end
Id of one-past-last surface of this boundary condition.

◆ flux()

template<MInt nDim, class SysEqn_ >
MFloat * DgBoundaryCondition< nDim, SysEqn_ >::flux ( const MInt  i)
inlineprotected

Definition at line 113 of file dgcartesianboundarycondition.h.

113{ return &surfaces().flux(i); }
SurfaceCollector & surfaces()
Return reference to surfaces.
MFloat & flux(const MInt srfcId)
Accessor for flux.

◆ getCellDataDlb() [1/2]

template<MInt nDim, class SysEqn_ >
virtual void DgBoundaryCondition< nDim, SysEqn_ >::getCellDataDlb ( const MInt   NotUseddataId,
MFloat *const   NotUseddata 
) const
inlinevirtual

Definition at line 91 of file dgcartesianboundarycondition.h.

91 {
92 TERMM(1, "Not implemented in derived class!");
93 }

◆ getCellDataDlb() [2/2]

template<MInt nDim, class SysEqn_ >
virtual void DgBoundaryCondition< nDim, SysEqn_ >::getCellDataDlb ( const MInt   NotUseddataId,
MInt *const   NotUseddata 
) const
inlinevirtual

Definition at line 94 of file dgcartesianboundarycondition.h.

94 {
95 TERMM(1, "Not implemented in derived class!");
96 }

◆ getElementByCellId()

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::getElementByCellId ( const MInt  cellId) const
inlineprotected

Definition at line 119 of file dgcartesianboundarycondition.h.

119{ return m_solver.m_elements.getElementByCellId(cellId); }
MInt getElementByCellId(const MInt cellId) const
Return element id for a given cell id (or -1 if not found).

◆ getHElementId()

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::getHElementId ( const MInt  elementId)
inlineprotected

Definition at line 134 of file dgcartesianboundarycondition.h.

134{ return m_solver.getHElementId(elementId); };
MInt getHElementId(const MInt elementId) const
Return h-element id of a given element (if it exists).

◆ getLocalNoNodes()

template<MInt nDim, class SysEqn_ >
virtual MInt DgBoundaryCondition< nDim, SysEqn_ >::getLocalNoNodes ( ) const
inlinevirtual

Reimplemented in DgBcAcousticPerturbRBC< nDim >.

Definition at line 69 of file dgcartesianboundarycondition.h.

69{ return 0; }

◆ getRestartVariable()

template<MInt nDim, class SysEqn_ >
virtual void DgBoundaryCondition< nDim, SysEqn_ >::getRestartVariable ( const MInt   NotUsedid,
MFloat *const   NotUseddata 
) const
inlinevirtual

Definition at line 80 of file dgcartesianboundarycondition.h.

80 {
81 TERMM(1, "Not implemented in derived class!");
82 }

◆ hasBcElement()

template<MInt nDim, class SysEqn_ >
virtual MBool DgBoundaryCondition< nDim, SysEqn_ >::hasBcElement ( const MInt   NotUsedelementId) const
inlinevirtual

Definition at line 85 of file dgcartesianboundarycondition.h.

85{ return false; }

◆ helements()

template<MInt nDim, class SysEqn_ >
HElementCollector & DgBoundaryCondition< nDim, SysEqn_ >::helements ( )
inlineprotected

Definition at line 122 of file dgcartesianboundarycondition.h.

122{ return m_solver.m_helements; }
HElementCollector m_helements

◆ id()

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::id ( ) const
inline

Definition at line 52 of file dgcartesianboundarycondition.h.

52{ return m_bcId; }

◆ init() [1/2]

template<MInt nDim, class SysEqn_ >
virtual void DgBoundaryCondition< nDim, SysEqn_ >::init ( )
inlineprivatevirtual

Init method that may be implemented by the derived classes. This is called by the init(MInt, MInt) method after the begin/end indices are saved and accessible by begin()/end(). If a derived class does not need to initialize anything, it just does not provide an implementation.

Reimplemented in DgBcAcousticPerturbRBC< nDim >, and DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >.

Definition at line 197 of file dgcartesianboundarycondition.h.

197{}

◆ init() [2/2]

template<MInt nDim, class SysEqn_ >
void DgBoundaryCondition< nDim, SysEqn_ >::init ( const MInt  begin_,
const MInt  end_ 
)
inline

Definition at line 36 of file dgcartesianboundarycondition.h.

36 {
37 m_begin = begin_;
38 m_end = end_;
39 init();
40 }

◆ integrationMethod()

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::integrationMethod ( ) const
inlineprotected

Definition at line 142 of file dgcartesianboundarycondition.h.

◆ interpolation()

template<MInt nDim, class SysEqn_ >
const DgInterpolation & DgBoundaryCondition< nDim, SysEqn_ >::interpolation ( const MInt  polyDeg,
const MInt  noNodes1D 
) const
inlineprotected

Definition at line 137 of file dgcartesianboundarycondition.h.

137 {
138 return m_solver.m_interpolation[polyDeg][noNodes1D];
139 }
std::vector< std::vector< DgInterpolation > > m_interpolation

◆ isMpiSurface()

template<MInt nDim, class SysEqn_ >
MBool DgBoundaryCondition< nDim, SysEqn_ >::isMpiSurface ( const MInt  id_) const
inlineprotected

Definition at line 128 of file dgcartesianboundarycondition.h.

128{ return m_solver.isMpiSurface(id_); };
MBool isMpiSurface(const MInt id) const
Return true if a surface is a MPI surface.

◆ isRestart()

template<MInt nDim, class SysEqn_ >
MBool DgBoundaryCondition< nDim, SysEqn_ >::isRestart ( ) const
inlineprotected

Definition at line 157 of file dgcartesianboundarycondition.h.

157{ return m_solver.m_restart; }
MBool m_restart
Definition: solver.h:97

◆ maxNoNodes1D()

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::maxNoNodes1D ( ) const
inlineprotected

Definition at line 151 of file dgcartesianboundarycondition.h.

◆ maxPolyDeg()

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::maxPolyDeg ( ) const
inlineprotected

Definition at line 148 of file dgcartesianboundarycondition.h.

◆ name()

◆ needHElementForCell()

template<MInt nDim, class SysEqn_ >
MBool DgBoundaryCondition< nDim, SysEqn_ >::needHElementForCell ( const MInt  cellId)
inlineprotected

Definition at line 131 of file dgcartesianboundarycondition.h.

131{ return m_solver.needHElementForCell(cellId); };
MBool needHElementForCell(const MInt cellId)
Return true if h-element is needed for cell, false otherwise.

◆ noBcElements()

template<MInt nDim, class SysEqn_ >
virtual MInt DgBoundaryCondition< nDim, SysEqn_ >::noBcElements ( ) const
inlinevirtual

Reimplemented in DgBcAcousticPerturbRBC< nDim >.

Definition at line 84 of file dgcartesianboundarycondition.h.

84{ return 0; }

◆ noCellDataDlb()

template<MInt nDim, class SysEqn_ >
virtual MInt DgBoundaryCondition< nDim, SysEqn_ >::noCellDataDlb ( ) const
inlinevirtual

Reimplemented in DgBcAcousticPerturbRBC< nDim >.

Definition at line 88 of file dgcartesianboundarycondition.h.

88{ return 0; }

◆ noRestartVars()

template<MInt nDim, class SysEqn_ >
virtual MInt DgBoundaryCondition< nDim, SysEqn_ >::noRestartVars ( ) const
inlinevirtual

Return number of restart variables that need to be stored/loaded. If a derived class does not need additional restart variables, it may omit an implementation.

Reimplemented in DgBcAcousticPerturbRBC< nDim >.

Definition at line 66 of file dgcartesianboundarycondition.h.

66{ return 0; }

◆ resetBuffer()

template<MInt nDim, class SysEqn_ >
void DgBoundaryCondition< nDim, SysEqn_ >::resetBuffer ( const MInt  totalSize,
MFloat *const  buffer 
)
inlineprotected

Definition at line 172 of file dgcartesianboundarycondition.h.

172{ m_solver.resetBuffer(totalSize, buffer); }
void resetBuffer(const MInt totalSize, MFloat *const buffer)
Reset the given buffer to zero.

◆ restartVarName()

template<MInt nDim, class SysEqn_ >
virtual MString DgBoundaryCondition< nDim, SysEqn_ >::restartVarName ( const MInt   NotUsedid) const
inlinevirtual

Definition at line 72 of file dgcartesianboundarycondition.h.

72{ TERMM(1, "Not implemented in derived class!"); }

◆ setCellDataDlb() [1/2]

template<MInt nDim, class SysEqn_ >
virtual void DgBoundaryCondition< nDim, SysEqn_ >::setCellDataDlb ( const MInt   NotUseddataId,
const MFloat *const   NotUseddata 
)
inlinevirtual

Definition at line 97 of file dgcartesianboundarycondition.h.

97 {
98 TERMM(1, "Not implemented in derived class!");
99 }

◆ setCellDataDlb() [2/2]

template<MInt nDim, class SysEqn_ >
virtual void DgBoundaryCondition< nDim, SysEqn_ >::setCellDataDlb ( const MInt   NotUseddataId,
const MInt *const   NotUseddata 
)
inlinevirtual

Definition at line 100 of file dgcartesianboundarycondition.h.

100 {
101 TERMM(1, "Not implemented in derived class!");
102 }

◆ setRestartVariable()

template<MInt nDim, class SysEqn_ >
virtual void DgBoundaryCondition< nDim, SysEqn_ >::setRestartVariable ( const MInt   NotUsedid,
const MFloat *const   NotUseddata 
)
inlinevirtual

Definition at line 75 of file dgcartesianboundarycondition.h.

75 {
76 TERMM(1, "Not implemented in derived class!");
77 }

◆ solver()

template<MInt nDim, class SysEqn_ >
SolverType & DgBoundaryCondition< nDim, SysEqn_ >::solver ( )
inlineprotected

Definition at line 107 of file dgcartesianboundarycondition.h.

107{ return m_solver; }

◆ subTimeStepRk()

template<MInt nDim, class SysEqn_ >
void DgBoundaryCondition< nDim, SysEqn_ >::subTimeStepRk ( const MFloat  dt_,
const MInt  stage,
const MInt  totalSize,
const MFloat *const  rhs,
MFloat *const  variables,
MFloat *const  timeIntStorage 
)
inlineprotected

Definition at line 160 of file dgcartesianboundarycondition.h.

161 {
162 m_solver.subTimeStepRk(dt_, stage, totalSize, rhs, variables, timeIntStorage);
163 }
void subTimeStepRk(const MFloat dt, const MInt stage, const MInt totalSize, const MFloat *const rhs, MFloat *const variables, MFloat *const timeIntStorage)
Perform one Runge-Kutta substep on the given elements.

◆ surfaces()

template<MInt nDim, class SysEqn_ >
SurfaceCollector & DgBoundaryCondition< nDim, SysEqn_ >::surfaces ( )
inlineprotected

Definition at line 125 of file dgcartesianboundarycondition.h.

125{ return m_solver.m_surfaces; }
SurfaceCollector m_surfaces

◆ sysEqn()

template<MInt nDim, class SysEqn_ >
SysEqn & DgBoundaryCondition< nDim, SysEqn_ >::sysEqn ( )
inlineprotected

Definition at line 110 of file dgcartesianboundarycondition.h.

110{ return m_solver.m_sysEqn; }

◆ timeIntegrationScheme()

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::timeIntegrationScheme ( ) const
inlineprotected

Definition at line 145 of file dgcartesianboundarycondition.h.

Member Data Documentation

◆ m_bcId

template<MInt nDim, class SysEqn_ >
const MInt DgBoundaryCondition< nDim, SysEqn_ >::m_bcId
private

Definition at line 204 of file dgcartesianboundarycondition.h.

◆ m_begin

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::m_begin
private

Definition at line 207 of file dgcartesianboundarycondition.h.

◆ m_end

template<MInt nDim, class SysEqn_ >
MInt DgBoundaryCondition< nDim, SysEqn_ >::m_end
private

Definition at line 210 of file dgcartesianboundarycondition.h.

◆ m_solver

template<MInt nDim, class SysEqn_ >
SolverType& DgBoundaryCondition< nDim, SysEqn_ >::m_solver
private

Definition at line 201 of file dgcartesianboundarycondition.h.


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