$extrastylesheet
fe_monomial.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/fe.h"
00022 #include "libmesh/elem.h"
00023 #include "libmesh/fe_interface.h"
00024 #include "libmesh/string_to_enum.h"
00025 
00026 namespace libMesh
00027 {
00028 
00029 // ------------------------------------------------------------
00030 // Monomials-specific implementations
00031 
00032 
00033 // Anonymous namespace for local helper functions
00034 namespace {
00035 
00036 void monomial_nodal_soln(const Elem* elem,
00037                          const Order order,
00038                          const std::vector<Number>& elem_soln,
00039                          std::vector<Number>&       nodal_soln,
00040                          const 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, MONOMIAL);
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         std::vector<Point> refspace_nodes;
00078         FEBase::get_refspace_nodes(elem_type,refspace_nodes);
00079         libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
00080 
00081         for (unsigned int n=0; n<n_nodes; n++)
00082           {
00083             libmesh_assert_equal_to (elem_soln.size(), n_sf);
00084 
00085             // Zero before summation
00086             nodal_soln[n] = 0;
00087 
00088             // u_i = Sum (alpha_i phi_i)
00089             for (unsigned int i=0; i<n_sf; i++)
00090               nodal_soln[n] += elem_soln[i] *
00091                 // FE<Dim,T>::shape(elem, order, i, mapped_point);
00092                 FEInterface::shape(Dim, fe_type, elem, i, refspace_nodes[n]);
00093           }
00094 
00095         return;
00096       } // default
00097     } // switch
00098 } // monomial_nodal_soln()
00099 
00100 
00101 
00102 
00103 unsigned int monomial_n_dofs(const ElemType t, const Order o)
00104 {
00105   switch (o)
00106     {
00107 
00108       // constant shape functions
00109       // no matter what shape there is only one DOF.
00110     case CONSTANT:
00111       return (t != INVALID_ELEM) ? 1 : 0;
00112 
00113 
00114       // Discontinuous linear shape functions
00115       // expressed in the monomials.
00116     case FIRST:
00117       {
00118         switch (t)
00119           {
00120           case NODEELEM:
00121             return 1;
00122 
00123           case EDGE2:
00124           case EDGE3:
00125           case EDGE4:
00126             return 2;
00127 
00128           case TRI3:
00129           case TRI6:
00130           case QUAD4:
00131           case QUAD8:
00132           case QUAD9:
00133             return 3;
00134 
00135           case TET4:
00136           case TET10:
00137           case HEX8:
00138           case HEX20:
00139           case HEX27:
00140           case PRISM6:
00141           case PRISM15:
00142           case PRISM18:
00143           case PYRAMID5:
00144           case PYRAMID13:
00145           case PYRAMID14:
00146             return 4;
00147 
00148           case INVALID_ELEM:
00149             return 0;
00150 
00151           default:
00152             libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
00153           }
00154       }
00155 
00156 
00157       // Discontinuous quadratic shape functions
00158       // expressed in the monomials.
00159     case SECOND:
00160       {
00161         switch (t)
00162           {
00163           case NODEELEM:
00164             return 1;
00165 
00166           case EDGE2:
00167           case EDGE3:
00168           case EDGE4:
00169             return 3;
00170 
00171           case TRI3:
00172           case TRI6:
00173           case QUAD4:
00174           case QUAD8:
00175           case QUAD9:
00176             return 6;
00177 
00178           case TET4:
00179           case TET10:
00180           case HEX8:
00181           case HEX20:
00182           case HEX27:
00183           case PRISM6:
00184           case PRISM15:
00185           case PRISM18:
00186           case PYRAMID5:
00187           case PYRAMID13:
00188           case PYRAMID14:
00189             return 10;
00190 
00191           case INVALID_ELEM:
00192             return 0;
00193 
00194           default:
00195             libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
00196           }
00197       }
00198 
00199 
00200       // Discontinuous cubic shape functions
00201       // expressed in the 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           case PYRAMID13:
00231           case PYRAMID14:
00232             return 20;
00233 
00234           case INVALID_ELEM:
00235             return 0;
00236 
00237           default:
00238             libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
00239           }
00240       }
00241 
00242 
00243       // Discontinuous quartic shape functions
00244       // expressed in the monomials.
00245     case FOURTH:
00246       {
00247         switch (t)
00248           {
00249           case NODEELEM:
00250             return 1;
00251 
00252           case EDGE2:
00253           case EDGE3:
00254             return 5;
00255 
00256           case TRI3:
00257           case TRI6:
00258           case QUAD4:
00259           case QUAD8:
00260           case QUAD9:
00261             return 15;
00262 
00263           case TET4:
00264           case TET10:
00265           case HEX8:
00266           case HEX20:
00267           case HEX27:
00268           case PRISM6:
00269           case PRISM15:
00270           case PRISM18:
00271           case PYRAMID5:
00272           case PYRAMID13:
00273           case PYRAMID14:
00274             return 35;
00275 
00276           case INVALID_ELEM:
00277             return 0;
00278 
00279           default:
00280             libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
00281           }
00282       }
00283 
00284 
00285     default:
00286       {
00287         const unsigned int order = static_cast<unsigned int>(o);
00288         switch (t)
00289           {
00290           case NODEELEM:
00291             return 1;
00292 
00293           case EDGE2:
00294           case EDGE3:
00295             return (order+1);
00296 
00297           case TRI3:
00298           case TRI6:
00299           case QUAD4:
00300           case QUAD8:
00301           case QUAD9:
00302             return (order+1)*(order+2)/2;
00303 
00304           case TET4:
00305           case TET10:
00306           case HEX8:
00307           case HEX20:
00308           case HEX27:
00309           case PRISM6:
00310           case PRISM15:
00311           case PRISM18:
00312           case PYRAMID5:
00313           case PYRAMID13:
00314           case PYRAMID14:
00315             return (order+1)*(order+2)*(order+3)/6;
00316 
00317           case INVALID_ELEM:
00318             return 0;
00319 
00320           default:
00321             libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
00322           }
00323       }
00324     }
00325 
00326   libmesh_error_msg("We'll never get here!");
00327   return 0;
00328 } // monomial_n_dofs()
00329 
00330 
00331 } // anonymous namespace
00332 
00333 
00334 
00335 
00336 
00337   // Do full-specialization for every dimension, instead
00338   // of explicit instantiation at the end of this file.
00339   // This could be macro-ified so that it fits on one line...
00340 template <>
00341 void FE<0,MONOMIAL>::nodal_soln(const Elem* elem,
00342                                 const Order order,
00343                                 const std::vector<Number>& elem_soln,
00344                                 std::vector<Number>& nodal_soln)
00345 { monomial_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/0); }
00346 
00347 template <>
00348 void FE<1,MONOMIAL>::nodal_soln(const Elem* elem,
00349                                 const Order order,
00350                                 const std::vector<Number>& elem_soln,
00351                                 std::vector<Number>& nodal_soln)
00352 { monomial_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/1); }
00353 
00354 template <>
00355 void FE<2,MONOMIAL>::nodal_soln(const Elem* elem,
00356                                 const Order order,
00357                                 const std::vector<Number>& elem_soln,
00358                                 std::vector<Number>& nodal_soln)
00359 { monomial_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/2); }
00360 
00361 template <>
00362 void FE<3,MONOMIAL>::nodal_soln(const Elem* elem,
00363                                 const Order order,
00364                                 const std::vector<Number>& elem_soln,
00365                                 std::vector<Number>& nodal_soln)
00366 { monomial_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/3); }
00367 
00368 
00369 // Full specialization of n_dofs() function for every dimension
00370 template <> unsigned int FE<0,MONOMIAL>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
00371 template <> unsigned int FE<1,MONOMIAL>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
00372 template <> unsigned int FE<2,MONOMIAL>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
00373 template <> unsigned int FE<3,MONOMIAL>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
00374 
00375 // Full specialization of n_dofs_at_node() function for every dimension.
00376 // Monomials have no dofs at nodes, only element dofs.
00377 template <> unsigned int FE<0,MONOMIAL>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
00378 template <> unsigned int FE<1,MONOMIAL>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
00379 template <> unsigned int FE<2,MONOMIAL>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
00380 template <> unsigned int FE<3,MONOMIAL>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
00381 
00382 // Full specialization of n_dofs_per_elem() function for every dimension.
00383 template <> unsigned int FE<0,MONOMIAL>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
00384 template <> unsigned int FE<1,MONOMIAL>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
00385 template <> unsigned int FE<2,MONOMIAL>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
00386 template <> unsigned int FE<3,MONOMIAL>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
00387 
00388 
00389 // Full specialization of get_continuity() function for every dimension.
00390 template <> FEContinuity FE<0,MONOMIAL>::get_continuity() const { return DISCONTINUOUS; }
00391 template <> FEContinuity FE<1,MONOMIAL>::get_continuity() const { return DISCONTINUOUS; }
00392 template <> FEContinuity FE<2,MONOMIAL>::get_continuity() const { return DISCONTINUOUS; }
00393 template <> FEContinuity FE<3,MONOMIAL>::get_continuity() const { return DISCONTINUOUS; }
00394 
00395 // Full specialization of is_hierarchic() function for every dimension.
00396 // The monomials are hierarchic!
00397 template <> bool FE<0,MONOMIAL>::is_hierarchic() const { return true; }
00398 template <> bool FE<1,MONOMIAL>::is_hierarchic() const { return true; }
00399 template <> bool FE<2,MONOMIAL>::is_hierarchic() const { return true; }
00400 template <> bool FE<3,MONOMIAL>::is_hierarchic() const { return true; }
00401 
00402 #ifdef LIBMESH_ENABLE_AMR
00403 
00404 // Full specialization of compute_constraints() function for 2D and
00405 // 3D only.  There are no constraints for discontinuous elements, so
00406 // we do nothing.
00407 template <> void FE<2,MONOMIAL>::compute_constraints (DofConstraints &, DofMap &, const unsigned int, const Elem*) {}
00408 template <> void FE<3,MONOMIAL>::compute_constraints (DofConstraints &, DofMap &, const unsigned int, const Elem*) {}
00409 
00410 #endif // #ifdef LIBMESH_ENABLE_AMR
00411 
00412 // Full specialization of shapes_need_reinit() function for every dimension.
00413 template <> bool FE<0,MONOMIAL>::shapes_need_reinit() const { return false; }
00414 template <> bool FE<1,MONOMIAL>::shapes_need_reinit() const { return false; }
00415 template <> bool FE<2,MONOMIAL>::shapes_need_reinit() const { return false; }
00416 template <> bool FE<3,MONOMIAL>::shapes_need_reinit() const { return false; }
00417 
00418 } // namespace libMesh