Nektar++
EmpiricalPressureArea.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: EmpiricalPressureArea.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: EmpiricalPressureArea class
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38namespace Nektar
39{
40
44 "Empirical pressure area relationship for the arterial system");
45
49 : PulseWavePressureArea(pVessel, pSession)
50{
51}
52
54 const NekDouble &A,
55 const NekDouble &A0,
56 const NekDouble &dAUdx,
57 const NekDouble &gamma,
58 const NekDouble &alpha)
59{
60 NekDouble kappa = 0.0;
61 GetKappa(kappa, A, A0, alpha);
62
63 P = m_PExt - beta * sqrt(A0) * log(kappa) / (2 * alpha) -
64 gamma * dAUdx / sqrt(A); // Viscoelasticity
65}
66
68 const NekDouble &A, const NekDouble &A0,
69 const NekDouble &alpha)
70{
71 NekDouble kappa = 0.0;
72 GetKappa(kappa, A, A0, alpha);
73
74 c = sqrt(beta * sqrt(A0) / (2 * m_rho * kappa)); // Elastic
75}
76
78 const NekDouble &beta, const NekDouble &A,
79 const NekDouble &A0, const NekDouble &alpha)
80{
81 NekDouble I = 0.0;
82 GetCharIntegral(I, beta, A, A0, alpha);
83
84 W1 = u + I; // Elastic and assumes u0 = 0
85}
86
88 const NekDouble &beta, const NekDouble &A,
89 const NekDouble &A0, const NekDouble &alpha)
90{
91 NekDouble I = 0.0;
92 GetCharIntegral(I, beta, A, A0, alpha);
93
94 W2 = u - I; // Elastic and assumes u0 = 0
95}
96
98 const NekDouble &W2,
99 const NekDouble &beta,
100 const NekDouble &A0,
101 const NekDouble &alpha)
102{
103 NekDouble xi = (W1 - W2) * sqrt(m_rho / (8 * beta * sqrt(A0)));
104
105 A = A0 * exp(xi * (2 - alpha * xi));
106}
107
109 const NekDouble &W2)
110{
111 u = (W1 + W2) / 2;
112}
113
115 const NekDouble &beta,
116 const NekDouble &A,
117 const NekDouble &A0,
118 const NekDouble &alpha)
119{
120 NekDouble kappa = 0.0;
121 GetKappa(kappa, A, A0, alpha);
122
123 I = sqrt(2 * beta * sqrt(A0) / m_rho) * (1 - sqrt(kappa)) / alpha;
124}
125
129 const Array<OneD, NekDouble> &A0, const Array<OneD, NekDouble> &alpha,
130 const std::string &type)
131{
132 // General formulation
133 if (type == "Bifurcation")
134 {
135 NekMatrix<NekDouble> J(6, 6);
137
138 for (int i = 0; i < 3; ++i)
139 {
140 GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
141 }
142
143 J.SetValue(0, 0, 1);
144 J.SetValue(0, 1, 0);
145 J.SetValue(0, 2, 0);
146 J.SetValue(0, 3, c[0] / Au[0]);
147 J.SetValue(0, 4, 0);
148 J.SetValue(0, 5, 0);
149
150 J.SetValue(1, 0, 0);
151 J.SetValue(1, 1, 1);
152 J.SetValue(1, 2, 0);
153 J.SetValue(1, 3, 0);
154 J.SetValue(1, 4, -c[1] / Au[1]);
155 J.SetValue(1, 5, 0);
156
157 J.SetValue(2, 0, 0);
158 J.SetValue(2, 1, 0);
159 J.SetValue(2, 2, 1);
160 J.SetValue(2, 3, 0);
161 J.SetValue(2, 4, 0);
162 J.SetValue(2, 5, -c[2] / Au[2]);
163
164 J.SetValue(3, 0, Au[0]);
165 J.SetValue(3, 1, -Au[1]);
166 J.SetValue(3, 2, -Au[2]);
167 J.SetValue(3, 3, uu[0]);
168 J.SetValue(3, 4, -uu[1]);
169 J.SetValue(3, 5, -uu[2]);
170
171 J.SetValue(4, 0, 2 * uu[0]);
172 J.SetValue(4, 1, -2 * uu[1]);
173 J.SetValue(4, 2, 0);
174 J.SetValue(4, 3, 2 * c[0] * c[0] / Au[0]);
175 J.SetValue(4, 4, -2 * c[1] * c[1] / Au[1]);
176 J.SetValue(4, 5, 0);
177
178 J.SetValue(5, 0, 2 * uu[0]);
179 J.SetValue(5, 1, 0);
180 J.SetValue(5, 2, -2 * uu[2]);
181 J.SetValue(5, 3, 2 * c[0] * c[0] / Au[0]);
182 J.SetValue(5, 4, 0);
183 J.SetValue(5, 5, -2 * c[2] * c[2] / Au[2]);
184
185 invJ = J;
186 invJ.Invert();
187 }
188 else if (type == "Merge")
189 {
190 NekMatrix<NekDouble> J(6, 6);
192
193 for (int i = 0; i < 3; ++i)
194 {
195 GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
196 }
197
198 J.SetValue(0, 0, 1);
199 J.SetValue(0, 1, 0);
200 J.SetValue(0, 2, 0);
201 J.SetValue(0, 3, -c[0] / Au[0]);
202 J.SetValue(0, 4, 0);
203 J.SetValue(0, 5, 0);
204
205 J.SetValue(1, 0, 0);
206 J.SetValue(1, 1, 1);
207 J.SetValue(1, 2, 0);
208 J.SetValue(1, 3, 0);
209 J.SetValue(1, 4, c[1] / Au[1]);
210 J.SetValue(1, 5, 0);
211
212 J.SetValue(2, 0, 0);
213 J.SetValue(2, 1, 0);
214 J.SetValue(2, 2, 1);
215 J.SetValue(2, 3, 0);
216 J.SetValue(2, 4, 0);
217 J.SetValue(2, 5, c[2] / Au[2]);
218
219 J.SetValue(3, 0, Au[0]);
220 J.SetValue(3, 1, -Au[1]);
221 J.SetValue(3, 2, -Au[2]);
222 J.SetValue(3, 3, uu[0]);
223 J.SetValue(3, 4, -uu[1]);
224 J.SetValue(3, 5, -uu[2]);
225
226 J.SetValue(4, 0, 2 * uu[0]);
227 J.SetValue(4, 1, -2 * uu[1]);
228 J.SetValue(4, 2, 0);
229 J.SetValue(4, 3, 2 * c[0] * c[0] / Au[0]);
230 J.SetValue(4, 4, -2 * c[1] * c[1] / Au[1]);
231 J.SetValue(4, 5, 0);
232
233 J.SetValue(5, 0, 2 * uu[0]);
234 J.SetValue(5, 1, 0);
235 J.SetValue(5, 2, -2 * uu[2]);
236 J.SetValue(5, 3, 2 * c[0] * c[0] / Au[0]);
237 J.SetValue(5, 4, 0);
238 J.SetValue(5, 5, -2 * c[2] * c[2] / Au[2]);
239
240 invJ = J;
241 invJ.Invert();
242 }
243 else if (type == "Interface")
244 {
245 NekMatrix<NekDouble> J(4, 4);
247
248 for (int i = 0; i < 2; ++i)
249 {
250 GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
251 }
252
253 J.SetValue(0, 0, 1);
254 J.SetValue(0, 1, 0);
255 J.SetValue(0, 2, c[0] / Au[0]);
256 J.SetValue(0, 3, 0);
257
258 J.SetValue(1, 0, 0);
259 J.SetValue(1, 1, 1);
260 J.SetValue(1, 2, 0);
261 J.SetValue(1, 3, -c[1] / Au[1]);
262
263 J.SetValue(2, 0, Au[0]);
264 J.SetValue(2, 1, -Au[1]);
265 J.SetValue(2, 2, uu[0]);
266 J.SetValue(2, 3, -uu[1]);
267
268 J.SetValue(3, 0, 2 * uu[0]);
269 J.SetValue(3, 1, -2 * uu[1]);
270 J.SetValue(3, 2, 2 * c[0] * c[0] / Au[0]);
271 J.SetValue(3, 3, -2 * c[1] * c[1] / Au[1]);
272
273 invJ = J;
274 invJ.Invert();
275 }
276}
277
279 const NekDouble &A0,
280 const NekDouble &alpha)
281{
282 kappa = 1 - alpha * log(A / A0);
283}
284
285} // namespace Nektar
void v_GetC(NekDouble &c, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5) override
void v_GetPressure(NekDouble &P, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &dAUdx, const NekDouble &gamma=0, const NekDouble &alpha=0.5) override
void GetKappa(NekDouble &kappa, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
EmpiricalPressureArea(Array< OneD, MultiRegions::ExpListSharedPtr > pVessel, const LibUtilities::SessionReaderSharedPtr pSession)
void v_GetW1(NekDouble &W1, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5) override
void v_GetJacobianInverse(NekMatrix< NekDouble > &invJ, const Array< OneD, NekDouble > &Au, const Array< OneD, NekDouble > &uu, const Array< OneD, NekDouble > &beta, const Array< OneD, NekDouble > &A0, const Array< OneD, NekDouble > &alpha, const std::string &type) override
void v_GetCharIntegral(NekDouble &I, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5) override
static PulseWavePressureAreaSharedPtr create(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession)
void v_GetUFromChars(NekDouble &u, const NekDouble &W1, const NekDouble &W2) override
void v_GetAFromChars(NekDouble &A, const NekDouble &W1, const NekDouble &W2, const NekDouble &beta, const NekDouble &A0, const NekDouble &alpha=0.5) override
void v_GetW2(NekDouble &W2, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5) override
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
void GetCharIntegral(NekDouble &I, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
void GetC(NekDouble &c, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
@ P
Monomial polynomials .
Definition: BasisType.h:62
PressureAreaFactory & GetPressureAreaFactory()
double NekDouble
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:294
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:285