$extrastylesheet
cell_tet.C
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 // C++ includes
00020 
00021 // Local includes
00022 #include "libmesh/cell_tet.h"
00023 #include "libmesh/cell_tet4.h"
00024 #include "libmesh/face_tri3.h"
00025 
00026 namespace libMesh
00027 {
00028 
00029 
00030 
00031 // ------------------------------------------------------------
00032 // Tet class member functions
00033 dof_id_type Tet::key (const unsigned int s) const
00034 {
00035   libmesh_assert_less (s, this->n_sides());
00036 
00037   switch (s)
00038     {
00039     case 0:
00040       return
00041         this->compute_key (this->node(0),
00042                            this->node(2),
00043                            this->node(1));
00044 
00045     case 1:
00046       return
00047         this->compute_key (this->node(0),
00048                            this->node(1),
00049                            this->node(3));
00050 
00051     case 2:
00052       return
00053         this->compute_key (this->node(1),
00054                            this->node(2),
00055                            this->node(3));
00056 
00057     case 3:
00058       return
00059         this->compute_key (this->node(2),
00060                            this->node(0),
00061                            this->node(3));
00062 
00063     default:
00064       libmesh_error_msg("Invalid side s = " << s);
00065     }
00066 
00067   libmesh_error_msg("We'll never get here!");
00068   return 0;
00069 }
00070 
00071 
00072 
00073 UniquePtr<Elem> Tet::side (const unsigned int i) const
00074 {
00075   libmesh_assert_less (i, this->n_sides());
00076 
00077   Elem* face = new Tri3;
00078 
00079   switch (i)
00080     {
00081     case 0:
00082       {
00083         face->set_node(0) = this->get_node(0);
00084         face->set_node(1) = this->get_node(2);
00085         face->set_node(2) = this->get_node(1);
00086         break;
00087       }
00088     case 1:
00089       {
00090         face->set_node(0) = this->get_node(0);
00091         face->set_node(1) = this->get_node(1);
00092         face->set_node(2) = this->get_node(3);
00093         break;
00094       }
00095     case 2:
00096       {
00097         face->set_node(0) = this->get_node(1);
00098         face->set_node(1) = this->get_node(2);
00099         face->set_node(2) = this->get_node(3);
00100         break;
00101       }
00102     case 3:
00103       {
00104         face->set_node(0) = this->get_node(2);
00105         face->set_node(1) = this->get_node(0);
00106         face->set_node(2) = this->get_node(3);
00107         break;
00108       }
00109     default:
00110       libmesh_error_msg("Invalid side i = " << i);
00111     }
00112 
00113   return UniquePtr<Elem>(face);
00114 }
00115 
00116 
00117 void Tet::select_diagonal (const Diagonal diag) const
00118 {
00119   libmesh_assert_equal_to (_diagonal_selection, INVALID_DIAG);
00120   _diagonal_selection = diag;
00121 }
00122 
00123 
00124 
00125 
00126 
00127 #ifdef LIBMESH_ENABLE_AMR
00128 
00129 
00130 bool Tet::is_child_on_side_helper(const unsigned int c,
00131                                   const unsigned int s,
00132                                   const unsigned int checked_nodes[][3]) const
00133 {
00134   libmesh_assert_less (c, this->n_children());
00135   libmesh_assert_less (s, this->n_sides());
00136 
00137   // For the 4 vertices, child c touches vertex c, so we can return
00138   // true if that vertex is on side s
00139   for (unsigned int i = 0; i != 3; ++i)
00140     if (Tet4::side_nodes_map[s][i] == c)
00141       return true;
00142 
00143   // If we are a "vertex child" and we didn't already return true,
00144   // we must not be on the side in question
00145   if (c < 4)
00146     return false;
00147 
00148   // For the 4 non-vertex children, the child ordering depends on the
00149   // diagonal selection.  We'll let the embedding matrix figure that
00150   // out: if this child has three nodes that don't depend on the
00151   // position of the node_facing_side[s], then we're on side s.  Which
00152   // three nodes those are depends on the subclass, so their responsibility
00153   // is to call this function with the proper check_nodes array
00154   const unsigned int node_facing_side[4] = {3, 2, 0, 1};
00155   const unsigned int n = node_facing_side[s];
00156 
00157   // Add up the absolute values of the entries of the embedding matrix for the
00158   // nodes opposite node n.  If it is equal to zero, then the child in question is
00159   // on side s, so return true.
00160   Real embedding_sum = 0.;
00161   for (unsigned i=0; i<3; ++i)
00162     embedding_sum += std::abs(this->embedding_matrix(c, checked_nodes[n][i], n));
00163 
00164   return ( std::abs(embedding_sum) < 1.e-3 );
00165 }
00166 
00167 #else
00168 
00169 bool Tet::is_child_on_side_helper(const unsigned int /*c*/,
00170                                   const unsigned int /*s*/,
00171                                   const unsigned int /*checked_nodes*/[][3]) const
00172 {
00173   libmesh_not_implemented();
00174   return false;
00175 }
00176 
00177 #endif //LIBMESH_ENABLE_AMR
00178 
00179 
00180 
00181 
00182 void Tet::choose_diagonal() const
00183 {
00184   // Check for uninitialized diagonal selection
00185   if (this->_diagonal_selection==INVALID_DIAG)
00186     {
00187       Real diag_01_23 = (this->point(0)+this->point(1)-this->point(2)-this->point(3)).size_sq();
00188       Real diag_02_13 = (this->point(0)-this->point(1)+this->point(2)-this->point(3)).size_sq();
00189       Real diag_03_12 = (this->point(0)-this->point(1)-this->point(2)+this->point(3)).size_sq();
00190 
00191       this->_diagonal_selection=DIAG_02_13;
00192 
00193       if (diag_01_23 < diag_02_13 || diag_03_12 < diag_02_13)
00194         {
00195           if (diag_01_23 < diag_03_12)
00196             this->_diagonal_selection=DIAG_01_23;
00197 
00198           else
00199             this->_diagonal_selection=DIAG_03_12;
00200         }
00201     }
00202 }
00203 
00204 
00205 
00206 bool Tet::is_edge_on_side(const unsigned int e,
00207                           const unsigned int s) const
00208 {
00209   libmesh_assert_less (e, this->n_edges());
00210   libmesh_assert_less (s, this->n_sides());
00211 
00212   return (is_node_on_side(Tet4::edge_nodes_map[e][0],s) &&
00213           is_node_on_side(Tet4::edge_nodes_map[e][1],s));
00214 }
00215 
00216 
00217 
00218 Real Tet::quality(const ElemQuality q) const
00219 {
00220   return Elem::quality(q); // Not implemented
00221 }
00222 
00223 
00224 
00225 
00226 std::pair<Real, Real> Tet::qual_bounds (const ElemQuality q) const
00227 {
00228   std::pair<Real, Real> bounds;
00229 
00230   switch (q)
00231     {
00232 
00233     case ASPECT_RATIO_BETA:
00234     case ASPECT_RATIO_GAMMA:
00235       bounds.first  = 1.;
00236       bounds.second = 3.;
00237       break;
00238 
00239     case SIZE:
00240     case SHAPE:
00241       bounds.first  = 0.2;
00242       bounds.second = 1.;
00243       break;
00244 
00245     case CONDITION:
00246       bounds.first  = 1.;
00247       bounds.second = 3.;
00248       break;
00249 
00250     case DISTORTION:
00251       bounds.first  = 0.6;
00252       bounds.second = 1.;
00253       break;
00254 
00255     case JACOBIAN:
00256       bounds.first  = 0.5;
00257       bounds.second = 1.414;
00258       break;
00259 
00260     default:
00261       libMesh::out << "Warning: Invalid quality measure chosen." << std::endl;
00262       bounds.first  = -1;
00263       bounds.second = -1;
00264     }
00265 
00266   return bounds;
00267 }
00268 
00269 } // namespace libMesh