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

#include <PowerPressureArea.h>

Inheritance diagram for Nektar::PowerPressureArea:
[legend]

Public Member Functions

 PowerPressureArea (Array< OneD, MultiRegions::ExpListSharedPtr > pVessel, const LibUtilities::SessionReaderSharedPtr pSession)
 
virtual ~PowerPressureArea ()
 
- 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

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)
 
virtual void v_GetC (NekDouble &c, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
 
virtual void v_GetW1 (NekDouble &W1, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
 
virtual void v_GetW2 (NekDouble &W2, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
 
virtual void v_GetAFromChars (NekDouble &A, const NekDouble &W1, const NekDouble &W2, const NekDouble &beta, const NekDouble &A0, const NekDouble &alpha=0.5)
 
virtual void v_GetUFromChars (NekDouble &u, const NekDouble &W1, const NekDouble &W2)
 
virtual void v_GetCharIntegral (NekDouble &I, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
 
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)
 
virtual void GetC0 (NekDouble &c0, const NekDouble &beta, const NekDouble &A0)
 
virtual void GetB (NekDouble &b, const NekDouble &c0)
 

Private Attributes

NekDouble P_Collapse
 

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 PowerPressureArea.h.

Constructor & Destructor Documentation

◆ PowerPressureArea()

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

Definition at line 47 of file PowerPressureArea.cpp.

50  : PulseWavePressureArea(pVessel, pSession)
51 {
52  m_session->LoadParameter("P_Collapse", P_Collapse, -13.3322); // -10mmHg converted to kg / (cm s^2)
53 }
LibUtilities::SessionReaderSharedPtr m_session
PulseWavePressureArea(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession)

References Nektar::PulseWavePressureArea::m_session, and P_Collapse.

◆ ~PowerPressureArea()

Nektar::PowerPressureArea::~PowerPressureArea ( )
virtual

Definition at line 55 of file PowerPressureArea.cpp.

56 {
57 }

Member Function Documentation

◆ create()

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

Definition at line 57 of file PowerPressureArea.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().

◆ GetB()

void Nektar::PowerPressureArea::GetB ( NekDouble b,
const NekDouble c0 
)
protectedvirtual

◆ GetC0()

void Nektar::PowerPressureArea::GetC0 ( NekDouble c0,
const NekDouble beta,
const NekDouble A0 
)
protectedvirtual

Definition at line 290 of file PowerPressureArea.cpp.

292 {
293  // Reference c0 from the beta law
294  c0 = sqrt(beta * sqrt(A0) / (2 * m_rho));
295 
296  /*
297  // Empirical approximation from Olufsen et al (1999)
298  NekDouble k1 = 3E3;
299  NekDouble k2 = -9;
300  NekDouble k3 = 337;
301  NekDouble PI = 3.14159265359;
302 
303  NekDouble R0 = sqrt(A0 / PI);
304 
305  c0 = sqrt(2 / (3 * m_rho) * (k1 * exp(k2 * R0) + k3));
306  */
307 }
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267

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

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

◆ v_GetAFromChars()

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

Implements Nektar::PulseWavePressureArea.

Definition at line 105 of file PowerPressureArea.cpp.

108 {
109  NekDouble c0 = 0.0;
110  GetC0(c0, beta, A0);
111 
112  NekDouble b = 0.0;
113  GetB(b, c0);
114 
115  A = A0 * pow(((b / (8 * c0)) * (W1 - W2) + 1), 4 / b);
116 }
virtual void GetC0(NekDouble &c0, const NekDouble &beta, const NekDouble &A0)
virtual void GetB(NekDouble &b, const NekDouble &c0)
double NekDouble

References GetB(), and GetC0().

◆ v_GetC()

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

Implements Nektar::PulseWavePressureArea.

Definition at line 73 of file PowerPressureArea.cpp.

75 {
76  NekDouble c0 = 0.0;
77  GetC0(c0, beta, A0);
78 
79  NekDouble b = 0.0;
80  GetB(b, c0);
81 
82  c = c0 * pow((A / A0), b / 4); // Elastic
83 }

References GetB(), and GetC0().

◆ v_GetCharIntegral()

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

Implements Nektar::PulseWavePressureArea.

Definition at line 124 of file PowerPressureArea.cpp.

126 {
127  NekDouble c = 0.0;
128  NekDouble c0 = 0.0;
129 
130  GetC0(c0, beta, A0);
131  GetC(c, beta, A, A0);
132 
133  NekDouble b = 0.0;
134  GetB(b, c0);
135 
136  I = (4 / b) * (c - c0);
137 }
void GetC(NekDouble &c, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)

References GetB(), Nektar::PulseWavePressureArea::GetC(), and GetC0().

◆ v_GetJacobianInverse()

void Nektar::PowerPressureArea::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 
)
protectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 139 of file PowerPressureArea.cpp.

143 {
144  // General formulation
145  if (type == "Bifurcation")
146  {
147  NekMatrix<NekDouble> J(6, 6);
148  Array<OneD, NekDouble> c(3);
149 
150  for (int i = 0; i < 3; ++i)
151  {
152  GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
153  }
154 
155  J.SetValue(0, 0, 1);
156  J.SetValue(0, 1, 0);
157  J.SetValue(0, 2, 0);
158  J.SetValue(0, 3, c[0] / Au[0]);
159  J.SetValue(0, 4, 0);
160  J.SetValue(0, 5, 0);
161 
162  J.SetValue(1, 0, 0);
163  J.SetValue(1, 1, 1);
164  J.SetValue(1, 2, 0);
165  J.SetValue(1, 3, 0);
166  J.SetValue(1, 4, -c[1] / Au[1]);
167  J.SetValue(1, 5, 0);
168 
169  J.SetValue(2, 0, 0);
170  J.SetValue(2, 1, 0);
171  J.SetValue(2, 2, 1);
172  J.SetValue(2, 3, 0);
173  J.SetValue(2, 4, 0);
174  J.SetValue(2, 5, -c[2] / Au[2]);
175 
176  J.SetValue(3, 0, Au[0]);
177  J.SetValue(3, 1, -Au[1]);
178  J.SetValue(3, 2, -Au[2]);
179  J.SetValue(3, 3, uu[0]);
180  J.SetValue(3, 4, -uu[1]);
181  J.SetValue(3, 5, -uu[2]);
182 
183  J.SetValue(4, 0, 2 * uu[0]);
184  J.SetValue(4, 1, -2 * uu[1]);
185  J.SetValue(4, 2, 0);
186  J.SetValue(4, 3, 2 * c[0] * c[0] / Au[0]);
187  J.SetValue(4, 4, -2 * c[1] * c[1] / Au[1]);
188  J.SetValue(4, 5, 0);
189 
190  J.SetValue(5, 0, 2 * uu[0]);
191  J.SetValue(5, 1, 0);
192  J.SetValue(5, 2, -2 * uu[2]);
193  J.SetValue(5, 3, 2 * c[0] * c[0] / Au[0]);
194  J.SetValue(5, 4, 0);
195  J.SetValue(5, 5, -2 * c[2] * c[2] / Au[2]);
196 
197  invJ = J;
198  invJ.Invert();
199  }
200  else if (type == "Merge")
201  {
202  NekMatrix<NekDouble> J(6, 6);
203  Array<OneD, NekDouble> c(3);
204 
205  for (int i = 0; i < 3; ++i)
206  {
207  GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
208  }
209 
210  J.SetValue(0, 0, 1);
211  J.SetValue(0, 1, 0);
212  J.SetValue(0, 2, 0);
213  J.SetValue(0, 3, -c[0] / Au[0]);
214  J.SetValue(0, 4, 0);
215  J.SetValue(0, 5, 0);
216 
217  J.SetValue(1, 0, 0);
218  J.SetValue(1, 1, 1);
219  J.SetValue(1, 2, 0);
220  J.SetValue(1, 3, 0);
221  J.SetValue(1, 4, c[1] / Au[1]);
222  J.SetValue(1, 5, 0);
223 
224  J.SetValue(2, 0, 0);
225  J.SetValue(2, 1, 0);
226  J.SetValue(2, 2, 1);
227  J.SetValue(2, 3, 0);
228  J.SetValue(2, 4, 0);
229  J.SetValue(2, 5, c[2] / Au[2]);
230 
231  J.SetValue(3, 0, Au[0]);
232  J.SetValue(3, 1, -Au[1]);
233  J.SetValue(3, 2, -Au[2]);
234  J.SetValue(3, 3, uu[0]);
235  J.SetValue(3, 4, -uu[1]);
236  J.SetValue(3, 5, -uu[2]);
237 
238  J.SetValue(4, 0, 2 * uu[0]);
239  J.SetValue(4, 1, -2 * uu[1]);
240  J.SetValue(4, 2, 0);
241  J.SetValue(4, 3, 2 * c[0] * c[0] / Au[0]);
242  J.SetValue(4, 4, -2 * c[1] * c[1] / Au[1]);
243  J.SetValue(4, 5, 0);
244 
245  J.SetValue(5, 0, 2 * uu[0]);
246  J.SetValue(5, 1, 0);
247  J.SetValue(5, 2, -2 * uu[2]);
248  J.SetValue(5, 3, 2 * c[0] * c[0] / Au[0]);
249  J.SetValue(5, 4, 0);
250  J.SetValue(5, 5, -2 * c[2] * c[2] / Au[2]);
251 
252  invJ = J;
253  invJ.Invert();
254  }
255  else if (type == "Junction")
256  {
257  NekMatrix<NekDouble> J(4, 4);
258  Array<OneD, NekDouble> c(2);
259 
260  for (int i = 0; i < 2; ++i)
261  {
262  GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
263  }
264 
265  J.SetValue(0, 0, 1);
266  J.SetValue(0, 1, 0);
267  J.SetValue(0, 2, c[0] / Au[0]);
268  J.SetValue(0, 3, 0);
269 
270  J.SetValue(1, 0, 0);
271  J.SetValue(1, 1, 1);
272  J.SetValue(1, 2, 0);
273  J.SetValue(1, 3, -c[1] / Au[1]);
274 
275  J.SetValue(2, 0, Au[0]);
276  J.SetValue(2, 1, -Au[1]);
277  J.SetValue(2, 2, uu[0]);
278  J.SetValue(2, 3, -uu[1]);
279 
280  J.SetValue(3, 0, 2 * uu[0]);
281  J.SetValue(3, 1, -2 * uu[1]);
282  J.SetValue(3, 2, 2 * c[0] * c[0] / Au[0]);
283  J.SetValue(3, 3, -2 * c[1] * c[1] / Au[1]);
284 
285  invJ = J;
286  invJ.Invert();
287  }
288 }

References Nektar::PulseWavePressureArea::GetC().

◆ v_GetPressure()

void Nektar::PowerPressureArea::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 
)
protectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 59 of file PowerPressureArea.cpp.

62 {
63  NekDouble c0 = 0.0;
64  GetC0(c0, beta, A0);
65 
66  NekDouble b = 0.0;
67  GetB(b, c0);
68 
69  P = m_PExt + (2 * m_rho * c0 * c0 / b) * (pow((A / A0), b / 2) - 1) // Power law by Smith/Canic/Mynard
70  - A0 * gamma * dAUdx / (A * sqrt(A)); // Viscoelasticity
71 }
P
Definition: main.py:133

References GetB(), GetC0(), Nektar::PulseWavePressureArea::m_PExt, Nektar::PulseWavePressureArea::m_rho, main::P, and tinysimd::sqrt().

◆ v_GetUFromChars()

void Nektar::PowerPressureArea::v_GetUFromChars ( NekDouble u,
const NekDouble W1,
const NekDouble W2 
)
protectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 118 of file PowerPressureArea.cpp.

120 {
121  u = (W1 + W2) / 2;
122 }

◆ v_GetW1()

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

Implements Nektar::PulseWavePressureArea.

Definition at line 85 of file PowerPressureArea.cpp.

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

References Nektar::PulseWavePressureArea::GetCharIntegral().

◆ v_GetW2()

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

Implements Nektar::PulseWavePressureArea.

Definition at line 95 of file PowerPressureArea.cpp.

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

References Nektar::PulseWavePressureArea::GetCharIntegral().

Member Data Documentation

◆ className

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

Definition at line 65 of file PowerPressureArea.h.

◆ P_Collapse

NekDouble Nektar::PowerPressureArea::P_Collapse
private

Definition at line 109 of file PowerPressureArea.h.

Referenced by GetB(), and PowerPressureArea().