fe_xyz.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_logging.h" 00022 #include "libmesh/fe.h" 00023 #include "libmesh/elem.h" 00024 #include "libmesh/fe_interface.h" 00025 00026 namespace libMesh 00027 { 00028 00029 // ------------------------------------------------------------ 00030 // XYZ-specific implementations 00031 00032 // Anonymous namespace for local helper functions 00033 namespace { 00034 00035 void xyz_nodal_soln(const Elem* elem, 00036 const Order order, 00037 const std::vector<Number>& elem_soln, 00038 std::vector<Number>& nodal_soln, 00039 unsigned Dim) 00040 { 00041 const unsigned int n_nodes = elem->n_nodes(); 00042 00043 const ElemType elem_type = elem->type(); 00044 00045 nodal_soln.resize(n_nodes); 00046 00047 const Order totalorder = static_cast<Order>(order + elem->p_level()); 00048 00049 switch (totalorder) 00050 { 00051 // Constant shape functions 00052 case CONSTANT: 00053 { 00054 libmesh_assert_equal_to (elem_soln.size(), 1); 00055 00056 const Number val = elem_soln[0]; 00057 00058 for (unsigned int n=0; n<n_nodes; n++) 00059 nodal_soln[n] = val; 00060 00061 return; 00062 } 00063 00064 00065 // For other orders do interpolation at the nodes 00066 // explicitly. 00067 default: 00068 { 00069 // FEType object to be passed to various FEInterface functions below. 00070 FEType fe_type(totalorder, XYZ); 00071 00072 const unsigned int n_sf = 00073 // FE<Dim,T>::n_shape_functions(elem_type, totalorder); 00074 FEInterface::n_shape_functions(Dim, fe_type, elem_type); 00075 00076 for (unsigned int n=0; n<n_nodes; n++) 00077 { 00078 libmesh_assert_equal_to (elem_soln.size(), n_sf); 00079 00080 // Zero before summation 00081 nodal_soln[n] = 0; 00082 00083 // u_i = Sum (alpha_i phi_i) 00084 for (unsigned int i=0; i<n_sf; i++) 00085 nodal_soln[n] += elem_soln[i] * 00086 // FE<Dim,T>::shape(elem, order, i, elem->point(n)); 00087 FEInterface::shape(Dim, fe_type, elem, i, elem->point(n)); 00088 } 00089 00090 return; 00091 } // default 00092 } // switch 00093 } // xyz_nodal_soln() 00094 00095 00096 00097 00098 00099 unsigned int xyz_n_dofs(const ElemType t, const Order o) 00100 { 00101 switch (o) 00102 { 00103 00104 // constant shape functions 00105 // no matter what shape there is only one DOF. 00106 case CONSTANT: 00107 return 1; 00108 00109 00110 // Discontinuous linear shape functions 00111 // expressed in the XYZ monomials. 00112 case FIRST: 00113 { 00114 switch (t) 00115 { 00116 case NODEELEM: 00117 return 1; 00118 00119 case EDGE2: 00120 case EDGE3: 00121 case EDGE4: 00122 return 2; 00123 00124 case TRI3: 00125 case TRI6: 00126 case QUAD4: 00127 case QUAD8: 00128 case QUAD9: 00129 return 3; 00130 00131 case TET4: 00132 case TET10: 00133 case HEX8: 00134 case HEX20: 00135 case HEX27: 00136 case PRISM6: 00137 case PRISM15: 00138 case PRISM18: 00139 case PYRAMID5: 00140 return 4; 00141 00142 default: 00143 { 00144 #ifdef DEBUG 00145 libMesh::err << "ERROR: Bad ElemType = " << t 00146 << " for " << o << "th order approximation!" 00147 << std::endl; 00148 #endif 00149 libmesh_error(); 00150 } 00151 } 00152 } 00153 00154 00155 // Discontinuous quadratic shape functions 00156 // expressed in the XYZ monomials. 00157 case SECOND: 00158 { 00159 switch (t) 00160 { 00161 case NODEELEM: 00162 return 1; 00163 00164 case EDGE2: 00165 case EDGE3: 00166 case EDGE4: 00167 return 3; 00168 00169 case TRI3: 00170 case TRI6: 00171 case QUAD4: 00172 case QUAD8: 00173 case QUAD9: 00174 return 6; 00175 00176 case TET4: 00177 case TET10: 00178 case HEX8: 00179 case HEX20: 00180 case HEX27: 00181 case PRISM6: 00182 case PRISM15: 00183 case PRISM18: 00184 case PYRAMID5: 00185 return 10; 00186 00187 default: 00188 { 00189 #ifdef DEBUG 00190 libMesh::err << "ERROR: Bad ElemType = " << t 00191 << " for " << o << "th order approximation!" 00192 << std::endl; 00193 #endif 00194 libmesh_error(); 00195 } 00196 } 00197 } 00198 00199 00200 // Discontinuous cubic shape functions 00201 // expressed in the XYZ monomials. 00202 case THIRD: 00203 { 00204 switch (t) 00205 { 00206 case NODEELEM: 00207 return 1; 00208 00209 case EDGE2: 00210 case EDGE3: 00211 case EDGE4: 00212 return 4; 00213 00214 case TRI3: 00215 case TRI6: 00216 case QUAD4: 00217 case QUAD8: 00218 case QUAD9: 00219 return 10; 00220 00221 case TET4: 00222 case TET10: 00223 case HEX8: 00224 case HEX20: 00225 case HEX27: 00226 case PRISM6: 00227 case PRISM15: 00228 case PRISM18: 00229 case PYRAMID5: 00230 return 20; 00231 00232 default: 00233 { 00234 #ifdef DEBUG 00235 libMesh::err << "ERROR: Bad ElemType = " << t 00236 << " for " << o << "th order approximation!" 00237 << std::endl; 00238 #endif 00239 libmesh_error(); 00240 } 00241 } 00242 } 00243 00244 00245 // Discontinuous quartic shape functions 00246 // expressed in the XYZ monomials. 00247 case FOURTH: 00248 { 00249 switch (t) 00250 { 00251 case NODEELEM: 00252 return 1; 00253 00254 case EDGE2: 00255 case EDGE3: 00256 return 5; 00257 00258 case TRI3: 00259 case TRI6: 00260 case QUAD4: 00261 case QUAD8: 00262 case QUAD9: 00263 return 15; 00264 00265 case TET4: 00266 case TET10: 00267 case HEX8: 00268 case HEX20: 00269 case HEX27: 00270 case PRISM6: 00271 case PRISM15: 00272 case PRISM18: 00273 case PYRAMID5: 00274 return 35; 00275 00276 default: 00277 { 00278 #ifdef DEBUG 00279 libMesh::err << "ERROR: Bad ElemType = " << t 00280 << " for " << o << "th order approximation!" 00281 << std::endl; 00282 #endif 00283 libmesh_error(); 00284 } 00285 } 00286 } 00287 00288 00289 default: 00290 { 00291 const unsigned int order = static_cast<unsigned int>(o); 00292 switch (t) 00293 { 00294 case NODEELEM: 00295 return 1; 00296 00297 case EDGE2: 00298 case EDGE3: 00299 return (order+1); 00300 00301 case TRI3: 00302 case TRI6: 00303 case QUAD4: 00304 case QUAD8: 00305 case QUAD9: 00306 return (order+1)*(order+2)/2; 00307 00308 case TET4: 00309 case TET10: 00310 case HEX8: 00311 case HEX20: 00312 case HEX27: 00313 case PRISM6: 00314 case PRISM15: 00315 case PRISM18: 00316 case PYRAMID5: 00317 return (order+1)*(order+2)*(order+3)/6; 00318 00319 default: 00320 { 00321 #ifdef DEBUG 00322 libMesh::err << "ERROR: Bad ElemType = " << t 00323 << " for " << o << "th order approximation!" 00324 << std::endl; 00325 #endif 00326 libmesh_error(); 00327 } 00328 } 00329 } 00330 } 00331 00332 libmesh_error(); 00333 00334 return 0; 00335 } 00336 00337 00338 00339 00340 unsigned int xyz_n_dofs_per_elem(const ElemType t, 00341 const Order o) 00342 { 00343 switch (o) 00344 { 00345 // constant shape functions always have 1 DOF per element 00346 case CONSTANT: 00347 return 1; 00348 00349 00350 // Discontinuous linear shape functions 00351 // expressed in the XYZ monomials. 00352 case FIRST: 00353 { 00354 switch (t) 00355 { 00356 case NODEELEM: 00357 return 1; 00358 00359 // 1D linears have 2 DOFs per element 00360 case EDGE2: 00361 case EDGE3: 00362 case EDGE4: 00363 return 2; 00364 00365 // 2D linears have 3 DOFs per element 00366 case TRI3: 00367 case TRI6: 00368 case QUAD4: 00369 case QUAD8: 00370 case QUAD9: 00371 return 3; 00372 00373 // 3D linears have 4 DOFs per element 00374 case TET4: 00375 case TET10: 00376 case HEX8: 00377 case HEX20: 00378 case HEX27: 00379 case PRISM6: 00380 case PRISM15: 00381 case PRISM18: 00382 case PYRAMID5: 00383 return 4; 00384 00385 default: 00386 { 00387 #ifdef DEBUG 00388 libMesh::err << "ERROR: Bad ElemType = " << t 00389 << " for " << o << "th order approximation!" 00390 << std::endl; 00391 #endif 00392 libmesh_error(); 00393 } 00394 } 00395 } 00396 00397 00398 // Discontinuous quadratic shape functions 00399 // expressed in the XYZ monomials. 00400 case SECOND: 00401 { 00402 switch (t) 00403 { 00404 case NODEELEM: 00405 return 1; 00406 00407 // 1D quadratics have 3 DOFs per element 00408 case EDGE2: 00409 case EDGE3: 00410 case EDGE4: 00411 return 3; 00412 00413 // 2D quadratics have 6 DOFs per element 00414 case TRI3: 00415 case TRI6: 00416 case QUAD4: 00417 case QUAD8: 00418 case QUAD9: 00419 return 6; 00420 00421 // 3D quadratics have 10 DOFs per element 00422 case TET4: 00423 case TET10: 00424 case HEX8: 00425 case HEX20: 00426 case HEX27: 00427 case PRISM6: 00428 case PRISM15: 00429 case PRISM18: 00430 case PYRAMID5: 00431 return 10; 00432 00433 default: 00434 { 00435 #ifdef DEBUG 00436 libMesh::err << "ERROR: Bad ElemType = " << t 00437 << " for " << o << "th order approximation!" 00438 << std::endl; 00439 #endif 00440 libmesh_error(); 00441 } 00442 } 00443 } 00444 00445 00446 // Discontinuous cubic shape functions 00447 // expressed in the XYZ monomials. 00448 case THIRD: 00449 { 00450 switch (t) 00451 { 00452 case NODEELEM: 00453 return 1; 00454 00455 case EDGE2: 00456 case EDGE3: 00457 case EDGE4: 00458 return 4; 00459 00460 case TRI3: 00461 case TRI6: 00462 case QUAD4: 00463 case QUAD8: 00464 case QUAD9: 00465 return 10; 00466 00467 case TET4: 00468 case TET10: 00469 case HEX8: 00470 case HEX20: 00471 case HEX27: 00472 case PRISM6: 00473 case PRISM15: 00474 case PRISM18: 00475 case PYRAMID5: 00476 return 20; 00477 00478 default: 00479 { 00480 #ifdef DEBUG 00481 libMesh::err << "ERROR: Bad ElemType = " << t 00482 << " for " << o << "th order approximation!" 00483 << std::endl; 00484 #endif 00485 libmesh_error(); 00486 } 00487 } 00488 } 00489 00490 00491 // Discontinuous quartic shape functions 00492 // expressed in the XYZ monomials. 00493 case FOURTH: 00494 { 00495 switch (t) 00496 { 00497 case NODEELEM: 00498 return 1; 00499 00500 case EDGE2: 00501 case EDGE3: 00502 case EDGE4: 00503 return 5; 00504 00505 case TRI3: 00506 case TRI6: 00507 case QUAD4: 00508 case QUAD8: 00509 case QUAD9: 00510 return 15; 00511 00512 case TET4: 00513 case TET10: 00514 case HEX8: 00515 case HEX20: 00516 case HEX27: 00517 case PRISM6: 00518 case PRISM15: 00519 case PRISM18: 00520 case PYRAMID5: 00521 return 35; 00522 00523 default: 00524 { 00525 #ifdef DEBUG 00526 libMesh::err << "ERROR: Bad ElemType = " << t 00527 << " for " << o << "th order approximation!" 00528 << std::endl; 00529 #endif 00530 libmesh_error(); 00531 } 00532 } 00533 } 00534 00535 default: 00536 { 00537 const unsigned int order = static_cast<unsigned int>(o); 00538 switch (t) 00539 { 00540 case NODEELEM: 00541 return 1; 00542 00543 case EDGE2: 00544 case EDGE3: 00545 return (order+1); 00546 00547 case TRI3: 00548 case TRI6: 00549 case QUAD4: 00550 case QUAD8: 00551 case QUAD9: 00552 return (order+1)*(order+2)/2; 00553 00554 case TET4: 00555 case TET10: 00556 case HEX8: 00557 case HEX20: 00558 case HEX27: 00559 case PRISM6: 00560 case PRISM15: 00561 case PRISM18: 00562 case PYRAMID5: 00563 return (order+1)*(order+2)*(order+3)/6; 00564 00565 default: 00566 { 00567 #ifdef DEBUG 00568 libMesh::err << "ERROR: Bad ElemType = " << t 00569 << " for " << o << "th order approximation!" 00570 << std::endl; 00571 #endif 00572 libmesh_error(); 00573 } 00574 } 00575 } 00576 return 0; 00577 } 00578 } 00579 00580 00581 } // anonymous namespace 00582 00583 00584 00585 00586 00587 00588 00589 template <unsigned int Dim> 00590 void FEXYZ<Dim>::init_shape_functions(const std::vector<Point>& qp, 00591 const Elem* elem) 00592 { 00593 libmesh_assert(elem); 00594 this->calculations_started = true; 00595 00596 // If the user forgot to request anything, we'll be safe and 00597 // calculate everything: 00598 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00599 if (!this->calculate_phi && !this->calculate_dphi && !this->calculate_d2phi) 00600 this->calculate_phi = this->calculate_dphi = this->calculate_d2phi = true; 00601 #else 00602 if (!this->calculate_phi && !this->calculate_dphi) 00603 this->calculate_phi = this->calculate_dphi = true; 00604 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 00605 00606 // Start logging the shape function initialization 00607 START_LOG("init_shape_functions()", "FE"); 00608 00609 00610 // The number of quadrature points. 00611 const std::size_t n_qp = qp.size(); 00612 00613 // Number of shape functions in the finite element approximation 00614 // space. 00615 const unsigned int n_approx_shape_functions = 00616 this->n_shape_functions(this->get_type(), 00617 this->get_order()); 00618 00619 // resize the vectors to hold current data 00620 // Phi are the shape functions used for the FE approximation 00621 { 00622 // (note: GCC 3.4.0 requires the use of this-> here) 00623 if (this->calculate_phi) 00624 this->phi.resize (n_approx_shape_functions); 00625 if (this->calculate_dphi) 00626 { 00627 this->dphi.resize (n_approx_shape_functions); 00628 this->dphidx.resize (n_approx_shape_functions); 00629 this->dphidy.resize (n_approx_shape_functions); 00630 this->dphidz.resize (n_approx_shape_functions); 00631 } 00632 00633 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00634 if (this->calculate_d2phi) 00635 { 00636 this->d2phi.resize (n_approx_shape_functions); 00637 this->d2phidx2.resize (n_approx_shape_functions); 00638 this->d2phidxdy.resize (n_approx_shape_functions); 00639 this->d2phidxdz.resize (n_approx_shape_functions); 00640 this->d2phidy2.resize (n_approx_shape_functions); 00641 this->d2phidydz.resize (n_approx_shape_functions); 00642 this->d2phidz2.resize (n_approx_shape_functions); 00643 this->d2phidxi2.resize (n_approx_shape_functions); 00644 } 00645 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00646 00647 for (unsigned int i=0; i<n_approx_shape_functions; i++) 00648 { 00649 if (this->calculate_phi) 00650 this->phi[i].resize (n_qp); 00651 if (this->calculate_dphi) 00652 { 00653 this->dphi[i].resize (n_qp); 00654 this->dphidx[i].resize (n_qp); 00655 this->dphidy[i].resize (n_qp); 00656 this->dphidz[i].resize (n_qp); 00657 } 00658 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00659 if (this->calculate_d2phi) 00660 { 00661 this->d2phi[i].resize (n_qp); 00662 this->d2phidx2[i].resize (n_qp); 00663 this->d2phidxdy[i].resize (n_qp); 00664 this->d2phidxdz[i].resize (n_qp); 00665 this->d2phidy2[i].resize (n_qp); 00666 this->d2phidydz[i].resize (n_qp); 00667 this->d2phidz2[i].resize (n_qp); 00668 } 00669 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00670 } 00671 } 00672 00673 00674 00675 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00676 //------------------------------------------------------------ 00677 // Initialize the data fields, which should only be used for infinite 00678 // elements, to some sensible values, so that using a FE with the 00679 // variational formulation of an InfFE, correct element matrices are 00680 // returned 00681 00682 { 00683 this->weight.resize (n_qp); 00684 this->dweight.resize (n_qp); 00685 this->dphase.resize (n_qp); 00686 00687 for (unsigned int p=0; p<n_qp; p++) 00688 { 00689 this->weight[p] = 1.; 00690 this->dweight[p].zero(); 00691 this->dphase[p].zero(); 00692 } 00693 00694 } 00695 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00696 00697 // Stop logging the shape function initialization 00698 STOP_LOG("init_shape_functions()", "FE"); 00699 } 00700 00701 00702 00703 00704 template <unsigned int Dim> 00705 void FEXYZ<Dim>::compute_shape_functions (const Elem* elem, const std::vector<Point>&) 00706 { 00707 libmesh_assert(elem); 00708 00709 //------------------------------------------------------------------------- 00710 // Compute the shape function values (and derivatives) 00711 // at the Quadrature points. Note that the actual values 00712 // have already been computed via init_shape_functions 00713 00714 // Start logging the shape function computation 00715 START_LOG("compute_shape_functions()", "FE"); 00716 00717 const std::vector<Point>& xyz_qp = this->get_xyz(); 00718 00719 // Compute the value of the derivative shape function i at quadrature point p 00720 switch (this->dim) 00721 { 00722 00723 case 1: 00724 { 00725 if (this->calculate_phi) 00726 for (unsigned int i=0; i<this->phi.size(); i++) 00727 for (unsigned int p=0; p<this->phi[i].size(); p++) 00728 this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz_qp[p]); 00729 if (this->calculate_dphi) 00730 for (unsigned int i=0; i<this->dphi.size(); i++) 00731 for (unsigned int p=0; p<this->dphi[i].size(); p++) 00732 { 00733 this->dphi[i][p](0) = 00734 this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]); 00735 00736 this->dphi[i][p](1) = this->dphidy[i][p] = 0.; 00737 this->dphi[i][p](2) = this->dphidz[i][p] = 0.; 00738 } 00739 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00740 if (this->calculate_d2phi) 00741 for (unsigned int i=0; i<this->d2phi.size(); i++) 00742 for (unsigned int p=0; p<this->d2phi[i].size(); p++) 00743 { 00744 this->d2phi[i][p](0,0) = 00745 this->d2phidx2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]); 00746 00747 #if LIBMESH_DIM>1 00748 this->d2phi[i][p](0,1) = this->d2phidxdy[i][p] = 00749 this->d2phi[i][p](1,0) = 0.; 00750 this->d2phi[i][p](1,1) = this->d2phidy2[i][p] = 0.; 00751 #if LIBMESH_DIM>2 00752 this->d2phi[i][p](0,2) = this->d2phidxdz[i][p] = 00753 this->d2phi[i][p](2,0) = 0.; 00754 this->d2phi[i][p](1,2) = this->d2phidydz[i][p] = 00755 this->d2phi[i][p](2,1) = 0.; 00756 this->d2phi[i][p](2,2) = this->d2phidz2[i][p] = 0.; 00757 #endif 00758 #endif 00759 } 00760 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00761 00762 // All done 00763 break; 00764 } 00765 00766 case 2: 00767 { 00768 if (this->calculate_phi) 00769 for (unsigned int i=0; i<this->phi.size(); i++) 00770 for (unsigned int p=0; p<this->phi[i].size(); p++) 00771 this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz_qp[p]); 00772 if (this->calculate_dphi) 00773 for (unsigned int i=0; i<this->dphi.size(); i++) 00774 for (unsigned int p=0; p<this->dphi[i].size(); p++) 00775 { 00776 this->dphi[i][p](0) = 00777 this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]); 00778 00779 this->dphi[i][p](1) = 00780 this->dphidy[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]); 00781 00782 #if LIBMESH_DIM == 3 00783 this->dphi[i][p](2) = // can only assign to the Z component if LIBMESH_DIM==3 00784 #endif 00785 this->dphidz[i][p] = 0.; 00786 } 00787 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00788 if (this->calculate_d2phi) 00789 for (unsigned int i=0; i<this->d2phi.size(); i++) 00790 for (unsigned int p=0; p<this->d2phi[i].size(); p++) 00791 { 00792 this->d2phi[i][p](0,0) = 00793 this->d2phidx2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]); 00794 00795 this->d2phi[i][p](0,1) = this->d2phidxdy[i][p] = 00796 this->d2phi[i][p](1,0) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]); 00797 this->d2phi[i][p](1,1) = 00798 this->d2phidy2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 2, xyz_qp[p]); 00799 #if LIBMESH_DIM>2 00800 this->d2phi[i][p](0,2) = this->d2phidxdz[i][p] = 00801 this->d2phi[i][p](2,0) = 0.; 00802 this->d2phi[i][p](1,2) = this->d2phidydz[i][p] = 00803 this->d2phi[i][p](2,1) = 0.; 00804 this->d2phi[i][p](2,2) = this->d2phidz2[i][p] = 0.; 00805 #endif 00806 } 00807 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00808 00809 // All done 00810 break; 00811 } 00812 00813 case 3: 00814 { 00815 if (this->calculate_phi) 00816 for (unsigned int i=0; i<this->phi.size(); i++) 00817 for (unsigned int p=0; p<this->phi[i].size(); p++) 00818 this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz_qp[p]); 00819 00820 if (this->calculate_dphi) 00821 for (unsigned int i=0; i<this->dphi.size(); i++) 00822 for (unsigned int p=0; p<this->dphi[i].size(); p++) 00823 { 00824 this->dphi[i][p](0) = 00825 this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]); 00826 00827 this->dphi[i][p](1) = 00828 this->dphidy[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]); 00829 00830 this->dphi[i][p](2) = 00831 this->dphidz[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 2, xyz_qp[p]); 00832 } 00833 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00834 if (this->calculate_d2phi) 00835 for (unsigned int i=0; i<this->d2phi.size(); i++) 00836 for (unsigned int p=0; p<this->d2phi[i].size(); p++) 00837 { 00838 this->d2phi[i][p](0,0) = 00839 this->d2phidx2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]); 00840 00841 this->d2phi[i][p](0,1) = this->d2phidxdy[i][p] = 00842 this->d2phi[i][p](1,0) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]); 00843 this->d2phi[i][p](1,1) = 00844 this->d2phidy2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 2, xyz_qp[p]); 00845 this->d2phi[i][p](0,2) = this->d2phidxdz[i][p] = 00846 this->d2phi[i][p](2,0) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 3, xyz_qp[p]); 00847 this->d2phi[i][p](1,2) = this->d2phidydz[i][p] = 00848 this->d2phi[i][p](2,1) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 4, xyz_qp[p]); 00849 this->d2phi[i][p](2,2) = this->d2phidz2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 5, xyz_qp[p]); 00850 } 00851 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00852 00853 // All done 00854 break; 00855 } 00856 00857 default: 00858 { 00859 libmesh_error(); 00860 } 00861 } 00862 00863 // Stop logging the shape function computation 00864 STOP_LOG("compute_shape_functions()", "FE"); 00865 } 00866 00867 00868 00869 00870 // Do full-specialization for every dimension, instead 00871 // of explicit instantiation at the end of this file. 00872 // This could be macro-ified so that it fits on one line... 00873 template <> 00874 void FE<0,XYZ>::nodal_soln(const Elem* elem, 00875 const Order order, 00876 const std::vector<Number>& elem_soln, 00877 std::vector<Number>& nodal_soln) 00878 { xyz_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/0); } 00879 00880 template <> 00881 void FE<1,XYZ>::nodal_soln(const Elem* elem, 00882 const Order order, 00883 const std::vector<Number>& elem_soln, 00884 std::vector<Number>& nodal_soln) 00885 { xyz_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/1); } 00886 00887 template <> 00888 void FE<2,XYZ>::nodal_soln(const Elem* elem, 00889 const Order order, 00890 const std::vector<Number>& elem_soln, 00891 std::vector<Number>& nodal_soln) 00892 { xyz_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/2); } 00893 00894 template <> 00895 void FE<3,XYZ>::nodal_soln(const Elem* elem, 00896 const Order order, 00897 const std::vector<Number>& elem_soln, 00898 std::vector<Number>& nodal_soln) 00899 { xyz_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/3); } 00900 00901 00902 00903 // Full specialization of n_dofs() function for every dimension 00904 template <> unsigned int FE<0,XYZ>::n_dofs(const ElemType t, const Order o) { return xyz_n_dofs(t, o); } 00905 template <> unsigned int FE<1,XYZ>::n_dofs(const ElemType t, const Order o) { return xyz_n_dofs(t, o); } 00906 template <> unsigned int FE<2,XYZ>::n_dofs(const ElemType t, const Order o) { return xyz_n_dofs(t, o); } 00907 template <> unsigned int FE<3,XYZ>::n_dofs(const ElemType t, const Order o) { return xyz_n_dofs(t, o); } 00908 00909 // Full specialization of n_dofs_at_node() function for every dimension. 00910 // XYZ FEMs have no dofs at nodes, only element dofs. 00911 template <> unsigned int FE<0,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; } 00912 template <> unsigned int FE<1,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; } 00913 template <> unsigned int FE<2,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; } 00914 template <> unsigned int FE<3,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; } 00915 00916 // Full specialization of n_dofs_per_elem() function for every dimension. 00917 template <> unsigned int FE<0,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return xyz_n_dofs_per_elem(t, o); } 00918 template <> unsigned int FE<1,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return xyz_n_dofs_per_elem(t, o); } 00919 template <> unsigned int FE<2,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return xyz_n_dofs_per_elem(t, o); } 00920 template <> unsigned int FE<3,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return xyz_n_dofs_per_elem(t, o); } 00921 00922 // Full specialization of get_continuity() function for every dimension. 00923 template <> FEContinuity FE<0,XYZ>::get_continuity() const { return DISCONTINUOUS; } 00924 template <> FEContinuity FE<1,XYZ>::get_continuity() const { return DISCONTINUOUS; } 00925 template <> FEContinuity FE<2,XYZ>::get_continuity() const { return DISCONTINUOUS; } 00926 template <> FEContinuity FE<3,XYZ>::get_continuity() const { return DISCONTINUOUS; } 00927 00928 // Full specialization of is_hierarchic() function for every dimension. 00929 // The XYZ shape functions are hierarchic! 00930 template <> bool FE<0,XYZ>::is_hierarchic() const { return true; } 00931 template <> bool FE<1,XYZ>::is_hierarchic() const { return true; } 00932 template <> bool FE<2,XYZ>::is_hierarchic() const { return true; } 00933 template <> bool FE<3,XYZ>::is_hierarchic() const { return true; } 00934 00935 #ifdef LIBMESH_ENABLE_AMR 00936 00937 // Full specialization of compute_constraints() function for 2D and 00938 // 3D only. There are no constraints for discontinuous elements, so 00939 // we do nothing. 00940 template <> void FE<2,XYZ>::compute_constraints (DofConstraints &, DofMap &, const unsigned int, const Elem*) {} 00941 template <> void FE<3,XYZ>::compute_constraints (DofConstraints &, DofMap &, const unsigned int, const Elem*) {} 00942 00943 #endif // #ifdef LIBMESH_ENABLE_AMR 00944 00945 // Full specialization of shapes_need_reinit() function for every dimension. 00946 template <> bool FE<0,XYZ>::shapes_need_reinit() const { return true; } 00947 template <> bool FE<1,XYZ>::shapes_need_reinit() const { return true; } 00948 template <> bool FE<2,XYZ>::shapes_need_reinit() const { return true; } 00949 template <> bool FE<3,XYZ>::shapes_need_reinit() const { return true; } 00950 00951 00952 // Explicit instantiations for non-static FEXYZ member functions. 00953 // These non-static member functions map more naturally to explicit 00954 // instantiations than the functions above: 00955 // 00956 // 1.) Since they are member functions, they rely on 00957 // private/protected member data, and therefore do not work well 00958 // with the "anonymous function call" model we've used above for 00959 // the specializations. 00960 // 00961 // 2.) There is (IMHO) less chance of the linker calling the 00962 // wrong version of one of these member functions, since there is 00963 // only one FEXYZ. 00964 template void FEXYZ<0>::init_shape_functions(const std::vector<Point>&, const Elem*); 00965 template void FEXYZ<1>::init_shape_functions(const std::vector<Point>&, const Elem*); 00966 template void FEXYZ<2>::init_shape_functions(const std::vector<Point>&, const Elem*); 00967 template void FEXYZ<3>::init_shape_functions(const std::vector<Point>&, const Elem*); 00968 00969 template void FEXYZ<0>::compute_shape_functions(const Elem*,const std::vector<Point>&); 00970 template void FEXYZ<1>::compute_shape_functions(const Elem*,const std::vector<Point>&); 00971 template void FEXYZ<2>::compute_shape_functions(const Elem*,const std::vector<Point>&); 00972 template void FEXYZ<3>::compute_shape_functions(const Elem*,const std::vector<Point>&); 00973 00974 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:46 UTC
Hosted By: