cell_prism6.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 00019 // C++ includes 00020 00021 // Local includes 00022 #include "libmesh/side.h" 00023 #include "libmesh/cell_prism6.h" 00024 #include "libmesh/edge_edge2.h" 00025 #include "libmesh/face_quad4.h" 00026 #include "libmesh/face_tri3.h" 00027 00028 namespace libMesh 00029 { 00030 00031 00032 00033 // ------------------------------------------------------------ 00034 // Prism6 class static member initializations 00035 const unsigned int Prism6::side_nodes_map[5][4] = 00036 { 00037 {0, 2, 1, 99}, // Side 0 00038 {0, 1, 4, 3}, // Side 1 00039 {1, 2, 5, 4}, // Side 2 00040 {2, 0, 3, 5}, // Side 3 00041 {3, 4, 5, 99} // Side 4 00042 }; 00043 00044 const unsigned int Prism6::side_elems_map[5][4] = 00045 { 00046 {0, 1, 2, 3}, // Side 0 00047 {0, 1, 4, 5}, // Side 1 00048 {1, 2, 5, 6}, // Side 2 00049 {0, 2, 4, 6}, // Side 3 00050 {4, 5, 6, 7} // Side 4 00051 }; 00052 00053 const unsigned int Prism6::edge_nodes_map[9][2] = 00054 { 00055 {0, 1}, // Side 0 00056 {1, 2}, // Side 1 00057 {0, 2}, // Side 2 00058 {0, 3}, // Side 3 00059 {1, 4}, // Side 4 00060 {2, 5}, // Side 5 00061 {3, 4}, // Side 6 00062 {4, 5}, // Side 7 00063 {3, 5} // Side 8 00064 }; 00065 00066 00067 // ------------------------------------------------------------ 00068 // Prism6 class member functions 00069 00070 bool Prism6::is_vertex(const unsigned int) const 00071 { 00072 return true; 00073 } 00074 00075 bool Prism6::is_edge(const unsigned int) const 00076 { 00077 return false; 00078 } 00079 00080 bool Prism6::is_face(const unsigned int) const 00081 { 00082 return false; 00083 } 00084 00085 bool Prism6::is_node_on_side(const unsigned int n, 00086 const unsigned int s) const 00087 { 00088 libmesh_assert_less (s, n_sides()); 00089 for (unsigned int i = 0; i != 4; ++i) 00090 if (side_nodes_map[s][i] == n) 00091 return true; 00092 return false; 00093 } 00094 00095 bool Prism6::is_node_on_edge(const unsigned int n, 00096 const unsigned int e) const 00097 { 00098 libmesh_assert_less (e, n_edges()); 00099 for (unsigned int i = 0; i != 2; ++i) 00100 if (edge_nodes_map[e][i] == n) 00101 return true; 00102 return false; 00103 } 00104 00105 00106 00107 bool Prism6::has_affine_map() const 00108 { 00109 // Make sure z edges are affine 00110 Point v = this->point(3) - this->point(0); 00111 if (!v.relative_fuzzy_equals(this->point(4) - this->point(1)) || 00112 !v.relative_fuzzy_equals(this->point(5) - this->point(2))) 00113 return false; 00114 return true; 00115 } 00116 00117 00118 00119 AutoPtr<Elem> Prism6::build_side (const unsigned int i, 00120 bool proxy) const 00121 { 00122 libmesh_assert_less (i, this->n_sides()); 00123 00124 if (proxy) 00125 { 00126 switch(i) 00127 { 00128 case 0: 00129 case 4: 00130 { 00131 AutoPtr<Elem> face(new Side<Tri3,Prism6>(this,i)); 00132 return face; 00133 } 00134 00135 case 1: 00136 case 2: 00137 case 3: 00138 { 00139 AutoPtr<Elem> face(new Side<Quad4,Prism6>(this,i)); 00140 return face; 00141 } 00142 00143 default: 00144 { 00145 libmesh_error(); 00146 } 00147 } 00148 } 00149 00150 else 00151 { 00152 switch (i) 00153 { 00154 case 0: // the triangular face at z=-1 00155 { 00156 AutoPtr<Elem> face(new Tri3); 00157 00158 face->set_node(0) = this->get_node(0); 00159 face->set_node(1) = this->get_node(2); 00160 face->set_node(2) = this->get_node(1); 00161 00162 return face; 00163 } 00164 case 1: // the quad face at y=0 00165 { 00166 AutoPtr<Elem> face(new Quad4); 00167 00168 face->set_node(0) = this->get_node(0); 00169 face->set_node(1) = this->get_node(1); 00170 face->set_node(2) = this->get_node(4); 00171 face->set_node(3) = this->get_node(3); 00172 00173 return face; 00174 } 00175 case 2: // the other quad face 00176 { 00177 AutoPtr<Elem> face(new Quad4); 00178 00179 face->set_node(0) = this->get_node(1); 00180 face->set_node(1) = this->get_node(2); 00181 face->set_node(2) = this->get_node(5); 00182 face->set_node(3) = this->get_node(4); 00183 00184 return face; 00185 } 00186 case 3: // the quad face at x=0 00187 { 00188 AutoPtr<Elem> face(new Quad4); 00189 00190 face->set_node(0) = this->get_node(2); 00191 face->set_node(1) = this->get_node(0); 00192 face->set_node(2) = this->get_node(3); 00193 face->set_node(3) = this->get_node(5); 00194 00195 return face; 00196 } 00197 case 4: // the triangular face at z=1 00198 { 00199 AutoPtr<Elem> face(new Tri3); 00200 00201 face->set_node(0) = this->get_node(3); 00202 face->set_node(1) = this->get_node(4); 00203 face->set_node(2) = this->get_node(5); 00204 00205 return face; 00206 } 00207 default: 00208 { 00209 libmesh_error(); 00210 } 00211 } 00212 } 00213 00214 // We'll never get here. 00215 libmesh_error(); 00216 AutoPtr<Elem> ap(NULL); return ap; 00217 } 00218 00219 00220 00221 AutoPtr<Elem> Prism6::build_edge (const unsigned int i) const 00222 { 00223 libmesh_assert_less (i, this->n_edges()); 00224 00225 return AutoPtr<Elem>(new SideEdge<Edge2,Prism6>(this,i)); 00226 } 00227 00228 00229 00230 void Prism6::connectivity(const unsigned int libmesh_dbg_var(sc), 00231 const IOPackage iop, 00232 std::vector<dof_id_type>& conn) const 00233 { 00234 libmesh_assert(_nodes); 00235 libmesh_assert_less (sc, this->n_sub_elem()); 00236 libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE); 00237 00238 switch (iop) 00239 { 00240 case TECPLOT: 00241 { 00242 conn.resize(8); 00243 conn[0] = this->node(0)+1; 00244 conn[1] = this->node(1)+1; 00245 conn[2] = this->node(2)+1; 00246 conn[3] = this->node(2)+1; 00247 conn[4] = this->node(3)+1; 00248 conn[5] = this->node(4)+1; 00249 conn[6] = this->node(5)+1; 00250 conn[7] = this->node(5)+1; 00251 return; 00252 } 00253 00254 case VTK: 00255 { 00256 conn.resize(6); 00257 conn[0] = this->node(0); 00258 conn[1] = this->node(2); 00259 conn[2] = this->node(1); 00260 conn[3] = this->node(3); 00261 conn[4] = this->node(5); 00262 conn[5] = this->node(4); 00263 return; 00264 } 00265 00266 default: 00267 libmesh_error(); 00268 } 00269 00270 libmesh_error(); 00271 } 00272 00273 00274 00275 #ifdef LIBMESH_ENABLE_AMR 00276 00277 const float Prism6::_embedding_matrix[8][6][6] = 00278 { 00279 // embedding matrix for child 0 00280 { 00281 // 0 1 2 3 4 5 00282 { 1.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0 00283 { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 1 00284 { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 2 00285 { 0.5, 0.0, 0.0, 0.5, 0.0, 0.0}, // 3 00286 { .25, .25, 0.0, .25, .25, 0.0}, // 4 00287 { .25, 0.0, .25, .25, 0.0, .25} // 5 00288 }, 00289 00290 // embedding matrix for child 1 00291 { 00292 // 0 1 2 3 4 5 00293 { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 0 00294 { 0.0, 1.0, 0.0, 0.0, 0.0, 0.0}, // 1 00295 { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 2 00296 { .25, .25, 0.0, .25, .25, 0.0}, // 3 00297 { 0.0, 0.5, 0.0, 0.0, 0.5, 0.0}, // 4 00298 { 0.0, .25, .25, 0.0, .25, .25} // 5 00299 }, 00300 00301 // embedding matrix for child 2 00302 { 00303 // 0 1 2 3 4 5 00304 { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 0 00305 { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 1 00306 { 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 2 00307 { .25, 0.0, .25, .25, 0.0, .25}, // 3 00308 { 0.0, .25, .25, 0.0, .25, .25}, // 4 00309 { 0.0, 0.0, 0.5, 0.0, 0.0, 0.5} // 5 00310 }, 00311 00312 // embedding matrix for child 3 00313 { 00314 // 0 1 2 3 4 5 00315 { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 0 00316 { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 1 00317 { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 2 00318 { .25, .25, 0.0, .25, .25, 0.0}, // 3 00319 { 0.0, .25, .25, 0.0, .25, .25}, // 4 00320 { .25, 0.0, .25, .25, 0.0, .25} // 5 00321 }, 00322 00323 // embedding matrix for child 4 00324 { 00325 // 0 1 2 3 4 5 00326 { 0.5, 0.0, 0.0, 0.5, 0.0, 0.0}, // 0 00327 { .25, .25, 0.0, .25, .25, 0.0}, // 1 00328 { .25, 0.0, .25, .25, 0.0, .25}, // 2 00329 { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 3 00330 { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 4 00331 { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5} // 5 00332 }, 00333 00334 // embedding matrix for child 5 00335 { 00336 // 0 1 2 3 4 5 00337 { .25, .25, 0.0, .25, .25, 0.0}, // 0 00338 { 0.0, 0.5, 0.0, 0.0, 0.5, 0.0}, // 1 00339 { 0.0, .25, .25, 0.0, .25, .25}, // 2 00340 { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 3 00341 { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 4 00342 { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5} // 5 00343 }, 00344 00345 // embedding matrix for child 6 00346 { 00347 // 0 1 2 3 4 5 00348 { .25, 0.0, .25, .25, 0.0, .25}, // 0 00349 { 0.0, .25, .25, 0.0, .25, .25}, // 1 00350 { 0.0, 0.0, 0.5, 0.0, 0.0, 0.5}, // 2 00351 { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5}, // 3 00352 { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5}, // 4 00353 { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0} // 5 00354 }, 00355 00356 // embedding matrix for child 7 00357 { 00358 // 0 1 2 3 4 5 00359 { .25, .25, 0.0, .25, .25, 0.0}, // 0 00360 { 0.0, .25, .25, 0.0, .25, .25}, // 1 00361 { .25, 0.0, .25, .25, 0.0, .25}, // 2 00362 { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 3 00363 { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5}, // 4 00364 { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5} // 5 00365 } 00366 }; 00367 00368 #endif 00369 00370 00371 00372 Real Prism6::volume () const 00373 { 00374 // The volume of the prism is computed by splitting 00375 // it into 2 tetrahedra and 3 pyramids with bilinear bases. 00376 // Then the volume formulae for the tetrahedron and pyramid 00377 // are applied and summed to obtain the prism's volume. 00378 00379 static const unsigned char sub_pyr[3][4] = 00380 { 00381 {0, 1, 4, 3}, 00382 {1, 2, 5, 4}, 00383 {0, 3, 5, 2} 00384 }; 00385 00386 static const unsigned char sub_tet[2][3] = 00387 { 00388 {0, 1, 2}, 00389 {5, 4, 3} 00390 }; 00391 00392 // The centroid is a convenient point to use 00393 // for the apex of all the pyramids. 00394 const Point R = this->centroid(); 00395 00396 // temporary storage for Nodes which form the base of the 00397 // subelements 00398 Node* base[4]; 00399 00400 // volume accumulation variable 00401 Real vol=0.; 00402 00403 // Add up the sub-pyramid volumes 00404 for (unsigned int n=0; n<3; ++n) 00405 { 00406 // Set the nodes of the pyramid base 00407 for (unsigned int i=0; i<4; ++i) 00408 base[i] = this->_nodes[sub_pyr[n][i]]; 00409 00410 // Compute diff vectors 00411 Point a ( *base[0] - R ); 00412 Point b ( *base[1] - *base[3] ); 00413 Point c ( *base[2] - *base[0] ); 00414 Point d ( *base[3] - *base[0] ); 00415 Point e ( *base[1] - *base[0] ); 00416 00417 // Compute pyramid volume 00418 Real sub_vol = (1./6.)*(a*(b.cross(c))) + (1./12.)*(c*(d.cross(e))); 00419 00420 libmesh_assert (sub_vol>0.); 00421 00422 vol += sub_vol; 00423 } 00424 00425 00426 // Add up the sub-tet volumes 00427 for (unsigned int n=0; n<2; ++n) 00428 { 00429 // Set the nodes of the pyramid base 00430 for (unsigned int i=0; i<3; ++i) 00431 base[i] = this->_nodes[sub_tet[n][i]]; 00432 00433 // The volume of a tetrahedron is 1/6 the box product formed 00434 // by its base and apex vectors 00435 Point a ( R - *base[0] ); 00436 00437 // b is the vector pointing from 0 to 1 00438 Point b ( *base[1] - *base[0] ); 00439 00440 // c is the vector pointing from 0 to 2 00441 Point c ( *base[2] - *base[0] ); 00442 00443 Real sub_vol = (1.0 / 6.0) * (a * (b.cross(c))); 00444 00445 libmesh_assert (sub_vol>0.); 00446 00447 vol += sub_vol; 00448 } 00449 00450 00451 // Done with all sub-volumes, so return 00452 return vol; 00453 } 00454 00455 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:45 UTC
Hosted By: