$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 // Local includes 00019 #include "libmesh/libmesh_config.h" 00020 00021 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00022 00023 // C++ includes 00024 00025 // Local includes cont'd 00026 #include "libmesh/cell_inf_prism12.h" 00027 #include "libmesh/edge_edge3.h" 00028 #include "libmesh/edge_inf_edge2.h" 00029 #include "libmesh/face_tri6.h" 00030 #include "libmesh/face_inf_quad6.h" 00031 #include "libmesh/side.h" 00032 00033 namespace libMesh 00034 { 00035 00036 00037 // ------------------------------------------------------------ 00038 // InfPrism12 class static member initializations 00039 const unsigned int InfPrism12::side_nodes_map[4][6] = 00040 { 00041 { 0, 1, 2, 6, 7, 8}, // Side 0 00042 { 0, 1, 3, 4, 6, 9}, // Side 1 00043 { 1, 2, 4, 5, 7, 10}, // Side 2 00044 { 2, 0, 5, 3, 8, 11} // Side 3 00045 }; 00046 00047 const unsigned int InfPrism12::edge_nodes_map[6][3] = 00048 { 00049 { 0, 1, 6}, // Side 0 00050 { 1, 2, 7}, // Side 1 00051 { 0, 2, 8}, // Side 2 00052 { 0, 3, 99}, // Side 3 00053 { 1, 4, 99}, // Side 4 00054 { 2, 5, 99} // Side 5 00055 }; 00056 00057 00058 // ------------------------------------------------------------ 00059 // InfPrism12 class member functions 00060 00061 bool InfPrism12::is_vertex(const unsigned int i) const 00062 { 00063 if (i < 3) 00064 return true; 00065 return false; 00066 } 00067 00068 bool InfPrism12::is_edge(const unsigned int i) const 00069 { 00070 if (i < 3) 00071 return false; 00072 if (i > 8) 00073 return false; 00074 return true; 00075 } 00076 00077 bool InfPrism12::is_face(const unsigned int i) const 00078 { 00079 if (i > 8) 00080 return true; 00081 return false; 00082 } 00083 00084 bool InfPrism12::is_node_on_side(const unsigned int n, 00085 const unsigned int s) const 00086 { 00087 libmesh_assert_less (s, n_sides()); 00088 for (unsigned int i = 0; i != 6; ++i) 00089 if (side_nodes_map[s][i] == n) 00090 return true; 00091 return false; 00092 } 00093 00094 bool InfPrism12::is_node_on_edge(const unsigned int n, 00095 const unsigned int e) const 00096 { 00097 libmesh_assert_less (e, n_edges()); 00098 for (unsigned int i = 0; i != 3; ++i) 00099 if (edge_nodes_map[e][i] == n) 00100 return true; 00101 return false; 00102 } 00103 00104 UniquePtr<Elem> InfPrism12::build_side (const unsigned int i, 00105 bool proxy) const 00106 { 00107 libmesh_assert_less (i, this->n_sides()); 00108 00109 if (proxy) 00110 { 00111 switch (i) 00112 { 00113 // base 00114 case 0: 00115 return UniquePtr<Elem>(new Side<Tri6,InfPrism12>(this,i)); 00116 00117 // ifem sides 00118 case 1: 00119 case 2: 00120 case 3: 00121 return UniquePtr<Elem>(new Side<InfQuad6,InfPrism12>(this,i)); 00122 00123 default: 00124 libmesh_error_msg("Invalid side i = " << i); 00125 } 00126 } 00127 00128 else 00129 { 00130 // Create NULL pointer to be initialized, returned later. 00131 Elem* face = NULL; 00132 00133 switch (i) 00134 { 00135 case 0: // the triangular face at z=-1, base face 00136 { 00137 face = new Tri6; 00138 00139 // Note that for this face element, the normal points inward 00140 face->set_node(0) = this->get_node(0); 00141 face->set_node(1) = this->get_node(1); 00142 face->set_node(2) = this->get_node(2); 00143 face->set_node(3) = this->get_node(6); 00144 face->set_node(4) = this->get_node(7); 00145 face->set_node(5) = this->get_node(8); 00146 00147 break; 00148 } 00149 00150 case 1: // the quad face at y=0 00151 { 00152 face = new InfQuad6; 00153 00154 face->set_node(0) = this->get_node(0); 00155 face->set_node(1) = this->get_node(1); 00156 face->set_node(2) = this->get_node(3); 00157 face->set_node(3) = this->get_node(4); 00158 face->set_node(4) = this->get_node(6); 00159 face->set_node(5) = this->get_node(9); 00160 00161 break; 00162 } 00163 00164 case 2: // the other quad face 00165 { 00166 face = new InfQuad6; 00167 00168 face->set_node(0) = this->get_node(1); 00169 face->set_node(1) = this->get_node(2); 00170 face->set_node(2) = this->get_node(4); 00171 face->set_node(3) = this->get_node(5); 00172 face->set_node(4) = this->get_node(7); 00173 face->set_node(5) = this->get_node(10); 00174 00175 break; 00176 } 00177 00178 case 3: // the quad face at x=0 00179 { 00180 face = new InfQuad6; 00181 00182 face->set_node(0) = this->get_node(2); 00183 face->set_node(1) = this->get_node(0); 00184 face->set_node(2) = this->get_node(5); 00185 face->set_node(3) = this->get_node(3); 00186 face->set_node(4) = this->get_node(8); 00187 face->set_node(5) = this->get_node(11); 00188 00189 break; 00190 } 00191 00192 default: 00193 libmesh_error_msg("Invalid side i = " << i); 00194 } 00195 00196 face->subdomain_id() = this->subdomain_id(); 00197 return UniquePtr<Elem>(face); 00198 } 00199 00200 libmesh_error_msg("We'll never get here!"); 00201 return UniquePtr<Elem>(); 00202 } 00203 00204 00205 UniquePtr<Elem> InfPrism12::build_edge (const unsigned int i) const 00206 { 00207 libmesh_assert_less (i, this->n_edges()); 00208 00209 if (i < 3) // base edges 00210 return UniquePtr<Elem>(new SideEdge<Edge3,InfPrism12>(this,i)); 00211 // infinite edges 00212 return UniquePtr<Elem>(new SideEdge<InfEdge2,InfPrism12>(this,i)); 00213 } 00214 00215 00216 void InfPrism12::connectivity(const unsigned int sc, 00217 const IOPackage iop, 00218 std::vector<dof_id_type>& conn) const 00219 { 00220 libmesh_assert(_nodes); 00221 libmesh_assert_less (sc, this->n_sub_elem()); 00222 libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE); 00223 00224 switch (iop) 00225 { 00226 case TECPLOT: 00227 { 00228 conn.resize(8); 00229 switch (sc) 00230 { 00231 case 0: 00232 00233 // guess this is a collapsed hex8 00234 conn[0] = this->node(0)+1; 00235 conn[1] = this->node(6)+1; 00236 conn[2] = this->node(8)+1; 00237 conn[3] = this->node(8)+1; 00238 conn[4] = this->node(3)+1; 00239 conn[5] = this->node(9)+1; 00240 conn[6] = this->node(11)+1; 00241 conn[7] = this->node(11)+1; 00242 00243 return; 00244 00245 case 1: 00246 00247 conn[0] = this->node(6)+1; 00248 conn[1] = this->node(7)+1; 00249 conn[2] = this->node(8)+1; 00250 conn[3] = this->node(8)+1; 00251 conn[4] = this->node(9)+1; 00252 conn[5] = this->node(10)+1; 00253 conn[6] = this->node(11)+1; 00254 conn[7] = this->node(11)+1; 00255 00256 return; 00257 00258 case 2: 00259 00260 conn[0] = this->node(6)+1; 00261 conn[1] = this->node(1)+1; 00262 conn[2] = this->node(7)+1; 00263 conn[3] = this->node(7)+1; 00264 conn[4] = this->node(9)+1; 00265 conn[5] = this->node(4)+1; 00266 conn[6] = this->node(10)+1; 00267 conn[7] = this->node(10)+1; 00268 00269 return; 00270 00271 case 3: 00272 00273 conn[0] = this->node(8)+1; 00274 conn[1] = this->node(7)+1; 00275 conn[2] = this->node(2)+1; 00276 conn[3] = this->node(2)+1; 00277 conn[4] = this->node(11)+1; 00278 conn[5] = this->node(10)+1; 00279 conn[6] = this->node(5)+1; 00280 conn[7] = this->node(5)+1; 00281 00282 return; 00283 00284 default: 00285 libmesh_error_msg("Invalid sc = " << sc); 00286 } 00287 } 00288 00289 default: 00290 libmesh_error_msg("Unsupported IO package " << iop); 00291 } 00292 } 00293 00294 00295 00296 00297 00298 unsigned short int InfPrism12::second_order_adjacent_vertex (const unsigned int n, 00299 const unsigned int v) const 00300 { 00301 libmesh_assert_greater_equal (n, this->n_vertices()); 00302 libmesh_assert_less (n, this->n_nodes()); 00303 libmesh_assert_less (v, 2); 00304 return _second_order_adjacent_vertices[n-this->n_vertices()][v]; 00305 } 00306 00307 00308 00309 const unsigned short int InfPrism12::_second_order_adjacent_vertices[6][2] = 00310 { 00311 { 0, 1}, // vertices adjacent to node 6 00312 { 1, 2}, // vertices adjacent to node 7 00313 { 0, 2}, // vertices adjacent to node 8 00314 00315 { 3, 4}, // vertices adjacent to node 9 00316 { 4, 5}, // vertices adjacent to node 10 00317 { 3, 5} // vertices adjacent to node 11 00318 }; 00319 00320 00321 00322 std::pair<unsigned short int, unsigned short int> 00323 InfPrism12::second_order_child_vertex (const unsigned int n) const 00324 { 00325 libmesh_assert_greater_equal (n, this->n_vertices()); 00326 libmesh_assert_less (n, this->n_nodes()); 00327 00328 return std::pair<unsigned short int, unsigned short int> 00329 (_second_order_vertex_child_number[n], 00330 _second_order_vertex_child_index[n]); 00331 } 00332 00333 00334 00335 const unsigned short int InfPrism12::_second_order_vertex_child_number[12] = 00336 { 00337 99,99,99,99,99,99, // Vertices 00338 0,1,0, // Edges 00339 0,1,0 // Faces 00340 }; 00341 00342 00343 00344 const unsigned short int InfPrism12::_second_order_vertex_child_index[12] = 00345 { 00346 99,99,99,99,99,99, // Vertices 00347 1,2,2, // Edges 00348 4,5,5 // Faces 00349 }; 00350 00351 00352 00353 #ifdef LIBMESH_ENABLE_AMR 00354 00355 const float InfPrism12::_embedding_matrix[4][12][12] = 00356 { 00357 // embedding matrix for child 0 00358 { 00359 // 0 1 2 3 4 5 6 7 8 9 10 11 th parent Node 00360 { 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0th child N. 00361 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 1 00362 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 2 00363 { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 3 00364 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 4 00365 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0}, // 5 00366 { 0.375, -0.125, 0.0, 0.0, 0.0, 0.0, 0.75, 0.0, 0.0, 0.0, 0.0, 0.0}, // 6 00367 { 0.0, -0.125, -0.125, 0.0, 0.0, 0.0, 0.5, 0.25, 0.5, 0.0, 0.0, 0.0}, // 7 00368 { 0.375, 0.0, -0.125, 0.0, 0.0, 0.0, 0.0, 0.0, 0.75, 0.0, 0.0, 0.0}, // 8 00369 { 0.0, 0.0, 0.0, 0.375, -0.125, 0.0, 0.0, 0.0, 0.0, 0.75, 0.0, 0.0}, // 9 00370 { 0.0, 0.0, 0.0, 0.0, -0.125, -0.125, 0.0, 0.0, 0.0, 0.5, 0.25, 0.5}, // 10 00371 { 0.0, 0.0, 0.0, 0.375, 0.0, -0.125, 0.0, 0.0, 0.0, 0.0, 0.0, 0.75} // 11 00372 }, 00373 00374 // embedding matrix for child 1 00375 { 00376 // 0 1 2 3 4 5 6 7 8 9 10 11 th parent Node 00377 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0th child N. 00378 { 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 1 00379 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0}, // 2 00380 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 3 00381 { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 4 00382 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 5 00383 { -0.125, 0.375, 0.0, 0.0, 0.0, 0.0, 0.75, 0.0, 0.0, 0.0, 0.0, 0.0}, // 6 00384 { 0.0, 0.375, -0.125, 0.0, 0.0, 0.0, 0.0, 0.75, 0.0, 0.0, 0.0, 0.0}, // 7 00385 { -0.125, 0.0, -0.125, 0.0, 0.0, 0.0, 0.5, 0.5, 0.25, 0.0, 0.0, 0.0}, // 8 00386 { 0.0, 0.0, 0.0, -0.125, 0.375, 0.0, 0.0, 0.0, 0.0, 0.75, 0.0, 0.0}, // 9 00387 { 0.0, 0.0, 0.0, 0.0, 0.375, -0.125, 0.0, 0.0, 0.0, 0.0, 0.75, 0.0}, // 10 00388 { 0.0, 0.0, 0.0, -0.125, 0.0, -0.125, 0.0, 0.0, 0.0, 0.5, 0.5, 0.25} // 11 00389 }, 00390 00391 // embedding matrix for child 2 00392 { 00393 // 0 1 2 3 4 5 6 7 8 9 10 11 th parent Node 00394 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 0th child N. 00395 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0}, // 1 00396 { 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 2 00397 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0}, // 3 00398 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 4 00399 { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 5 00400 { -0.125, -0.125, 0.0, 0.0, 0.0, 0.0, 0.25, 0.5, 0.5, 0.0, 0.0, 0.0}, // 6 00401 { 0.0, -0.125, 0.375, 0.0, 0.0, 0.0, 0.0, 0.75, 0.0, 0.0, 0.0, 0.0}, // 7 00402 { -0.125, 0.0, 0.375, 0.0, 0.0, 0.0, 0.0, 0.0, 0.75, 0.0, 0.0, 0.0}, // 8 00403 { 0.0, 0.0, 0.0, -0.125, -0.125, 0.0, 0.0, 0.0, 0.0, 0.25, 0.5, 0.5}, // 9 00404 { 0.0, 0.0, 0.0, 0.0, -0.125, 0.375, 0.0, 0.0, 0.0, 0.0, 0.75, 0.0}, // 10 00405 { 0.0, 0.0, 0.0, -0.125, 0.0, 0.375, 0.0, 0.0, 0.0, 0.0, 0.0, 0.75} // 11 00406 }, 00407 00408 // embedding matrix for child 3 00409 { 00410 // 0 1 2 3 4 5 6 7 8 9 10 11 th parent Node 00411 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0th child N. 00412 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0}, // 1 00413 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 2 00414 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 3 00415 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 4 00416 { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0}, // 5 00417 { -0.125, 0.0, -0.125, 0.0, 0.0, 0.0, 0.5, 0.5, 0.25, 0.0, 0.0, 0.0}, // 6 00418 { -0.125, -0.125, 0.0, 0.0, 0.0, 0.0, 0.25, 0.5, 0.5, 0.0, 0.0, 0.0}, // 7 00419 { 0.0, -0.125, -0.125, 0.0, 0.0, 0.0, 0.5, 0.25, 0.5, 0.0, 0.0, 0.0}, // 8 00420 { 0.0, 0.0, 0.0, -0.125, 0.0, -0.125, 0.0, 0.0, 0.0, 0.5, 0.5, 0.25}, // 9 00421 { 0.0, 0.0, 0.0, -0.125, -0.125, 0.0, 0.0, 0.0, 0.0, 0.25, 0.5, 0.5}, // 10 00422 { 0.0, 0.0, 0.0, 0.0, -0.125, -0.125, 0.0, 0.0, 0.0, 0.5, 0.25, 0.5} // 11 00423 } 00424 00425 }; 00426 00427 00428 00429 00430 #endif 00431 00432 } // namespace libMesh 00433 00434 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS