![]() |
DamBreak++ Wave Simulator 0.3
DamBreak++ Simulation Framework
|
Algorithm based on the Godunov-type Scheme (explicit). We use a two-time steps integrator belongs to Runge-Kutta family and the HLL algorithm for convective flux evaluation. Order is achieved by MUSCL reconstruction procedure to extrapolate state variables at cell face j+1/2 (piecewise-constant state variable on each cell). This algorithm solve at cell-face a Riemann problem (system of conservation law coupled with 2-state tied by a discontinuity define a Riemann Problem). More...
#include <dbpp_RK2HLLScheme.h>
Public Member Functions | |
| Rk2HllScheme () | |
| default ctor | |
| Rk2HllScheme (dbpp::PhysicalSystem *aPhysys) | |
| ~Rk2HllScheme () | |
| destructor | |
| void | setSweRhsAlgorithm (SweRhsAlgorithm *aRhsAlgo) override |
| Set. | |
| SweRhsAlgorithm * | getSweRhsAlgorithm () const override final |
| Right-hand-side disretization algorithm in use. | |
| SweRhsAlgorithm * | getSweRhsAlgorithm () override final |
| Right-hand-side disretization algorithm in use. | |
| void | setSourceTerm (dbpp::SrcNumericalTreatment *aSrcTreatment) override final |
| Set source term discretization. | |
| void | setFaceFluxAlgorithm () override final |
| Set numerical flux algorithm. | |
| Uh | updatedValues () override final |
| update nodal values | |
| Public Member Functions inherited from dbpp::GodunovTypeScheme | |
| GodunovTypeScheme ()=default | |
| virtual | ~GodunovTypeScheme ()=default |
| destructor | |
| void | mainLoop (const std::shared_ptr< dbpp::FiniteVolumeDiscretization > &aGblDiscr, Sfx::SfxTimePrm &aTime) override final |
| update nodal value | |
| virtual dbpp::FluxTensor | calculFF (const dbpp::cellFace &aCellFace, const Sfx::StateVector &aUL, const Sfx::StateVector &aUR) |
| Compute the numerical flux (solve Riemann problem at each cell face) default implementation is HLL (Harten-Lax-van Leer algorithm) | |
| bool | isGodunovType () const override final |
| Approximate Riemann solver. | |
| Public Member Functions inherited from dbpp::NumericalMethod | |
| virtual void | initialize (const std::shared_ptr< FiniteVolumeDiscretization > &aGlbDiscr) |
| Initialization. | |
| virtual void | mainLoop (const std::shared_ptr< FiniteVolumeDiscretization > &aGlbDiscr, Sfx::SfxTimePrm &aTime)=0 |
| Main loop to update all nodal values. | |
| virtual bool | isFluxDiffSplitting () const |
| Numerical method type. | |
| virtual bool | isSemiDiscreteMethod () const |
| Numerical method type. | |
| virtual bool | isTimeDependent () const |
| System time dependency. | |
Protected Member Functions | |
| void | advance (const std::shared_ptr< dbpp::FiniteVolumeDiscretization > &aGblDiscr, const double aDt) override final |
| time stepping algorithm | |
| Protected Member Functions inherited from dbpp::GodunovTypeScheme | |
| GodunovTypeScheme (const GodunovTypeScheme &aOther)=default | |
| implicit conversion only by subclass | |
| GodunovTypeScheme (GodunovTypeScheme &&aOther)=default | |
Private Member Functions | |
| FluxTensorMap | computeFaceFlux (const dbpp::Omega &aOmega, const std::vector< float64 > &aVecA, const std::vector< float64 > &aVecQ) |
| Compute numerical face flux. | |
| FluxTensorMap | computeFaceFlux (const dbpp::Omega &aOmega, std::vector< float64 > &&aVecA, std::vector< float64 > &&aVecQ) |
| Compute numerical face flux (move semantic) | |
Private Attributes | |
| NujicIntegrator | m_twoStepsIntegrator |
| Physical boundary condition (deprecated) | |
| dbpp::PhysicalSystem * | m_physys |
| dbpp::SweRhsAlgorithm * | m_rhsAlgo |
Algorithm based on the Godunov-type Scheme (explicit). We use a two-time steps integrator belongs to Runge-Kutta family and the HLL algorithm for convective flux evaluation. Order is achieved by MUSCL reconstruction procedure to extrapolate state variables at cell face j+1/2 (piecewise-constant state variable on each cell). This algorithm solve at cell-face a Riemann problem (system of conservation law coupled with 2-state tied by a discontinuity define a Riemann Problem).
| dbpp::Rk2HllScheme::Rk2HllScheme | ( | ) |
default ctor
| dbpp::Rk2HllScheme::Rk2HllScheme | ( | dbpp::PhysicalSystem * | aPhysys | ) |
| aPhysys |
| dbpp::Rk2HllScheme::~Rk2HllScheme | ( | ) |
destructor
|
finaloverrideprotectedvirtual |
time stepping algorithm
| aU | state vector (state variables) |
| aDt | time step |
Implements dbpp::GodunovTypeScheme.
|
private |
Compute numerical face flux.
| aOmega | computational domain |
| aVecA | vector first state variable |
| aVecQ | vector second d state variable |
|
private |
Compute numerical face flux (move semantic)
| aOmega | computational domain |
| aVecA | vector first state variable |
| aVecQ | vector second d state variable |
|
inlinefinaloverridevirtual |
Right-hand-side disretization algorithm in use.
Implements dbpp::NumericalMethod.
|
inlinefinaloverridevirtual |
Right-hand-side disretization algorithm in use.
Implements dbpp::NumericalMethod.
|
inlinefinaloverridevirtual |
Set numerical flux algorithm.
Reimplemented from dbpp::GodunovTypeScheme.
|
inlinefinaloverridevirtual |
Set source term discretization.
| aSrcTreatment |
Reimplemented from dbpp::GodunovTypeScheme.
|
inlineoverridevirtual |
|
finaloverridevirtual |
|
private |
|
private |
|
private |
Physical boundary condition (deprecated)