DamBreak++ Wave Simulator 0.3
DamBreak++ Simulation Framework
Loading...
Searching...
No Matches
Sfx Namespace Reference

Classes

class  DamBreakData
 Hold data for validation. Two type: EMcneil, Hudson (DamBreak config) More...
class  DamBreakExSol
class  DamBreakInitialCnd
 Sole responsibility is to provide wave profile at initial start up based on DamBreakData. More...
class  DbgLogger
 Helper utility that save result to a file to be used for debugging and visualizing. More...
class  EMcNeilBndCnd
 Calculate the characteristic information to set boundary node values (A,Q,H) More...
class  ExplicitOperator
 Some prototype. More...
class  FiniteDifferenceDerivative
 First derivative at second and fourth order. provide a stencil concept to evaluate derivative without spurious oscillations (ENO: Essentially Non-Oscillatory. More...
class  LDeltaOperator
 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...
class  Logger
 Class to output a message log file. Singleton base class manage the creation/deletion of the logger. That the instance function returns a pointer to a static variable and thus is declared static. Only the class function Instance can call the constructor. Public access to the constructor is denied. The constructor, copy constructor and assignment operator are all private to ensure that the programmer using the singleton class can only create a single instance of the class using only the Instance() function. The life of the singleton instantiation is for the duration of the application. More...
class  PhysicalAlgorithm
 Abstract class provide an interface with services to implement physical based algorithm to solve st-Venant equations in a conservative form (divergence). Numerical integrate the numerical model in time ... advance ... to be completed. More...
class  PhysicalConfigure
 Responsible to set initial configuration of the physical system (e.g. initial condition). More...
class  PhysicalMeasurement
 Physical measuremt taken on physical system (e.g. energy, friction, ...). More...
class  SectionGeometryHandler
class  Simulation
 Bean that represents a simulation in the framework. Many of the attributes of the simulation bean are defined in the Simulation.properties file by the user of the framework. The Controller class uses the simulation bean to control the creation of the correct classes for the simulation. More...
class  SimulationMgr
 Global instance of the simulation. Hold the different parameters of run (simulation parameters). Use the CRTP (Curious Recursive Template Pattern) known as static polymorphism) More...
class  StVenant1D
 Shallow-water (one-dimension) is a system of conservation laws. Conserved quantities (state variables) are the mass and the momentum. It's the differential form of Newton law (f=ma) but for a mass of fluid in a movement under gravity. Forces that act on the mass of fluid are friction, gravity and pressure. More...
class  StVenant1DTerms
 One-dimensional St-Venant equation spatial terms. More...

Typedefs

using mapFluxSchemeName = std::map<std::string, CalculFF>
 Alias (map name and function)

Functions

constexpr float64 computeVelocity (float64 aU1, float64 aU2)
 "Model of" binary function (STL function objects).
constexpr float64 WettedPerimeter (float64 aArea, float64 aWidth=1.)
 Function object (callable) that perform calculation cross-section ( wetted area, hydraulic radius, width at free surface.
constexpr float64 HydRadius (float64 aArea, float64 aWidth=1.)
 Function specific to a rectangular channel of width B Perimeter established according to a flow area A.
constexpr float64 WettedArea (float64 aDepth, float64 aWidth=1.)
constexpr float64 SectionWidth (float64 aWidth=1.)
std::valarray< float64 > cell_face_average (std::ranges::view auto aRng)
 Compute average of a variable at cell face.
template<typename Range1, typename Range2>
std::common_type_t< Range1, Range2 > computeLimitedGradient (const Range1 &aRng1, const Range2 &aRng2)
 Compute gradient by using a limiter function.
template<typename Range, typename FuncLimiter>
requires std::same_as<Range, std::vector<float64>>
auto computeDU (Range aRng, FuncLimiter &&aFuncLimiter)
 Compute gradient over cell by applying a slope limiter function. Range aRng Range of values (computational grid)
template<typename Range, typename FuncLimiter>
decltype(auto) computeDU (Range &aRng, FuncLimiter &&aFuncLimiter)
auto computeDU_v (std::ranges::view auto aView)
 state variables gradient over each cell (1st order)
template<typename Range>
auto cell_face_average (const Range &aRng)
float64 wrapPi (float64 aTheta)
 "Wrap" an angle in range -pi...pi by adding the correct multiple of 2 pi
float64 safeAcos (float64 x)
 "Safe" inverse trig functions
constexpr float64 degToRad (float64 aDeg)
 Convert between degrees and radians.
constexpr float64 radToDeg (float64 aRad)
void sinCos (float *aReturnSin, float *aReturnCos, float aTheta)
 Compute the sin and cosine of an angle. On some platforms, if we know that we need both values, it can be computed faster than computing the two values seperately.
template<typename Range>
auto StVenant1D_Incomplete_Flux (Range &aU1, Range &aU2)
template<std::ranges::input_range Rng>
requires std::same_as<Rng,std::vector<double>>
auto UpwindDerivatr1st (Rng aRng)
 Implementation of the upwind derivative scheme using generic programming.
void HLL_Scheme (DBSArrayType &aFF1, DBSArrayType &aFF2, const DBSArrayType &aU1, const DBSArrayType &aU2)
 Harten-Lax-Leer Algorithm ( belongs to Approximate Riemann Solver family)
void Nujic_SchemeI (DBSArrayType &aFF1, DBSArrayType &aFF2, const DBSArrayType &aU1, const DBSArrayType &aU2)
 Nujic(1995) ENO (Essentially Non-Oscillatory) extrapolation scheme (belongs to flux-splitting family)
void Nujic_SchemeII (DBSArrayType &aFF1, DBSArrayType &aFF2, const DBSArrayType &aU1, const DBSArrayType &aU2)
 Nujic (1995) simplified version of the "Approximate Roe Solver" using LxF (second-order accurate)
void LXF_Scheme (DBSArrayType &aFF1, DBSArrayType &aFF2, const DBSArrayType &aU1, const DBSArrayType &aU2)
 First-order Lax-Friedrich scheme ( belongs to ... family)

Variables

int Ns =0
int Nv =0
int Nst =0
int Nsu =0
int testing =0
void(* algorithm [])(DBSArrayType &aFF1, DBSArrayType &aFF2, const DBSArrayType &aU1, const DBSArrayType &aU2) = { HLL_Scheme, Nujic_SchemeI, Nujic_SchemeII, LXF_Scheme}
const char algorithmName [][20]

Typedef Documentation

◆ mapFluxSchemeName

using Sfx::mapFluxSchemeName = std::map<std::string, CalculFF>

Alias (map name and function)

Function Documentation

◆ cell_face_average() [1/2]

template<typename Range>
auto Sfx::cell_face_average ( const Range & aRng)

◆ cell_face_average() [2/2]

std::valarray< float64 > Sfx::cell_face_average ( std::ranges::view auto aRng)

Compute average of a variable at cell face.

Parameters
aRngrange concept (cheap to copy and move)
Returns
average value at each cell face

◆ computeDU() [1/2]

template<typename Range, typename FuncLimiter>
decltype(auto) Sfx::computeDU ( Range & aRng,
FuncLimiter && aFuncLimiter )
@brief Limiter function gradient

Range aRng pass by ref since const ref is part of the type FuncLimiter aFuncLimiter slope limiter function

◆ computeDU() [2/2]

template<typename Range, typename FuncLimiter>
requires std::same_as<Range, std::vector<float64>>
auto Sfx::computeDU ( Range aRng,
FuncLimiter && aFuncLimiter )

Compute gradient over cell by applying a slope limiter function. Range aRng Range of values (computational grid)

Parameters
FuncLimiteraFuncLimiter slope limiter function that belongs to TVD-family limit spurious oscillations near sharp peaked profile (select minimum slope)

◆ computeDU_v()

auto Sfx::computeDU_v ( std::ranges::view auto aView)

state variables gradient over each cell (1st order)

Template Parameters
Rangecontainer of values as view
Parameters
aRngnodal values (state variables)
Returns
derivative at nodal point Usage (adaptor convert range to view, cheap to copy) computeDU_v(std::ranges::all(rng)) -> ref_view (lvalue) computeDU_v(std::ranges::all(move(rng))) -> owning_view (rvalue)

◆ computeLimitedGradient()

template<typename Range1, typename Range2>
std::common_type_t< Range1, Range2 > Sfx::computeLimitedGradient ( const Range1 & aRng1,
const Range2 & aRng2 )

Compute gradient by using a limiter function.

Template Parameters
Range1range of values
Range2range of values
Sametype range
Parameters
aRng1first component value range
aRng2second component value range
Returns
pair of value range

◆ computeVelocity()

float64 Sfx::computeVelocity ( float64 aU1,
float64 aU2 )
inlineconstexpr

"Model of" binary function (STL function objects).

◆ degToRad()

float64 Sfx::degToRad ( float64 aDeg)
constexpr

Convert between degrees and radians.

Parameters
aDegdegree to convert
Returns
angle in radians

◆ HLL_Scheme()

void Sfx::HLL_Scheme ( DBSArrayType & aFF1,
DBSArrayType & aFF2,
const DBSArrayType & aU1,
const DBSArrayType & aU2 )

Harten-Lax-Leer Algorithm ( belongs to Approximate Riemann Solver family)

Parameters
aFF1Numerical flux at cell face j+1/2 (first state variable)
aFF2Numerical flux at cell face j+1/2 (second state variable)
aU1First state variable at global nodes
aU2Second state variable at global nodes

◆ HydRadius()

float64 Sfx::HydRadius ( float64 aArea,
float64 aWidth = 1. )
inlineconstexpr

Function specific to a rectangular channel of width B Perimeter established according to a flow area A.

◆ LXF_Scheme()

void Sfx::LXF_Scheme ( DBSArrayType & aFF1,
DBSArrayType & aFF2,
const DBSArrayType & aU1,
const DBSArrayType & aU2 )

First-order Lax-Friedrich scheme ( belongs to ... family)

Parameters
aFF1Numerical flux at cell face j+1/2 (first state variable)
aFF2Numerical flux at cell face j+1/2 (second state variable)
aU1First state variable at global nodes
aU2Second state variable at global nodes

◆ Nujic_SchemeI()

void Sfx::Nujic_SchemeI ( DBSArrayType & aFF1,
DBSArrayType & aFF2,
const DBSArrayType & aU1,
const DBSArrayType & aU2 )

Nujic(1995) ENO (Essentially Non-Oscillatory) extrapolation scheme (belongs to flux-splitting family)

Parameters
aFF1Numerical flux at cell face j+1/2 (first state variable)
aFF2Numerical flux at cell face j+1/2 (second state variable)
aU1First state variable at global nodes
aU2Second state variable at global nodes

◆ Nujic_SchemeII()

void Sfx::Nujic_SchemeII ( DBSArrayType & aFF1,
DBSArrayType & aFF2,
const DBSArrayType & aU1,
const DBSArrayType & aU2 )

Nujic (1995) simplified version of the "Approximate Roe Solver" using LxF (second-order accurate)

Parameters
aFF1Numerical flux at cell face j+1/2 (first state variable)
aFF2Numerical flux at cell face j+1/2 (second state variable)
aU1First state variable at global nodes
aU2Second state variable at global nodes

◆ radToDeg()

float64 Sfx::radToDeg ( float64 aRad)
constexpr

◆ safeAcos()

float64 Sfx::safeAcos ( float64 x)
inline

"Safe" inverse trig functions

Parameters
x
Returns
inverse

Same as acos(x), but if x is out of range, it is "clamped" to the nearest valid value. The value returned is in range 0...pi, the same as the standard C acos() function

◆ SectionWidth()

float64 Sfx::SectionWidth ( float64 aWidth = 1.)
inlineconstexpr

◆ sinCos()

void Sfx::sinCos ( float * aReturnSin,
float * aReturnCos,
float aTheta )
inline

Compute the sin and cosine of an angle. On some platforms, if we know that we need both values, it can be computed faster than computing the two values seperately.

Parameters
aReturnSin
aReturnCos
aTheta

◆ StVenant1D_Incomplete_Flux()

template<typename Range>
auto Sfx::StVenant1D_Incomplete_Flux ( Range & aU1,
Range & aU2 )

◆ UpwindDerivatr1st()

template<std::ranges::input_range Rng>
requires std::same_as<Rng,std::vector<double>>
auto Sfx::UpwindDerivatr1st ( Rng aRng)

Implementation of the upwind derivative scheme using generic programming.

Parameters
aRngcontainer of nodal values
Returns
derivative values

◆ WettedArea()

float64 Sfx::WettedArea ( float64 aDepth,
float64 aWidth = 1. )
inlineconstexpr

◆ WettedPerimeter()

float64 Sfx::WettedPerimeter ( float64 aArea,
float64 aWidth = 1. )
inlineconstexpr

Function object (callable) that perform calculation cross-section ( wetted area, hydraulic radius, width at free surface.

◆ wrapPi()

float64 Sfx::wrapPi ( float64 aTheta)
inline

"Wrap" an angle in range -pi...pi by adding the correct multiple of 2 pi

Parameters
thetaangle to wwrap
Returns
wrapped angle

Variable Documentation

◆ algorithm

void(* Sfx::algorithm[])(DBSArrayType &aFF1, DBSArrayType &aFF2, const DBSArrayType &aU1, const DBSArrayType &aU2) ( DBSArrayType & aFF1,
DBSArrayType & aFF2,
const DBSArrayType & aU1,
const DBSArrayType & aU2 ) = { HLL_Scheme, Nujic_SchemeI, Nujic_SchemeII, LXF_Scheme}

◆ algorithmName

const char Sfx::algorithmName[][20]
Initial value:
=
{ "HLL_Scheme", "Nujic_SchemeI", "Nujic_SchemeII", "LXF_Scheme"}

◆ Ns

int Sfx::Ns =0

◆ Nst

int Sfx::Nst =0

◆ Nsu

int Sfx::Nsu =0

◆ Nv

int Sfx::Nv =0

◆ testing

int Sfx::testing =0