cell_tet4.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_CELL_TET4_H
21 #define LIBMESH_CELL_TET4_H
22 
23 // Local includes
24 #include "libmesh/cell_tet.h"
25 
26 // C++ includes
27 #include <cstddef>
28 
29 namespace libMesh
30 {
31 
32 
33 
34 
54 // ------------------------------------------------------------
55 // Tet4 class definition
56 class Tet4 : public Tet
57 {
58 public:
59 
63  explicit
64  Tet4 (Elem* p=NULL);
65 
69  ElemType type () const { return TET4; }
70 
74  unsigned int n_nodes() const { return 4; }
75 
79  unsigned int n_sub_elem() const { return 1; }
80 
84  virtual bool is_vertex(const unsigned int i) const;
85 
89  virtual bool is_edge(const unsigned int i) const;
90 
94  virtual bool is_face(const unsigned int i) const;
95 
96  /*
97  * @returns true iff the specified (local) node number is on the
98  * specified side
99  */
100  virtual bool is_node_on_side(const unsigned int n,
101  const unsigned int s) const;
102 
103  /*
104  * @returns true iff the specified (local) node number is on the
105  * specified edge
106  */
107  virtual bool is_node_on_edge(const unsigned int n,
108  const unsigned int e) const;
109 
110  /*
111  * @returns true iff the specified child is on the
112  * specified side
113  */
114  virtual bool is_child_on_side(const unsigned int c,
115  const unsigned int s) const;
116 
117  /*
118  * @returns true iff the element map is definitely affine within
119  * numerical tolerances
120  */
121  virtual bool has_affine_map () const { return true; }
122 
127  virtual bool is_linear () const { return true; }
128 
132  Order default_order() const { return FIRST; }
133 
138  AutoPtr<Elem> build_side (const unsigned int i,
139  bool proxy) const;
140 
145  AutoPtr<Elem> build_edge (const unsigned int i) const;
146 
147  virtual void connectivity(const unsigned int sc,
148  const IOPackage iop,
149  std::vector<dof_id_type>& conn) const;
150 
155  static const unsigned int side_nodes_map[4][3];
156 
161  static const unsigned int edge_nodes_map[6][2];
162 
167  virtual Real volume () const;
168 
177  std::pair<Real, Real> min_and_max_angle() const;
178 
179 protected:
180 
185 
186 
187 
188 #ifdef LIBMESH_ENABLE_AMR
189 
193  float embedding_matrix (const unsigned int i,
194  const unsigned int j,
195  const unsigned int k) const;
196 
201  static const float _embedding_matrix[8][4][4];
202 
203 // public:
204 //
205 // /**
206 // * Allows the user to reselect the diagonal after refinement. This
207 // * function may only be called directly after the element is refined
208 // * for the first time (and before the \p EquationSystems::reinit()
209 // * is called). It will destroy and re-create the children if
210 // * necessary.
211 // */
212 // void reselect_diagonal (const Diagonal diag);
213 //
214 // /**
215 // * Reselects the diagonal after refinement to be the optimal one.
216 // * This makes sense if the user has moved some grid points, so that
217 // * the former optimal choice is no longer optimal. Also, the user
218 // * may exclude one diagonal from this selection by giving it as
219 // * argument. In this case, the more optimal one of the remaining
220 // * two diagonals is chosen.
221 // */
222 // void reselect_optimal_diagonal (const Diagonal exclude_this=INVALID_DIAG);
223 
224 #endif
225 
226 };
227 
228 
229 
230 // ------------------------------------------------------------
231 // Tet4 class member functions
232 inline
234  Tet(Tet4::n_nodes(), p, _nodelinks_data)
235 {
236 }
237 
238 } // namespace libMesh
239 
240 
241 #endif // LIBMESH_CELL_TET4_H

Site Created By: libMesh Developers
Last modified: February 07 2014 16:57:04 UTC

Hosted By:
SourceForge.net Logo