Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | List of all members
Nektar::QInflow Class Reference

A global linear system. More...

#include <QInflow.h>

Inheritance diagram for Nektar::QInflow:
Inheritance graph
[legend]
Collaboration diagram for Nektar::QInflow:
Collaboration graph
[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, 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)
 Creates an instance of this class. More...
 

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, 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 &Au, NekDouble &uu)
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::PulseWaveBoundary
Array< OneD,
MultiRegions::ExpListSharedPtr
m_vessels
 
LibUtilities::SessionReaderSharedPtr m_session
 
PulseWavePressureAreaSharedPtr m_pressureArea
 
NekDouble m_pext
 
NekDouble m_pout
 
NekDouble m_rho
 

Detailed Description

A global linear system.

Definition at line 51 of file QInflow.h.

Constructor & Destructor Documentation

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

Definition at line 52 of file QInflow.cpp.

55  : PulseWaveBoundary(pVessel,pSession,pressureArea)
56  {
57  }
PulseWaveBoundary(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession, PulseWavePressureAreaSharedPtr &pressureArea)
Nektar::QInflow::~QInflow ( )
virtual

Definition at line 62 of file QInflow.cpp.

63  {
64 
65  }

Member Function Documentation

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

Creates an instance of this class.

Definition at line 55 of file QInflow.h.

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

58  {
59  return MemoryManager<QInflow>::AllocateSharedPtr(pVessel,pSession,pressureArea);
60  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void Nektar::QInflow::Q_RiemannSolver ( NekDouble  Q,
NekDouble  A_r,
NekDouble  u_r,
NekDouble  A_0,
NekDouble  beta,
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 114 of file QInflow.cpp.

References Nektar::PulseWaveBoundary::m_rho.

Referenced by v_DoBoundary().

116  {
117  NekDouble W2 = 0.0;
118  NekDouble A_calc = 0.0;
119  NekDouble fa = 0.0;
120  NekDouble dfa = 0.0;
121  NekDouble delta_A_calc = 0.0;
122  NekDouble rho = m_rho;
123 
124  int proceed = 1;
125  int iter = 0;
126  int MAX_ITER = 200;
127 
128  // Tolerances for the algorithm
129  NekDouble Tol = 1.0e-10;
130 
131  // Riemann invariant \f$W_2(Ar,ur)\f$
132  W2 = u_r - 4*sqrt(beta/(2*rho))*sqrt(sqrt(A_r));
133 
134  // Newton Iteration (Area only)
135  A_calc = A_r;
136  while ((proceed) && (iter < MAX_ITER))
137  {
138  iter =iter+1;
139 
140  fa = Q - W2*A_calc - A_calc*4*sqrt(beta/(2*rho))*sqrt(sqrt(A_calc));
141  dfa = -W2 - 5*sqrt(beta/(2*rho))*sqrt(sqrt(A_calc));
142  delta_A_calc = fa/dfa;
143  A_calc = A_calc - delta_A_calc;
144 
145  if (sqrt(delta_A_calc*delta_A_calc) < Tol)
146  proceed = 0;
147  }
148 
149  // Obtain u_u and A_u
150  uu = W2+4*sqrt(beta/(2*rho))*sqrt(sqrt(A_calc));
151  Au = A_calc;
152  }
double NekDouble
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,
const NekDouble  time,
int  omega,
int  offset,
int  n 
)
protectedvirtual

Implements Nektar::PulseWaveBoundary.

Definition at line 67 of file QInflow.cpp.

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

75  {
76  NekDouble Q, A_u, u_u;
77  NekDouble A_r, u_r, A_l, u_l;
78  Array<OneD, MultiRegions::ExpListSharedPtr> vessel(2);
79 
80  vessel[0] = m_vessels[2*omega];
81  vessel[1] = m_vessels[2*omega+1];
82 
83  vessel[0]->EvaluateBoundaryConditions(time);
84 
85  // Note: The Q value is contained in A in the
86  // inputfile, the value in u has to be 1.0
87  //ASSERTL0((vessel[1]->UpdateBndCondExpansion(n))->UpdatePhys()[0] == 1.0, "For the Q-inflow BC the value in u must be 1.0");
88 
89  // Get the values of all variables needed for the Riemann problem
90  Q = (vessel[0]->UpdateBndCondExpansion(n))->GetCoeffs()[0];
91 
92  A_r = inarray[0][offset];
93  u_r = inarray[1][offset];
94 
95  // Call the Q-inflow Riemann solver
96  Q_RiemannSolver(Q,A_r,u_r,A_0[omega][0],beta[omega][0],A_u,u_u);
97 
98  // Set the boundary conditions to prescribe
99  A_l=A_r;
100  u_l=2*u_u-u_r;
101 
102  // Store the updated values in the boundary condition
103  (vessel[0]->UpdateBndCondExpansion(n))->UpdatePhys()[0] = A_l;
104  (vessel[1]->UpdateBndCondExpansion(n))->UpdatePhys()[0] = u_l;
105  }
void Q_RiemannSolver(NekDouble Q, NekDouble A_r, NekDouble u_r, NekDouble A_0, NekDouble beta, NekDouble &Au, NekDouble &uu)
Definition: QInflow.cpp:114
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
double NekDouble

Member Data Documentation

std::string Nektar::QInflow::className
static
Initial value:
"Q-inflow",
"Inflow boundary condition")

Name of class.

Definition at line 63 of file QInflow.h.