inf_fe.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 00020 // Local includes 00021 #include "libmesh/libmesh_config.h" 00022 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00023 #include "libmesh/inf_fe.h" 00024 #include "libmesh/quadrature_gauss.h" 00025 #include "libmesh/elem.h" 00026 #include "libmesh/libmesh_logging.h" 00027 00028 namespace libMesh 00029 { 00030 00031 00032 00033 // ------------------------------------------------------------ 00034 // InfFE class members 00035 00036 00037 00038 // Constructor 00039 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00040 InfFE<Dim,T_radial,T_map>::InfFE (const FEType& fet) : 00041 FEBase (Dim, fet), 00042 00043 _n_total_approx_sf (0), 00044 _n_total_qp (0), 00045 00046 base_qrule (NULL), 00047 radial_qrule (NULL), 00048 base_elem (NULL), 00049 base_fe (NULL), 00050 00051 // initialize the current_fe_type to all the same 00052 // values as \p fet (since the FE families and coordinate 00053 // map type should @e not change), but use an invalid order 00054 // for the radial part (since this is the only order 00055 // that may change!). 00056 // the data structures like \p phi etc are not initialized 00057 // through the constructor, but throught reinit() 00058 current_fe_type ( FEType(fet.order, 00059 fet.family, 00060 INVALID_ORDER, 00061 fet.radial_family, 00062 fet.inf_map) ) 00063 00064 { 00065 // Sanity checks 00066 libmesh_assert_equal_to (T_radial, fe_type.radial_family); 00067 libmesh_assert_equal_to (T_map, fe_type.inf_map); 00068 00069 // build the base_fe object, handle the AutoPtr 00070 if (Dim != 1) 00071 { 00072 AutoPtr<FEBase> ap_fb(FEBase::build(Dim-1, fet)); 00073 base_fe = ap_fb.release(); 00074 } 00075 } 00076 00077 00078 00079 00080 // Destructor 00081 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00082 InfFE<Dim,T_radial,T_map>::~InfFE () 00083 { 00084 // delete pointers, if necessary 00085 if (base_qrule != NULL) 00086 { 00087 delete base_qrule; 00088 base_qrule = NULL; 00089 } 00090 00091 if (radial_qrule != NULL) 00092 { 00093 delete radial_qrule; 00094 radial_qrule = NULL; 00095 } 00096 00097 if (base_elem != NULL) 00098 { 00099 delete base_elem; 00100 base_elem = NULL; 00101 } 00102 00103 if (base_fe != NULL) 00104 { 00105 delete base_fe; 00106 base_fe = NULL; 00107 } 00108 } 00109 00110 00111 00112 00113 00114 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00115 void InfFE<Dim,T_radial,T_map>:: attach_quadrature_rule (QBase* q) 00116 { 00117 libmesh_assert(q); 00118 libmesh_assert(base_fe); 00119 00120 const Order base_int_order = q->get_order(); 00121 const Order radial_int_order = static_cast<Order>(2 * (static_cast<unsigned int>(fe_type.radial_order) + 1) +2); 00122 const unsigned int qrule_dim = q->get_dim(); 00123 00124 if (Dim != 1) 00125 { 00126 // build a Dim-1 quadrature rule of the type that we received 00127 AutoPtr<QBase> apq( QBase::build(q->type(), qrule_dim-1, base_int_order) ); 00128 base_qrule = apq.release(); 00129 base_fe->attach_quadrature_rule(base_qrule); 00130 } 00131 00132 // in radial direction, always use Gauss quadrature 00133 radial_qrule = new QGauss(1, radial_int_order); 00134 00135 // currently not used. But maybe helpful to store the QBase* 00136 // with which we initialized our own quadrature rules 00137 qrule = q; 00138 } 00139 00140 00141 00142 00143 00144 template <unsigned int Dim, FEFamily T_radial, InfMapType T_base> 00145 void InfFE<Dim,T_radial,T_base>::update_base_elem (const Elem* inf_elem) 00146 { 00147 if (base_elem != NULL) 00148 delete base_elem; 00149 base_elem = Base::build_elem(inf_elem); 00150 } 00151 00152 00153 00154 00155 00156 00157 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00158 void InfFE<Dim,T_radial,T_map>::reinit(const Elem* inf_elem, 00159 const std::vector<Point>* const pts, 00160 const std::vector<Real>* const weights) 00161 { 00162 libmesh_assert(base_fe); 00163 libmesh_assert(base_fe->qrule); 00164 libmesh_assert_equal_to (base_fe->qrule, base_qrule); 00165 libmesh_assert(radial_qrule); 00166 libmesh_assert(inf_elem); 00167 00168 if (pts == NULL) 00169 { 00170 bool init_shape_functions_required = false; 00171 00172 // ----------------------------------------------------------------- 00173 // init the radial data fields only when the radial order changes 00174 if (current_fe_type.radial_order != fe_type.radial_order) 00175 { 00176 current_fe_type.radial_order = fe_type.radial_order; 00177 00178 // Watch out: this call to QBase->init() only works for 00179 // current_fe_type = const! To allow variable Order, 00180 // the init() of QBase has to be modified... 00181 radial_qrule->init(EDGE2); 00182 00183 // initialize the radial shape functions 00184 this->init_radial_shape_functions(inf_elem); 00185 00186 init_shape_functions_required=true; 00187 } 00188 00189 00190 bool update_base_elem_required=true; 00191 00192 // ----------------------------------------------------------------- 00193 // update the type in accordance to the current cell 00194 // and reinit if the cell type has changed or (as in 00195 // the case of the hierarchics) the shape functions 00196 // depend on the particular element and need a reinit 00197 if ( ( Dim != 1) && 00198 ( (this->get_type() != inf_elem->type()) || 00199 (base_fe->shapes_need_reinit()) ) ) 00200 { 00201 // store the new element type, update base_elem 00202 // here. Through \p update_base_elem_required, 00203 // remember whether it has to be updated (see below). 00204 elem_type = inf_elem->type(); 00205 this->update_base_elem(inf_elem); 00206 update_base_elem_required=false; 00207 00208 // initialize the base quadrature rule for the new element 00209 base_qrule->init(base_elem->type()); 00210 00211 // initialize the shape functions in the base 00212 base_fe->init_base_shape_functions(base_fe->qrule->get_points(), 00213 base_elem); 00214 00215 init_shape_functions_required=true; 00216 } 00217 00218 00219 // when either the radial or base part change, 00220 // we have to init the whole fields 00221 if (init_shape_functions_required) 00222 this->init_shape_functions (inf_elem); 00223 00224 // computing the distance only works when we have the current 00225 // base_elem stored. This happens when fe_type is const, 00226 // the inf_elem->type remains the same. Then we have to 00227 // update the base elem _here_. 00228 if (update_base_elem_required) 00229 this->update_base_elem(inf_elem); 00230 00231 // compute dist (depends on geometry, therefore has to be updated for 00232 // each and every new element), throw radial and base part together 00233 this->combine_base_radial (inf_elem); 00234 00235 this->_fe_map->compute_map (this->dim,_total_qrule_weights, inf_elem); 00236 00237 // Compute the shape functions and the derivatives 00238 // at all quadrature points. 00239 this->compute_shape_functions (inf_elem,base_fe->qrule->get_points()); 00240 } 00241 00242 else // if pts != NULL 00243 { 00244 // update the elem_type 00245 elem_type = inf_elem->type(); 00246 00247 // init radial shapes 00248 this->init_radial_shape_functions(inf_elem); 00249 00250 // update the base 00251 this->update_base_elem(inf_elem); 00252 00253 // the finite element on the ifem base 00254 { 00255 AutoPtr<FEBase> ap_fb(FEBase::build(Dim-1, this->fe_type)); 00256 if (base_fe != NULL) 00257 delete base_fe; 00258 base_fe = ap_fb.release(); 00259 } 00260 00261 // inite base shapes 00262 base_fe->init_base_shape_functions(*pts, 00263 base_elem); 00264 00265 this->init_shape_functions (inf_elem); 00266 00267 // combine the base and radial shapes 00268 this->combine_base_radial (inf_elem); 00269 00270 // weights 00271 if (weights != NULL) 00272 { 00273 this->_fe_map->compute_map (this->dim, *weights, inf_elem); 00274 } 00275 else 00276 { 00277 std::vector<Real> dummy_weights (pts->size(), 1.); 00278 this->_fe_map->compute_map (this->dim, dummy_weights, inf_elem); 00279 } 00280 00281 // finally compute the ifem shapes 00282 this->compute_shape_functions (inf_elem,*pts); 00283 } 00284 00285 } 00286 00287 00288 00289 00290 00291 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00292 void InfFE<Dim,T_radial,T_map>::init_radial_shape_functions(const Elem* libmesh_dbg_var(inf_elem)) 00293 { 00294 libmesh_assert(radial_qrule); 00295 libmesh_assert(inf_elem); 00296 00297 00301 START_LOG("init_radial_shape_functions()", "InfFE"); 00302 00303 00304 // ----------------------------------------------------------------- 00305 // initialize most of the things related to mapping 00306 00307 // The order to use in the radial map (currently independent of the element type) 00308 const Order radial_mapping_order (Radial::mapping_order()); 00309 const unsigned int n_radial_mapping_shape_functions (Radial::n_dofs(radial_mapping_order)); 00310 00311 00312 00313 // ----------------------------------------------------------------- 00314 // initialize most of the things related to physical approximation 00315 00316 const Order radial_approx_order (fe_type.radial_order); 00317 const unsigned int n_radial_approx_shape_functions (Radial::n_dofs(radial_approx_order)); 00318 00319 const unsigned int n_radial_qp = radial_qrule->n_points(); 00320 const std::vector<Point>& radial_qp = radial_qrule->get_points(); 00321 00322 00323 00324 // ----------------------------------------------------------------- 00325 // resize the radial data fields 00326 00327 mode.resize (n_radial_approx_shape_functions); // the radial polynomials (eval) 00328 dmodedv.resize (n_radial_approx_shape_functions); 00329 00330 som.resize (n_radial_qp); // the (1-v)/2 weight 00331 dsomdv.resize (n_radial_qp); 00332 00333 radial_map.resize (n_radial_mapping_shape_functions); // the radial map 00334 dradialdv_map.resize (n_radial_mapping_shape_functions); 00335 00336 00337 for (unsigned int i=0; i<n_radial_mapping_shape_functions; i++) 00338 { 00339 radial_map[i].resize (n_radial_qp); 00340 dradialdv_map[i].resize (n_radial_qp); 00341 } 00342 00343 00344 for (unsigned int i=0; i<n_radial_approx_shape_functions; i++) 00345 { 00346 mode[i].resize (n_radial_qp); 00347 dmodedv[i].resize (n_radial_qp); 00348 } 00349 00350 00351 // compute scalar values at radial quadrature points 00352 for (unsigned int p=0; p<n_radial_qp; p++) 00353 { 00354 som[p] = Radial::decay (radial_qp[p](0)); 00355 dsomdv[p] = Radial::decay_deriv (radial_qp[p](0)); 00356 } 00357 00358 00359 // evaluate the mode shapes in radial direction at radial quadrature points 00360 for (unsigned int i=0; i<n_radial_approx_shape_functions; i++) 00361 for (unsigned int p=0; p<n_radial_qp; p++) 00362 { 00363 mode[i][p] = InfFE<Dim,T_radial,T_map>::eval (radial_qp[p](0), radial_approx_order, i); 00364 dmodedv[i][p] = InfFE<Dim,T_radial,T_map>::eval_deriv (radial_qp[p](0), radial_approx_order, i); 00365 } 00366 00367 00368 // evaluate the mapping functions in radial direction at radial quadrature points 00369 for (unsigned int i=0; i<n_radial_mapping_shape_functions; i++) 00370 for (unsigned int p=0; p<n_radial_qp; p++) 00371 { 00372 radial_map[i][p] = InfFE<Dim,INFINITE_MAP,T_map>::eval (radial_qp[p](0), radial_mapping_order, i); 00373 dradialdv_map[i][p] = InfFE<Dim,INFINITE_MAP,T_map>::eval_deriv (radial_qp[p](0), radial_mapping_order, i); 00374 } 00375 00379 STOP_LOG("init_radial_shape_functions()", "InfFE"); 00380 00381 } 00382 00383 00384 00385 00386 00387 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00388 void InfFE<Dim,T_radial,T_map>::init_shape_functions(const Elem* inf_elem) 00389 { 00390 libmesh_assert(inf_elem); 00391 00392 00393 // Start logging the radial shape function initialization 00394 START_LOG("init_shape_functions()", "InfFE"); 00395 00396 00397 // ----------------------------------------------------------------- 00398 // fast access to some const int's for the radial data 00399 const unsigned int n_radial_mapping_sf = 00400 libmesh_cast_int<unsigned int>(radial_map.size()); 00401 const unsigned int n_radial_approx_sf = 00402 libmesh_cast_int<unsigned int>(mode.size()); 00403 const unsigned int n_radial_qp = 00404 libmesh_cast_int<unsigned int>(som.size()); 00405 00406 00407 // ----------------------------------------------------------------- 00408 // initialize most of the things related to mapping 00409 00410 // The element type and order to use in the base map 00411 const Order base_mapping_order ( base_elem->default_order() ); 00412 const ElemType base_mapping_elem_type ( base_elem->type() ); 00413 00414 // the number of base shape functions used to construct the map 00415 // (Lagrange shape functions are used for mapping in the base) 00416 unsigned int n_base_mapping_shape_functions = Base::n_base_mapping_sf(base_mapping_elem_type, 00417 base_mapping_order); 00418 00419 const unsigned int n_total_mapping_shape_functions = 00420 n_radial_mapping_sf * n_base_mapping_shape_functions; 00421 00422 00423 00424 // ----------------------------------------------------------------- 00425 // initialize most of the things related to physical approximation 00426 00427 unsigned int n_base_approx_shape_functions; 00428 if (Dim > 1) 00429 n_base_approx_shape_functions = base_fe->n_shape_functions(); 00430 else 00431 n_base_approx_shape_functions = 1; 00432 00433 00434 const unsigned int n_total_approx_shape_functions = 00435 n_radial_approx_sf * n_base_approx_shape_functions; 00436 00437 // update class member field 00438 _n_total_approx_sf = n_total_approx_shape_functions; 00439 00440 00441 // The number of the base quadrature points. 00442 const unsigned int n_base_qp = base_qrule->n_points(); 00443 00444 // The total number of quadrature points. 00445 const unsigned int n_total_qp = n_radial_qp * n_base_qp; 00446 00447 00448 // update class member field 00449 _n_total_qp = n_total_qp; 00450 00451 00452 00453 // ----------------------------------------------------------------- 00454 // initialize the node and shape numbering maps 00455 { 00456 // these vectors work as follows: the i-th entry stores 00457 // the associated base/radial node number 00458 _radial_node_index.resize (n_total_mapping_shape_functions); 00459 _base_node_index.resize (n_total_mapping_shape_functions); 00460 00461 // similar for the shapes: the i-th entry stores 00462 // the associated base/radial shape number 00463 _radial_shape_index.resize (n_total_approx_shape_functions); 00464 _base_shape_index.resize (n_total_approx_shape_functions); 00465 00466 const ElemType inf_elem_type (inf_elem->type()); 00467 00468 // fill the node index map 00469 for (unsigned int n=0; n<n_total_mapping_shape_functions; n++) 00470 { 00471 compute_node_indices (inf_elem_type, 00472 n, 00473 _base_node_index[n], 00474 _radial_node_index[n]); 00475 libmesh_assert_less (_base_node_index[n], n_base_mapping_shape_functions); 00476 libmesh_assert_less (_radial_node_index[n], n_radial_mapping_sf); 00477 } 00478 00479 // fill the shape index map 00480 for (unsigned int n=0; n<n_total_approx_shape_functions; n++) 00481 { 00482 compute_shape_indices (this->fe_type, 00483 inf_elem_type, 00484 n, 00485 _base_shape_index[n], 00486 _radial_shape_index[n]); 00487 libmesh_assert_less (_base_shape_index[n], n_base_approx_shape_functions); 00488 libmesh_assert_less (_radial_shape_index[n], n_radial_approx_sf); 00489 } 00490 } 00491 00492 00493 00494 00495 00496 // ----------------------------------------------------------------- 00497 // resize the base data fields 00498 dist.resize(n_base_mapping_shape_functions); 00499 00500 00501 00502 // ----------------------------------------------------------------- 00503 // resize the total data fields 00504 00505 // the phase term varies with xi, eta and zeta(v): store it for _all_ qp 00506 // 00507 // when computing the phase, we need the base approximations 00508 // therefore, initialize the phase here, but evaluate it 00509 // in combine_base_radial(). 00510 // 00511 // the weight, though, is only needed at the radial quadrature points, n_radial_qp. 00512 // but for a uniform interface to the protected data fields 00513 // the weight data field (which are accessible from the outside) are expanded to n_total_qp. 00514 weight.resize (n_total_qp); 00515 dweightdv.resize (n_total_qp); 00516 dweight.resize (n_total_qp); 00517 00518 dphase.resize (n_total_qp); 00519 dphasedxi.resize (n_total_qp); 00520 dphasedeta.resize (n_total_qp); 00521 dphasedzeta.resize (n_total_qp); 00522 00523 // this vector contains the integration weights for the combined quadrature rule 00524 _total_qrule_weights.resize(n_total_qp); 00525 00526 00527 // ----------------------------------------------------------------- 00528 // InfFE's data fields phi, dphi, dphidx, phi_map etc hold the _total_ 00529 // shape and mapping functions, respectively 00530 { 00531 phi.resize (n_total_approx_shape_functions); 00532 dphi.resize (n_total_approx_shape_functions); 00533 dphidx.resize (n_total_approx_shape_functions); 00534 dphidy.resize (n_total_approx_shape_functions); 00535 dphidz.resize (n_total_approx_shape_functions); 00536 dphidxi.resize (n_total_approx_shape_functions); 00537 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00538 libmesh_do_once(libMesh::err << "Second derivatives for Infinite elements" 00539 << " are not yet implemented!" 00540 << std::endl); 00541 00542 d2phi.resize (n_total_approx_shape_functions); 00543 d2phidx2.resize (n_total_approx_shape_functions); 00544 d2phidxdy.resize (n_total_approx_shape_functions); 00545 d2phidxdz.resize (n_total_approx_shape_functions); 00546 d2phidy2.resize (n_total_approx_shape_functions); 00547 d2phidydz.resize (n_total_approx_shape_functions); 00548 d2phidz2.resize (n_total_approx_shape_functions); 00549 d2phidxi2.resize (n_total_approx_shape_functions); 00550 00551 if (Dim > 1) 00552 { 00553 d2phidxideta.resize (n_total_approx_shape_functions); 00554 d2phideta2.resize (n_total_approx_shape_functions); 00555 } 00556 00557 if (Dim > 2) 00558 { 00559 d2phidetadzeta.resize (n_total_approx_shape_functions); 00560 d2phidxidzeta.resize (n_total_approx_shape_functions); 00561 d2phidzeta2.resize (n_total_approx_shape_functions); 00562 } 00563 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00564 00565 if (Dim > 1) 00566 dphideta.resize (n_total_approx_shape_functions); 00567 00568 if (Dim == 3) 00569 dphidzeta.resize (n_total_approx_shape_functions); 00570 00571 00572 00573 std::vector<std::vector<Real> >& phi_map = this->_fe_map->get_phi_map(); 00574 std::vector<std::vector<Real> >& dphidxi_map = this->_fe_map->get_dphidxi_map(); 00575 00576 phi_map.resize (n_total_mapping_shape_functions); 00577 dphidxi_map.resize (n_total_mapping_shape_functions); 00578 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00579 std::vector<std::vector<Real> >& d2phidxi2_map = this->_fe_map->get_d2phidxi2_map(); 00580 d2phidxi2_map.resize (n_total_mapping_shape_functions); 00581 00582 if (Dim > 1) 00583 { 00584 std::vector<std::vector<Real> >& d2phidxideta_map = this->_fe_map->get_d2phidxideta_map(); 00585 std::vector<std::vector<Real> >& d2phideta2_map = this->_fe_map->get_d2phideta2_map(); 00586 d2phidxideta_map.resize (n_total_mapping_shape_functions); 00587 d2phideta2_map.resize (n_total_mapping_shape_functions); 00588 } 00589 00590 if (Dim == 3) 00591 { 00592 std::vector<std::vector<Real> >& d2phidxidzeta_map = this->_fe_map->get_d2phidxidzeta_map(); 00593 std::vector<std::vector<Real> >& d2phidetadzeta_map = this->_fe_map->get_d2phidetadzeta_map(); 00594 std::vector<std::vector<Real> >& d2phidzeta2_map = this->_fe_map->get_d2phidzeta2_map(); 00595 d2phidxidzeta_map.resize (n_total_mapping_shape_functions); 00596 d2phidetadzeta_map.resize (n_total_mapping_shape_functions); 00597 d2phidzeta2_map.resize (n_total_mapping_shape_functions); 00598 } 00599 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00600 00601 if (Dim > 1) 00602 { 00603 std::vector<std::vector<Real> >& dphideta_map = this->_fe_map->get_dphideta_map(); 00604 dphideta_map.resize (n_total_mapping_shape_functions); 00605 } 00606 00607 if (Dim == 3) 00608 { 00609 std::vector<std::vector<Real> >& dphidzeta_map = this->_fe_map->get_dphidzeta_map(); 00610 dphidzeta_map.resize (n_total_mapping_shape_functions); 00611 } 00612 } 00613 00614 00615 00616 // ----------------------------------------------------------------- 00617 // collect all the for loops, where inner vectors are 00618 // resized to the appropriate number of quadrature points 00619 { 00620 for (unsigned int i=0; i<n_total_approx_shape_functions; i++) 00621 { 00622 phi[i].resize (n_total_qp); 00623 dphi[i].resize (n_total_qp); 00624 dphidx[i].resize (n_total_qp); 00625 dphidy[i].resize (n_total_qp); 00626 dphidz[i].resize (n_total_qp); 00627 dphidxi[i].resize (n_total_qp); 00628 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00629 d2phi[i].resize (n_total_qp); 00630 d2phidx2[i].resize (n_total_qp); 00631 d2phidxdy[i].resize (n_total_qp); 00632 d2phidxdz[i].resize (n_total_qp); 00633 d2phidy2[i].resize (n_total_qp); 00634 d2phidydz[i].resize (n_total_qp); 00635 d2phidy2[i].resize (n_total_qp); 00636 d2phidxi2[i].resize (n_total_qp); 00637 00638 if (Dim > 1) 00639 { 00640 d2phidxideta[i].resize (n_total_qp); 00641 d2phideta2[i].resize (n_total_qp); 00642 } 00643 if (Dim > 2) 00644 { 00645 d2phidxidzeta[i].resize (n_total_qp); 00646 d2phidetadzeta[i].resize (n_total_qp); 00647 d2phidzeta2[i].resize (n_total_qp); 00648 } 00649 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00650 00651 if (Dim > 1) 00652 dphideta[i].resize (n_total_qp); 00653 00654 if (Dim == 3) 00655 dphidzeta[i].resize (n_total_qp); 00656 00657 } 00658 00659 for (unsigned int i=0; i<n_total_mapping_shape_functions; i++) 00660 { 00661 std::vector<std::vector<Real> >& phi_map = this->_fe_map->get_phi_map(); 00662 std::vector<std::vector<Real> >& dphidxi_map = this->_fe_map->get_dphidxi_map(); 00663 phi_map[i].resize (n_total_qp); 00664 dphidxi_map[i].resize (n_total_qp); 00665 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00666 std::vector<std::vector<Real> >& d2phidxi2_map = this->_fe_map->get_d2phidxi2_map(); 00667 d2phidxi2_map[i].resize (n_total_qp); 00668 if (Dim > 1) 00669 { 00670 std::vector<std::vector<Real> >& d2phidxideta_map = this->_fe_map->get_d2phidxideta_map(); 00671 std::vector<std::vector<Real> >& d2phideta2_map = this->_fe_map->get_d2phideta2_map(); 00672 d2phidxideta_map[i].resize (n_total_qp); 00673 d2phideta2_map[i].resize (n_total_qp); 00674 } 00675 00676 if (Dim > 2) 00677 { 00678 std::vector<std::vector<Real> >& d2phidxidzeta_map = this->_fe_map->get_d2phidxidzeta_map(); 00679 std::vector<std::vector<Real> >& d2phidetadzeta_map = this->_fe_map->get_d2phidetadzeta_map(); 00680 std::vector<std::vector<Real> >& d2phidzeta2_map = this->_fe_map->get_d2phidzeta2_map(); 00681 d2phidxidzeta_map[i].resize (n_total_qp); 00682 d2phidetadzeta_map[i].resize (n_total_qp); 00683 d2phidzeta2_map[i].resize (n_total_qp); 00684 } 00685 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00686 00687 if (Dim > 1) 00688 { 00689 std::vector<std::vector<Real> >& dphideta_map = this->_fe_map->get_dphideta_map(); 00690 dphideta_map[i].resize (n_total_qp); 00691 } 00692 00693 if (Dim == 3) 00694 { 00695 std::vector<std::vector<Real> >& dphidzeta_map = this->_fe_map->get_dphidzeta_map(); 00696 dphidzeta_map[i].resize (n_total_qp); 00697 } 00698 } 00699 } 00700 00701 00702 00703 { 00704 // ----------------------------------------------------------------- 00705 // (a) compute scalar values at _all_ quadrature points -- for uniform 00706 // access from the outside to these fields 00707 // (b) form a std::vector<Real> which contains the appropriate weights 00708 // of the combined quadrature rule! 00709 const std::vector<Point>& radial_qp = radial_qrule->get_points(); 00710 libmesh_assert_equal_to (radial_qp.size(), n_radial_qp); 00711 00712 const std::vector<Real>& radial_qw = radial_qrule->get_weights(); 00713 const std::vector<Real>& base_qw = base_qrule->get_weights(); 00714 libmesh_assert_equal_to (radial_qw.size(), n_radial_qp); 00715 libmesh_assert_equal_to (base_qw.size(), n_base_qp); 00716 00717 for (unsigned int rp=0; rp<n_radial_qp; rp++) 00718 for (unsigned int bp=0; bp<n_base_qp; bp++) 00719 { 00720 weight [ bp+rp*n_base_qp ] = Radial::D (radial_qp[rp](0)); 00721 dweightdv[ bp+rp*n_base_qp ] = Radial::D_deriv (radial_qp[rp](0)); 00722 00723 _total_qrule_weights[ bp+rp*n_base_qp ] = radial_qw[rp] * base_qw[bp]; 00724 } 00725 } 00726 00727 00731 STOP_LOG("init_shape_functions()", "InfFE"); 00732 00733 } 00734 00735 00736 00737 00738 00739 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00740 void InfFE<Dim,T_radial,T_map>::combine_base_radial(const Elem* inf_elem) 00741 { 00742 libmesh_assert(inf_elem); 00743 // at least check whether the base element type is correct. 00744 // otherwise this version of computing dist would give problems 00745 libmesh_assert_equal_to (base_elem->type(), Base::get_elem_type(inf_elem->type())); 00746 00747 00751 START_LOG("combine_base_radial()", "InfFE"); 00752 00753 00754 // zero the phase, since it is to be summed up 00755 std::fill (dphasedxi.begin(), dphasedxi.end(), 0.); 00756 std::fill (dphasedeta.begin(), dphasedeta.end(), 0.); 00757 std::fill (dphasedzeta.begin(), dphasedzeta.end(), 0.); 00758 00759 00760 const unsigned int n_base_mapping_sf = 00761 libmesh_cast_int<unsigned int>(dist.size()); 00762 const Point origin = inf_elem->origin(); 00763 00764 // for each new infinite element, compute the radial distances 00765 for (unsigned int n=0; n<n_base_mapping_sf; n++) 00766 dist[n] = Point(base_elem->point(n) - origin).size(); 00767 00768 00769 switch (Dim) 00770 { 00771 00772 //------------------------------------------------------------ 00773 // 1D 00774 case 1: 00775 { 00776 libmesh_not_implemented(); 00777 break; 00778 } 00779 00780 00781 00782 //------------------------------------------------------------ 00783 // 2D 00784 case 2: 00785 { 00786 libmesh_not_implemented(); 00787 break; 00788 } 00789 00790 00791 00792 //------------------------------------------------------------ 00793 // 3D 00794 case 3: 00795 { 00796 // fast access to the approximation and mapping shapes of base_fe 00797 const std::vector<std::vector<Real> >& S = base_fe->phi; 00798 const std::vector<std::vector<Real> >& Ss = base_fe->dphidxi; 00799 const std::vector<std::vector<Real> >& St = base_fe->dphideta; 00800 const std::vector<std::vector<Real> >& S_map = (base_fe->get_fe_map()).get_phi_map(); 00801 const std::vector<std::vector<Real> >& Ss_map = (base_fe->get_fe_map()).get_dphidxi_map(); 00802 const std::vector<std::vector<Real> >& St_map = (base_fe->get_fe_map()).get_dphideta_map(); 00803 00804 const unsigned int n_radial_qp = radial_qrule->n_points(); 00805 const unsigned int n_base_qp = base_qrule-> n_points(); 00806 00807 const unsigned int n_total_mapping_sf = 00808 libmesh_cast_int<unsigned int>(radial_map.size()) * n_base_mapping_sf; 00809 00810 const unsigned int n_total_approx_sf = Radial::n_dofs(fe_type.radial_order) * base_fe->n_shape_functions(); 00811 00812 00813 // compute the phase term derivatives 00814 { 00815 unsigned int tp=0; 00816 for (unsigned int rp=0; rp<n_radial_qp; rp++) // over radial qp's 00817 for (unsigned int bp=0; bp<n_base_qp; bp++) // over base qp's 00818 { 00819 // sum over all base shapes, to get the average distance 00820 for (unsigned int i=0; i<n_base_mapping_sf; i++) 00821 { 00822 dphasedxi[tp] += Ss_map[i][bp] * dist[i] * radial_map [1][rp]; 00823 dphasedeta[tp] += St_map[i][bp] * dist[i] * radial_map [1][rp]; 00824 dphasedzeta[tp] += S_map [i][bp] * dist[i] * dradialdv_map[1][rp]; 00825 } 00826 00827 tp++; 00828 00829 } // loop radial and base qp's 00830 00831 } 00832 00833 libmesh_assert_equal_to (phi.size(), n_total_approx_sf); 00834 libmesh_assert_equal_to (dphidxi.size(), n_total_approx_sf); 00835 libmesh_assert_equal_to (dphideta.size(), n_total_approx_sf); 00836 libmesh_assert_equal_to (dphidzeta.size(), n_total_approx_sf); 00837 00838 // compute the overall approximation shape functions, 00839 // pick the appropriate radial and base shapes through using 00840 // _base_shape_index and _radial_shape_index 00841 for (unsigned int rp=0; rp<n_radial_qp; rp++) // over radial qp's 00842 for (unsigned int bp=0; bp<n_base_qp; bp++) // over base qp's 00843 for (unsigned int ti=0; ti<n_total_approx_sf; ti++) // over _all_ approx_sf 00844 { 00845 // let the index vectors take care of selecting the appropriate base/radial shape 00846 const unsigned int bi = _base_shape_index [ti]; 00847 const unsigned int ri = _radial_shape_index[ti]; 00848 phi [ti][bp+rp*n_base_qp] = S [bi][bp] * mode[ri][rp] * som[rp]; 00849 dphidxi [ti][bp+rp*n_base_qp] = Ss[bi][bp] * mode[ri][rp] * som[rp]; 00850 dphideta [ti][bp+rp*n_base_qp] = St[bi][bp] * mode[ri][rp] * som[rp]; 00851 dphidzeta[ti][bp+rp*n_base_qp] = S [bi][bp] 00852 * (dmodedv[ri][rp] * som[rp] + mode[ri][rp] * dsomdv[rp]); 00853 } 00854 00855 std::vector<std::vector<Real> >& phi_map = this->_fe_map->get_phi_map(); 00856 std::vector<std::vector<Real> >& dphidxi_map = this->_fe_map->get_dphidxi_map(); 00857 std::vector<std::vector<Real> >& dphideta_map = this->_fe_map->get_dphideta_map(); 00858 std::vector<std::vector<Real> >& dphidzeta_map = this->_fe_map->get_dphidzeta_map(); 00859 00860 libmesh_assert_equal_to (phi_map.size(), n_total_mapping_sf); 00861 libmesh_assert_equal_to (dphidxi_map.size(), n_total_mapping_sf); 00862 libmesh_assert_equal_to (dphideta_map.size(), n_total_mapping_sf); 00863 libmesh_assert_equal_to (dphidzeta_map.size(), n_total_mapping_sf); 00864 00865 // compute the overall mapping functions, 00866 // pick the appropriate radial and base entries through using 00867 // _base_node_index and _radial_node_index 00868 for (unsigned int rp=0; rp<n_radial_qp; rp++) // over radial qp's 00869 for (unsigned int bp=0; bp<n_base_qp; bp++) // over base qp's 00870 for (unsigned int ti=0; ti<n_total_mapping_sf; ti++) // over all mapping shapes 00871 { 00872 // let the index vectors take care of selecting the appropriate base/radial mapping shape 00873 const unsigned int bi = _base_node_index [ti]; 00874 const unsigned int ri = _radial_node_index[ti]; 00875 phi_map [ti][bp+rp*n_base_qp] = S_map [bi][bp] * radial_map [ri][rp]; 00876 dphidxi_map [ti][bp+rp*n_base_qp] = Ss_map[bi][bp] * radial_map [ri][rp]; 00877 dphideta_map [ti][bp+rp*n_base_qp] = St_map[bi][bp] * radial_map [ri][rp]; 00878 dphidzeta_map[ti][bp+rp*n_base_qp] = S_map [bi][bp] * dradialdv_map[ri][rp]; 00879 } 00880 00881 00882 break; 00883 } 00884 00885 00886 default: 00887 libmesh_error(); 00888 } 00889 00890 00894 STOP_LOG("combine_base_radial()", "InfFE"); 00895 00896 } 00897 00898 00899 00900 00901 00902 00903 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00904 void InfFE<Dim,T_radial,T_map>::compute_shape_functions(const Elem*, const std::vector<Point>&) 00905 { 00906 libmesh_assert(radial_qrule); 00907 00908 00909 00910 // Start logging the overall computation of shape functions 00911 START_LOG("compute_shape_functions()", "InfFE"); 00912 00913 00914 const unsigned int n_total_qp = _n_total_qp; 00915 00916 00917 //------------------------------------------------------------------------- 00918 // Compute the shape function values (and derivatives) 00919 // at the Quadrature points. Note that the actual values 00920 // have already been computed via init_shape_functions 00921 00922 // Compute the value of the derivative shape function i at quadrature point p 00923 switch (dim) 00924 { 00925 00926 case 1: 00927 { 00928 libmesh_not_implemented(); 00929 break; 00930 } 00931 00932 case 2: 00933 { 00934 libmesh_not_implemented(); 00935 break; 00936 } 00937 00938 case 3: 00939 { 00940 const std::vector<Real>& dxidx_map = this->_fe_map->get_dxidx(); 00941 const std::vector<Real>& dxidy_map = this->_fe_map->get_dxidy(); 00942 const std::vector<Real>& dxidz_map = this->_fe_map->get_dxidz(); 00943 00944 const std::vector<Real>& detadx_map = this->_fe_map->get_detadx(); 00945 const std::vector<Real>& detady_map = this->_fe_map->get_detady(); 00946 const std::vector<Real>& detadz_map = this->_fe_map->get_detadz(); 00947 00948 const std::vector<Real>& dzetadx_map = this->_fe_map->get_dzetadx(); 00949 const std::vector<Real>& dzetady_map = this->_fe_map->get_dzetady(); 00950 const std::vector<Real>& dzetadz_map = this->_fe_map->get_dzetadz(); 00951 00952 // These are _all_ shape functions of this infinite element 00953 for (unsigned int i=0; i<phi.size(); i++) 00954 for (unsigned int p=0; p<n_total_qp; p++) 00955 { 00956 // dphi/dx = (dphi/dxi)*(dxi/dx) + (dphi/deta)*(deta/dx) + (dphi/dzeta)*(dzeta/dx); 00957 dphi[i][p](0) = 00958 dphidx[i][p] = (dphidxi[i][p]*dxidx_map[p] + 00959 dphideta[i][p]*detadx_map[p] + 00960 dphidzeta[i][p]*dzetadx_map[p]); 00961 00962 // dphi/dy = (dphi/dxi)*(dxi/dy) + (dphi/deta)*(deta/dy) + (dphi/dzeta)*(dzeta/dy); 00963 dphi[i][p](1) = 00964 dphidy[i][p] = (dphidxi[i][p]*dxidy_map[p] + 00965 dphideta[i][p]*detady_map[p] + 00966 dphidzeta[i][p]*dzetady_map[p]); 00967 00968 // dphi/dz = (dphi/dxi)*(dxi/dz) + (dphi/deta)*(deta/dz) + (dphi/dzeta)*(dzeta/dz); 00969 dphi[i][p](2) = 00970 dphidz[i][p] = (dphidxi[i][p]*dxidz_map[p] + 00971 dphideta[i][p]*detadz_map[p] + 00972 dphidzeta[i][p]*dzetadz_map[p]); 00973 } 00974 00975 00976 // This is the derivative of the phase term of this infinite element 00977 for (unsigned int p=0; p<n_total_qp; p++) 00978 { 00979 // the derivative of the phase term 00980 dphase[p](0) = (dphasedxi[p] * dxidx_map[p] + 00981 dphasedeta[p] * detadx_map[p] + 00982 dphasedzeta[p] * dzetadx_map[p]); 00983 00984 dphase[p](1) = (dphasedxi[p] * dxidy_map[p] + 00985 dphasedeta[p] * detady_map[p] + 00986 dphasedzeta[p] * dzetady_map[p]); 00987 00988 dphase[p](2) = (dphasedxi[p] * dxidz_map[p] + 00989 dphasedeta[p] * detadz_map[p] + 00990 dphasedzeta[p] * dzetadz_map[p]); 00991 00992 // the derivative of the radial weight - varies only in radial direction, 00993 // therefore dweightdxi = dweightdeta = 0. 00994 dweight[p](0) = dweightdv[p] * dzetadx_map[p]; 00995 00996 dweight[p](1) = dweightdv[p] * dzetady_map[p]; 00997 00998 dweight[p](2) = dweightdv[p] * dzetadz_map[p]; 00999 01000 } 01001 01002 break; 01003 } 01004 01005 01006 01007 default: 01008 { 01009 libmesh_error(); 01010 } 01011 } 01012 01013 01014 01015 // Stop logging the overall computation of shape functions 01016 STOP_LOG("compute_shape_functions()", "InfFE"); 01017 01018 } 01019 01020 01021 01022 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 01023 bool InfFE<Dim,T_radial,T_map>::shapes_need_reinit() const 01024 { 01025 return false; 01026 } 01027 01028 } // namespace libMesh 01029 01030 01031 //-------------------------------------------------------------- 01032 // Explicit instantiations 01033 #include "libmesh/inf_fe_instantiate_1D.h" 01034 #include "libmesh/inf_fe_instantiate_2D.h" 01035 #include "libmesh/inf_fe_instantiate_3D.h" 01036 01037 01038 01039 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 01040
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: