oomph::OneDimLagrange Namespace Reference

Functions

template<unsigned NNODE_1D>
void shape (const double &s, double *Psi)
 
template<unsigned NNODE_1D>
void dshape (const double &s, double *DPsi)
 
template<unsigned NNODE_1D>
void d2shape (const double &s, double *DPsi)
 
template<>
void shape< 2 > (const double &s, double *Psi)
 1D shape functions specialised to linear order (2 Nodes) More...
 
template<>
void dshape< 2 > (const double &s, double *DPsi)
 Derivatives of 1D shape functions specialised to linear order (2 Nodes) More...
 
template<>
void d2shape< 2 > (const double &s, double *DPsi)
 
template<>
void shape< 3 > (const double &s, double *Psi)
 1D shape functions specialised to quadratic order (3 Nodes) More...
 
template<>
void dshape< 3 > (const double &s, double *DPsi)
 
template<>
void d2shape< 3 > (const double &s, double *DPsi)
 
template<>
void shape< 4 > (const double &s, double *Psi)
 1D shape functions specialised to cubic order (4 Nodes) More...
 
template<>
void dshape< 4 > (const double &s, double *DPsi)
 Derivatives of 1D shape functions specialised to cubic order (4 Nodes) More...
 
template<>
void d2shape< 4 > (const double &s, double *DPsi)
 

Function Documentation

◆ d2shape()

template<unsigned NNODE_1D>
void oomph::OneDimLagrange::d2shape ( const double s,
double DPsi 
)

Definition for second derivatives of 1D Lagrange shape functions. The value of all the shape function derivatives at the local coordinate s are returned in the array DPsi.

594  {
595  std::ostringstream error_stream;
596  error_stream << "One dimensional Lagrange shape function "
597  << "second derivatives "
598  << "have not been defined "
599  << "for " << NNODE_1D << " nodes." << std::endl;
600  throw OomphLibError(
601  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
602  }
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ d2shape< 2 >()

template<>
void oomph::OneDimLagrange::d2shape< 2 > ( const double s,
double DPsi 
)
inline

Second Derivatives of 1D shape functions, specialised to linear order (2 Nodes)

626  {
627  DPsi[0] = 0.0;
628  DPsi[1] = 0.0;
629  }

◆ d2shape< 3 >()

template<>
void oomph::OneDimLagrange::d2shape< 3 > ( const double s,
double DPsi 
)
inline

Second Derivatives of 1D shape functions, specialised to quadratic order (3 Nodes)

657  {
658  DPsi[0] = 1.0;
659  DPsi[1] = -2.0;
660  DPsi[2] = 1.0;
661  }

◆ d2shape< 4 >()

template<>
void oomph::OneDimLagrange::d2shape< 4 > ( const double s,
double DPsi 
)
inline

Second Derivatives of 1D shape functions specialised to cubic order (4 Nodes)

701  {
702  // Output from Maple (modified by ALH, CHECK IT)
703  double t1 = 2.0 * s;
704  double t2 = 0.16875E1 * t1;
705  double t5 = 0.50625E1 * t1;
706  DPsi[0] = -t2 + 0.1125E1;
707  DPsi[1] = t5 - 0.1125E1;
708  DPsi[2] = -t5 - 0.1125E1;
709  DPsi[3] = t2 + 0.1125E1;
710  }
RealScalar s
Definition: level1_cplx_impl.h:130

References s.

◆ dshape()

template<unsigned NNODE_1D>
void oomph::OneDimLagrange::dshape ( const double s,
double DPsi 
)

Definition for derivatives of 1D Lagrange shape functions. The value of all the shape function derivatives at the local coordinate s are returned in the array DPsi.

579  {
580  std::ostringstream error_stream;
581  error_stream << "One dimensional Lagrange shape function derivatives "
582  << "have not been defined "
583  << "for " << NNODE_1D << " nodes." << std::endl;
584  throw OomphLibError(
585  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
586  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ dshape< 2 >()

◆ dshape< 3 >()

template<>
void oomph::OneDimLagrange::dshape< 3 > ( const double s,
double DPsi 
)
inline

◆ dshape< 4 >()

template<>
void oomph::OneDimLagrange::dshape< 4 > ( const double s,
double DPsi 
)
inline

Derivatives of 1D shape functions specialised to cubic order (4 Nodes)

685  {
686  // Output from Maple
687  double t1 = s * s;
688  double t2 = 0.16875E1 * t1;
689  double t3 = 0.1125E1 * s;
690  double t5 = 0.50625E1 * t1;
691  DPsi[0] = -t2 + t3 + 0.625E-1;
692  DPsi[1] = t5 - t3 - 0.16875E1;
693  DPsi[2] = -t5 - t3 + 0.16875E1;
694  DPsi[3] = t2 + t3 - 0.625E-1;
695  }

References s.

◆ shape()

template<unsigned NNODE_1D>
void oomph::OneDimLagrange::shape ( const double s,
double Psi 
)

Definition for 1D Lagrange shape functions. The value of all the shape functions at the local coordinate s are returned in the array Psi.

565  {
566  std::ostringstream error_stream;
567  error_stream << "One dimensional Lagrange shape functions "
568  << "have not been defined "
569  << "for " << NNODE_1D << " nodes." << std::endl;
570  throw OomphLibError(
571  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
572  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by TwoNodeGeometricElement::dshape_local(), oomph::TElementShape< 1, 2 >::dshape_local(), oomph::TElementShape< 1, 3 >::dshape_local(), oomph::TElementShape< 1, 4 >::dshape_local(), oomph::TElementShape< 2, 2 >::dshape_local(), oomph::TElementShape< 2, 3 >::dshape_local(), oomph::TElementShape< 2, 4 >::dshape_local(), oomph::TBubbleEnrichedElementShape< 2, 3 >::dshape_local(), oomph::TElementShape< 3, 2 >::dshape_local(), oomph::TElementShape< 3, 3 >::dshape_local(), oomph::TBubbleEnrichedElementShape< 3, 3 >::dshape_local(), oomph::HelmholtzPointSourceElement< ELEMENT >::fill_in_point_source_contribution_to_residuals(), oomph::PMLHelmholtzPointSourceElement< ELEMENT >::fill_in_point_source_contribution_to_residuals(), oomph::ProjectableLinearHeatAndElasticityElement< LINEAR_HEAT_AND_ELAST_ELEMENT >::get_field(), oomph::ProjectableTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::get_field(), oomph::ProjectableTaylorHoodSpaceTimeElement< TAYLOR_HOOD_ELEMENT >::get_field(), SpineGravityTractionElement< ELEMENT >::get_flux(), oomph::ElementWithMortaringStatusAtNodes< ELEMENT >::get_interpolated_mortared_status(), oomph::VorticitySmootherElement< ELEMENT >::get_interpolated_values(), RefineableModalPoissonEquations< DIM >::get_interpolated_values(), VolumeCoupling::getProjectionMatrix(), oomph::SpineAxisymmetricMarangoniSurfactantFluidInterfaceElement< ELEMENT >::interpolated_C(), oomph::SpineLineMarangoniSurfactantFluidInterfaceElement< ELEMENT >::interpolated_C(), ElasticAxisymmetricSolubleSurfactantTransportInterfaceElement< ELEMENT >::interpolated_C_bulk(), SSPorousChannelEquations::interpolated_f(), oomph::NavierStokesElementWithSingularity< BASIC_NAVIER_STOKES_ELEMENT >::interpolated_p_nst(), oomph::DummyBrickElement::interpolated_s_tet(), oomph::SpineLineMarangoniSurfactantFluidInterfaceElement< ELEMENT >::interpolated_T(), ComplexHarmonicEquations::interpolated_u(), HarmonicEquations::interpolated_u(), OrrSommerfeldEquations< DIM >::interpolated_u(), oomph::interpolated_u_cons_axisym_adv_diff(), oomph::FluxTransportEquations< DIM >::interpolated_u_flux_transport(), oomph::NavierStokesElementWithSingularity< BASIC_NAVIER_STOKES_ELEMENT >::interpolated_u_nst(), oomph::NavierStokesElementWithSingularity< BASIC_NAVIER_STOKES_ELEMENT >::interpolated_u_nst_fe_only(), PoissonElementWithSingularity< BASIC_POISSON_ELEMENT >::interpolated_u_poisson(), PoissonElementWithSingularity< BASIC_POISSON_ELEMENT >::interpolated_u_poisson_fe_only(), ComplexHarmonicEquations::interpolated_w(), oomph::RefineableQElement< 2 >::interpolated_zeta_on_edge(), oomph::RefineableQElement< 3 >::interpolated_zeta_on_face(), DistanceWrapper< ELEMENT >::output(), oomph::QSpectralElement< 1, NNODE_1D >::output(), oomph::VolumeCoupledElement< ELEMENT >::output(), oomph::VorticitySmootherElement< ELEMENT >::output(), oomph::PoroelasticityEquations< DIM >::output(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::output_veloc(), oomph::NavierStokesEquations< DIM >::output_veloc(), oomph::HeatedPenetratorFluxElement< ELEMENT >::penetrator_temperature(), oomph::TElement< 1, NNODE_1D >::shape(), oomph::TElement< 2, NNODE_1D >::shape(), oomph::TElement< 3, NNODE_1D >::shape(), oomph::TBubbleEnrichedElement< DIM, 3 >::shape(), oomph::ClampedHermiteShellBoundaryConditionElement::shape(), oomph::HeatedPenetratorFluxElement< ELEMENT >::shape_and_test(), oomph::SurfaceContactElementBase< ELEMENT >::shape_p(), oomph::SpineLineMarangoniFluidInterfaceElement< ELEMENT >::sigma(), oomph::FixedVolumeSpineLineMarangoniFluidInterfaceElement< ELEMENT >::sigma(), oomph::SpineAxisymmetricMarangoniSurfactantFluidInterfaceElement< ELEMENT >::sigma(), oomph::SpineLineMarangoniSurfactantFluidInterfaceElement< ELEMENT >::sigma(), and oomph::VorticitySmootherElement< ELEMENT >::vorticity_and_its_derivs().

◆ shape< 2 >()

template<>
void oomph::OneDimLagrange::shape< 2 > ( const double s,
double Psi 
)
inline

1D shape functions specialised to linear order (2 Nodes)

609  {
610  Psi[0] = 0.5 * (1.0 - s);
611  Psi[1] = 0.5 * (1.0 + s);
612  }

References s.

Referenced by oomph::GeneralisedNewtonianQTaylorHoodElement< DIM >::dpshape_and_dptest_eulerian_nst(), oomph::QTaylorHoodElement< DIM >::dpshape_and_dptest_eulerian_nst(), oomph::QTaylorHoodSpaceTimeElement< DIM >::dpshape_and_dptest_eulerian_nst(), oomph::QTaylorHoodSpaceTimeElement< DIM >::dpshape_eulerian(), oomph::QTaylorHoodMixedOrderSpaceTimeElement< DIM >::dpshape_eulerian(), oomph::QTaylorHoodSpaceTimeElement< DIM >::dptest_eulerian(), oomph::QTaylorHoodMixedOrderSpaceTimeElement< DIM >::dptest_eulerian(), oomph::AxisymmetricQTaylorHoodElement::pshape_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricQTaylorHoodElement::pshape_axi_nst(), oomph::LinearisedAxisymmetricQTaylorHoodElement::pshape_lin_axi_nst(), oomph::LinearisedAxisymmetricQTaylorHoodElement::pshape_linearised_axi_nst(), oomph::LinearisedQTaylorHoodElement::pshape_linearised_nst(), oomph::GeneralisedNewtonianQTaylorHoodElement< DIM >::pshape_nst(), oomph::QTaylorHoodElement< DIM >::pshape_nst(), oomph::QTaylorHoodSpaceTimeElement< DIM >::pshape_nst(), oomph::QTaylorHoodMixedOrderSpaceTimeElement< DIM >::pshape_nst(), oomph::PolarTaylorHoodElement::pshape_pnst(), oomph::QSphericalTaylorHoodElement::pshape_spherical_nst(), oomph::QTaylorHoodSpaceTimeElement< DIM >::ptest_nst(), oomph::QTaylorHoodMixedOrderSpaceTimeElement< DIM >::ptest_nst(), and oomph::QPVDElementWithContinuousPressure< DIM >::solid_pshape().

◆ shape< 3 >()

◆ shape< 4 >()

template<>
void oomph::OneDimLagrange::shape< 4 > ( const double s,
double Psi 
)
inline

1D shape functions specialised to cubic order (4 Nodes)

666  {
667  // Output from Maple
668  double t1 = s * s;
669  double t2 = t1 * s;
670  double t3 = 0.5625 * t2;
671  double t4 = 0.5625 * t1;
672  double t5 = 0.625E-1 * s;
673  double t7 = 0.16875E1 * t2;
674  double t8 = 0.16875E1 * s;
675  Psi[0] = -t3 + t4 + t5 - 0.625E-1;
676  Psi[1] = t7 - t4 - t8 + 0.5625;
677  Psi[2] = -t7 - t4 + t8 + 0.5625;
678  Psi[3] = t3 + t4 - t5 - 0.625E-1;
679  }

References s.