inf_fe_static.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 // 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 AutoPtr<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 AutoPtr<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 #ifdef DEBUG 00298 default: 00299 libmesh_error(); 00300 #endif 00301 } 00302 00303 00304 00305 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00306 00307 // assumption on time-harmonic behavior 00308 const short int sign (-1); 00309 00310 // the wave number 00311 const Real wavenumber = 2. * libMesh::pi * data.frequency / data.speed; 00312 00313 // the exponent for time-harmonic behavior 00314 const Real exponent = sign /* +1. or -1. */ 00315 * wavenumber /* k */ 00316 * interpolated_dist /* together with next line: */ 00317 * InfFE<Dim,INFINITE_MAP,T_map>::eval(v, radial_mapping_order, 1); /* phase(s,t,v) */ 00318 00319 const Number time_harmonic = Number(cos(exponent), sin(exponent)); /* e^(sign*i*k*phase(s,t,v)) */ 00320 00321 /* 00322 * compute \p shape for all dof in the element 00323 */ 00324 if (Dim > 1) 00325 { 00326 const unsigned int n_dof = n_dofs (fet, inf_elem->type()); 00327 data.shape.resize(n_dof); 00328 00329 for (unsigned int i=0; i<n_dof; i++) 00330 { 00331 // compute base and radial shape indices 00332 unsigned int i_base, i_radial; 00333 compute_shape_indices(fet, inf_elem->type(), i, i_base, i_radial); 00334 00335 data.shape[i] = (InfFE<Dim,T_radial,T_map>::Radial::decay(v) /* (1.-v)/2. in 3D */ 00336 * FEInterface::shape(Dim-1, fet, base_el.get(), i_base, p) /* S_n(s,t) */ 00337 * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)) /* L_n(v) */ 00338 * time_harmonic; /* e^(sign*i*k*phase(s,t,v) */ 00339 } 00340 } 00341 else 00342 { 00343 libMesh::err << "compute_data() for 1-dimensional InfFE not implemented." << std::endl; 00344 libmesh_error(); 00345 } 00346 00347 #else 00348 00349 const Real speed = data.speed; 00350 00351 /* 00352 * This is quite weird: the phase is actually 00353 * a measure how @e advanced the pressure is that 00354 * we compute. In other words: the further away 00355 * the node \p data.p is, the further we look into 00356 * the future... 00357 */ 00358 data.phase = interpolated_dist /* phase(s,t,v)/c */ 00359 * InfFE<Dim,INFINITE_MAP,T_map>::eval(v, radial_mapping_order, 1) / speed; 00360 00361 if (Dim > 1) 00362 { 00363 const unsigned int n_dof = n_dofs (fet, inf_elem->type()); 00364 data.shape.resize(n_dof); 00365 00366 for (unsigned int i=0; i<n_dof; i++) 00367 { 00368 // compute base and radial shape indices 00369 unsigned int i_base, i_radial; 00370 compute_shape_indices(fet, inf_elem->type(), i, i_base, i_radial); 00371 00372 data.shape[i] = InfFE<Dim,T_radial,T_map>::Radial::decay(v) /* (1.-v)/2. in 3D */ 00373 * FEInterface::shape(Dim-1, fet, base_el.get(), i_base, p) /* S_n(s,t) */ 00374 * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial); /* L_n(v) */ 00375 } 00376 } 00377 else 00378 { 00379 libMesh::err << "compute_data() for 1-dimensional InfFE not implemented." << std::endl; 00380 libmesh_error(); 00381 } 00382 00383 #endif 00384 } 00385 00386 00387 00388 00389 00390 00391 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00392 void InfFE<Dim,T_radial,T_map>::compute_node_indices (const ElemType inf_elem_type, 00393 const unsigned int outer_node_index, 00394 unsigned int& base_node, 00395 unsigned int& radial_node) 00396 { 00397 switch (inf_elem_type) 00398 { 00399 case INFEDGE2: 00400 { 00401 libmesh_assert_less (outer_node_index, 2); 00402 base_node = 0; 00403 radial_node = outer_node_index; 00404 return; 00405 } 00406 00407 00408 // linear base approximation, easy to determine 00409 case INFQUAD4: 00410 { 00411 libmesh_assert_less (outer_node_index, 4); 00412 base_node = outer_node_index % 2; 00413 radial_node = outer_node_index / 2; 00414 return; 00415 } 00416 00417 case INFPRISM6: 00418 { 00419 libmesh_assert_less (outer_node_index, 6); 00420 base_node = outer_node_index % 3; 00421 radial_node = outer_node_index / 3; 00422 return; 00423 } 00424 00425 case INFHEX8: 00426 { 00427 libmesh_assert_less (outer_node_index, 8); 00428 base_node = outer_node_index % 4; 00429 radial_node = outer_node_index / 4; 00430 return; 00431 } 00432 00433 00434 // higher order base approximation, more work necessary 00435 case INFQUAD6: 00436 { 00437 switch (outer_node_index) 00438 { 00439 case 0: 00440 case 1: 00441 { 00442 radial_node = 0; 00443 base_node = outer_node_index; 00444 return; 00445 } 00446 00447 case 2: 00448 case 3: 00449 { 00450 radial_node = 1; 00451 base_node = outer_node_index-2; 00452 return; 00453 } 00454 00455 case 4: 00456 { 00457 radial_node = 0; 00458 base_node = 2; 00459 return; 00460 } 00461 00462 case 5: 00463 { 00464 radial_node = 1; 00465 base_node = 2; 00466 return; 00467 } 00468 00469 default: 00470 { 00471 libmesh_error(); 00472 return; 00473 } 00474 } 00475 } 00476 00477 00478 case INFHEX16: 00479 case INFHEX18: 00480 { 00481 switch (outer_node_index) 00482 { 00483 case 0: 00484 case 1: 00485 case 2: 00486 case 3: 00487 { 00488 radial_node = 0; 00489 base_node = outer_node_index; 00490 return; 00491 } 00492 00493 case 4: 00494 case 5: 00495 case 6: 00496 case 7: 00497 { 00498 radial_node = 1; 00499 base_node = outer_node_index-4; 00500 return; 00501 } 00502 00503 case 8: 00504 case 9: 00505 case 10: 00506 case 11: 00507 { 00508 radial_node = 0; 00509 base_node = outer_node_index-4; 00510 return; 00511 } 00512 00513 case 12: 00514 case 13: 00515 case 14: 00516 case 15: 00517 { 00518 radial_node = 1; 00519 base_node = outer_node_index-8; 00520 return; 00521 } 00522 00523 case 16: 00524 { 00525 libmesh_assert_equal_to (inf_elem_type, INFHEX18); 00526 radial_node = 0; 00527 base_node = 8; 00528 return; 00529 } 00530 00531 case 17: 00532 { 00533 libmesh_assert_equal_to (inf_elem_type, INFHEX18); 00534 radial_node = 1; 00535 base_node = 8; 00536 return; 00537 } 00538 00539 default: 00540 { 00541 libmesh_error(); 00542 return; 00543 } 00544 } 00545 } 00546 00547 00548 case INFPRISM12: 00549 { 00550 switch (outer_node_index) 00551 { 00552 case 0: 00553 case 1: 00554 case 2: 00555 { 00556 radial_node = 0; 00557 base_node = outer_node_index; 00558 return; 00559 } 00560 00561 case 3: 00562 case 4: 00563 case 5: 00564 { 00565 radial_node = 1; 00566 base_node = outer_node_index-3; 00567 return; 00568 } 00569 00570 case 6: 00571 case 7: 00572 case 8: 00573 { 00574 radial_node = 0; 00575 base_node = outer_node_index-3; 00576 return; 00577 } 00578 00579 case 9: 00580 case 10: 00581 case 11: 00582 { 00583 radial_node = 1; 00584 base_node = outer_node_index-6; 00585 return; 00586 } 00587 00588 default: 00589 { 00590 libmesh_error(); 00591 return; 00592 } 00593 } 00594 } 00595 00596 00597 default: 00598 { 00599 libMesh::err << "ERROR: Bad infinite element type=" << inf_elem_type 00600 << ", node=" << outer_node_index << std::endl; 00601 libmesh_error(); 00602 return; 00603 } 00604 } 00605 } 00606 00607 00608 00609 00610 00611 00612 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00613 void InfFE<Dim,T_radial,T_map>::compute_node_indices_fast (const ElemType inf_elem_type, 00614 const unsigned int outer_node_index, 00615 unsigned int& base_node, 00616 unsigned int& radial_node) 00617 { 00618 libmesh_assert_not_equal_to (inf_elem_type, INVALID_ELEM); 00619 00620 static std::vector<unsigned int> _static_base_node_index; 00621 static std::vector<unsigned int> _static_radial_node_index; 00622 00623 /* 00624 * fast counterpart to compute_node_indices(), uses local static buffers 00625 * to store the index maps. The class member 00626 * \p _compute_node_indices_fast_current_elem_type remembers 00627 * the current element type. 00628 * 00629 * Note that there exist non-static members storing the 00630 * same data. However, you never know what element type 00631 * is currently used by the \p InfFE object, and what 00632 * request is currently directed to the static \p InfFE 00633 * members (which use \p compute_node_indices_fast()). 00634 * So separate these. 00635 * 00636 * check whether the work for this elemtype has already 00637 * been done. If so, use this index. Otherwise, refresh 00638 * the buffer to this element type. 00639 */ 00640 if (inf_elem_type==_compute_node_indices_fast_current_elem_type) 00641 { 00642 base_node = _static_base_node_index [outer_node_index]; 00643 radial_node = _static_radial_node_index[outer_node_index]; 00644 return; 00645 } 00646 else 00647 { 00648 // store the map for _all_ nodes for this element type 00649 _compute_node_indices_fast_current_elem_type = inf_elem_type; 00650 00651 unsigned int n_nodes = libMesh::invalid_uint; 00652 00653 switch (inf_elem_type) 00654 { 00655 case INFEDGE2: 00656 { 00657 n_nodes = 2; 00658 break; 00659 } 00660 case INFQUAD4: 00661 { 00662 n_nodes = 4; 00663 break; 00664 } 00665 case INFQUAD6: 00666 { 00667 n_nodes = 6; 00668 break; 00669 } 00670 case INFHEX8: 00671 { 00672 n_nodes = 8; 00673 break; 00674 } 00675 case INFHEX16: 00676 { 00677 n_nodes = 16; 00678 break; 00679 } 00680 case INFHEX18: 00681 { 00682 n_nodes = 18; 00683 break; 00684 } 00685 case INFPRISM6: 00686 { 00687 n_nodes = 6; 00688 break; 00689 } 00690 case INFPRISM12: 00691 { 00692 n_nodes = 12; 00693 break; 00694 } 00695 default: 00696 { 00697 libMesh::err << "ERROR: Bad infinite element type=" << inf_elem_type 00698 << ", node=" << outer_node_index << std::endl; 00699 libmesh_error(); 00700 break; 00701 } 00702 } 00703 00704 00705 _static_base_node_index.resize (n_nodes); 00706 _static_radial_node_index.resize(n_nodes); 00707 00708 for (unsigned int n=0; n<n_nodes; n++) 00709 compute_node_indices (inf_elem_type, 00710 n, 00711 _static_base_node_index [outer_node_index], 00712 _static_radial_node_index[outer_node_index]); 00713 00714 // and return for the specified node 00715 base_node = _static_base_node_index [outer_node_index]; 00716 radial_node = _static_radial_node_index[outer_node_index]; 00717 return; 00718 } 00719 } 00720 00721 00722 00723 00724 00725 00726 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00727 void InfFE<Dim,T_radial,T_map>::compute_shape_indices (const FEType& fet, 00728 const ElemType inf_elem_type, 00729 const unsigned int i, 00730 unsigned int& base_shape, 00731 unsigned int& radial_shape) 00732 { 00733 00734 /* 00735 * An example is provided: the numbers in comments refer to 00736 * a fictitious InfHex18. The numbers are chosen as exemplary 00737 * values. There is currently no base approximation that 00738 * requires this many dof's at nodes, sides, faces and in the element. 00739 * 00740 * the order of the shape functions is heavily related with the 00741 * order the dofs are assigned in \p DofMap::distributed_dofs(). 00742 * Due to the infinite elements with higher-order base approximation, 00743 * some more effort is necessary. 00744 * 00745 * numbering scheme: 00746 * 1. all vertices in the base, assign node->n_comp() dofs to each vertex 00747 * 2. all vertices further out: innermost loop: radial shapes, 00748 * then the base approximation shapes 00749 * 3. all side nodes in the base, assign node->n_comp() dofs to each side node 00750 * 4. all side nodes further out: innermost loop: radial shapes, 00751 * then the base approximation shapes 00752 * 5. (all) face nodes in the base, assign node->n_comp() dofs to each face node 00753 * 6. (all) face nodes further out: innermost loop: radial shapes, 00754 * then the base approximation shapes 00755 * 7. element-associated dof in the base 00756 * 8. element-associated dof further out 00757 */ 00758 00759 const unsigned int radial_order = static_cast<unsigned int>(fet.radial_order); // 4 00760 const unsigned int radial_order_p_one = radial_order+1; // 5 00761 00762 const ElemType base_elem_type (Base::get_elem_type(inf_elem_type)); // QUAD9 00763 00764 // assume that the number of dof is the same for all vertices 00765 unsigned int n_base_vertices = libMesh::invalid_uint; // 4 00766 const unsigned int n_base_vertex_dof = FEInterface::n_dofs_at_node (Dim-1, fet, base_elem_type, 0);// 2 00767 00768 unsigned int n_base_side_nodes = libMesh::invalid_uint; // 4 00769 unsigned int n_base_side_dof = libMesh::invalid_uint; // 3 00770 00771 unsigned int n_base_face_nodes = libMesh::invalid_uint; // 1 00772 unsigned int n_base_face_dof = libMesh::invalid_uint; // 5 00773 00774 const unsigned int n_base_elem_dof = FEInterface::n_dofs_per_elem (Dim-1, fet, base_elem_type);// 9 00775 00776 00777 switch (inf_elem_type) 00778 { 00779 case INFEDGE2: 00780 { 00781 n_base_vertices = 1; 00782 n_base_side_nodes = 0; 00783 n_base_face_nodes = 0; 00784 n_base_side_dof = 0; 00785 n_base_face_dof = 0; 00786 break; 00787 } 00788 00789 case INFQUAD4: 00790 { 00791 n_base_vertices = 2; 00792 n_base_side_nodes = 0; 00793 n_base_face_nodes = 0; 00794 n_base_side_dof = 0; 00795 n_base_face_dof = 0; 00796 break; 00797 } 00798 00799 case INFQUAD6: 00800 { 00801 n_base_vertices = 2; 00802 n_base_side_nodes = 1; 00803 n_base_face_nodes = 0; 00804 n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices); 00805 n_base_face_dof = 0; 00806 break; 00807 } 00808 00809 case INFHEX8: 00810 { 00811 n_base_vertices = 4; 00812 n_base_side_nodes = 0; 00813 n_base_face_nodes = 0; 00814 n_base_side_dof = 0; 00815 n_base_face_dof = 0; 00816 break; 00817 } 00818 00819 case INFHEX16: 00820 { 00821 n_base_vertices = 4; 00822 n_base_side_nodes = 4; 00823 n_base_face_nodes = 0; 00824 n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices); 00825 n_base_face_dof = 0; 00826 break; 00827 } 00828 00829 case INFHEX18: 00830 { 00831 n_base_vertices = 4; 00832 n_base_side_nodes = 4; 00833 n_base_face_nodes = 1; 00834 n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices); 00835 n_base_face_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, 8); 00836 break; 00837 } 00838 00839 00840 case INFPRISM6: 00841 { 00842 n_base_vertices = 3; 00843 n_base_side_nodes = 0; 00844 n_base_face_nodes = 0; 00845 n_base_side_dof = 0; 00846 n_base_face_dof = 0; 00847 break; 00848 } 00849 00850 case INFPRISM12: 00851 { 00852 n_base_vertices = 3; 00853 n_base_side_nodes = 3; 00854 n_base_face_nodes = 0; 00855 n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices); 00856 n_base_face_dof = 0; 00857 break; 00858 } 00859 00860 default: 00861 libmesh_error(); 00862 } 00863 00864 00865 { 00866 // these are the limits describing the intervals where the shape function lies 00867 const unsigned int n_dof_at_base_vertices = n_base_vertices*n_base_vertex_dof; // 8 00868 const unsigned int n_dof_at_all_vertices = n_dof_at_base_vertices*radial_order_p_one; // 40 00869 00870 const unsigned int n_dof_at_base_sides = n_base_side_nodes*n_base_side_dof; // 12 00871 const unsigned int n_dof_at_all_sides = n_dof_at_base_sides*radial_order_p_one; // 60 00872 00873 const unsigned int n_dof_at_base_face = n_base_face_nodes*n_base_face_dof; // 5 00874 const unsigned int n_dof_at_all_faces = n_dof_at_base_face*radial_order_p_one; // 25 00875 00876 00877 // start locating the shape function 00878 if (i < n_dof_at_base_vertices) // range of i: 0..7 00879 { 00880 // belongs to vertex in the base 00881 radial_shape = 0; 00882 base_shape = i; 00883 } 00884 00885 else if (i < n_dof_at_all_vertices) // range of i: 8..39 00886 { 00887 /* belongs to vertex in the outer shell 00888 * 00889 * subtract the number of dof already counted, 00890 * so that i_offset contains only the offset for the base 00891 */ 00892 const unsigned int i_offset = i - n_dof_at_base_vertices; // 0..31 00893 00894 // first the radial dof are counted, then the base dof 00895 radial_shape = (i_offset % radial_order) + 1; 00896 base_shape = i_offset / radial_order; 00897 } 00898 00899 else if (i < n_dof_at_all_vertices+n_dof_at_base_sides) // range of i: 40..51 00900 { 00901 // belongs to base, is a side node 00902 radial_shape = 0; 00903 base_shape = i - radial_order * n_dof_at_base_vertices; // 8..19 00904 } 00905 00906 else if (i < n_dof_at_all_vertices+n_dof_at_all_sides) // range of i: 52..99 00907 { 00908 // belongs to side node in the outer shell 00909 const unsigned int i_offset = i - (n_dof_at_all_vertices 00910 + n_dof_at_base_sides); // 0..47 00911 radial_shape = (i_offset % radial_order) + 1; 00912 base_shape = (i_offset / radial_order) + n_dof_at_base_vertices; 00913 } 00914 00915 else if (i < n_dof_at_all_vertices+n_dof_at_all_sides+n_dof_at_base_face) // range of i: 100..104 00916 { 00917 // belongs to the node in the base face 00918 radial_shape = 0; 00919 base_shape = i - radial_order*(n_dof_at_base_vertices 00920 + n_dof_at_base_sides); // 20..24 00921 } 00922 00923 else if (i < n_dof_at_all_vertices+n_dof_at_all_sides+n_dof_at_all_faces) // range of i: 105..124 00924 { 00925 // belongs to the node in the outer face 00926 const unsigned int i_offset = i - (n_dof_at_all_vertices 00927 + n_dof_at_all_sides 00928 + n_dof_at_base_face); // 0..19 00929 radial_shape = (i_offset % radial_order) + 1; 00930 base_shape = (i_offset / radial_order) + n_dof_at_base_vertices + n_dof_at_base_sides; 00931 } 00932 00933 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 00934 { 00935 // belongs to the base and is an element associated shape 00936 radial_shape = 0; 00937 base_shape = i - (n_dof_at_all_vertices 00938 + n_dof_at_all_sides 00939 + n_dof_at_all_faces); // 0..8 00940 } 00941 00942 else // range of i: 134..169 00943 { 00944 libmesh_assert_less (i, n_dofs(fet, inf_elem_type)); 00945 // belongs to the outer shell and is an element associated shape 00946 const unsigned int i_offset = i - (n_dof_at_all_vertices 00947 + n_dof_at_all_sides 00948 + n_dof_at_all_faces 00949 + n_base_elem_dof); // 0..19 00950 radial_shape = (i_offset % radial_order) + 1; 00951 base_shape = (i_offset / radial_order) + n_dof_at_base_vertices + n_dof_at_base_sides + n_dof_at_base_face; 00952 } 00953 } 00954 00955 return; 00956 } 00957 00958 00959 00960 //-------------------------------------------------------------- 00961 // Explicit instantiations 00962 //#include "libmesh/inf_fe_instantiate_1D.h" 00963 //#include "libmesh/inf_fe_instantiate_2D.h" 00964 //#include "libmesh/inf_fe_instantiate_3D.h" 00965 00966 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,unsigned int,n_dofs(const FEType&,const ElemType)); 00967 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,unsigned int,n_dofs(const FEType&,const ElemType)); 00968 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,unsigned int,n_dofs(const FEType&,const ElemType)); 00969 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,unsigned int,n_dofs_per_elem(const FEType&,const ElemType)); 00970 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,unsigned int,n_dofs_per_elem(const FEType&,const ElemType)); 00971 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,unsigned int,n_dofs_per_elem(const FEType&,const ElemType)); 00972 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,unsigned int,n_dofs_at_node(const FEType&,const ElemType,const unsigned int)); 00973 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,unsigned int,n_dofs_at_node(const FEType&,const ElemType,const unsigned int)); 00974 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,unsigned int,n_dofs_at_node(const FEType&,const ElemType,const unsigned int)); 00975 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,compute_shape_indices(const FEType&,const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00976 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,compute_shape_indices(const FEType&,const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00977 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,compute_shape_indices(const FEType&,const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00978 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,compute_node_indices(const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00979 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,compute_node_indices(const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00980 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,compute_node_indices(const ElemType,const unsigned int,unsigned int&,unsigned int&)); 00981 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,Real,shape(const FEType&,const Elem*,const unsigned int,const Point& p)); 00982 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,Real,shape(const FEType&,const Elem*,const unsigned int,const Point& p)); 00983 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,Real,shape(const FEType&,const Elem*,const unsigned int,const Point& p)); 00984 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,Real,shape(const FEType&,const ElemType,const unsigned int,const Point&)); 00985 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,Real,shape(const FEType&,const ElemType,const unsigned int,const Point&)); 00986 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,Real,shape(const FEType&,const ElemType,const unsigned int,const Point&)); 00987 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,compute_data(const FEType&,const Elem*,FEComputeData&)); 00988 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,compute_data(const FEType&,const Elem*,FEComputeData&)); 00989 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,compute_data(const FEType&,const Elem*,FEComputeData&)); 00990 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,void,nodal_soln(const FEType&,const Elem*,const std::vector<Number>&,std::vector<Number>&)); 00991 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,void,nodal_soln(const FEType&,const Elem*,const std::vector<Number>&,std::vector<Number>&)); 00992 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,void,nodal_soln(const FEType&,const Elem*,const std::vector<Number>&,std::vector<Number>&)); 00993 00994 } // namespace libMesh 00995 00996 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00997
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: