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

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 53 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 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 101 of file TemplatePressureArea.cpp.

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

References Nektar::LibUtilities::beta, 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 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 74 of file TemplatePressureArea.cpp.

77 {
78  // PWV here
79 }

◆ v_GetCharIntegral()

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

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

References Nektar::LibUtilities::beta, and 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 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 184 of file TemplatePressureArea.cpp.

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

References Nektar::LibUtilities::beta, and 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 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 59 of file TemplatePressureArea.cpp.

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

◆ v_GetUFromChars()

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

Implements Nektar::PulseWavePressureArea.

Definition at line 144 of file TemplatePressureArea.cpp.

146 {
147  u = (W1 + W2) / 2; // Necessarily the case for all tube laws
148 }

◆ 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 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 81 of file TemplatePressureArea.cpp.

84 {
85  NekDouble I = 0.0;
86  GetCharIntegral(I, beta, A, A0);
87 
88  W1 = u + I; // Elastic and assumes u0 = 0
89 }

References Nektar::LibUtilities::beta, and 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 
)
overrideprotectedvirtual

Implements Nektar::PulseWavePressureArea.

Definition at line 91 of file TemplatePressureArea.cpp.

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

References Nektar::LibUtilities::beta, and 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:198
static PulseWavePressureAreaSharedPtr create(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession)
PressureAreaFactory & GetPressureAreaFactory()

Definition at line 61 of file TemplatePressureArea.h.