$extrastylesheet
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/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 UniquePtr<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 return UniquePtr<Elem>(new Side<Edge2,Tri3>(this,i)); 00116 00117 else 00118 { 00119 Elem* edge = new Edge2; 00120 edge->subdomain_id() = this->subdomain_id(); 00121 00122 switch (i) 00123 { 00124 case 0: 00125 { 00126 edge->set_node(0) = this->get_node(0); 00127 edge->set_node(1) = this->get_node(1); 00128 break; 00129 } 00130 case 1: 00131 { 00132 edge->set_node(0) = this->get_node(1); 00133 edge->set_node(1) = this->get_node(2); 00134 break; 00135 } 00136 case 2: 00137 { 00138 edge->set_node(0) = this->get_node(2); 00139 edge->set_node(1) = this->get_node(0); 00140 break; 00141 } 00142 default: 00143 libmesh_error_msg("Invalid side i = " << i); 00144 } 00145 00146 return UniquePtr<Elem>(edge); 00147 } 00148 00149 libmesh_error_msg("We'll never get here!"); 00150 return UniquePtr<Elem>(); 00151 } 00152 00153 00154 void Tri3::connectivity(const unsigned int libmesh_dbg_var(sf), 00155 const IOPackage iop, 00156 std::vector<dof_id_type>& conn) const 00157 { 00158 libmesh_assert_less (sf, this->n_sub_elem()); 00159 libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE); 00160 00161 switch (iop) 00162 { 00163 case TECPLOT: 00164 { 00165 conn.resize(4); 00166 conn[0] = this->node(0)+1; 00167 conn[1] = this->node(1)+1; 00168 conn[2] = this->node(2)+1; 00169 conn[3] = this->node(2)+1; 00170 return; 00171 } 00172 00173 case VTK: 00174 { 00175 conn.resize(3); 00176 conn[0] = this->node(0); 00177 conn[1] = this->node(1); 00178 conn[2] = this->node(2); 00179 return; 00180 } 00181 00182 default: 00183 libmesh_error_msg("Unsupported IO package " << iop); 00184 } 00185 } 00186 00187 00188 00189 00190 00191 00192 Real Tri3::volume () const 00193 { 00194 // 3-node triangles have the following formula for computing the area 00195 Point v10 ( *(this->get_node(1)) - *(this->get_node(0)) ); 00196 00197 Point v20 ( *(this->get_node(2)) - *(this->get_node(0)) ); 00198 00199 return 0.5 * (v10.cross(v20)).size() ; 00200 } 00201 00202 00203 00204 std::pair<Real, Real> Tri3::min_and_max_angle() const 00205 { 00206 Point v10 ( this->point(1) - this->point(0) ); 00207 Point v20 ( this->point(2) - this->point(0) ); 00208 Point v21 ( this->point(2) - this->point(1) ); 00209 00210 const Real 00211 len_10=v10.size(), 00212 len_20=v20.size(), 00213 len_21=v21.size() 00214 ; 00215 00216 const Real 00217 theta0=std::acos(( v10*v20)/len_10/len_20), 00218 theta1=std::acos((-v10*v21)/len_10/len_21), 00219 theta2=libMesh::pi - theta0 - theta1 00220 ; 00221 00222 libmesh_assert_greater (theta0, 0.); 00223 libmesh_assert_greater (theta1, 0.); 00224 libmesh_assert_greater (theta2, 0.); 00225 00226 return std::make_pair(std::min(theta0, std::min(theta1,theta2)), 00227 std::max(theta0, std::max(theta1,theta2))); 00228 } 00229 00230 } // namespace libMesh