face_tri3.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/side.h" 00022 #include "libmesh/edge_edge2.h" 00023 #include "libmesh/face_tri3.h" 00024 00025 namespace libMesh 00026 { 00027 00028 00029 00030 // ------------------------------------------------------------ 00031 // Tri3 class static member initializations 00032 const unsigned int Tri3::side_nodes_map[3][2] = 00033 { 00034 {0, 1}, // Side 0 00035 {1, 2}, // Side 1 00036 {2, 0} // Side 2 00037 }; 00038 00039 00040 #ifdef LIBMESH_ENABLE_AMR 00041 00042 const float Tri3::_embedding_matrix[4][3][3] = 00043 { 00044 // embedding matrix for child 0 00045 { 00046 // 0 1 2 00047 {1.0, 0.0, 0.0}, // 0 00048 {0.5, 0.5, 0.0}, // 1 00049 {0.5, 0.0, 0.5} // 2 00050 }, 00051 00052 // embedding matrix for child 1 00053 { 00054 // 0 1 2 00055 {0.5, 0.5, 0.0}, // 0 00056 {0.0, 1.0, 0.0}, // 1 00057 {0.0, 0.5, 0.5} // 2 00058 }, 00059 00060 // embedding matrix for child 2 00061 { 00062 // 0 1 2 00063 {0.5, 0.0, 0.5}, // 0 00064 {0.0, 0.5, 0.5}, // 1 00065 {0.0, 0.0, 1.0} // 2 00066 }, 00067 00068 // embedding matrix for child 3 00069 { 00070 // 0 1 2 00071 {0.5, 0.5, 0.0}, // 0 00072 {0.0, 0.5, 0.5}, // 1 00073 {0.5, 0.0, 0.5} // 2 00074 } 00075 }; 00076 00077 #endif 00078 00079 00080 00081 // ------------------------------------------------------------ 00082 // Tri3 class member functions 00083 00084 bool Tri3::is_vertex(const unsigned int) const 00085 { 00086 return true; 00087 } 00088 00089 bool Tri3::is_edge(const unsigned int) const 00090 { 00091 return false; 00092 } 00093 00094 bool Tri3::is_face(const unsigned int) const 00095 { 00096 return false; 00097 } 00098 00099 bool Tri3::is_node_on_side(const unsigned int n, 00100 const unsigned int s) const 00101 { 00102 libmesh_assert_less (s, n_sides()); 00103 for (unsigned int i = 0; i != 2; ++i) 00104 if (side_nodes_map[s][i] == n) 00105 return true; 00106 return false; 00107 } 00108 00109 AutoPtr<Elem> Tri3::build_side (const unsigned int i, 00110 bool proxy) const 00111 { 00112 libmesh_assert_less (i, this->n_sides()); 00113 00114 if (proxy) 00115 { 00116 AutoPtr<Elem> ap(new Side<Edge2,Tri3>(this,i)); 00117 return ap; 00118 } 00119 00120 else 00121 { 00122 Edge2* edge = new Edge2; 00123 00124 switch (i) 00125 { 00126 case 0: 00127 { 00128 edge->set_node(0) = this->get_node(0); 00129 edge->set_node(1) = this->get_node(1); 00130 00131 AutoPtr<Elem> ap(edge); return ap; 00132 } 00133 case 1: 00134 { 00135 edge->set_node(0) = this->get_node(1); 00136 edge->set_node(1) = this->get_node(2); 00137 00138 AutoPtr<Elem> ap(edge); return ap; 00139 } 00140 case 2: 00141 { 00142 edge->set_node(0) = this->get_node(2); 00143 edge->set_node(1) = this->get_node(0); 00144 00145 AutoPtr<Elem> ap(edge); return ap; 00146 } 00147 default: 00148 { 00149 libmesh_error(); 00150 } 00151 } 00152 } 00153 00154 // We will never get here... Look at the code above. 00155 libmesh_error(); 00156 AutoPtr<Elem> ap(NULL); return ap; 00157 } 00158 00159 00160 void Tri3::connectivity(const unsigned int libmesh_dbg_var(sf), 00161 const IOPackage iop, 00162 std::vector<dof_id_type>& conn) const 00163 { 00164 libmesh_assert_less (sf, this->n_sub_elem()); 00165 libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE); 00166 00167 switch (iop) 00168 { 00169 case TECPLOT: 00170 { 00171 conn.resize(4); 00172 conn[0] = this->node(0)+1; 00173 conn[1] = this->node(1)+1; 00174 conn[2] = this->node(2)+1; 00175 conn[3] = this->node(2)+1; 00176 return; 00177 } 00178 00179 case VTK: 00180 { 00181 conn.resize(3); 00182 conn[0] = this->node(0); 00183 conn[1] = this->node(1); 00184 conn[2] = this->node(2); 00185 return; 00186 } 00187 00188 default: 00189 libmesh_error(); 00190 } 00191 00192 libmesh_error(); 00193 } 00194 00195 00196 00197 00198 00199 00200 Real Tri3::volume () const 00201 { 00202 // 3-node triangles have the following formula for computing the area 00203 Point v10 ( *(this->get_node(1)) - *(this->get_node(0)) ); 00204 00205 Point v20 ( *(this->get_node(2)) - *(this->get_node(0)) ); 00206 00207 return 0.5 * (v10.cross(v20)).size() ; 00208 } 00209 00210 00211 00212 std::pair<Real, Real> Tri3::min_and_max_angle() const 00213 { 00214 Point v10 ( this->point(1) - this->point(0) ); 00215 Point v20 ( this->point(2) - this->point(0) ); 00216 Point v21 ( this->point(2) - this->point(1) ); 00217 00218 const Real 00219 len_10=v10.size(), 00220 len_20=v20.size(), 00221 len_21=v21.size() 00222 ; 00223 00224 const Real 00225 theta0=std::acos(( v10*v20)/len_10/len_20), 00226 theta1=std::acos((-v10*v21)/len_10/len_21), 00227 theta2=libMesh::pi - theta0 - theta1 00228 ; 00229 00230 libmesh_assert_greater (theta0, 0.); 00231 libmesh_assert_greater (theta1, 0.); 00232 libmesh_assert_greater (theta2, 0.); 00233 00234 return std::make_pair(std::min(theta0, std::min(theta1,theta2)), 00235 std::max(theta0, std::max(theta1,theta2))); 00236 } 00237 00238 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: