MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn > Class Template Referencefinal

#include <dgcartesianbcacousticperturbstraightductexit.h>

Inheritance diagram for DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >:
[legend]
Collaboration diagram for DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >:
[legend]

Public Types

using Base = DgBoundaryCondition< nDim, SysEqn >
 
- Public Types inherited from DgBoundaryCondition< nDim, SysEqn >
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

 DgBcAcousticPerturbStraightDuctExit (SolverType &solver_, MInt bcId)
 
MString name () const override
 Returns name of boundary condition. More...
 
void init () override
 
void apply (const MFloat time) override
 Apply method to apply boundary condition. More...
 
void applyAtSurface (const MInt surfaceId, const MFloat time)
 
- Public Member Functions inherited from DgBoundaryCondition< nDim, SysEqn >
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))
 

Private Attributes

DgBcAcousticPerturbSolidWall< nDim, SysEqnm_wallBc
 
MFloat m_ductPosition [2] {}
 

Additional Inherited Members

- Protected Member Functions inherited from DgBoundaryCondition< nDim, SysEqn >
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...
 
void calcVolumeIntegral (const MInt noElements, ElementCollector &elem, F &fluxFct)
 
void resetBuffer (const MInt totalSize, MFloat *const buffer)
 
void applyJacobian (const MInt noElements, ElementCollector &elem)
 
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)
 
void calcRegularSurfaceFlux (const MInt begin_, const MInt end_, SurfaceCollector &surf, F &riemannFct)
 

Detailed Description

template<MInt nDim, class SysEqn>
class DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >

Definition at line 20 of file dgcartesianbcacousticperturbstraightductexit.h.

Member Typedef Documentation

◆ Base

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

Constructor & Destructor Documentation

◆ DgBcAcousticPerturbStraightDuctExit()

template<MInt nDim, class SysEqn >
DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >::DgBcAcousticPerturbStraightDuctExit ( SolverType solver_,
MInt  bcId 
)
inline

Definition at line 33 of file dgcartesianbcacousticperturbstraightductexit.h.

33: Base(solver_, bcId), m_wallBc(solver_, bcId) {}
DgBcAcousticPerturbSolidWall< nDim, SysEqn > m_wallBc

Member Function Documentation

◆ apply()

template<MInt nDim, class SysEqn >
void DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >::apply ( const MFloat  time)
inlineoverridevirtual

Implements DgBoundaryCondition< nDim, SysEqn >.

Definition at line 38 of file dgcartesianbcacousticperturbstraightductexit.h.

38 {
40 }
void applyAtSurface(const MInt surfaceId, const MFloat time)
MInt begin() const
Return index of first surface.
MInt end() const
Return index of one-past-last surface.
void loop(Functor &&fun, Class *object, const IdType begin, const IdType end, Args... args)

◆ applyAtSurface()

template<MInt nDim, class SysEqn >
void DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >::applyAtSurface ( const MInt  surfaceId,
const MFloat  time 
)

Definition at line 81 of file dgcartesianbcacousticperturbstraightductexit.h.

81 {
82 const MFloat* const coordinates = &surfaces().coords(surfaceId, 0);
83 const MInt sOrientation = surfaces().orientation(surfaceId);
84
85 const MFloat pos = coordinates[(sOrientation + 1) % 2];
86 const MFloat lb = m_ductPosition[0];
87 const MFloat ub = m_ductPosition[1];
88
89 // Check position of surface center
90 if(lb > pos || pos > ub) {
91 // Apply wall boundary condition
92 m_wallBc.applyAtSurface(surfaceId, time);
93 } else {
94 // Calculate flux with the inner state of a surface and boundary state
95 const MInt noNodes1D = surfaces().noNodes1D(surfaceId);
96
97 const MFloat* surfaceStateL = &surfaces().variables(surfaceId, 0);
98 const MFloat* surfaceStateR = &surfaces().variables(surfaceId, 1);
99 const MFloat* outerState = (surfaces().nghbrElementIds(surfaceId, 0) == -1) ? surfaceStateL : surfaceStateR;
100 MFloatTensor oS(const_cast<MFloat*>(outerState), noNodes1D, SysEqn::noVars());
101
102 // Check the orientation of the Duct-exit and set the outer state
103 if(sOrientation == 1) {
104 for(MInt n = 0; n < noNodes1D; n++) {
105 // for y-orientation x-velocity is zero
106 oS(n, 0) = 0.0;
107 oS(n, 1) = sin(2 * PI * time);
108 oS(n, 2) = sin(2 * PI * time);
109 }
110 } else {
111 for(MInt n = 0.0; n < noNodes1D; n++) {
112 // for x orientation y-velocity is zero
113 oS(n, 0) = sin(2 * PI * time);
114 oS(n, 1) = 0;
115 oS(n, 2) = sin(2 * PI * time);
116 }
117 }
118
119 // Call regular Riemann solver
120 const MFloat* nodeVarsL = &surfaces().nodeVars(surfaceId, 0);
121 const MFloat* nodeVarsR = &surfaces().nodeVars(surfaceId, 1);
122 const MFloat* innerNodeVars = (surfaces().nghbrElementIds(surfaceId, 0) == -1) ? nodeVarsR : nodeVarsL;
123 sysEqn().calcRiemann(innerNodeVars, innerNodeVars, surfaceStateL, surfaceStateR, noNodes1D, sOrientation,
124 flux(surfaceId));
125 }
126}
SurfaceCollector & surfaces()
Return reference to surfaces.
SysEqn & sysEqn()
Return reference to SysEqn object.
MFloat * flux(const MInt i)
Return pointer to surface flux.
MInt & noNodes1D(const MInt srfcId)
Accessor for number of nodes 1D.
MFloat & nodeVars(const MInt srfcId, const MInt side)
Accessor for node variables.
MInt & nghbrElementIds(const MInt srfcId, const MInt side)
Accessor for neighbor element ids.
MFloat & variables(const MInt srfcId, const MInt side)
Accessor for variables.
MInt & orientation(const MInt srfcId)
Accessor for orientation.
MFloat & coords(const MInt srfcId, const MInt dir)
Accessor for coordinates.
int32_t MInt
Definition: maiatypes.h:62
double MFloat
Definition: maiatypes.h:52

◆ init()

template<MInt nDim, class SysEqn >
void DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >::init
overridevirtual

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 from DgBoundaryCondition< nDim, SysEqn >.

Definition at line 62 of file dgcartesianbcacousticperturbstraightductexit.h.

◆ name()

template<MInt nDim, class SysEqn >
MString DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >::name ( ) const
inlineoverridevirtual

Implements DgBoundaryCondition< nDim, SysEqn >.

Definition at line 34 of file dgcartesianbcacousticperturbstraightductexit.h.

34{ return "Straight duct"; }

Member Data Documentation

◆ m_ductPosition

template<MInt nDim, class SysEqn >
MFloat DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >::m_ductPosition[2] {}
private

◆ m_wallBc

template<MInt nDim, class SysEqn >
DgBcAcousticPerturbSolidWall<nDim, SysEqn> DgBcAcousticPerturbStraightDuctExit< nDim, SysEqn >::m_wallBc
private

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