Nektar++
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 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 
44 std::string QInflow::className = GetBoundaryFactory().RegisterCreatorFunction(
45  "Q-inflow", QInflow::create, "Inflow boundary condition");
46 
49  PulseWavePressureAreaSharedPtr pressureArea)
50  : PulseWaveBoundary(pVessel, pSession, pressureArea)
51 {
52 }
53 
55 {
56 }
57 
59  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
62  Array<OneD, Array<OneD, NekDouble>> &alpha, const NekDouble time, int omega,
63  int offset, int n)
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 
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 }
102 
103 /**
104  * Q-inflow Riemann solver for pulse wave propagation. This Riemann solver is
105  * called by DoOdeProjection in case of a Q-inflow boundary condition. It is
106  * based on the conservation of mass and total pressure and on the
107  * characteristic information. For further details see "Pulse wave propagation
108  * in the human vascular system", section 3.4.1 Returns the upwinded quantities
109  * \f$(A_u,u_u)\f$ and stores them into the boundary values
110  */
112  NekDouble A_0, NekDouble beta, NekDouble alpha,
113  NekDouble &Au, NekDouble &uu)
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 }
159 
160 } // namespace Nektar
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
PulseWavePressureAreaSharedPtr m_pressureArea
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
virtual ~QInflow()
Definition: QInflow.cpp:54
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) override
Definition: QInflow.cpp:58
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:61
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< PulseWavePressureArea > PulseWavePressureAreaSharedPtr
BoundaryFactory & GetBoundaryFactory()
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294