cell_tet4.h
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 00020 #ifndef LIBMESH_CELL_TET4_H 00021 #define LIBMESH_CELL_TET4_H 00022 00023 // Local includes 00024 #include "libmesh/cell_tet.h" 00025 00026 // C++ includes 00027 #include <cstddef> 00028 00029 namespace libMesh 00030 { 00031 00032 00033 00034 00054 // ------------------------------------------------------------ 00055 // Tet4 class definition 00056 class Tet4 : public Tet 00057 { 00058 public: 00059 00063 explicit 00064 Tet4 (Elem* p=NULL); 00065 00069 ElemType type () const { return TET4; } 00070 00074 unsigned int n_nodes() const { return 4; } 00075 00079 unsigned int n_sub_elem() const { return 1; } 00080 00084 virtual bool is_vertex(const unsigned int i) const; 00085 00089 virtual bool is_edge(const unsigned int i) const; 00090 00094 virtual bool is_face(const unsigned int i) const; 00095 00096 /* 00097 * @returns true iff the specified (local) node number is on the 00098 * specified side 00099 */ 00100 virtual bool is_node_on_side(const unsigned int n, 00101 const unsigned int s) const; 00102 00103 /* 00104 * @returns true iff the specified (local) node number is on the 00105 * specified edge 00106 */ 00107 virtual bool is_node_on_edge(const unsigned int n, 00108 const unsigned int e) const; 00109 00110 /* 00111 * @returns true iff the specified child is on the 00112 * specified side 00113 */ 00114 virtual bool is_child_on_side(const unsigned int c, 00115 const unsigned int s) const; 00116 00117 /* 00118 * @returns true iff the element map is definitely affine within 00119 * numerical tolerances 00120 */ 00121 virtual bool has_affine_map () const { return true; } 00122 00127 virtual bool is_linear () const { return true; } 00128 00132 Order default_order() const { return FIRST; } 00133 00138 AutoPtr<Elem> build_side (const unsigned int i, 00139 bool proxy) const; 00140 00145 AutoPtr<Elem> build_edge (const unsigned int i) const; 00146 00147 virtual void connectivity(const unsigned int sc, 00148 const IOPackage iop, 00149 std::vector<dof_id_type>& conn) const; 00150 00155 static const unsigned int side_nodes_map[4][3]; 00156 00161 static const unsigned int edge_nodes_map[6][2]; 00162 00167 virtual Real volume () const; 00168 00177 std::pair<Real, Real> min_and_max_angle() const; 00178 00179 protected: 00180 00184 Node* _nodelinks_data[4]; 00185 00186 00187 00188 #ifdef LIBMESH_ENABLE_AMR 00189 00193 float embedding_matrix (const unsigned int i, 00194 const unsigned int j, 00195 const unsigned int k) const; 00196 00201 static const float _embedding_matrix[8][4][4]; 00202 00203 // public: 00204 // 00205 // /** 00206 // * Allows the user to reselect the diagonal after refinement. This 00207 // * function may only be called directly after the element is refined 00208 // * for the first time (and before the \p EquationSystems::reinit() 00209 // * is called). It will destroy and re-create the children if 00210 // * necessary. 00211 // */ 00212 // void reselect_diagonal (const Diagonal diag); 00213 // 00214 // /** 00215 // * Reselects the diagonal after refinement to be the optimal one. 00216 // * This makes sense if the user has moved some grid points, so that 00217 // * the former optimal choice is no longer optimal. Also, the user 00218 // * may exclude one diagonal from this selection by giving it as 00219 // * argument. In this case, the more optimal one of the remaining 00220 // * two diagonals is chosen. 00221 // */ 00222 // void reselect_optimal_diagonal (const Diagonal exclude_this=INVALID_DIAG); 00223 00224 #endif 00225 00226 }; 00227 00228 00229 00230 // ------------------------------------------------------------ 00231 // Tet4 class member functions 00232 inline 00233 Tet4::Tet4(Elem* p) : 00234 Tet(Tet4::n_nodes(), p, _nodelinks_data) 00235 { 00236 } 00237 00238 } // namespace libMesh 00239 00240 00241 #endif // LIBMESH_CELL_TET4_H
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:45 UTC
Hosted By: