cell_tet4.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_tet4.h" 00024 #include "libmesh/edge_edge2.h" 00025 #include "libmesh/face_tri3.h" 00026 00027 namespace libMesh 00028 { 00029 00030 00031 00032 // ------------------------------------------------------------ 00033 // Tet4 class static member initializations 00034 const unsigned int Tet4::side_nodes_map[4][3] = 00035 { 00036 {0, 2, 1}, // Side 0 00037 {0, 1, 3}, // Side 1 00038 {1, 2, 3}, // Side 2 00039 {2, 0, 3} // Side 3 00040 }; 00041 00042 const unsigned int Tet4::edge_nodes_map[6][2] = 00043 { 00044 {0, 1}, // Side 0 00045 {1, 2}, // Side 1 00046 {0, 2}, // Side 2 00047 {0, 3}, // Side 3 00048 {1, 3}, // Side 4 00049 {2, 3} // Side 5 00050 }; 00051 00052 00053 // ------------------------------------------------------------ 00054 // Tet4 class member functions 00055 00056 bool Tet4::is_vertex(const unsigned int) const 00057 { 00058 return true; 00059 } 00060 00061 bool Tet4::is_edge(const unsigned int) const 00062 { 00063 return false; 00064 } 00065 00066 bool Tet4::is_face(const unsigned int) const 00067 { 00068 return false; 00069 } 00070 00071 bool Tet4::is_node_on_edge(const unsigned int n, 00072 const unsigned int e) const 00073 { 00074 libmesh_assert_less (e, n_edges()); 00075 for (unsigned int i = 0; i != 2; ++i) 00076 if (edge_nodes_map[e][i] == n) 00077 return true; 00078 return false; 00079 } 00080 00081 00082 00083 00084 #ifdef LIBMESH_ENABLE_AMR 00085 00086 // This function only works if LIBMESH_ENABLE_AMR... 00087 bool Tet4::is_child_on_side(const unsigned int c, 00088 const unsigned int s) const 00089 { 00090 // OK, for the Tet4, this is pretty obvious... it is sets of nodes 00091 // not equal to the current node. But if we want this algorithm to 00092 // be generic and work for Tet10 also it helps to do it this way. 00093 const unsigned int nodes_opposite[4][3] = 00094 { 00095 {1,2,3}, // nodes opposite node 0 00096 {0,2,3}, // nodes opposite node 1 00097 {0,1,3}, // nodes opposite node 2 00098 {0,1,2} // nodes opposite node 3 00099 }; 00100 00101 // Call the base class helper function 00102 return Tet::is_child_on_side_helper(c, s, nodes_opposite); 00103 } 00104 00105 #else 00106 00107 bool Tet4::is_child_on_side(const unsigned int /*c*/, 00108 const unsigned int /*s*/) const 00109 { 00110 libmesh_not_implemented(); 00111 return false; 00112 } 00113 00114 #endif //LIBMESH_ENABLE_AMR 00115 00116 00117 00118 00119 bool Tet4::is_node_on_side(const unsigned int n, 00120 const unsigned int s) const 00121 { 00122 libmesh_assert_less (s, n_sides()); 00123 for (unsigned int i = 0; i != 3; ++i) 00124 if (side_nodes_map[s][i] == n) 00125 return true; 00126 return false; 00127 } 00128 00129 AutoPtr<Elem> Tet4::build_side (const unsigned int i, 00130 bool proxy) const 00131 { 00132 libmesh_assert_less (i, this->n_sides()); 00133 00134 if (proxy) 00135 { 00136 AutoPtr<Elem> ap(new Side<Tri3,Tet4>(this,i)); 00137 return ap; 00138 } 00139 00140 else 00141 { 00142 AutoPtr<Elem> face(new Tri3); 00143 00144 switch (i) 00145 { 00146 case 0: 00147 { 00148 face->set_node(0) = this->get_node(0); 00149 face->set_node(1) = this->get_node(2); 00150 face->set_node(2) = this->get_node(1); 00151 00152 return face; 00153 } 00154 case 1: 00155 { 00156 face->set_node(0) = this->get_node(0); 00157 face->set_node(1) = this->get_node(1); 00158 face->set_node(2) = this->get_node(3); 00159 00160 return face; 00161 } 00162 case 2: 00163 { 00164 face->set_node(0) = this->get_node(1); 00165 face->set_node(1) = this->get_node(2); 00166 face->set_node(2) = this->get_node(3); 00167 00168 return face; 00169 } 00170 case 3: 00171 { 00172 face->set_node(0) = this->get_node(2); 00173 face->set_node(1) = this->get_node(0); 00174 face->set_node(2) = this->get_node(3); 00175 00176 return face; 00177 } 00178 default: 00179 { 00180 libmesh_error(); 00181 } 00182 } 00183 } 00184 00185 // We'll never get here. 00186 libmesh_error(); 00187 AutoPtr<Elem> ap(NULL); return ap; 00188 } 00189 00190 00191 AutoPtr<Elem> Tet4::build_edge (const unsigned int i) const 00192 { 00193 libmesh_assert_less (i, this->n_edges()); 00194 00195 return AutoPtr<Elem>(new SideEdge<Edge2,Tet4>(this,i)); 00196 } 00197 00198 00199 void Tet4::connectivity(const unsigned int libmesh_dbg_var(sc), 00200 const IOPackage iop, 00201 std::vector<dof_id_type>& conn) const 00202 { 00203 libmesh_assert(_nodes); 00204 libmesh_assert_less (sc, this->n_sub_elem()); 00205 libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE); 00206 00207 00208 switch (iop) 00209 { 00210 case TECPLOT: 00211 { 00212 conn.resize(8); 00213 conn[0] = this->node(0)+1; 00214 conn[1] = this->node(1)+1; 00215 conn[2] = this->node(2)+1; 00216 conn[3] = this->node(2)+1; 00217 conn[4] = this->node(3)+1; 00218 conn[5] = this->node(3)+1; 00219 conn[6] = this->node(3)+1; 00220 conn[7] = this->node(3)+1; 00221 return; 00222 } 00223 00224 case VTK: 00225 { 00226 conn.resize(4); 00227 conn[0] = this->node(0); 00228 conn[1] = this->node(1); 00229 conn[2] = this->node(2); 00230 conn[3] = this->node(3); 00231 return; 00232 } 00233 00234 default: 00235 libmesh_error(); 00236 } 00237 00238 libmesh_error(); 00239 } 00240 00241 00242 00243 #ifdef LIBMESH_ENABLE_AMR 00244 00245 const float Tet4::_embedding_matrix[8][4][4] = 00246 { 00247 // embedding matrix for child 0 00248 { 00249 // 0 1 2 3 00250 {1.0, 0.0, 0.0, 0.0}, // 0 00251 {0.5, 0.5, 0.0, 0.0}, // 1 00252 {0.5, 0.0, 0.5, 0.0}, // 2 00253 {0.5, 0.0, 0.0, 0.5} // 3 00254 }, 00255 00256 // embedding matrix for child 1 00257 { 00258 // 0 1 2 3 00259 {0.5, 0.5, 0.0, 0.0}, // 0 00260 {0.0, 1.0, 0.0, 0.0}, // 1 00261 {0.0, 0.5, 0.5, 0.0}, // 2 00262 {0.0, 0.5, 0.0, 0.5} // 3 00263 }, 00264 00265 // embedding matrix for child 2 00266 { 00267 // 0 1 2 3 00268 {0.5, 0.0, 0.5, 0.0}, // 0 00269 {0.0, 0.5, 0.5, 0.0}, // 1 00270 {0.0, 0.0, 1.0, 0.0}, // 2 00271 {0.0, 0.0, 0.5, 0.5} // 3 00272 }, 00273 00274 // embedding matrix for child 3 00275 { 00276 // 0 1 2 3 00277 {0.5, 0.0, 0.0, 0.5}, // 0 00278 {0.0, 0.5, 0.0, 0.5}, // 1 00279 {0.0, 0.0, 0.5, 0.5}, // 2 00280 {0.0, 0.0, 0.0, 1.0} // 3 00281 }, 00282 00283 // embedding matrix for child 4 00284 { 00285 // 0 1 2 3 00286 {0.5, 0.5, 0.0, 0.0}, // 0 00287 {0.0, 0.5, 0.0, 0.5}, // 1 00288 {0.5, 0.0, 0.5, 0.0}, // 2 00289 {0.5, 0.0, 0.0, 0.5} // 3 00290 }, 00291 00292 // embedding matrix for child 5 00293 { 00294 // 0 1 2 3 00295 {0.5, 0.5, 0.0, 0.0}, // 0 00296 {0.0, 0.5, 0.5, 0.0}, // 1 00297 {0.5, 0.0, 0.5, 0.0}, // 2 00298 {0.0, 0.5, 0.0, 0.5} // 3 00299 }, 00300 00301 // embedding matrix for child 6 00302 { 00303 // 0 1 2 3 00304 {0.5, 0.0, 0.5, 0.0}, // 0 00305 {0.0, 0.5, 0.5, 0.0}, // 1 00306 {0.0, 0.0, 0.5, 0.5}, // 2 00307 {0.0, 0.5, 0.0, 0.5} // 3 00308 }, 00309 00310 // embedding matrix for child 7 00311 { 00312 // 0 1 2 3 00313 {0.5, 0.0, 0.5, 0.0}, // 0 00314 {0.0, 0.5, 0.0, 0.5}, // 1 00315 {0.0, 0.0, 0.5, 0.5}, // 2 00316 {0.5, 0.0, 0.0, 0.5} // 3 00317 } 00318 }; 00319 00320 #endif // #ifdef LIBMESH_ENABLE_AMR 00321 00322 00323 00324 00325 00326 Real Tet4::volume () const 00327 { 00328 // The volume of a tetrahedron is 1/6 the box product formed 00329 // by its base and apex vectors 00330 Point a ( *this->get_node(3) - *this->get_node(0) ); 00331 00332 // b is the vector pointing from 0 to 1 00333 Point b ( *this->get_node(1) - *this->get_node(0) ); 00334 00335 // c is the vector pointing from 0 to 2 00336 Point c ( *this->get_node(2) - *this->get_node(0) ); 00337 00338 return (1.0 / 6.0) * (a * (b.cross(c))); 00339 } 00340 00341 00342 00343 00344 std::pair<Real, Real> Tet4::min_and_max_angle() const 00345 { 00346 Point n[4]; 00347 00348 // Compute the outward normal vectors on each face 00349 n[0] = (this->point(2) - this->point(0)).cross(this->point(1) - this->point(0)); 00350 n[1] = (this->point(1) - this->point(0)).cross(this->point(3) - this->point(0)); 00351 n[2] = (this->point(2) - this->point(1)).cross(this->point(3) - this->point(1)); 00352 n[3] = (this->point(0) - this->point(2)).cross(this->point(3) - this->point(2)); 00353 00354 Real dihedral_angles[6]; // 01, 02, 03, 12, 13, 23 00355 00356 // Compute dihedral angles 00357 for (unsigned int k=0,i=0; i<4; ++i) 00358 for (unsigned int j=i+1; j<4; ++j,k+=1) 00359 dihedral_angles[k] = std::acos(n[i]*n[j] / n[i].size() / n[j].size()); // return value is between 0 and PI 00360 00361 // Return max/min dihedral angles 00362 return std::make_pair(*std::min_element(dihedral_angles, dihedral_angles+6), 00363 *std::max_element(dihedral_angles, dihedral_angles+6)); 00364 00365 } 00366 00367 00368 00369 #ifdef LIBMESH_ENABLE_AMR 00370 float Tet4::embedding_matrix (const unsigned int i, 00371 const unsigned int j, 00372 const unsigned int k) const 00373 { 00374 // Choose an optimal diagonal, if one has not already been selected 00375 this->choose_diagonal(); 00376 00377 // Permuted j and k indices 00378 unsigned int 00379 jp=j, 00380 kp=k; 00381 00382 if ((i>3) && (this->_diagonal_selection!=DIAG_02_13)) 00383 { 00384 // Just the enum value cast to an unsigned int... 00385 const unsigned ds = static_cast<unsigned int>(this->_diagonal_selection); 00386 00387 // Permute j, k: 00388 // ds==1 ds==2 00389 // 0 -> 1 0 -> 2 00390 // 1 -> 2 1 -> 0 00391 // 2 -> 0 2 -> 1 00392 if (jp != 3) 00393 jp = (jp+ds)%3; 00394 00395 if (kp != 3) 00396 kp = (kp+ds)%3; 00397 } 00398 00399 // Debugging 00400 // libMesh::err << "Selected diagonal " << _diagonal_selection << std::endl; 00401 // libMesh::err << "j=" << j << std::endl; 00402 // libMesh::err << "k=" << k << std::endl; 00403 // libMesh::err << "jp=" << jp << std::endl; 00404 // libMesh::err << "kp=" << kp << std::endl; 00405 00406 // Call embedding matrx with permuted indices 00407 return this->_embedding_matrix[i][jp][kp]; 00408 } 00409 00410 00411 00412 00413 00414 00415 // void Tet4::reselect_diagonal (const Diagonal diag) 00416 // { 00417 // /* Make sure that the element has just been refined. */ 00418 // libmesh_assert(_children); 00419 // libmesh_assert_equal_to (n_children(), 8); 00420 // libmesh_assert_equal_to (_children[0]->refinement_flag(), JUST_REFINED); 00421 // libmesh_assert_equal_to (_children[1]->refinement_flag(), JUST_REFINED); 00422 // libmesh_assert_equal_to (_children[2]->refinement_flag(), JUST_REFINED); 00423 // libmesh_assert_equal_to (_children[3]->refinement_flag(), JUST_REFINED); 00424 // libmesh_assert_equal_to (_children[4]->refinement_flag(), JUST_REFINED); 00425 // libmesh_assert_equal_to (_children[5]->refinement_flag(), JUST_REFINED); 00426 // libmesh_assert_equal_to (_children[6]->refinement_flag(), JUST_REFINED); 00427 // libmesh_assert_equal_to (_children[7]->refinement_flag(), JUST_REFINED); 00428 // 00429 // /* Check whether anything has to be changed. */ 00430 // if (_diagonal_selection!=diag) 00431 // { 00432 // /* Set new diagonal selection. */ 00433 // _diagonal_selection = diag; 00434 // 00435 // /* The first four children do not have to be changed. For the 00436 // others, only the nodes have to be changed. Note that we have 00437 // to keep track of the nodes ourselves since there is no \p 00438 // MeshRefinement object with a valid \p _new_nodes_map 00439 // available. */ 00440 // for (unsigned int c=4; c<this->n_children(); c++) 00441 // { 00442 // Elem *child = this->child(c); 00443 // for (unsigned int nc=0; nc<child->n_nodes(); nc++) 00444 // { 00445 // /* Unassign the current node. */ 00446 // child->set_node(nc) = NULL; 00447 // 00448 // /* We have to find the correct new node now. We know 00449 // that it exists somewhere. We make use of the fact 00450 // that the embedding matrix for these children consists 00451 // of entries 0.0 and 0.5 only. Also, we make use of 00452 // the properties of the embedding matrix for the first 00453 // (unchanged) four children, which allow us to use a 00454 // simple mechanism to find the required node. */ 00455 // 00456 // 00457 // unsigned int first_05_in_embedding_matrix = libMesh::invalid_uint; 00458 // for (unsigned int n=0; n<this->n_nodes(); n++) 00459 // { 00460 // if (this->embedding_matrix(c,nc,n) != 0.0) 00461 // { 00462 // /* It must be 0.5 then. Check whether it's the 00463 // first or second time that we get a 0.5 00464 // value. */ 00465 // if (first_05_in_embedding_matrix==libMesh::invalid_uint) 00466 // { 00467 // /* First time, so just remeber this position. */ 00468 // first_05_in_embedding_matrix = n; 00469 // } 00470 // else 00471 // { 00472 // /* Second time, so we know now which node to 00473 // use. */ 00474 // child->set_node(nc) = this->child(n)->get_node(first_05_in_embedding_matrix); 00475 // } 00476 // 00477 // } 00478 // } 00479 // 00480 // /* Make sure that a node has been found. */ 00481 // libmesh_assert(child->get_node(nc)); 00482 // } 00483 // } 00484 // } 00485 // } 00486 00487 00488 00489 // void Tet4::reselect_optimal_diagonal (const Diagonal exclude_this) 00490 // { 00491 // Real diag_01_23 = (this->point(0)+this->point(1)-this->point(2)-this->point(3)).size_sq(); 00492 // Real diag_02_13 = (this->point(0)-this->point(1)+this->point(2)-this->point(3)).size_sq(); 00493 // Real diag_03_12 = (this->point(0)-this->point(1)-this->point(2)+this->point(3)).size_sq(); 00494 // 00495 // Diagonal use_this = INVALID_DIAG; 00496 // switch (exclude_this) 00497 // { 00498 // case DIAG_01_23: 00499 // use_this = DIAG_02_13; 00500 // if (diag_03_12 < diag_02_13) 00501 // { 00502 // use_this = DIAG_03_12; 00503 // } 00504 // break; 00505 // 00506 // case DIAG_02_13: 00507 // use_this = DIAG_03_12; 00508 // if (diag_01_23 < diag_03_12) 00509 // { 00510 // use_this = DIAG_01_23; 00511 // } 00512 // break; 00513 // 00514 // case DIAG_03_12: 00515 // use_this = DIAG_02_13; 00516 // if (diag_01_23 < diag_02_13) 00517 // { 00518 // use_this = DIAG_01_23; 00519 // } 00520 // break; 00521 // 00522 // default: 00523 // use_this = DIAG_02_13; 00524 // if (diag_01_23 < diag_02_13 || diag_03_12 < diag_02_13) 00525 // { 00526 // if (diag_01_23 < diag_03_12) 00527 // { 00528 // use_this = DIAG_01_23; 00529 // } 00530 // else 00531 // { 00532 // use_this = DIAG_03_12; 00533 // } 00534 // } 00535 // break; 00536 // } 00537 // 00538 // reselect_diagonal (use_this); 00539 // } 00540 #endif // #ifdef LIBMESH_ENABLE_AMR 00541 00542 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:45 UTC
Hosted By: