![]() |
DamBreak++ Wave Simulator 0.3
DamBreak++ Simulation Framework
|
Physics algorithm based on Nujic paper (1995), but use HLL (Harten-Lax-van Leer) for convective flux solver (E. McNeil). More...
#include <dbpp_RhsFluxSrc.h>
Public Member Functions | |
RhsHLLFluxSrc (std::string aName, const PhysicalSystem *aPhysys, bool aUsePressure=false) | |
Constructor from physical system. | |
RhsHLLFluxSrc (std::string aName, HLLFluxAlgorithm *aFluxAlgo, SrcNumericalTreatment *aSrcDiscr, const PhysicalBoundaryCnd &aPhybc) | |
flux algorithm and source discretization | |
void | setPhysicalBoundaryCnd (const PhysicalBoundaryCnd &aPhyBnd) override final |
set physical boundary condition | |
PhysicalBoundaryCnd | getPhysicalBoundaryCnd () const override final |
physical boundary condition in use | |
void | setFluxAlgorithm (FluxAlgorithm *aFluxAlgo, bool useIncompleteFlux=true) override |
Set Numerical flux algorithm. | |
void | setSourceTermDiscr (SrcNumericalTreatment *aSrcDiscr, bool useManningFormula=true) override |
set numerical source treatment (discretization) | |
bool | useReconstruction () const noexcept override final |
variable reconstruction is use | |
bool | isFrictionLess () const noexcept override final |
friction coefficient is applied | |
bool | useSourceTerms () const noexcept override final |
source terms discretization is used | |
SweRhsData | calculate (const Sfx::StateVectorField &aU, const std::shared_ptr< FiniteVolumeDiscretization > &aGblDiscr) override |
calculate rhs terms (apply discretization at each term) | |
SweRhsData | calculate (const Sfx::StateVectorField &aU) override final |
calculate rhs terms (apply discretization at each term) | |
Public Member Functions inherited from dbpp::SweRhsAlgorithm | |
virtual | ~SweRhsAlgorithm ()=default |
Dtor (disable move semantic) If base class has no members, not supporting move semantic has no effect. | |
virtual std::string | name () const |
Class name. | |
virtual void | setPtr2FuncAlgo (CalculFF aPtr2Func) |
Flux algorithm (Pointer to function) | |
virtual void | setPressureTermDiscr () |
Set pressure term discretization. | |
virtual bool | usePressureTerm () const noexcept |
Used pressure terms. | |
virtual bool | useManningFormula () const noexcept |
Used friction formula. | |
virtual bool | useIncompleteFlux () const noexcept |
Used convective. | |
virtual bool | usePtr2FuncLegacy () const noexcept |
Used legacy algorithm. |
Private Attributes | |
std::string | m_name {"RhsHLLFluxSrc"} |
bool | m_usePressure |
PhyBCNdlConstraint | m_leftEnd |
PhyBCNdlConstraint | m_rightEnd |
const PhysicalSystem * | m_phySystem |
PhysicalBoundaryCnd | m_phyBC |
FluxAlgorithm * | m_fluxAlgo {nullptr} |
SrcNumericalTreatment * | m_srcDiscr { nullptr } |
Additional Inherited Members | |
Protected Member Functions inherited from dbpp::SweRhsAlgorithm | |
SweRhsAlgorithm & | operator= (const SweRhsAlgorithm &)=delete |
disable assign operator to avoid slicing problem | |
SweRhsAlgorithm & | operator= (SweRhsAlgorithm &&)=delete |
disable assign operator to avoid slicing problem |
Physics algorithm based on Nujic paper (1995), but use HLL (Harten-Lax-van Leer) for convective flux solver (E. McNeil).
dbpp::RhsHLLFluxSrc::RhsHLLFluxSrc | ( | std::string | aName, |
const PhysicalSystem * | aPhysys, | ||
bool | aUsePressure = false ) |
Constructor from physical system.
aName | class name |
aPhysys | physical system |
aUsePressure | flag to set pressure term |
dbpp::RhsHLLFluxSrc::RhsHLLFluxSrc | ( | std::string | aName, |
HLLFluxAlgorithm * | aFluxAlgo, | ||
SrcNumericalTreatment * | aSrcDiscr, | ||
const PhysicalBoundaryCnd & | aPhybc ) |
flux algorithm and source discretization
aName | string description |
aFluxAlgo | |
aSrcDiscr | |
aPhybc |
|
finaloverridevirtual |
calculate rhs terms (apply discretization at each term)
aU | state variables vector field |
Implements dbpp::SweRhsAlgorithm.
|
overridevirtual |
calculate rhs terms (apply discretization at each term)
aU | state vector of state variables |
aGblDiscr | finite volume discretization |
Implements dbpp::SweRhsAlgorithm.
|
inlinefinaloverridevirtual |
physical boundary condition in use
Implements dbpp::SweRhsAlgorithm.
|
inlinefinaloverridevirtualnoexcept |
|
inlineoverridevirtual |
Set Numerical flux algorithm.
aFluxAlgo | flux algorithm |
useIncompleteFlux | flux type take account ppressure or not |
Reimplemented from dbpp::SweRhsAlgorithm.
|
inlinefinaloverridevirtual |
set physical boundary condition
aPhyBnd | physical boundary condition |
Implements dbpp::SweRhsAlgorithm.
|
inlineoverridevirtual |
set numerical source treatment (discretization)
Reimplemented from dbpp::SweRhsAlgorithm.
|
inlinefinaloverridevirtualnoexcept |
|
inlinefinaloverridevirtualnoexcept |
|
private |
Convective flux algorithm
|
private |
far left domain node b.c.
|
private |
name
|
private |
physical boundary condition
|
private |
physical system
|
private |
far right domain node b.c.
|
private |
Source numerical discretization
|
private |
pressure term flag