DamBreak++ Wave Simulator 0.3
DamBreak++ Simulation Framework
Loading...
Searching...
No Matches
dbpp_HLLegacyFluxAlgorithm.h
Go to the documentation of this file.
1#pragma once
2
3// Nov19 DamBreak includes
4#include "dbpp_FluxTensor.h"
6#include "dbpp_HLLSolver1D.h"
8#include "SimulationConfig.h" // pointer-2-function
9// SfxBase19 include
10#include "include/Sfx_cellFaceVariables.h"
11
12// forward declarations
13namespace dbpp {
14 class cellFace;
16}
17
18namespace dbpp
19{
20 using valarray64_t = std::valarray<float64>;
21 using ff12pair = std::pair<valarray64_t, valarray64_t>;
22
27 {
28 private:
29 std::string m_ptr2fName;
30 int32 m_varOrder;
32 protected:
37 FluxTensorMap calculFF( const std::valarray<float64>& aU1, const std::valarray<float64>& aU2) const;
38 public:
44 Ptr2FLegacyFluxAlgorithm(std::string aFluxAlgoName, CalculFF aPtr2Func);
50 FluxTensor calculFF( const Sfx::cellFaceVariables& aCellFaceVar) override final;
60 calculFF( const Sfx::scalarField1D& U1, const Sfx::scalarField1D& U2,
61 const Omega& aDomain, const PhysicalBoundaryCnd& aPhysbc) override final;
69 FluxTensor calculFF( const Sfx::StateVector& aUL, const Sfx::StateVector& aUR,
70 const cellFace& aCellFace) override final;
78 FluxTensorMap calculFF( const Sfx::scalarField1D& U1, const Sfx::scalarField1D& U2,
79 const PhysicalBoundaryCnd& aPhysbc) override final;
80 // MUSCL as default
81 // eReconstrType getReconstrType() const { return m_reconstrType; }
86 int32 getReconstrVarOrder() const noexcept { return m_varOrder; }
91 bool useReconstr() const noexcept override { return true; }
96 void setReconstrType(eReconstrType aRecnstrType) noexcept override {/*no implemnentation*/ } //ON/OFF??
101 virtual eReconstrType getReconstrType() const noexcept override { return eReconstrType::None; }
106 bool usePhysicalInCompleteFlux() const noexcept override { return true; }
111 // void setNujicAlphaCoeff(float64 aAlphaCoeff) noexcept { m_alphaCoeff = aAlphaCoeff; }
116 // float64 getNujicAlphaCoeff() const noexcept { return m_alphaCoeff; }
121 std::string PtrToFuncName() const noexcept { return m_ptr2fName; }
122 };
123} // End of namespace
void(*)(DBSArrayType &aFF1, DBSArrayType &aFF2, const DBSArrayType &aU1, const DBSArrayType &aU2) CalculFF
PointerToFunction(numerical flux computation)
Definition SimulationConfig.h:32
Abstract base class for numerical flux algorithm.
Definition dbpp_FluxAlgorithm.h:24
Global Domain (part of global discretization) list of elements and geomeric nodes used by numerical m...
Definition dbpp_Omega.h:19
Physical boundary condition (computational domain) based on characteristic equation ....
Definition dbpp_PhysicalBoundaryCnd.h:22
void setReconstrType(eReconstrType aRecnstrType) noexcept override
set reconstruction type
Definition dbpp_HLLegacyFluxAlgorithm.h:96
std::string PtrToFuncName() const noexcept
Nujic coefficient.
Definition dbpp_HLLegacyFluxAlgorithm.h:121
std::string m_ptr2fName
Definition dbpp_HLLegacyFluxAlgorithm.h:29
bool useReconstr() const noexcept override
Reconstruction of state variables at cell face.
Definition dbpp_HLLegacyFluxAlgorithm.h:91
Ptr2FLegacyFluxAlgorithm(std::string aFluxAlgoName, CalculFF aPtr2Func)
support legacy code (flux algorithm scheme pointer-to-function).
Definition dbpp_HLLegacyFluxAlgorithm.cpp:12
int32 getReconstrVarOrder() const noexcept
set reconstruction variable order
Definition dbpp_HLLegacyFluxAlgorithm.h:86
virtual eReconstrType getReconstrType() const noexcept override
MUSCL as default.
Definition dbpp_HLLegacyFluxAlgorithm.h:101
CalculFF m_ptr2func
Definition dbpp_HLLegacyFluxAlgorithm.h:31
int32 m_varOrder
Definition dbpp_HLLegacyFluxAlgorithm.h:30
bool usePhysicalInCompleteFlux() const noexcept override
Physical flux type.
Definition dbpp_HLLegacyFluxAlgorithm.h:106
FluxTensorMap calculFF(const std::valarray< float64 > &aU1, const std::valarray< float64 > &aU2) const
Legacy code support of pointer-to-function flux algorithm. Reconstr procedure of state variables use ...
Definition dbpp_HLLegacyFluxAlgorithm.cpp:20
Cell face in the finite volume discretization Usage: caculFF(const cellFace& aFace) compute the numer...
Definition dbpp_CellFace.h:29
Definition DamBreakProb.h:15
std::valarray< float64 > valarray64_t
Definition dbpp_HLLegacyFluxAlgorithm.h:20
eReconstrType
Definition dbpp_Enumerations.h:33
@ None
Definition dbpp_Enumerations.h:34
std::pair< valarray64_t, valarray64_t > ff12pair
Definition dbpp_HLLegacyFluxAlgorithm.h:21
Flux tensor field (aggregate initialization).
Definition dbpp_FluxTensor.h:21
Map cell face and flux values.
Definition dbpp_FluxTensor.h:52