$extrastylesheet
face_tri.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 // C++ includes
00019 
00020 // Local includes
00021 #include "libmesh/face_tri.h"
00022 #include "libmesh/edge_edge2.h"
00023 
00024 namespace libMesh
00025 {
00026 
00027 
00028 
00029 
00030 
00031 // ------------------------------------------------------------
00032 // Tri class member functions
00033 dof_id_type Tri::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(1));
00043 
00044     case 1:
00045       return
00046         this->compute_key (this->node(1),
00047                            this->node(2));
00048     case 2:
00049       return
00050         this->compute_key (this->node(2),
00051                            this->node(0));
00052 
00053     default:
00054       libmesh_error_msg("Invalid side s = " << s);
00055     }
00056 
00057   libmesh_error_msg("We'll never get here!");
00058   return 0;
00059 }
00060 
00061 
00062 
00063 UniquePtr<Elem> Tri::side (const unsigned int i) const
00064 {
00065   libmesh_assert_less (i, this->n_sides());
00066 
00067   Elem* edge = new Edge2;
00068 
00069   switch (i)
00070     {
00071     case 0:
00072       {
00073         edge->set_node(0) = this->get_node(0);
00074         edge->set_node(1) = this->get_node(1);
00075         break;
00076       }
00077     case 1:
00078       {
00079         edge->set_node(0) = this->get_node(1);
00080         edge->set_node(1) = this->get_node(2);
00081         break;
00082       }
00083     case 2:
00084       {
00085         edge->set_node(0) = this->get_node(2);
00086         edge->set_node(1) = this->get_node(0);
00087         break;
00088       }
00089     default:
00090       libmesh_error_msg("Invalid side i = " << i);
00091     }
00092 
00093   return UniquePtr<Elem>(edge);
00094 }
00095 
00096 
00097 
00098 bool Tri::is_child_on_side(const unsigned int c,
00099                            const unsigned int s) const
00100 {
00101   libmesh_assert_less (c, this->n_children());
00102   libmesh_assert_less (s, this->n_sides());
00103 
00104   return (c == s || c == (s+1)%3);
00105 }
00106 
00107 
00108 
00109 Real Tri::quality (const ElemQuality q) const
00110 {
00111   switch (q)
00112     {
00113 
00117     case DISTORTION:
00118     case STRETCH:
00119       {
00120         const Node* p1 = this->get_node(0);
00121         const Node* p2 = this->get_node(1);
00122         const Node* p3 = this->get_node(2);
00123 
00124         Point v1 = (*p2) - (*p1);
00125         Point v2 = (*p3) - (*p1);
00126         Point v3 = (*p3) - (*p2);
00127         const Real l1 = v1.size();
00128         const Real l2 = v2.size();
00129         const Real l3 = v3.size();
00130 
00131         // if one length is 0, quality is quite bad!
00132         if ((l1 <=0.) || (l2 <= 0.) || (l3 <= 0.))
00133           return 0.;
00134 
00135         const Real s1 = std::sin(std::acos(v1*v2/l1/l2)/2.);
00136         v1 *= -1;
00137         const Real s2 = std::sin(std::acos(v1*v3/l1/l3)/2.);
00138         const Real s3 = std::sin(std::acos(v2*v3/l2/l3)/2.);
00139 
00140         return 8. * s1 * s2 * s3;
00141 
00142       }
00143     default:
00144       return Elem::quality(q);
00145     }
00146 
00152   return Elem::quality(q);
00153 
00154 }
00155 
00156 
00157 
00158 
00159 
00160 
00161 std::pair<Real, Real> Tri::qual_bounds (const ElemQuality q) const
00162 {
00163   std::pair<Real, Real> bounds;
00164 
00165   switch (q)
00166     {
00167 
00168     case MAX_ANGLE:
00169       bounds.first  = 60.;
00170       bounds.second = 90.;
00171       break;
00172 
00173     case MIN_ANGLE:
00174       bounds.first  = 30.;
00175       bounds.second = 60.;
00176       break;
00177 
00178     case CONDITION:
00179       bounds.first  = 1.;
00180       bounds.second = 1.3;
00181       break;
00182 
00183     case JACOBIAN:
00184       bounds.first  = 0.5;
00185       bounds.second = 1.155;
00186       break;
00187 
00188     case SIZE:
00189     case SHAPE:
00190       bounds.first  = 0.25;
00191       bounds.second = 1.;
00192       break;
00193 
00194     case DISTORTION:
00195       bounds.first  = 0.6;
00196       bounds.second = 1.;
00197       break;
00198 
00199     default:
00200       libMesh::out << "Warning: Invalid quality measure chosen." << std::endl;
00201       bounds.first  = -1;
00202       bounds.second = -1;
00203     }
00204 
00205   return bounds;
00206 }
00207 
00208 } // namespace libMesh