MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
dgcartesianboundarycondition.h
Go to the documentation of this file.
1// Copyright (C) 2024 The m-AIA AUTHORS
2//
3// This file is part of m-AIA (https://git.rwth-aachen.de/aia/m-AIA/m-AIA)
4//
5// SPDX-License-Identifier: LGPL-3.0-only
6
7#ifndef DGBOUNDARYCONDITION_H_
8#define DGBOUNDARYCONDITION_H_
9
10#include "INCLUDE/maiatypes.h"
14#include "dgcartesiansolver.h"
16
17template <MInt nDim, class SysEqn>
19
20template <MInt nDim, class SysEqn_>
22 // Typedefs
23 public:
24 using SysEqn = SysEqn_;
29
30 // Interface methods
31 public:
34
36 void init(const MInt begin_, const MInt end_) {
37 m_begin = begin_;
38 m_end = end_;
39 init();
40 }
41
43 virtual void apply(const MFloat time) = 0;
44
46 virtual MString name() const = 0;
47
49 DgBoundaryCondition(SolverType& solver_, MInt bcId) : m_solver(solver_), m_bcId(bcId) {}
50
52 MInt id() const { return m_bcId; }
53
55 MInt begin() const { return m_begin; }
56
58 MInt end() const { return m_end; }
59
61 MInt count() const { return end() - begin(); }
62
66 virtual MInt noRestartVars() const { return 0; }
67
69 virtual MInt getLocalNoNodes() const { return 0; }
70
72 virtual MString restartVarName(const MInt NotUsed(id)) const { TERMM(1, "Not implemented in derived class!"); }
73
75 virtual void setRestartVariable(const MInt NotUsed(id), const MFloat* const NotUsed(data)) {
76 TERMM(1, "Not implemented in derived class!");
77 }
78
80 virtual void getRestartVariable(const MInt NotUsed(id), MFloat* const NotUsed(data)) const {
81 TERMM(1, "Not implemented in derived class!");
82 }
83
84 virtual MInt noBcElements() const { return 0; }
85 virtual MBool hasBcElement(const MInt NotUsed(elementId)) const { return false; }
86
87 // DLB methods
88 virtual MInt noCellDataDlb() const { return 0; }
89 virtual MInt cellDataTypeDlb(const MInt NotUsed(dataId)) const { return -1; }
90 virtual MInt cellDataSizeDlb(const MInt NotUsed(dataId), const MInt NotUsed(cellId)) const { return 0; }
91 virtual void getCellDataDlb(const MInt NotUsed(dataId), MFloat* const NotUsed(data)) const {
92 TERMM(1, "Not implemented in derived class!");
93 }
94 virtual void getCellDataDlb(const MInt NotUsed(dataId), MInt* const NotUsed(data)) const {
95 TERMM(1, "Not implemented in derived class!");
96 }
97 virtual void setCellDataDlb(const MInt NotUsed(dataId), const MFloat* const NotUsed(data)) {
98 TERMM(1, "Not implemented in derived class!");
99 }
100 virtual void setCellDataDlb(const MInt NotUsed(dataId), const MInt* const NotUsed(data)) {
101 TERMM(1, "Not implemented in derived class!");
102 }
103
104 // Methods for derived classes
105 protected:
108
111
113 MFloat* flux(const MInt i) { return &surfaces().flux(i); }
114
117
119 MInt getElementByCellId(const MInt cellId) const { return m_solver.m_elements.getElementByCellId(cellId); }
120
123
126
128 MBool isMpiSurface(const MInt id_) const { return m_solver.isMpiSurface(id_); };
129
131 MBool needHElementForCell(const MInt cellId) { return m_solver.needHElementForCell(cellId); };
132
134 MInt getHElementId(const MInt elementId) { return m_solver.getHElementId(elementId); };
135
137 const DgInterpolation& interpolation(const MInt polyDeg, const MInt noNodes1D) const {
138 return m_solver.m_interpolation[polyDeg][noNodes1D];
139 }
140
143
146
149
152
154 MFloat dt() const { return m_solver.m_dt; }
155
157 MBool isRestart() const { return m_solver.m_restart; }
158
160 void subTimeStepRk(const MFloat dt_, const MInt stage, const MInt totalSize, const MFloat* const rhs,
161 MFloat* const variables, MFloat* const timeIntStorage) {
162 m_solver.subTimeStepRk(dt_, stage, totalSize, rhs, variables, timeIntStorage);
163 }
164
165 // Access to DG-operator methods
166
167 template <class F>
168 void calcVolumeIntegral(const MInt noElements, ElementCollector& elem, F& fluxFct) {
169 m_solver.calcVolumeIntegral(noElements, elem, fluxFct);
170 }
171
172 void resetBuffer(const MInt totalSize, MFloat* const buffer) { m_solver.resetBuffer(totalSize, buffer); }
173
174 void applyJacobian(const MInt noElements, ElementCollector& elem) { m_solver.applyJacobian(noElements, elem); }
175
176 template <class F>
177 void calcSourceTerms(const MFloat t, const MInt noElements, ElementCollector& elem, F& sourceFct) {
178 m_solver.calcSourceTerms(t, noElements, elem, sourceFct);
179 }
180
181 void calcSurfaceIntegral(const MInt begin_, const MInt end_, ElementCollector& elem, SurfaceCollector& surf,
182 HElementCollector& helem, const MInt noHElements) {
183 m_solver.calcSurfaceIntegral(begin_, end_, elem, surf, helem, noHElements);
184 }
185
186 template <class F>
187 void calcRegularSurfaceFlux(const MInt begin_, const MInt end_, SurfaceCollector& surf, F& riemannFct) {
188 m_solver.calcRegularSurfaceFlux(begin_, end_, surf, riemannFct);
189 }
190
191 // Methods that should not be called from anywhere else
192 private:
197 virtual void init() {}
198
199 private:
202
205
208
211};
212
213
214namespace maia {
215namespace dg {
216namespace bc {
217
218template <class Functor, class Class, class IdType, class... Args>
219void loop(Functor&& fun, Class* object, const IdType begin, const IdType end, Args... args) {
220#ifdef _OPENMP
221#pragma omp parallel for
222#endif
223 for(IdType i = begin; i < end; i++) {
224 (object->*fun)(i, std::forward<Args>(args)...);
225 }
226}
227
228} // namespace bc
229} // namespace dg
230} // namespace maia
231
232#endif // DGBOUNDARYCONDITION_H_
const DgInterpolation & interpolation(const MInt polyDeg, const MInt noNodes1D) const
Return interpolation.
virtual void apply(const MFloat time)=0
Apply method to apply boundary condition.
void calcSourceTerms(const MFloat t, const MInt noElements, ElementCollector &elem, F &sourceFct)
MInt id() const
Return boundary condition if of this boundary condition.
MInt m_end
Id of one-past-last surface of this boundary condition.
virtual void setCellDataDlb(const MInt NotUsed(dataId), const MInt *const NotUsed(data))
MInt m_begin
Id of first surface of this boundary condition.
virtual MInt noCellDataDlb() const
virtual void getRestartVariable(const MInt NotUsed(id), MFloat *const NotUsed(data)) const
Copy restart variable data from boundary condition class to pointer.
SurfaceCollector & surfaces()
Return reference to surfaces.
virtual void setCellDataDlb(const MInt NotUsed(dataId), const MFloat *const NotUsed(data))
virtual MInt noBcElements() const
MInt timeIntegrationScheme() const
Return time integration scheme.
MBool needHElementForCell(const MInt cellId)
Return if h-element is needed for given cell.
ElementCollector & elements()
Return reference to elements.
virtual MString name() const =0
Returns name of boundary condition.
virtual void setRestartVariable(const MInt NotUsed(id), const MFloat *const NotUsed(data))
Copy restart variable data from pointer to boundary condition class.
MBool isRestart() const
Return if a restart is performed.
virtual ~DgBoundaryCondition()
Destructor must be virtual.
MInt getElementByCellId(const MInt cellId) const
Return element id corresponding to given cell id.
void calcSurfaceIntegral(const MInt begin_, const MInt end_, ElementCollector &elem, SurfaceCollector &surf, HElementCollector &helem, const MInt noHElements)
virtual void getCellDataDlb(const MInt NotUsed(dataId), MInt *const NotUsed(data)) const
void applyJacobian(const MInt noElements, ElementCollector &elem)
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.
SysEqn & sysEqn()
Return reference to SysEqn object.
SolverType & solver()
Return reference to solver.
HElementCollector & helements()
Return reference to h-elements.
MInt getHElementId(const MInt elementId)
Return h-element id for an element.
MInt maxPolyDeg() const
Return maximum polynomial degree.
MBool isMpiSurface(const MInt id_) const
Return true if surface is a MPI surface.
MFloat dt() const
Return current time step size.
virtual MInt getLocalNoNodes() const
Return local number of nodes.
void calcVolumeIntegral(const MInt noElements, ElementCollector &elem, F &fluxFct)
virtual MInt cellDataSizeDlb(const MInt NotUsed(dataId), const MInt NotUsed(cellId)) const
virtual void getCellDataDlb(const MInt NotUsed(dataId), MFloat *const NotUsed(data)) const
MInt begin() const
Return index of first surface.
SolverType & m_solver
Store a reference to the solver.
MInt count() const
Return number of boundary surfaces.
void resetBuffer(const MInt totalSize, MFloat *const buffer)
MInt end() const
Return index of one-past-last surface.
virtual MBool hasBcElement(const MInt NotUsed(elementId)) const
void init(const MInt begin_, const MInt end_)
Init method to initialize boundary condition for range of surfaces.
MInt integrationMethod() const
Return integration method.
const MInt m_bcId
The boundary condition id of this boundary condition.
DgBoundaryCondition(SolverType &solver_, MInt bcId)
Constructor saves arguments to member variables.
void calcRegularSurfaceFlux(const MInt begin_, const MInt end_, SurfaceCollector &surf, F &riemannFct)
virtual MString restartVarName(const MInt NotUsed(id)) const
Return name of restart variable.
MInt maxNoNodes1D() const
Return maximum number of nodes.
virtual MInt noRestartVars() const
MFloat * flux(const MInt i)
Return pointer to surface flux.
virtual MInt cellDataTypeDlb(const MInt NotUsed(dataId)) const
MBool needHElementForCell(const MInt cellId)
Return true if h-element is needed for cell, false otherwise.
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.
void resetBuffer(const MInt totalSize, MFloat *const buffer)
Reset the given buffer to zero.
void calcSurfaceIntegral()
Calculate the surface integral for all faces of element and update dU/dt.
void applyJacobian()
Adds the negative of the inverse Jacobian to the time derivative.
MInt getHElementId(const MInt elementId) const
Return h-element id of a given element (if it exists).
std::vector< std::vector< DgInterpolation > > m_interpolation
MBool isMpiSurface(const MInt id) const
Return true if a surface is a MPI surface.
void calcSourceTerms(MFloat t)
Calculates the source terms for each node and adds them to the time derivative of the conservative va...
void calcVolumeIntegral()
Calculate the volume integral for all elements and update m_rightHandSide.
SurfaceCollector m_surfaces
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.
ElementCollector m_elements
HElementCollector m_helements
Class stores precalculated values for interpolation & integration on the reference interval [-1,...
MBool m_restart
Definition: solver.h:97
Class that represents DG element collector.
MInt getElementByCellId(const MInt cellId) const
Return element id for a given cell id (or -1 if not found).
Class that represents DG element collector.
Class that represents DG element collector.
MFloat & flux(const MInt srfcId)
Accessor for flux.
int32_t MInt
Definition: maiatypes.h:62
std::basic_string< char > MString
Definition: maiatypes.h:55
double MFloat
Definition: maiatypes.h:52
bool MBool
Definition: maiatypes.h:58
void loop(Functor &&fun, Class *object, const IdType begin, const IdType end, Args... args)
Namespace for auxiliary functions/classes.