face_tri.C
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 // 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 00054 00055 // We will never get here... Look at the code above. 00056 libmesh_error(); 00057 return 0; 00058 } 00059 00060 00061 00062 AutoPtr<Elem> Tri::side (const unsigned int i) const 00063 { 00064 libmesh_assert_less (i, this->n_sides()); 00065 00066 Elem* edge = new Edge2; 00067 00068 switch (i) 00069 { 00070 case 0: 00071 { 00072 edge->set_node(0) = this->get_node(0); 00073 edge->set_node(1) = this->get_node(1); 00074 00075 AutoPtr<Elem> ap_edge(edge); 00076 return ap_edge; 00077 } 00078 case 1: 00079 { 00080 edge->set_node(0) = this->get_node(1); 00081 edge->set_node(1) = this->get_node(2); 00082 00083 AutoPtr<Elem> ap_edge(edge); 00084 return ap_edge; 00085 } 00086 case 2: 00087 { 00088 edge->set_node(0) = this->get_node(2); 00089 edge->set_node(1) = this->get_node(0); 00090 00091 AutoPtr<Elem> ap_edge(edge); 00092 return ap_edge; 00093 } 00094 default: 00095 { 00096 libmesh_error(); 00097 } 00098 } 00099 00100 00101 // We will never get here... Look at the code above. 00102 libmesh_error(); 00103 AutoPtr<Elem> ap_edge(edge); 00104 return ap_edge; 00105 } 00106 00107 00108 00109 bool Tri::is_child_on_side(const unsigned int c, 00110 const unsigned int s) const 00111 { 00112 libmesh_assert_less (c, this->n_children()); 00113 libmesh_assert_less (s, this->n_sides()); 00114 00115 return (c == s || c == (s+1)%3); 00116 } 00117 00118 00119 00120 Real Tri::quality (const ElemQuality q) const 00121 { 00122 switch (q) 00123 { 00124 00128 case DISTORTION: 00129 case STRETCH: 00130 { 00131 const Node* p1 = this->get_node(0); 00132 const Node* p2 = this->get_node(1); 00133 const Node* p3 = this->get_node(2); 00134 00135 Point v1 = (*p2) - (*p1); 00136 Point v2 = (*p3) - (*p1); 00137 Point v3 = (*p3) - (*p2); 00138 const Real l1 = v1.size(); 00139 const Real l2 = v2.size(); 00140 const Real l3 = v3.size(); 00141 00142 // if one length is 0, quality is quite bad! 00143 if ((l1 <=0.) || (l2 <= 0.) || (l3 <= 0.)) 00144 return 0.; 00145 00146 const Real s1 = std::sin(std::acos(v1*v2/l1/l2)/2.); 00147 v1 *= -1; 00148 const Real s2 = std::sin(std::acos(v1*v3/l1/l3)/2.); 00149 const Real s3 = std::sin(std::acos(v2*v3/l2/l3)/2.); 00150 00151 return 8. * s1 * s2 * s3; 00152 00153 } 00154 default: 00155 return Elem::quality(q); 00156 } 00157 00163 return Elem::quality(q); 00164 00165 } 00166 00167 00168 00169 00170 00171 00172 std::pair<Real, Real> Tri::qual_bounds (const ElemQuality q) const 00173 { 00174 std::pair<Real, Real> bounds; 00175 00176 switch (q) 00177 { 00178 00179 case MAX_ANGLE: 00180 bounds.first = 60.; 00181 bounds.second = 90.; 00182 break; 00183 00184 case MIN_ANGLE: 00185 bounds.first = 30.; 00186 bounds.second = 60.; 00187 break; 00188 00189 case CONDITION: 00190 bounds.first = 1.; 00191 bounds.second = 1.3; 00192 break; 00193 00194 case JACOBIAN: 00195 bounds.first = 0.5; 00196 bounds.second = 1.155; 00197 break; 00198 00199 case SIZE: 00200 case SHAPE: 00201 bounds.first = 0.25; 00202 bounds.second = 1.; 00203 break; 00204 00205 case DISTORTION: 00206 bounds.first = 0.6; 00207 bounds.second = 1.; 00208 break; 00209 00210 default: 00211 libMesh::out << "Warning: Invalid quality measure chosen." << std::endl; 00212 bounds.first = -1; 00213 bounds.second = -1; 00214 } 00215 00216 return bounds; 00217 } 00218 00219 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: