$extrastylesheet
fe_xyz_boundary.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 // C++ includes
00021 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
00022 #include <cmath> // for std::sqrt
00023 
00024 
00025 // Local includes
00026 #include "libmesh/libmesh_common.h"
00027 #include "libmesh/fe.h"
00028 #include "libmesh/quadrature.h"
00029 #include "libmesh/elem.h"
00030 #include "libmesh/libmesh_logging.h"
00031 
00032 namespace libMesh
00033 {
00034 
00035 
00036 
00037 
00038 //-------------------------------------------------------
00039 // Full specialization for 0D when this is a useless method
00040 template <>
00041 void FEXYZ<0>::reinit(const Elem*,
00042                       const unsigned int,
00043                       const Real,
00044                       const std::vector<Point>* const,
00045                       const std::vector<Real>* const)
00046 {
00047   libmesh_error_msg("ERROR: This method only makes sense for 2D/3D elements!");
00048 }
00049 
00050 
00051 //-------------------------------------------------------
00052 // Full specialization for 1D when this is a useless method
00053 template <>
00054 void FEXYZ<1>::reinit(const Elem*,
00055                       const unsigned int,
00056                       const Real,
00057                       const std::vector<Point>* const,
00058                       const std::vector<Real>* const)
00059 {
00060   libmesh_error_msg("ERROR: This method only makes sense for 2D/3D elements!");
00061 }
00062 
00063 
00064 //-------------------------------------------------------
00065 // Method for 2D, 3D
00066 template <unsigned int Dim>
00067 void FEXYZ<Dim>::reinit(const Elem* elem,
00068                         const unsigned int s,
00069                         const Real,
00070                         const std::vector<Point>* const pts,
00071                         const std::vector<Real>* const weights)
00072 {
00073   libmesh_assert(elem);
00074   libmesh_assert (this->qrule != NULL || pts != NULL);
00075   // We don't do this for 1D elements!
00076   libmesh_assert_not_equal_to (Dim, 1);
00077 
00078   // Build the side of interest
00079   const UniquePtr<Elem> side(elem->build_side(s));
00080 
00081   // Initialize the shape functions at the user-specified
00082   // points
00083   if (pts != NULL)
00084     {
00085       // We can't get away without recomputing shape functions next
00086       // time
00087       this->shapes_on_quadrature = false;
00088 
00089       // Set the element type
00090       this->elem_type = elem->type();
00091 
00092       // Initialize the face shape functions
00093       this->_fe_map->template init_face_shape_functions<Dim>(*pts,  side.get());
00094       if (weights != NULL)
00095         {
00096           this->compute_face_values (elem, side.get(), *weights);
00097         }
00098       else
00099         {
00100           std::vector<Real> dummy_weights (pts->size(), 1.);
00101           // Compute data on the face for integration
00102           this->compute_face_values (elem, side.get(), dummy_weights);
00103         }
00104     }
00105   else
00106     {
00107       // initialize quadrature rule
00108       this->qrule->init(side->type(), elem->p_level());
00109 
00110       {
00111         // Set the element type
00112         this->elem_type = elem->type();
00113 
00114         // Initialize the face shape functions
00115         this->_fe_map->template init_face_shape_functions<Dim>(this->qrule->get_points(),  side.get());
00116       }
00117       // We can't get away without recomputing shape functions next
00118       // time
00119       this->shapes_on_quadrature = false;
00120       // Compute data on the face for integration
00121       this->compute_face_values (elem, side.get(), this->qrule->get_weights());
00122     }
00123 }
00124 
00125 
00126 
00127 template <unsigned int Dim>
00128 void FEXYZ<Dim>::compute_face_values(const Elem* elem,
00129                                      const Elem* side,
00130                                      const std::vector<Real>& qw)
00131 {
00132   libmesh_assert(elem);
00133   libmesh_assert(side);
00134 
00135   START_LOG("compute_face_values()", "FEXYZ");
00136 
00137   // The number of quadrature points.
00138   const std::size_t n_qp = qw.size();
00139 
00140   // Number of shape functions in the finite element approximation
00141   // space.
00142   const unsigned int n_approx_shape_functions =
00143     this->n_shape_functions(this->get_type(),
00144                             this->get_order());
00145 
00146   // Resize the shape functions and their gradients
00147   this->phi.resize    (n_approx_shape_functions);
00148   this->dphi.resize   (n_approx_shape_functions);
00149   this->dphidx.resize (n_approx_shape_functions);
00150   this->dphidy.resize (n_approx_shape_functions);
00151   this->dphidz.resize (n_approx_shape_functions);
00152 
00153   for (unsigned int i=0; i<n_approx_shape_functions; i++)
00154     {
00155       this->phi[i].resize    (n_qp);
00156       this->dphi[i].resize   (n_qp);
00157       this->dphidx[i].resize (n_qp);
00158       this->dphidy[i].resize (n_qp);
00159       this->dphidz[i].resize (n_qp);
00160     }
00161 
00162   this->_fe_map->compute_face_map(this->dim, qw, side);
00163 
00164   const std::vector<libMesh::Point>& xyz = this->_fe_map->get_xyz();
00165 
00166   switch (this->dim)
00167     {
00168       // A 2D finite element living in either 2D or 3D space.
00169       // This means the boundary is a 1D finite element, i.e.
00170       // and EDGE2 or EDGE3.
00171     case 2:
00172       {
00173         // compute the shape function values & gradients
00174         for (unsigned int i=0; i<n_approx_shape_functions; i++)
00175           for (std::size_t p=0; p<n_qp; p++)
00176             {
00177               this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz[p]);
00178 
00179               this->dphi[i][p](0) =
00180                 this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz[p]);
00181 
00182               this->dphi[i][p](1) =
00183                 this->dphidy[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 1, xyz[p]);
00184 
00185 #if LIBMESH_DIM == 3
00186               this->dphi[i][p](2) = // can only assign to the Z component if LIBMESH_DIM==3
00187 #endif
00188                 this->dphidz[i][p] = 0.;
00189             }
00190 
00191         // done computing face values
00192         break;
00193       }
00194 
00195       // A 3D finite element living in 3D space.
00196     case 3:
00197       {
00198         // compute the shape function values & gradients
00199         for (unsigned int i=0; i<n_approx_shape_functions; i++)
00200           for (std::size_t p=0; p<n_qp; p++)
00201             {
00202               this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz[p]);
00203 
00204               this->dphi[i][p](0) =
00205                 this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz[p]);
00206 
00207               this->dphi[i][p](1) =
00208                 this->dphidy[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 1, xyz[p]);
00209 
00210               this->dphi[i][p](2) =
00211                 this->dphidz[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 2, xyz[p]);
00212             }
00213 
00214         // done computing face values
00215         break;
00216       }
00217 
00218     default:
00219       libmesh_error_msg("Invalid dim " << this->dim);
00220     }
00221 
00222   STOP_LOG("compute_face_values()", "FEXYZ");
00223 }
00224 
00225 
00226 
00227 
00228 //--------------------------------------------------------------
00229 // Explicit instantiations (doesn't make sense in 1D!) using fe_macro.h's macro
00230 
00231 template class FEXYZ<2>;
00232 template class FEXYZ<3>;
00233 
00234 
00235 } // namespace libMesh