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
38using namespace std;
39
40namespace Nektar
41{
42
46 "Empirical pressure area relationship for the arterial system");
47
51 : PulseWavePressureArea(pVessel, pSession)
52{
53}
54
56{
57}
58
60 const NekDouble &A,
61 const NekDouble &A0,
62 const NekDouble &dAUdx,
63 const NekDouble &gamma,
64 const NekDouble &alpha)
65{
66 NekDouble kappa = 0.0;
67 GetKappa(kappa, A, A0, alpha);
68
69 P = m_PExt - beta * sqrt(A0) * log(kappa) / (2 * alpha) -
70 gamma * dAUdx / sqrt(A); // Viscoelasticity
71}
72
74 const NekDouble &A, const NekDouble &A0,
75 const NekDouble &alpha)
76{
77 NekDouble kappa = 0.0;
78 GetKappa(kappa, A, A0, alpha);
79
80 c = sqrt(beta * sqrt(A0) / (2 * m_rho * kappa)); // Elastic
81}
82
84 const NekDouble &beta, const NekDouble &A,
85 const NekDouble &A0, const NekDouble &alpha)
86{
87 NekDouble I = 0.0;
88 GetCharIntegral(I, beta, A, A0, alpha);
89
90 W1 = u + I; // Elastic and assumes u0 = 0
91}
92
94 const NekDouble &beta, const NekDouble &A,
95 const NekDouble &A0, const NekDouble &alpha)
96{
97 NekDouble I = 0.0;
98 GetCharIntegral(I, beta, A, A0, alpha);
99
100 W2 = u - I; // Elastic and assumes u0 = 0
101}
102
104 const NekDouble &W2,
105 const NekDouble &beta,
106 const NekDouble &A0,
107 const NekDouble &alpha)
108{
109 NekDouble xi = (W1 - W2) * sqrt(m_rho / (8 * beta * sqrt(A0)));
110
111 A = A0 * exp(xi * (2 - alpha * xi));
112}
113
115 const NekDouble &W2)
116{
117 u = (W1 + W2) / 2;
118}
119
121 const NekDouble &beta,
122 const NekDouble &A,
123 const NekDouble &A0,
124 const NekDouble &alpha)
125{
126 NekDouble kappa = 0.0;
127 GetKappa(kappa, A, A0, alpha);
128
129 I = sqrt(2 * beta * sqrt(A0) / m_rho) * (1 - sqrt(kappa)) / alpha;
130}
131
135 const Array<OneD, NekDouble> &A0, const Array<OneD, NekDouble> &alpha,
136 const std::string &type)
137{
138 // General formulation
139 if (type == "Bifurcation")
140 {
141 NekMatrix<NekDouble> J(6, 6);
143
144 for (int i = 0; i < 3; ++i)
145 {
146 GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
147 }
148
149 J.SetValue(0, 0, 1);
150 J.SetValue(0, 1, 0);
151 J.SetValue(0, 2, 0);
152 J.SetValue(0, 3, c[0] / Au[0]);
153 J.SetValue(0, 4, 0);
154 J.SetValue(0, 5, 0);
155
156 J.SetValue(1, 0, 0);
157 J.SetValue(1, 1, 1);
158 J.SetValue(1, 2, 0);
159 J.SetValue(1, 3, 0);
160 J.SetValue(1, 4, -c[1] / Au[1]);
161 J.SetValue(1, 5, 0);
162
163 J.SetValue(2, 0, 0);
164 J.SetValue(2, 1, 0);
165 J.SetValue(2, 2, 1);
166 J.SetValue(2, 3, 0);
167 J.SetValue(2, 4, 0);
168 J.SetValue(2, 5, -c[2] / Au[2]);
169
170 J.SetValue(3, 0, Au[0]);
171 J.SetValue(3, 1, -Au[1]);
172 J.SetValue(3, 2, -Au[2]);
173 J.SetValue(3, 3, uu[0]);
174 J.SetValue(3, 4, -uu[1]);
175 J.SetValue(3, 5, -uu[2]);
176
177 J.SetValue(4, 0, 2 * uu[0]);
178 J.SetValue(4, 1, -2 * uu[1]);
179 J.SetValue(4, 2, 0);
180 J.SetValue(4, 3, 2 * c[0] * c[0] / Au[0]);
181 J.SetValue(4, 4, -2 * c[1] * c[1] / Au[1]);
182 J.SetValue(4, 5, 0);
183
184 J.SetValue(5, 0, 2 * uu[0]);
185 J.SetValue(5, 1, 0);
186 J.SetValue(5, 2, -2 * uu[2]);
187 J.SetValue(5, 3, 2 * c[0] * c[0] / Au[0]);
188 J.SetValue(5, 4, 0);
189 J.SetValue(5, 5, -2 * c[2] * c[2] / Au[2]);
190
191 invJ = J;
192 invJ.Invert();
193 }
194 else if (type == "Merge")
195 {
196 NekMatrix<NekDouble> J(6, 6);
198
199 for (int i = 0; i < 3; ++i)
200 {
201 GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
202 }
203
204 J.SetValue(0, 0, 1);
205 J.SetValue(0, 1, 0);
206 J.SetValue(0, 2, 0);
207 J.SetValue(0, 3, -c[0] / Au[0]);
208 J.SetValue(0, 4, 0);
209 J.SetValue(0, 5, 0);
210
211 J.SetValue(1, 0, 0);
212 J.SetValue(1, 1, 1);
213 J.SetValue(1, 2, 0);
214 J.SetValue(1, 3, 0);
215 J.SetValue(1, 4, c[1] / Au[1]);
216 J.SetValue(1, 5, 0);
217
218 J.SetValue(2, 0, 0);
219 J.SetValue(2, 1, 0);
220 J.SetValue(2, 2, 1);
221 J.SetValue(2, 3, 0);
222 J.SetValue(2, 4, 0);
223 J.SetValue(2, 5, c[2] / Au[2]);
224
225 J.SetValue(3, 0, Au[0]);
226 J.SetValue(3, 1, -Au[1]);
227 J.SetValue(3, 2, -Au[2]);
228 J.SetValue(3, 3, uu[0]);
229 J.SetValue(3, 4, -uu[1]);
230 J.SetValue(3, 5, -uu[2]);
231
232 J.SetValue(4, 0, 2 * uu[0]);
233 J.SetValue(4, 1, -2 * uu[1]);
234 J.SetValue(4, 2, 0);
235 J.SetValue(4, 3, 2 * c[0] * c[0] / Au[0]);
236 J.SetValue(4, 4, -2 * c[1] * c[1] / Au[1]);
237 J.SetValue(4, 5, 0);
238
239 J.SetValue(5, 0, 2 * uu[0]);
240 J.SetValue(5, 1, 0);
241 J.SetValue(5, 2, -2 * uu[2]);
242 J.SetValue(5, 3, 2 * c[0] * c[0] / Au[0]);
243 J.SetValue(5, 4, 0);
244 J.SetValue(5, 5, -2 * c[2] * c[2] / Au[2]);
245
246 invJ = J;
247 invJ.Invert();
248 }
249 else if (type == "Interface")
250 {
251 NekMatrix<NekDouble> J(4, 4);
253
254 for (int i = 0; i < 2; ++i)
255 {
256 GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
257 }
258
259 J.SetValue(0, 0, 1);
260 J.SetValue(0, 1, 0);
261 J.SetValue(0, 2, c[0] / Au[0]);
262 J.SetValue(0, 3, 0);
263
264 J.SetValue(1, 0, 0);
265 J.SetValue(1, 1, 1);
266 J.SetValue(1, 2, 0);
267 J.SetValue(1, 3, -c[1] / Au[1]);
268
269 J.SetValue(2, 0, Au[0]);
270 J.SetValue(2, 1, -Au[1]);
271 J.SetValue(2, 2, uu[0]);
272 J.SetValue(2, 3, -uu[1]);
273
274 J.SetValue(3, 0, 2 * uu[0]);
275 J.SetValue(3, 1, -2 * uu[1]);
276 J.SetValue(3, 2, 2 * c[0] * c[0] / Au[0]);
277 J.SetValue(3, 3, -2 * c[1] * c[1] / Au[1]);
278
279 invJ = J;
280 invJ.Invert();
281 }
282}
283
285 const NekDouble &A0,
286 const NekDouble &alpha)
287{
288 kappa = 1 - alpha * log(A / A0);
289}
290
291} // 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.
Definition: NekFactory.hpp:197
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:303
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294