$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 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 UniquePtr<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 return UniquePtr<Elem>(new Side<Tri3,Prism6>(this,i)); 00131 00132 case 1: 00133 case 2: 00134 case 3: 00135 return UniquePtr<Elem>(new Side<Quad4,Prism6>(this,i)); 00136 00137 default: 00138 libmesh_error_msg("Invalid side i = " << i); 00139 } 00140 } 00141 00142 else 00143 { 00144 // Create NULL pointer to be initialized, returned later. 00145 Elem* face = NULL; 00146 00147 switch (i) 00148 { 00149 case 0: // the triangular face at z=-1 00150 { 00151 face = new Tri3; 00152 00153 face->set_node(0) = this->get_node(0); 00154 face->set_node(1) = this->get_node(2); 00155 face->set_node(2) = this->get_node(1); 00156 00157 break; 00158 } 00159 case 1: // the quad face at y=0 00160 { 00161 face = new Quad4; 00162 00163 face->set_node(0) = this->get_node(0); 00164 face->set_node(1) = this->get_node(1); 00165 face->set_node(2) = this->get_node(4); 00166 face->set_node(3) = this->get_node(3); 00167 00168 break; 00169 } 00170 case 2: // the other quad face 00171 { 00172 face = new Quad4; 00173 00174 face->set_node(0) = this->get_node(1); 00175 face->set_node(1) = this->get_node(2); 00176 face->set_node(2) = this->get_node(5); 00177 face->set_node(3) = this->get_node(4); 00178 00179 break; 00180 } 00181 case 3: // the quad face at x=0 00182 { 00183 face = new Quad4; 00184 00185 face->set_node(0) = this->get_node(2); 00186 face->set_node(1) = this->get_node(0); 00187 face->set_node(2) = this->get_node(3); 00188 face->set_node(3) = this->get_node(5); 00189 00190 break; 00191 } 00192 case 4: // the triangular face at z=1 00193 { 00194 face = new Tri3; 00195 00196 face->set_node(0) = this->get_node(3); 00197 face->set_node(1) = this->get_node(4); 00198 face->set_node(2) = this->get_node(5); 00199 00200 break; 00201 } 00202 default: 00203 libmesh_error_msg("Invalid side i = " << i); 00204 } 00205 00206 face->subdomain_id() = this->subdomain_id(); 00207 return UniquePtr<Elem>(face); 00208 } 00209 00210 libmesh_error_msg("We'll never get here!"); 00211 return UniquePtr<Elem>(); 00212 } 00213 00214 00215 00216 UniquePtr<Elem> Prism6::build_edge (const unsigned int i) const 00217 { 00218 libmesh_assert_less (i, this->n_edges()); 00219 00220 return UniquePtr<Elem>(new SideEdge<Edge2,Prism6>(this,i)); 00221 } 00222 00223 00224 00225 void Prism6::connectivity(const unsigned int libmesh_dbg_var(sc), 00226 const IOPackage iop, 00227 std::vector<dof_id_type>& conn) const 00228 { 00229 libmesh_assert(_nodes); 00230 libmesh_assert_less (sc, this->n_sub_elem()); 00231 libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE); 00232 00233 switch (iop) 00234 { 00235 case TECPLOT: 00236 { 00237 conn.resize(8); 00238 conn[0] = this->node(0)+1; 00239 conn[1] = this->node(1)+1; 00240 conn[2] = this->node(2)+1; 00241 conn[3] = this->node(2)+1; 00242 conn[4] = this->node(3)+1; 00243 conn[5] = this->node(4)+1; 00244 conn[6] = this->node(5)+1; 00245 conn[7] = this->node(5)+1; 00246 return; 00247 } 00248 00249 case VTK: 00250 { 00251 conn.resize(6); 00252 conn[0] = this->node(0); 00253 conn[1] = this->node(2); 00254 conn[2] = this->node(1); 00255 conn[3] = this->node(3); 00256 conn[4] = this->node(5); 00257 conn[5] = this->node(4); 00258 return; 00259 } 00260 00261 default: 00262 libmesh_error_msg("Unsupported IO package " << iop); 00263 } 00264 } 00265 00266 00267 00268 #ifdef LIBMESH_ENABLE_AMR 00269 00270 const float Prism6::_embedding_matrix[8][6][6] = 00271 { 00272 // embedding matrix for child 0 00273 { 00274 // 0 1 2 3 4 5 00275 { 1.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0 00276 { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 1 00277 { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 2 00278 { 0.5, 0.0, 0.0, 0.5, 0.0, 0.0}, // 3 00279 { .25, .25, 0.0, .25, .25, 0.0}, // 4 00280 { .25, 0.0, .25, .25, 0.0, .25} // 5 00281 }, 00282 00283 // embedding matrix for child 1 00284 { 00285 // 0 1 2 3 4 5 00286 { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 0 00287 { 0.0, 1.0, 0.0, 0.0, 0.0, 0.0}, // 1 00288 { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 2 00289 { .25, .25, 0.0, .25, .25, 0.0}, // 3 00290 { 0.0, 0.5, 0.0, 0.0, 0.5, 0.0}, // 4 00291 { 0.0, .25, .25, 0.0, .25, .25} // 5 00292 }, 00293 00294 // embedding matrix for child 2 00295 { 00296 // 0 1 2 3 4 5 00297 { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 0 00298 { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 1 00299 { 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 2 00300 { .25, 0.0, .25, .25, 0.0, .25}, // 3 00301 { 0.0, .25, .25, 0.0, .25, .25}, // 4 00302 { 0.0, 0.0, 0.5, 0.0, 0.0, 0.5} // 5 00303 }, 00304 00305 // embedding matrix for child 3 00306 { 00307 // 0 1 2 3 4 5 00308 { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 0 00309 { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0}, // 1 00310 { 0.5, 0.0, 0.5, 0.0, 0.0, 0.0}, // 2 00311 { .25, .25, 0.0, .25, .25, 0.0}, // 3 00312 { 0.0, .25, .25, 0.0, .25, .25}, // 4 00313 { .25, 0.0, .25, .25, 0.0, .25} // 5 00314 }, 00315 00316 // embedding matrix for child 4 00317 { 00318 // 0 1 2 3 4 5 00319 { 0.5, 0.0, 0.0, 0.5, 0.0, 0.0}, // 0 00320 { .25, .25, 0.0, .25, .25, 0.0}, // 1 00321 { .25, 0.0, .25, .25, 0.0, .25}, // 2 00322 { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 3 00323 { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 4 00324 { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5} // 5 00325 }, 00326 00327 // embedding matrix for child 5 00328 { 00329 // 0 1 2 3 4 5 00330 { .25, .25, 0.0, .25, .25, 0.0}, // 0 00331 { 0.0, 0.5, 0.0, 0.0, 0.5, 0.0}, // 1 00332 { 0.0, .25, .25, 0.0, .25, .25}, // 2 00333 { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 3 00334 { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 4 00335 { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5} // 5 00336 }, 00337 00338 // embedding matrix for child 6 00339 { 00340 // 0 1 2 3 4 5 00341 { .25, 0.0, .25, .25, 0.0, .25}, // 0 00342 { 0.0, .25, .25, 0.0, .25, .25}, // 1 00343 { 0.0, 0.0, 0.5, 0.0, 0.0, 0.5}, // 2 00344 { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5}, // 3 00345 { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5}, // 4 00346 { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0} // 5 00347 }, 00348 00349 // embedding matrix for child 7 00350 { 00351 // 0 1 2 3 4 5 00352 { .25, .25, 0.0, .25, .25, 0.0}, // 0 00353 { 0.0, .25, .25, 0.0, .25, .25}, // 1 00354 { .25, 0.0, .25, .25, 0.0, .25}, // 2 00355 { 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 3 00356 { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5}, // 4 00357 { 0.0, 0.0, 0.0, 0.5, 0.0, 0.5} // 5 00358 } 00359 }; 00360 00361 #endif 00362 00363 00364 00365 Real Prism6::volume () const 00366 { 00367 // The volume of the prism is computed by splitting 00368 // it into 2 tetrahedra and 3 pyramids with bilinear bases. 00369 // Then the volume formulae for the tetrahedron and pyramid 00370 // are applied and summed to obtain the prism's volume. 00371 00372 static const unsigned char sub_pyr[3][4] = 00373 { 00374 {0, 1, 4, 3}, 00375 {1, 2, 5, 4}, 00376 {0, 3, 5, 2} 00377 }; 00378 00379 static const unsigned char sub_tet[2][3] = 00380 { 00381 {0, 1, 2}, 00382 {5, 4, 3} 00383 }; 00384 00385 // The centroid is a convenient point to use 00386 // for the apex of all the pyramids. 00387 const Point R = this->centroid(); 00388 00389 // temporary storage for Nodes which form the base of the 00390 // subelements 00391 Node* base[4]; 00392 00393 // volume accumulation variable 00394 Real vol=0.; 00395 00396 // Add up the sub-pyramid volumes 00397 for (unsigned int n=0; n<3; ++n) 00398 { 00399 // Set the nodes of the pyramid base 00400 for (unsigned int i=0; i<4; ++i) 00401 base[i] = this->_nodes[sub_pyr[n][i]]; 00402 00403 // Compute diff vectors 00404 Point a ( *base[0] - R ); 00405 Point b ( *base[1] - *base[3] ); 00406 Point c ( *base[2] - *base[0] ); 00407 Point d ( *base[3] - *base[0] ); 00408 Point e ( *base[1] - *base[0] ); 00409 00410 // Compute pyramid volume 00411 Real sub_vol = (1./6.)*(a*(b.cross(c))) + (1./12.)*(c*(d.cross(e))); 00412 00413 libmesh_assert (sub_vol>0.); 00414 00415 vol += sub_vol; 00416 } 00417 00418 00419 // Add up the sub-tet volumes 00420 for (unsigned int n=0; n<2; ++n) 00421 { 00422 // Set the nodes of the pyramid base 00423 for (unsigned int i=0; i<3; ++i) 00424 base[i] = this->_nodes[sub_tet[n][i]]; 00425 00426 // The volume of a tetrahedron is 1/6 the box product formed 00427 // by its base and apex vectors 00428 Point a ( R - *base[0] ); 00429 00430 // b is the vector pointing from 0 to 1 00431 Point b ( *base[1] - *base[0] ); 00432 00433 // c is the vector pointing from 0 to 2 00434 Point c ( *base[2] - *base[0] ); 00435 00436 Real sub_vol = (1.0 / 6.0) * (a * (b.cross(c))); 00437 00438 libmesh_assert (sub_vol>0.); 00439 00440 vol += sub_vol; 00441 } 00442 00443 00444 // Done with all sub-volumes, so return 00445 return vol; 00446 } 00447 00448 } // namespace libMesh