Nektar++
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QInflow.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: QInflow.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28
29// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
30// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
31// DEALINGS IN THE SOFTWARE.
32//
33// Description: QInflow class
34//
35///////////////////////////////////////////////////////////////////////////////
36
38
39namespace Nektar
40{
41
43 "Q-inflow", QInflow::create, "Inflow boundary condition");
44
48 : PulseWaveBoundary(pVessel, pSession, pressureArea)
49{
50}
51
53 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
56 Array<OneD, Array<OneD, NekDouble>> &alpha, const NekDouble time, int omega,
57 int offset, int n)
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
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}
96
97/**
98 * Q-inflow Riemann solver for pulse wave propagation. This Riemann solver is
99 * called by DoOdeProjection in case of a Q-inflow boundary condition. It is
100 * based on the conservation of mass and total pressure and on the
101 * characteristic information. For further details see "Pulse wave propagation
102 * in the human vascular system", section 3.4.1 Returns the upwinded quantities
103 * \f$(A_u,u_u)\f$ and stores them into the boundary values
104 */
106 NekDouble A_0, NekDouble beta, NekDouble alpha,
107 NekDouble &Au, NekDouble &uu)
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}
153
154} // namespace Nektar
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
PulseWavePressureAreaSharedPtr m_pressureArea
static std::string className
Name of class.
Definition: QInflow.h:68
static PulseWaveBoundarySharedPtr create(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession, PulseWavePressureAreaSharedPtr &pressureArea)
Definition: QInflow.h:58
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
Definition: QInflow.cpp:52
QInflow(Array< OneD, MultiRegions::ExpListSharedPtr > pVessel, const LibUtilities::SessionReaderSharedPtr pSession, PulseWavePressureAreaSharedPtr pressureArea)
Definition: QInflow.cpp:45
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
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
std::shared_ptr< PulseWavePressureArea > PulseWavePressureAreaSharedPtr
BoundaryFactory & GetBoundaryFactory()
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:285