fe_abstract.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_FE_ABSTRACT_H
21 #define LIBMESH_FE_ABSTRACT_H
22 
23 // Local includes
25 #include "libmesh/point.h"
26 #include "libmesh/vector_value.h"
27 #include "libmesh/fe_type.h"
28 #include "libmesh/auto_ptr.h" // deprecated
29 #include "libmesh/fe_map.h"
30 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
31 #include "libmesh/tensor_value.h"
32 #endif
33 
34 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS
35 namespace libMesh
36 {
37 enum ElemType : int;
38 }
39 #else
40 #include "libmesh/enum_elem_type.h"
41 #endif
42 
43 // C++ includes
44 #include <cstddef>
45 #include <vector>
46 #include <memory>
47 
48 namespace libMesh
49 {
50 
51 
52 // forward declarations
53 template <typename T> class DenseMatrix;
54 template <typename T> class DenseVector;
55 class BoundaryInfo;
56 class DofConstraints;
57 class DofMap;
58 class Elem;
59 class MeshBase;
60 template <typename T> class NumericVector;
61 class QBase;
62 
63 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
64 class NodeConstraints;
65 #endif
66 
67 #ifdef LIBMESH_ENABLE_PERIODIC
68 class PeriodicBoundaries;
69 class PointLocatorBase;
70 #endif
71 
72 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
73 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
74 class InfFE;
75 #endif
76 
77 
78 
101 class FEAbstract : public ReferenceCountedObject<FEAbstract>
102 {
103 protected:
104 
110  FEAbstract (const unsigned int dim,
111  const FEType & fet);
112 
113 public:
114 
118  virtual ~FEAbstract();
119 
126  static std::unique_ptr<FEAbstract> build (const unsigned int dim,
127  const FEType & type);
128 
143  virtual void reinit (const Elem * elem,
144  const std::vector<Point> * const pts = nullptr,
145  const std::vector<Real> * const weights = nullptr) = 0;
146 
155  virtual void reinit (const Elem * elem,
156  const unsigned int side,
157  const Real tolerance = TOLERANCE,
158  const std::vector<Point> * const pts = nullptr,
159  const std::vector<Real> * const weights = nullptr) = 0;
160 
169  virtual void edge_reinit (const Elem * elem,
170  const unsigned int edge,
171  const Real tolerance = TOLERANCE,
172  const std::vector<Point> * pts = nullptr,
173  const std::vector<Real> * weights = nullptr) = 0;
174 
179  virtual void side_map (const Elem * elem,
180  const Elem * side,
181  const unsigned int s,
182  const std::vector<Point> & reference_side_points,
183  std::vector<Point> & reference_points) = 0;
184 
192  static bool on_reference_element(const Point & p,
193  const ElemType t,
194  const Real eps = TOLERANCE);
199  static void get_refspace_nodes(const ElemType t,
200  std::vector<Point> & nodes);
201 
202 
203 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
204 
208  static void compute_node_constraints (NodeConstraints & constraints,
209  const Elem * elem);
210 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
211 
212 #ifdef LIBMESH_ENABLE_PERIODIC
213 
214 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
215 
219  static void compute_periodic_node_constraints (NodeConstraints & constraints,
220  const PeriodicBoundaries & boundaries,
221  const MeshBase & mesh,
222  const PointLocatorBase * point_locator,
223  const Elem * elem);
224 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
225 
226 #endif // LIBMESH_ENABLE_PERIODIC
227 
231  unsigned int get_dim() const
232  { return dim; }
233 
238  const std::vector<Point> & get_xyz() const
239  { return this->_fe_map->get_xyz(); }
240 
245  const std::vector<Real> & get_JxW() const
246  { return this->_fe_map->get_JxW(); }
247 
252  const std::vector<RealGradient> & get_dxyzdxi() const
253  { return this->_fe_map->get_dxyzdxi(); }
254 
259  const std::vector<RealGradient> & get_dxyzdeta() const
260  { return this->_fe_map->get_dxyzdeta(); }
261 
266  const std::vector<RealGradient> & get_dxyzdzeta() const
267  { return _fe_map->get_dxyzdzeta(); }
268 
272  const std::vector<RealGradient> & get_d2xyzdxi2() const
273  { return this->_fe_map->get_d2xyzdxi2(); }
274 
278  const std::vector<RealGradient> & get_d2xyzdeta2() const
279  { return this->_fe_map->get_d2xyzdeta2(); }
280 
281 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
282 
286  const std::vector<RealGradient> & get_d2xyzdzeta2() const
287  { return this->_fe_map->get_d2xyzdzeta2(); }
288 
289 #endif
290 
294  const std::vector<RealGradient> & get_d2xyzdxideta() const
295  { return this->_fe_map->get_d2xyzdxideta(); }
296 
297 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
298 
302  const std::vector<RealGradient> & get_d2xyzdxidzeta() const
303  { return this->_fe_map->get_d2xyzdxidzeta(); }
304 
308  const std::vector<RealGradient> & get_d2xyzdetadzeta() const
309  { return this->_fe_map->get_d2xyzdetadzeta(); }
310 
311 #endif
312 
317  const std::vector<Real> & get_dxidx() const
318  { return this->_fe_map->get_dxidx(); }
319 
324  const std::vector<Real> & get_dxidy() const
325  { return this->_fe_map->get_dxidy(); }
326 
331  const std::vector<Real> & get_dxidz() const
332  { return this->_fe_map->get_dxidz(); }
333 
338  const std::vector<Real> & get_detadx() const
339  { return this->_fe_map->get_detadx(); }
340 
345  const std::vector<Real> & get_detady() const
346  { return this->_fe_map->get_detady(); }
347 
352  const std::vector<Real> & get_detadz() const
353  { return this->_fe_map->get_detadz(); }
354 
359  const std::vector<Real> & get_dzetadx() const
360  { return this->_fe_map->get_dzetadx(); }
361 
366  const std::vector<Real> & get_dzetady() const
367  { return this->_fe_map->get_dzetady(); }
368 
373  const std::vector<Real> & get_dzetadz() const
374  { return this->_fe_map->get_dzetadz(); }
375 
379  const std::vector<std::vector<Point>> & get_tangents() const
380  { return this->_fe_map->get_tangents(); }
381 
385  const std::vector<Point> & get_normals() const
386  { return this->_fe_map->get_normals(); }
387 
391  const std::vector<Real> & get_curvatures() const
392  { return this->_fe_map->get_curvatures();}
393 
398  virtual void attach_quadrature_rule (QBase * q) = 0;
399 
405  virtual unsigned int n_shape_functions () const = 0;
406 
411  virtual unsigned int n_quadrature_points () const = 0;
412 
418  ElemType get_type() const { return elem_type; }
419 
424  unsigned int get_p_level() const { return _p_level; }
425 
429  FEType get_fe_type() const { return fe_type; }
430 
434  Order get_order() const { return static_cast<Order>(fe_type.order + _p_level); }
435 
439  void set_fe_order(int new_order) { fe_type.order = new_order; }
440 
444  virtual FEContinuity get_continuity() const = 0;
445 
450  virtual bool is_hierarchic() const = 0;
451 
455  FEFamily get_family() const { return fe_type.family; }
456 
460  const FEMap & get_fe_map() const { return *_fe_map.get(); }
461 
465  void print_JxW(std::ostream & os) const;
466 
472  virtual void print_phi(std::ostream & os) const =0;
473 
479  virtual void print_dphi(std::ostream & os) const =0;
480 
481 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
482 
488  virtual void print_d2phi(std::ostream & os) const =0;
489 
490 #endif
491 
496  void print_xyz(std::ostream & os) const;
497 
501  void print_info(std::ostream & os) const;
502 
506  friend std::ostream & operator << (std::ostream & os, const FEAbstract & fe);
507 
508 
509 protected:
510 
523  virtual void compute_shape_functions(const Elem *, const std::vector<Point> & ) =0;
524 
525  std::unique_ptr<FEMap> _fe_map;
526 
527 
531  const unsigned int dim;
532 
537  mutable bool calculations_started;
538 
542  mutable bool calculate_phi;
543 
547  mutable bool calculate_dphi;
548 
552  mutable bool calculate_d2phi;
553 
557  mutable bool calculate_curl_phi;
558 
562  mutable bool calculate_div_phi;
563 
567  mutable bool calculate_dphiref;
568 
569 
576 
582 
587  unsigned int _p_level;
588 
593 
599 
606  virtual bool shapes_need_reinit() const = 0;
607 
608 };
609 
610 } // namespace libMesh
611 
612 #endif // LIBMESH_FE_ABSTRACT_H
Order get_order() const
Definition: fe_abstract.h:434
const std::vector< Real > & get_curvatures() const
Definition: fe_abstract.h:391
Manages the family, order, etc. parameters for a given FE.
Definition: fe_type.h:179
FEFamily family
Definition: fe_type.h:204
const std::vector< Real > & get_detady() const
Definition: fe_abstract.h:345
virtual bool is_hierarchic() const =0
virtual void print_phi(std::ostream &os) const =0
unsigned int get_p_level() const
Definition: fe_abstract.h:424
unsigned int _p_level
Definition: fe_abstract.h:587
const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:238
const std::vector< Real > & get_dzetadz() const
Definition: fe_abstract.h:373
const std::vector< Real > & get_detadx() const
Definition: fe_abstract.h:338
static void compute_node_constraints(NodeConstraints &constraints, const Elem *elem)
Definition: fe_abstract.C:820
unsigned int get_dim() const
Definition: fe_abstract.h:231
ElemType get_type() const
Definition: fe_abstract.h:418
const std::vector< std::vector< Point > > & get_tangents() const
Definition: fe_abstract.h:379
Maps between boundary ids and PeriodicBoundaryBase objects.
const std::vector< RealGradient > & get_d2xyzdetadzeta() const
Definition: fe_abstract.h:308
const std::vector< RealGradient > & get_dxyzdzeta() const
Definition: fe_abstract.h:266
const std::vector< Real > & get_dxidy() const
Definition: fe_abstract.h:324
const std::vector< RealGradient > & get_d2xyzdeta2() const
Definition: fe_abstract.h:278
unsigned short int side
Definition: xdr_io.C:50
const std::vector< RealGradient > & get_dxyzdxi() const
Definition: fe_abstract.h:252
The base class for all geometric element types.
Definition: elem.h:100
MeshBase & mesh
virtual unsigned int n_quadrature_points() const =0
FEAbstract(const unsigned int dim, const FEType &fet)
Definition: fe_abstract.C:46
OrderWrapper order
Definition: fe_type.h:198
static const Real TOLERANCE
static bool on_reference_element(const Point &p, const ElemType t, const Real eps=TOLERANCE)
Definition: fe_abstract.C:581
Base class for Mesh.
Definition: mesh_base.h:77
virtual FEContinuity get_continuity() const =0
FEType get_fe_type() const
Definition: fe_abstract.h:429
virtual unsigned int n_shape_functions() const =0
virtual void print_dphi(std::ostream &os) const =0
const std::vector< Point > & get_normals() const
Definition: fe_abstract.h:385
const std::vector< Real > & get_dzetady() const
Definition: fe_abstract.h:366
const unsigned int dim
Definition: fe_abstract.h:531
void print_xyz(std::ostream &os) const
Definition: fe_abstract.C:787
const std::vector< RealGradient > & get_d2xyzdzeta2() const
Definition: fe_abstract.h:286
const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:245
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
Definition: fe_abstract.C:283
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr)=0
virtual bool shapes_need_reinit() const =0
const std::vector< RealGradient > & get_d2xyzdxideta() const
Definition: fe_abstract.h:294
const std::vector< Real > & get_dxidz() const
Definition: fe_abstract.h:331
void print_JxW(std::ostream &os) const
Definition: fe_abstract.C:780
void print_info(std::ostream &os) const
Definition: fe_abstract.C:793
virtual void edge_reinit(const Elem *elem, const unsigned int edge, const Real tolerance=TOLERANCE, const std::vector< Point > *pts=nullptr, const std::vector< Real > *weights=nullptr)=0
virtual void compute_shape_functions(const Elem *, const std::vector< Point > &)=0
const std::vector< RealGradient > & get_dxyzdeta() const
Definition: fe_abstract.h:259
const std::vector< Real > & get_detadz() const
Definition: fe_abstract.h:352
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static std::unique_ptr< FEAbstract > build(const unsigned int dim, const FEType &type)
Definition: fe_abstract.C:71
const std::vector< Real > & get_dzetadx() const
Definition: fe_abstract.h:359
FEFamily get_family() const
Definition: fe_abstract.h:455
friend std::ostream & operator<<(std::ostream &os, const FEAbstract &fe)
Definition: fe_abstract.C:809
virtual ~FEAbstract()
Definition: fe_abstract.C:66
void set_fe_order(int new_order)
Definition: fe_abstract.h:439
const std::vector< Real > & get_dxidx() const
Definition: fe_abstract.h:317
virtual void print_d2phi(std::ostream &os) const =0
const FEMap & get_fe_map() const
Definition: fe_abstract.h:460
std::unique_ptr< FEMap > _fe_map
Definition: fe_abstract.h:525
const std::vector< RealGradient > & get_d2xyzdxidzeta() const
Definition: fe_abstract.h:302
Computes finite element mapping function values, gradients, etc.
Definition: fe_map.h:48
A geometric point in (x,y,z) space.
Definition: point.h:38
virtual void attach_quadrature_rule(QBase *q)=0
Base class for all quadrature families and orders.
Definition: quadrature.h:62
const std::vector< RealGradient > & get_d2xyzdxi2() const
Definition: fe_abstract.h:272
static void compute_periodic_node_constraints(NodeConstraints &constraints, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const Elem *elem)
Definition: fe_abstract.C:965
virtual void side_map(const Elem *elem, const Elem *side, const unsigned int s, const std::vector< Point > &reference_side_points, std::vector< Point > &reference_points)=0