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

#include <TemplatePressureArea.h>

Inheritance diagram for Nektar::TemplatePressureArea:
[legend]

Public Member Functions

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

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 50 of file TemplatePressureArea.h.

Constructor & Destructor Documentation

◆ TemplatePressureArea()

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

Definition at line 48 of file TemplatePressureArea.cpp.

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

◆ ~TemplatePressureArea()

Nektar::TemplatePressureArea::~TemplatePressureArea ( )
virtual

Definition at line 55 of file TemplatePressureArea.cpp.

56 {
57 }

Member Function Documentation

◆ create()

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

Definition at line 54 of file TemplatePressureArea.h.

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

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

◆ v_GetAFromChars()

void Nektar::TemplatePressureArea::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 97 of file TemplatePressureArea.cpp.

100 {
101  /*
102  If an anayltical formula exists for this, then it should be used here. If
103  not, a numerical method is used.
104  */
105 
106  // Iterative Newton-Raphson
107  NekDouble I = 0.0;
108  NekDouble c = 0.0;
109  NekDouble FA = 0.0;
110  NekDouble dFdA = 0.0;
111  NekDouble dA = 0.0;
112 
113  int proceed = 1;
114  int iter = 0;
115  int MAX_ITER = 200;
116  NekDouble Tol = 1.0E-10;
117 
118  while ((proceed) && (iter < MAX_ITER))
119  {
120  iter += 1;
121 
122  GetCharIntegral(I, beta, A, A0, alpha);
123  GetC(c, beta, A, A0, alpha);
124 
125  FA = I - ((W1 - W2) / 2);
126  dFdA = c / A;
127 
128  dA = FA / dFdA;
129  A -= dA;
130 
131  if (sqrt(dA * dA) < Tol)
132  {
133  proceed = 0;
134  }
135  }
136 
137 }
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)
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267

References Nektar::PulseWavePressureArea::GetC(), Nektar::PulseWavePressureArea::GetCharIntegral(), and tinysimd::sqrt().

◆ v_GetC()

void Nektar::TemplatePressureArea::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 71 of file TemplatePressureArea.cpp.

73 {
74  // PWV here
75 }

◆ v_GetCharIntegral()

void Nektar::TemplatePressureArea::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 145 of file TemplatePressureArea.cpp.

147 {
148  /*
149  If an anayltical formula exists for this, then it should be used here. If
150  not, a numerical method is used.
151  */
152 
153  // Numerical integration via the trapezium rule
154  // f(x) = c / A, x1 = A0, x2 = A
155 
156  int n = 500;
157  NekDouble h = (A - A0) / n;
158  NekDouble A_n = A0;
159  NekDouble c_n = 0.0;
160 
161  for (int i = 1; i < n; ++i)
162  {
163  A_n += h;
164  GetC(c_n, beta, A_n, A0, alpha);
165  I += c_n / A_n;
166  }
167 
168  NekDouble c = 0.0;
169  NekDouble c0 = 0.0;
170  GetC(c, beta, A, A0, alpha);
171  GetC(c0, beta, A0, A0, alpha);
172 
173  I += ((c / A) + (c0 / A0)) / 2;
174  I *= h;
175 }

References Nektar::PulseWavePressureArea::GetC().

◆ v_GetJacobianInverse()

void Nektar::TemplatePressureArea::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 177 of file TemplatePressureArea.cpp.

181 {
182  /*
183  This is a general method that will work for all tube laws. Simulations will
184  run faster if specific analytical formulae are derived as in the case of the
185  beta law, but it's not necessary.
186  */
187 
188  // General formulation
189  if (type == "Bifurcation")
190  {
191  NekMatrix<NekDouble> J(6, 6);
192  Array<OneD, NekDouble> c(3);
193 
194  for (int i = 0; i < 3; ++i)
195  {
196  GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
197  }
198 
199  J.SetValue(0, 0, 1);
200  J.SetValue(0, 1, 0);
201  J.SetValue(0, 2, 0);
202  J.SetValue(0, 3, c[0] / Au[0]);
203  J.SetValue(0, 4, 0);
204  J.SetValue(0, 5, 0);
205 
206  J.SetValue(1, 0, 0);
207  J.SetValue(1, 1, 1);
208  J.SetValue(1, 2, 0);
209  J.SetValue(1, 3, 0);
210  J.SetValue(1, 4, -c[1] / Au[1]);
211  J.SetValue(1, 5, 0);
212 
213  J.SetValue(2, 0, 0);
214  J.SetValue(2, 1, 0);
215  J.SetValue(2, 2, 1);
216  J.SetValue(2, 3, 0);
217  J.SetValue(2, 4, 0);
218  J.SetValue(2, 5, -c[2] / Au[2]);
219 
220  J.SetValue(3, 0, Au[0]);
221  J.SetValue(3, 1, -Au[1]);
222  J.SetValue(3, 2, -Au[2]);
223  J.SetValue(3, 3, uu[0]);
224  J.SetValue(3, 4, -uu[1]);
225  J.SetValue(3, 5, -uu[2]);
226 
227  J.SetValue(4, 0, 2 * uu[0]);
228  J.SetValue(4, 1, -2 * uu[1]);
229  J.SetValue(4, 2, 0);
230  J.SetValue(4, 3, 2 * c[0] * c[0] / Au[0]);
231  J.SetValue(4, 4, -2 * c[1] * c[1] / Au[1]);
232  J.SetValue(4, 5, 0);
233 
234  J.SetValue(5, 0, 2 * uu[0]);
235  J.SetValue(5, 1, 0);
236  J.SetValue(5, 2, -2 * uu[2]);
237  J.SetValue(5, 3, 2 * c[0] * c[0] / Au[0]);
238  J.SetValue(5, 4, 0);
239  J.SetValue(5, 5, -2 * c[2] * c[2] / Au[2]);
240 
241  invJ = J;
242  invJ.Invert();
243  }
244  else if (type == "Merge")
245  {
246  NekMatrix<NekDouble> J(6, 6);
247  Array<OneD, NekDouble> c(3);
248 
249  for (int i = 0; i < 3; ++i)
250  {
251  GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
252  }
253 
254  J.SetValue(0, 0, 1);
255  J.SetValue(0, 1, 0);
256  J.SetValue(0, 2, 0);
257  J.SetValue(0, 3, -c[0] / Au[0]);
258  J.SetValue(0, 4, 0);
259  J.SetValue(0, 5, 0);
260 
261  J.SetValue(1, 0, 0);
262  J.SetValue(1, 1, 1);
263  J.SetValue(1, 2, 0);
264  J.SetValue(1, 3, 0);
265  J.SetValue(1, 4, c[1] / Au[1]);
266  J.SetValue(1, 5, 0);
267 
268  J.SetValue(2, 0, 0);
269  J.SetValue(2, 1, 0);
270  J.SetValue(2, 2, 1);
271  J.SetValue(2, 3, 0);
272  J.SetValue(2, 4, 0);
273  J.SetValue(2, 5, c[2] / Au[2]);
274 
275  J.SetValue(3, 0, Au[0]);
276  J.SetValue(3, 1, -Au[1]);
277  J.SetValue(3, 2, -Au[2]);
278  J.SetValue(3, 3, uu[0]);
279  J.SetValue(3, 4, -uu[1]);
280  J.SetValue(3, 5, -uu[2]);
281 
282  J.SetValue(4, 0, 2 * uu[0]);
283  J.SetValue(4, 1, -2 * uu[1]);
284  J.SetValue(4, 2, 0);
285  J.SetValue(4, 3, 2 * c[0] * c[0] / Au[0]);
286  J.SetValue(4, 4, -2 * c[1] * c[1] / Au[1]);
287  J.SetValue(4, 5, 0);
288 
289  J.SetValue(5, 0, 2 * uu[0]);
290  J.SetValue(5, 1, 0);
291  J.SetValue(5, 2, -2 * uu[2]);
292  J.SetValue(5, 3, 2 * c[0] * c[0] / Au[0]);
293  J.SetValue(5, 4, 0);
294  J.SetValue(5, 5, -2 * c[2] * c[2] / Au[2]);
295 
296  invJ = J;
297  invJ.Invert();
298  }
299  else if (type == "Junction")
300  {
301  NekMatrix<NekDouble> J(4, 4);
302  Array<OneD, NekDouble> c(2);
303 
304  for (int i = 0; i < 2; ++i)
305  {
306  GetC(c[i], beta[i], Au[i], A0[i], alpha[i]);
307  }
308 
309  J.SetValue(0, 0, 1);
310  J.SetValue(0, 1, 0);
311  J.SetValue(0, 2, c[0] / Au[0]);
312  J.SetValue(0, 3, 0);
313 
314  J.SetValue(1, 0, 0);
315  J.SetValue(1, 1, 1);
316  J.SetValue(1, 2, 0);
317  J.SetValue(1, 3, -c[1] / Au[1]);
318 
319  J.SetValue(2, 0, Au[0]);
320  J.SetValue(2, 1, -Au[1]);
321  J.SetValue(2, 2, uu[0]);
322  J.SetValue(2, 3, -uu[1]);
323 
324  J.SetValue(3, 0, 2 * uu[0]);
325  J.SetValue(3, 1, -2 * uu[1]);
326  J.SetValue(3, 2, 2 * c[0] * c[0] / Au[0]);
327  J.SetValue(3, 3, -2 * c[1] * c[1] / Au[1]);
328 
329  invJ = J;
330  invJ.Invert();
331  }
332 }

References Nektar::PulseWavePressureArea::GetC().

◆ v_GetPressure()

void Nektar::TemplatePressureArea::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 TemplatePressureArea.cpp.

62 {
63  // Relation here
64 
65  /*
66  NOTE: to add a new relation, you also need to add the name of the .h and
67  .cpp files to the CMakeLists.txt file in the PulseWaveSolver folder.
68  */
69 }

◆ v_GetUFromChars()

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

Implements Nektar::PulseWavePressureArea.

Definition at line 139 of file TemplatePressureArea.cpp.

141 {
142  u = (W1 + W2) / 2; // Necessarily the case for all tube laws
143 }

◆ v_GetW1()

void Nektar::TemplatePressureArea::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 77 of file TemplatePressureArea.cpp.

80 {
81  NekDouble I = 0.0;
82  GetCharIntegral(I, beta, A, A0);
83 
84  W1 = u + I; // Elastic and assumes u0 = 0
85 }

References Nektar::PulseWavePressureArea::GetCharIntegral().

◆ v_GetW2()

void Nektar::TemplatePressureArea::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 87 of file TemplatePressureArea.cpp.

90 {
91  NekDouble I = 0.0;
92  GetCharIntegral(I, beta, A, A0);
93 
94  W2 = u - I; // Elastic and assumes u0 = 0
95 }

References Nektar::PulseWavePressureArea::GetCharIntegral().

Member Data Documentation

◆ className

std::string Nektar::TemplatePressureArea::className
static
Initial value:
=
"Template to add a new Pressure-Area relation")
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 61 of file TemplatePressureArea.h.