libMesh::FEXYZMap Class Reference
#include <fe_xyz_map.h>

Public Member Functions | |
| FEXYZMap () | |
| virtual | ~FEXYZMap () |
| virtual void | compute_face_map (int dim, const std::vector< Real > &qw, const Elem *side) |
| template<unsigned int Dim> | |
| void | init_reference_to_physical_map (const std::vector< Point > &qp, const Elem *elem) |
| void | compute_single_point_map (const unsigned int dim, const std::vector< Real > &qw, const Elem *elem, unsigned int p) |
| virtual void | compute_affine_map (const unsigned int dim, const std::vector< Real > &qw, const Elem *elem) |
| virtual void | compute_map (const unsigned int dim, const std::vector< Real > &qw, const Elem *elem) |
| void | compute_edge_map (int dim, const std::vector< Real > &qw, const Elem *side) |
| template<unsigned int Dim> | |
| void | init_face_shape_functions (const std::vector< Point > &qp, const Elem *side) |
| template<unsigned int Dim> | |
| void | init_edge_shape_functions (const std::vector< Point > &qp, const Elem *edge) |
| const std::vector< Point > & | get_xyz () const |
| const std::vector< Real > & | get_jacobian () const |
| const std::vector< Real > & | get_JxW () const |
| std::vector< Real > & | get_JxW () |
| const std::vector< RealGradient > & | get_dxyzdxi () const |
| const std::vector< RealGradient > & | get_dxyzdeta () const |
| const std::vector< RealGradient > & | get_dxyzdzeta () const |
| const std::vector< RealGradient > & | get_d2xyzdxi2 () const |
| const std::vector< RealGradient > & | get_d2xyzdeta2 () const |
| const std::vector< RealGradient > & | get_d2xyzdzeta2 () const |
| const std::vector< RealGradient > & | get_d2xyzdxideta () const |
| const std::vector< RealGradient > & | get_d2xyzdxidzeta () const |
| const std::vector< RealGradient > & | get_d2xyzdetadzeta () const |
| const std::vector< Real > & | get_dxidx () const |
| const std::vector< Real > & | get_dxidy () const |
| const std::vector< Real > & | get_dxidz () const |
| const std::vector< Real > & | get_detadx () const |
| const std::vector< Real > & | get_detady () const |
| const std::vector< Real > & | get_detadz () const |
| const std::vector< Real > & | get_dzetadx () const |
| const std::vector< Real > & | get_dzetady () const |
| const std::vector< Real > & | get_dzetadz () const |
| const std::vector< std::vector < Real > > & | get_psi () const |
| std::vector< std::vector< Real > > & | get_psi () |
| const std::vector< std::vector < Real > > & | get_phi_map () const |
| std::vector< std::vector< Real > > & | get_phi_map () |
| const std::vector< std::vector < Real > > & | get_dphidxi_map () const |
| std::vector< std::vector< Real > > & | get_dphidxi_map () |
| const std::vector< std::vector < Real > > & | get_dphideta_map () const |
| std::vector< std::vector< Real > > & | get_dphideta_map () |
| const std::vector< std::vector < Real > > & | get_dphidzeta_map () const |
| std::vector< std::vector< Real > > & | get_dphidzeta_map () |
| const std::vector< std::vector < Point > > & | get_tangents () const |
| const std::vector< Point > & | get_normals () const |
| const std::vector< Real > & | get_curvatures () const |
| void | print_JxW (std::ostream &os) const |
| void | print_xyz (std::ostream &os) const |
| std::vector< std::vector< Real > > & | get_dpsidxi () |
| std::vector< std::vector< Real > > & | get_dpsideta () |
| std::vector< std::vector< Real > > & | get_d2psidxi2 () |
| std::vector< std::vector< Real > > & | get_d2psidxideta () |
| std::vector< std::vector< Real > > & | get_d2psideta2 () |
| std::vector< std::vector< Real > > & | get_d2phidxi2_map () |
| std::vector< std::vector< Real > > & | get_d2phidxideta_map () |
| std::vector< std::vector< Real > > & | get_d2phidxidzeta_map () |
| std::vector< std::vector< Real > > & | get_d2phideta2_map () |
| std::vector< std::vector< Real > > & | get_d2phidetadzeta_map () |
| std::vector< std::vector< Real > > & | get_d2phidzeta2_map () |
Static Public Member Functions | |
| static AutoPtr< FEMap > | build (FEType fe_type) |
Protected Member Functions | |
| void | resize_quadrature_map_vectors (const unsigned int dim, unsigned int n_qp) |
| Real | dxdxi_map (const unsigned int p) const |
| Real | dydxi_map (const unsigned int p) const |
| Real | dzdxi_map (const unsigned int p) const |
| Real | dxdeta_map (const unsigned int p) const |
| Real | dydeta_map (const unsigned int p) const |
| Real | dzdeta_map (const unsigned int p) const |
| Real | dxdzeta_map (const unsigned int p) const |
| Real | dydzeta_map (const unsigned int p) const |
| Real | dzdzeta_map (const unsigned int p) const |
Protected Attributes | |
| std::vector< Point > | xyz |
| std::vector< RealGradient > | dxyzdxi_map |
| std::vector< RealGradient > | dxyzdeta_map |
| std::vector< RealGradient > | dxyzdzeta_map |
| std::vector< RealGradient > | d2xyzdxi2_map |
| std::vector< RealGradient > | d2xyzdxideta_map |
| std::vector< RealGradient > | d2xyzdeta2_map |
| std::vector< RealGradient > | d2xyzdxidzeta_map |
| std::vector< RealGradient > | d2xyzdetadzeta_map |
| std::vector< RealGradient > | d2xyzdzeta2_map |
| std::vector< Real > | dxidx_map |
| std::vector< Real > | dxidy_map |
| std::vector< Real > | dxidz_map |
| std::vector< Real > | detadx_map |
| std::vector< Real > | detady_map |
| std::vector< Real > | detadz_map |
| std::vector< Real > | dzetadx_map |
| std::vector< Real > | dzetady_map |
| std::vector< Real > | dzetadz_map |
| std::vector< std::vector< Real > > | phi_map |
| std::vector< std::vector< Real > > | dphidxi_map |
| std::vector< std::vector< Real > > | dphideta_map |
| std::vector< std::vector< Real > > | dphidzeta_map |
| std::vector< std::vector< Real > > | d2phidxi2_map |
| std::vector< std::vector< Real > > | d2phidxideta_map |
| std::vector< std::vector< Real > > | d2phidxidzeta_map |
| std::vector< std::vector< Real > > | d2phideta2_map |
| std::vector< std::vector< Real > > | d2phidetadzeta_map |
| std::vector< std::vector< Real > > | d2phidzeta2_map |
| std::vector< std::vector< Real > > | psi_map |
| std::vector< std::vector< Real > > | dpsidxi_map |
| std::vector< std::vector< Real > > | dpsideta_map |
| std::vector< std::vector< Real > > | d2psidxi2_map |
| std::vector< std::vector< Real > > | d2psidxideta_map |
| std::vector< std::vector< Real > > | d2psideta2_map |
| std::vector< std::vector< Point > > | tangents |
| std::vector< Point > | normals |
| std::vector< Real > | curvatures |
| std::vector< Real > | jac |
| std::vector< Real > | JxW |
Detailed Description
Definition at line 30 of file fe_xyz_map.h.
Constructor & Destructor Documentation
| libMesh::FEXYZMap::FEXYZMap | ( | ) | [inline] |
Definition at line 34 of file fe_xyz_map.h.
00035 : FEMap(){};
| virtual libMesh::FEXYZMap::~FEXYZMap | ( | ) | [inline, virtual] |
Definition at line 37 of file fe_xyz_map.h.
Member Function Documentation
Definition at line 36 of file fe_map.C.
References libMesh::FEType::family, and libMeshEnums::XYZ.
00037 { 00038 switch( fe_type.family ) 00039 { 00040 case XYZ: 00041 { 00042 AutoPtr<FEMap> ap( new FEXYZMap ); 00043 return ap; 00044 } 00045 default: 00046 { 00047 AutoPtr<FEMap> ap( new FEMap ); 00048 return ap; 00049 } 00050 } 00051 00052 //Shouldn't ever get here 00053 libmesh_error(); 00054 return AutoPtr<FEMap>(); 00055 }
| void libMesh::FEMap::compute_affine_map | ( | const unsigned int | dim, | |
| const std::vector< Real > & | qw, | |||
| const Elem * | elem | |||
| ) | [virtual, inherited] |
Compute the jacobian and some other additional data fields. Takes the integration weights as input, along with a pointer to the element. The element is assumed to have a constant Jacobian
Definition at line 704 of file fe_map.C.
References libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::d2xyzdeta2_map, libMesh::FEMap::d2xyzdetadzeta_map, libMesh::FEMap::d2xyzdxi2_map, libMesh::FEMap::d2xyzdxideta_map, libMesh::FEMap::d2xyzdxidzeta_map, libMesh::FEMap::d2xyzdzeta2_map, libMesh::FEMap::detadx_map, libMesh::FEMap::detady_map, libMesh::FEMap::detadz_map, libMesh::FEMap::dxidx_map, libMesh::FEMap::dxidy_map, libMesh::FEMap::dxidz_map, libMesh::FEMap::dxyzdeta_map, libMesh::FEMap::dxyzdxi_map, libMesh::FEMap::dxyzdzeta_map, libMesh::FEMap::dzetadx_map, libMesh::FEMap::dzetady_map, libMesh::FEMap::dzetadz_map, libMesh::FEMap::jac, libMesh::FEMap::JxW, libMesh::FEMap::phi_map, libMesh::Elem::point(), libMesh::FEMap::resize_quadrature_map_vectors(), and libMesh::FEMap::xyz.
Referenced by libMesh::FEMap::compute_map().
00707 { 00708 // Start logging the map computation. 00709 START_LOG("compute_affine_map()", "FEMap"); 00710 00711 libmesh_assert(elem); 00712 00713 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qw.size()); 00714 00715 // Resize the vectors to hold data at the quadrature points 00716 this->resize_quadrature_map_vectors(dim, n_qp); 00717 00718 // Compute map at quadrature point 0 00719 this->compute_single_point_map(dim, qw, elem, 0); 00720 00721 // Compute xyz at all other quadrature points 00722 for (unsigned int p=1; p<n_qp; p++) 00723 { 00724 xyz[p].zero(); 00725 for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes 00726 xyz[p].add_scaled (elem->point(i), phi_map[i][p] ); 00727 } 00728 00729 // Copy other map data from quadrature point 0 00730 for (unsigned int p=1; p<n_qp; p++) // for each extra quadrature point 00731 { 00732 dxyzdxi_map[p] = dxyzdxi_map[0]; 00733 dxidx_map[p] = dxidx_map[0]; 00734 dxidy_map[p] = dxidy_map[0]; 00735 dxidz_map[p] = dxidz_map[0]; 00736 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00737 // The map should be affine, so second derivatives are zero 00738 d2xyzdxi2_map[p] = 0.; 00739 #endif 00740 if (dim > 1) 00741 { 00742 dxyzdeta_map[p] = dxyzdeta_map[0]; 00743 detadx_map[p] = detadx_map[0]; 00744 detady_map[p] = detady_map[0]; 00745 detadz_map[p] = detadz_map[0]; 00746 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00747 d2xyzdxideta_map[p] = 0.; 00748 d2xyzdeta2_map[p] = 0.; 00749 #endif 00750 if (dim > 2) 00751 { 00752 dxyzdzeta_map[p] = dxyzdzeta_map[0]; 00753 dzetadx_map[p] = dzetadx_map[0]; 00754 dzetady_map[p] = dzetady_map[0]; 00755 dzetadz_map[p] = dzetadz_map[0]; 00756 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00757 d2xyzdxidzeta_map[p] = 0.; 00758 d2xyzdetadzeta_map[p] = 0.; 00759 d2xyzdzeta2_map[p] = 0.; 00760 #endif 00761 } 00762 } 00763 jac[p] = jac[0]; 00764 JxW[p] = JxW[0] / qw[0] * qw[p]; 00765 } 00766 00767 STOP_LOG("compute_affine_map()", "FEMap"); 00768 }
| void libMesh::FEMap::compute_edge_map | ( | int | dim, | |
| const std::vector< Real > & | qw, | |||
| const Elem * | side | |||
| ) | [inherited] |
Same as before, but for an edge. Useful for some projections.
Definition at line 805 of file fe_boundary.C.
References libMesh::FEMap::compute_face_map(), libMesh::FEMap::curvatures, libMesh::FEMap::d2psidxi2_map, libMesh::FEMap::d2xyzdeta2_map, libMesh::FEMap::d2xyzdxi2_map, libMesh::FEMap::d2xyzdxideta_map, libMesh::FEMap::dpsidxi_map, libMesh::FEMap::dxdxi_map(), libMesh::FEMap::dxyzdeta_map, libMesh::FEMap::dxyzdxi_map, libMesh::FEMap::dydxi_map(), libMesh::FEMap::dzdxi_map(), libMesh::FEMap::JxW, libMesh::FEMap::normals, libMesh::Elem::point(), libMesh::FEMap::psi_map, libMesh::Real, libMesh::FEMap::tangents, and libMesh::FEMap::xyz.
00808 { 00809 libmesh_assert(edge); 00810 00811 if (dim == 2) 00812 { 00813 // A 2D finite element living in either 2D or 3D space. 00814 // The edges here are the sides of the element, so the 00815 // (misnamed) compute_face_map function does what we want 00816 this->compute_face_map(dim, qw, edge); 00817 return; 00818 } 00819 00820 libmesh_assert_equal_to (dim, 3); // 1D is unnecessary and currently unsupported 00821 00822 START_LOG("compute_edge_map()", "FEMap"); 00823 00824 // The number of quadrature points. 00825 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qw.size()); 00826 00827 // Resize the vectors to hold data at the quadrature points 00828 this->xyz.resize(n_qp); 00829 this->dxyzdxi_map.resize(n_qp); 00830 this->dxyzdeta_map.resize(n_qp); 00831 this->d2xyzdxi2_map.resize(n_qp); 00832 this->d2xyzdxideta_map.resize(n_qp); 00833 this->d2xyzdeta2_map.resize(n_qp); 00834 this->tangents.resize(n_qp); 00835 this->normals.resize(n_qp); 00836 this->curvatures.resize(n_qp); 00837 00838 this->JxW.resize(n_qp); 00839 00840 // Clear the entities that will be summed 00841 for (unsigned int p=0; p<n_qp; p++) 00842 { 00843 this->tangents[p].resize(1); 00844 this->xyz[p].zero(); 00845 this->dxyzdxi_map[p].zero(); 00846 this->dxyzdeta_map[p].zero(); 00847 this->d2xyzdxi2_map[p].zero(); 00848 this->d2xyzdxideta_map[p].zero(); 00849 this->d2xyzdeta2_map[p].zero(); 00850 } 00851 00852 // compute x, dxdxi at the quadrature points 00853 for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes 00854 { 00855 const Point& edge_point = edge->point(i); 00856 00857 for (unsigned int p=0; p<n_qp; p++) // for each quadrature point... 00858 { 00859 this->xyz[p].add_scaled (edge_point, this->psi_map[i][p]); 00860 this->dxyzdxi_map[p].add_scaled (edge_point, this->dpsidxi_map[i][p]); 00861 this->d2xyzdxi2_map[p].add_scaled (edge_point, this->d2psidxi2_map[i][p]); 00862 } 00863 } 00864 00865 // Compute the tangents at the quadrature point 00866 // FIXME: normals (plural!) and curvatures are uncalculated 00867 for (unsigned int p=0; p<n_qp; p++) 00868 { 00869 const Point n = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]); 00870 this->tangents[p][0] = this->dxyzdxi_map[p].unit(); 00871 00872 // compute the jacobian at the quadrature points 00873 const Real the_jac = std::sqrt(this->dxdxi_map(p)*this->dxdxi_map(p) + 00874 this->dydxi_map(p)*this->dydxi_map(p) + 00875 this->dzdxi_map(p)*this->dzdxi_map(p)); 00876 00877 libmesh_assert_greater (the_jac, 0.); 00878 00879 this->JxW[p] = the_jac*qw[p]; 00880 } 00881 00882 STOP_LOG("compute_edge_map()", "FEMap"); 00883 }
| void libMesh::FEXYZMap::compute_face_map | ( | int | dim, | |
| const std::vector< Real > & | qw, | |||
| const Elem * | side | |||
| ) | [virtual] |
Special implementation for XYZ finite elements
Reimplemented from libMesh::FEMap.
Definition at line 22 of file fe_xyz_map.C.
References libMesh::TypeVector< T >::cross(), libMesh::FEMap::curvatures, libMesh::FEMap::d2psideta2_map, libMesh::FEMap::d2psidxi2_map, libMesh::FEMap::d2psidxideta_map, libMesh::FEMap::d2xyzdeta2_map, libMesh::FEMap::d2xyzdxi2_map, libMesh::FEMap::d2xyzdxideta_map, libMesh::FEMap::dpsideta_map, libMesh::FEMap::dpsidxi_map, libMesh::FEMap::dxdeta_map(), libMesh::FEMap::dxdxi_map(), libMesh::FEMap::dxyzdeta_map, libMesh::FEMap::dxyzdxi_map, libMesh::FEMap::dydeta_map(), libMesh::FEMap::dydxi_map(), libMesh::FEMap::dzdeta_map(), libMesh::FEMap::dzdxi_map(), libMesh::FEMap::JxW, libMesh::FEMap::normals, libMesh::Elem::point(), libMesh::FEMap::psi_map, libMesh::Real, libMesh::FEMap::tangents, libMesh::TypeVector< T >::unit(), and libMesh::FEMap::xyz.
00023 { 00024 libmesh_assert(side); 00025 00026 START_LOG("compute_face_map()", "FEXYZMap"); 00027 00028 // The number of quadrature points. 00029 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qw.size()); 00030 00031 switch(dim) 00032 { 00033 case 2: 00034 { 00035 00036 // Resize the vectors to hold data at the quadrature points 00037 { 00038 this->xyz.resize(n_qp); 00039 this->dxyzdxi_map.resize(n_qp); 00040 this->d2xyzdxi2_map.resize(n_qp); 00041 this->tangents.resize(n_qp); 00042 this->normals.resize(n_qp); 00043 this->curvatures.resize(n_qp); 00044 00045 this->JxW.resize(n_qp); 00046 } 00047 00048 // Clear the entities that will be summed 00049 // Compute the tangent & normal at the quadrature point 00050 for (unsigned int p=0; p<n_qp; p++) 00051 { 00052 this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D 00053 this->xyz[p].zero(); 00054 this->dxyzdxi_map[p].zero(); 00055 this->d2xyzdxi2_map[p].zero(); 00056 } 00057 00058 // compute x, dxdxi at the quadrature points 00059 for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes 00060 { 00061 const Point& side_point = side->point(i); 00062 00063 for (unsigned int p=0; p<n_qp; p++) // for each quadrature point... 00064 { 00065 this->xyz[p].add_scaled (side_point, this->psi_map[i][p]); 00066 this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]); 00067 this->d2xyzdxi2_map[p].add_scaled(side_point, this->d2psidxi2_map[i][p]); 00068 } 00069 } 00070 00071 // Compute the tangent & normal at the quadrature point 00072 for (unsigned int p=0; p<n_qp; p++) 00073 { 00074 const Point n(this->dxyzdxi_map[p](1), -this->dxyzdxi_map[p](0), 0.); 00075 00076 this->normals[p] = n.unit(); 00077 this->tangents[p][0] = this->dxyzdxi_map[p].unit(); 00078 #if LIBMESH_DIM == 3 // Only good in 3D space 00079 this->tangents[p][1] = this->dxyzdxi_map[p].cross(n).unit(); 00080 #endif 00081 // The curvature is computed via the familiar Frenet formula: 00082 // curvature = [d^2(x) / d (xi)^2] dot [normal] 00083 // For a reference, see: 00084 // F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill, p. 310 00085 // 00086 // Note: The sign convention here is different from the 00087 // 3D case. Concave-upward curves (smiles) have a positive 00088 // curvature. Concave-downward curves (frowns) have a 00089 // negative curvature. Be sure to take that into account! 00090 const Real numerator = this->d2xyzdxi2_map[p] * this->normals[p]; 00091 const Real denominator = this->dxyzdxi_map[p].size_sq(); 00092 libmesh_assert_not_equal_to (denominator, 0); 00093 this->curvatures[p] = numerator / denominator; 00094 } 00095 00096 // compute the jacobian at the quadrature points 00097 for (unsigned int p=0; p<n_qp; p++) 00098 { 00099 const Real the_jac = std::sqrt(this->dxdxi_map(p)*this->dxdxi_map(p) + 00100 this->dydxi_map(p)*this->dydxi_map(p)); 00101 00102 libmesh_assert_greater (the_jac, 0.); 00103 00104 this->JxW[p] = the_jac*qw[p]; 00105 } 00106 00107 break; 00108 } 00109 00110 case 3: 00111 { 00112 // Resize the vectors to hold data at the quadrature points 00113 { 00114 this->xyz.resize(n_qp); 00115 this->dxyzdxi_map.resize(n_qp); 00116 this->dxyzdeta_map.resize(n_qp); 00117 this->d2xyzdxi2_map.resize(n_qp); 00118 this->d2xyzdxideta_map.resize(n_qp); 00119 this->d2xyzdeta2_map.resize(n_qp); 00120 this->tangents.resize(n_qp); 00121 this->normals.resize(n_qp); 00122 this->curvatures.resize(n_qp); 00123 00124 this->JxW.resize(n_qp); 00125 } 00126 00127 // Clear the entities that will be summed 00128 for (unsigned int p=0; p<n_qp; p++) 00129 { 00130 this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D 00131 this->xyz[p].zero(); 00132 this->dxyzdxi_map[p].zero(); 00133 this->dxyzdeta_map[p].zero(); 00134 this->d2xyzdxi2_map[p].zero(); 00135 this->d2xyzdxideta_map[p].zero(); 00136 this->d2xyzdeta2_map[p].zero(); 00137 } 00138 00139 // compute x, dxdxi at the quadrature points 00140 for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes 00141 { 00142 const Point& side_point = side->point(i); 00143 00144 for (unsigned int p=0; p<n_qp; p++) // for each quadrature point... 00145 { 00146 this->xyz[p].add_scaled (side_point, this->psi_map[i][p]); 00147 this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]); 00148 this->dxyzdeta_map[p].add_scaled (side_point, this->dpsideta_map[i][p]); 00149 this->d2xyzdxi2_map[p].add_scaled (side_point, this->d2psidxi2_map[i][p]); 00150 this->d2xyzdxideta_map[p].add_scaled(side_point, this->d2psidxideta_map[i][p]); 00151 this->d2xyzdeta2_map[p].add_scaled (side_point, this->d2psideta2_map[i][p]); 00152 } 00153 } 00154 00155 // Compute the tangents, normal, and curvature at the quadrature point 00156 for (unsigned int p=0; p<n_qp; p++) 00157 { 00158 const Point n = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]); 00159 this->normals[p] = n.unit(); 00160 this->tangents[p][0] = this->dxyzdxi_map[p].unit(); 00161 this->tangents[p][1] = n.cross(this->dxyzdxi_map[p]).unit(); 00162 00163 // Compute curvature using the typical nomenclature 00164 // of the first and second fundamental forms. 00165 // For reference, see: 00166 // 1) http://mathworld.wolfram.com/MeanCurvature.html 00167 // (note -- they are using inward normal) 00168 // 2) F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill 00169 const Real L = -this->d2xyzdxi2_map[p] * this->normals[p]; 00170 const Real M = -this->d2xyzdxideta_map[p] * this->normals[p]; 00171 const Real N = -this->d2xyzdeta2_map[p] * this->normals[p]; 00172 const Real E = this->dxyzdxi_map[p].size_sq(); 00173 const Real F = this->dxyzdxi_map[p] * this->dxyzdeta_map[p]; 00174 const Real G = this->dxyzdeta_map[p].size_sq(); 00175 00176 const Real numerator = E*N -2.*F*M + G*L; 00177 const Real denominator = E*G - F*F; 00178 libmesh_assert_not_equal_to (denominator, 0.); 00179 this->curvatures[p] = 0.5*numerator/denominator; 00180 } 00181 00182 // compute the jacobian at the quadrature points, see 00183 // http://sp81.msi.umn.edu:999/fluent/fidap/help/theory/thtoc.htm 00184 for (unsigned int p=0; p<n_qp; p++) 00185 { 00186 const Real g11 = (this->dxdxi_map(p)*this->dxdxi_map(p) + 00187 this->dydxi_map(p)*this->dydxi_map(p) + 00188 this->dzdxi_map(p)*this->dzdxi_map(p)); 00189 00190 const Real g12 = (this->dxdxi_map(p)*this->dxdeta_map(p) + 00191 this->dydxi_map(p)*this->dydeta_map(p) + 00192 this->dzdxi_map(p)*this->dzdeta_map(p)); 00193 00194 const Real g21 = g12; 00195 00196 const Real g22 = (this->dxdeta_map(p)*this->dxdeta_map(p) + 00197 this->dydeta_map(p)*this->dydeta_map(p) + 00198 this->dzdeta_map(p)*this->dzdeta_map(p)); 00199 00200 00201 const Real the_jac = std::sqrt(g11*g22 - g12*g21); 00202 00203 libmesh_assert_greater (the_jac, 0.); 00204 00205 this->JxW[p] = the_jac*qw[p]; 00206 } 00207 00208 break; 00209 } 00210 default: 00211 libmesh_error(); 00212 00213 } // switch(dim) 00214 00215 STOP_LOG("compute_face_map()", "FEXYZMap"); 00216 00217 return; 00218 }
| void libMesh::FEMap::compute_map | ( | const unsigned int | dim, | |
| const std::vector< Real > & | qw, | |||
| const Elem * | elem | |||
| ) | [virtual, inherited] |
Compute the jacobian and some other additional data fields. Takes the integration weights as input, along with a pointer to the element.
Definition at line 772 of file fe_map.C.
References libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::Elem::has_affine_map(), and libMesh::FEMap::resize_quadrature_map_vectors().
00775 { 00776 if (elem->has_affine_map()) 00777 { 00778 compute_affine_map(dim, qw, elem); 00779 return; 00780 } 00781 00782 // Start logging the map computation. 00783 START_LOG("compute_map()", "FEMap"); 00784 00785 libmesh_assert(elem); 00786 00787 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qw.size()); 00788 00789 // Resize the vectors to hold data at the quadrature points 00790 this->resize_quadrature_map_vectors(dim, n_qp); 00791 00792 // Compute map at all quadrature points 00793 for (unsigned int p=0; p!=n_qp; p++) 00794 this->compute_single_point_map(dim, qw, elem, p); 00795 00796 // Stop logging the map computation. 00797 STOP_LOG("compute_map()", "FEMap"); 00798 }
| void libMesh::FEMap::compute_single_point_map | ( | const unsigned int | dim, | |
| const std::vector< Real > & | qw, | |||
| const Elem * | elem, | |||
| unsigned int | p | |||
| ) | [inherited] |
Compute the jacobian and some other additional data fields at the single point with index p.
Definition at line 308 of file fe_map.C.
References libMesh::FEMap::d2phideta2_map, libMesh::FEMap::d2phidetadzeta_map, libMesh::FEMap::d2phidxi2_map, libMesh::FEMap::d2phidxideta_map, libMesh::FEMap::d2phidxidzeta_map, libMesh::FEMap::d2phidzeta2_map, libMesh::FEMap::d2xyzdeta2_map, libMesh::FEMap::d2xyzdetadzeta_map, libMesh::FEMap::d2xyzdxi2_map, libMesh::FEMap::d2xyzdxideta_map, libMesh::FEMap::d2xyzdxidzeta_map, libMesh::FEMap::d2xyzdzeta2_map, libMesh::FEMap::detadx_map, libMesh::FEMap::detady_map, libMesh::FEMap::detadz_map, libMesh::FEMap::dphideta_map, libMesh::FEMap::dphidxi_map, libMesh::FEMap::dphidzeta_map, libMesh::FEMap::dxdeta_map(), libMesh::FEMap::dxdxi_map(), libMesh::FEMap::dxdzeta_map(), libMesh::FEMap::dxidx_map, libMesh::FEMap::dxidy_map, libMesh::FEMap::dxidz_map, libMesh::FEMap::dxyzdeta_map, libMesh::FEMap::dxyzdxi_map, libMesh::FEMap::dxyzdzeta_map, libMesh::FEMap::dydeta_map(), libMesh::FEMap::dydxi_map(), libMesh::FEMap::dydzeta_map(), libMesh::FEMap::dzdeta_map(), libMesh::FEMap::dzdxi_map(), libMesh::FEMap::dzdzeta_map(), libMesh::FEMap::dzetadx_map, libMesh::FEMap::dzetady_map, libMesh::FEMap::dzetadz_map, libMesh::err, libMesh::DofObject::id(), libMesh::FEMap::jac, libMesh::FEMap::JxW, libMesh::FEMap::phi_map, libMesh::Elem::point(), libMesh::Real, and libMesh::FEMap::xyz.
Referenced by libMesh::FEMap::compute_affine_map(), and libMesh::FEMap::compute_map().
00312 { 00313 libmesh_assert(elem); 00314 00315 switch (dim) 00316 { 00317 //-------------------------------------------------------------------- 00318 // 0D 00319 case 0: 00320 { 00321 xyz[p] = elem->point(0); 00322 jac[p] = 1.0; 00323 JxW[p] = qw[p]; 00324 break; 00325 } 00326 00327 //-------------------------------------------------------------------- 00328 // 1D 00329 case 1: 00330 { 00331 // Clear the entities that will be summed 00332 xyz[p].zero(); 00333 dxyzdxi_map[p].zero(); 00334 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00335 d2xyzdxi2_map[p].zero(); 00336 #endif 00337 00338 // compute x, dx, d2x at the quadrature point 00339 for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes 00340 { 00341 // Reference to the point, helps eliminate 00342 // exessive temporaries in the inner loop 00343 const Point& elem_point = elem->point(i); 00344 00345 xyz[p].add_scaled (elem_point, phi_map[i][p] ); 00346 dxyzdxi_map[p].add_scaled (elem_point, dphidxi_map[i][p]); 00347 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00348 d2xyzdxi2_map[p].add_scaled(elem_point, d2phidxi2_map[i][p]); 00349 #endif 00350 } 00351 00352 // Compute the jacobian 00353 // 00354 // 1D elements can live in 2D or 3D space. 00355 // The transformation matrix from local->global 00356 // coordinates is 00357 // 00358 // T = | dx/dxi | 00359 // | dy/dxi | 00360 // | dz/dxi | 00361 // 00362 // The generalized determinant of T (from the 00363 // so-called "normal" eqns.) is 00364 // jac = "det(T)" = sqrt(det(T'T)) 00365 // 00366 // where T'= transpose of T, so 00367 // 00368 // jac = sqrt( (dx/dxi)^2 + (dy/dxi)^2 + (dz/dxi)^2 ) 00369 jac[p] = dxyzdxi_map[p].size(); 00370 00371 if (jac[p] <= 0.) 00372 { 00373 libMesh::err << "ERROR: negative Jacobian: " 00374 << jac[p] 00375 << " in element " 00376 << elem->id() 00377 << std::endl; 00378 libmesh_error(); 00379 } 00380 00381 // The inverse Jacobian entries also come from the 00382 // generalized inverse of T (see also the 2D element 00383 // living in 3D code). 00384 const Real jacm2 = 1./jac[p]/jac[p]; 00385 dxidx_map[p] = jacm2*dxdxi_map(p); 00386 #if LIBMESH_DIM > 1 00387 dxidy_map[p] = jacm2*dydxi_map(p); 00388 #endif 00389 #if LIBMESH_DIM > 2 00390 dxidz_map[p] = jacm2*dzdxi_map(p); 00391 #endif 00392 00393 JxW[p] = jac[p]*qw[p]; 00394 00395 // done computing the map 00396 break; 00397 } 00398 00399 00400 //-------------------------------------------------------------------- 00401 // 2D 00402 case 2: 00403 { 00404 //------------------------------------------------------------------ 00405 // Compute the (x,y) values at the quadrature points, 00406 // the Jacobian at the quadrature points 00407 00408 xyz[p].zero(); 00409 00410 dxyzdxi_map[p].zero(); 00411 dxyzdeta_map[p].zero(); 00412 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00413 d2xyzdxi2_map[p].zero(); 00414 d2xyzdxideta_map[p].zero(); 00415 d2xyzdeta2_map[p].zero(); 00416 #endif 00417 00418 00419 // compute (x,y) at the quadrature points, derivatives once 00420 for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes 00421 { 00422 // Reference to the point, helps eliminate 00423 // exessive temporaries in the inner loop 00424 const Point& elem_point = elem->point(i); 00425 00426 xyz[p].add_scaled (elem_point, phi_map[i][p] ); 00427 00428 dxyzdxi_map[p].add_scaled (elem_point, dphidxi_map[i][p] ); 00429 dxyzdeta_map[p].add_scaled (elem_point, dphideta_map[i][p]); 00430 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00431 d2xyzdxi2_map[p].add_scaled (elem_point, d2phidxi2_map[i][p]); 00432 d2xyzdxideta_map[p].add_scaled (elem_point, d2phidxideta_map[i][p]); 00433 d2xyzdeta2_map[p].add_scaled (elem_point, d2phideta2_map[i][p]); 00434 #endif 00435 } 00436 00437 // compute the jacobian once 00438 const Real dx_dxi = dxdxi_map(p), 00439 dx_deta = dxdeta_map(p), 00440 dy_dxi = dydxi_map(p), 00441 dy_deta = dydeta_map(p); 00442 00443 #if LIBMESH_DIM == 2 00444 // Compute the Jacobian. This assumes the 2D face 00445 // lives in 2D space 00446 // 00447 // Symbolically, the matrix determinant is 00448 // 00449 // | dx/dxi dx/deta | 00450 // jac = | dy/dxi dy/deta | 00451 // 00452 // jac = dx/dxi*dy/deta - dx/deta*dy/dxi 00453 jac[p] = (dx_dxi*dy_deta - dx_deta*dy_dxi); 00454 00455 if (jac[p] <= 0.) 00456 { 00457 libMesh::err << "ERROR: negative Jacobian: " 00458 << jac[p] 00459 << " in element " 00460 << elem->id() 00461 << std::endl; 00462 libmesh_error(); 00463 } 00464 00465 JxW[p] = jac[p]*qw[p]; 00466 00467 // Compute the shape function derivatives wrt x,y at the 00468 // quadrature points 00469 const Real inv_jac = 1./jac[p]; 00470 00471 dxidx_map[p] = dy_deta*inv_jac; //dxi/dx = (1/J)*dy/deta 00472 dxidy_map[p] = -dx_deta*inv_jac; //dxi/dy = -(1/J)*dx/deta 00473 detadx_map[p] = -dy_dxi* inv_jac; //deta/dx = -(1/J)*dy/dxi 00474 detady_map[p] = dx_dxi* inv_jac; //deta/dy = (1/J)*dx/dxi 00475 00476 dxidz_map[p] = detadz_map[p] = 0.; 00477 #else 00478 00479 const Real dz_dxi = dzdxi_map(p), 00480 dz_deta = dzdeta_map(p); 00481 00482 // Compute the Jacobian. This assumes a 2D face in 00483 // 3D space. 00484 // 00485 // The transformation matrix T from local to global 00486 // coordinates is 00487 // 00488 // | dx/dxi dx/deta | 00489 // T = | dy/dxi dy/deta | 00490 // | dz/dxi dz/deta | 00491 // note det(T' T) = det(T')det(T) = det(T)det(T) 00492 // so det(T) = std::sqrt(det(T' T)) 00493 // 00494 //---------------------------------------------- 00495 // Notes: 00496 // 00497 // dX = R dXi -> R'dX = R'R dXi 00498 // (R^-1)dX = dXi [(R'R)^-1 R']dX = dXi 00499 // 00500 // so R^-1 = (R'R)^-1 R' 00501 // 00502 // and R^-1 R = (R'R)^-1 R'R = I. 00503 // 00504 const Real g11 = (dx_dxi*dx_dxi + 00505 dy_dxi*dy_dxi + 00506 dz_dxi*dz_dxi); 00507 00508 const Real g12 = (dx_dxi*dx_deta + 00509 dy_dxi*dy_deta + 00510 dz_dxi*dz_deta); 00511 00512 const Real g21 = g12; 00513 00514 const Real g22 = (dx_deta*dx_deta + 00515 dy_deta*dy_deta + 00516 dz_deta*dz_deta); 00517 00518 const Real det = (g11*g22 - g12*g21); 00519 00520 if (det <= 0.) 00521 { 00522 libMesh::err << "ERROR: negative Jacobian! " 00523 << " in element " 00524 << elem->id() 00525 << std::endl; 00526 libmesh_error(); 00527 } 00528 00529 const Real inv_det = 1./det; 00530 jac[p] = std::sqrt(det); 00531 00532 JxW[p] = jac[p]*qw[p]; 00533 00534 const Real g11inv = g22*inv_det; 00535 const Real g12inv = -g12*inv_det; 00536 const Real g21inv = -g21*inv_det; 00537 const Real g22inv = g11*inv_det; 00538 00539 dxidx_map[p] = g11inv*dx_dxi + g12inv*dx_deta; 00540 dxidy_map[p] = g11inv*dy_dxi + g12inv*dy_deta; 00541 dxidz_map[p] = g11inv*dz_dxi + g12inv*dz_deta; 00542 00543 detadx_map[p] = g21inv*dx_dxi + g22inv*dx_deta; 00544 detady_map[p] = g21inv*dy_dxi + g22inv*dy_deta; 00545 detadz_map[p] = g21inv*dz_dxi + g22inv*dz_deta; 00546 00547 #endif 00548 // done computing the map 00549 break; 00550 } 00551 00552 00553 00554 //-------------------------------------------------------------------- 00555 // 3D 00556 case 3: 00557 { 00558 //------------------------------------------------------------------ 00559 // Compute the (x,y,z) values at the quadrature points, 00560 // the Jacobian at the quadrature point 00561 00562 // Clear the entities that will be summed 00563 xyz[p].zero (); 00564 dxyzdxi_map[p].zero (); 00565 dxyzdeta_map[p].zero (); 00566 dxyzdzeta_map[p].zero (); 00567 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00568 d2xyzdxi2_map[p].zero(); 00569 d2xyzdxideta_map[p].zero(); 00570 d2xyzdxidzeta_map[p].zero(); 00571 d2xyzdeta2_map[p].zero(); 00572 d2xyzdetadzeta_map[p].zero(); 00573 d2xyzdzeta2_map[p].zero(); 00574 #endif 00575 00576 00577 // compute (x,y,z) at the quadrature points, 00578 // dxdxi, dydxi, dzdxi, 00579 // dxdeta, dydeta, dzdeta, 00580 // dxdzeta, dydzeta, dzdzeta all once 00581 for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes 00582 { 00583 // Reference to the point, helps eliminate 00584 // exessive temporaries in the inner loop 00585 const Point& elem_point = elem->point(i); 00586 00587 xyz[p].add_scaled (elem_point, phi_map[i][p] ); 00588 dxyzdxi_map[p].add_scaled (elem_point, dphidxi_map[i][p] ); 00589 dxyzdeta_map[p].add_scaled (elem_point, dphideta_map[i][p] ); 00590 dxyzdzeta_map[p].add_scaled (elem_point, dphidzeta_map[i][p]); 00591 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00592 d2xyzdxi2_map[p].add_scaled (elem_point, 00593 d2phidxi2_map[i][p]); 00594 d2xyzdxideta_map[p].add_scaled (elem_point, 00595 d2phidxideta_map[i][p]); 00596 d2xyzdxidzeta_map[p].add_scaled (elem_point, 00597 d2phidxidzeta_map[i][p]); 00598 d2xyzdeta2_map[p].add_scaled (elem_point, 00599 d2phideta2_map[i][p]); 00600 d2xyzdetadzeta_map[p].add_scaled (elem_point, 00601 d2phidetadzeta_map[i][p]); 00602 d2xyzdzeta2_map[p].add_scaled (elem_point, 00603 d2phidzeta2_map[i][p]); 00604 #endif 00605 } 00606 00607 // compute the jacobian 00608 const Real 00609 dx_dxi = dxdxi_map(p), dy_dxi = dydxi_map(p), dz_dxi = dzdxi_map(p), 00610 dx_deta = dxdeta_map(p), dy_deta = dydeta_map(p), dz_deta = dzdeta_map(p), 00611 dx_dzeta = dxdzeta_map(p), dy_dzeta = dydzeta_map(p), dz_dzeta = dzdzeta_map(p); 00612 00613 // Symbolically, the matrix determinant is 00614 // 00615 // | dx/dxi dy/dxi dz/dxi | 00616 // jac = | dx/deta dy/deta dz/deta | 00617 // | dx/dzeta dy/dzeta dz/dzeta | 00618 // 00619 // jac = dx/dxi*(dy/deta*dz/dzeta - dz/deta*dy/dzeta) + 00620 // dy/dxi*(dz/deta*dx/dzeta - dx/deta*dz/dzeta) + 00621 // dz/dxi*(dx/deta*dy/dzeta - dy/deta*dx/dzeta) 00622 00623 jac[p] = (dx_dxi*(dy_deta*dz_dzeta - dz_deta*dy_dzeta) + 00624 dy_dxi*(dz_deta*dx_dzeta - dx_deta*dz_dzeta) + 00625 dz_dxi*(dx_deta*dy_dzeta - dy_deta*dx_dzeta)); 00626 00627 if (jac[p] <= 0.) 00628 { 00629 libMesh::err << "ERROR: negative Jacobian: " 00630 << jac[p] 00631 << " in element " 00632 << elem->id() 00633 << std::endl; 00634 libmesh_error(); 00635 } 00636 00637 JxW[p] = jac[p]*qw[p]; 00638 00639 // Compute the shape function derivatives wrt x,y at the 00640 // quadrature points 00641 const Real inv_jac = 1./jac[p]; 00642 00643 dxidx_map[p] = (dy_deta*dz_dzeta - dz_deta*dy_dzeta)*inv_jac; 00644 dxidy_map[p] = (dz_deta*dx_dzeta - dx_deta*dz_dzeta)*inv_jac; 00645 dxidz_map[p] = (dx_deta*dy_dzeta - dy_deta*dx_dzeta)*inv_jac; 00646 00647 detadx_map[p] = (dz_dxi*dy_dzeta - dy_dxi*dz_dzeta )*inv_jac; 00648 detady_map[p] = (dx_dxi*dz_dzeta - dz_dxi*dx_dzeta )*inv_jac; 00649 detadz_map[p] = (dy_dxi*dx_dzeta - dx_dxi*dy_dzeta )*inv_jac; 00650 00651 dzetadx_map[p] = (dy_dxi*dz_deta - dz_dxi*dy_deta )*inv_jac; 00652 dzetady_map[p] = (dz_dxi*dx_deta - dx_dxi*dz_deta )*inv_jac; 00653 dzetadz_map[p] = (dx_dxi*dy_deta - dy_dxi*dx_deta )*inv_jac; 00654 00655 // done computing the map 00656 break; 00657 } 00658 00659 default: 00660 libmesh_error(); 00661 } 00662 }
| Real libMesh::FEMap::dxdeta_map | ( | const unsigned int | p | ) | const [inline, protected, inherited] |
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the x value of the pth entry of the dxzydeta_map.
Definition at line 456 of file fe_map.h.
References libMesh::FEMap::dxyzdeta_map.
Referenced by compute_face_map(), libMesh::FEMap::compute_face_map(), and libMesh::FEMap::compute_single_point_map().
00456 { return dxyzdeta_map[p](0); }
| Real libMesh::FEMap::dxdxi_map | ( | const unsigned int | p | ) | const [inline, protected, inherited] |
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the x value of the pth entry of the dxzydxi_map.
Definition at line 435 of file fe_map.h.
References libMesh::FEMap::dxyzdxi_map.
Referenced by libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), and libMesh::FEMap::compute_single_point_map().
00435 { return dxyzdxi_map[p](0); }
| Real libMesh::FEMap::dxdzeta_map | ( | const unsigned int | p | ) | const [inline, protected, inherited] |
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the x value of the pth entry of the dxzydzeta_map.
Definition at line 477 of file fe_map.h.
References libMesh::FEMap::dxyzdzeta_map.
Referenced by libMesh::FEMap::compute_single_point_map().
00477 { return dxyzdzeta_map[p](0); }
| Real libMesh::FEMap::dydeta_map | ( | const unsigned int | p | ) | const [inline, protected, inherited] |
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the y value of the pth entry of the dxzydeta_map.
Definition at line 463 of file fe_map.h.
References libMesh::FEMap::dxyzdeta_map.
Referenced by compute_face_map(), libMesh::FEMap::compute_face_map(), and libMesh::FEMap::compute_single_point_map().
00463 { return dxyzdeta_map[p](1); }
| Real libMesh::FEMap::dydxi_map | ( | const unsigned int | p | ) | const [inline, protected, inherited] |
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the y value of the pth entry of the dxzydxi_map.
Definition at line 442 of file fe_map.h.
References libMesh::FEMap::dxyzdxi_map.
Referenced by libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), and libMesh::FEMap::compute_single_point_map().
00442 { return dxyzdxi_map[p](1); }
| Real libMesh::FEMap::dydzeta_map | ( | const unsigned int | p | ) | const [inline, protected, inherited] |
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the y value of the pth entry of the dxzydzeta_map.
Definition at line 484 of file fe_map.h.
References libMesh::FEMap::dxyzdzeta_map.
Referenced by libMesh::FEMap::compute_single_point_map().
00484 { return dxyzdzeta_map[p](1); }
| Real libMesh::FEMap::dzdeta_map | ( | const unsigned int | p | ) | const [inline, protected, inherited] |
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the z value of the pth entry of the dxzydeta_map.
Definition at line 470 of file fe_map.h.
References libMesh::FEMap::dxyzdeta_map.
Referenced by compute_face_map(), libMesh::FEMap::compute_face_map(), and libMesh::FEMap::compute_single_point_map().
00470 { return dxyzdeta_map[p](2); }
| Real libMesh::FEMap::dzdxi_map | ( | const unsigned int | p | ) | const [inline, protected, inherited] |
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the z value of the pth entry of the dxzydxi_map.
Definition at line 449 of file fe_map.h.
References libMesh::FEMap::dxyzdxi_map.
Referenced by libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), and libMesh::FEMap::compute_single_point_map().
00449 { return dxyzdxi_map[p](2); }
| Real libMesh::FEMap::dzdzeta_map | ( | const unsigned int | p | ) | const [inline, protected, inherited] |
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected. Returns the z value of the pth entry of the dxzydzeta_map.
Definition at line 491 of file fe_map.h.
References libMesh::FEMap::dxyzdzeta_map.
Referenced by libMesh::FEMap::compute_single_point_map().
00491 { return dxyzdzeta_map[p](2); }
| const std::vector<Real>& libMesh::FEMap::get_curvatures | ( | ) | const [inline, inherited] |
- Returns:
- the curvatures for use in face integration.
Definition at line 298 of file fe_map.h.
References libMesh::FEMap::curvatures.
00299 { return curvatures;}
| std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phideta2_map | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map 2nd derivative
Definition at line 398 of file fe_map.h.
References libMesh::FEMap::d2phideta2_map.
00399 { return d2phideta2_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phidetadzeta_map | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map 2nd derivative
Definition at line 404 of file fe_map.h.
References libMesh::FEMap::d2phidetadzeta_map.
00405 { return d2phidetadzeta_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phidxi2_map | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map 2nd derivative
Definition at line 380 of file fe_map.h.
References libMesh::FEMap::d2phidxi2_map.
00381 { return d2phidxi2_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phidxideta_map | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map 2nd derivative
Definition at line 386 of file fe_map.h.
References libMesh::FEMap::d2phidxideta_map.
00387 { return d2phidxideta_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phidxidzeta_map | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map 2nd derivative
Definition at line 392 of file fe_map.h.
References libMesh::FEMap::d2phidxidzeta_map.
00393 { return d2phidxidzeta_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_d2phidzeta2_map | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map 2nd derivative
Definition at line 410 of file fe_map.h.
References libMesh::FEMap::d2phidzeta2_map.
00411 { return d2phidzeta2_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_d2psideta2 | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map 2nd derivative for the side/edge
Definition at line 349 of file fe_map.h.
References libMesh::FEMap::d2psideta2_map.
00350 { return d2psideta2_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_d2psidxi2 | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map 2nd derivative for the side/edge
Definition at line 337 of file fe_map.h.
References libMesh::FEMap::d2psidxi2_map.
00338 { return d2psidxi2_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_d2psidxideta | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map 2nd derivative for the side/edge
Definition at line 343 of file fe_map.h.
References libMesh::FEMap::d2psidxideta_map.
00344 { return d2psidxideta_map; }
| const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdeta2 | ( | ) | const [inline, inherited] |
- Returns:
- the second partial derivatives in eta.
Definition at line 155 of file fe_map.h.
References libMesh::FEMap::d2xyzdeta2_map.
00156 { return d2xyzdeta2_map; }
| const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdetadzeta | ( | ) | const [inline, inherited] |
- Returns:
- the second partial derivatives in eta-zeta.
Definition at line 185 of file fe_map.h.
References libMesh::FEMap::d2xyzdetadzeta_map.
00186 { return d2xyzdetadzeta_map; }
| const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdxi2 | ( | ) | const [inline, inherited] |
- Returns:
- the second partial derivatives in xi.
Definition at line 149 of file fe_map.h.
References libMesh::FEMap::d2xyzdxi2_map.
00150 { return d2xyzdxi2_map; }
| const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdxideta | ( | ) | const [inline, inherited] |
- Returns:
- the second partial derivatives in xi-eta.
Definition at line 171 of file fe_map.h.
References libMesh::FEMap::d2xyzdxideta_map.
00172 { return d2xyzdxideta_map; }
| const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdxidzeta | ( | ) | const [inline, inherited] |
- Returns:
- the second partial derivatives in xi-zeta.
Definition at line 179 of file fe_map.h.
References libMesh::FEMap::d2xyzdxidzeta_map.
00180 { return d2xyzdxidzeta_map; }
| const std::vector<RealGradient>& libMesh::FEMap::get_d2xyzdzeta2 | ( | ) | const [inline, inherited] |
- Returns:
- the second partial derivatives in zeta.
Definition at line 163 of file fe_map.h.
References libMesh::FEMap::d2xyzdzeta2_map.
00164 { return d2xyzdzeta2_map; }
| const std::vector<Real>& libMesh::FEMap::get_detadx | ( | ) | const [inline, inherited] |
- Returns:
- the deta/dx entry in the transformation matrix from physical to local coordinates.
Definition at line 215 of file fe_map.h.
References libMesh::FEMap::detadx_map.
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
00216 { return detadx_map; }
| const std::vector<Real>& libMesh::FEMap::get_detady | ( | ) | const [inline, inherited] |
- Returns:
- the deta/dy entry in the transformation matrix from physical to local coordinates.
Definition at line 222 of file fe_map.h.
References libMesh::FEMap::detady_map.
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
00223 { return detady_map; }
| const std::vector<Real>& libMesh::FEMap::get_detadz | ( | ) | const [inline, inherited] |
- Returns:
- the deta/dz entry in the transformation matrix from physical to local coordinates.
Definition at line 229 of file fe_map.h.
References libMesh::FEMap::detadz_map.
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
00230 { return detadz_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_dphideta_map | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map derivative
Definition at line 367 of file fe_map.h.
References libMesh::FEMap::dphideta_map.
00368 { return dphideta_map; }
| const std::vector<std::vector<Real> >& libMesh::FEMap::get_dphideta_map | ( | ) | const [inline, inherited] |
- Returns:
- the reference to physical map derivative
Definition at line 274 of file fe_map.h.
References libMesh::FEMap::dphideta_map.
00275 { return dphideta_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_dphidxi_map | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map derivative
Definition at line 361 of file fe_map.h.
References libMesh::FEMap::dphidxi_map.
00362 { return dphidxi_map; }
| const std::vector<std::vector<Real> >& libMesh::FEMap::get_dphidxi_map | ( | ) | const [inline, inherited] |
- Returns:
- the reference to physical map derivative
Definition at line 268 of file fe_map.h.
References libMesh::FEMap::dphidxi_map.
00269 { return dphidxi_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_dphidzeta_map | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map derivative
Definition at line 373 of file fe_map.h.
References libMesh::FEMap::dphidzeta_map.
00374 { return dphidzeta_map; }
| const std::vector<std::vector<Real> >& libMesh::FEMap::get_dphidzeta_map | ( | ) | const [inline, inherited] |
- Returns:
- the reference to physical map derivative
Definition at line 280 of file fe_map.h.
References libMesh::FEMap::dphidzeta_map.
00281 { return dphidzeta_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_dpsideta | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map derivative for the side/edge
Definition at line 331 of file fe_map.h.
References libMesh::FEMap::dpsideta_map.
00332 { return dpsideta_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_dpsidxi | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map derivative for the side/edge
Definition at line 325 of file fe_map.h.
References libMesh::FEMap::dpsidxi_map.
00326 { return dpsidxi_map; }
| const std::vector<Real>& libMesh::FEMap::get_dxidx | ( | ) | const [inline, inherited] |
- Returns:
- the dxi/dx entry in the transformation matrix from physical to local coordinates.
Definition at line 194 of file fe_map.h.
References libMesh::FEMap::dxidx_map.
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
00195 { return dxidx_map; }
| const std::vector<Real>& libMesh::FEMap::get_dxidy | ( | ) | const [inline, inherited] |
- Returns:
- the dxi/dy entry in the transformation matrix from physical to local coordinates.
Definition at line 201 of file fe_map.h.
References libMesh::FEMap::dxidy_map.
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
00202 { return dxidy_map; }
| const std::vector<Real>& libMesh::FEMap::get_dxidz | ( | ) | const [inline, inherited] |
- Returns:
- the dxi/dz entry in the transformation matrix from physical to local coordinates.
Definition at line 208 of file fe_map.h.
References libMesh::FEMap::dxidz_map.
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
00209 { return dxidz_map; }
| const std::vector<RealGradient>& libMesh::FEMap::get_dxyzdeta | ( | ) | const [inline, inherited] |
- Returns:
- the element tangents in eta-direction at the quadrature points.
Definition at line 136 of file fe_map.h.
References libMesh::FEMap::dxyzdeta_map.
Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl().
00137 { return dxyzdeta_map; }
| const std::vector<RealGradient>& libMesh::FEMap::get_dxyzdxi | ( | ) | const [inline, inherited] |
- Returns:
- the element tangents in xi-direction at the quadrature points.
Definition at line 129 of file fe_map.h.
References libMesh::FEMap::dxyzdxi_map.
Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl().
00130 { return dxyzdxi_map; }
| const std::vector<RealGradient>& libMesh::FEMap::get_dxyzdzeta | ( | ) | const [inline, inherited] |
- Returns:
- the element tangents in zeta-direction at the quadrature points.
Definition at line 143 of file fe_map.h.
References libMesh::FEMap::dxyzdzeta_map.
Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl().
00144 { return dxyzdzeta_map; }
| const std::vector<Real>& libMesh::FEMap::get_dzetadx | ( | ) | const [inline, inherited] |
- Returns:
- the dzeta/dx entry in the transformation matrix from physical to local coordinates.
Definition at line 236 of file fe_map.h.
References libMesh::FEMap::dzetadx_map.
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
00237 { return dzetadx_map; }
| const std::vector<Real>& libMesh::FEMap::get_dzetady | ( | ) | const [inline, inherited] |
- Returns:
- the dzeta/dy entry in the transformation matrix from physical to local coordinates.
Definition at line 243 of file fe_map.h.
References libMesh::FEMap::dzetady_map.
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
00244 { return dzetady_map; }
| const std::vector<Real>& libMesh::FEMap::get_dzetadz | ( | ) | const [inline, inherited] |
- Returns:
- the dzeta/dz entry in the transformation matrix from physical to local coordinates.
Definition at line 250 of file fe_map.h.
References libMesh::FEMap::dzetadz_map.
Referenced by libMesh::H1FETransformation< OutputShape >::map_curl(), libMesh::H1FETransformation< OutputShape >::map_d2phi(), libMesh::H1FETransformation< OutputShape >::map_div(), libMesh::H1FETransformation< OutputShape >::map_dphi(), and libMesh::HCurlFETransformation< OutputShape >::map_phi().
00251 { return dzetadz_map; }
| const std::vector<Real>& libMesh::FEMap::get_jacobian | ( | ) | const [inline, inherited] |
- Returns:
- the element Jacobian for each quadrature point.
Definition at line 115 of file fe_map.h.
References libMesh::FEMap::jac.
Referenced by libMesh::HCurlFETransformation< OutputShape >::map_curl().
00116 { return jac; }
| std::vector<Real>& libMesh::FEMap::get_JxW | ( | ) | [inline, inherited] |
- Returns:
- writable reference to the element Jacobian times the quadrature weight for each quadrature point.
Definition at line 420 of file fe_map.h.
References libMesh::FEMap::JxW.
00421 { return JxW; }
| const std::vector<Real>& libMesh::FEMap::get_JxW | ( | ) | const [inline, inherited] |
- Returns:
- the element Jacobian times the quadrature weight for each quadrature point.
Definition at line 122 of file fe_map.h.
References libMesh::FEMap::JxW.
00123 { return JxW; }
| const std::vector<Point>& libMesh::FEMap::get_normals | ( | ) | const [inline, inherited] |
- Returns:
- the normal vectors for face integration.
Definition at line 292 of file fe_map.h.
References libMesh::FEMap::normals.
00293 { return normals; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_phi_map | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map for the element
Definition at line 355 of file fe_map.h.
References libMesh::FEMap::phi_map.
00356 { return phi_map; }
| const std::vector<std::vector<Real> >& libMesh::FEMap::get_phi_map | ( | ) | const [inline, inherited] |
- Returns:
- the reference to physical map for the element
Definition at line 262 of file fe_map.h.
References libMesh::FEMap::phi_map.
00263 { return phi_map; }
| std::vector<std::vector<Real> >& libMesh::FEMap::get_psi | ( | ) | [inline, inherited] |
- Returns:
- the reference to physical map for the side/edge
Definition at line 319 of file fe_map.h.
References libMesh::FEMap::psi_map.
00320 { return psi_map; }
| const std::vector<std::vector<Real> >& libMesh::FEMap::get_psi | ( | ) | const [inline, inherited] |
- Returns:
- the reference to physical map for the side/edge
Definition at line 256 of file fe_map.h.
References libMesh::FEMap::psi_map.
00257 { return psi_map; }
| const std::vector<std::vector<Point> >& libMesh::FEMap::get_tangents | ( | ) | const [inline, inherited] |
- Returns:
- the tangent vectors for face integration.
Definition at line 286 of file fe_map.h.
References libMesh::FEMap::tangents.
00287 { return tangents; }
| const std::vector<Point>& libMesh::FEMap::get_xyz | ( | ) | const [inline, inherited] |
- Returns:
- the
xyzspatial locations of the quadrature points on the element.
Definition at line 109 of file fe_map.h.
References libMesh::FEMap::xyz.
00110 { return xyz; }
| void libMesh::FEMap::init_edge_shape_functions | ( | const std::vector< Point > & | qp, | |
| const Elem * | edge | |||
| ) | [inline, inherited] |
Same as before, but for an edge. This is used for some projection operators.
Start logging the shape function initialization
Stop logging the shape function initialization
Definition at line 473 of file fe_boundary.C.
References libMesh::FEMap::d2psidxi2_map, libMesh::Elem::default_order(), libMesh::FEMap::dpsidxi_map, libMesh::FEMap::psi_map, and libMesh::Elem::type().
00475 { 00476 libmesh_assert(edge); 00477 00481 START_LOG("init_edge_shape_functions()", "FEMap"); 00482 00483 // The element type and order to use in 00484 // the map 00485 const Order mapping_order (edge->default_order()); 00486 const ElemType mapping_elem_type (edge->type()); 00487 00488 // The number of quadrature points. 00489 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qp.size()); 00490 00491 const unsigned int n_mapping_shape_functions = 00492 FE<Dim,LAGRANGE>::n_shape_functions (mapping_elem_type, 00493 mapping_order); 00494 00495 // resize the vectors to hold current data 00496 // Psi are the shape functions used for the FE mapping 00497 this->psi_map.resize (n_mapping_shape_functions); 00498 this->dpsidxi_map.resize (n_mapping_shape_functions); 00499 this->d2psidxi2_map.resize (n_mapping_shape_functions); 00500 00501 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00502 { 00503 // Allocate space to store the values of the shape functions 00504 // and their first and second derivatives at the quadrature points. 00505 this->psi_map[i].resize (n_qp); 00506 this->dpsidxi_map[i].resize (n_qp); 00507 this->d2psidxi2_map[i].resize (n_qp); 00508 00509 // Compute the value of shape function i, and its first and 00510 // second derivatives at quadrature point p 00511 // (Lagrange shape functions are used for the mapping) 00512 for (unsigned int p=0; p<n_qp; p++) 00513 { 00514 this->psi_map[i][p] = FE<1,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00515 this->dpsidxi_map[i][p] = FE<1,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]); 00516 this->d2psidxi2_map[i][p] = FE<1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 0, qp[p]); 00517 } 00518 } 00519 00523 STOP_LOG("init_edge_shape_functions()", "FEMap"); 00524 }
| void libMesh::FEMap::init_face_shape_functions | ( | const std::vector< Point > & | qp, | |
| const Elem * | side | |||
| ) | [inline, inherited] |
Initalizes the reference to physical element map for a side. This is used for boundary integration.
Start logging the shape function initialization
Stop logging the shape function initialization
Definition at line 384 of file fe_boundary.C.
References libMesh::FEMap::d2psideta2_map, libMesh::FEMap::d2psidxi2_map, libMesh::FEMap::d2psidxideta_map, libMesh::Elem::default_order(), libMesh::FEMap::dpsideta_map, libMesh::FEMap::dpsidxi_map, libMesh::FEMap::psi_map, and libMesh::Elem::type().
00386 { 00387 libmesh_assert(side); 00388 00392 START_LOG("init_face_shape_functions()", "FEMap"); 00393 00394 // The element type and order to use in 00395 // the map 00396 const Order mapping_order (side->default_order()); 00397 const ElemType mapping_elem_type (side->type()); 00398 00399 // The number of quadrature points. 00400 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qp.size()); 00401 00402 const unsigned int n_mapping_shape_functions = 00403 FE<Dim,LAGRANGE>::n_shape_functions (mapping_elem_type, 00404 mapping_order); 00405 00406 // resize the vectors to hold current data 00407 // Psi are the shape functions used for the FE mapping 00408 this->psi_map.resize (n_mapping_shape_functions); 00409 00410 if (Dim > 1) 00411 { 00412 this->dpsidxi_map.resize (n_mapping_shape_functions); 00413 this->d2psidxi2_map.resize (n_mapping_shape_functions); 00414 } 00415 00416 if (Dim == 3) 00417 { 00418 this->dpsideta_map.resize (n_mapping_shape_functions); 00419 this->d2psidxideta_map.resize (n_mapping_shape_functions); 00420 this->d2psideta2_map.resize (n_mapping_shape_functions); 00421 } 00422 00423 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00424 { 00425 // Allocate space to store the values of the shape functions 00426 // and their first and second derivatives at the quadrature points. 00427 this->psi_map[i].resize (n_qp); 00428 if (Dim > 1) 00429 { 00430 this->dpsidxi_map[i].resize (n_qp); 00431 this->d2psidxi2_map[i].resize (n_qp); 00432 } 00433 if (Dim == 3) 00434 { 00435 this->dpsideta_map[i].resize (n_qp); 00436 this->d2psidxideta_map[i].resize (n_qp); 00437 this->d2psideta2_map[i].resize (n_qp); 00438 } 00439 00440 // Compute the value of shape function i, and its first and 00441 // second derivatives at quadrature point p 00442 // (Lagrange shape functions are used for the mapping) 00443 for (unsigned int p=0; p<n_qp; p++) 00444 { 00445 this->psi_map[i][p] = FE<Dim-1,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00446 if (Dim > 1) 00447 { 00448 this->dpsidxi_map[i][p] = FE<Dim-1,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]); 00449 this->d2psidxi2_map[i][p] = FE<Dim-1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 0, qp[p]); 00450 } 00451 // libMesh::out << "this->d2psidxi2_map["<<i<<"][p]=" << d2psidxi2_map[i][p] << std::endl; 00452 00453 // If we are in 3D, then our sides are 2D faces. 00454 // For the second derivatives, we must also compute the cross 00455 // derivative d^2() / dxi deta 00456 if (Dim == 3) 00457 { 00458 this->dpsideta_map[i][p] = FE<Dim-1,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 1, qp[p]); 00459 this->d2psidxideta_map[i][p] = FE<Dim-1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 1, qp[p]); 00460 this->d2psideta2_map[i][p] = FE<Dim-1,LAGRANGE>::shape_second_deriv(mapping_elem_type, mapping_order, i, 2, qp[p]); 00461 } 00462 } 00463 } 00464 00465 00469 STOP_LOG("init_face_shape_functions()", "FEMap"); 00470 }
| template void libMesh::FEMap::init_reference_to_physical_map< 3 > | ( | const std::vector< Point > & | qp, | |
| const Elem * | elem | |||
| ) | [inline, inherited] |
Definition at line 58 of file fe_map.C.
References libMesh::FEMap::d2phideta2_map, libMesh::FEMap::d2phidetadzeta_map, libMesh::FEMap::d2phidxi2_map, libMesh::FEMap::d2phidxideta_map, libMesh::FEMap::d2phidxidzeta_map, libMesh::FEMap::d2phidzeta2_map, libMesh::Elem::default_order(), libMesh::FEMap::dphideta_map, libMesh::FEMap::dphidxi_map, libMesh::FEMap::dphidzeta_map, libMesh::Elem::is_linear(), libMesh::FEMap::phi_map, and libMesh::Elem::type().
00060 { 00061 // Start logging the reference->physical map initialization 00062 START_LOG("init_reference_to_physical_map()", "FEMap"); 00063 00064 // The number of quadrature points. 00065 const std::size_t n_qp = qp.size(); 00066 00067 // The element type and order to use in 00068 // the map 00069 const Order mapping_order (elem->default_order()); 00070 const ElemType mapping_elem_type (elem->type()); 00071 00072 // Number of shape functions used to construt the map 00073 // (Lagrange shape functions are used for mapping) 00074 const unsigned int n_mapping_shape_functions = 00075 FE<Dim,LAGRANGE>::n_shape_functions (mapping_elem_type, 00076 mapping_order); 00077 00078 this->phi_map.resize (n_mapping_shape_functions); 00079 if (Dim > 0) 00080 { 00081 this->dphidxi_map.resize (n_mapping_shape_functions); 00082 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00083 this->d2phidxi2_map.resize (n_mapping_shape_functions); 00084 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00085 } 00086 00087 if (Dim > 1) 00088 { 00089 this->dphideta_map.resize (n_mapping_shape_functions); 00090 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00091 this->d2phidxideta_map.resize (n_mapping_shape_functions); 00092 this->d2phideta2_map.resize (n_mapping_shape_functions); 00093 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00094 } 00095 00096 if (Dim > 2) 00097 { 00098 this->dphidzeta_map.resize (n_mapping_shape_functions); 00099 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00100 this->d2phidxidzeta_map.resize (n_mapping_shape_functions); 00101 this->d2phidetadzeta_map.resize (n_mapping_shape_functions); 00102 this->d2phidzeta2_map.resize (n_mapping_shape_functions); 00103 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00104 } 00105 00106 00107 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00108 { 00109 this->phi_map[i].resize (n_qp); 00110 if (Dim > 0) 00111 { 00112 this->dphidxi_map[i].resize (n_qp); 00113 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00114 this->d2phidxi2_map[i].resize (n_qp); 00115 if (Dim > 1) 00116 { 00117 this->d2phidxideta_map[i].resize (n_qp); 00118 this->d2phideta2_map[i].resize (n_qp); 00119 } 00120 if (Dim > 2) 00121 { 00122 this->d2phidxidzeta_map[i].resize (n_qp); 00123 this->d2phidetadzeta_map[i].resize (n_qp); 00124 this->d2phidzeta2_map[i].resize (n_qp); 00125 } 00126 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00127 00128 if (Dim > 1) 00129 this->dphideta_map[i].resize (n_qp); 00130 00131 if (Dim > 2) 00132 this->dphidzeta_map[i].resize (n_qp); 00133 } 00134 } 00135 00136 // Optimize for the *linear* geometric elements case: 00137 bool is_linear = elem->is_linear(); 00138 00139 switch (Dim) 00140 { 00141 00142 //------------------------------------------------------------ 00143 // 0D 00144 case 0: 00145 { 00146 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00147 for (std::size_t p=0; p<n_qp; p++) 00148 this->phi_map[i][p] = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00149 00150 break; 00151 } 00152 00153 //------------------------------------------------------------ 00154 // 1D 00155 case 1: 00156 { 00157 // Compute the value of the mapping shape function i at quadrature point p 00158 // (Lagrange shape functions are used for mapping) 00159 if (is_linear) 00160 { 00161 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00162 { 00163 this->phi_map[i][0] = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[0]); 00164 this->dphidxi_map[i][0] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]); 00165 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00166 this->d2phidxi2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]); 00167 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00168 for (std::size_t p=1; p<n_qp; p++) 00169 { 00170 this->phi_map[i][p] = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00171 this->dphidxi_map[i][p] = this->dphidxi_map[i][0]; 00172 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00173 this->d2phidxi2_map[i][p] = this->d2phidxi2_map[i][0]; 00174 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00175 } 00176 } 00177 } 00178 else 00179 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00180 for (std::size_t p=0; p<n_qp; p++) 00181 { 00182 this->phi_map[i][p] = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00183 this->dphidxi_map[i][p] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]); 00184 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00185 this->d2phidxi2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]); 00186 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00187 } 00188 00189 break; 00190 } 00191 //------------------------------------------------------------ 00192 // 2D 00193 case 2: 00194 { 00195 // Compute the value of the mapping shape function i at quadrature point p 00196 // (Lagrange shape functions are used for mapping) 00197 if (is_linear) 00198 { 00199 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00200 { 00201 this->phi_map[i][0] = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[0]); 00202 this->dphidxi_map[i][0] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]); 00203 this->dphideta_map[i][0] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 1, qp[0]); 00204 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00205 this->d2phidxi2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]); 00206 this->d2phidxideta_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 1, qp[0]); 00207 this->d2phideta2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 2, qp[0]); 00208 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00209 for (std::size_t p=1; p<n_qp; p++) 00210 { 00211 this->phi_map[i][p] = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00212 this->dphidxi_map[i][p] = this->dphidxi_map[i][0]; 00213 this->dphideta_map[i][p] = this->dphideta_map[i][0]; 00214 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00215 this->d2phidxi2_map[i][p] = this->d2phidxi2_map[i][0]; 00216 this->d2phidxideta_map[i][p] = this->d2phidxideta_map[i][0]; 00217 this->d2phideta2_map[i][p] = this->d2phideta2_map[i][0]; 00218 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00219 } 00220 } 00221 } 00222 else 00223 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00224 for (std::size_t p=0; p<n_qp; p++) 00225 { 00226 this->phi_map[i][p] = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00227 this->dphidxi_map[i][p] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]); 00228 this->dphideta_map[i][p] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 1, qp[p]); 00229 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00230 this->d2phidxi2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]); 00231 this->d2phidxideta_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 1, qp[p]); 00232 this->d2phideta2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 2, qp[p]); 00233 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00234 } 00235 00236 break; 00237 } 00238 00239 //------------------------------------------------------------ 00240 // 3D 00241 case 3: 00242 { 00243 // Compute the value of the mapping shape function i at quadrature point p 00244 // (Lagrange shape functions are used for mapping) 00245 if (is_linear) 00246 { 00247 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00248 { 00249 this->phi_map[i][0] = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[0]); 00250 this->dphidxi_map[i][0] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]); 00251 this->dphideta_map[i][0] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 1, qp[0]); 00252 this->dphidzeta_map[i][0] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 2, qp[0]); 00253 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00254 this->d2phidxi2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[0]); 00255 this->d2phidxideta_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 1, qp[0]); 00256 this->d2phideta2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 2, qp[0]); 00257 this->d2phidxideta_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 3, qp[0]); 00258 this->d2phidetadzeta_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 4, qp[0]); 00259 this->d2phidzeta2_map[i][0] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 5, qp[0]); 00260 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00261 for (std::size_t p=1; p<n_qp; p++) 00262 { 00263 this->phi_map[i][p] = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00264 this->dphidxi_map[i][p] = this->dphidxi_map[i][0]; 00265 this->dphideta_map[i][p] = this->dphideta_map[i][0]; 00266 this->dphidzeta_map[i][p] = this->dphidzeta_map[i][0]; 00267 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00268 this->d2phidxi2_map[i][p] = this->d2phidxi2_map[i][0]; 00269 this->d2phidxideta_map[i][p] = this->d2phidxideta_map[i][0]; 00270 this->d2phideta2_map[i][p] = this->d2phideta2_map[i][0]; 00271 this->d2phidxideta_map[i][p] = this->d2phidxideta_map[i][0]; 00272 this->d2phidetadzeta_map[i][p] = this->d2phidetadzeta_map[i][0]; 00273 this->d2phidzeta2_map[i][p] = this->d2phidzeta2_map[i][0]; 00274 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00275 } 00276 } 00277 } 00278 else 00279 for (unsigned int i=0; i<n_mapping_shape_functions; i++) 00280 for (std::size_t p=0; p<n_qp; p++) 00281 { 00282 this->phi_map[i][p] = FE<Dim,LAGRANGE>::shape (mapping_elem_type, mapping_order, i, qp[p]); 00283 this->dphidxi_map[i][p] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]); 00284 this->dphideta_map[i][p] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 1, qp[p]); 00285 this->dphidzeta_map[i][p] = FE<Dim,LAGRANGE>::shape_deriv (mapping_elem_type, mapping_order, i, 2, qp[p]); 00286 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00287 this->d2phidxi2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 0, qp[p]); 00288 this->d2phidxideta_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 1, qp[p]); 00289 this->d2phideta2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 2, qp[p]); 00290 this->d2phidxideta_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 3, qp[p]); 00291 this->d2phidetadzeta_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 4, qp[p]); 00292 this->d2phidzeta2_map[i][p] = FE<Dim,LAGRANGE>::shape_second_deriv (mapping_elem_type, mapping_order, i, 5, qp[p]); 00293 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00294 } 00295 00296 break; 00297 } 00298 00299 default: 00300 libmesh_error(); 00301 } 00302 00303 // Stop logging the reference->physical map initialization 00304 STOP_LOG("init_reference_to_physical_map()", "FEMap"); 00305 return; 00306 }
| void libMesh::FEMap::print_JxW | ( | std::ostream & | os | ) | const [inherited] |
Prints the Jacobian times the weight for each quadrature point.
Definition at line 801 of file fe_map.C.
References libMesh::FEMap::JxW.
| void libMesh::FEMap::print_xyz | ( | std::ostream & | os | ) | const [inherited] |
Prints the spatial location of each quadrature point (on the physical element).
Definition at line 809 of file fe_map.C.
References libMesh::FEMap::xyz.
| void libMesh::FEMap::resize_quadrature_map_vectors | ( | const unsigned int | dim, | |
| unsigned int | n_qp | |||
| ) | [protected, inherited] |
A utility function for use by compute_*_map
Definition at line 665 of file fe_map.C.
References libMesh::FEMap::d2xyzdeta2_map, libMesh::FEMap::d2xyzdetadzeta_map, libMesh::FEMap::d2xyzdxi2_map, libMesh::FEMap::d2xyzdxideta_map, libMesh::FEMap::d2xyzdxidzeta_map, libMesh::FEMap::d2xyzdzeta2_map, libMesh::FEMap::detadx_map, libMesh::FEMap::detady_map, libMesh::FEMap::detadz_map, libMesh::FEMap::dxidx_map, libMesh::FEMap::dxidy_map, libMesh::FEMap::dxidz_map, libMesh::FEMap::dxyzdeta_map, libMesh::FEMap::dxyzdxi_map, libMesh::FEMap::dxyzdzeta_map, libMesh::FEMap::dzetadx_map, libMesh::FEMap::dzetady_map, libMesh::FEMap::dzetadz_map, libMesh::FEMap::jac, libMesh::FEMap::JxW, and libMesh::FEMap::xyz.
Referenced by libMesh::FEMap::compute_affine_map(), and libMesh::FEMap::compute_map().
00666 { 00667 // Resize the vectors to hold data at the quadrature points 00668 xyz.resize(n_qp); 00669 dxyzdxi_map.resize(n_qp); 00670 dxidx_map.resize(n_qp); 00671 dxidy_map.resize(n_qp); // 1D element may live in 2D ... 00672 dxidz_map.resize(n_qp); // ... or 3D 00673 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00674 d2xyzdxi2_map.resize(n_qp); 00675 #endif 00676 if (dim > 1) 00677 { 00678 dxyzdeta_map.resize(n_qp); 00679 detadx_map.resize(n_qp); 00680 detady_map.resize(n_qp); 00681 detadz_map.resize(n_qp); 00682 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00683 d2xyzdxideta_map.resize(n_qp); 00684 d2xyzdeta2_map.resize(n_qp); 00685 #endif 00686 if (dim > 2) 00687 { 00688 dxyzdzeta_map.resize (n_qp); 00689 dzetadx_map.resize (n_qp); 00690 dzetady_map.resize (n_qp); 00691 dzetadz_map.resize (n_qp); 00692 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00693 d2xyzdxidzeta_map.resize(n_qp); 00694 d2xyzdetadzeta_map.resize(n_qp); 00695 d2xyzdzeta2_map.resize(n_qp); 00696 #endif 00697 } 00698 } 00699 00700 jac.resize(n_qp); 00701 JxW.resize(n_qp); 00702 }
Member Data Documentation
std::vector<Real> libMesh::FEMap::curvatures [protected, inherited] |
The mean curvature (= one half the sum of the principal curvatures) on the boundary at the quadrature points. The mean curvature is a scalar value.
Definition at line 719 of file fe_map.h.
Referenced by libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), and libMesh::FEMap::get_curvatures().
std::vector<std::vector<Real> > libMesh::FEMap::d2phideta2_map [protected, inherited] |
Map for the second derivative, d^2(phi)/d(eta)^2.
Definition at line 652 of file fe_map.h.
Referenced by libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2phideta2_map(), and libMesh::FEMap::init_reference_to_physical_map().
std::vector<std::vector<Real> > libMesh::FEMap::d2phidetadzeta_map [protected, inherited] |
Map for the second derivative, d^2(phi)/d(eta)d(zeta).
Definition at line 657 of file fe_map.h.
Referenced by libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2phidetadzeta_map(), and libMesh::FEMap::init_reference_to_physical_map().
std::vector<std::vector<Real> > libMesh::FEMap::d2phidxi2_map [protected, inherited] |
Map for the second derivative, d^2(phi)/d(xi)^2.
Definition at line 637 of file fe_map.h.
Referenced by libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2phidxi2_map(), and libMesh::FEMap::init_reference_to_physical_map().
std::vector<std::vector<Real> > libMesh::FEMap::d2phidxideta_map [protected, inherited] |
Map for the second derivative, d^2(phi)/d(xi)d(eta).
Definition at line 642 of file fe_map.h.
Referenced by libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2phidxideta_map(), and libMesh::FEMap::init_reference_to_physical_map().
std::vector<std::vector<Real> > libMesh::FEMap::d2phidxidzeta_map [protected, inherited] |
Map for the second derivative, d^2(phi)/d(xi)d(zeta).
Definition at line 647 of file fe_map.h.
Referenced by libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2phidxidzeta_map(), and libMesh::FEMap::init_reference_to_physical_map().
std::vector<std::vector<Real> > libMesh::FEMap::d2phidzeta2_map [protected, inherited] |
Map for the second derivative, d^2(phi)/d(zeta)^2.
Definition at line 662 of file fe_map.h.
Referenced by libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2phidzeta2_map(), and libMesh::FEMap::init_reference_to_physical_map().
std::vector<std::vector<Real> > libMesh::FEMap::d2psideta2_map [protected, inherited] |
Map for the second derivatives (in eta) of the side shape functions. Useful for computing the curvature at the quadrature points.
Definition at line 702 of file fe_map.h.
Referenced by compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::get_d2psideta2(), and libMesh::FEMap::init_face_shape_functions().
std::vector<std::vector<Real> > libMesh::FEMap::d2psidxi2_map [protected, inherited] |
Map for the second derivatives (in xi) of the side shape functions. Useful for computing the curvature at the quadrature points.
Definition at line 688 of file fe_map.h.
Referenced by libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::get_d2psidxi2(), libMesh::FEMap::init_edge_shape_functions(), and libMesh::FEMap::init_face_shape_functions().
std::vector<std::vector<Real> > libMesh::FEMap::d2psidxideta_map [protected, inherited] |
Map for the second (cross) derivatives in xi, eta of the side shape functions. Useful for computing the curvature at the quadrature points.
Definition at line 695 of file fe_map.h.
Referenced by compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::get_d2psidxideta(), and libMesh::FEMap::init_face_shape_functions().
std::vector<RealGradient> libMesh::FEMap::d2xyzdeta2_map [protected, inherited] |
Vector of second partial derivatives in eta: d^2(x)/d(eta)^2
Definition at line 532 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2xyzdeta2(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<RealGradient> libMesh::FEMap::d2xyzdetadzeta_map [protected, inherited] |
Vector of mixed second partial derivatives in eta-zeta: d^2(x)/d(eta)d(zeta) d^2(y)/d(eta)d(zeta) d^2(z)/d(eta)d(zeta)
Definition at line 546 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2xyzdetadzeta(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<RealGradient> libMesh::FEMap::d2xyzdxi2_map [protected, inherited] |
Vector of second partial derivatives in xi: d^2(x)/d(xi)^2, d^2(y)/d(xi)^2, d^2(z)/d(xi)^2
Definition at line 520 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2xyzdxi2(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<RealGradient> libMesh::FEMap::d2xyzdxideta_map [protected, inherited] |
Vector of mixed second partial derivatives in xi-eta: d^2(x)/d(xi)d(eta) d^2(y)/d(xi)d(eta) d^2(z)/d(xi)d(eta)
Definition at line 526 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2xyzdxideta(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<RealGradient> libMesh::FEMap::d2xyzdxidzeta_map [protected, inherited] |
Vector of second partial derivatives in xi-zeta: d^2(x)/d(xi)d(zeta), d^2(y)/d(xi)d(zeta), d^2(z)/d(xi)d(zeta)
Definition at line 540 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2xyzdxidzeta(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<RealGradient> libMesh::FEMap::d2xyzdzeta2_map [protected, inherited] |
Vector of second partial derivatives in zeta: d^2(x)/d(zeta)^2
Definition at line 552 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_d2xyzdzeta2(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Real> libMesh::FEMap::detadx_map [protected, inherited] |
Map for partial derivatives: d(eta)/d(x). Needed for the Jacobian.
Definition at line 579 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_detadx(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Real> libMesh::FEMap::detady_map [protected, inherited] |
Map for partial derivatives: d(eta)/d(y). Needed for the Jacobian.
Definition at line 585 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_detady(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Real> libMesh::FEMap::detadz_map [protected, inherited] |
Map for partial derivatives: d(eta)/d(z). Needed for the Jacobian.
Definition at line 591 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_detadz(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<std::vector<Real> > libMesh::FEMap::dphideta_map [protected, inherited] |
Map for the derivative, d(phi)/d(eta).
Definition at line 625 of file fe_map.h.
Referenced by libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_dphideta_map(), and libMesh::FEMap::init_reference_to_physical_map().
std::vector<std::vector<Real> > libMesh::FEMap::dphidxi_map [protected, inherited] |
Map for the derivative, d(phi)/d(xi).
Definition at line 620 of file fe_map.h.
Referenced by libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_dphidxi_map(), and libMesh::FEMap::init_reference_to_physical_map().
std::vector<std::vector<Real> > libMesh::FEMap::dphidzeta_map [protected, inherited] |
Map for the derivative, d(phi)/d(zeta).
Definition at line 630 of file fe_map.h.
Referenced by libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_dphidzeta_map(), and libMesh::FEMap::init_reference_to_physical_map().
std::vector<std::vector<Real> > libMesh::FEMap::dpsideta_map [protected, inherited] |
Map for the derivative of the side function, d(psi)/d(eta).
Definition at line 681 of file fe_map.h.
Referenced by compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::get_dpsideta(), and libMesh::FEMap::init_face_shape_functions().
std::vector<std::vector<Real> > libMesh::FEMap::dpsidxi_map [protected, inherited] |
Map for the derivative of the side functions, d(psi)/d(xi).
Definition at line 675 of file fe_map.h.
Referenced by libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::get_dpsidxi(), libMesh::FEMap::init_edge_shape_functions(), and libMesh::FEMap::init_face_shape_functions().
std::vector<Real> libMesh::FEMap::dxidx_map [protected, inherited] |
Map for partial derivatives: d(xi)/d(x). Needed for the Jacobian.
Definition at line 560 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_dxidx(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Real> libMesh::FEMap::dxidy_map [protected, inherited] |
Map for partial derivatives: d(xi)/d(y). Needed for the Jacobian.
Definition at line 566 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_dxidy(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Real> libMesh::FEMap::dxidz_map [protected, inherited] |
Map for partial derivatives: d(xi)/d(z). Needed for the Jacobian.
Definition at line 572 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_dxidz(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<RealGradient> libMesh::FEMap::dxyzdeta_map [protected, inherited] |
Vector of parital derivatives: d(x)/d(eta), d(y)/d(eta), d(z)/d(eta)
Definition at line 508 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::dxdeta_map(), libMesh::FEMap::dydeta_map(), libMesh::FEMap::dzdeta_map(), libMesh::FEMap::get_dxyzdeta(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<RealGradient> libMesh::FEMap::dxyzdxi_map [protected, inherited] |
Vector of parital derivatives: d(x)/d(xi), d(y)/d(xi), d(z)/d(xi)
Definition at line 502 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::dxdxi_map(), libMesh::FEMap::dydxi_map(), libMesh::FEMap::dzdxi_map(), libMesh::FEMap::get_dxyzdxi(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<RealGradient> libMesh::FEMap::dxyzdzeta_map [protected, inherited] |
Vector of parital derivatives: d(x)/d(zeta), d(y)/d(zeta), d(z)/d(zeta)
Definition at line 514 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::dxdzeta_map(), libMesh::FEMap::dydzeta_map(), libMesh::FEMap::dzdzeta_map(), libMesh::FEMap::get_dxyzdzeta(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Real> libMesh::FEMap::dzetadx_map [protected, inherited] |
Map for partial derivatives: d(zeta)/d(x). Needed for the Jacobian.
Definition at line 598 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_dzetadx(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Real> libMesh::FEMap::dzetady_map [protected, inherited] |
Map for partial derivatives: d(zeta)/d(y). Needed for the Jacobian.
Definition at line 604 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_dzetady(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Real> libMesh::FEMap::dzetadz_map [protected, inherited] |
Map for partial derivatives: d(zeta)/d(z). Needed for the Jacobian.
Definition at line 610 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_dzetadz(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Real> libMesh::FEMap::jac [protected, inherited] |
Jacobian values at quadrature points
Definition at line 724 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_jacobian(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Real> libMesh::FEMap::JxW [protected, inherited] |
Jacobian*Weight values at quadrature points
Definition at line 729 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_JxW(), libMesh::FEMap::print_JxW(), and libMesh::FEMap::resize_quadrature_map_vectors().
std::vector<Point> libMesh::FEMap::normals [protected, inherited] |
Normal vectors on boundary at quadrature points
Definition at line 712 of file fe_map.h.
Referenced by libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), and libMesh::FEMap::get_normals().
std::vector<std::vector<Real> > libMesh::FEMap::phi_map [protected, inherited] |
Map for the shape function phi.
Definition at line 615 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_phi_map(), and libMesh::FEMap::init_reference_to_physical_map().
std::vector<std::vector<Real> > libMesh::FEMap::psi_map [protected, inherited] |
Map for the side shape functions, psi.
Definition at line 669 of file fe_map.h.
Referenced by libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::get_psi(), libMesh::FEMap::init_edge_shape_functions(), and libMesh::FEMap::init_face_shape_functions().
std::vector<std::vector<Point> > libMesh::FEMap::tangents [protected, inherited] |
Tangent vectors on boundary at quadrature points.
Definition at line 707 of file fe_map.h.
Referenced by libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), and libMesh::FEMap::get_tangents().
std::vector<Point> libMesh::FEMap::xyz [protected, inherited] |
The spatial locations of the quadrature points
Definition at line 496 of file fe_map.h.
Referenced by libMesh::FEMap::compute_affine_map(), libMesh::FEMap::compute_edge_map(), compute_face_map(), libMesh::FEMap::compute_face_map(), libMesh::FEMap::compute_single_point_map(), libMesh::FEMap::get_xyz(), libMesh::FEMap::print_xyz(), and libMesh::FEMap::resize_quadrature_map_vectors().
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:25 UTC
Hosted By: