$extrastylesheet
cell_tet4.h
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 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 
00053 class Tet4 : public Tet
00054 {
00055 public:
00056 
00060   explicit
00061   Tet4  (Elem* p=NULL);
00062 
00066   ElemType type () const { return TET4; }
00067 
00071   unsigned int n_nodes() const { return 4; }
00072 
00076   unsigned int n_sub_elem() const { return 1; }
00077 
00081   virtual bool is_vertex(const unsigned int i) const;
00082 
00086   virtual bool is_edge(const unsigned int i) const;
00087 
00091   virtual bool is_face(const unsigned int i) const;
00092 
00093   /*
00094    * @returns true iff the specified (local) node number is on the
00095    * specified side
00096    */
00097   virtual bool is_node_on_side(const unsigned int n,
00098                                const unsigned int s) const;
00099 
00100   /*
00101    * @returns true iff the specified (local) node number is on the
00102    * specified edge
00103    */
00104   virtual bool is_node_on_edge(const unsigned int n,
00105                                const unsigned int e) const;
00106 
00107   /*
00108    * @returns true iff the specified child is on the
00109    * specified side
00110    */
00111   virtual bool is_child_on_side(const unsigned int c,
00112                                 const unsigned int s) const;
00113 
00114   /*
00115    * @returns true iff the element map is definitely affine within
00116    * numerical tolerances
00117    */
00118   virtual bool has_affine_map () const { return true; }
00119 
00124   virtual bool is_linear () const { return true; }
00125 
00129   Order default_order() const { return FIRST; }
00130 
00135   UniquePtr<Elem> build_side (const unsigned int i,
00136                               bool proxy) const;
00137 
00142   UniquePtr<Elem> build_edge (const unsigned int i) const;
00143 
00144   virtual void connectivity(const unsigned int sc,
00145                             const IOPackage iop,
00146                             std::vector<dof_id_type>& conn) const;
00147 
00152   static const unsigned int side_nodes_map[4][3];
00153 
00158   static const unsigned int edge_nodes_map[6][2];
00159 
00164   virtual Real volume () const;
00165 
00174   std::pair<Real, Real> min_and_max_angle() const;
00175 
00176 protected:
00177 
00181   Node* _nodelinks_data[4];
00182 
00183 
00184 
00185 #ifdef LIBMESH_ENABLE_AMR
00186 
00190   float embedding_matrix (const unsigned int i,
00191                           const unsigned int j,
00192                           const unsigned int k) const;
00193 
00198   static const float _embedding_matrix[8][4][4];
00199 
00200   // public:
00201   //
00202   //  /**
00203   //   * Allows the user to reselect the diagonal after refinement.  This
00204   //   * function may only be called directly after the element is refined
00205   //   * for the first time (and before the \p EquationSystems::reinit()
00206   //   * is called).  It will destroy and re-create the children if
00207   //   * necessary.
00208   //   */
00209   //  void reselect_diagonal (const Diagonal diag);
00210   //
00211   //  /**
00212   //   * Reselects the diagonal after refinement to be the optimal one.
00213   //   * This makes sense if the user has moved some grid points, so that
00214   //   * the former optimal choice is no longer optimal.  Also, the user
00215   //   * may exclude one diagonal from this selection by giving it as
00216   //   * argument.  In this case, the more optimal one of the remaining
00217   //   * two diagonals is chosen.
00218   //   */
00219   //  void reselect_optimal_diagonal (const Diagonal exclude_this=INVALID_DIAG);
00220 
00221 #endif
00222 
00223 };
00224 
00225 
00226 
00227 // ------------------------------------------------------------
00228 // Tet4 class member functions
00229 inline
00230 Tet4::Tet4(Elem* p) :
00231   Tet(Tet4::n_nodes(), p, _nodelinks_data)
00232 {
00233 }
00234 
00235 } // namespace libMesh
00236 
00237 
00238 #endif // LIBMESH_CELL_TET4_H