$extrastylesheet
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 #ifndef LIBMESH_FE_BASE_H 00021 #define LIBMESH_FE_BASE_H 00022 00023 // Local includes 00024 #include "libmesh/libmesh_common.h" 00025 #include "libmesh/auto_ptr.h" 00026 #include "libmesh/compare_types.h" 00027 #include "libmesh/enum_elem_type.h" 00028 #include "libmesh/fe_abstract.h" 00029 #include "libmesh/fe_transformation_base.h" 00030 #include "libmesh/fe_type.h" 00031 #include "libmesh/point.h" 00032 #include "libmesh/reference_counted_object.h" 00033 #include "libmesh/tensor_tools.h" 00034 #include "libmesh/type_n_tensor.h" 00035 #include "libmesh/vector_value.h" 00036 00037 // C++ includes 00038 #include <cstddef> 00039 #include <vector> 00040 00041 00042 namespace libMesh 00043 { 00044 00045 00046 // forward declarations 00047 template <typename T> class DenseMatrix; 00048 template <typename T> class DenseVector; 00049 class BoundaryInfo; 00050 class DofConstraints; 00051 class DofMap; 00052 class Elem; 00053 class MeshBase; 00054 template <typename T> class NumericVector; 00055 class QBase; 00056 template <typename T> class FETransformationBase; 00057 00058 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS 00059 class NodeConstraints; 00060 #endif 00061 00062 #ifdef LIBMESH_ENABLE_PERIODIC 00063 class PeriodicBoundaries; 00064 class PointLocatorBase; 00065 #endif 00066 00067 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00068 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00069 class InfFE; 00070 #endif 00071 00099 // ------------------------------------------------------------ 00100 // FEBase class definition 00101 template <typename OutputType> 00102 class FEGenericBase : public FEAbstract 00103 { 00104 protected: 00105 00111 FEGenericBase (const unsigned int dim, 00112 const FEType& fet); 00113 00114 public: 00115 00119 virtual ~FEGenericBase(); 00120 00129 static UniquePtr<FEGenericBase> build (const unsigned int dim, 00130 const FEType& type); 00131 00136 typedef OutputType OutputShape; 00137 typedef typename TensorTools::IncrementRank<OutputShape>::type OutputGradient; 00138 typedef typename TensorTools::IncrementRank<OutputGradient>::type OutputTensor; 00139 typedef typename TensorTools::DecrementRank<OutputShape>::type OutputDivergence; 00140 typedef typename TensorTools::MakeNumber<OutputShape>::type OutputNumber; 00141 typedef typename TensorTools::IncrementRank<OutputNumber>::type OutputNumberGradient; 00142 typedef typename TensorTools::IncrementRank<OutputNumberGradient>::type OutputNumberTensor; 00143 typedef typename TensorTools::DecrementRank<OutputNumber>::type OutputNumberDivergence; 00144 00145 00146 00147 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00148 00157 static UniquePtr<FEGenericBase> build_InfFE (const unsigned int dim, 00158 const FEType& type); 00159 00160 #endif 00161 00162 #ifdef LIBMESH_ENABLE_AMR 00163 00170 static void compute_proj_constraints (DofConstraints &constraints, 00171 DofMap &dof_map, 00172 const unsigned int variable_number, 00173 const Elem* elem); 00174 00180 static void coarsened_dof_values(const NumericVector<Number> &global_vector, 00181 const DofMap &dof_map, 00182 const Elem *coarse_elem, 00183 DenseVector<Number> &coarse_dofs, 00184 const unsigned int var, 00185 const bool use_old_dof_indices = false); 00186 00187 #endif // #ifdef LIBMESH_ENABLE_AMR 00188 00189 #ifdef LIBMESH_ENABLE_PERIODIC 00190 00196 static void compute_periodic_constraints (DofConstraints &constraints, 00197 DofMap &dof_map, 00198 const PeriodicBoundaries &boundaries, 00199 const MeshBase& mesh, 00200 const PointLocatorBase* point_locator, 00201 const unsigned int variable_number, 00202 const Elem* elem); 00203 00204 #endif // LIBMESH_ENABLE_PERIODIC 00205 00210 const std::vector<std::vector<OutputShape> >& get_phi() const 00211 { libmesh_assert(!calculations_started || calculate_phi); 00212 calculate_phi = true; return phi; } 00213 00218 const std::vector<std::vector<OutputGradient> >& get_dphi() const 00219 { libmesh_assert(!calculations_started || calculate_dphi); 00220 calculate_dphi = calculate_dphiref = true; return dphi; } 00221 00226 const std::vector<std::vector<OutputShape> >& get_curl_phi() const 00227 { libmesh_assert(!calculations_started || calculate_curl_phi); 00228 calculate_curl_phi = calculate_dphiref = true; return curl_phi; } 00229 00234 const std::vector<std::vector<OutputDivergence> >& get_div_phi() const 00235 { libmesh_assert(!calculations_started || calculate_div_phi); 00236 calculate_div_phi = calculate_dphiref = true; return div_phi; } 00237 00242 const std::vector<std::vector<OutputShape> >& get_dphidx() const 00243 { libmesh_assert(!calculations_started || calculate_dphi); 00244 calculate_dphi = calculate_dphiref = true; return dphidx; } 00245 00250 const std::vector<std::vector<OutputShape> >& get_dphidy() const 00251 { libmesh_assert(!calculations_started || calculate_dphi); 00252 calculate_dphi = calculate_dphiref = true; return dphidy; } 00253 00258 const std::vector<std::vector<OutputShape> >& get_dphidz() const 00259 { libmesh_assert(!calculations_started || calculate_dphi); 00260 calculate_dphi = calculate_dphiref = true; return dphidz; } 00261 00266 const std::vector<std::vector<OutputShape> >& get_dphidxi() const 00267 { libmesh_assert(!calculations_started || calculate_dphiref); 00268 calculate_dphiref = true; return dphidxi; } 00269 00274 const std::vector<std::vector<OutputShape> >& get_dphideta() const 00275 { libmesh_assert(!calculations_started || calculate_dphiref); 00276 calculate_dphiref = true; return dphideta; } 00277 00282 const std::vector<std::vector<OutputShape> >& get_dphidzeta() const 00283 { libmesh_assert(!calculations_started || calculate_dphiref); 00284 calculate_dphiref = true; return dphidzeta; } 00285 00286 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00287 00292 const std::vector<std::vector<OutputTensor> >& get_d2phi() const 00293 { libmesh_assert(!calculations_started || calculate_d2phi); 00294 calculate_d2phi = calculate_dphiref = true; return d2phi; } 00295 00300 const std::vector<std::vector<OutputShape> >& get_d2phidx2() const 00301 { libmesh_assert(!calculations_started || calculate_d2phi); 00302 calculate_d2phi = calculate_dphiref = true; return d2phidx2; } 00303 00308 const std::vector<std::vector<OutputShape> >& get_d2phidxdy() const 00309 { libmesh_assert(!calculations_started || calculate_d2phi); 00310 calculate_d2phi = calculate_dphiref = true; return d2phidxdy; } 00311 00316 const std::vector<std::vector<OutputShape> >& get_d2phidxdz() const 00317 { libmesh_assert(!calculations_started || calculate_d2phi); 00318 calculate_d2phi = calculate_dphiref = true; return d2phidxdz; } 00319 00324 const std::vector<std::vector<OutputShape> >& get_d2phidy2() const 00325 { libmesh_assert(!calculations_started || calculate_d2phi); 00326 calculate_d2phi = calculate_dphiref = true; return d2phidy2; } 00327 00332 const std::vector<std::vector<OutputShape> >& get_d2phidydz() const 00333 { libmesh_assert(!calculations_started || calculate_d2phi); 00334 calculate_d2phi = calculate_dphiref = true; return d2phidydz; } 00335 00340 const std::vector<std::vector<OutputShape> >& get_d2phidz2() const 00341 { libmesh_assert(!calculations_started || calculate_d2phi); 00342 calculate_d2phi = calculate_dphiref = true; return d2phidz2; } 00343 00348 const std::vector<std::vector<OutputShape> >& get_d2phidxi2() const 00349 { libmesh_assert(!calculations_started || calculate_d2phi); 00350 calculate_d2phi = calculate_dphiref = true; return d2phidxi2; } 00351 00356 const std::vector<std::vector<OutputShape> >& get_d2phidxideta() const 00357 { libmesh_assert(!calculations_started || calculate_d2phi); 00358 calculate_d2phi = calculate_dphiref = true; return d2phidxideta; } 00359 00364 const std::vector<std::vector<OutputShape> >& get_d2phidxidzeta() const 00365 { libmesh_assert(!calculations_started || calculate_d2phi); 00366 calculate_d2phi = calculate_dphiref = true; return d2phidxidzeta; } 00367 00372 const std::vector<std::vector<OutputShape> >& get_d2phideta2() const 00373 { libmesh_assert(!calculations_started || calculate_d2phi); 00374 calculate_d2phi = calculate_dphiref = true; return d2phideta2; } 00375 00380 const std::vector<std::vector<OutputShape> >& get_d2phidetadzeta() const 00381 { libmesh_assert(!calculations_started || calculate_d2phi); 00382 calculate_d2phi = calculate_dphiref = true; return d2phidetadzeta; } 00383 00388 const std::vector<std::vector<OutputShape> >& get_d2phidzeta2() const 00389 { libmesh_assert(!calculations_started || calculate_d2phi); 00390 calculate_d2phi = calculate_dphiref = true; return d2phidzeta2; } 00391 00392 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES 00393 00394 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00395 00406 const std::vector<OutputGradient>& get_dphase() const 00407 { return dphase; } 00408 00409 00422 const std::vector<Real>& get_Sobolev_weight() const 00423 { return weight; } 00424 00430 const std::vector<RealGradient>& get_Sobolev_dweight() const 00431 { return dweight; } 00432 00433 #endif 00434 00435 00439 void print_phi(std::ostream& os) const; 00440 00445 void print_dphi(std::ostream& os) const; 00446 00447 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00448 00453 void print_d2phi(std::ostream& os) const; 00454 00455 #endif 00456 00457 00458 protected: 00459 00460 00461 00462 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00463 00469 virtual void init_base_shape_functions(const std::vector<Point>& qp, 00470 const Elem* e) = 0; 00471 00472 #endif 00473 00484 virtual void compute_shape_functions(const Elem* elem, const std::vector<Point>& qp); 00485 00490 UniquePtr<FETransformationBase<OutputType> > _fe_trans; 00491 00495 std::vector<std::vector<OutputShape> > phi; 00496 00500 std::vector<std::vector<OutputGradient> > dphi; 00501 00505 std::vector<std::vector<OutputShape> > curl_phi; 00506 00510 std::vector<std::vector<OutputDivergence> > div_phi; 00511 00515 std::vector<std::vector<OutputShape> > dphidxi; 00516 00520 std::vector<std::vector<OutputShape> > dphideta; 00521 00525 std::vector<std::vector<OutputShape> > dphidzeta; 00526 00530 std::vector<std::vector<OutputShape> > dphidx; 00531 00535 std::vector<std::vector<OutputShape> > dphidy; 00536 00540 std::vector<std::vector<OutputShape> > dphidz; 00541 00542 00543 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00544 00548 std::vector<std::vector<OutputTensor> > d2phi; 00549 00553 std::vector<std::vector<OutputShape> > d2phidxi2; 00554 00558 std::vector<std::vector<OutputShape> > d2phidxideta; 00559 00563 std::vector<std::vector<OutputShape> > d2phidxidzeta; 00564 00568 std::vector<std::vector<OutputShape> > d2phideta2; 00569 00573 std::vector<std::vector<OutputShape> > d2phidetadzeta; 00574 00578 std::vector<std::vector<OutputShape> > d2phidzeta2; 00579 00583 std::vector<std::vector<OutputShape> > d2phidx2; 00584 00588 std::vector<std::vector<OutputShape> > d2phidxdy; 00589 00593 std::vector<std::vector<OutputShape> > d2phidxdz; 00594 00598 std::vector<std::vector<OutputShape> > d2phidy2; 00599 00603 std::vector<std::vector<OutputShape> > d2phidydz; 00604 00608 std::vector<std::vector<OutputShape> > d2phidz2; 00609 00610 #endif 00611 00612 00613 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00614 00615 //-------------------------------------------------------------- 00616 /* protected members for infinite elements, which are accessed 00617 * from the outside through some inline functions 00618 */ 00619 00620 00626 std::vector<OutputGradient> dphase; 00627 00633 std::vector<RealGradient> dweight; 00634 00640 std::vector<Real> weight; 00641 00642 #endif 00643 00644 private: 00645 00646 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00647 00653 template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map> 00654 friend class InfFE; 00655 00656 #endif 00657 00658 00659 }; 00660 00661 00662 // Typedefs for convenience and backwards compatibility 00663 typedef FEGenericBase<Real> FEBase; 00664 typedef FEGenericBase<RealGradient> FEVectorBase; 00665 00666 00667 00668 00669 // ------------------------------------------------------------ 00670 // FEGenericBase class inline members 00671 template <typename OutputType> 00672 inline 00673 FEGenericBase<OutputType>::FEGenericBase(const unsigned int d, 00674 const FEType& fet) : 00675 FEAbstract(d,fet), 00676 _fe_trans( FETransformationBase<OutputType>::build(fet) ), 00677 phi(), 00678 dphi(), 00679 curl_phi(), 00680 div_phi(), 00681 dphidxi(), 00682 dphideta(), 00683 dphidzeta(), 00684 dphidx(), 00685 dphidy(), 00686 dphidz() 00687 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00688 ,d2phi(), 00689 d2phidxi2(), 00690 d2phidxideta(), 00691 d2phidxidzeta(), 00692 d2phideta2(), 00693 d2phidetadzeta(), 00694 d2phidzeta2(), 00695 d2phidx2(), 00696 d2phidxdy(), 00697 d2phidxdz(), 00698 d2phidy2(), 00699 d2phidydz(), 00700 d2phidz2() 00701 #endif 00702 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00703 ,dphase(), 00704 dweight(), 00705 weight() 00706 #endif 00707 { 00708 } 00709 00710 00711 00712 template <typename OutputType> 00713 inline 00714 FEGenericBase<OutputType>::~FEGenericBase() 00715 { 00716 } 00717 00718 } // namespace libMesh 00719 00720 #endif // LIBMESH_FE_BASE_H