$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 // Local includes 00020 #include "libmesh/libmesh_config.h" 00021 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00022 #include "libmesh/inf_fe.h" 00023 #include "libmesh/inf_fe_macro.h" 00024 #include "libmesh/fe.h" 00025 #include "libmesh/fe_interface.h" 00026 #include "libmesh/fe_compute_data.h" 00027 #include "libmesh/elem.h" 00028 00029 namespace libMesh 00030 { 00031 00032 00033 // ------------------------------------------------------------ 00034 // InfFE class static member initialization 00035 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00036 ElemType InfFE<Dim,T_radial,T_map>::_compute_node_indices_fast_current_elem_type = INVALID_ELEM; 00037 00038 #ifdef DEBUG 00039 00040 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00041 bool InfFE<Dim,T_radial,T_map>::_warned_for_nodal_soln = false; 00042 00043 00044 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00045 bool InfFE<Dim,T_radial,T_map>::_warned_for_shape = false; 00046 00047 #endif 00048 00049 00050 00051 00052 // ------------------------------------------------------------ 00053 // InfFE static class members 00054 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00055 unsigned int InfFE<Dim,T_radial,T_map>::n_dofs (const FEType& fet, 00056 const ElemType inf_elem_type) 00057 { 00058 const ElemType base_et (Base::get_elem_type(inf_elem_type)); 00059 00060 if (Dim > 1) 00061 return FEInterface::n_dofs(Dim-1, fet, base_et) * Radial::n_dofs(fet.radial_order); 00062 else 00063 return Radial::n_dofs(fet.radial_order); 00064 } 00065 00066 00067 00068 00069 00070 00071 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00072 unsigned int InfFE<Dim,T_radial,T_map>::n_dofs_at_node (const FEType& fet, 00073 const ElemType inf_elem_type, 00074 const unsigned int n) 00075 { 00076 const ElemType base_et (Base::get_elem_type(inf_elem_type)); 00077 00078 unsigned int n_base, n_radial; 00079 compute_node_indices(inf_elem_type, n, n_base, n_radial); 00080 00081 // libMesh::out << "elem_type=" << inf_elem_type 00082 // << ", fet.radial_order=" << fet.radial_order 00083 // << ", n=" << n 00084 // << ", n_radial=" << n_radial 00085 // << ", n_base=" << n_base 00086 // << std::endl; 00087 00088 if (Dim > 1) 00089 return FEInterface::n_dofs_at_node(Dim-1, fet, base_et, n_base) 00090 * Radial::n_dofs_at_node(fet.radial_order, n_radial); 00091 else 00092 return Radial::n_dofs_at_node(fet.radial_order, n_radial); 00093 } 00094 00095 00096 00097 00098 00099 00100 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00101 unsigned int InfFE<Dim,T_radial,T_map>::n_dofs_per_elem (const FEType& fet, 00102 const ElemType inf_elem_type) 00103 { 00104 const ElemType base_et (Base::get_elem_type(inf_elem_type)); 00105 00106 if (Dim > 1) 00107 return FEInterface::n_dofs_per_elem(Dim-1, fet, base_et) 00108 * Radial::n_dofs_per_elem(fet.radial_order); 00109 else 00110 return Radial::n_dofs_per_elem(fet.radial_order); 00111 } 00112 00113 00114 00115 00116 00117 00118 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00119 void InfFE<Dim,T_radial,T_map>::nodal_soln(const FEType& /* fet */, 00120 const Elem* /* elem */, 00121 const std::vector<Number>& /* elem_soln */, 00122 std::vector<Number>& nodal_soln) 00123 { 00124 #ifdef DEBUG 00125 if (!_warned_for_nodal_soln) 00126 { 00127 libMesh::err << "WARNING: nodal_soln(...) does _not_ work for infinite elements." << std::endl 00128 << " Will return an empty nodal solution. Use " << std::endl 00129 << " InfFE<Dim,T_radial,T_map>::compute_data(..) instead!" << std::endl; 00130 _warned_for_nodal_soln = true; 00131 } 00132 #endif 00133 00134 /* 00135 * In the base the infinite element couples to 00136 * conventional finite elements. To not destroy 00137 * the values there, clear \p nodal_soln. This 00138 * indicates to the user of \p nodal_soln to 00139 * not use this result. 00140 */ 00141 nodal_soln.clear(); 00142 libmesh_assert (nodal_soln.empty()); 00143 return; 00144 } 00145 00146 00147 00148 00149 00150 00151 00152 00153 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00154 Real InfFE<Dim,T_radial,T_map>::shape(const FEType& fet, 00155 const ElemType inf_elem_type, 00156 const unsigned int i, 00157 const Point& p) 00158 { 00159 libmesh_assert_not_equal_to (Dim, 0); 00160 00161 #ifdef DEBUG 00162 // this makes only sense when used for mapping 00163 if ((T_radial != INFINITE_MAP) && !_warned_for_shape) 00164 { 00165 libMesh::err << "WARNING: InfFE<Dim,T_radial,T_map>::shape(...) does _not_" << std::endl 00166 << " return the correct trial function! Use " << std::endl 00167 << " InfFE<Dim,T_radial,T_map>::compute_data(..) instead!" 00168 << std::endl; 00169 _warned_for_shape = true; 00170 } 00171 #endif 00172 00173 const ElemType base_et (Base::get_elem_type(inf_elem_type)); 00174 const Order o_radial (fet.radial_order); 00175 const Real v (p(Dim-1)); 00176 00177 unsigned int i_base, i_radial; 00178 compute_shape_indices(fet, inf_elem_type, i, i_base, i_radial); 00179 00180 //TODO:[SP/DD] exp(ikr) is still missing here! 00181 if (Dim > 1) 00182 return FEInterface::shape(Dim-1, fet, base_et, i_base, p) 00183 * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial) 00184 * InfFE<Dim,T_radial,T_map>::Radial::decay(v); 00185 else 00186 return InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial) 00187 * InfFE<Dim,T_radial,T_map>::Radial::decay(v); 00188 } 00189 00190 00191 00192 00193 00194 00195 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00196 Real InfFE<Dim,T_radial,T_map>::shape(const FEType& fet, 00197 const Elem* inf_elem, 00198 const unsigned int i, 00199 const Point& p) 00200 { 00201 libmesh_assert(inf_elem); 00202 libmesh_assert_not_equal_to (Dim, 0); 00203 00204 #ifdef DEBUG 00205 // this makes only sense when used for mapping 00206 if ((T_radial != INFINITE_MAP) && !_warned_for_shape) 00207 { 00208 libMesh::err << "WARNING: InfFE<Dim,T_radial,T_map>::shape(...) does _not_" << std::endl 00209 << " return the correct trial function! Use " << std::endl 00210 << " InfFE<Dim,T_radial,T_map>::compute_data(..) instead!" 00211 << std::endl; 00212 _warned_for_shape = true; 00213 } 00214 #endif 00215 00216 const Order o_radial (fet.radial_order); 00217 const Real v (p(Dim-1)); 00218 UniquePtr<Elem> base_el (inf_elem->build_side(0)); 00219 00220 unsigned int i_base, i_radial; 00221 compute_shape_indices(fet, inf_elem->type(), i, i_base, i_radial); 00222 00223 if (Dim > 1) 00224 return FEInterface::shape(Dim-1, fet, base_el.get(), i_base, p) 00225 * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial) 00226 * InfFE<Dim,T_radial,T_map>::Radial::decay(v); 00227 else 00228 return InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial) 00229 * InfFE<Dim,T_radial,T_map>::Radial::decay(v); 00230 } 00231 00232 00233 00234 00235 00236 00237 00238 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00239 void InfFE<Dim,T_radial,T_map>::compute_data(const FEType& fet, 00240 const Elem* inf_elem, 00241 FEComputeData& data) 00242 { 00243 libmesh_assert(inf_elem); 00244 libmesh_assert_not_equal_to (Dim, 0); 00245 00246 const Order o_radial (fet.radial_order); 00247 const Order radial_mapping_order (Radial::mapping_order()); 00248 const Point& p (data.p); 00249 const Real v (p(Dim-1)); 00250 UniquePtr<Elem> base_el (inf_elem->build_side(0)); 00251 00252 /* 00253 * compute \p interpolated_dist containing the mapping-interpolated 00254 * distance of the base point to the origin. This is the same 00255 * for all shape functions. Set \p interpolated_dist to 0, it 00256 * is added to. 00257 */ 00258 Real interpolated_dist = 0.; 00259 switch (Dim) 00260 { 00261 case 1: 00262 { 00263 libmesh_assert_equal_to (inf_elem->type(), INFEDGE2); 00264 interpolated_dist = Point(inf_elem->point(0) - inf_elem->point(1)).size(); 00265 break; 00266 } 00267 00268 case 2: 00269 { 00270 const unsigned int n_base_nodes = base_el->n_nodes(); 00271 00272 const Point origin = inf_elem->origin(); 00273 const Order base_mapping_order (base_el->default_order()); 00274 const ElemType base_mapping_elem_type (base_el->type()); 00275 00276 // interpolate the base nodes' distances 00277 for (unsigned int n=0; n<n_base_nodes; n++) 00278 interpolated_dist += Point(base_el->point(n) - origin).size() 00279 * FE<1,LAGRANGE>::shape (base_mapping_elem_type, base_mapping_order, n, p); 00280 break; 00281 } 00282 00283 case 3: 00284 { 00285 const unsigned int n_base_nodes = base_el->n_nodes(); 00286 00287 const Point origin = inf_elem->origin(); 00288 const Order base_mapping_order (base_el->default_order()); 00289 const ElemType base_mapping_elem_type (base_el->type()); 00290 00291 // interpolate the base nodes' distances 00292 for (unsigned int n=0; n<n_base_nodes; n++) 00293 interpolated_dist += Point(base_el->point(n) - origin).size() 00294 * FE<2,LAGRANGE>::shape (base_mapping_elem_type, base_mapping_order, n, p); 00295 break; 00296 } 00297 00298 default: 00299 libmesh_error_msg("Unknown Dim = " << Dim); 00300 } 00301 00302 00303 00304 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00305 00306 // assumption on time-harmonic behavior 00307 const short int sign (-1); 00308 00309 // the wave number 00310 const Real wavenumber = 2. * libMesh::pi * data.frequency / data.speed; 00311 00312 // the exponent for time-harmonic behavior 00313 const Real exponent = sign /* +1. or -1. */ 00314 * wavenumber /* k */ 00315 * interpolated_dist /* together with next line: */ 00316 * InfFE<Dim,INFINITE_MAP,T_map>::eval(v, radial_mapping_order, 1); /* phase(s,t,v) */ 00317 00318 const Number time_harmonic = Number(cos(exponent), sin(exponent)); /* e^(sign*i*k*phase(s,t,v)) */ 00319 00320 /* 00321 * compute \p shape for all dof in the element 00322 */ 00323 if (Dim > 1) 00324 { 00325 const unsigned int n_dof = n_dofs (fet, inf_elem->type()); 00326 data.shape.resize(n_dof); 00327 00328 for (unsigned int i=0; i<n_dof; i++) 00329 { 00330 // compute base and radial shape indices 00331 unsigned int i_base, i_radial; 00332 compute_shape_indices(fet, inf_elem->type(), i, i_base, i_radial); 00333 00334 data.shape[i] = (InfFE<Dim,T_radial,T_map>::Radial::decay(v) /* (1.-v)/2. in 3D */ 00335 * FEInterface::shape(Dim-1, fet, base_el.get(), i_base, p) /* S_n(s,t) */ 00336 * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)) /* L_n(v) */ 00337 * time_harmonic; /* e^(sign*i*k*phase(s,t,v) */ 00338 } 00339 } 00340 else 00341 libmesh_error_msg("compute_data() for 1-dimensional InfFE not implemented."); 00342 00343 #else 00344 00345 const Real speed = data.speed; 00346 00347 /* 00348 * This is quite weird: the phase is actually 00349 * a measure how @e advanced the pressure is that 00350 * we compute. In other words: the further away 00351 * the node \p data.p is, the further we look into 00352 * the future... 00353 */ 00354 data.phase = interpolated_dist /* phase(s,t,v)/c */ 00355 * InfFE<Dim,INFINITE_MAP,T_map>::eval(v, radial_mapping_order, 1) / speed; 00356 00357 if (Dim > 1) 00358 { 00359 const unsigned int n_dof = n_dofs (fet, inf_elem->type()); 00360 data.shape.resize(n_dof); 00361 00362 for (unsigned int i=0; i<n_dof; i++) 00363 { 00364 // compute base and radial shape indices 00365 unsigned int i_base, i_radial; 00366 compute_shape_indices(fet, inf_elem->type(), i, i_base, i_radial); 00367 00368 data.shape[i] = InfFE<Dim,T_radial,T_map>::Radial::decay(v) /* (1.-v)/2. in 3D */ 00369 * FEInterface::shape(Dim-1, fet, base_el.get(), i_base, p) /* S_n(s,t) */ 00370 * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial); /* L_n(v) */ 00371 } 00372 } 00373 else 00374 libmesh_error_msg("compute_data() for 1-dimensional InfFE not implemented."); 00375 00376 #endif 00377 } 00378 00379 00380 00381 00382 00383 00384 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00385 void InfFE<Dim,T_radial,T_map>::compute_node_indices (const ElemType inf_elem_type, 00386 const unsigned int outer_node_index, 00387 unsigned int& base_node, 00388 unsigned int& radial_node) 00389 { 00390 switch (inf_elem_type) 00391 { 00392 case INFEDGE2: 00393 { 00394 libmesh_assert_less (outer_node_index, 2); 00395 base_node = 0; 00396 radial_node = outer_node_index; 00397 return; 00398 } 00399 00400 00401 // linear base approximation, easy to determine 00402 case INFQUAD4: 00403 { 00404 libmesh_assert_less (outer_node_index, 4); 00405 base_node = outer_node_index % 2; 00406 radial_node = outer_node_index / 2; 00407 return; 00408 } 00409 00410 case INFPRISM6: 00411 { 00412 libmesh_assert_less (outer_node_index, 6); 00413 base_node = outer_node_index % 3; 00414 radial_node = outer_node_index / 3; 00415 return; 00416 } 00417 00418 case INFHEX8: 00419 { 00420 libmesh_assert_less (outer_node_index, 8); 00421 base_node = outer_node_index % 4; 00422 radial_node = outer_node_index / 4; 00423 return; 00424 } 00425 00426 00427 // higher order base approximation, more work necessary 00428 case INFQUAD6: 00429 { 00430 switch (outer_node_index) 00431 { 00432 case 0: 00433 case 1: 00434 { 00435 radial_node = 0; 00436 base_node = outer_node_index; 00437 return; 00438 } 00439 00440 case 2: 00441 case 3: 00442 { 00443 radial_node = 1; 00444 base_node = outer_node_index-2; 00445 return; 00446 } 00447 00448 case 4: 00449 { 00450 radial_node = 0; 00451 base_node = 2; 00452 return; 00453 } 00454 00455 case 5: 00456 { 00457 radial_node = 1; 00458 base_node = 2; 00459 return; 00460 } 00461 00462 default: 00463 libmesh_error_msg("Unrecognized outer_node_index = " << outer_node_index); 00464 } 00465 } 00466 00467 00468 case INFHEX16: 00469 case INFHEX18: 00470 { 00471 switch (outer_node_index) 00472 { 00473 case 0: 00474 case 1: 00475 case 2: 00476 case 3: 00477 { 00478 radial_node = 0; 00479 base_node = outer_node_index; 00480 return; 00481 } 00482 00483 case 4: 00484 case 5: 00485 case 6: 00486 case 7: 00487 { 00488 radial_node = 1; 00489 base_node = outer_node_index-4; 00490 return; 00491 } 00492 00493 case 8: 00494 case 9: 00495 case 10: 00496 case 11: 00497 { 00498 radial_node = 0; 00499 base_node = outer_node_index-4; 00500 return; 00501 } 00502 00503 case 12: 00504 case 13: 00505 case 14: 00506 case 15: 00507 { 00508 radial_node = 1; 00509 base_node = outer_node_index-8; 00510 return; 00511 } 00512 00513 case 16: 00514 { 00515 libmesh_assert_equal_to (inf_elem_type, INFHEX18); 00516 radial_node = 0; 00517 base_node = 8; 00518 return; 00519 } 00520 00521 case 17: 00522 { 00523 libmesh_assert_equal_to (inf_elem_type, INFHEX18); 00524 radial_node = 1; 00525 base_node = 8; 00526 return; 00527 } 00528 00529 default: 00530 libmesh_error_msg("Unrecognized outer_node_index = " << outer_node_index); 00531 } 00532 } 00533 00534 00535 case INFPRISM12: 00536 { 00537 switch (outer_node_index) 00538 { 00539 case 0: 00540 case 1: 00541 case 2: 00542 { 00543 radial_node = 0; 00544 base_node = outer_node_index; 00545 return; 00546 } 00547 00548 case 3: 00549 case 4: 00550 case 5: 00551 { 00552 radial_node = 1; 00553 base_node = outer_node_index-3; 00554 return; 00555 } 00556 00557 case 6: 00558 case 7: 00559 case 8: 00560 { 00561 radial_node = 0; 00562 base_node = outer_node_index-3; 00563 return; 00564 } 00565 00566 case 9: 00567 case 10: 00568 case 11: 00569 { 00570 radial_node = 1; 00571 base_node = outer_node_index-6; 00572 return; 00573 } 00574 00575 default: 00576 libmesh_error_msg("Unrecognized outer_node_index = " << outer_node_index); 00577 } 00578 } 00579 00580 00581 default: 00582 libmesh_error_msg("ERROR: Bad infinite element type=" << inf_elem_type << ", node=" << outer_node_index); 00583 } 00584 } 00585 00586 00587 00588 00589 00590 00591 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00592 void InfFE<Dim,T_radial,T_map>::compute_node_indices_fast (const ElemType inf_elem_type, 00593 const unsigned int outer_node_index, 00594 unsigned int& base_node, 00595 unsigned int& radial_node) 00596 { 00597 libmesh_assert_not_equal_to (inf_elem_type, INVALID_ELEM); 00598 00599 static std::vector<unsigned int> _static_base_node_index; 00600 static std::vector<unsigned int> _static_radial_node_index; 00601 00602 /* 00603 * fast counterpart to compute_node_indices(), uses local static buffers 00604 * to store the index maps. The class member 00605 * \p _compute_node_indices_fast_current_elem_type remembers 00606 * the current element type. 00607 * 00608 * Note that there exist non-static members storing the 00609 * same data. However, you never know what element type 00610 * is currently used by the \p InfFE object, and what 00611 * request is currently directed to the static \p InfFE 00612 * members (which use \p compute_node_indices_fast()). 00613 * So separate these. 00614 * 00615 * check whether the work for this elemtype has already 00616 * been done. If so, use this index. Otherwise, refresh 00617 * the buffer to this element type. 00618 */ 00619 if (inf_elem_type==_compute_node_indices_fast_current_elem_type) 00620 { 00621 base_node = _static_base_node_index [outer_node_index]; 00622 radial_node = _static_radial_node_index[outer_node_index]; 00623 return; 00624 } 00625 else 00626 { 00627 // store the map for _all_ nodes for this element type 00628 _compute_node_indices_fast_current_elem_type = inf_elem_type; 00629 00630 unsigned int n_nodes = libMesh::invalid_uint; 00631 00632 switch (inf_elem_type) 00633 { 00634 case INFEDGE2: 00635 { 00636 n_nodes = 2; 00637 break; 00638 } 00639 case INFQUAD4: 00640 { 00641 n_nodes = 4; 00642 break; 00643 } 00644 case INFQUAD6: 00645 { 00646 n_nodes = 6; 00647 break; 00648 } 00649 case INFHEX8: 00650 { 00651 n_nodes = 8; 00652 break; 00653 } 00654 case INFHEX16: 00655 { 00656 n_nodes = 16; 00657 break; 00658 } 00659 case INFHEX18: 00660 { 00661 n_nodes = 18; 00662 break; 00663 } 00664 case INFPRISM6: 00665 { 00666 n_nodes = 6; 00667 break; 00668 } 00669 case INFPRISM12: 00670 { 00671 n_nodes = 12; 00672 break; 00673 } 00674 default: 00675 libmesh_error_msg("ERROR: Bad infinite element type=" << inf_elem_type << ", node=" << outer_node_index); 00676 } 00677 00678 00679 _static_base_node_index.resize (n_nodes); 00680 _static_radial_node_index.resize(n_nodes); 00681 00682 for (unsigned int n=0; n<n_nodes; n++) 00683 compute_node_indices (inf_elem_type, 00684 n, 00685 _static_base_node_index [outer_node_index], 00686 _static_radial_node_index[outer_node_index]); 00687 00688 // and return for the specified node 00689 base_node = _static_base_node_index [outer_node_index]; 00690 radial_node = _static_radial_node_index[outer_node_index]; 00691 return; 00692 } 00693 } 00694 00695 00696 00697 00698 00699 00700 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00701 void InfFE<Dim,T_radial,T_map>::compute_shape_indices (const FEType& fet, 00702 const ElemType inf_elem_type, 00703 const unsigned int i, 00704 unsigned int& base_shape, 00705 unsigned int& radial_shape) 00706 { 00707 00708 /* 00709 * An example is provided: the numbers in comments refer to 00710 * a fictitious InfHex18. The numbers are chosen as exemplary 00711 * values. There is currently no base approximation that 00712 * requires this many dof's at nodes, sides, faces and in the element. 00713 * 00714 * the order of the shape functions is heavily related with the 00715 * order the dofs are assigned in \p DofMap::distributed_dofs(). 00716 * Due to the infinite elements with higher-order base approximation, 00717 * some more effort is necessary. 00718 * 00719 * numbering scheme: 00720 * 1. all vertices in the base, assign node->n_comp() dofs to each vertex 00721 * 2. all vertices further out: innermost loop: radial shapes, 00722 * then the base approximation shapes 00723 * 3. all side nodes in the base, assign node->n_comp() dofs to each side node 00724 * 4. all side nodes further out: innermost loop: radial shapes, 00725 * then the base approximation shapes 00726 * 5. (all) face nodes in the base, assign node->n_comp() dofs to each face node 00727 * 6. (all) face nodes further out: innermost loop: radial shapes, 00728 * then the base approximation shapes 00729 * 7. element-associated dof in the base 00730 * 8. element-associated dof further out 00731 */ 00732 00733 const unsigned int radial_order = static_cast<unsigned int>(fet.radial_order); // 4 00734 const unsigned int radial_order_p_one = radial_order+1; // 5 00735 00736 const ElemType base_elem_type (Base::get_elem_type(inf_elem_type)); // QUAD9 00737 00738 // assume that the number of dof is the same for all vertices 00739 unsigned int n_base_vertices = libMesh::invalid_uint; // 4 00740 const unsigned int n_base_vertex_dof = FEInterface::n_dofs_at_node (Dim-1, fet, base_elem_type, 0);// 2 00741 00742 unsigned int n_base_side_nodes = libMesh::invalid_uint; // 4 00743 unsigned int n_base_side_dof = libMesh::invalid_uint; // 3 00744 00745 unsigned int n_base_face_nodes = libMesh::invalid_uint; // 1 00746 unsigned int n_base_face_dof = libMesh::invalid_uint; // 5 00747 00748 const unsigned int n_base_elem_dof = FEInterface::n_dofs_per_elem (Dim-1, fet, base_elem_type);// 9 00749 00750 00751 switch (inf_elem_type) 00752 { 00753 case INFEDGE2: 00754 { 00755 n_base_vertices = 1; 00756 n_base_side_nodes = 0; 00757 n_base_face_nodes = 0; 00758 n_base_side_dof = 0; 00759 n_base_face_dof = 0; 00760 break; 00761 } 00762 00763 case INFQUAD4: 00764 { 00765 n_base_vertices = 2; 00766 n_base_side_nodes = 0; 00767 n_base_face_nodes = 0; 00768 n_base_side_dof = 0; 00769 n_base_face_dof = 0; 00770 break; 00771 } 00772 00773 case INFQUAD6: 00774 { 00775 n_base_vertices = 2; 00776 n_base_side_nodes = 1; 00777 n_base_face_nodes = 0; 00778 n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices); 00779 n_base_face_dof = 0; 00780 break; 00781 } 00782 00783 case INFHEX8: 00784 { 00785 n_base_vertices = 4; 00786 n_base_side_nodes = 0; 00787 n_base_face_nodes = 0; 00788 n_base_side_dof = 0; 00789 n_base_face_dof = 0; 00790 break; 00791 } 00792 00793 case INFHEX16: 00794 { 00795 n_base_vertices = 4; 00796 n_base_side_nodes = 4; 00797 n_base_face_nodes = 0; 00798 n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices); 00799 n_base_face_dof = 0; 00800 break; 00801 } 00802 00803 case INFHEX18: 00804 { 00805 n_base_vertices = 4; 00806 n_base_side_nodes = 4; 00807 n_base_face_nodes = 1; 00808 n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices); 00809 n_base_face_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, 8); 00810 break; 00811 } 00812 00813 00814 case INFPRISM6: 00815 { 00816 n_base_vertices = 3; 00817 n_base_side_nodes = 0; 00818 n_base_face_nodes = 0; 00819 n_base_side_dof = 0; 00820 n_base_face_dof = 0; 00821 break; 00822 } 00823 00824 case INFPRISM12: 00825 { 00826 n_base_vertices = 3; 00827 n_base_side_nodes = 3; 00828 n_base_face_nodes = 0; 00829 n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices); 00830 n_base_face_dof = 0; 00831 break; 00832 } 00833 00834 default: 00835 libmesh_error_msg("Unrecognized inf_elem_type = " << inf_elem_type); 00836 } 00837 00838 00839 { 00840 // these are the limits describing the intervals where the shape function lies 00841 const unsigned int n_dof_at_base_vertices = n_base_vertices*n_base_vertex_dof; // 8 00842 const unsigned int n_dof_at_all_vertices = n_dof_at_base_vertices*radial_order_p_one; // 40 00843 00844 const unsigned int n_dof_at_base_sides = n_base_side_nodes*n_base_side_dof; // 12 00845 const unsigned int n_dof_at_all_sides = n_dof_at_base_sides*radial_order_p_one; // 60 00846 00847 const unsigned int n_dof_at_base_face = n_base_face_nodes*n_base_face_dof; // 5 00848 const unsigned int n_dof_at_all_faces = n_dof_at_base_face*radial_order_p_one; // 25 00849 00850 00851 // start locating the shape function 00852 if (i < n_dof_at_base_vertices) // range of i: 0..7 00853 { 00854 // belongs to vertex in the base 00855 radial_shape = 0; 00856 base_shape = i; 00857 } 00858 00859 else if (i < n_dof_at_all_vertices) // range of i: 8..39 00860 { 00861 /* belongs to vertex in the outer shell 00862 * 00863 * subtract the number of dof already counted, 00864 * so that i_offset contains only the offset for the base 00865 */ 00866 const unsigned int i_offset = i - n_dof_at_base_vertices; // 0..31 00867 00868 // first the radial dof are counted, then the base dof 00869 radial_shape = (i_offset % radial_order) + 1; 00870 base_shape = i_offset / radial_order; 00871 } 00872 00873 else if (i < n_dof_at_all_vertices+n_dof_at_base_sides) // range of i: 40..51 00874 { 00875 // belongs to base, is a side node 00876 radial_shape = 0; 00877 base_shape = i - radial_order * n_dof_at_base_vertices; // 8..19 00878 } 00879 00880 else if (i < n_dof_at_all_vertices+n_dof_at_all_sides) // range of i: 52..99 00881 { 00882 // belongs to side node in the outer shell 00883 const unsigned int i_offset = i - (n_dof_at_all_vertices 00884 + n_dof_at_base_sides); // 0..47 00885 radial_shape = (i_offset % radial_order) + 1; 00886 base_shape = (i_offset / radial_order) + n_dof_at_base_vertices; 00887 } 00888 00889 else if (i < n_dof_at_all_vertices+n_dof_at_all_sides+n_dof_at_base_face) // range of i: 100..104 00890 { 00891 // belongs to the node in the base face 00892 radial_shape = 0; 00893 base_shape = i - radial_order*(n_dof_at_base_vertices 00894 + n_dof_at_base_sides); // 20..24 00895 } 00896 00897 else if (i < n_dof_at_all_vertices+n_dof_at_all_sides+n_dof_at_all_faces) // range of i: 105..124 00898 { 00899 // belongs to the node in the outer face 00900 const unsigned int i_offset = i - (n_dof_at_all_vertices 00901 + n_dof_at_all_sides 00902 + n_dof_at_base_face); // 0..19 00903 radial_shape = (i_offset % radial_order) + 1; 00904 base_shape = (i_offset / radial_order) + n_dof_at_base_vertices + n_dof_at_base_sides; 00905 } 00906 00907 else if (i < n_dof_at_all_vertices+n_dof_at_all_sides+n_dof_at_all_faces+n_base_elem_dof) // range of i: 125..133 00908 { 00909 // belongs to the base and is an element associated shape 00910 radial_shape = 0; 00911 base_shape = i - (n_dof_at_all_vertices 00912 + n_dof_at_all_sides 00913 + n_dof_at_all_faces); // 0..8 00914 } 00915 00916 else // range of i: 134..169 00917 { 00918 libmesh_assert_less (i, n_dofs(fet, inf_elem_type)); 00919 // belongs to the outer shell and is an element associated shape 00920 const unsigned int i_offset = i - (n_dof_at_all_vertices 00921 + n_dof_at_all_sides 00922 + n_dof_at_all_faces 00923 + n_base_elem_dof); // 0..19 00924 radial_shape = (i_offset % radial_order) + 1; 00925 base_shape = (i_offset / radial_order) + n_dof_at_base_vertices + n_dof_at_base_sides + n_dof_at_base_face; 00926 } 00927 } 00928 00929 return; 00930 } 00931 00932 00933 00934 //-------------------------------------------------------------- 00935 // Explicit instantiations 00936 //#include "libmesh/inf_fe_instantiate_1D.h" 00937 //#include "libmesh/inf_fe_instantiate_2D.h" 00938 //#include "libmesh/inf_fe_instantiate_3D.h" 00939 00940 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,unsigned int,n_dofs(const FEType&,const ElemType)); 00941 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,unsigned int,n_dofs(const FEType&,const ElemType)); 00942 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,unsigned int,n_dofs(const FEType&,const ElemType)); 00943 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,unsigned int,n_dofs_per_elem(const FEType&,const ElemType)); 00944 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,unsigned int,n_dofs_per_elem(const FEType&,const ElemType)); 00945 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,unsigned int,n_dofs_per_elem(const FEType&,const ElemType)); 00946 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,unsigned int,n_dofs_at_node(const FEType&,const ElemType,const unsigned int)); 00947 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,unsigned int,n_dofs_at_node(const FEType&,const ElemType,const unsigned int)); 00948 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,unsigned int,n_dofs_at_node(const FEType&,const ElemType,const unsigned int)); 00949 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,compute_shape_indices(const FEType&,const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00950 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,compute_shape_indices(const FEType&,const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00951 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,compute_shape_indices(const FEType&,const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00952 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,compute_node_indices(const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00953 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,compute_node_indices(const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00954 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,compute_node_indices(const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00955 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,Real,shape(const FEType&,const Elem*,const unsigned int,const Point& p)); 00956 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,Real,shape(const FEType&,const Elem*,const unsigned int,const Point& p)); 00957 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,Real,shape(const FEType&,const Elem*,const unsigned int,const Point& p)); 00958 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,Real,shape(const FEType&,const ElemType,const unsigned int,const Point&)); 00959 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,Real,shape(const FEType&,const ElemType,const unsigned int,const Point&)); 00960 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,Real,shape(const FEType&,const ElemType,const unsigned int,const Point&)); 00961 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,compute_data(const FEType&,const Elem*,FEComputeData&)); 00962 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,compute_data(const FEType&,const Elem*,FEComputeData&)); 00963 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,compute_data(const FEType&,const Elem*,FEComputeData&)); 00964 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,nodal_soln(const FEType&,const Elem*,const std::vector<Number>&,std::vector<Number>&)); 00965 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,nodal_soln(const FEType&,const Elem*,const std::vector<Number>&,std::vector<Number>&)); 00966 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,nodal_soln(const FEType&,const Elem*,const std::vector<Number>&,std::vector<Number>&)); 00967 00968 } // namespace libMesh 00969 00970 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS