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/elem.h" 00022 #include "libmesh/fe.h" 00023 #include "libmesh/fe_interface.h" 00024 #include "libmesh/fe_macro.h" 00025 #include "libmesh/libmesh_logging.h" 00026 #include "libmesh/quadrature.h" 00027 #include "libmesh/tensor_value.h" 00028 00029 namespace libMesh 00030 { 00031 00032 00033 // ------------------------------------------------------------ 00034 // FE class members 00035 template <unsigned int Dim, FEFamily T> 00036 unsigned int FE<Dim,T>::n_shape_functions () const 00037 { 00038 return FE<Dim,T>::n_dofs (this->elem_type, 00039 static_cast<Order>(this->fe_type.order + this->_p_level)); 00040 } 00041 00042 00043 template <unsigned int Dim, FEFamily T> 00044 void FE<Dim,T>::attach_quadrature_rule (QBase* q) 00045 { 00046 libmesh_assert(q); 00047 this->qrule = q; 00048 // make sure we don't cache results from a previous quadrature rule 00049 this->elem_type = INVALID_ELEM; 00050 return; 00051 } 00052 00053 00054 template <unsigned int Dim, FEFamily T> 00055 unsigned int FE<Dim,T>::n_quadrature_points () const 00056 { 00057 libmesh_assert(this->qrule); 00058 return this->qrule->n_points(); 00059 } 00060 00061 00062 template <unsigned int Dim, FEFamily T> 00063 void FE<Dim,T>::dofs_on_side(const Elem* const elem, 00064 const Order o, 00065 unsigned int s, 00066 std::vector<unsigned int>& di) 00067 { 00068 libmesh_assert(elem); 00069 libmesh_assert_less (s, elem->n_sides()); 00070 00071 di.clear(); 00072 unsigned int nodenum = 0; 00073 const unsigned int n_nodes = elem->n_nodes(); 00074 for (unsigned int n = 0; n != n_nodes; ++n) 00075 { 00076 const unsigned int n_dofs = n_dofs_at_node(elem->type(), 00077 static_cast<Order>(o + elem->p_level()), n); 00078 if (elem->is_node_on_side(n, s)) 00079 for (unsigned int i = 0; i != n_dofs; ++i) 00080 di.push_back(nodenum++); 00081 else 00082 nodenum += n_dofs; 00083 } 00084 } 00085 00086 00087 00088 template <unsigned int Dim, FEFamily T> 00089 void FE<Dim,T>::dofs_on_edge(const Elem* const elem, 00090 const Order o, 00091 unsigned int e, 00092 std::vector<unsigned int>& di) 00093 { 00094 libmesh_assert(elem); 00095 libmesh_assert_less (e, elem->n_edges()); 00096 00097 di.clear(); 00098 unsigned int nodenum = 0; 00099 const unsigned int n_nodes = elem->n_nodes(); 00100 for (unsigned int n = 0; n != n_nodes; ++n) 00101 { 00102 const unsigned int n_dofs = n_dofs_at_node(elem->type(), 00103 static_cast<Order>(o + elem->p_level()), n); 00104 if (elem->is_node_on_edge(n, e)) 00105 for (unsigned int i = 0; i != n_dofs; ++i) 00106 di.push_back(nodenum++); 00107 else 00108 nodenum += n_dofs; 00109 } 00110 } 00111 00112 00113 00114 template <unsigned int Dim, FEFamily T> 00115 void FE<Dim,T>::reinit(const Elem* elem, 00116 const std::vector<Point>* const pts, 00117 const std::vector<Real>* const weights) 00118 { 00119 libmesh_assert(elem); 00120 00121 // Try to avoid calling init_shape_functions 00122 // even when shapes_need_reinit 00123 bool cached_nodes_still_fit = false; 00124 00125 // Initialize the shape functions at the user-specified 00126 // points 00127 if (pts != NULL) 00128 { 00129 // Set the type and p level for this element 00130 this->elem_type = elem->type(); 00131 this->_p_level = elem->p_level(); 00132 00133 // Initialize the shape functions 00134 this->_fe_map->template init_reference_to_physical_map<Dim>(*pts, elem); 00135 this->init_shape_functions (*pts, elem); 00136 00137 // The shape functions do not correspond to the qrule 00138 this->shapes_on_quadrature = false; 00139 } 00140 00141 // If there are no user specified points, we use the 00142 // quadrature rule 00143 00144 // update the type in accordance to the current cell 00145 // and reinit if the cell type has changed or (as in 00146 // the case of the hierarchics) the shape functions need 00147 // reinit, since they depend on the particular element shape 00148 else 00149 { 00150 libmesh_assert(this->qrule); 00151 this->qrule->init(elem->type(), elem->p_level()); 00152 00153 if(this->qrule->shapes_need_reinit()) 00154 this->shapes_on_quadrature = false; 00155 00156 if (this->elem_type != elem->type() || 00157 this->_p_level != elem->p_level() || 00158 !this->shapes_on_quadrature) 00159 { 00160 // Set the type and p level for this element 00161 this->elem_type = elem->type(); 00162 this->_p_level = elem->p_level(); 00163 // Initialize the shape functions 00164 this->_fe_map->template init_reference_to_physical_map<Dim>(this->qrule->get_points(), elem); 00165 this->init_shape_functions (this->qrule->get_points(), elem); 00166 00167 if (this->shapes_need_reinit()) 00168 { 00169 cached_nodes.resize(elem->n_nodes()); 00170 for (unsigned int n = 0; n != elem->n_nodes(); ++n) 00171 { 00172 cached_nodes[n] = elem->point(n); 00173 } 00174 } 00175 } 00176 else 00177 { 00178 // libmesh_assert_greater (elem->n_nodes(), 1); 00179 00180 cached_nodes_still_fit = true; 00181 if (cached_nodes.size() != elem->n_nodes()) 00182 cached_nodes_still_fit = false; 00183 else 00184 for (unsigned int n = 1; n < elem->n_nodes(); ++n) 00185 { 00186 if (!(elem->point(n) - elem->point(0)).relative_fuzzy_equals( 00187 (cached_nodes[n] - cached_nodes[0]), 1e-13)) 00188 { 00189 cached_nodes_still_fit = false; 00190 break; 00191 } 00192 } 00193 00194 if (this->shapes_need_reinit() && !cached_nodes_still_fit) 00195 { 00196 this->_fe_map->template init_reference_to_physical_map<Dim>(this->qrule->get_points(), elem); 00197 this->init_shape_functions (this->qrule->get_points(), elem); 00198 cached_nodes.resize(elem->n_nodes()); 00199 for (unsigned int n = 0; n != elem->n_nodes(); ++n) 00200 cached_nodes[n] = elem->point(n); 00201 } 00202 } 00203 00204 // The shape functions correspond to the qrule 00205 this->shapes_on_quadrature = true; 00206 } 00207 00208 // Compute the map for this element. In the future we can specify 00209 // different types of maps 00210 if (pts != NULL) 00211 { 00212 if (weights != NULL) 00213 { 00214 this->_fe_map->compute_map (this->dim,*weights, elem); 00215 } 00216 else 00217 { 00218 std::vector<Real> dummy_weights (pts->size(), 1.); 00219 this->_fe_map->compute_map (this->dim,dummy_weights, elem); 00220 } 00221 } 00222 else 00223 { 00224 this->_fe_map->compute_map (this->dim,this->qrule->get_weights(), elem); 00225 } 00226 00227 // Compute the shape functions and the derivatives at all of the 00228 // quadrature points. 00229 if (!cached_nodes_still_fit) 00230 { 00231 if (pts != NULL) 00232 this->compute_shape_functions (elem,*pts); 00233 else 00234 this->compute_shape_functions(elem,this->qrule->get_points()); 00235 } 00236 } 00237 00238 00239 00240 template <unsigned int Dim, FEFamily T> 00241 void FE<Dim,T>::init_shape_functions(const std::vector<Point>& qp, 00242 const Elem* elem) 00243 { 00244 libmesh_assert(elem); 00245 this->calculations_started = true; 00246 00247 // If the user forgot to request anything, we'll be safe and 00248 // calculate everything: 00249 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00250 if (!this->calculate_phi && !this->calculate_dphi && !this->calculate_d2phi 00251 && !this->calculate_curl_phi && !this->calculate_div_phi) 00252 { 00253 this->calculate_phi = this->calculate_dphi = this->calculate_d2phi = this->calculate_dphiref = true; 00254 if( FEInterface::field_type(T) == TYPE_VECTOR ) 00255 { 00256 this->calculate_curl_phi = true; 00257 this->calculate_div_phi = true; 00258 } 00259 } 00260 #else 00261 if (!this->calculate_phi && !this->calculate_dphi && !this->calculate_curl_phi && !this->calculate_div_phi) 00262 { 00263 this->calculate_phi = this->calculate_dphi = this->calculate_dphiref = true; 00264 if( FEInterface::field_type(T) == TYPE_VECTOR ) 00265 { 00266 this->calculate_curl_phi = true; 00267 this->calculate_div_phi = true; 00268 } 00269 } 00270 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 00271 00272 // Start logging the shape function initialization 00273 START_LOG("init_shape_functions()", "FE"); 00274 00275 00276 // The number of quadrature points. 00277 const unsigned int n_qp = libmesh_cast_int<unsigned int>(qp.size()); 00278 00279 // Number of shape functions in the finite element approximation 00280 // space. 00281 const unsigned int n_approx_shape_functions = 00282 this->n_shape_functions(this->get_type(), 00283 this->get_order()); 00284 00285 // resize the vectors to hold current data 00286 // Phi are the shape functions used for the FE approximation 00287 // Phi_map are the shape functions used for the FE mapping 00288 if (this->calculate_phi) 00289 this->phi.resize (n_approx_shape_functions); 00290 00291 if (this->calculate_dphi) 00292 { 00293 this->dphi.resize (n_approx_shape_functions); 00294 this->dphidx.resize (n_approx_shape_functions); 00295 this->dphidy.resize (n_approx_shape_functions); 00296 this->dphidz.resize (n_approx_shape_functions); 00297 } 00298 00299 if(this->calculate_dphiref) 00300 { 00301 if (Dim > 0) 00302 this->dphidxi.resize (n_approx_shape_functions); 00303 00304 if (Dim > 1) 00305 this->dphideta.resize (n_approx_shape_functions); 00306 00307 if (Dim > 2) 00308 this->dphidzeta.resize (n_approx_shape_functions); 00309 } 00310 00311 if( this->calculate_curl_phi && (FEInterface::field_type(T) == TYPE_VECTOR) ) 00312 this->curl_phi.resize(n_approx_shape_functions); 00313 00314 if( this->calculate_div_phi && (FEInterface::field_type(T) == TYPE_VECTOR) ) 00315 this->div_phi.resize(n_approx_shape_functions); 00316 00317 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00318 if (this->calculate_d2phi) 00319 { 00320 this->d2phi.resize (n_approx_shape_functions); 00321 this->d2phidx2.resize (n_approx_shape_functions); 00322 this->d2phidxdy.resize (n_approx_shape_functions); 00323 this->d2phidxdz.resize (n_approx_shape_functions); 00324 this->d2phidy2.resize (n_approx_shape_functions); 00325 this->d2phidydz.resize (n_approx_shape_functions); 00326 this->d2phidz2.resize (n_approx_shape_functions); 00327 00328 if (Dim > 0) 00329 this->d2phidxi2.resize (n_approx_shape_functions); 00330 00331 if (Dim > 1) 00332 { 00333 this->d2phidxideta.resize (n_approx_shape_functions); 00334 this->d2phideta2.resize (n_approx_shape_functions); 00335 } 00336 if (Dim > 2) 00337 { 00338 this->d2phidxidzeta.resize (n_approx_shape_functions); 00339 this->d2phidetadzeta.resize (n_approx_shape_functions); 00340 this->d2phidzeta2.resize (n_approx_shape_functions); 00341 } 00342 } 00343 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00344 00345 for (unsigned int i=0; i<n_approx_shape_functions; i++) 00346 { 00347 if (this->calculate_phi) 00348 this->phi[i].resize (n_qp); 00349 if (this->calculate_dphi) 00350 { 00351 this->dphi[i].resize (n_qp); 00352 this->dphidx[i].resize (n_qp); 00353 this->dphidy[i].resize (n_qp); 00354 this->dphidz[i].resize (n_qp); 00355 } 00356 00357 if(this->calculate_dphiref) 00358 { 00359 if (Dim > 0) 00360 this->dphidxi[i].resize(n_qp); 00361 00362 if (Dim > 1) 00363 this->dphideta[i].resize(n_qp); 00364 00365 if (Dim > 2) 00366 this->dphidzeta[i].resize(n_qp); 00367 } 00368 00369 if(this->calculate_curl_phi && (FEInterface::field_type(T) == TYPE_VECTOR) ) 00370 this->curl_phi[i].resize(n_qp); 00371 00372 if(this->calculate_div_phi && (FEInterface::field_type(T) == TYPE_VECTOR) ) 00373 this->div_phi[i].resize(n_qp); 00374 00375 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00376 if (this->calculate_d2phi) 00377 { 00378 this->d2phi[i].resize (n_qp); 00379 this->d2phidx2[i].resize (n_qp); 00380 this->d2phidxdy[i].resize (n_qp); 00381 this->d2phidxdz[i].resize (n_qp); 00382 this->d2phidy2[i].resize (n_qp); 00383 this->d2phidydz[i].resize (n_qp); 00384 this->d2phidz2[i].resize (n_qp); 00385 if (Dim > 0) 00386 this->d2phidxi2[i].resize (n_qp); 00387 if (Dim > 1) 00388 { 00389 this->d2phidxideta[i].resize (n_qp); 00390 this->d2phideta2[i].resize (n_qp); 00391 } 00392 if (Dim > 2) 00393 { 00394 this->d2phidxidzeta[i].resize (n_qp); 00395 this->d2phidetadzeta[i].resize (n_qp); 00396 this->d2phidzeta2[i].resize (n_qp); 00397 } 00398 } 00399 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00400 } 00401 00402 00403 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00404 //------------------------------------------------------------ 00405 // Initialize the data fields, which should only be used for infinite 00406 // elements, to some sensible values, so that using a FE with the 00407 // variational formulation of an InfFE, correct element matrices are 00408 // returned 00409 00410 { 00411 this->weight.resize (n_qp); 00412 this->dweight.resize (n_qp); 00413 this->dphase.resize (n_qp); 00414 00415 for (unsigned int p=0; p<n_qp; p++) 00416 { 00417 this->weight[p] = 1.; 00418 this->dweight[p].zero(); 00419 this->dphase[p].zero(); 00420 } 00421 00422 } 00423 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00424 00425 switch (Dim) 00426 { 00427 00428 //------------------------------------------------------------ 00429 // 0D 00430 case 0: 00431 { 00432 break; 00433 } 00434 00435 //------------------------------------------------------------ 00436 // 1D 00437 case 1: 00438 { 00439 // Compute the value of the approximation shape function i at quadrature point p 00440 if (this->calculate_dphiref) 00441 for (unsigned int i=0; i<n_approx_shape_functions; i++) 00442 for (unsigned int p=0; p<n_qp; p++) 00443 this->dphidxi[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 0, qp[p]); 00444 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00445 if (this->calculate_d2phi) 00446 for (unsigned int i=0; i<n_approx_shape_functions; i++) 00447 for (unsigned int p=0; p<n_qp; p++) 00448 this->d2phidxi2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 0, qp[p]); 00449 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00450 00451 break; 00452 } 00453 00454 00455 00456 //------------------------------------------------------------ 00457 // 2D 00458 case 2: 00459 { 00460 // Compute the value of the approximation shape function i at quadrature point p 00461 if (this->calculate_dphiref) 00462 for (unsigned int i=0; i<n_approx_shape_functions; i++) 00463 for (unsigned int p=0; p<n_qp; p++) 00464 { 00465 this->dphidxi[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 0, qp[p]); 00466 this->dphideta[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 1, qp[p]); 00467 } 00468 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00469 if (this->calculate_d2phi) 00470 for (unsigned int i=0; i<n_approx_shape_functions; i++) 00471 for (unsigned int p=0; p<n_qp; p++) 00472 { 00473 this->d2phidxi2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 0, qp[p]); 00474 this->d2phidxideta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 1, qp[p]); 00475 this->d2phideta2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 2, qp[p]); 00476 } 00477 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00478 00479 00480 break; 00481 } 00482 00483 00484 00485 //------------------------------------------------------------ 00486 // 3D 00487 case 3: 00488 { 00489 // Compute the value of the approximation shape function i at quadrature point p 00490 if (this->calculate_dphiref) 00491 for (unsigned int i=0; i<n_approx_shape_functions; i++) 00492 for (unsigned int p=0; p<n_qp; p++) 00493 { 00494 this->dphidxi[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 0, qp[p]); 00495 this->dphideta[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 1, qp[p]); 00496 this->dphidzeta[i][p] = FE<Dim,T>::shape_deriv (elem, this->fe_type.order, i, 2, qp[p]); 00497 } 00498 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00499 if (this->calculate_d2phi) 00500 for (unsigned int i=0; i<n_approx_shape_functions; i++) 00501 for (unsigned int p=0; p<n_qp; p++) 00502 { 00503 this->d2phidxi2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 0, qp[p]); 00504 this->d2phidxideta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 1, qp[p]); 00505 this->d2phideta2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 2, qp[p]); 00506 this->d2phidxidzeta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 3, qp[p]); 00507 this->d2phidetadzeta[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 4, qp[p]); 00508 this->d2phidzeta2[i][p] = FE<Dim,T>::shape_second_deriv (elem, this->fe_type.order, i, 5, qp[p]); 00509 } 00510 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00511 00512 break; 00513 } 00514 00515 00516 default: 00517 libmesh_error(); 00518 } 00519 00520 // Stop logging the shape function initialization 00521 STOP_LOG("init_shape_functions()", "FE"); 00522 } 00523 00524 00525 00526 00527 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00528 00529 template <unsigned int Dim, FEFamily T> 00530 void FE<Dim,T>::init_base_shape_functions(const std::vector<Point>& qp, 00531 const Elem* e) 00532 { 00533 // I don't understand infinite elements well enough to risk 00534 // calculating too little. :-( RHS 00535 this->calculate_phi = this->calculate_dphi = this->calculate_d2phi = true; 00536 00537 this->elem_type = e->type(); 00538 this->_fe_map->template init_reference_to_physical_map<Dim>(qp, e); 00539 init_shape_functions(qp, e); 00540 } 00541 00542 #endif // LIBMESH_ENABLE_INFINITE_ELEMENTS 00543 00544 00545 00546 //-------------------------------------------------------------- 00547 // Explicit instantiations using macro from fe_macro.h 00548 00549 INSTANTIATE_FE(0); 00550 00551 INSTANTIATE_FE(1); 00552 00553 INSTANTIATE_FE(2); 00554 00555 INSTANTIATE_FE(3); 00556 00557 00558 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: