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