$extrastylesheet
fe_hermite_shape_2D.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 // C++ inlcludes
00020 
00021 // Local includes
00022 #include "libmesh/fe.h"
00023 #include "libmesh/elem.h"
00024 #include "libmesh/number_lookups.h"
00025 
00026 namespace
00027 {
00028 using namespace libMesh;
00029 
00030 // Compute the static coefficients for an element
00031 void hermite_compute_coefs(const Elem* elem,  std::vector<std::vector<Real> > & dxdxi
00032 #ifdef DEBUG
00033                            , std::vector<Real> & dxdeta, std::vector<Real> & dydxi
00034 #endif
00035                            )
00036 {
00037   const Order mapping_order        (elem->default_order());
00038   const ElemType mapping_elem_type (elem->type());
00039   const int n_mapping_shape_functions =
00040     FE<2,LAGRANGE>::n_shape_functions(mapping_elem_type,
00041                                       mapping_order);
00042 
00043   std::vector<Point> dofpt;
00044   dofpt.push_back(Point(-1,-1));
00045   dofpt.push_back(Point(1,1));
00046 
00047   for (int p = 0; p != 2; ++p)
00048     {
00049       dxdxi[0][p] = 0;
00050       dxdxi[1][p] = 0;
00051 #ifdef DEBUG
00052       dxdeta[p] = 0;
00053       dydxi[p] = 0;
00054 #endif
00055       for (int i = 0; i != n_mapping_shape_functions; ++i)
00056         {
00057           const Real ddxi = FE<2,LAGRANGE>::shape_deriv
00058             (mapping_elem_type, mapping_order, i, 0, dofpt[p]);
00059           const Real ddeta = FE<2,LAGRANGE>::shape_deriv
00060             (mapping_elem_type, mapping_order, i, 1, dofpt[p]);
00061 
00062           dxdxi[0][p] += elem->point(i)(0) * ddxi;
00063           dxdxi[1][p] += elem->point(i)(1) * ddeta;
00064           // dxdeta and dydxi should be 0!
00065 #ifdef DEBUG
00066           dxdeta[p] += elem->point(i)(0) * ddeta;
00067           dydxi[p] += elem->point(i)(1) * ddxi;
00068 #endif
00069         }
00070       // No singular elements!
00071       libmesh_assert(dxdxi[0][p]);
00072       libmesh_assert(dxdxi[1][p]);
00073       // No non-rectilinear or non-axis-aligned elements!
00074 #ifdef DEBUG
00075       libmesh_assert_less (std::abs(dxdeta[p]), 1e-9);
00076       libmesh_assert_less (std::abs(dydxi[p]), 1e-9);
00077 #endif
00078     }
00079 }
00080 
00081 
00082 
00083 Real hermite_bases_2D
00084 (std::vector<unsigned int> &bases1D,
00085  const std::vector<std::vector<Real> > &dxdxi,
00086  const Order &o,
00087  unsigned int i)
00088 {
00089   bases1D.clear();
00090   bases1D.resize(2,0);
00091   Real coef = 1.0;
00092 
00093   unsigned int e = o-3;
00094 
00095   // Nodes
00096   if (i < 16)
00097     {
00098       switch (i / 4)
00099         {
00100         case 0:
00101           break;
00102         case 1:
00103           bases1D[0] = 1;
00104           break;
00105         case 2:
00106           bases1D[0] = 1;
00107           bases1D[1] = 1;
00108           break;
00109         case 3:
00110           bases1D[1] = 1;
00111           break;
00112         default:
00113           libmesh_error_msg("Invalid basis node = " << i/4);
00114         }
00115 
00116       unsigned int basisnum = i%4;
00117       switch (basisnum)
00118         {
00119         case 0: // DoF = value at node
00120           coef = 1.0;
00121           break;
00122         case 1: // DoF = x derivative at node
00123           coef = dxdxi[0][bases1D[0]];
00124           bases1D[0] += 2; break;
00125         case 2: // DoF = y derivative at node
00126           coef = dxdxi[1][bases1D[1]];
00127           bases1D[1] += 2; break;
00128         case 3: // DoF = xy derivative at node
00129           coef = dxdxi[0][bases1D[0]] * dxdxi[1][bases1D[1]];
00130           bases1D[0] += 2; bases1D[1] += 2; break;
00131         default:
00132           libmesh_error_msg("Invalid basisnum = " << basisnum);
00133         }
00134     }
00135   // Edges
00136   else if (i < 16 + 4*2*e)
00137     {
00138       unsigned int basisnum = (i - 16) % (2*e);
00139       switch ((i - 16) / (2*e))
00140         {
00141         case 0:
00142           bases1D[0] = basisnum/2 + 4;
00143           bases1D[1] = basisnum%2*2;
00144           if (basisnum%2)
00145             coef = dxdxi[1][0];
00146           break;
00147         case 1:
00148           bases1D[0] = basisnum%2*2 + 1;
00149           bases1D[1] = basisnum/2 + 4;
00150           if (basisnum%2)
00151             coef = dxdxi[0][1];
00152           break;
00153         case 2:
00154           bases1D[0] = basisnum/2 + 4;
00155           bases1D[1] = basisnum%2*2 + 1;
00156           if (basisnum%2)
00157             coef = dxdxi[1][1];
00158           break;
00159         case 3:
00160           bases1D[0] = basisnum%2*2;
00161           bases1D[1] = basisnum/2 + 4;
00162           if (basisnum%2)
00163             coef = dxdxi[0][0];
00164           break;
00165         default:
00166           libmesh_error_msg("Invalid basisnum = " << basisnum);
00167         }
00168     }
00169   // Interior
00170   else
00171     {
00172       unsigned int basisnum = i - 16 - 8*e;
00173       bases1D[0] = square_number_row[basisnum]+4;
00174       bases1D[1] = square_number_column[basisnum]+4;
00175     }
00176 
00177   // No singular elements
00178   libmesh_assert(coef);
00179   return coef;
00180 }
00181 
00182 } // end anonymous namespace
00183 
00184 
00185 namespace libMesh
00186 {
00187 
00188 
00189 template <>
00190 Real FE<2,HERMITE>::shape(const ElemType,
00191                           const Order,
00192                           const unsigned int,
00193                           const Point&)
00194 {
00195   libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom.");
00196   return 0.;
00197 }
00198 
00199 
00200 
00201 template <>
00202 Real FE<2,HERMITE>::shape(const Elem* elem,
00203                           const Order order,
00204                           const unsigned int i,
00205                           const Point& p)
00206 {
00207   libmesh_assert(elem);
00208 
00209   std::vector<std::vector<Real> > dxdxi(2, std::vector<Real>(2, 0));
00210 
00211 #ifdef DEBUG
00212   std::vector<Real> dxdeta(2), dydxi(2);
00213 #endif
00214 
00215   hermite_compute_coefs(elem,dxdxi
00216 #ifdef DEBUG
00217                         ,dxdeta,dydxi
00218 #endif
00219                         );
00220 
00221   const ElemType type = elem->type();
00222 
00223   const Order totalorder = static_cast<Order>(order + elem->p_level());
00224 
00225   switch (type)
00226     {
00227     case QUAD4:
00228       libmesh_assert_less (totalorder, 4);
00229     case QUAD8:
00230     case QUAD9:
00231       {
00232         libmesh_assert_less (i, (totalorder+1u)*(totalorder+1u));
00233 
00234         std::vector<unsigned int> bases1D;
00235 
00236         Real coef = hermite_bases_2D(bases1D, dxdxi, totalorder, i);
00237 
00238         return coef * FEHermite<1>::hermite_raw_shape(bases1D[0],p(0)) *
00239           FEHermite<1>::hermite_raw_shape(bases1D[1],p(1));
00240       }
00241     default:
00242       libmesh_error_msg("ERROR: Unsupported element type = " << type);
00243     }
00244 
00245   libmesh_error_msg("We'll never get here!");
00246   return 0.;
00247 }
00248 
00249 
00250 
00251 template <>
00252 Real FE<2,HERMITE>::shape_deriv(const ElemType,
00253                                 const Order,
00254                                 const unsigned int,
00255                                 const unsigned int,
00256                                 const Point&)
00257 {
00258   libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom.");
00259   return 0.;
00260 }
00261 
00262 
00263 
00264 template <>
00265 Real FE<2,HERMITE>::shape_deriv(const Elem* elem,
00266                                 const Order order,
00267                                 const unsigned int i,
00268                                 const unsigned int j,
00269                                 const Point& p)
00270 {
00271   libmesh_assert(elem);
00272   libmesh_assert (j == 0 || j == 1);
00273 
00274   std::vector<std::vector<Real> > dxdxi(2, std::vector<Real>(2, 0));
00275 
00276 #ifdef DEBUG
00277   std::vector<Real> dxdeta(2), dydxi(2);
00278 #endif
00279 
00280   hermite_compute_coefs(elem,dxdxi
00281 #ifdef DEBUG
00282                         ,dxdeta,dydxi
00283 #endif
00284                         );
00285 
00286   const ElemType type = elem->type();
00287 
00288   const Order totalorder = static_cast<Order>(order + elem->p_level());
00289 
00290   switch (type)
00291     {
00292     case QUAD4:
00293       libmesh_assert_less (totalorder, 4);
00294     case QUAD8:
00295     case QUAD9:
00296       {
00297         libmesh_assert_less (i, (totalorder+1u)*(totalorder+1u));
00298 
00299         std::vector<unsigned int> bases1D;
00300 
00301         Real coef = hermite_bases_2D(bases1D, dxdxi, totalorder, i);
00302 
00303         switch (j)
00304           {
00305           case 0:
00306             return coef *
00307               FEHermite<1>::hermite_raw_shape_deriv(bases1D[0],p(0)) *
00308               FEHermite<1>::hermite_raw_shape(bases1D[1],p(1));
00309           case 1:
00310             return coef *
00311               FEHermite<1>::hermite_raw_shape(bases1D[0],p(0)) *
00312               FEHermite<1>::hermite_raw_shape_deriv(bases1D[1],p(1));
00313           default:
00314             libmesh_error_msg("Invalid derivative index j = " << j);
00315           }
00316       }
00317     default:
00318       libmesh_error_msg("ERROR: Unsupported element type = " << type);
00319     }
00320 
00321   libmesh_error_msg("We'll never get here!");
00322   return 0.;
00323 }
00324 
00325 
00326 
00327 template <>
00328 Real FE<2,HERMITE>::shape_second_deriv(const Elem* elem,
00329                                        const Order order,
00330                                        const unsigned int i,
00331                                        const unsigned int j,
00332                                        const Point& p)
00333 {
00334   libmesh_assert(elem);
00335   libmesh_assert (j == 0 || j == 1 || j == 2);
00336 
00337   std::vector<std::vector<Real> > dxdxi(2, std::vector<Real>(2, 0));
00338 
00339 #ifdef DEBUG
00340   std::vector<Real> dxdeta(2), dydxi(2);
00341 #endif
00342 
00343   hermite_compute_coefs(elem,dxdxi
00344 #ifdef DEBUG
00345                         ,dxdeta,dydxi
00346 #endif
00347                         );
00348 
00349   const ElemType type = elem->type();
00350 
00351   const Order totalorder = static_cast<Order>(order + elem->p_level());
00352 
00353   switch (type)
00354     {
00355     case QUAD4:
00356       libmesh_assert_less (totalorder, 4);
00357     case QUAD8:
00358     case QUAD9:
00359       {
00360         libmesh_assert_less (i, (totalorder+1u)*(totalorder+1u));
00361 
00362         std::vector<unsigned int> bases1D;
00363 
00364         Real coef = hermite_bases_2D(bases1D, dxdxi, totalorder, i);
00365 
00366         switch (j)
00367           {
00368           case 0:
00369             return coef *
00370               FEHermite<1>::hermite_raw_shape_second_deriv(bases1D[0],p(0)) *
00371               FEHermite<1>::hermite_raw_shape(bases1D[1],p(1));
00372           case 1:
00373             return coef *
00374               FEHermite<1>::hermite_raw_shape_deriv(bases1D[0],p(0)) *
00375               FEHermite<1>::hermite_raw_shape_deriv(bases1D[1],p(1));
00376           case 2:
00377             return coef *
00378               FEHermite<1>::hermite_raw_shape(bases1D[0],p(0)) *
00379               FEHermite<1>::hermite_raw_shape_second_deriv(bases1D[1],p(1));
00380           default:
00381             libmesh_error_msg("Invalid derivative index j = " << j);
00382           }
00383       }
00384     default:
00385       libmesh_error_msg("ERROR: Unsupported element type = " << type);
00386     }
00387 
00388   libmesh_error_msg("We'll never get here!");
00389   return 0.;
00390 }
00391 
00392 } // namespace libMesh