MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
dgcartesianbcacousticperturbcbc.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 DGBOUNDARYCONDITIONACOUSTICPERTURBCBC_H_
8#define DGBOUNDARYCONDITIONACOUSTICPERTURBCBC_H_
9
11
12template <MInt nDim>
13class DgBcAcousticPerturbCBC final : public DgBoundaryCondition<nDim, DgSysEqnAcousticPerturb<nDim>> {
14 // Typedefs
15 public:
17 using Base::begin;
18 using Base::end;
19 using Base::flux;
20 using Base::id;
21 using Base::surfaces;
22 using Base::sysEqn;
23 using typename Base::SolverType;
24 using typename Base::SysEqn;
25
26 // Methods
27 DgBcAcousticPerturbCBC(SolverType& solver_, MInt bcId) : Base(solver_, bcId) {}
28 MString name() const override { return "CBC"; }
29
30 void apply(const MFloat time) override {
32 }
33
34 void applyAtSurface(const MInt surfaceId, const MFloat NotUsed(time));
35};
36
37
38template <MInt nDim>
39void DgBcAcousticPerturbCBC<nDim>::applyAtSurface(const MInt surfaceId, const MFloat NotUsed(time)) {
40 // Storing the surface-information and define the inner state
41 MFloat* nodeVarsL = &surfaces().nodeVars(surfaceId, 0);
42 MFloat* nodeVarsR = &surfaces().nodeVars(surfaceId, 1);
43 MFloat* stateL = &surfaces().variables(surfaceId, 0);
44 MFloat* stateR = &surfaces().variables(surfaceId, 1);
45 // TODO labels:DG,totest This might not work for all cases
46 MFloat* innerNodeVars = (surfaces().internalSideId(surfaceId) != 1) ? nodeVarsL : nodeVarsR;
47 MFloat* boundaryNodeVars = (surfaces().internalSideId(surfaceId) != 1) ? nodeVarsR : nodeVarsL;
48 MFloat* innerState = (surfaces().nghbrElementIds(surfaceId, 0) == -1) ? stateR : stateL;
49 MFloat* boundaryState = (surfaces().nghbrElementIds(surfaceId, 0) == -1) ? stateL : stateR;
50
51 // Collect data for flux-computation
52 const MInt dirId = surfaces().orientation(surfaceId);
53 const MInt noNodes1D = surfaces().noNodes1D(surfaceId);
54 const MInt noNodesXD = surfaces().noNodesXD(surfaceId);
55 MFloat* f = flux(surfaceId);
56
57 switch(id()) {
58 // CBC OutFlow
59 case 301: {
60 std::fill_n(&boundaryState[0], noNodesXD * Base::SysEqn::noVars(), 0.0);
61 std::fill_n(&boundaryNodeVars[0], noNodesXD * Base::SysEqn::noNodeVars(), 0.0);
62 // Set default node vars of boundary state that are not zero
63 for(MInt i = 0; i < noNodesXD; i++) {
64 MFloat* const nodeVars = &boundaryNodeVars[i * Base::SysEqn::noNodeVars()];
65 // Set default density
66 nodeVars[Base::SysEqn::CV::RHO0] = 1.0;
67 // Set default speed of sound
68 nodeVars[Base::SysEqn::CV::C0] = 1.0;
69 }
70 sysEqn().calcRiemannRoe(nodeVarsL, nodeVarsR, stateL, stateR, noNodes1D, dirId, f);
71 break;
72 }
73
74 // "periodic" extension of domain for 1D test
75 case 304: {
76 sysEqn().calcRiemannRoe(innerNodeVars, innerNodeVars, innerState, innerState, noNodes1D, dirId, f);
77 break;
78 }
79
80 default: {
81 TERMM(1, "Bad boundary condition id for CBC boundary conditions");
82 break;
83 }
84 }
85}
86
87#endif // DGBOUNDARYCONDITIONACOUSTICPERTURBCBC_H_
DgBcAcousticPerturbCBC(SolverType &solver_, MInt bcId)
void apply(const MFloat time) override
Apply method to apply boundary condition.
MString name() const override
Returns name of boundary condition.
void applyAtSurface(const MInt surfaceId, const MFloat NotUsed(time))
MInt id() const
Return boundary condition if of this boundary condition.
SurfaceCollector & surfaces()
Return reference to surfaces.
MInt end() const
Return index of one-past-last surface.
MFloat * flux(const MInt i)
Return pointer to surface flux.
int32_t MInt
Definition: maiatypes.h:62
std::basic_string< char > MString
Definition: maiatypes.h:55
double MFloat
Definition: maiatypes.h:52
void loop(Functor &&fun, Class *object, const IdType begin, const IdType end, Args... args)