oomph::QuarterPipeDomain Class Reference

Domain representing a quarter pipe. More...

#include <quarter_pipe_domain.h>

+ Inheritance diagram for oomph::QuarterPipeDomain:

Public Types

typedef double(* AxialSpacingFctPt) (const double &xi)
 

Public Member Functions

 QuarterPipeDomain (const unsigned &ntheta, const unsigned &nr, const unsigned &nz, const double &rmin, const double &rmax, const double &length)
 
 QuarterPipeDomain (const QuarterPipeDomain &)=delete
 Broken copy constructor. More...
 
void operator= (const QuarterPipeDomain &)=delete
 Broken assignment operator. More...
 
 ~QuarterPipeDomain ()
 Destructor: More...
 
AxialSpacingFctPtaxial_spacing_fct_pt ()
 
double axial_spacing_fct (const double &xi)
 
void macro_element_boundary (const unsigned &t, const unsigned &i_macro, const unsigned &i_direct, const Vector< double > &s, Vector< double > &f)
 
- Public Member Functions inherited from oomph::Domain
 Domain ()
 Constructor. More...
 
 Domain (const Domain &)=delete
 Broken copy constructor. More...
 
void operator= (const Domain &)=delete
 Broken assignment operator. More...
 
virtual ~Domain ()
 
MacroElementmacro_element_pt (const unsigned &i)
 Access to i-th macro element. More...
 
unsigned nmacro_element ()
 Number of macro elements in domain. More...
 
void output (const std::string &filename, const unsigned &nplot)
 Output macro elements. More...
 
void output (std::ostream &outfile, const unsigned &nplot)
 Output macro elements. More...
 
virtual void macro_element_boundary (const double &t, const unsigned &i_macro, const unsigned &i_direct, const Vector< double > &s, Vector< double > &f)
 
void macro_element_boundary (const unsigned &i_macro, const unsigned &i_direct, const Vector< double > &s, Vector< double > &f)
 
void output_macro_element_boundaries (const std::string &filename, const unsigned &nplot)
 Output all macro element boundaries as tecplot zones. More...
 
void output_macro_element_boundaries (std::ostream &outfile, const unsigned &nplot)
 Output all macro element boundaries as tecplot zones. More...
 
virtual void dmacro_element_boundary (const unsigned &t, const unsigned &i_macro, const unsigned &i_direct, const Vector< double > &s, Vector< double > &f)
 
virtual void dmacro_element_boundary (const double &t, const unsigned &i_macro, const unsigned &i_direct, const Vector< double > &s, Vector< double > &f)
 
void dmacro_element_boundary (const unsigned &i_macro, const unsigned &i_direct, const Vector< double > &s, Vector< double > &f)
 
virtual void d2macro_element_boundary (const unsigned &t, const unsigned &i_macro, const unsigned &i_direct, const Vector< double > &s, Vector< double > &f)
 
virtual void d2macro_element_boundary (const double &t, const unsigned &i_macro, const unsigned &i_direct, const Vector< double > &s, Vector< double > &f)
 
void d2macro_element_boundary (const unsigned &i_macro, const unsigned &i_direct, const Vector< double > &s, Vector< double > &f)
 

Private Member Functions

void r_U (const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
 Boundary of macro element zeta \( \in [-1,1]x[-1,1] \). More...
 
void r_L (const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
 Boundary of macro element zeta \( \in [-1,1]x[-1,1] \). More...
 
void r_D (const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
 Boundary of macro element zeta \( \in [-1,1]x[-1,1] \). More...
 
void r_R (const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
 Boundary of macro element zeta \( \in [-1,1]x[-1,1] \). More...
 
void r_F (const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
 Boundary of macro element zeta \( \in [-1,1]x[-1,1] \). More...
 
void r_B (const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
 Boundary of macro element zeta \( \in [-1,1]x[-1,1] \). More...
 

Static Private Member Functions

static double default_axial_spacing_fct (const double &xi)
 

Private Attributes

unsigned Ntheta
 Number of elements azimuthal direction. More...
 
unsigned Nr
 Number of elements radial direction. More...
 
unsigned Nz
 Number of elements axial direction. More...
 
double Rmin
 Inner radius. More...
 
double Rmax
 Outer radius. More...
 
double Length
 Length. More...
 
GeomObjectOuter_boundary_cross_section_pt
 
GeomObjectInner_boundary_cross_section_pt
 
AxialSpacingFctPt Axial_spacing_fct_pt
 

Additional Inherited Members

- Protected Attributes inherited from oomph::Domain
Vector< MacroElement * > Macro_element_pt
 Vector of pointers to macro elements. More...
 

Detailed Description

Domain representing a quarter pipe.

Member Typedef Documentation

◆ AxialSpacingFctPt

typedef double(* oomph::QuarterPipeDomain::AxialSpacingFctPt) (const double &xi)

Typedef for function pointer for function that implements axial spacing of macro elements

Constructor & Destructor Documentation

◆ QuarterPipeDomain() [1/2]

oomph::QuarterPipeDomain::QuarterPipeDomain ( const unsigned ntheta,
const unsigned nr,
const unsigned nz,
const double rmin,
const double rmax,
const double length 
)
inline

Constructor: Pass number of elements in various directions, the inner and outer radius and the length of the tube

51  : Ntheta(ntheta),
52  Nr(nr),
53  Nz(nz),
54  Rmin(rmin),
55  Rmax(rmax),
56  Length(length),
58  {
59  // Number of macroelements
60  unsigned nmacro = nr * ntheta * nz;
61 
62  // Resize
63  Macro_element_pt.resize(nmacro);
64 
65  // Create macro elements
66  for (unsigned i = 0; i < nmacro; i++)
67  {
68  Macro_element_pt[i] = new QMacroElement<3>(this, i);
69  }
70 
71  // Make geom object representing the outer and inner boundaries of
72  // the cross section
75  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Vector< MacroElement * > Macro_element_pt
Vector of pointers to macro elements.
Definition: domain.h:301
double Length
Length.
Definition: quarter_pipe_domain.h:135
unsigned Nr
Number of elements radial direction.
Definition: quarter_pipe_domain.h:123
GeomObject * Inner_boundary_cross_section_pt
Definition: quarter_pipe_domain.h:143
static double default_axial_spacing_fct(const double &xi)
Definition: quarter_pipe_domain.h:151
double Rmax
Outer radius.
Definition: quarter_pipe_domain.h:132
unsigned Nz
Number of elements axial direction.
Definition: quarter_pipe_domain.h:126
double Rmin
Inner radius.
Definition: quarter_pipe_domain.h:129
AxialSpacingFctPt Axial_spacing_fct_pt
Definition: quarter_pipe_domain.h:147
GeomObject * Outer_boundary_cross_section_pt
Definition: quarter_pipe_domain.h:139
unsigned Ntheta
Number of elements azimuthal direction.
Definition: quarter_pipe_domain.h:120
double rmax
Definition: time_harmonic_fourier_decomposed_linear_elasticity/cylinder/cylinder.cc:70
double rmin
Definition: time_harmonic_fourier_decomposed_linear_elasticity/cylinder/cylinder.cc:68
const unsigned nz
Definition: ConstraintElementsUnitTest.cpp:32

References i, Inner_boundary_cross_section_pt, oomph::Domain::Macro_element_pt, Mesh_Parameters::nz, Outer_boundary_cross_section_pt, Global_Parameters::rmax, and Global_Parameters::rmin.

◆ QuarterPipeDomain() [2/2]

oomph::QuarterPipeDomain::QuarterPipeDomain ( const QuarterPipeDomain )
delete

Broken copy constructor.

◆ ~QuarterPipeDomain()

oomph::QuarterPipeDomain::~QuarterPipeDomain ( )
inline

Destructor:

85  {
86  // Note: macro elements are cleaned up in base class...
89  }

References Inner_boundary_cross_section_pt, and Outer_boundary_cross_section_pt.

Member Function Documentation

◆ axial_spacing_fct()

double oomph::QuarterPipeDomain::axial_spacing_fct ( const double xi)
inline

Function that implements axial spacing of macro elements

105  {
106  return Axial_spacing_fct_pt(xi);
107  }

References Axial_spacing_fct_pt.

Referenced by macro_element_boundary().

◆ axial_spacing_fct_pt()

AxialSpacingFctPt& oomph::QuarterPipeDomain::axial_spacing_fct_pt ( )
inline

Function pointer for function that implements axial spacing of macro elements

98  {
99  return Axial_spacing_fct_pt;
100  }

References Axial_spacing_fct_pt.

◆ default_axial_spacing_fct()

static double oomph::QuarterPipeDomain::default_axial_spacing_fct ( const double xi)
inlinestaticprivate

Default for function that implements axial spacing of macro elements

152  {
153  return xi;
154  }

◆ macro_element_boundary()

void oomph::QuarterPipeDomain::macro_element_boundary ( const unsigned t,
const unsigned imacro,
const unsigned idirect,
const Vector< double > &  s,
Vector< double > &  f 
)
virtual

Vector representation of the i_macro-th macro element boundary i_direct (U/D/L/R/F/B) at time level t (t=0: present; t>0: previous): f(s).

Vector representation of the imacro-th macro element boundary idirect (U/D/L/R/F/B) at time level t: f(s)

Implements oomph::Domain.

235  {
236  using namespace OcTreeNames;
237 
238  const double pi = MathematicalConstants::Pi;
239 
240  // Match the elements number with the position
241  unsigned num_z = imacro / (Nr * Ntheta);
242  unsigned num_y = (imacro % (Nr * Ntheta)) / Ntheta;
243  unsigned num_x = imacro % Ntheta;
244 
245  // Define the extreme coordinates
246 
247  // radial direction
248  double rmin = Rmin + (Rmax - Rmin) * double(num_y) / double(Nr);
249  double rmax = Rmin + (Rmax - Rmin) * double(num_y + 1) / double(Nr);
250 
251  // theta direction
252  double thetamin = (pi / 2.0) * (1.0 - double(num_x + 1) / double(Ntheta));
253  double thetamax = (pi / 2.0) * (1.0 - double(num_x) / double(Ntheta));
254 
255  // zdirection (tube)
256  double zmin = double(num_z) * Length / double(Nz);
257  double zmax = double(num_z + 1) * Length / double(Nz);
258 
259 
260  // Which direction?
261  if (idirect == U)
262  {
263  r_U(t, s, f, rmin, rmax, thetamin, thetamax, zmin, zmax);
264  }
265  else if (idirect == D)
266  {
267  r_D(t, s, f, rmin, rmax, thetamin, thetamax, zmin, zmax);
268  }
269  else if (idirect == L)
270  {
271  r_L(t, s, f, rmin, rmax, thetamin, thetamax, zmin, zmax);
272  }
273  else if (idirect == R)
274  {
275  r_R(t, s, f, rmin, rmax, thetamin, thetamax, zmin, zmax);
276  }
277  else if (idirect == F)
278  {
279  r_F(t, s, f, rmin, rmax, thetamin, thetamax, zmin, zmax);
280  }
281  else if (idirect == B)
282  {
283  r_B(t, s, f, rmin, rmax, thetamin, thetamax, zmin, zmax);
284  }
285  else
286  {
287  std::ostringstream error_stream;
288  error_stream << "idirect is " << idirect << " not one of U, D, L, R, F, B"
289  << std::endl;
290 
291  throw OomphLibError(
292  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
293  }
294 
295  // Now redistribute points in the axial direction
296  double z_frac = f[2] / Length;
297  f[2] = Length * axial_spacing_fct(z_frac);
298  }
dominoes D
Definition: Domino.cpp:55
MatrixXd L
Definition: LLT_example.cpp:6
@ R
Definition: StatisticsVector.h:21
Definition: matrices.h:74
double axial_spacing_fct(const double &xi)
Definition: quarter_pipe_domain.h:104
void r_L(const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
Boundary of macro element zeta .
Definition: quarter_pipe_domain.h:304
void r_U(const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
Boundary of macro element zeta .
Definition: quarter_pipe_domain.h:421
void r_F(const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
Boundary of macro element zeta .
Definition: quarter_pipe_domain.h:455
void r_R(const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
Boundary of macro element zeta .
Definition: quarter_pipe_domain.h:345
void r_B(const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
Boundary of macro element zeta .
Definition: quarter_pipe_domain.h:490
void r_D(const unsigned &t, const Vector< double > &zeta, Vector< double > &f, const double &rmin, const double &rmax, const double &thetamin, const double &thetamax, const double &zmin, const double &zmax)
Boundary of macro element zeta .
Definition: quarter_pipe_domain.h:387
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
RealScalar s
Definition: level1_cplx_impl.h:130
double zmax
Definition: time_harmonic_fourier_decomposed_linear_elasticity/cylinder/cylinder.cc:71
double zmin
Definition: time_harmonic_fourier_decomposed_linear_elasticity/cylinder/cylinder.cc:69
double U
Swimming speed.
Definition: two_d_variable_diff_adapt.cc:53
const Mdouble pi
Definition: ExtendedMath.h:23
const double Pi
50 digits from maple
Definition: oomph_utilities.h:157
@ F
Definition: octree.h:74
t
Definition: plotPSD.py:36
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References axial_spacing_fct(), D, f(), oomph::OcTreeNames::F, L, Length, Nr, Ntheta, Nz, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, constants::pi, oomph::MathematicalConstants::Pi, R, r_B(), r_D(), r_F(), r_L(), r_R(), r_U(), Global_Parameters::rmax, Rmax, Global_Parameters::rmin, Rmin, s, plotPSD::t, RachelsAdvectionDiffusion::U, Global_Parameters::zmax, and Global_Parameters::zmin.

◆ operator=()

void oomph::QuarterPipeDomain::operator= ( const QuarterPipeDomain )
delete

Broken assignment operator.

◆ r_B()

void oomph::QuarterPipeDomain::r_B ( const unsigned t,
const Vector< double > &  zeta,
Vector< double > &  f,
const double rmin,
const double rmax,
const double thetamin,
const double thetamax,
const double zmin,
const double zmax 
)
private

Boundary of macro element zeta \( \in [-1,1]x[-1,1] \).

Back face of a macro element \( s \in [-1,1]*[-1,1] \).

499  {
500  Vector<double> x(1);
501  x[0] = thetamax + (0.5 * (s[0] + 1.0)) * (thetamin - thetamax);
502 
503  // Point on outer wall
504  Vector<double> r_outer(2);
506 
507  // Point on inner wall
508  Vector<double> r_inner(2);
510 
511  // Get layer
512  double rad = rmin + (0.5 * (s[1] + 1.0)) * (rmax - rmin);
513  for (unsigned i = 0; i < 2; i++)
514  {
515  f[i] =
516  r_inner[i] + (r_outer[i] - r_inner[i]) * (rad - Rmin) / (Rmax - Rmin);
517  }
518  f[2] = zmin;
519  }
virtual void position(const Vector< double > &zeta, Vector< double > &r) const =0
Parametrised position on object at current time: r(zeta).
list x
Definition: plotDoE.py:28

References f(), i, Inner_boundary_cross_section_pt, Outer_boundary_cross_section_pt, oomph::GeomObject::position(), Global_Parameters::rmax, Rmax, Global_Parameters::rmin, Rmin, s, plotPSD::t, plotDoE::x, and Global_Parameters::zmin.

Referenced by macro_element_boundary().

◆ r_D()

void oomph::QuarterPipeDomain::r_D ( const unsigned t,
const Vector< double > &  zeta,
Vector< double > &  f,
const double rmin,
const double rmax,
const double thetamin,
const double thetamax,
const double zmin,
const double zmax 
)
private

Boundary of macro element zeta \( \in [-1,1]x[-1,1] \).

Left face of a macro element \(s \in [-1,1]*[-1,1] \).

396  {
397  Vector<double> x(1);
398  x[0] = thetamax + (0.5 * (s[0] + 1.0)) * (thetamin - thetamax);
399 
400  // Point on outer wall
401  Vector<double> r_outer(2);
403 
404  // Point on inner wall
405  Vector<double> r_inner(2);
407 
408  // Get layer
409  for (unsigned i = 0; i < 2; i++)
410  {
411  f[i] =
412  r_inner[i] + (r_outer[i] - r_inner[i]) * (rmin - Rmin) / (Rmax - Rmin);
413  }
414  f[2] = zmin + (zmax - zmin) * (0.5 * (s[1] + 1.0));
415  }

References f(), i, Inner_boundary_cross_section_pt, Outer_boundary_cross_section_pt, oomph::GeomObject::position(), Rmax, Global_Parameters::rmin, Rmin, s, plotPSD::t, plotDoE::x, Global_Parameters::zmax, and Global_Parameters::zmin.

Referenced by macro_element_boundary().

◆ r_F()

void oomph::QuarterPipeDomain::r_F ( const unsigned t,
const Vector< double > &  zeta,
Vector< double > &  f,
const double rmin,
const double rmax,
const double thetamin,
const double thetamax,
const double zmin,
const double zmax 
)
private

Boundary of macro element zeta \( \in [-1,1]x[-1,1] \).

Front face of a macro element \( s \in [-1,1]*[-1,1] \).

464  {
465  Vector<double> x(1);
466  x[0] = thetamax + (0.5 * (s[0] + 1.0)) * (thetamin - thetamax);
467 
468  // Point on outer wall
469  Vector<double> r_outer(2);
471 
472  // Point on inner wall
473  Vector<double> r_inner(2);
475 
476  // Get layer
477  double rad = rmin + (0.5 * (s[1] + 1.0)) * (rmax - rmin);
478  for (unsigned i = 0; i < 2; i++)
479  {
480  f[i] =
481  r_inner[i] + (r_outer[i] - r_inner[i]) * (rad - Rmin) / (Rmax - Rmin);
482  }
483  f[2] = zmax;
484  }

References f(), i, Inner_boundary_cross_section_pt, Outer_boundary_cross_section_pt, oomph::GeomObject::position(), Global_Parameters::rmax, Rmax, Global_Parameters::rmin, Rmin, s, plotPSD::t, plotDoE::x, and Global_Parameters::zmax.

Referenced by macro_element_boundary().

◆ r_L()

void oomph::QuarterPipeDomain::r_L ( const unsigned t,
const Vector< double > &  zeta,
Vector< double > &  f,
const double rmin,
const double rmax,
const double thetamin,
const double thetamax,
const double zmin,
const double zmax 
)
private

Boundary of macro element zeta \( \in [-1,1]x[-1,1] \).

Left face of a macro element \( s \in [-1,1]*[-1,1] \).

313  {
314  Vector<double> x(1);
315  x[0] = thetamax;
316 
317  // Point on outer wall
318  Vector<double> r_outer(2);
320 
321  // Point on inner wall
322  Vector<double> r_inner(2);
324 
325  // Get layer boundaries
326  Vector<double> r_top(2);
327  Vector<double> r_bot(2);
328  for (unsigned i = 0; i < 2; i++)
329  {
330  r_top[i] =
331  r_inner[i] + (r_outer[i] - r_inner[i]) * (rmax - Rmin) / (Rmax - Rmin);
332  r_bot[i] =
333  r_inner[i] + (r_outer[i] - r_inner[i]) * (rmin - Rmin) / (Rmax - Rmin);
334  }
335 
336  // Compute coordinates
337  f[0] = r_bot[0] + (0.5 * (s[0] + 1.0)) * (r_top[0] - r_bot[0]);
338  f[1] = r_bot[1] + (0.5 * (s[0] + 1.0)) * (r_top[1] - r_bot[1]);
339  f[2] = zmin + (zmax - zmin) * (0.5 * (s[1] + 1.0));
340  }

References f(), i, Inner_boundary_cross_section_pt, Outer_boundary_cross_section_pt, oomph::GeomObject::position(), Global_Parameters::rmax, Rmax, Global_Parameters::rmin, Rmin, s, plotPSD::t, plotDoE::x, Global_Parameters::zmax, and Global_Parameters::zmin.

Referenced by macro_element_boundary().

◆ r_R()

void oomph::QuarterPipeDomain::r_R ( const unsigned t,
const Vector< double > &  zeta,
Vector< double > &  f,
const double rmin,
const double rmax,
const double thetamin,
const double thetamax,
const double zmin,
const double zmax 
)
private

Boundary of macro element zeta \( \in [-1,1]x[-1,1] \).

Right face of a macro element \( s \in [-1,1]*[-1,1] \).

354  {
355  Vector<double> x(1);
356  x[0] = thetamin;
357 
358  // Point on outer wall
359  Vector<double> r_outer(2);
361 
362  // Point on inner wall
363  Vector<double> r_inner(2);
365 
366  // Get layer boundaries
367  Vector<double> r_top(2);
368  Vector<double> r_bot(2);
369  for (unsigned i = 0; i < 2; i++)
370  {
371  r_top[i] =
372  r_inner[i] + (r_outer[i] - r_inner[i]) * (rmax - Rmin) / (Rmax - Rmin);
373  r_bot[i] =
374  r_inner[i] + (r_outer[i] - r_inner[i]) * (rmin - Rmin) / (Rmax - Rmin);
375  }
376 
377  // Compute coordinates
378  f[0] = r_bot[0] + (0.5 * (s[0] + 1.0)) * (r_top[0] - r_bot[0]);
379  f[1] = r_bot[1] + (0.5 * (s[0] + 1.0)) * (r_top[1] - r_bot[1]);
380  f[2] = zmin + (zmax - zmin) * (0.5 * (s[1] + 1.0));
381  }

References f(), i, Inner_boundary_cross_section_pt, Outer_boundary_cross_section_pt, oomph::GeomObject::position(), Global_Parameters::rmax, Rmax, Global_Parameters::rmin, Rmin, s, plotPSD::t, plotDoE::x, Global_Parameters::zmax, and Global_Parameters::zmin.

Referenced by macro_element_boundary().

◆ r_U()

void oomph::QuarterPipeDomain::r_U ( const unsigned t,
const Vector< double > &  zeta,
Vector< double > &  f,
const double rmin,
const double rmax,
const double thetamin,
const double thetamax,
const double zmin,
const double zmax 
)
private

Boundary of macro element zeta \( \in [-1,1]x[-1,1] \).

Right face of a macro element \( s \in [-1,1]*[-1,1] \).

430  {
431  Vector<double> x(1);
432  x[0] = thetamax + (0.5 * (s[0] + 1.0)) * (thetamin - thetamax);
433 
434  // Point on outer wall
435  Vector<double> r_outer(2);
437 
438  // Point on inner wall
439  Vector<double> r_inner(2);
441 
442  // Get layer
443  for (unsigned i = 0; i < 2; i++)
444  {
445  f[i] =
446  r_inner[i] + (r_outer[i] - r_inner[i]) * (rmax - Rmin) / (Rmax - Rmin);
447  }
448  f[2] = zmin + (zmax - zmin) * (0.5 * (s[1] + 1.0));
449  }

References f(), i, Inner_boundary_cross_section_pt, Outer_boundary_cross_section_pt, oomph::GeomObject::position(), Global_Parameters::rmax, Rmax, Rmin, s, plotPSD::t, plotDoE::x, Global_Parameters::zmax, and Global_Parameters::zmin.

Referenced by macro_element_boundary().

Member Data Documentation

◆ Axial_spacing_fct_pt

AxialSpacingFctPt oomph::QuarterPipeDomain::Axial_spacing_fct_pt
private

Function pointer for function that implements axial spacing of macro elements

Referenced by axial_spacing_fct(), and axial_spacing_fct_pt().

◆ Inner_boundary_cross_section_pt

GeomObject* oomph::QuarterPipeDomain::Inner_boundary_cross_section_pt
private

Geom object representing the inner boundary of the cross section

Referenced by QuarterPipeDomain(), r_B(), r_D(), r_F(), r_L(), r_R(), r_U(), and ~QuarterPipeDomain().

◆ Length

double oomph::QuarterPipeDomain::Length
private

Length.

Referenced by macro_element_boundary().

◆ Nr

unsigned oomph::QuarterPipeDomain::Nr
private

Number of elements radial direction.

Referenced by macro_element_boundary().

◆ Ntheta

unsigned oomph::QuarterPipeDomain::Ntheta
private

Number of elements azimuthal direction.

Referenced by macro_element_boundary().

◆ Nz

unsigned oomph::QuarterPipeDomain::Nz
private

Number of elements axial direction.

Referenced by macro_element_boundary().

◆ Outer_boundary_cross_section_pt

GeomObject* oomph::QuarterPipeDomain::Outer_boundary_cross_section_pt
private

Geom object representing the outer boundary of the cross section

Referenced by QuarterPipeDomain(), r_B(), r_D(), r_F(), r_L(), r_R(), r_U(), and ~QuarterPipeDomain().

◆ Rmax

double oomph::QuarterPipeDomain::Rmax
private

Outer radius.

Referenced by macro_element_boundary(), r_B(), r_D(), r_F(), r_L(), r_R(), and r_U().

◆ Rmin

double oomph::QuarterPipeDomain::Rmin
private

Inner radius.

Referenced by macro_element_boundary(), r_B(), r_D(), r_F(), r_L(), r_R(), and r_U().


The documentation for this class was generated from the following file: