QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D > Class Template Reference

#include <gen_axisym_advection_diffusion_elements.h>

+ Inheritance diagram for QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >:

Public Member Functions

 QGeneralisedAxisymAdvectionDiffusionElement ()
 
 QGeneralisedAxisymAdvectionDiffusionElement (const QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D > &dummy)=delete
 Broken copy constructor. More...
 
unsigned required_nvalue (const unsigned &n) const
 Broken assignment operator. More...
 
void output (std::ostream &outfile)
 
void output (std::ostream &outfile, const unsigned &n_plot)
 
void output (FILE *file_pt)
 
void output (FILE *file_pt, const unsigned &n_plot)
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 

Protected Member Functions

double dshape_and_dtest_eulerian_cons_axisym_adv_diff (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
 
double dshape_and_dtest_eulerian_at_knot_cons_axisym_adv_diff (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
 

Static Private Attributes

static const unsigned Initial_Nvalue = 1
 

Detailed Description

template<unsigned NNODE_1D>
class QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >

//////////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////////// QGeneralisedAxisymAdvectionDiffusionElement elements are linear/quadrilateral/brick-shaped Advection Diffusion elements with isoparametric interpolation for the function.

Constructor & Destructor Documentation

◆ QGeneralisedAxisymAdvectionDiffusionElement() [1/2]

Constructor: Call constructors for QElement and Advection Diffusion equations

◆ QGeneralisedAxisymAdvectionDiffusionElement() [2/2]

Broken copy constructor.

Member Function Documentation

◆ dshape_and_dtest_eulerian_at_knot_cons_axisym_adv_diff()

template<unsigned NNODE_1D>
double QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >::dshape_and_dtest_eulerian_at_knot_cons_axisym_adv_diff ( const unsigned ipt,
Shape &  psi,
DShape &  dpsidx,
Shape &  test,
DShape &  dtestdx 
) const
inlineprotected

Shape, test functions & derivs. w.r.t. to global coords. at integration point ipt. Return Jacobian.

Define the shape functions and test functions and derivatives w.r.t. global coordinates and return Jacobian of mapping.

Galerkin: Test functions = shape functions

813  {
814  // Call the geometrical shape functions and derivatives
815  double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
816 
817  // Set the test functions equal to the shape functions (pointer copy)
818  test = psi;
819  dtestdx = dpsidx;
820 
821  // Return the jacobian
822  return J;
823  }
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
Definition: indexed_view.cpp:20

References J.

◆ dshape_and_dtest_eulerian_cons_axisym_adv_diff()

template<unsigned NNODE_1D>
double QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >::dshape_and_dtest_eulerian_cons_axisym_adv_diff ( const Vector< double > &  s,
Shape &  psi,
DShape &  dpsidx,
Shape &  test,
DShape &  dtestdx 
) const
inlineprotected

Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.

Define the shape functions and test functions and derivatives w.r.t. global coordinates and return Jacobian of mapping.

Galerkin: Test functions = shape functions

779  {
780  // Call the geometrical shape functions and derivatives
781  double J = this->dshape_eulerian(s, psi, dpsidx);
782 
783  // Loop over the test functions and derivatives and set them equal to the
784  // shape functions
785  for (unsigned i = 0; i < NNODE_1D; i++)
786  {
787  test[i] = psi[i];
788  for (unsigned j = 0; j < 2; j++)
789  {
790  dtestdx(i, j) = dpsidx(i, j);
791  }
792  }
793 
794  // Return the jacobian
795  return J;
796  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
RealScalar s
Definition: level1_cplx_impl.h:130
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References i, J, j, and s.

◆ output() [1/4]

template<unsigned NNODE_1D>
void QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >::output ( FILE *  file_pt)
inline

C-style output function: r,z,u

708  {
710  }
void output(FILE *file_pt, const unsigned &n_plot)

References oomph::output().

◆ output() [2/4]

template<unsigned NNODE_1D>
void QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >::output ( FILE *  file_pt,
const unsigned n_plot 
)
inline

C-style output function: r,z,u at n_plot^2 plot points

715  {
717  }

References oomph::output().

◆ output() [3/4]

template<unsigned NNODE_1D>
void QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >::output ( std::ostream &  outfile)
inline

Output function: r,z,u

References oomph::output().

◆ output() [4/4]

template<unsigned NNODE_1D>
void QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >::output ( std::ostream &  outfile,
const unsigned n_plot 
)
inline

Output function: r,z,u at n_plot^2 plot points

700  {
702  }

References oomph::output().

◆ output_fct() [1/2]

template<unsigned NNODE_1D>
void QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >::output_fct ( std::ostream &  outfile,
const unsigned n_plot,
const double time,
FiniteElement::UnsteadyExactSolutionFctPt  exact_soln_pt 
)
inline

Output function for a time-dependent exact solution. r,z,u_exact at n_plot^2 plot points (Calls the steady version)

737  {
739  outfile, n_plot, time, exact_soln_pt);
740  }
virtual void output_fct(std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: gen_axisym_advection_diffusion_elements.h:191

References oomph::output_fct().

◆ output_fct() [2/2]

template<unsigned NNODE_1D>
void QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >::output_fct ( std::ostream &  outfile,
const unsigned n_plot,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt 
)
inline

Output function for an exact solution: r,z,u_exact at n_plot^2 plot points

724  {
726  outfile, n_plot, exact_soln_pt);
727  }

References oomph::output_fct().

◆ required_nvalue()

template<unsigned NNODE_1D>
unsigned QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >::required_nvalue ( const unsigned n) const
inline

Broken assignment operator.

Required # of ‘values’ (pinned or dofs) at node n

686  {
687  return Initial_Nvalue;
688  }
static const unsigned Initial_Nvalue
Definition: gen_axisym_advection_diffusion_elements.h:664

Member Data Documentation

◆ Initial_Nvalue

template<unsigned NNODE_1D>
const unsigned QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D >::Initial_Nvalue = 1
staticprivate

Static array of ints to hold number of variables at nodes: Initial_Nvalue[n]


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