$extrastylesheet
#include <h1_fe_transformation.h>

Public Member Functions | |
| H1FETransformation () | |
| virtual | ~H1FETransformation () |
| virtual void | map_phi (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< OutputShape > &, std::vector< std::vector< OutputShape > > &) const |
| virtual void | map_dphi (const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputGradient > > &dphi, std::vector< std::vector< OutputShape > > &dphidx, std::vector< std::vector< OutputShape > > &dphidy, std::vector< std::vector< OutputShape > > &dphidz) const |
| virtual void | map_d2phi (const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputTensor > > &d2phi, std::vector< std::vector< OutputShape > > &d2phidx2, std::vector< std::vector< OutputShape > > &d2phidxdy, std::vector< std::vector< OutputShape > > &d2phidxdz, std::vector< std::vector< OutputShape > > &d2phidy2, std::vector< std::vector< OutputShape > > &d2phidydz, std::vector< std::vector< OutputShape > > &d2phidz2) const |
| virtual void | map_curl (const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< OutputShape > > &curl_phi) const |
| virtual void | map_div (const unsigned int dim, const Elem *const elem, const std::vector< Point > &qp, const FEGenericBase< OutputShape > &fe, std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputDivergence > > &div_phi) const |
| template<> | |
| void | map_curl (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< Real > &, std::vector< std::vector< Real > > &) const |
| template<> | |
| void | map_curl (const unsigned int dim, const Elem *const, const std::vector< Point > &, const FEGenericBase< RealGradient > &fe, std::vector< std::vector< RealGradient > > &curl_phi) const |
| template<> | |
| void | map_div (const unsigned int, const Elem *const, const std::vector< Point > &, const FEGenericBase< Real > &, std::vector< std::vector< FEGenericBase< Real >::OutputDivergence > > &) const |
| template<> | |
| void | map_div (const unsigned int dim, const Elem *const, const std::vector< Point > &, const FEGenericBase< RealGradient > &fe, std::vector< std::vector< FEGenericBase< RealGradient >::OutputDivergence > > &div_phi) const |
Static Public Member Functions | |
| static UniquePtr < FETransformationBase < OutputShape > > | build (const FEType &type) |
This class handles the computation of the shape functions in the physical domain for H1 conforming elements. This class assumes the FEGenericBase object has been initialized in the reference domain (i.e. init_shape_functions has been called).
Definition at line 35 of file h1_fe_transformation.h.
| libMesh::H1FETransformation< OutputShape >::H1FETransformation | ( | ) | [inline] |
Definition at line 39 of file h1_fe_transformation.h.
: FETransformationBase<OutputShape>(){}
| virtual libMesh::H1FETransformation< OutputShape >::~H1FETransformation | ( | ) | [inline, virtual] |
Definition at line 42 of file h1_fe_transformation.h.
{}
| UniquePtr< FETransformationBase< OutputShape > > libMesh::FETransformationBase< OutputShape >::build | ( | const FEType & | type | ) | [static, inherited] |
Builds an FETransformation object based on the finite element type
Definition at line 26 of file fe_transformation_base.C.
References libMesh::BERNSTEIN, libMesh::CLOUGH, libMesh::FEType::family, libMesh::HERMITE, libMesh::HIERARCHIC, libMesh::JACOBI_20_00, libMesh::JACOBI_30_00, libMesh::L2_HIERARCHIC, libMesh::L2_LAGRANGE, libMesh::LAGRANGE, libMesh::LAGRANGE_VEC, libMesh::MONOMIAL, libMesh::NEDELEC_ONE, libMesh::SCALAR, libMesh::SUBDIVISION, libMesh::SZABAB, and libMesh::XYZ.
{
switch (fe_type.family)
{
// H1 Conforming Elements
case LAGRANGE:
case HIERARCHIC:
case BERNSTEIN:
case SZABAB:
case CLOUGH: // PB: Really H2
case HERMITE: // PB: Really H2
case SUBDIVISION:
case LAGRANGE_VEC:
case MONOMIAL: // PB: Shouldn't this be L2 conforming?
case XYZ: // PB: Shouldn't this be L2 conforming?
case L2_HIERARCHIC: // PB: Shouldn't this be L2 conforming?
case L2_LAGRANGE: // PB: Shouldn't this be L2 conforming?
case JACOBI_20_00: // PB: For infinite elements...
case JACOBI_30_00: // PB: For infinite elements...
return UniquePtr<FETransformationBase<OutputShape> >(new H1FETransformation<OutputShape>);
// HCurl Conforming Elements
case NEDELEC_ONE:
return UniquePtr<FETransformationBase<OutputShape> >(new HCurlFETransformation<OutputShape>);
// HDiv Conforming Elements
// L2 Conforming Elements
// Other...
case SCALAR:
// Should never need this for SCALARs
return UniquePtr<FETransformationBase<OutputShape> >(new H1FETransformation<OutputShape>);
default:
libmesh_error_msg("Unknown family = " << fe_type.family);
}
libmesh_error_msg("We'll never get here!");
return UniquePtr<FETransformationBase<OutputShape> >();
}
| virtual void libMesh::H1FETransformation< OutputShape >::map_curl | ( | const unsigned int | dim, |
| const Elem *const | elem, | ||
| const std::vector< Point > & | qp, | ||
| const FEGenericBase< OutputShape > & | fe, | ||
| std::vector< std::vector< OutputShape > > & | curl_phi | ||
| ) | const [virtual] |
Evaluates the shape function curl in physical coordinates based on H1 conforming finite element transformation.
Implements libMesh::FETransformationBase< OutputShape >.
| void libMesh::H1FETransformation< Real >::map_curl | ( | const unsigned int | , |
| const Elem * | const, | ||
| const std::vector< Point > & | , | ||
| const FEGenericBase< Real > & | , | ||
| std::vector< std::vector< Real > > & | |||
| ) | const |
Definition at line 449 of file h1_fe_transformation.C.
{
libmesh_error_msg("Computing the curl of a shape function only \nmakes sense for vector-valued elements.");
}
| void libMesh::H1FETransformation< RealGradient >::map_curl | ( | const unsigned int | dim, |
| const Elem * | const, | ||
| const std::vector< Point > & | , | ||
| const FEGenericBase< RealGradient > & | fe, | ||
| std::vector< std::vector< RealGradient > > & | curl_phi | ||
| ) | const |
Definition at line 459 of file h1_fe_transformation.C.
References libMesh::FEMap::get_detadx(), libMesh::FEMap::get_detady(), libMesh::FEMap::get_detadz(), libMesh::FEGenericBase< OutputType >::get_dphideta(), libMesh::FEGenericBase< OutputType >::get_dphidxi(), libMesh::FEGenericBase< OutputType >::get_dphidzeta(), libMesh::FEMap::get_dxidx(), libMesh::FEMap::get_dxidy(), libMesh::FEMap::get_dxidz(), libMesh::FEMap::get_dzetadx(), libMesh::FEMap::get_dzetady(), libMesh::FEMap::get_dzetadz(), libMesh::FEAbstract::get_fe_map(), and libMesh::Real.
{
switch (dim)
{
case 0:
case 1:
libmesh_error_msg("The curl only make sense in 2D and 3D");
case 2:
{
const std::vector<std::vector<RealGradient> >& dphidxi = fe.get_dphidxi();
const std::vector<std::vector<RealGradient> >& dphideta = fe.get_dphideta();
const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
#if LIBMESH_DIM > 2
const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();
#endif
const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
#if LIBMESH_DIM > 2
const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz();
#endif
/*
For 2D elements in 3D space:
curl = ( -dphi_y/dz, dphi_x/dz, dphi_y/dx - dphi_x/dy )
*/
for (unsigned int i=0; i<curl_phi.size(); i++)
for (unsigned int p=0; p<curl_phi[i].size(); p++)
{
Real dphiy_dx = (dphidxi[i][p].slice(1))*dxidx_map[p]
+ (dphideta[i][p].slice(1))*detadx_map[p];
Real dphix_dy = (dphidxi[i][p].slice(0))*dxidy_map[p]
+ (dphideta[i][p].slice(0))*detady_map[p];
curl_phi[i][p].slice(2) = dphiy_dx - dphix_dy;
#if LIBMESH_DIM > 2
Real dphiy_dz = (dphidxi[i][p].slice(1))*dxidz_map[p]
+ (dphideta[i][p].slice(1))*detadz_map[p];
Real dphix_dz = (dphidxi[i][p].slice(0))*dxidz_map[p]
+ (dphideta[i][p].slice(0))*detadz_map[p];
curl_phi[i][p].slice(0) = -dphiy_dz;
curl_phi[i][p].slice(1) = dphix_dz;
#endif
}
break;
}
case 3:
{
const std::vector<std::vector<RealGradient> >& dphidxi = fe.get_dphidxi();
const std::vector<std::vector<RealGradient> >& dphideta = fe.get_dphideta();
const std::vector<std::vector<RealGradient> >& dphidzeta = fe.get_dphidzeta();
const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();
const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz();
const std::vector<Real>& dzetadx_map = fe.get_fe_map().get_dzetadx();
const std::vector<Real>& dzetady_map = fe.get_fe_map().get_dzetady();
const std::vector<Real>& dzetadz_map = fe.get_fe_map().get_dzetadz();
/*
In 3D: curl = ( dphi_z/dy - dphi_y/dz, dphi_x/dz - dphi_z/dx, dphi_y/dx - dphi_x/dy )
*/
for (unsigned int i=0; i<curl_phi.size(); i++)
for (unsigned int p=0; p<curl_phi[i].size(); p++)
{
Real dphiz_dy = (dphidxi[i][p].slice(2))*dxidy_map[p]
+ (dphideta[i][p].slice(2))*detady_map[p]
+ (dphidzeta[i][p].slice(2))*dzetady_map[p];
Real dphiy_dz = (dphidxi[i][p].slice(1))*dxidz_map[p]
+ (dphideta[i][p].slice(1))*detadz_map[p]
+ (dphidzeta[i][p].slice(1))*dzetadz_map[p];
Real dphix_dz = (dphidxi[i][p].slice(0))*dxidz_map[p]
+ (dphideta[i][p].slice(0))*detadz_map[p]
+ (dphidzeta[i][p].slice(0))*dzetadz_map[p];
Real dphiz_dx = (dphidxi[i][p].slice(2))*dxidx_map[p]
+ (dphideta[i][p].slice(2))*detadx_map[p]
+ (dphidzeta[i][p].slice(2))*dzetadx_map[p];
Real dphiy_dx = (dphidxi[i][p].slice(1))*dxidx_map[p]
+ (dphideta[i][p].slice(1))*detadx_map[p]
+ (dphidzeta[i][p].slice(1))*dzetadx_map[p];
Real dphix_dy = (dphidxi[i][p].slice(0))*dxidy_map[p]
+ (dphideta[i][p].slice(0))*detady_map[p]
+ (dphidzeta[i][p].slice(0))*dzetady_map[p];
curl_phi[i][p].slice(0) = dphiz_dy - dphiy_dz;
curl_phi[i][p].slice(1) = dphix_dz - dphiz_dx;
curl_phi[i][p].slice(2) = dphiy_dx - dphix_dy;
}
break;
}
default:
libmesh_error_msg("Invalid dim = " << dim);
}
}
| void libMesh::H1FETransformation< OutputShape >::map_d2phi | ( | const unsigned int | dim, |
| const Elem *const | elem, | ||
| const std::vector< Point > & | qp, | ||
| const FEGenericBase< OutputShape > & | fe, | ||
| std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputTensor > > & | d2phi, | ||
| std::vector< std::vector< OutputShape > > & | d2phidx2, | ||
| std::vector< std::vector< OutputShape > > & | d2phidxdy, | ||
| std::vector< std::vector< OutputShape > > & | d2phidxdz, | ||
| std::vector< std::vector< OutputShape > > & | d2phidy2, | ||
| std::vector< std::vector< OutputShape > > & | d2phidydz, | ||
| std::vector< std::vector< OutputShape > > & | d2phidz2 | ||
| ) | const [virtual] |
Evaluates shape function Hessians in physical coordinates based on H1 conforming finite element transformation. FIXME: These are currently not calculated correctly for non-affine elements. The second derivative terms of the FEMap are not implemented.
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 215 of file h1_fe_transformation.C.
References libMesh::err, libMesh::FEGenericBase< OutputType >::get_d2phideta2(), libMesh::FEGenericBase< OutputType >::get_d2phidetadzeta(), libMesh::FEGenericBase< OutputType >::get_d2phidxi2(), libMesh::FEGenericBase< OutputType >::get_d2phidxideta(), libMesh::FEGenericBase< OutputType >::get_d2phidxidzeta(), libMesh::FEGenericBase< OutputType >::get_d2phidzeta2(), libMesh::FEMap::get_detadx(), libMesh::FEMap::get_detady(), libMesh::FEMap::get_detadz(), libMesh::FEMap::get_dxidx(), libMesh::FEMap::get_dxidy(), libMesh::FEMap::get_dxidz(), libMesh::FEMap::get_dzetadx(), libMesh::FEMap::get_dzetady(), libMesh::FEMap::get_dzetadz(), libMesh::FEAbstract::get_fe_map(), and libMesh::Elem::has_affine_map().
{
libmesh_do_once(
if (!elem->has_affine_map())
{
libMesh::err << "WARNING: Second derivatives are not currently "
<< "correctly calculated on non-affine elements!"
<< std::endl;
}
);
switch (dim)
{
case 0: // No derivatives in 0D
{
for (unsigned int i=0; i<d2phi.size(); i++)
for (unsigned int p=0; p<d2phi[i].size(); p++)
{
d2phi[i][p] = 0.;
}
break;
}
case 1:
{
const std::vector<std::vector<OutputShape> >& d2phidxi2 = fe.get_d2phidxi2();
const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
#if LIBMESH_DIM>1
const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
#endif
#if LIBMESH_DIM>2
const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();
#endif
for (unsigned int i=0; i<d2phi.size(); i++)
for (unsigned int p=0; p<d2phi[i].size(); p++)
{
d2phi[i][p].slice(0).slice(0) = d2phidx2[i][p] =
d2phidxi2[i][p]*dxidx_map[p]*dxidx_map[p];
#if LIBMESH_DIM>1
d2phi[i][p].slice(0).slice(1) =
d2phi[i][p].slice(1).slice(0) = d2phidxdy[i][p] =
d2phidxi2[i][p]*dxidx_map[p]*dxidy_map[p];
d2phi[i][p].slice(1).slice(1) = d2phidy2[i][p] =
d2phidxi2[i][p]*dxidy_map[p]*dxidy_map[p];
#endif
#if LIBMESH_DIM>2
d2phi[i][p].slice(0).slice(2) =
d2phi[i][p].slice(2).slice(0) = d2phidxdz[i][p] =
d2phidxi2[i][p]*dxidx_map[p]*dxidz_map[p];
d2phi[i][p].slice(1).slice(2) =
d2phi[i][p].slice(2).slice(1) = d2phidydz[i][p] =
d2phidxi2[i][p]*dxidy_map[p]*dxidz_map[p];
d2phi[i][p].slice(2).slice(2) = d2phidz2[i][p] =
d2phidxi2[i][p]*dxidz_map[p]*dxidz_map[p];
#endif
}
break;
}
case 2:
{
const std::vector<std::vector<OutputShape> >& d2phidxi2 = fe.get_d2phidxi2();
const std::vector<std::vector<OutputShape> >& d2phidxideta = fe.get_d2phidxideta();
const std::vector<std::vector<OutputShape> >& d2phideta2 = fe.get_d2phideta2();
const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
#if LIBMESH_DIM > 2
const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();
#endif
const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
#if LIBMESH_DIM > 2
const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz();
#endif
for (unsigned int i=0; i<d2phi.size(); i++)
for (unsigned int p=0; p<d2phi[i].size(); p++)
{
d2phi[i][p].slice(0).slice(0) = d2phidx2[i][p] =
d2phidxi2[i][p]*dxidx_map[p]*dxidx_map[p] +
2*d2phidxideta[i][p]*dxidx_map[p]*detadx_map[p] +
d2phideta2[i][p]*detadx_map[p]*detadx_map[p];
d2phi[i][p].slice(0).slice(1) =
d2phi[i][p].slice(1).slice(0) = d2phidxdy[i][p] =
d2phidxi2[i][p]*dxidx_map[p]*dxidy_map[p] +
d2phidxideta[i][p]*dxidx_map[p]*detady_map[p] +
d2phideta2[i][p]*detadx_map[p]*detady_map[p] +
d2phidxideta[i][p]*detadx_map[p]*dxidy_map[p];
d2phi[i][p].slice(1).slice(1) = d2phidy2[i][p] =
d2phidxi2[i][p]*dxidy_map[p]*dxidy_map[p] +
2*d2phidxideta[i][p]*dxidy_map[p]*detady_map[p] +
d2phideta2[i][p]*detady_map[p]*detady_map[p];
#if LIBMESH_DIM > 2
d2phi[i][p].slice(0).slice(2) =
d2phi[i][p].slice(2).slice(0) = d2phidxdz[i][p] =
d2phidxi2[i][p]*dxidx_map[p]*dxidz_map[p] +
d2phidxideta[i][p]*dxidx_map[p]*detadz_map[p] +
d2phideta2[i][p]*detadx_map[p]*detadz_map[p] +
d2phidxideta[i][p]*detadx_map[p]*dxidz_map[p];
d2phi[i][p].slice(1).slice(2) =
d2phi[i][p].slice(2).slice(1) = d2phidydz[i][p] =
d2phidxi2[i][p]*dxidy_map[p]*dxidz_map[p] +
d2phidxideta[i][p]*dxidy_map[p]*detadz_map[p] +
d2phideta2[i][p]*detady_map[p]*detadz_map[p] +
d2phidxideta[i][p]*detady_map[p]*dxidz_map[p];
d2phi[i][p].slice(2).slice(2) = d2phidz2[i][p] =
d2phidxi2[i][p]*dxidz_map[p]*dxidz_map[p] +
2*d2phidxideta[i][p]*dxidz_map[p]*detadz_map[p] +
d2phideta2[i][p]*detadz_map[p]*detadz_map[p];
#endif
}
break;
}
case 3:
{
const std::vector<std::vector<OutputShape> >& d2phidxi2 = fe.get_d2phidxi2();
const std::vector<std::vector<OutputShape> >& d2phidxideta = fe.get_d2phidxideta();
const std::vector<std::vector<OutputShape> >& d2phideta2 = fe.get_d2phideta2();
const std::vector<std::vector<OutputShape> >& d2phidxidzeta = fe.get_d2phidxidzeta();
const std::vector<std::vector<OutputShape> >& d2phidetadzeta = fe.get_d2phidetadzeta();
const std::vector<std::vector<OutputShape> >& d2phidzeta2 = fe.get_d2phidzeta2();
const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();
const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz();
const std::vector<Real>& dzetadx_map = fe.get_fe_map().get_dzetadx();
const std::vector<Real>& dzetady_map = fe.get_fe_map().get_dzetady();
const std::vector<Real>& dzetadz_map = fe.get_fe_map().get_dzetadz();
for (unsigned int i=0; i<d2phi.size(); i++)
for (unsigned int p=0; p<d2phi[i].size(); p++)
{
d2phi[i][p].slice(0).slice(0) = d2phidx2[i][p] =
d2phidxi2[i][p]*dxidx_map[p]*dxidx_map[p] +
2*d2phidxideta[i][p]*dxidx_map[p]*detadx_map[p] +
2*d2phidxidzeta[i][p]*dxidx_map[p]*dzetadx_map[p] +
2*d2phidetadzeta[i][p]*detadx_map[p]*dzetadx_map[p] +
d2phideta2[i][p]*detadx_map[p]*detadx_map[p] +
d2phidzeta2[i][p]*dzetadx_map[p]*dzetadx_map[p];
d2phi[i][p].slice(0).slice(1) =
d2phi[i][p].slice(1).slice(0) = d2phidxdy[i][p] =
d2phidxi2[i][p]*dxidx_map[p]*dxidy_map[p] +
d2phidxideta[i][p]*dxidx_map[p]*detady_map[p] +
d2phidxidzeta[i][p]*dxidx_map[p]*dzetady_map[p] +
d2phideta2[i][p]*detadx_map[p]*detady_map[p] +
d2phidxideta[i][p]*detadx_map[p]*dxidy_map[p] +
d2phidetadzeta[i][p]*detadx_map[p]*dzetady_map[p] +
d2phidzeta2[i][p]*dzetadx_map[p]*dzetady_map[p] +
d2phidxidzeta[i][p]*dzetadx_map[p]*dxidy_map[p] +
d2phidetadzeta[i][p]*dzetadx_map[p]*detady_map[p];
d2phi[i][p].slice(0).slice(2) =
d2phi[i][p].slice(2).slice(0) = d2phidy2[i][p] =
d2phidxi2[i][p]*dxidx_map[p]*dxidz_map[p] +
d2phidxideta[i][p]*dxidx_map[p]*detadz_map[p] +
d2phidxidzeta[i][p]*dxidx_map[p]*dzetadz_map[p] +
d2phideta2[i][p]*detadx_map[p]*detadz_map[p] +
d2phidxideta[i][p]*detadx_map[p]*dxidz_map[p] +
d2phidetadzeta[i][p]*detadx_map[p]*dzetadz_map[p] +
d2phidzeta2[i][p]*dzetadx_map[p]*dzetadz_map[p] +
d2phidxidzeta[i][p]*dzetadx_map[p]*dxidz_map[p] +
d2phidetadzeta[i][p]*dzetadx_map[p]*detadz_map[p];
d2phi[i][p].slice(1).slice(1) = d2phidxdz[i][p] =
d2phidxi2[i][p]*dxidy_map[p]*dxidy_map[p] +
2*d2phidxideta[i][p]*dxidy_map[p]*detady_map[p] +
2*d2phidxidzeta[i][p]*dxidy_map[p]*dzetady_map[p] +
2*d2phidetadzeta[i][p]*detady_map[p]*dzetady_map[p] +
d2phideta2[i][p]*detady_map[p]*detady_map[p] +
d2phidzeta2[i][p]*dzetady_map[p]*dzetady_map[p];
d2phi[i][p].slice(1).slice(2) =
d2phi[i][p].slice(2).slice(1) = d2phidydz[i][p] =
d2phidxi2[i][p]*dxidy_map[p]*dxidz_map[p] +
d2phidxideta[i][p]*dxidy_map[p]*detadz_map[p] +
d2phidxidzeta[i][p]*dxidy_map[p]*dzetadz_map[p] +
d2phideta2[i][p]*detady_map[p]*detadz_map[p] +
d2phidxideta[i][p]*detady_map[p]*dxidz_map[p] +
d2phidetadzeta[i][p]*detady_map[p]*dzetadz_map[p] +
d2phidzeta2[i][p]*dzetady_map[p]*dzetadz_map[p] +
d2phidxidzeta[i][p]*dzetady_map[p]*dxidz_map[p] +
d2phidetadzeta[i][p]*dzetady_map[p]*detadz_map[p];
d2phi[i][p].slice(2).slice(2) = d2phidz2[i][p] =
d2phidxi2[i][p]*dxidz_map[p]*dxidz_map[p] +
2*d2phidxideta[i][p]*dxidz_map[p]*detadz_map[p] +
2*d2phidxidzeta[i][p]*dxidz_map[p]*dzetadz_map[p] +
2*d2phidetadzeta[i][p]*detadz_map[p]*dzetadz_map[p] +
d2phideta2[i][p]*detadz_map[p]*detadz_map[p] +
d2phidzeta2[i][p]*dzetadz_map[p]*dzetadz_map[p];
}
break;
}
default:
libmesh_error_msg("Invalid dim = " << dim);
} // switch(dim)
}
| virtual void libMesh::H1FETransformation< OutputShape >::map_div | ( | const unsigned int | dim, |
| const Elem *const | elem, | ||
| const std::vector< Point > & | qp, | ||
| const FEGenericBase< OutputShape > & | fe, | ||
| std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputDivergence > > & | div_phi | ||
| ) | const [virtual] |
Evaluates the shape function divergence in physical coordinates based on H1 conforming finite element transformation.
Implements libMesh::FETransformationBase< OutputShape >.
| void libMesh::H1FETransformation< Real >::map_div | ( | const unsigned int | , |
| const Elem * | const, | ||
| const std::vector< Point > & | , | ||
| const FEGenericBase< Real > & | , | ||
| std::vector< std::vector< FEGenericBase< Real >::OutputDivergence > > & | |||
| ) | const |
Definition at line 584 of file h1_fe_transformation.C.
{
libmesh_error_msg("Computing the divergence of a shape function only \nmakes sense for vector-valued elements.");
}
| void libMesh::H1FETransformation< RealGradient >::map_div | ( | const unsigned int | dim, |
| const Elem * | const, | ||
| const std::vector< Point > & | , | ||
| const FEGenericBase< RealGradient > & | fe, | ||
| std::vector< std::vector< FEGenericBase< RealGradient >::OutputDivergence > > & | div_phi | ||
| ) | const |
Definition at line 596 of file h1_fe_transformation.C.
References libMesh::FEMap::get_detadx(), libMesh::FEMap::get_detady(), libMesh::FEMap::get_detadz(), libMesh::FEGenericBase< OutputType >::get_dphideta(), libMesh::FEGenericBase< OutputType >::get_dphidxi(), libMesh::FEGenericBase< OutputType >::get_dphidzeta(), libMesh::FEMap::get_dxidx(), libMesh::FEMap::get_dxidy(), libMesh::FEMap::get_dxidz(), libMesh::FEMap::get_dzetadx(), libMesh::FEMap::get_dzetady(), libMesh::FEMap::get_dzetadz(), libMesh::FEAbstract::get_fe_map(), and libMesh::Real.
{
switch (dim)
{
case 0:
case 1:
libmesh_error_msg("The divergence only make sense in 2D and 3D");
case 2:
{
const std::vector<std::vector<RealGradient> >& dphidxi = fe.get_dphidxi();
const std::vector<std::vector<RealGradient> >& dphideta = fe.get_dphideta();
const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
/*
In 2D: div = dphi_x/dx + dphi_y/dy
*/
for (unsigned int i=0; i<div_phi.size(); i++)
for (unsigned int p=0; p<div_phi[i].size(); p++)
{
Real dphix_dx = (dphidxi[i][p].slice(0))*dxidx_map[p]
+ (dphideta[i][p].slice(0))*detadx_map[p];
Real dphiy_dy = (dphidxi[i][p].slice(1))*dxidy_map[p]
+ (dphideta[i][p].slice(1))*detady_map[p];
div_phi[i][p] = dphix_dx + dphiy_dy;
}
break;
}
case 3:
{
const std::vector<std::vector<RealGradient> >& dphidxi = fe.get_dphidxi();
const std::vector<std::vector<RealGradient> >& dphideta = fe.get_dphideta();
const std::vector<std::vector<RealGradient> >& dphidzeta = fe.get_dphidzeta();
const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();
const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz();
const std::vector<Real>& dzetadx_map = fe.get_fe_map().get_dzetadx();
const std::vector<Real>& dzetady_map = fe.get_fe_map().get_dzetady();
const std::vector<Real>& dzetadz_map = fe.get_fe_map().get_dzetadz();
/*
In 3D: div = dphi_x/dx + dphi_y/dy + dphi_z/dz
*/
for (unsigned int i=0; i<div_phi.size(); i++)
for (unsigned int p=0; p<div_phi[i].size(); p++)
{
Real dphix_dx = (dphidxi[i][p].slice(0))*dxidx_map[p]
+ (dphideta[i][p].slice(0))*detadx_map[p]
+ (dphidzeta[i][p].slice(0))*dzetadx_map[p];
Real dphiy_dy = (dphidxi[i][p].slice(1))*dxidy_map[p]
+ (dphideta[i][p].slice(1))*detady_map[p]
+ (dphidzeta[i][p].slice(1))*dzetady_map[p];
Real dphiz_dz = (dphidxi[i][p].slice(2))*dxidz_map[p]
+ (dphideta[i][p].slice(2))*detadz_map[p]
+ (dphidzeta[i][p].slice(2))*dzetadz_map[p];
div_phi[i][p] = dphix_dx + dphiy_dy + dphiz_dz;
}
break;
}
default: \
libmesh_error_msg("Invalid dim = " << dim); \
} // switch(dim)
return;
}
| void libMesh::H1FETransformation< OutputShape >::map_dphi | ( | const unsigned int | dim, |
| const Elem *const | elem, | ||
| const std::vector< Point > & | qp, | ||
| const FEGenericBase< OutputShape > & | fe, | ||
| std::vector< std::vector< typename FEGenericBase< OutputShape >::OutputGradient > > & | dphi, | ||
| std::vector< std::vector< OutputShape > > & | dphidx, | ||
| std::vector< std::vector< OutputShape > > & | dphidy, | ||
| std::vector< std::vector< OutputShape > > & | dphidz | ||
| ) | const [virtual] |
Evaluates shape function gradients in physical coordinates for H1 conforming elements. dphi/dx = dphi/dxi * dxi/dx, etc.
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 81 of file h1_fe_transformation.C.
References libMesh::FEMap::get_detadx(), libMesh::FEMap::get_detady(), libMesh::FEMap::get_detadz(), libMesh::FEGenericBase< OutputType >::get_dphideta(), libMesh::FEGenericBase< OutputType >::get_dphidxi(), libMesh::FEGenericBase< OutputType >::get_dphidzeta(), libMesh::FEMap::get_dxidx(), libMesh::FEMap::get_dxidy(), libMesh::FEMap::get_dxidz(), libMesh::FEMap::get_dzetadx(), libMesh::FEMap::get_dzetady(), libMesh::FEMap::get_dzetadz(), and libMesh::FEAbstract::get_fe_map().
{
switch(dim)
{
case 0: // No derivatives in 0D
{
for (unsigned int i=0; i<dphi.size(); i++)
for (unsigned int p=0; p<dphi[i].size(); p++)
{
dphi[i][p] = 0.;
}
break;
}
case 1:
{
const std::vector<std::vector<OutputShape> >& dphidxi = fe.get_dphidxi();
const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
#if LIBMESH_DIM>1
const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
#endif
#if LIBMESH_DIM>2
const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();
#endif
for (unsigned int i=0; i<dphi.size(); i++)
for (unsigned int p=0; p<dphi[i].size(); p++)
{
// dphi/dx = (dphi/dxi)*(dxi/dx)
dphi[i][p].slice(0) = dphidx[i][p] = dphidxi[i][p]*dxidx_map[p];
#if LIBMESH_DIM>1
dphi[i][p].slice(1) = dphidy[i][p] = dphidxi[i][p]*dxidy_map[p];
#endif
#if LIBMESH_DIM>2
dphi[i][p].slice(2) = dphidz[i][p] = dphidxi[i][p]*dxidz_map[p];
#endif
}
break;
}
case 2:
{
const std::vector<std::vector<OutputShape> >& dphidxi = fe.get_dphidxi();
const std::vector<std::vector<OutputShape> >& dphideta = fe.get_dphideta();
const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
#if LIBMESH_DIM > 2
const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();
#endif
const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
#if LIBMESH_DIM > 2
const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz();
#endif
for (unsigned int i=0; i<dphi.size(); i++)
for (unsigned int p=0; p<dphi[i].size(); p++)
{
// dphi/dx = (dphi/dxi)*(dxi/dx) + (dphi/deta)*(deta/dx)
dphi[i][p].slice(0) = dphidx[i][p] = (dphidxi[i][p]*dxidx_map[p] +
dphideta[i][p]*detadx_map[p]);
// dphi/dy = (dphi/dxi)*(dxi/dy) + (dphi/deta)*(deta/dy)
dphi[i][p].slice(1) = dphidy[i][p] = (dphidxi[i][p]*dxidy_map[p] +
dphideta[i][p]*detady_map[p]);
#if LIBMESH_DIM > 2
// dphi/dz = (dphi/dxi)*(dxi/dz) + (dphi/deta)*(deta/dz)
dphi[i][p].slice(2) = dphidz[i][p] = (dphidxi[i][p]*dxidz_map[p] +
dphideta[i][p]*detadz_map[p]);
#endif
}
break;
}
case 3:
{
const std::vector<std::vector<OutputShape> >& dphidxi = fe.get_dphidxi();
const std::vector<std::vector<OutputShape> >& dphideta = fe.get_dphideta();
const std::vector<std::vector<OutputShape> >& dphidzeta = fe.get_dphidzeta();
const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx();
const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy();
const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz();
const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx();
const std::vector<Real>& detady_map = fe.get_fe_map().get_detady();
const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz();
const std::vector<Real>& dzetadx_map = fe.get_fe_map().get_dzetadx();
const std::vector<Real>& dzetady_map = fe.get_fe_map().get_dzetady();
const std::vector<Real>& dzetadz_map = fe.get_fe_map().get_dzetadz();
for (unsigned int i=0; i<dphi.size(); i++)
for (unsigned int p=0; p<dphi[i].size(); p++)
{
// dphi/dx = (dphi/dxi)*(dxi/dx) + (dphi/deta)*(deta/dx) + (dphi/dzeta)*(dzeta/dx);
dphi[i][p].slice(0) = dphidx[i][p] = (dphidxi[i][p]*dxidx_map[p] +
dphideta[i][p]*detadx_map[p] +
dphidzeta[i][p]*dzetadx_map[p]);
// dphi/dy = (dphi/dxi)*(dxi/dy) + (dphi/deta)*(deta/dy) + (dphi/dzeta)*(dzeta/dy);
dphi[i][p].slice(1) = dphidy[i][p] = (dphidxi[i][p]*dxidy_map[p] +
dphideta[i][p]*detady_map[p] +
dphidzeta[i][p]*dzetady_map[p]);
// dphi/dz = (dphi/dxi)*(dxi/dz) + (dphi/deta)*(deta/dz) + (dphi/dzeta)*(dzeta/dz);
dphi[i][p].slice(2) = dphidz[i][p] = (dphidxi[i][p]*dxidz_map[p] +
dphideta[i][p]*detadz_map[p] +
dphidzeta[i][p]*dzetadz_map[p]);
}
break;
}
default:
libmesh_error_msg("Invalid dim = " << dim);
} // switch(dim)
}
| void libMesh::H1FETransformation< OutputShape >::map_phi | ( | const unsigned int | dim, |
| const Elem * const | elem, | ||
| const std::vector< Point > & | qp, | ||
| const FEGenericBase< OutputShape > & | fe, | ||
| std::vector< std::vector< OutputShape > > & | phi | ||
| ) | const [virtual] |
Evaluates shape functions in physical coordinates for H1 conforming elements. In this case
Implements libMesh::FETransformationBase< OutputShape >.
Definition at line 25 of file h1_fe_transformation.C.
References libMesh::FEAbstract::get_fe_type().
{
switch(dim)
{
case 0:
{
for (unsigned int i=0; i<phi.size(); i++)
{
libmesh_assert_equal_to ( qp.size(), phi[i].size() );
for (unsigned int p=0; p<phi[i].size(); p++)
FEInterface::shape<OutputShape>(0, fe.get_fe_type(), elem, i, qp[p], phi[i][p]);
}
break;
}
case 1:
{
for (unsigned int i=0; i<phi.size(); i++)
{
libmesh_assert_equal_to ( qp.size(), phi[i].size() );
for (unsigned int p=0; p<phi[i].size(); p++)
FEInterface::shape<OutputShape>(1, fe.get_fe_type(), elem, i, qp[p], phi[i][p]);
}
break;
}
case 2:
{
for (unsigned int i=0; i<phi.size(); i++)
{
libmesh_assert_equal_to ( qp.size(), phi[i].size() );
for (unsigned int p=0; p<phi[i].size(); p++)
FEInterface::shape<OutputShape>(2, fe.get_fe_type(), elem, i, qp[p], phi[i][p]);
}
break;
}
case 3:
{
for (unsigned int i=0; i<phi.size(); i++)
{
libmesh_assert_equal_to ( qp.size(), phi[i].size() );
for (unsigned int p=0; p<phi[i].size(); p++)
FEInterface::shape<OutputShape>(3, fe.get_fe_type(), elem, i, qp[p], phi[i][p]);
}
break;
}
default:
libmesh_error_msg("Invalid dim = " << dim);
}
}