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 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: QInflow class
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
38 using namespace std;
39 
40 namespace Nektar
41 {
42 
43 std::string QInflow::className = GetBoundaryFactory().RegisterCreatorFunction(
44  "Q-inflow", QInflow::create, "Inflow boundary condition");
45 
48  PulseWavePressureAreaSharedPtr pressureArea)
49  : PulseWaveBoundary(pVessel, pSession, pressureArea)
50 {
51 }
52 
54 {
55 }
56 
58  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
62  const NekDouble time, int omega, int offset, int n)
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 
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 }
100 
101 /**
102  * Q-inflow Riemann solver for pulse wave propagation. This Riemann solver is
103  * called by DoOdeProjection in case of a Q-inflow boundary condition. It is
104  * based on the conservation of mass and total pressure and on the
105  * characteristic information. For further details see "Pulse wave propagation
106  * in the human vascular system", section 3.4.1 Returns the upwinded quantities
107  * \f$(A_u,u_u)\f$ and stores them into the boundary values
108  */
110  NekDouble A_0, NekDouble beta, NekDouble alpha,
111  NekDouble &Au, NekDouble &uu)
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 }
157 
158 } // namespace Nektar
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
PulseWavePressureAreaSharedPtr m_pressureArea
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)
Definition: QInflow.cpp:57
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
virtual ~QInflow()
Definition: QInflow.cpp:53
std::shared_ptr< SessionReader > SessionReaderSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< PulseWavePressureArea > PulseWavePressureAreaSharedPtr
BoundaryFactory & GetBoundaryFactory()
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267