DamBreak++ Wave Simulator 0.3
DamBreak++ Simulation Framework
Loading...
Searching...
No Matches
dbpp_HLLSolver1D.h
Go to the documentation of this file.
1#pragma once
2
3// SfxBase19 include
4#include "include/Sfx_DefineTypes.h"
5// DamBreak app include
6#include "dbpp_FluxTensor.h"
7
8namespace dbpp { class RiemannProblem; }
9namespace Sfx { class cellFaceVariables; }
10namespace dbpp
11{
21 {
22 public:
26 HLLSolver1D() = default;
27 //HLLSolver1D( const CellFaceVariables& aCellFvar); // not sure
32 HLLSolver1D( const HLLSolver1D& aOther) = delete;
38 HLLSolver1D& operator= (const HLLSolver1D& aOther) = delete;
44 FluxTensor solve(const RiemannProblem& aRprob); //override
49 std::pair<double/*SL*/, double/*SR*/> getShockSpeed() const
50 {
51 return {0.,0./* m_shockSpeed.m_SL, m_shockSpeed.m_SR*/ };
52 }
53#if 0 // another interesting implementation (will be back later)
59 virtual fluxcomp solve(const RiemannProblem& aRprob)
60 {
61 // retrieve conservative system
62 const auto w_stVen1D = aRprob.getMthEquations();
63
64 // compute state variables physical flux at cell face
65 const auto FL1 = aRprob.getUL().Q();
66 const auto FR1 = aRprob.getUR().Q();
67
68 // physical flux (math model in use)
69 const auto FL2 = w_stVen1D->flux(Sfx::StateVector{ aRprob.getUL().A(), aRprob.getUL().Q() }); //FL2
70 const auto FR2 = w_stVen1D->flux(Sfx::StateVector{ aRprob.getUR().A(), aRprob.getUR().Q() }); //FR2
71
72 // compute flux parameters according to HLL HLL algorithm
73 emcil::HLLParameters w_hllprms;
74 // structured binding C++17 (tuple-like API)
75 const auto& [CR, CL, CS, uL, uR, uS] =
76 w_hllprms.setParameters(aRprob.getUL(), aRprob.getUR());
77
78 // shock speed on both side of the shock
79 const auto SL = std::min(uL - CL, uS - CS);
80 const auto SR = std::max(uR + CR, uS + CS);
81 // ...
82 HLLFlux1D w_hllFlux; // call shock speed
83 w_hllFlux.setShockSpeed(std::make_pair(SL, SR));
84 w_hllFlux.setphysicalFlux(std::make_tuple(FL1, FR1, FL2, FR2));
85
86 // [U1L,U2L] and [U1R,U2R]
88 Sfx::cellFaceVariables w_cellFaceVar{ static_cast<unsigned>(aRprob.getFaceIdx()), Sfx::StateVector{aRprob.getUL()},
89 Sfx::StateVector{aRprob.getUR()} };
90
91 // calculate numerical cell face flux
92 auto [FF1, FF2] = w_hllFlux.calculFF(w_cellFaceVar);
93
94 // Design note
95 // maybe return a struct??
96 return { FF1,FF2 }; // numeric face flux
97 }
98#endif // 0
99
100 private:
111
112 protected:
118 virtual std::pair<double, double> computeShockSpeed( const Sfx::cellFaceVariables& aFace); // override
124 virtual hllParams computeParams( const Sfx::cellFaceVariables& aFace); // override
125 };
126}//End of namespace
HLLSolver1D & operator=(const HLLSolver1D &aOther)=delete
denied copy/assignment
virtual hllParams computeParams(const Sfx::cellFaceVariables &aFace)
Compute HLL algorithm parameters.
Definition dbpp_HLLSolver1D.cpp:34
HLLSolver1D()=default
default ctor
HLLSolver1D(const HLLSolver1D &aOther)=delete
denied copy/assignment
virtual std::pair< double, double > computeShockSpeed(const Sfx::cellFaceVariables &aFace)
Compute shock speed values at cell faces.
Definition dbpp_HLLSolver1D.cpp:15
std::pair< double, double > getShockSpeed() const
Shock speed.
Definition dbpp_HLLSolver1D.h:49
FluxTensor solve(const RiemannProblem &aRprob)
solution of the Riemann problem at cell face j+1/2
Definition dbpp_HLLSolver1D.cpp:68
System of conservation Law with an initial discontinuity in it. Two states connected by a left and ri...
Definition dbpp_RiemannProblem.h:15
Sfx::StateVector getUR() const
Right state vector.
Definition dbpp_RiemannProblem.h:47
Sfx::SCLEquation * getMthEquations() const
System of Conservation LAw.
Definition dbpp_RiemannProblem.h:52
Sfx::StateVector getUL() const
Left state vector.
Definition dbpp_RiemannProblem.h:42
auto getFaceIdx() const
Face index.
Definition dbpp_RiemannProblem.h:62
Definition HydUtils.h:15
Definition DamBreakProb.h:15
double float64
Definition dbpp_LDeltaOperator.h:12
Flux tensor field (aggregate initialization).
Definition dbpp_FluxTensor.h:21
aggregate that hold hll algorithm parameters
Definition dbpp_HLLSolver1D.h:105
float64 uL
Definition dbpp_HLLSolver1D.h:106
float64 uR
Definition dbpp_HLLSolver1D.h:107
float64 CL
Definition dbpp_HLLSolver1D.h:108
float64 CR
Definition dbpp_HLLSolver1D.h:109