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

#include <QInflow.h>

Inheritance diagram for Nektar::QInflow:
[legend]

Public Member Functions

 QInflow (Array< OneD, MultiRegions::ExpListSharedPtr > pVessel, const LibUtilities::SessionReaderSharedPtr pSession, PulseWavePressureAreaSharedPtr pressureArea)
 
 ~QInflow () override
 
- Public Member Functions inherited from Nektar::PulseWaveBoundary
 PulseWaveBoundary (Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession, PulseWavePressureAreaSharedPtr &pressureArea)
 
virtual ~PulseWaveBoundary ()
 
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)
 

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

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)
 

Additional Inherited Members

- 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 
)

Definition at line 47 of file QInflow.cpp.

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

◆ ~QInflow()

Nektar::QInflow::~QInflow ( )
override

Definition at line 54 of file QInflow.cpp.

55{
56}

Member Function Documentation

◆ create()

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

Definition at line 56 of file QInflow.h.

60 {
61 return MemoryManager<QInflow>::AllocateSharedPtr(pVessel, pSession,
62 pressureArea);
63 }
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 111 of file QInflow.cpp.

114{
115 NekDouble W2 = 0.0;
116 NekDouble A_calc = 0.0;
117 NekDouble FA = 0.0;
118 NekDouble dFdA = 0.0;
119 NekDouble delta_A_calc = 0.0;
120 NekDouble I = 0.0;
121 NekDouble c = 0.0;
122
123 int proceed = 1;
124 int iter = 0;
125 int MAX_ITER = 200;
126
127 // Tolerances for the algorithm
128 NekDouble Tol = 1.0e-10;
129
130 // Riemann invariant \f$W_2(Ar,ur)\f$
131 m_pressureArea->GetW2(W2, u_r, beta, A_r, A_0, alpha);
132
133 // Newton Iteration (Area only)
134 A_calc = A_r;
135 while ((proceed) && (iter < MAX_ITER))
136 {
137 iter += 1;
138
139 m_pressureArea->GetCharIntegral(I, beta, A_calc, A_0, alpha);
140 m_pressureArea->GetC(c, beta, A_calc, A_0, alpha);
141
142 FA = Q - A_calc * (W2 + I);
143 dFdA = -W2 - I - c;
144 delta_A_calc = FA / dFdA;
145 A_calc -= delta_A_calc;
146
147 if (sqrt(delta_A_calc * delta_A_calc) < Tol)
148 {
149 proceed = 0;
150 }
151 }
152
153 m_pressureArea->GetCharIntegral(I, beta, A_calc, A_0, alpha);
154
155 // Obtain u_u and A_u
156 uu = W2 + I;
157 Au = A_calc;
158}
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:294

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 58 of file QInflow.cpp.

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

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

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.
Definition: NekFactory.hpp:197
static PulseWaveBoundarySharedPtr create(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession, PulseWavePressureAreaSharedPtr &pressureArea)
Definition: QInflow.h:56
BoundaryFactory & GetBoundaryFactory()

Name of class.

Definition at line 66 of file QInflow.h.