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

Solid (slip) wall boundary condition. More...

#include <dgcartesianbcacousticperturbsolidwall.h>

Inheritance diagram for DgBcAcousticPerturbSolidWall< nDim, SysEqn, slipWall >:
[legend]
Collaboration diagram for DgBcAcousticPerturbSolidWall< nDim, SysEqn, slipWall >:
[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

 DgBcAcousticPerturbSolidWall (SolverType &solver_, MInt bcId)
 
MString name () const override
 Returns name of boundary condition. More...
 
void apply (const MFloat time) override
 Apply method to apply boundary condition. More...
 
void applyAtSurface (const MInt surfaceId, const MFloat NotUsed(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))
 

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, MBool slipWall = false>
class DgBcAcousticPerturbSolidWall< nDim, SysEqn, slipWall >
Author
Lev Liberson, Hsun-Jen Cheng
Date
2014

Definition at line 17 of file dgcartesianbcacousticperturbsolidwall.h.

Member Typedef Documentation

◆ Base

template<MInt nDim, class SysEqn , MBool slipWall = false>
using DgBcAcousticPerturbSolidWall< nDim, SysEqn, slipWall >::Base = DgBoundaryCondition<nDim, SysEqn>

Definition at line 20 of file dgcartesianbcacousticperturbsolidwall.h.

Constructor & Destructor Documentation

◆ DgBcAcousticPerturbSolidWall()

template<MInt nDim, class SysEqn , MBool slipWall = false>
DgBcAcousticPerturbSolidWall< nDim, SysEqn, slipWall >::DgBcAcousticPerturbSolidWall ( SolverType solver_,
MInt  bcId 
)
inline

Definition at line 29 of file dgcartesianbcacousticperturbsolidwall.h.

29: Base(solver_, bcId) {}
DgBoundaryCondition< nDim, SysEqn > Base

Member Function Documentation

◆ apply()

template<MInt nDim, class SysEqn , MBool slipWall = false>
void DgBcAcousticPerturbSolidWall< nDim, SysEqn, slipWall >::apply ( const MFloat  time)
inlineoverridevirtual

Implements DgBoundaryCondition< nDim, SysEqn >.

Definition at line 34 of file dgcartesianbcacousticperturbsolidwall.h.

34 {
36 }
void applyAtSurface(const MInt surfaceId, const MFloat NotUsed(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 , MBool slipWall>
void DgBcAcousticPerturbSolidWall< nDim, SysEqn, slipWall >::applyAtSurface ( const MInt  surfaceId,
const MFloat   NotUsedtime 
)

Definition at line 43 of file dgcartesianbcacousticperturbsolidwall.h.

44 {
45 // Storing the surface-information and define the inner state
46 MFloat* stateL = &surfaces().variables(surfaceId, 0);
47 MFloat* stateR = &surfaces().variables(surfaceId, 1);
48 MFloat* innerState = (surfaces().internalSideId(surfaceId) == 0) ? stateL : stateR;
49 MFloat* nodeVarsL = &surfaces().nodeVars(surfaceId, 0);
50 MFloat* nodeVarsR = &surfaces().nodeVars(surfaceId, 1);
51 MFloat* innerNodeVars = (surfaces().nghbrElementIds(surfaceId, 0) == -1) ? nodeVarsR : nodeVarsL;
52
53 // Collect data for flux-computation
54 const MInt noNodes1D = surfaces().noNodes1D(surfaceId);
55 const MInt noNodes1D3 = (nDim == 3) ? noNodes1D : 1;
56 const MInt dirId = surfaces().orientation(surfaceId);
57
58 MFloatTensor variables(innerState, noNodes1D, noNodes1D3, SysEqn::noVars());
59 MFloatTensor nodeVars(innerNodeVars, noNodes1D, noNodes1D3, SysEqn::noNodeVars());
60 MFloatTensor f(flux(surfaceId), noNodes1D, noNodes1D3, SysEqn::noVars());
61
62 // Computing the flux at all nodes considering zero-velocities at the wall
63 //
64 // Flux-components for boundary-surfaces (no-slip wall version):
65 // u' v' w' p'
66 // Fx = (p', 0, 0, 0)
67 // Fy = (0, p', 0, 0)
68 // Fz = (0, 0, p', 0)
69 //
70 // Reference: M. Bauer, Airframe noise prediction using a discontinuous
71 // Galerkin method, PhD thesis, p. 18, 2011
72 //
73 // Additionally, if slipWall == true, additional velocity terms are added to
74 // account for non-zero mean flow in the wall-parallel direction.
75 //
76 // Flux-components for boundary-surfaces (slip-wall version):
77 // u' v' w' p'
78 // Fx = (p' + v'v0 + w'w0, 0, 0, 0)
79 // Fy = (0, p' + u'u0 + w'w0, 0, 0)
80 // Fz = (0, 0, p' + u'u0 + v'v0, 0)
81
82 // Set the flux to zero
83 f.set(0.0);
84
85 // Set the x- or y- or z-component of the velocity to p', depending on the
86 // surface orientation
87 if(slipWall) {
88 // If it is a slip wall, add pressure and additional velocity terms
89 for(MInt i = 0; i < noNodes1D; i++) {
90 for(MInt j = 0; j < noNodes1D3; j++) {
91 // This following line is just so that std::inner_product can be used
92 variables(i, j, dirId) = 0.0;
93
94 // Set flux normal to wall (see above for equation)
95 f(i, j, dirId) =
96 std::inner_product(&variables(i, j, 0), &variables(i, j, nDim), &nodeVars(i, j, 0), variables(i, j, nDim));
97 }
98 }
99 } else {
100 // If it is a no-slip wall, just add the pressure term and implicity assume
101 // that all velocity compontents are zero at the wall
102 for(MInt i = 0; i < noNodes1D; i++) {
103 for(MInt j = 0; j < noNodes1D3; j++) {
104 // Set flux normal to wall (see above for equation)
105 f(i, j, dirId) = variables(i, j, nDim);
106 }
107 }
108 }
109}
SurfaceCollector & surfaces()
Return reference to surfaces.
MFloat * flux(const MInt i)
Return pointer to surface flux.
MInt & internalSideId(const MInt srfcId)
Accessor for internal side id.
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.
int32_t MInt
Definition: maiatypes.h:62
double MFloat
Definition: maiatypes.h:52

◆ name()

template<MInt nDim, class SysEqn , MBool slipWall = false>
MString DgBcAcousticPerturbSolidWall< nDim, SysEqn, slipWall >::name ( ) const
inlineoverridevirtual

Implements DgBoundaryCondition< nDim, SysEqn >.

Definition at line 30 of file dgcartesianbcacousticperturbsolidwall.h.

30 {
31 return slipWall ? "solid wall (mean flow slip wall)" : "solid wall (mean flow no-slip wall)";
32 }

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