face_quad.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_quad.h" 00022 #include "libmesh/edge_edge2.h" 00023 00024 namespace libMesh 00025 { 00026 00027 00028 // ------------------------------------------------------------ 00029 // Quad class member functions 00030 dof_id_type Quad::key (const unsigned int s) const 00031 { 00032 libmesh_assert_less (s, this->n_sides()); 00033 00034 switch (s) 00035 { 00036 case 0: 00037 return 00038 this->compute_key (this->node(0), 00039 this->node(1)); 00040 00041 case 1: 00042 return 00043 this->compute_key (this->node(1), 00044 this->node(2)); 00045 00046 case 2: 00047 return 00048 this->compute_key (this->node(2), 00049 this->node(3)); 00050 00051 case 3: 00052 return 00053 this->compute_key (this->node(3), 00054 this->node(0)); 00055 } 00056 00057 00058 // We will never get here... Look at the code above. 00059 libmesh_error(); 00060 return 0; 00061 } 00062 00063 00064 00065 AutoPtr<Elem> Quad::side (const unsigned int i) const 00066 { 00067 libmesh_assert_less (i, this->n_sides()); 00068 00069 Elem* edge = new Edge2; 00070 00071 switch (i) 00072 { 00073 case 0: 00074 { 00075 edge->set_node(0) = this->get_node(0); 00076 edge->set_node(1) = this->get_node(1); 00077 00078 AutoPtr<Elem> ap_edge(edge); 00079 return ap_edge; 00080 } 00081 case 1: 00082 { 00083 edge->set_node(0) = this->get_node(1); 00084 edge->set_node(1) = this->get_node(2); 00085 00086 AutoPtr<Elem> ap_edge(edge); 00087 return ap_edge; 00088 } 00089 case 2: 00090 { 00091 edge->set_node(0) = this->get_node(2); 00092 edge->set_node(1) = this->get_node(3); 00093 00094 AutoPtr<Elem> ap_edge(edge); 00095 return ap_edge; 00096 } 00097 case 3: 00098 { 00099 edge->set_node(0) = this->get_node(3); 00100 edge->set_node(1) = this->get_node(0); 00101 00102 AutoPtr<Elem> ap_edge(edge); 00103 return ap_edge; 00104 } 00105 default: 00106 { 00107 libmesh_error(); 00108 } 00109 } 00110 00111 00112 // We will never get here... Look at the code above. 00113 libmesh_error(); 00114 AutoPtr<Elem> ap_edge(edge); 00115 return ap_edge; 00116 } 00117 00118 00119 00120 bool Quad::is_child_on_side(const unsigned int c, 00121 const unsigned int s) const 00122 { 00123 libmesh_assert_less (c, this->n_children()); 00124 libmesh_assert_less (s, this->n_sides()); 00125 00126 // A quad's children and nodes don't share the same ordering: 00127 // child 2 and 3 are swapped; 00128 unsigned int n = (c < 2) ? c : 5-c; 00129 return (n == s || n == (s+1)%4); 00130 } 00131 00132 00133 00134 unsigned int Quad::opposite_side(const unsigned int side_in) const 00135 { 00136 libmesh_assert_less (side_in, 4); 00137 00138 return (side_in + 2) % 4; 00139 } 00140 00141 00142 00143 unsigned int Quad::opposite_node(const unsigned int node_in, 00144 const unsigned int side_in) const 00145 { 00146 libmesh_assert_less (node_in, 8); 00147 libmesh_assert_less (node_in, this->n_nodes()); 00148 libmesh_assert_less (side_in, this->n_sides()); 00149 libmesh_assert(this->is_node_on_side(node_in, side_in)); 00150 00151 //unsigned int opposite; 00152 00153 static const unsigned char side02_nodes_map[] = 00154 {3, 2, 1, 0, 6, 255, 4, 255}; 00155 static const unsigned char side13_nodes_map[] = 00156 {1, 0, 3, 2, 255, 7, 255, 5}; 00157 00158 switch (side_in) 00159 { 00160 case 0: 00161 case 2: 00162 return side02_nodes_map[node_in]; 00163 break; 00164 case 1: 00165 case 3: 00166 return side13_nodes_map[node_in]; 00167 break; 00168 } 00169 00170 libmesh_error(); 00171 return 255; 00172 } 00173 00174 00175 Real Quad::quality (const ElemQuality q) const 00176 { 00177 switch (q) 00178 { 00179 00184 case DISTORTION: 00185 case DIAGONAL: 00186 case STRETCH: 00187 { 00188 // Diagonal between node 0 and node 2 00189 const Real d02 = this->length(0,2); 00190 00191 // Diagonal between node 1 and node 3 00192 const Real d13 = this->length(1,3); 00193 00194 // Find the biggest and smallest diagonals 00195 if ( (d02 > 0.) && (d13 >0.) ) 00196 if (d02 < d13) return d02 / d13; 00197 else return d13 / d02; 00198 else 00199 return 0.; 00200 break; 00201 } 00202 00203 default: 00204 return Elem::quality(q); 00205 } 00206 00212 return Elem::quality(q); 00213 } 00214 00215 00216 00217 00218 00219 00220 std::pair<Real, Real> Quad::qual_bounds (const ElemQuality q) const 00221 { 00222 std::pair<Real, Real> bounds; 00223 00224 switch (q) 00225 { 00226 00227 case ASPECT_RATIO: 00228 bounds.first = 1.; 00229 bounds.second = 4.; 00230 break; 00231 00232 case SKEW: 00233 bounds.first = 0.; 00234 bounds.second = 0.5; 00235 break; 00236 00237 case TAPER: 00238 bounds.first = 0.; 00239 bounds.second = 0.7; 00240 break; 00241 00242 case WARP: 00243 bounds.first = 0.9; 00244 bounds.second = 1.; 00245 break; 00246 00247 case STRETCH: 00248 bounds.first = 0.25; 00249 bounds.second = 1.; 00250 break; 00251 00252 case MIN_ANGLE: 00253 bounds.first = 45.; 00254 bounds.second = 90.; 00255 break; 00256 00257 case MAX_ANGLE: 00258 bounds.first = 90.; 00259 bounds.second = 135.; 00260 break; 00261 00262 case CONDITION: 00263 bounds.first = 1.; 00264 bounds.second = 4.; 00265 break; 00266 00267 case JACOBIAN: 00268 bounds.first = 0.5; 00269 bounds.second = 1.; 00270 break; 00271 00272 case SHEAR: 00273 case SHAPE: 00274 case SIZE: 00275 bounds.first = 0.3; 00276 bounds.second = 1.; 00277 break; 00278 00279 case DISTORTION: 00280 bounds.first = 0.6; 00281 bounds.second = 1.; 00282 break; 00283 00284 default: 00285 libMesh::out << "Warning: Invalid quality measure chosen." << std::endl; 00286 bounds.first = -1; 00287 bounds.second = -1; 00288 } 00289 00290 return bounds; 00291 } 00292 00293 00294 00295 00296 const unsigned short int Quad::_second_order_adjacent_vertices[4][2] = 00297 { 00298 {0, 1}, // vertices adjacent to node 4 00299 {1, 2}, // vertices adjacent to node 5 00300 {2, 3}, // vertices adjacent to node 6 00301 {0, 3} // vertices adjacent to node 7 00302 }; 00303 00304 00305 00306 const unsigned short int Quad::_second_order_vertex_child_number[9] = 00307 { 00308 99,99,99,99, // Vertices 00309 0,1,2,0, // Edges 00310 0 // Interior 00311 }; 00312 00313 00314 00315 const unsigned short int Quad::_second_order_vertex_child_index[9] = 00316 { 00317 99,99,99,99, // Vertices 00318 1,2,3,3, // Edges 00319 2 // Interior 00320 }; 00321 00322 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: