$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_ABSTRACT_H 00021 #define LIBMESH_FE_ABSTRACT_H 00022 00023 // Local includes 00024 #include "libmesh/reference_counted_object.h" 00025 #include "libmesh/point.h" 00026 #include "libmesh/vector_value.h" 00027 #include "libmesh/enum_elem_type.h" 00028 #include "libmesh/fe_type.h" 00029 #include "libmesh/auto_ptr.h" 00030 #include "libmesh/fe_map.h" 00031 00032 // C++ includes 00033 #include <cstddef> 00034 #include <vector> 00035 00036 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00037 #include "libmesh/tensor_value.h" 00038 #endif 00039 00040 namespace libMesh 00041 { 00042 00043 00044 // forward declarations 00045 template <typename T> class DenseMatrix; 00046 template <typename T> class DenseVector; 00047 class BoundaryInfo; 00048 class DofConstraints; 00049 class DofMap; 00050 class Elem; 00051 class MeshBase; 00052 template <typename T> class NumericVector; 00053 class QBase; 00054 00055 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS 00056 class NodeConstraints; 00057 #endif 00058 00059 #ifdef LIBMESH_ENABLE_PERIODIC 00060 class PeriodicBoundaries; 00061 class PointLocatorBase; 00062 #endif 00063 00064 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 00065 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map> 00066 class InfFE; 00067 #endif 00068 00069 00070 00095 // ------------------------------------------------------------ 00096 // FEAbstract class definition 00097 class FEAbstract : public ReferenceCountedObject<FEAbstract> 00098 { 00099 protected: 00100 00106 FEAbstract (const unsigned int dim, 00107 const FEType& fet); 00108 00109 public: 00110 00114 virtual ~FEAbstract(); 00115 00121 static UniquePtr<FEAbstract> build (const unsigned int dim, 00122 const FEType& type); 00123 00138 virtual void reinit (const Elem* elem, 00139 const std::vector<Point>* const pts = NULL, 00140 const std::vector<Real>* const weights = NULL) = 0; 00141 00150 virtual void reinit (const Elem* elem, 00151 const unsigned int side, 00152 const Real tolerance = TOLERANCE, 00153 const std::vector<Point>* const pts = NULL, 00154 const std::vector<Real>* const weights = NULL) = 0; 00155 00164 virtual void edge_reinit (const Elem* elem, 00165 const unsigned int edge, 00166 const Real tolerance = TOLERANCE, 00167 const std::vector<Point>* pts = NULL, 00168 const std::vector<Real>* weights = NULL) = 0; 00169 00174 virtual void side_map (const Elem* elem, 00175 const Elem* side, 00176 const unsigned int s, 00177 const std::vector<Point>& reference_side_points, 00178 std::vector<Point>& reference_points) = 0; 00179 00187 static bool on_reference_element(const Point& p, 00188 const ElemType t, 00189 const Real eps = TOLERANCE); 00194 static void get_refspace_nodes(const ElemType t, 00195 std::vector<Point>& nodes); 00196 00197 00198 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS 00199 00203 static void compute_node_constraints (NodeConstraints &constraints, 00204 const Elem* elem); 00205 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS 00206 00207 #ifdef LIBMESH_ENABLE_PERIODIC 00208 00209 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS 00210 00214 static void compute_periodic_node_constraints (NodeConstraints &constraints, 00215 const PeriodicBoundaries &boundaries, 00216 const MeshBase& mesh, 00217 const PointLocatorBase* point_locator, 00218 const Elem* elem); 00219 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS 00220 00221 #endif // LIBMESH_ENABLE_PERIODIC 00222 00227 const std::vector<Point>& get_xyz() const 00228 { return this->_fe_map->get_xyz(); } 00229 00234 const std::vector<Real>& get_JxW() const 00235 { return this->_fe_map->get_JxW(); } 00236 00241 const std::vector<RealGradient>& get_dxyzdxi() const 00242 { return this->_fe_map->get_dxyzdxi(); } 00243 00248 const std::vector<RealGradient>& get_dxyzdeta() const 00249 { return this->_fe_map->get_dxyzdeta(); } 00250 00255 const std::vector<RealGradient>& get_dxyzdzeta() const 00256 { return _fe_map->get_dxyzdzeta(); } 00257 00261 const std::vector<RealGradient>& get_d2xyzdxi2() const 00262 { return this->_fe_map->get_d2xyzdxi2(); } 00263 00267 const std::vector<RealGradient>& get_d2xyzdeta2() const 00268 { return this->_fe_map->get_d2xyzdeta2(); } 00269 00270 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00271 00275 const std::vector<RealGradient>& get_d2xyzdzeta2() const 00276 { return this->_fe_map->get_d2xyzdzeta2(); } 00277 00278 #endif 00279 00283 const std::vector<RealGradient>& get_d2xyzdxideta() const 00284 { return this->_fe_map->get_d2xyzdxideta(); } 00285 00286 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00287 00291 const std::vector<RealGradient>& get_d2xyzdxidzeta() const 00292 { return this->_fe_map->get_d2xyzdxidzeta(); } 00293 00297 const std::vector<RealGradient>& get_d2xyzdetadzeta() const 00298 { return this->_fe_map->get_d2xyzdetadzeta(); } 00299 00300 #endif 00301 00306 const std::vector<Real>& get_dxidx() const 00307 { return this->_fe_map->get_dxidx(); } 00308 00313 const std::vector<Real>& get_dxidy() const 00314 { return this->_fe_map->get_dxidy(); } 00315 00320 const std::vector<Real>& get_dxidz() const 00321 { return this->_fe_map->get_dxidz(); } 00322 00327 const std::vector<Real>& get_detadx() const 00328 { return this->_fe_map->get_detadx(); } 00329 00334 const std::vector<Real>& get_detady() const 00335 { return this->_fe_map->get_detady(); } 00336 00341 const std::vector<Real>& get_detadz() const 00342 { return this->_fe_map->get_detadz(); } 00343 00348 const std::vector<Real>& get_dzetadx() const 00349 { return this->_fe_map->get_dzetadx(); } 00350 00355 const std::vector<Real>& get_dzetady() const 00356 { return this->_fe_map->get_dzetady(); } 00357 00362 const std::vector<Real>& get_dzetadz() const 00363 { return this->_fe_map->get_dzetadz(); } 00364 00368 const std::vector<std::vector<Point> >& get_tangents() const 00369 { return this->_fe_map->get_tangents(); } 00370 00374 const std::vector<Point>& get_normals() const 00375 { return this->_fe_map->get_normals(); } 00376 00380 const std::vector<Real>& get_curvatures() const 00381 { return this->_fe_map->get_curvatures();} 00382 00387 virtual void attach_quadrature_rule (QBase* q) = 0; 00388 00394 virtual unsigned int n_shape_functions () const = 0; 00395 00400 virtual unsigned int n_quadrature_points () const = 0; 00401 00407 ElemType get_type() const { return elem_type; } 00408 00413 unsigned int get_p_level() const { return _p_level; } 00414 00418 FEType get_fe_type() const { return fe_type; } 00419 00423 Order get_order() const { return static_cast<Order>(fe_type.order + _p_level); } 00424 00428 virtual FEContinuity get_continuity() const = 0; 00429 00434 virtual bool is_hierarchic() const = 0; 00435 00439 FEFamily get_family() const { return fe_type.family; } 00440 00444 const FEMap& get_fe_map() const { return *_fe_map.get(); } 00445 00449 void print_JxW(std::ostream& os) const; 00450 00456 virtual void print_phi(std::ostream& os) const =0; 00457 00463 virtual void print_dphi(std::ostream& os) const =0; 00464 00465 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 00466 00472 virtual void print_d2phi(std::ostream& os) const =0; 00473 00474 #endif 00475 00480 void print_xyz(std::ostream& os) const; 00481 00485 void print_info(std::ostream& os) const; 00486 00490 friend std::ostream& operator << (std::ostream& os, const FEAbstract& fe); 00491 00492 00493 protected: 00494 00507 virtual void compute_shape_functions(const Elem*, const std::vector<Point>& ) =0; 00508 00509 UniquePtr<FEMap> _fe_map; 00510 00511 00515 const unsigned int dim; 00516 00521 mutable bool calculations_started; 00522 00526 mutable bool calculate_phi; 00527 00531 mutable bool calculate_dphi; 00532 00536 mutable bool calculate_d2phi; 00537 00541 mutable bool calculate_curl_phi; 00542 00546 mutable bool calculate_div_phi; 00547 00551 mutable bool calculate_dphiref; 00552 00553 00558 const FEType fe_type; 00559 00564 ElemType elem_type; 00565 00570 unsigned int _p_level; 00571 00575 QBase* qrule; 00576 00581 bool shapes_on_quadrature; 00582 00589 virtual bool shapes_need_reinit() const = 0; 00590 00591 }; 00592 00593 00594 00595 00596 // ------------------------------------------------------------ 00597 // FEAbstract class inline members 00598 inline 00599 FEAbstract::FEAbstract(const unsigned int d, 00600 const FEType& fet) : 00601 _fe_map( FEMap::build(fet) ), 00602 dim(d), 00603 calculations_started(false), 00604 calculate_phi(false), 00605 calculate_dphi(false), 00606 calculate_d2phi(false), 00607 calculate_curl_phi(false), 00608 calculate_div_phi(false), 00609 calculate_dphiref(false), 00610 fe_type(fet), 00611 elem_type(INVALID_ELEM), 00612 _p_level(0), 00613 qrule(NULL), 00614 shapes_on_quadrature(false) 00615 { 00616 } 00617 00618 00619 inline 00620 FEAbstract::~FEAbstract() 00621 { 00622 } 00623 00624 } // namespace libMesh 00625 00626 #endif // LIBMESH_FE_ABSTRACT_H