Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Friends | List of all members
Nektar::QInflow Class Reference

#include <QInflow.h>

Inheritance diagram for Nektar::QInflow:
[legend]

Static Public Member Functions

static PulseWaveBoundarySharedPtr create (Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession, PulseWavePressureAreaSharedPtr &pressureArea)
 

Static Public Attributes

static std::string className
 Name of class. More...
 

Protected Member Functions

 QInflow (Array< OneD, MultiRegions::ExpListSharedPtr > pVessel, const LibUtilities::SessionReaderSharedPtr pSession, PulseWavePressureAreaSharedPtr pressureArea)
 
 ~QInflow () override=default
 
void v_DoBoundary (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &A_0, Array< OneD, Array< OneD, NekDouble > > &beta, Array< OneD, Array< OneD, NekDouble > > &alpha, const NekDouble time, int omega, int offset, int n) override
 
virtual void v_DoBoundary (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &A_0, Array< OneD, Array< OneD, NekDouble > > &beta, Array< OneD, Array< OneD, NekDouble > > &alpha, const NekDouble time, int omega, int offset, int n)=0
 

Private Member Functions

void Q_RiemannSolver (NekDouble Q, NekDouble A_r, NekDouble u_r, NekDouble A_0, NekDouble beta, NekDouble alpha, NekDouble &Au, NekDouble &uu)
 

Friends

class MemoryManager< QInflow >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::PulseWaveBoundary
 PulseWaveBoundary (Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession, PulseWavePressureAreaSharedPtr &pressureArea)
 
virtual ~PulseWaveBoundary ()=default
 
void DoBoundary (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &A_0, Array< OneD, Array< OneD, NekDouble > > &beta, Array< OneD, Array< OneD, NekDouble > > &alpha, const NekDouble time, int omega, int offset, int n)
 
- Protected Attributes inherited from Nektar::PulseWaveBoundary
Array< OneD, MultiRegions::ExpListSharedPtrm_vessels
 
LibUtilities::SessionReaderSharedPtr m_session
 
PulseWavePressureAreaSharedPtr m_pressureArea
 
NekDouble m_pext
 
NekDouble m_pout
 
NekDouble m_rho
 

Detailed Description

Definition at line 52 of file QInflow.h.

Constructor & Destructor Documentation

◆ QInflow()

Nektar::QInflow::QInflow ( Array< OneD, MultiRegions::ExpListSharedPtr pVessel,
const LibUtilities::SessionReaderSharedPtr  pSession,
PulseWavePressureAreaSharedPtr  pressureArea 
)
protected

Definition at line 45 of file QInflow.cpp.

48 : PulseWaveBoundary(pVessel, pSession, pressureArea)
49{
50}
PulseWaveBoundary(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession, PulseWavePressureAreaSharedPtr &pressureArea)

◆ ~QInflow()

Nektar::QInflow::~QInflow ( )
overrideprotecteddefault

Member Function Documentation

◆ create()

static PulseWaveBoundarySharedPtr Nektar::QInflow::create ( Array< OneD, MultiRegions::ExpListSharedPtr > &  pVessel,
const LibUtilities::SessionReaderSharedPtr pSession,
PulseWavePressureAreaSharedPtr pressureArea 
)
inlinestatic

Definition at line 58 of file QInflow.h.

62 {
63 return MemoryManager<QInflow>::AllocateSharedPtr(pVessel, pSession,
64 pressureArea);
65 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ Q_RiemannSolver()

void Nektar::QInflow::Q_RiemannSolver ( NekDouble  Q,
NekDouble  A_r,
NekDouble  u_r,
NekDouble  A_0,
NekDouble  beta,
NekDouble  alpha,
NekDouble Au,
NekDouble uu 
)
private

Q-inflow Riemann solver for pulse wave propagation. This Riemann solver is called by DoOdeProjection in case of a Q-inflow boundary condition. It is based on the conservation of mass and total pressure and on the characteristic information. For further details see "Pulse wave propagation in the human vascular system", section 3.4.1 Returns the upwinded quantities \((A_u,u_u)\) and stores them into the boundary values

Definition at line 105 of file QInflow.cpp.

108{
109 NekDouble W2 = 0.0;
110 NekDouble A_calc = 0.0;
111 NekDouble FA = 0.0;
112 NekDouble dFdA = 0.0;
113 NekDouble delta_A_calc = 0.0;
114 NekDouble I = 0.0;
115 NekDouble c = 0.0;
116
117 int proceed = 1;
118 int iter = 0;
119 int MAX_ITER = 200;
120
121 // Tolerances for the algorithm
122 NekDouble Tol = 1.0e-10;
123
124 // Riemann invariant \f$W_2(Ar,ur)\f$
125 m_pressureArea->GetW2(W2, u_r, beta, A_r, A_0, alpha);
126
127 // Newton Iteration (Area only)
128 A_calc = A_r;
129 while ((proceed) && (iter < MAX_ITER))
130 {
131 iter += 1;
132
133 m_pressureArea->GetCharIntegral(I, beta, A_calc, A_0, alpha);
134 m_pressureArea->GetC(c, beta, A_calc, A_0, alpha);
135
136 FA = Q - A_calc * (W2 + I);
137 dFdA = -W2 - I - c;
138 delta_A_calc = FA / dFdA;
139 A_calc -= delta_A_calc;
140
141 if (sqrt(delta_A_calc * delta_A_calc) < Tol)
142 {
143 proceed = 0;
144 }
145 }
146
147 m_pressureArea->GetCharIntegral(I, beta, A_calc, A_0, alpha);
148
149 // Obtain u_u and A_u
150 uu = W2 + I;
151 Au = A_calc;
152}
PulseWavePressureAreaSharedPtr m_pressureArea
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:285

References Nektar::LibUtilities::beta, Nektar::PulseWaveBoundary::m_pressureArea, and tinysimd::sqrt().

Referenced by v_DoBoundary().

◆ v_DoBoundary()

void Nektar::QInflow::v_DoBoundary ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  A_0,
Array< OneD, Array< OneD, NekDouble > > &  beta,
Array< OneD, Array< OneD, NekDouble > > &  alpha,
const NekDouble  time,
int  omega,
int  offset,
int  n 
)
overrideprotectedvirtual

Implements Nektar::PulseWaveBoundary.

Definition at line 52 of file QInflow.cpp.

58{
59 NekDouble Q;
60 NekDouble A_u;
61 NekDouble u_u;
62 NekDouble A_r;
63 NekDouble u_r;
64 NekDouble A_l;
65 NekDouble u_l;
66
67 Array<OneD, MultiRegions::ExpListSharedPtr> vessel(2);
68
69 // Pointers to the domains
70 vessel[0] = m_vessels[2 * omega];
71 vessel[1] = m_vessels[2 * omega + 1];
72
73 // Evaluates the time-dependent Q
74 vessel[0]->EvaluateBoundaryConditions(time);
75
76 // Q is contained as A in the input file
77 Q = (vessel[0]->UpdateBndCondExpansion(n))->GetCoeffs()[0];
78
79 // Initial conditions in the inlet vessel
80 A_r = inarray[0][offset];
81 u_r = inarray[1][offset];
82
83 // Call the Q-inflow Riemann solver
84 Q_RiemannSolver(Q, A_r, u_r, A_0[omega][0], beta[omega][0], alpha[omega][0],
85 A_u, u_u);
86
87 /* Fix the boundary conditions in the virtual region to ensure
88 upwind state matches the boundary condition at the next time step */
89 A_l = A_r;
90 u_l = 2 * u_u - u_r;
91
92 // Store the updated values in the boundary condition
93 (vessel[0]->UpdateBndCondExpansion(n))->UpdatePhys()[0] = A_l;
94 (vessel[1]->UpdateBndCondExpansion(n))->UpdatePhys()[0] = u_l;
95}
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
void Q_RiemannSolver(NekDouble Q, NekDouble A_r, NekDouble u_r, NekDouble A_0, NekDouble beta, NekDouble alpha, NekDouble &Au, NekDouble &uu)
Definition: QInflow.cpp:105

References Nektar::LibUtilities::beta, Nektar::PulseWaveBoundary::m_vessels, and Q_RiemannSolver().

Friends And Related Function Documentation

◆ MemoryManager< QInflow >

friend class MemoryManager< QInflow >
friend

Definition at line 49 of file QInflow.h.

Member Data Documentation

◆ className

std::string Nektar::QInflow::className
static
Initial value:
"Q-inflow", QInflow::create, "Inflow boundary condition")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static PulseWaveBoundarySharedPtr create(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession, PulseWavePressureAreaSharedPtr &pressureArea)
Definition: QInflow.h:58
BoundaryFactory & GetBoundaryFactory()

Name of class.

Definition at line 68 of file QInflow.h.