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)
 
virtual ~QInflow ()
 
- 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

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)
 

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

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

◆ ~QInflow()

Nektar::QInflow::~QInflow ( )
virtual

Definition at line 53 of file QInflow.cpp.

54 {
55 }

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.

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

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

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

Implements Nektar::PulseWaveBoundary.

Definition at line 57 of file QInflow.cpp.

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

References 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:200
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 65 of file QInflow.h.