cell_tet10.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_tet10.h" 00024 #include "libmesh/edge_edge3.h" 00025 #include "libmesh/face_tri6.h" 00026 00027 namespace libMesh 00028 { 00029 00030 00031 00032 // ------------------------------------------------------------ 00033 // Tet10 class static member initializations 00034 const unsigned int Tet10::side_nodes_map[4][6] = 00035 { 00036 {0, 2, 1, 6, 5, 4}, // Side 0 00037 {0, 1, 3, 4, 8, 7}, // Side 1 00038 {1, 2, 3, 5, 9, 8}, // Side 2 00039 {2, 0, 3, 6, 7, 9} // Side 3 00040 }; 00041 00042 const unsigned int Tet10::edge_nodes_map[6][3] = 00043 { 00044 {0, 1, 4}, // Side 0 00045 {1, 2, 5}, // Side 1 00046 {0, 2, 6}, // Side 2 00047 {0, 3, 7}, // Side 3 00048 {1, 3, 8}, // Side 4 00049 {2, 3, 9} // Side 5 00050 }; 00051 00052 00053 00054 // ------------------------------------------------------------ 00055 // Tet10 class member functions 00056 00057 bool Tet10::is_vertex(const unsigned int i) const 00058 { 00059 if (i < 4) 00060 return true; 00061 return false; 00062 } 00063 00064 bool Tet10::is_edge(const unsigned int i) const 00065 { 00066 if (i < 4) 00067 return false; 00068 return true; 00069 } 00070 00071 bool Tet10::is_face(const unsigned int) const 00072 { 00073 return false; 00074 } 00075 00076 bool Tet10::is_node_on_side(const unsigned int n, 00077 const unsigned int s) const 00078 { 00079 libmesh_assert_less (s, n_sides()); 00080 for (unsigned int i = 0; i != 6; ++i) 00081 if (side_nodes_map[s][i] == n) 00082 return true; 00083 return false; 00084 } 00085 00086 bool Tet10::is_node_on_edge(const unsigned int n, 00087 const unsigned int e) const 00088 { 00089 libmesh_assert_less (e, n_edges()); 00090 for (unsigned int i = 0; i != 3; ++i) 00091 if (edge_nodes_map[e][i] == n) 00092 return true; 00093 return false; 00094 } 00095 00096 00097 #ifdef LIBMESH_ENABLE_AMR 00098 00099 // This function only works if LIBMESH_ENABLE_AMR... 00100 bool Tet10::is_child_on_side(const unsigned int c, 00101 const unsigned int s) const 00102 { 00103 // Table of local IDs for the midege nodes on the side opposite a given node. 00104 // See the ASCII art in the header file for this class to confirm this. 00105 const unsigned int midedge_nodes_opposite[4][3] = 00106 { 00107 {5,8,9}, // midedge nodes opposite node 0 00108 {6,7,9}, // midedge nodes opposite node 1 00109 {4,7,8}, // midedge nodes opposite node 2 00110 {4,5,6} // midedge nodes opposite node 3 00111 }; 00112 00113 // Call the base class helper function 00114 return Tet::is_child_on_side_helper(c, s, midedge_nodes_opposite); 00115 } 00116 00117 #else 00118 00119 bool Tet10::is_child_on_side(const unsigned int /*c*/, 00120 const unsigned int /*s*/) const 00121 { 00122 libmesh_not_implemented(); 00123 return false; 00124 } 00125 00126 #endif //LIBMESH_ENABLE_AMR 00127 00128 00129 00130 bool Tet10::has_affine_map() const 00131 { 00132 // Make sure edges are straight 00133 if (!this->point(4).relative_fuzzy_equals 00134 ((this->point(0) + this->point(1))/2)) 00135 return false; 00136 if (!this->point(5).relative_fuzzy_equals 00137 ((this->point(1) + this->point(2))/2)) 00138 return false; 00139 if (!this->point(6).relative_fuzzy_equals 00140 ((this->point(2) + this->point(0))/2)) 00141 return false; 00142 if (!this->point(7).relative_fuzzy_equals 00143 ((this->point(3) + this->point(0))/2)) 00144 return false; 00145 if (!this->point(8).relative_fuzzy_equals 00146 ((this->point(3) + this->point(1))/2)) 00147 return false; 00148 if (!this->point(9).relative_fuzzy_equals 00149 ((this->point(3) + this->point(2))/2)) 00150 return false; 00151 return true; 00152 } 00153 00154 00155 00156 AutoPtr<Elem> Tet10::build_side (const unsigned int i, 00157 bool proxy) const 00158 { 00159 libmesh_assert_less (i, this->n_sides()); 00160 00161 if (proxy) 00162 { 00163 AutoPtr<Elem> ap(new Side<Tri6,Tet10>(this,i)); 00164 return ap; 00165 } 00166 00167 else 00168 { 00169 AutoPtr<Elem> face(new Tri6); 00170 00171 switch (i) 00172 { 00173 case 0: 00174 { 00175 face->set_node(0) = this->get_node(0); 00176 face->set_node(1) = this->get_node(2); 00177 face->set_node(2) = this->get_node(1); 00178 face->set_node(3) = this->get_node(6); 00179 face->set_node(4) = this->get_node(5); 00180 face->set_node(5) = this->get_node(4); 00181 00182 return face; 00183 } 00184 case 1: 00185 { 00186 face->set_node(0) = this->get_node(0); 00187 face->set_node(1) = this->get_node(1); 00188 face->set_node(2) = this->get_node(3); 00189 face->set_node(3) = this->get_node(4); 00190 face->set_node(4) = this->get_node(8); 00191 face->set_node(5) = this->get_node(7); 00192 00193 return face; 00194 } 00195 case 2: 00196 { 00197 face->set_node(0) = this->get_node(1); 00198 face->set_node(1) = this->get_node(2); 00199 face->set_node(2) = this->get_node(3); 00200 face->set_node(3) = this->get_node(5); 00201 face->set_node(4) = this->get_node(9); 00202 face->set_node(5) = this->get_node(8); 00203 00204 return face; 00205 } 00206 case 3: 00207 { 00208 face->set_node(0) = this->get_node(2); 00209 face->set_node(1) = this->get_node(0); 00210 face->set_node(2) = this->get_node(3); 00211 face->set_node(3) = this->get_node(6); 00212 face->set_node(4) = this->get_node(7); 00213 face->set_node(5) = this->get_node(9); 00214 00215 return face; 00216 } 00217 default: 00218 { 00219 libmesh_error(); 00220 } 00221 } 00222 } 00223 00224 00225 // We'll never get here. 00226 libmesh_error(); 00227 AutoPtr<Elem> ap(NULL); return ap; 00228 } 00229 00230 00231 00232 AutoPtr<Elem> Tet10::build_edge (const unsigned int i) const 00233 { 00234 libmesh_assert_less (i, this->n_edges()); 00235 00236 return AutoPtr<Elem>(new SideEdge<Edge3,Tet10>(this,i)); 00237 } 00238 00239 00240 00241 void Tet10::connectivity(const unsigned int sc, 00242 const IOPackage iop, 00243 std::vector<dof_id_type>& conn) const 00244 { 00245 libmesh_assert(_nodes); 00246 libmesh_assert_less (sc, this->n_sub_elem()); 00247 libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE); 00248 00249 switch (iop) 00250 { 00251 case TECPLOT: 00252 { 00253 conn.resize(8); 00254 switch (sc) 00255 { 00256 00257 00258 // Linear sub-tet 0 00259 case 0: 00260 00261 conn[0] = this->node(0)+1; 00262 conn[1] = this->node(4)+1; 00263 conn[2] = this->node(6)+1; 00264 conn[3] = this->node(6)+1; 00265 conn[4] = this->node(7)+1; 00266 conn[5] = this->node(7)+1; 00267 conn[6] = this->node(7)+1; 00268 conn[7] = this->node(7)+1; 00269 00270 return; 00271 00272 // Linear sub-tet 1 00273 case 1: 00274 00275 conn[0] = this->node(4)+1; 00276 conn[1] = this->node(1)+1; 00277 conn[2] = this->node(5)+1; 00278 conn[3] = this->node(5)+1; 00279 conn[4] = this->node(8)+1; 00280 conn[5] = this->node(8)+1; 00281 conn[6] = this->node(8)+1; 00282 conn[7] = this->node(8)+1; 00283 00284 return; 00285 00286 // Linear sub-tet 2 00287 case 2: 00288 00289 conn[0] = this->node(5)+1; 00290 conn[1] = this->node(2)+1; 00291 conn[2] = this->node(6)+1; 00292 conn[3] = this->node(6)+1; 00293 conn[4] = this->node(9)+1; 00294 conn[5] = this->node(9)+1; 00295 conn[6] = this->node(9)+1; 00296 conn[7] = this->node(9)+1; 00297 00298 return; 00299 00300 // Linear sub-tet 3 00301 case 3: 00302 00303 conn[0] = this->node(7)+1; 00304 conn[1] = this->node(8)+1; 00305 conn[2] = this->node(9)+1; 00306 conn[3] = this->node(9)+1; 00307 conn[4] = this->node(3)+1; 00308 conn[5] = this->node(3)+1; 00309 conn[6] = this->node(3)+1; 00310 conn[7] = this->node(3)+1; 00311 00312 return; 00313 00314 // Linear sub-tet 4 00315 case 4: 00316 00317 conn[0] = this->node(4)+1; 00318 conn[1] = this->node(8)+1; 00319 conn[2] = this->node(6)+1; 00320 conn[3] = this->node(6)+1; 00321 conn[4] = this->node(7)+1; 00322 conn[5] = this->node(7)+1; 00323 conn[6] = this->node(7)+1; 00324 conn[7] = this->node(7)+1; 00325 00326 return; 00327 00328 // Linear sub-tet 5 00329 case 5: 00330 00331 conn[0] = this->node(4)+1; 00332 conn[1] = this->node(5)+1; 00333 conn[2] = this->node(6)+1; 00334 conn[3] = this->node(6)+1; 00335 conn[4] = this->node(8)+1; 00336 conn[5] = this->node(8)+1; 00337 conn[6] = this->node(8)+1; 00338 conn[7] = this->node(8)+1; 00339 00340 return; 00341 00342 // Linear sub-tet 6 00343 case 6: 00344 00345 conn[0] = this->node(5)+1; 00346 conn[1] = this->node(9)+1; 00347 conn[2] = this->node(6)+1; 00348 conn[3] = this->node(6)+1; 00349 conn[4] = this->node(8)+1; 00350 conn[5] = this->node(8)+1; 00351 conn[6] = this->node(8)+1; 00352 conn[7] = this->node(8)+1; 00353 00354 return; 00355 00356 // Linear sub-tet 7 00357 case 7: 00358 00359 conn[0] = this->node(7)+1; 00360 conn[1] = this->node(6)+1; 00361 conn[2] = this->node(9)+1; 00362 conn[3] = this->node(9)+1; 00363 conn[4] = this->node(8)+1; 00364 conn[5] = this->node(8)+1; 00365 conn[6] = this->node(8)+1; 00366 conn[7] = this->node(8)+1; 00367 00368 return; 00369 00370 00371 default: 00372 00373 libmesh_error(); 00374 } 00375 } 00376 00377 case VTK: 00378 { 00379 conn.resize(10); 00380 conn[0] = this->node(0); 00381 conn[1] = this->node(1); 00382 conn[2] = this->node(2); 00383 conn[3] = this->node(3); 00384 conn[4] = this->node(4); 00385 conn[5] = this->node(5); 00386 conn[6] = this->node(6); 00387 conn[7] = this->node(7); 00388 conn[8] = this->node(8); 00389 conn[9] = this->node(9); 00390 return; 00391 /* 00392 conn.resize(4); 00393 switch (sc) 00394 { 00395 // Linear sub-tet 0 00396 case 0: 00397 00398 conn[0] = this->node(0); 00399 conn[1] = this->node(4); 00400 conn[2] = this->node(6); 00401 conn[3] = this->node(7); 00402 00403 return; 00404 00405 // Linear sub-tet 1 00406 case 1: 00407 00408 conn[0] = this->node(4); 00409 conn[1] = this->node(1); 00410 conn[2] = this->node(5); 00411 conn[3] = this->node(8); 00412 00413 return; 00414 00415 // Linear sub-tet 2 00416 case 2: 00417 00418 conn[0] = this->node(5); 00419 conn[1] = this->node(2); 00420 conn[2] = this->node(6); 00421 conn[3] = this->node(9); 00422 00423 return; 00424 00425 // Linear sub-tet 3 00426 case 3: 00427 00428 conn[0] = this->node(7); 00429 conn[1] = this->node(8); 00430 conn[2] = this->node(9); 00431 conn[3] = this->node(3); 00432 00433 return; 00434 00435 // Linear sub-tet 4 00436 case 4: 00437 00438 conn[0] = this->node(4); 00439 conn[1] = this->node(8); 00440 conn[2] = this->node(6); 00441 conn[3] = this->node(7); 00442 00443 return; 00444 00445 // Linear sub-tet 5 00446 case 5: 00447 00448 conn[0] = this->node(4); 00449 conn[1] = this->node(5); 00450 conn[2] = this->node(6); 00451 conn[3] = this->node(8); 00452 00453 return; 00454 00455 // Linear sub-tet 6 00456 case 6: 00457 00458 conn[0] = this->node(5); 00459 conn[1] = this->node(9); 00460 conn[2] = this->node(6); 00461 conn[3] = this->node(8); 00462 00463 return; 00464 00465 // Linear sub-tet 7 00466 case 7: 00467 00468 conn[0] = this->node(7); 00469 conn[1] = this->node(6); 00470 conn[2] = this->node(9); 00471 conn[3] = this->node(8); 00472 00473 return; 00474 00475 00476 default: 00477 00478 libmesh_error(); 00479 } 00480 */ 00481 } 00482 00483 default: 00484 libmesh_error(); 00485 } 00486 00487 libmesh_error(); 00488 } 00489 00490 00491 00492 const unsigned short int Tet10::_second_order_vertex_child_number[10] = 00493 { 00494 99,99,99,99, // Vertices 00495 0,1,0,0,1,2 // Edges 00496 }; 00497 00498 00499 00500 const unsigned short int Tet10::_second_order_vertex_child_index[10] = 00501 { 00502 99,99,99,99, // Vertices 00503 1,2,2,3,3,3 // Edges 00504 }; 00505 00506 00507 00508 std::pair<unsigned short int, unsigned short int> 00509 Tet10::second_order_child_vertex (const unsigned int n) const 00510 { 00511 libmesh_assert_greater_equal (n, this->n_vertices()); 00512 libmesh_assert_less (n, this->n_nodes()); 00513 return std::pair<unsigned short int, unsigned short int> 00514 (_second_order_vertex_child_number[n], 00515 _second_order_vertex_child_index[n]); 00516 } 00517 00518 00519 00520 unsigned short int Tet10::second_order_adjacent_vertex (const unsigned int n, 00521 const unsigned int v) const 00522 { 00523 libmesh_assert_greater_equal (n, this->n_vertices()); 00524 libmesh_assert_less (n, this->n_nodes()); 00525 libmesh_assert_less (v, 2); 00526 return _second_order_adjacent_vertices[n-this->n_vertices()][v]; 00527 } 00528 00529 00530 00531 const unsigned short int Tet10::_second_order_adjacent_vertices[6][2] = 00532 { 00533 {0, 1}, // vertices adjacent to node 4 00534 {1, 2}, // vertices adjacent to node 5 00535 {0, 2}, // vertices adjacent to node 6 00536 {0, 3}, // vertices adjacent to node 7 00537 {1, 3}, // vertices adjacent to node 8 00538 {2, 3} // vertices adjacent to node 9 00539 }; 00540 00541 00542 00543 00544 00545 #ifdef LIBMESH_ENABLE_AMR 00546 00547 const float Tet10::_embedding_matrix[8][10][10] = 00548 { 00549 // embedding matrix for child 0 00550 { 00551 // 0 1 2 3 4 5 6 7 8 9 00552 { 1., 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 0 00553 { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 1 00554 { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 2 00555 { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0.}, // 3 00556 { 0.375,-0.125, 0., 0., 0.75, 0., 0., 0., 0., 0.}, // 4 00557 { 0.,-0.125,-0.125, 0., 0.5, 0.25, 0.5, 0., 0., 0.}, // 5 00558 { 0.375, 0.,-0.125, 0., 0., 0., 0.75, 0., 0., 0.}, // 6 00559 { 0.375, 0., 0.,-0.125, 0., 0., 0., 0.75, 0., 0.}, // 7 00560 { 0.,-0.125, 0.,-0.125, 0.5, 0., 0., 0.5, 0.25, 0.}, // 8 00561 { 0., 0.,-0.125,-0.125, 0., 0., 0.5, 0.5, 0., 0.25} // 9 00562 }, 00563 00564 // embedding matrix for child 1 00565 { 00566 // 0 1 2 3 4 5 6 7 8 9 00567 { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 0 00568 { 0., 1., 0., 0., 0., 0., 0., 0., 0., 0.}, // 1 00569 { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 2 00570 { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 3 00571 {-0.125, 0.375, 0., 0., 0.75, 0., 0., 0., 0., 0.}, // 4 00572 { 0., 0.375,-0.125, 0., 0., 0.75, 0., 0., 0., 0.}, // 5 00573 {-0.125, 0.,-0.125, 0., 0.5, 0.5, 0.25, 0., 0., 0.}, // 6 00574 {-0.125, 0., 0.,-0.125, 0.5, 0., 0., 0.25, 0.5, 0.}, // 7 00575 { 0., 0.375, 0.,-0.125, 0., 0., 0., 0., 0.75, 0.}, // 8 00576 { 0., 0.,-0.125,-0.125, 0., 0.5, 0., 0., 0.5, 0.25} // 9 00577 }, 00578 00579 // embedding matrix for child 2 00580 { 00581 // 0 1 2 3 4 5 6 7 8 9 00582 { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 0 00583 { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 1 00584 { 0., 0., 1., 0., 0., 0., 0., 0., 0., 0.}, // 2 00585 { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1.}, // 3 00586 {-0.125,-0.125, 0., 0., 0.25, 0.5, 0.5, 0., 0., 0.}, // 4 00587 { 0.,-0.125, 0.375, 0., 0., 0.75, 0., 0., 0., 0.}, // 5 00588 {-0.125, 0., 0.375, 0., 0., 0., 0.75, 0., 0., 0.}, // 6 00589 {-0.125, 0., 0.,-0.125, 0., 0., 0.5, 0.25, 0., 0.5}, // 7 00590 { 0.,-0.125, 0.,-0.125, 0., 0.5, 0., 0., 0.25, 0.5}, // 8 00591 { 0., 0., 0.375,-0.125, 0., 0., 0., 0., 0., 0.75} // 9 00592 }, 00593 00594 // embedding matrix for child 3 00595 { 00596 // 0 1 2 3 4 5 6 7 8 9 00597 { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0.}, // 0 00598 { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 1 00599 { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1.}, // 2 00600 { 0., 0., 0., 1., 0., 0., 0., 0., 0., 0.}, // 3 00601 {-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5, 0.}, // 4 00602 { 0.,-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5}, // 5 00603 {-0.125, 0.,-0.125, 0., 0., 0., 0.25, 0.5, 0., 0.5}, // 6 00604 {-0.125, 0., 0., 0.375, 0., 0., 0., 0.75, 0., 0.}, // 7 00605 { 0.,-0.125, 0., 0.375, 0., 0., 0., 0., 0.75, 0.}, // 8 00606 { 0., 0.,-0.125, 0.375, 0., 0., 0., 0., 0., 0.75} // 9 00607 }, 00608 00609 // embedding matrix for child 4 00610 { 00611 // 0 1 2 3 4 5 6 7 8 9 00612 { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 0 00613 { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 1 00614 { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 2 00615 { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0.}, // 3 00616 {-0.125, 0., 0.,-0.125, 0.5, 0., 0., 0.25, 0.5, 0.}, // 4 00617 {-0.125,-0.125,-0.125,-0.125, 0.25, 0.25, 0.25, 0.25, 0.25, 0.25}, // 5 00618 { 0.,-0.125,-0.125, 0., 0.5, 0.25, 0.5, 0., 0., 0.}, // 6 00619 { 0.,-0.125, 0.,-0.125, 0.5, 0., 0., 0.5, 0.25, 0.}, // 7 00620 {-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5, 0.}, // 8 00621 { 0., 0.,-0.125,-0.125, 0., 0., 0.5, 0.5, 0., 0.25} // 9 00622 }, 00623 00624 // embedding matrix for child 5 00625 { 00626 // 0 1 2 3 4 5 6 7 8 9 00627 { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 0 00628 { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 1 00629 { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 2 00630 { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 3 00631 {-0.125, 0.,-0.125, 0., 0.5, 0.5, 0.25, 0., 0., 0.}, // 4 00632 {-0.125,-0.125, 0., 0., 0.25, 0.5, 0.5, 0., 0., 0.}, // 5 00633 { 0.,-0.125,-0.125, 0., 0.5, 0.25, 0.5, 0., 0., 0.}, // 6 00634 {-0.125, 0., 0.,-0.125, 0.5, 0., 0., 0.25, 0.5, 0.}, // 7 00635 { 0., 0.,-0.125,-0.125, 0., 0.5, 0., 0., 0.5, 0.25}, // 8 00636 {-0.125,-0.125,-0.125,-0.125, 0.25, 0.25, 0.25, 0.25, 0.25, 0.25} // 9 00637 }, 00638 00639 // embedding matrix for child 6 00640 { 00641 // 0 1 2 3 4 5 6 7 8 9 00642 { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 0 00643 { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 1 00644 { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1.}, // 2 00645 { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 3 00646 {-0.125,-0.125, 0., 0., 0.25, 0.5, 0.5, 0., 0., 0.}, // 4 00647 { 0.,-0.125, 0.,-0.125, 0., 0.5, 0., 0., 0.25, 0.5}, // 5 00648 {-0.125, 0., 0.,-0.125, 0., 0., 0.5, 0.25, 0., 0.5}, // 6 00649 {-0.125,-0.125,-0.125,-0.125, 0.25, 0.25, 0.25, 0.25, 0.25, 0.25}, // 7 00650 { 0., 0.,-0.125,-0.125, 0., 0.5, 0., 0., 0.5, 0.25}, // 8 00651 { 0.,-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5} // 9 00652 }, 00653 00654 // embedding matrix for child 7 00655 { 00656 // 0 1 2 3 4 5 6 7 8 9 00657 { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 0 00658 { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 1 00659 { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1.}, // 2 00660 { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0.}, // 3 00661 {-0.125,-0.125,-0.125,-0.125, 0.25, 0.25, 0.25, 0.25, 0.25, 0.25}, // 4 00662 { 0.,-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5}, // 5 00663 {-0.125, 0., 0.,-0.125, 0., 0., 0.5, 0.25, 0., 0.5}, // 6 00664 { 0., 0.,-0.125,-0.125, 0., 0., 0.5, 0.5, 0., 0.25}, // 7 00665 {-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5, 0.}, // 8 00666 {-0.125, 0.,-0.125, 0., 0., 0., 0.25, 0.5, 0., 0.5} // 9 00667 } 00668 }; 00669 00670 00671 00672 00673 00674 00675 float Tet10::embedding_matrix (const unsigned int i, 00676 const unsigned int j, 00677 const unsigned int k) const 00678 { 00679 // Choose an optimal diagonal, if one has not already been selected 00680 this->choose_diagonal(); 00681 00682 // Permuted j and k indices 00683 unsigned int 00684 jp=j, 00685 kp=k; 00686 00687 if ((i>3) && (this->_diagonal_selection!=DIAG_02_13)) 00688 { 00689 // Just the enum value cast to an unsigned int... 00690 const unsigned ds = static_cast<unsigned int>(this->_diagonal_selection); // == 1 or 2 00691 00692 // Instead of doing a lot of arithmetic, use these 00693 // straightforward arrays for the permutations. Note that 3 -> 00694 // 3, and the first array consists of "forward" permutations of 00695 // the sets {0,1,2}, {4,5,6}, and {7,8,9} while the second array 00696 // consists of "reverse" permutations of the same sets. 00697 const unsigned int perms[2][10] = 00698 { 00699 {1, 2, 0, 3, 5, 6, 4, 8, 9, 7}, 00700 {2, 0, 1, 3, 6, 4, 5, 9, 7, 8} 00701 }; 00702 00703 // Permute j 00704 jp = perms[ds-1][j]; 00705 // if (jp<3) 00706 // jp = (jp+ds)%3; 00707 // else if (jp>3) 00708 // jp = (jp-1+ds)%3 + 1 + 3*((jp-1)/3); 00709 00710 // Permute k 00711 kp = perms[ds-1][k]; 00712 // if (kp<3) 00713 // kp = (kp+ds)%3; 00714 // else if (kp>3) 00715 // kp = (kp-1+ds)%3 + 1 + 3*((kp-1)/3); 00716 } 00717 00718 // Debugging: 00719 // libMesh::err << "Selected diagonal " << _diagonal_selection << std::endl; 00720 // libMesh::err << "j=" << j << std::endl; 00721 // libMesh::err << "k=" << k << std::endl; 00722 // libMesh::err << "jp=" << jp << std::endl; 00723 // libMesh::err << "kp=" << kp << std::endl; 00724 00725 // Call embedding matrx with permuted indices 00726 return this->_embedding_matrix[i][jp][kp]; 00727 } 00728 00729 #endif // #ifdef LIBMESH_ENABLE_AMR 00730 00731 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:45 UTC
Hosted By: