Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::EmpiricalPressureArea Class Reference

#include <EmpiricalPressureArea.h>

Inheritance diagram for Nektar::EmpiricalPressureArea:
[legend]

Public Member Functions

 EmpiricalPressureArea (Array< OneD, MultiRegions::ExpListSharedPtr > pVessel, const LibUtilities::SessionReaderSharedPtr pSession)
 
 ~EmpiricalPressureArea () override
 
- Public Member Functions inherited from Nektar::PulseWavePressureArea
 PulseWavePressureArea (Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual ~PulseWavePressureArea ()
 
void GetPressure (NekDouble &P, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &dAUdx, const NekDouble &gamma=0, const NekDouble &alpha=0.5)
 
void GetC (NekDouble &c, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
 
void GetW1 (NekDouble &W1, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
 
void GetW2 (NekDouble &W2, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
 
void GetAFromChars (NekDouble &A, const NekDouble &W1, const NekDouble &W2, const NekDouble &beta, const NekDouble &A0, const NekDouble &alpha=0.5)
 
void GetUFromChars (NekDouble &u, const NekDouble &W1, const NekDouble &W2)
 
void GetCharIntegral (NekDouble &I, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
 
void 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)
 

Static Public Member Functions

static PulseWavePressureAreaSharedPtr create (Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession)
 

Static Public Attributes

static std::string className
 

Protected Member Functions

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 v_GetC (NekDouble &c, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5) override
 
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_GetW2 (NekDouble &W2, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5) 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_GetUFromChars (NekDouble &u, const NekDouble &W1, const NekDouble &W2) override
 
void v_GetCharIntegral (NekDouble &I, 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 GetKappa (NekDouble &kappa, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
 
virtual 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)=0
 
virtual void v_GetC (NekDouble &c, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)=0
 
virtual void v_GetW1 (NekDouble &W1, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)=0
 
virtual void v_GetW2 (NekDouble &W2, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)=0
 
virtual void v_GetAFromChars (NekDouble &A, const NekDouble &W1, const NekDouble &W2, const NekDouble &beta, const NekDouble &A0, const NekDouble &alpha=0.5)=0
 
virtual void v_GetUFromChars (NekDouble &u, const NekDouble &W1, const NekDouble &W2)=0
 
virtual void v_GetCharIntegral (NekDouble &I, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)=0
 
virtual 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)=0
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::PulseWavePressureArea
Array< OneD, MultiRegions::ExpListSharedPtrm_vessels
 
LibUtilities::SessionReaderSharedPtr m_session
 
NekDouble m_PExt
 
NekDouble m_rho
 

Detailed Description

Definition at line 52 of file EmpiricalPressureArea.h.

Constructor & Destructor Documentation

◆ EmpiricalPressureArea()

Nektar::EmpiricalPressureArea::EmpiricalPressureArea ( Array< OneD, MultiRegions::ExpListSharedPtr pVessel,
const LibUtilities::SessionReaderSharedPtr  pSession 
)

Definition at line 48 of file EmpiricalPressureArea.cpp.

51 : PulseWavePressureArea(pVessel, pSession)
52{
53}
PulseWavePressureArea(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession)

◆ ~EmpiricalPressureArea()

Nektar::EmpiricalPressureArea::~EmpiricalPressureArea ( )
override

Definition at line 55 of file EmpiricalPressureArea.cpp.

56{
57}

Member Function Documentation

◆ create()

static PulseWavePressureAreaSharedPtr Nektar::EmpiricalPressureArea::create ( Array< OneD, MultiRegions::ExpListSharedPtr > &  pVessel,
const LibUtilities::SessionReaderSharedPtr pSession 
)
inlinestatic

Definition at line 56 of file EmpiricalPressureArea.h.

59 {
61 pVessel, pSession);
62 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ GetKappa()

void Nektar::EmpiricalPressureArea::GetKappa ( NekDouble kappa,
const NekDouble A,
const NekDouble A0,
const NekDouble alpha = 0.5 
)
protected

Definition at line 284 of file EmpiricalPressureArea.cpp.

287{
288 kappa = 1 - alpha * log(A / A0);
289}
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:303

References tinysimd::log().

Referenced by v_GetC(), v_GetCharIntegral(), and v_GetPressure().

◆ v_GetAFromChars()

void Nektar::EmpiricalPressureArea::v_GetAFromChars ( NekDouble A,
const NekDouble W1,
const NekDouble W2,
const NekDouble beta,
const NekDouble A0,
const NekDouble alpha = 0.5 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 103 of file EmpiricalPressureArea.cpp.

108{
109 NekDouble xi = (W1 - W2) * sqrt(m_rho / (8 * beta * sqrt(A0)));
110
111 A = A0 * exp(xi * (2 - alpha * xi));
112}
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::LibUtilities::beta, Nektar::PulseWavePressureArea::m_rho, and tinysimd::sqrt().

◆ v_GetC()

void Nektar::EmpiricalPressureArea::v_GetC ( NekDouble c,
const NekDouble beta,
const NekDouble A,
const NekDouble A0,
const NekDouble alpha = 0.5 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 73 of file EmpiricalPressureArea.cpp.

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}
void GetKappa(NekDouble &kappa, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)

References Nektar::LibUtilities::beta, GetKappa(), Nektar::PulseWavePressureArea::m_rho, and tinysimd::sqrt().

◆ v_GetCharIntegral()

void Nektar::EmpiricalPressureArea::v_GetCharIntegral ( NekDouble I,
const NekDouble beta,
const NekDouble A,
const NekDouble A0,
const NekDouble alpha = 0.5 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 120 of file EmpiricalPressureArea.cpp.

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}

References Nektar::LibUtilities::beta, GetKappa(), Nektar::PulseWavePressureArea::m_rho, and tinysimd::sqrt().

◆ v_GetJacobianInverse()

void Nektar::EmpiricalPressureArea::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 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 132 of file EmpiricalPressureArea.cpp.

137{
138 // General formulation
139 if (type == "Bifurcation")
140 {
141 NekMatrix<NekDouble> J(6, 6);
142 Array<OneD, NekDouble> c(3);
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);
197 Array<OneD, NekDouble> c(3);
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);
252 Array<OneD, NekDouble> c(2);
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}
void GetC(NekDouble &c, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)

References Nektar::LibUtilities::beta, and Nektar::PulseWavePressureArea::GetC().

◆ v_GetPressure()

void Nektar::EmpiricalPressureArea::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 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 59 of file EmpiricalPressureArea.cpp.

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}
@ P
Monomial polynomials .
Definition: BasisType.h:62

References Nektar::LibUtilities::beta, GetKappa(), tinysimd::log(), Nektar::PulseWavePressureArea::m_PExt, Nektar::LibUtilities::P, and tinysimd::sqrt().

◆ v_GetUFromChars()

void Nektar::EmpiricalPressureArea::v_GetUFromChars ( NekDouble u,
const NekDouble W1,
const NekDouble W2 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 114 of file EmpiricalPressureArea.cpp.

116{
117 u = (W1 + W2) / 2;
118}

◆ v_GetW1()

void Nektar::EmpiricalPressureArea::v_GetW1 ( NekDouble W1,
const NekDouble u,
const NekDouble beta,
const NekDouble A,
const NekDouble A0,
const NekDouble alpha = 0.5 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 83 of file EmpiricalPressureArea.cpp.

86{
87 NekDouble I = 0.0;
88 GetCharIntegral(I, beta, A, A0, alpha);
89
90 W1 = u + I; // Elastic and assumes u0 = 0
91}
void GetCharIntegral(NekDouble &I, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)

References Nektar::LibUtilities::beta, and Nektar::PulseWavePressureArea::GetCharIntegral().

◆ v_GetW2()

void Nektar::EmpiricalPressureArea::v_GetW2 ( NekDouble W2,
const NekDouble u,
const NekDouble beta,
const NekDouble A,
const NekDouble A0,
const NekDouble alpha = 0.5 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 93 of file EmpiricalPressureArea.cpp.

96{
97 NekDouble I = 0.0;
98 GetCharIntegral(I, beta, A, A0, alpha);
99
100 W2 = u - I; // Elastic and assumes u0 = 0
101}

References Nektar::LibUtilities::beta, and Nektar::PulseWavePressureArea::GetCharIntegral().

Member Data Documentation

◆ className

std::string Nektar::EmpiricalPressureArea::className
static
Initial value:
=
"Empirical pressure area relationship for the arterial system")
static PulseWavePressureAreaSharedPtr create(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
PressureAreaFactory & GetPressureAreaFactory()

Definition at line 65 of file EmpiricalPressureArea.h.