![]() |
DamBreak++ Wave Simulator 0.3
DamBreak++ Simulation Framework
|
Responsible to evaluate the spatial terms according to spatial discretization. The HOperator provides type and method to compute numerically the different terms of the right-hand-side of the equation. It also provide numerical container used for fast floating-point operation. More...
#include <Sfx_LDeltaOperator.h>
Classes | |
struct | SWE_RHS |
Data structure RHS terms. More... |
Public Types | |
enum class | eNumArrayType { std_valarray =0 , std_vector =1 , swe_NumArray =2 , ublas_vector =3 } |
Numerical array type for fast floating point operation. More... | |
enum class | eSourceTermType { Hderivative =0 , Zderivative =1 } |
enum class | eSourceDxType { central_2ndorder =0 , central_4thorder =1 } |
Source term derivative order. More... | |
enum class | eDxSchemeType { central_2ndorder = 0 , central_4thorder = 1 } |
Derivative order. More... | |
enum class | eDerivBCType { noncentered_2ndorder = 0 , emcil_2nd_order = 1 } |
Derivative boundary condition (order) More... | |
enum class | eFluxType { convective =0 , pressure =1 , conv_pressure =2 } |
flux discretization More... | |
enum class | eFluxDxType { upwind_1st = 0 , upwind_2nd = 1 } |
flux derivative type More... | |
using | ListSectFlowPtr = std::shared_ptr<Testvs19::ListSectionsFlow<dbpp::SectionFlow>> |
alias (List of section pointer) | |
using | uvec64 = boost::numeric::ublas::vector<float64> |
using | CalculFF |
Pointer-to-function (numerical flux algorithm) | |
using | sourcediscr = eSourceDxType |
using | sourceTerm = eSourceTermType |
using | dx_f = eFluxDxType |
Public Member Functions | |
LDeltaOperator () | |
default ctor | |
LDeltaOperator (dbpp::SweRhsAlgorithm *aRhsDiscr) | |
Ctor from rhs algorithm. | |
LDeltaOperator (const LDeltaOperator &aOther)=delete | |
copy ctor | |
LDeltaOperator & | operator= (const LDeltaOperator &aOther)=delete |
assign operator | |
LDeltaOperator (CalculFF aFluxDiscr) | |
Set flux algorithm. | |
~LDeltaOperator () | |
dtor | |
void * | applyTo (std::shared_ptr< Sfx::FieldLattice > &aU1, std::shared_ptr< Sfx::FieldLattice > &aU2) override |
Apply rhs discretization to each term. | |
void * | applyTo (const Sfx::FieldLattice &aU1, const Sfx::FieldLattice &aU2) |
Apply rhs discretization to each term. | |
void | solveFor () override |
apply rhs discretization | |
eNumArrayType | arrayType () const noexcept |
Getter. | |
void | setArrayType (const eNumArrayType aType) noexcept |
Set type of numerical array. | |
bool | haspressureTerm () const noexcept |
taking account pressure term in numerical model | |
bool | useFriction () const noexcept |
friction ON/OFF | |
void | setTraitementS2 (const eSourceTermType aSdx) noexcept |
Set type of treatment to evaluate source term (bed slope) | |
void | setDxSourceTerm (const sourcediscr aSdx) noexcept |
Set derivative order of bed slope term. | |
void | setTraitementF (const dx_f aDxF) noexcept |
numerical flux evaluation | |
void | setDerivativType (eDxSchemeType aDerivType) noexcept |
eDxSchemeType | getderivativeType () const noexcept |
void | setD1xBCAtBothEnd (eDerivBCType aBCtype) noexcept |
eDerivBCType | getD1xBCtype () const noexcept |
Getter. | |
void | setSrcTermRhsData (const ListSectFlowPtr &aSrcData) noexcept |
Some data to compute rhs. | |
void | setFluxAlgorithm (CalculFF aConvF) noexcept |
Set the convective numerical algorithm (ptr-2-function) | |
CalculFF | getFluxAlgorithm () noexcept |
Numerical flux algorithm. |
Static Public Attributes | |
static constexpr int32 | DIM = 100 |
Responsible to evaluate the spatial terms according to spatial discretization. The HOperator provides type and method to compute numerically the different terms of the right-hand-side of the equation. It also provide numerical container used for fast floating-point operation.
Pointer-to-function (numerical flux algorithm)
using Sfx::LDeltaOperator::dx_f = eFluxDxType |
using Sfx::LDeltaOperator::ListSectFlowPtr = std::shared_ptr<Testvs19::ListSectionsFlow<dbpp::SectionFlow>> |
alias (List of section pointer)
using Sfx::LDeltaOperator::uvec64 = boost::numeric::ublas::vector<float64> |
|
strong |
|
strong |
|
strong |
|
strong |
flux discretization
Enumerator | |
---|---|
convective | Harten-Lax-Levy |
pressure | Nujic ENO flux |
conv_pressure | Nujic Roe simplified |
|
strong |
|
strong |
Source term derivative order.
Enumerator | |
---|---|
central_2ndorder | second-order derivative |
central_4thorder | fourth-order derivative |
|
strong |
Sfx::LDeltaOperator::LDeltaOperator | ( | ) |
default ctor
< flux discretization scheme
Sfx::LDeltaOperator::LDeltaOperator | ( | dbpp::SweRhsAlgorithm * | aRhsDiscr | ) |
Ctor from rhs algorithm.
Rhs | term discretization algorithm |
< flux discretization scheme
|
delete |
copy ctor
object | to copy from |
Sfx::LDeltaOperator::LDeltaOperator | ( | CalculFF | aFluxDiscr | ) |
Set flux algorithm.
function | pointer |
Sfx::LDeltaOperator::~LDeltaOperator | ( | ) |
dtor
void * Sfx::LDeltaOperator::applyTo | ( | const Sfx::FieldLattice & | aU1, |
const Sfx::FieldLattice & | aU2 ) |
Apply rhs discretization to each term.
aU1 | first state variable |
aU2 | second state variable |
|
override |
Apply rhs discretization to each term.
aU1 | first state variable |
aU2 | second state variable |
|
inlinenoexcept |
Getter.
|
inlinenoexcept |
Getter.
|
inlinenoexcept |
|
inlinenoexcept |
Numerical flux algorithm.
|
inlinenoexcept |
taking account pressure term in numerical model
|
delete |
assign operator
object | to assign from |
|
inlinenoexcept |
Set type of numerical array.
aType | type |
|
inlinenoexcept |
aBCtype |
|
inlinenoexcept |
aDerivType |
|
inlinenoexcept |
Set derivative order of bed slope term.
aSdx | derivative order type |
|
inlinenoexcept |
Set the convective numerical algorithm (ptr-2-function)
aConvF | pointer-to-function type |
|
inlinenoexcept |
|
inlinenoexcept |
numerical flux evaluation
aDxF | flux type (...) |
|
inlinenoexcept |
Set type of treatment to evaluate source term (bed slope)
aSdx | source type (H/Z) variables |
|
inlineoverride |
apply rhs discretization
|
inlinenoexcept |
|
staticconstexpr |
comptational domain
|
private |
numerical array for fast floating point op
|
private |
boundary condition
|
private |
ptr-2-function flux algorithm flux discretization (convective)()
|
private |
finite difference scheme type
|
private |
numerical treatment of physical flux
|
private |
flux discretization scheme
|
private |
rhs terms discretization
|
private |
numerical treatment of source terms
|
private |
source term to be considered (H or Z Nujic algorithm)
|
private |
List section flow
|
private |
right-hand-side values