libMesh
fe_base.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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_BASE_H
21 #define LIBMESH_FE_BASE_H
22 
23 // Local includes
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/compare_types.h"
26 #include "libmesh/fe_abstract.h"
27 #include "libmesh/fe_transformation_base.h"
28 #include "libmesh/point.h"
29 #include "libmesh/reference_counted_object.h"
30 #include "libmesh/tensor_tools.h"
31 #include "libmesh/type_n_tensor.h"
32 #include "libmesh/vector_value.h"
33 
34 // C++ includes
35 #include <cstddef>
36 #include <vector>
37 #include <memory>
38 
39 namespace libMesh
40 {
41 
42 
43 // forward declarations
44 template <typename T> class DenseMatrix;
45 template <typename T> class DenseVector;
46 class BoundaryInfo;
47 class DofConstraints;
48 class DofMap;
49 class Elem;
50 class MeshBase;
51 template <typename T> class NumericVector;
52 class QBase;
53 template <typename T> class FETransformationBase;
54 class FEType;
55 
56 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
57 class NodeConstraints;
58 #endif
59 
60 #ifdef LIBMESH_ENABLE_PERIODIC
61 class PeriodicBoundaries;
62 class PointLocatorBase;
63 #endif
64 
65 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
66 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
67 class InfFE;
68 #endif
69 
83 template <typename OutputType>
84 class FEGenericBase : public FEAbstract
85 {
86 protected:
87 
93  FEGenericBase (const unsigned int dim,
94  const FEType & fet);
95 
96 public:
97 
101  virtual ~FEGenericBase();
102 
111  static std::unique_ptr<FEGenericBase> build (const unsigned int dim,
112  const FEType & type);
113 
118  typedef OutputType OutputShape;
126 
127 
128 
129 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
130 
139  static std::unique_ptr<FEGenericBase> build_InfFE (const unsigned int dim,
140  const FEType & type);
141 
142 #endif
143 
144 #ifdef LIBMESH_ENABLE_AMR
145 
152  static void compute_proj_constraints (DofConstraints & constraints,
153  DofMap & dof_map,
154  const unsigned int variable_number,
155  const Elem * elem);
156 
164  static void coarsened_dof_values(const NumericVector<Number> & global_vector,
165  const DofMap & dof_map,
166  const Elem * coarse_elem,
167  DenseVector<Number> & coarse_dofs,
168  const unsigned int var,
169  const bool use_old_dof_indices = false);
170 
177  static void coarsened_dof_values(const NumericVector<Number> & global_vector,
178  const DofMap & dof_map,
179  const Elem * coarse_elem,
180  DenseVector<Number> & coarse_dofs,
181  const bool use_old_dof_indices = false);
182 
183 #endif // #ifdef LIBMESH_ENABLE_AMR
184 
185 #ifdef LIBMESH_ENABLE_PERIODIC
186 
192  static void compute_periodic_constraints (DofConstraints & constraints,
193  DofMap & dof_map,
194  const PeriodicBoundaries & boundaries,
195  const MeshBase & mesh,
196  const PointLocatorBase * point_locator,
197  const unsigned int variable_number,
198  const Elem * elem);
199 
200 #endif // LIBMESH_ENABLE_PERIODIC
201 
206  const std::vector<std::vector<OutputShape>> & get_phi() const
208  calculate_phi = true; return phi; }
209 
214  const std::vector<std::vector<OutputGradient>> & get_dphi() const
216  calculate_dphi = calculate_dphiref = true; return dphi; }
217 
222  const std::vector<std::vector<OutputShape>> & get_curl_phi() const
224  calculate_curl_phi = calculate_dphiref = true; return curl_phi; }
225 
230  const std::vector<std::vector<OutputDivergence>> & get_div_phi() const
232  calculate_div_phi = calculate_dphiref = true; return div_phi; }
233 
238  const std::vector<std::vector<OutputShape>> & get_dphidx() const
240  calculate_dphi = calculate_dphiref = true; return dphidx; }
241 
246  const std::vector<std::vector<OutputShape>> & get_dphidy() const
248  calculate_dphi = calculate_dphiref = true; return dphidy; }
249 
254  const std::vector<std::vector<OutputShape>> & get_dphidz() const
256  calculate_dphi = calculate_dphiref = true; return dphidz; }
257 
262  const std::vector<std::vector<OutputShape>> & get_dphidxi() const
264  calculate_dphiref = true; return dphidxi; }
265 
270  const std::vector<std::vector<OutputShape>> & get_dphideta() const
272  calculate_dphiref = true; return dphideta; }
273 
278  const std::vector<std::vector<OutputShape>> & get_dphidzeta() const
280  calculate_dphiref = true; return dphidzeta; }
281 
282 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
283 
288  const std::vector<std::vector<OutputTensor>> & get_d2phi() const
290  calculate_d2phi = calculate_dphiref = true; return d2phi; }
291 
296  const std::vector<std::vector<OutputShape>> & get_d2phidx2() const
298  calculate_d2phi = calculate_dphiref = true; return d2phidx2; }
299 
304  const std::vector<std::vector<OutputShape>> & get_d2phidxdy() const
306  calculate_d2phi = calculate_dphiref = true; return d2phidxdy; }
307 
312  const std::vector<std::vector<OutputShape>> & get_d2phidxdz() const
314  calculate_d2phi = calculate_dphiref = true; return d2phidxdz; }
315 
320  const std::vector<std::vector<OutputShape>> & get_d2phidy2() const
322  calculate_d2phi = calculate_dphiref = true; return d2phidy2; }
323 
328  const std::vector<std::vector<OutputShape>> & get_d2phidydz() const
330  calculate_d2phi = calculate_dphiref = true; return d2phidydz; }
331 
336  const std::vector<std::vector<OutputShape>> & get_d2phidz2() const
338  calculate_d2phi = calculate_dphiref = true; return d2phidz2; }
339 
344  const std::vector<std::vector<OutputShape>> & get_d2phidxi2() const
346  calculate_d2phi = calculate_dphiref = true; return d2phidxi2; }
347 
352  const std::vector<std::vector<OutputShape>> & get_d2phidxideta() const
355 
360  const std::vector<std::vector<OutputShape>> & get_d2phidxidzeta() const
363 
368  const std::vector<std::vector<OutputShape>> & get_d2phideta2() const
370  calculate_d2phi = calculate_dphiref = true; return d2phideta2; }
371 
376  const std::vector<std::vector<OutputShape>> & get_d2phidetadzeta() const
379 
384  const std::vector<std::vector<OutputShape>> & get_d2phidzeta2() const
386  calculate_d2phi = calculate_dphiref = true; return d2phidzeta2; }
387 
388 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
389 
390 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
391 
402  const std::vector<OutputGradient> & get_dphase() const
403  { return dphase; }
404 
405 
418  const std::vector<Real> & get_Sobolev_weight() const
419  { return weight; }
420 
426  const std::vector<RealGradient> & get_Sobolev_dweight() const
427  { return dweight; }
428 
429 #endif
430 
431 
435  void print_phi(std::ostream & os) const;
436 
441  void print_dphi(std::ostream & os) const;
442 
443 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
444 
449  void print_d2phi(std::ostream & os) const;
450 
451 #endif
452 
453 
454 protected:
455 
456 
457 
458 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
459 
465  virtual void init_base_shape_functions(const std::vector<Point> & qp,
466  const Elem * e) = 0;
467 
468 #endif
469 
474  void determine_calculations();
475 
486  virtual void compute_shape_functions(const Elem * elem, const std::vector<Point> & qp);
487 
492  std::unique_ptr<FETransformationBase<OutputType>> _fe_trans;
493 
497  std::vector<std::vector<OutputShape>> phi;
498 
502  std::vector<std::vector<OutputGradient>> dphi;
503 
507  std::vector<std::vector<OutputShape>> curl_phi;
508 
512  std::vector<std::vector<OutputDivergence>> div_phi;
513 
517  std::vector<std::vector<OutputShape>> dphidxi;
518 
522  std::vector<std::vector<OutputShape>> dphideta;
523 
527  std::vector<std::vector<OutputShape>> dphidzeta;
528 
532  std::vector<std::vector<OutputShape>> dphidx;
533 
537  std::vector<std::vector<OutputShape>> dphidy;
538 
542  std::vector<std::vector<OutputShape>> dphidz;
543 
544 
545 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
546 
550  std::vector<std::vector<OutputTensor>> d2phi;
551 
555  std::vector<std::vector<OutputShape>> d2phidxi2;
556 
560  std::vector<std::vector<OutputShape>> d2phidxideta;
561 
565  std::vector<std::vector<OutputShape>> d2phidxidzeta;
566 
570  std::vector<std::vector<OutputShape>> d2phideta2;
571 
575  std::vector<std::vector<OutputShape>> d2phidetadzeta;
576 
580  std::vector<std::vector<OutputShape>> d2phidzeta2;
581 
585  std::vector<std::vector<OutputShape>> d2phidx2;
586 
590  std::vector<std::vector<OutputShape>> d2phidxdy;
591 
595  std::vector<std::vector<OutputShape>> d2phidxdz;
596 
600  std::vector<std::vector<OutputShape>> d2phidy2;
601 
605  std::vector<std::vector<OutputShape>> d2phidydz;
606 
610  std::vector<std::vector<OutputShape>> d2phidz2;
611 
612 #endif
613 
614 
615 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
616 
617  //--------------------------------------------------------------
618  /* protected members for infinite elements, which are accessed
619  * from the outside through some inline functions
620  */
621 
622 
628  std::vector<OutputGradient> dphase;
629 
635  std::vector<RealGradient> dweight;
636 
642  std::vector<Real> weight;
643 
644 #endif
645 
646 private:
647 
648 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
649 
655  template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
656  friend class InfFE;
657 
658 #endif
659 
660 
661 };
662 
663 
664 // Typedefs for convenience and backwards compatibility
667 
668 
669 
670 
671 // ------------------------------------------------------------
672 // FEGenericBase class inline members
673 template <typename OutputType>
674 inline
676  const FEType & fet) :
677  FEAbstract(d,fet),
678  _fe_trans( FETransformationBase<OutputType>::build(fet) ),
679  phi(),
680  dphi(),
681  curl_phi(),
682  div_phi(),
683  dphidxi(),
684  dphideta(),
685  dphidzeta(),
686  dphidx(),
687  dphidy(),
688  dphidz()
689 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
690  ,d2phi(),
691  d2phidxi2(),
692  d2phidxideta(),
693  d2phidxidzeta(),
694  d2phideta2(),
695  d2phidetadzeta(),
696  d2phidzeta2(),
697  d2phidx2(),
698  d2phidxdy(),
699  d2phidxdz(),
700  d2phidy2(),
701  d2phidydz(),
702  d2phidz2()
703 #endif
704 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
705  ,dphase(),
706  dweight(),
707  weight()
708 #endif
709 {
710 }
711 
712 
713 
714 template <typename OutputType>
715 inline
717 {
718 }
719 
720 } // namespace libMesh
721 
722 #endif // LIBMESH_FE_BASE_H
libMesh::FEGenericBase::print_d2phi
void print_d2phi(std::ostream &os) const
Prints the value of each shape function's second derivatives at each quadrature point.
Definition: fe_base.C:807
libMesh::DofConstraints
The constraint matrix storage format.
Definition: dof_map.h:105
libMesh::FEGenericBase::d2phidxdz
std::vector< std::vector< OutputShape > > d2phidxdz
Shape function second derivatives in the x-z direction.
Definition: fe_base.h:595
libMesh::FEGenericBase::get_d2phideta2
const std::vector< std::vector< OutputShape > > & get_d2phideta2() const
Definition: fe_base.h:368
libMesh::PeriodicBoundaries
We're using a class instead of a typedef to allow forward declarations and future flexibility.
Definition: periodic_boundaries.h:51
libMesh::FEGenericBase::get_d2phidzeta2
const std::vector< std::vector< OutputShape > > & get_d2phidzeta2() const
Definition: fe_base.h:384
libMesh::FEGenericBase::compute_periodic_constraints
static void compute_periodic_constraints(DofConstraints &constraints, DofMap &dof_map, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for meshes with periodic boundary conditions) correspon...
Definition: fe_base.C:1692
libMesh::FEGenericBase::_fe_trans
std::unique_ptr< FETransformationBase< OutputType > > _fe_trans
Object that handles computing shape function values, gradients, etc in the physical domain.
Definition: fe_base.h:492
libMesh::FEAbstract::calculate_dphi
bool calculate_dphi
Should we calculate shape function gradients?
Definition: fe_abstract.h:551
libMesh::FEGenericBase::get_d2phidxidzeta
const std::vector< std::vector< OutputShape > > & get_d2phidxidzeta() const
Definition: fe_base.h:360
libMesh::FEGenericBase::d2phideta2
std::vector< std::vector< OutputShape > > d2phideta2
Shape function second derivatives in the eta direction.
Definition: fe_base.h:570
libMesh::FEGenericBase::curl_phi
std::vector< std::vector< OutputShape > > curl_phi
Shape function curl values.
Definition: fe_base.h:507
libMesh::FEGenericBase::FEGenericBase
FEGenericBase(const unsigned int dim, const FEType &fet)
Constructor.
Definition: fe_base.h:675
libMesh::FEGenericBase::OutputShape
OutputType OutputShape
Convenient typedefs for gradients of output, hessians of output, and potentially-complex-valued versi...
Definition: fe_base.h:118
libMesh::FEGenericBase::~FEGenericBase
virtual ~FEGenericBase()
Destructor.
libMesh::FEGenericBase::dphidzeta
std::vector< std::vector< OutputShape > > dphidzeta
Shape function derivatives in the zeta direction.
Definition: fe_base.h:527
libMesh::FEGenericBase::dphideta
std::vector< std::vector< OutputShape > > dphideta
Shape function derivatives in the eta direction.
Definition: fe_base.h:522
libMesh::FEGenericBase::get_dphidxi
const std::vector< std::vector< OutputShape > > & get_dphidxi() const
Definition: fe_base.h:262
libMesh::FEGenericBase::get_dphase
const std::vector< OutputGradient > & get_dphase() const
Definition: fe_base.h:402
libMesh::FEGenericBase::dphidz
std::vector< std::vector< OutputShape > > dphidz
Shape function derivatives in the z direction.
Definition: fe_base.h:542
libMesh::InfFE
A specific instantiation of the FEBase class.
Definition: fe.h:40
libMesh::FEGenericBase::get_dphidx
const std::vector< std::vector< OutputShape > > & get_dphidx() const
Definition: fe_base.h:238
libMesh::FEGenericBase::div_phi
std::vector< std::vector< OutputDivergence > > div_phi
Shape function divergence values.
Definition: fe_base.h:512
libMesh::FEGenericBase::dphidy
std::vector< std::vector< OutputShape > > dphidy
Shape function derivatives in the y direction.
Definition: fe_base.h:537
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::FEAbstract
This class forms the foundation from which generic finite elements may be derived.
Definition: fe_abstract.h:100
libMesh::FEGenericBase
This class forms the foundation from which generic finite elements may be derived.
Definition: exact_error_estimator.h:39
libMesh::FEAbstract::calculate_dphiref
bool calculate_dphiref
Should we calculate reference shape function gradients?
Definition: fe_abstract.h:577
libMesh::FEAbstract::calculations_started
bool calculations_started
Have calculations with this object already been started? Then all get_* functions should already have...
Definition: fe_abstract.h:536
libMesh::FEGenericBase::get_d2phi
const std::vector< std::vector< OutputTensor > > & get_d2phi() const
Definition: fe_base.h:288
libMesh::FEAbstract::calculate_phi
bool calculate_phi
Should we calculate shape functions?
Definition: fe_abstract.h:546
libMesh::FEVectorBase
FEGenericBase< RealGradient > FEVectorBase
Definition: fe_base.h:666
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::FEGenericBase::d2phidxidzeta
std::vector< std::vector< OutputShape > > d2phidxidzeta
Shape function second derivatives in the xi-zeta direction.
Definition: fe_base.h:565
libMesh::FEGenericBase::get_d2phidy2
const std::vector< std::vector< OutputShape > > & get_d2phidy2() const
Definition: fe_base.h:320
libMesh::FEGenericBase::get_dphi
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:214
libMesh::FEAbstract::calculate_curl_phi
bool calculate_curl_phi
Should we calculate shape function curls?
Definition: fe_abstract.h:567
libMesh::FEGenericBase::get_div_phi
const std::vector< std::vector< OutputDivergence > > & get_div_phi() const
Definition: fe_base.h:230
libMesh::TensorTools::DecrementRank::type
T type
Definition: tensor_tools.h:158
libMesh::FEGenericBase::phi
std::vector< std::vector< OutputShape > > phi
Shape function values.
Definition: fe_base.h:497
libMesh::FEGenericBase::build_InfFE
static std::unique_ptr< FEGenericBase > build_InfFE(const unsigned int dim, const FEType &type)
Builds a specific infinite element type.
libMesh::FEGenericBase::d2phidetadzeta
std::vector< std::vector< OutputShape > > d2phidetadzeta
Shape function second derivatives in the eta-zeta direction.
Definition: fe_base.h:575
libMesh::FEGenericBase::d2phidz2
std::vector< std::vector< OutputShape > > d2phidz2
Shape function second derivatives in the z direction.
Definition: fe_base.h:610
libMesh::VectorValue
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
Definition: exact_solution.h:54
libMesh::FEGenericBase::get_dphideta
const std::vector< std::vector< OutputShape > > & get_dphideta() const
Definition: fe_base.h:270
libMesh::NumericVector< Number >
libMesh::FEGenericBase::dphidxi
std::vector< std::vector< OutputShape > > dphidxi
Shape function derivatives in the xi direction.
Definition: fe_base.h:517
libMesh::FEGenericBase::d2phidzeta2
std::vector< std::vector< OutputShape > > d2phidzeta2
Shape function second derivatives in the zeta direction.
Definition: fe_base.h:580
libMesh::FEGenericBase::dweight
std::vector< RealGradient > dweight
Used for certain infinite element families: the global derivative of the additional radial weight ,...
Definition: fe_base.h:635
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::FEGenericBase::get_d2phidxideta
const std::vector< std::vector< OutputShape > > & get_d2phidxideta() const
Definition: fe_base.h:352
libMesh::FEGenericBase::d2phidxideta
std::vector< std::vector< OutputShape > > d2phidxideta
Shape function second derivatives in the xi-eta direction.
Definition: fe_base.h:560
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::FEBase
FEGenericBase< Real > FEBase
Definition: exact_error_estimator.h:39
libMesh::FEGenericBase::d2phi
std::vector< std::vector< OutputTensor > > d2phi
Shape function second derivative values.
Definition: fe_base.h:550
libMesh::FEGenericBase::d2phidydz
std::vector< std::vector< OutputShape > > d2phidydz
Shape function second derivatives in the y-z direction.
Definition: fe_base.h:605
libMesh::FEAbstract::calculate_div_phi
bool calculate_div_phi
Should we calculate shape function divergences?
Definition: fe_abstract.h:572
libMesh::FEGenericBase::build
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libMesh::FEGenericBase::get_d2phidxdz
const std::vector< std::vector< OutputShape > > & get_d2phidxdz() const
Definition: fe_base.h:312
libMesh::FEGenericBase::OutputDivergence
TensorTools::DecrementRank< OutputShape >::type OutputDivergence
Definition: fe_base.h:121
libMesh::FEGenericBase::dphi
std::vector< std::vector< OutputGradient > > dphi
Shape function derivative values.
Definition: fe_base.h:502
libMesh::FEGenericBase::d2phidxi2
std::vector< std::vector< OutputShape > > d2phidxi2
Shape function second derivatives in the xi direction.
Definition: fe_base.h:555
libMesh::FEGenericBase::print_dphi
void print_dphi(std::ostream &os) const
Prints the value of each shape function's derivative at each quadrature point.
Definition: fe_base.C:746
libMesh::FEAbstract::dim
const unsigned int dim
The dimensionality of the object.
Definition: fe_abstract.h:530
libMesh::FEGenericBase::coarsened_dof_values
static void coarsened_dof_values(const NumericVector< Number > &global_vector, const DofMap &dof_map, const Elem *coarse_elem, DenseVector< Number > &coarse_dofs, const unsigned int var, const bool use_old_dof_indices=false)
Creates a local projection on coarse_elem, based on the DoF values in global_vector for it's children...
Definition: fe_base.C:822
libMesh::FEGenericBase::OutputNumber
TensorTools::MakeNumber< OutputShape >::type OutputNumber
Definition: fe_base.h:122
libMesh::FEGenericBase::get_d2phidxi2
const std::vector< std::vector< OutputShape > > & get_d2phidxi2() const
Definition: fe_base.h:344
libMesh::FEGenericBase::get_Sobolev_dweight
const std::vector< RealGradient > & get_Sobolev_dweight() const
Definition: fe_base.h:426
libMesh::FEGenericBase::compute_shape_functions
virtual void compute_shape_functions(const Elem *elem, const std::vector< Point > &qp)
After having updated the jacobian and the transformation from local to global coordinates in FEAbstra...
Definition: fe_base.C:697
libMesh::FEGenericBase::OutputNumberDivergence
TensorTools::DecrementRank< OutputNumber >::type OutputNumberDivergence
Definition: fe_base.h:125
libMesh::FEGenericBase::get_d2phidydz
const std::vector< std::vector< OutputShape > > & get_d2phidydz() const
Definition: fe_base.h:328
libMesh::NodeConstraints
The Node constraint storage format.
Definition: dof_map.h:153
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::FEGenericBase::get_dphidy
const std::vector< std::vector< OutputShape > > & get_dphidy() const
Definition: fe_base.h:246
libMesh::FEGenericBase::OutputGradient
TensorTools::IncrementRank< OutputShape >::type OutputGradient
Definition: fe_base.h:119
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::FEGenericBase::get_curl_phi
const std::vector< std::vector< OutputShape > > & get_curl_phi() const
Definition: fe_base.h:222
libMesh::FEGenericBase::compute_proj_constraints
static void compute_proj_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
Definition: fe_base.C:1396
libMesh::FEGenericBase::print_phi
void print_phi(std::ostream &os) const
Prints the value of each shape function at each quadrature point.
Definition: fe_base.C:735
libMesh::FEGenericBase::dphase
std::vector< OutputGradient > dphase
Used for certain infinite element families: the first derivatives of the phase term in global coordin...
Definition: fe_base.h:628
libMesh::FEGenericBase::get_d2phidxdy
const std::vector< std::vector< OutputShape > > & get_d2phidxdy() const
Definition: fe_base.h:304
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::FEGenericBase::dphidx
std::vector< std::vector< OutputShape > > dphidx
Shape function derivatives in the x direction.
Definition: fe_base.h:532
libMesh::FEGenericBase::d2phidxdy
std::vector< std::vector< OutputShape > > d2phidxdy
Shape function second derivatives in the x-y direction.
Definition: fe_base.h:590
libMesh::FEGenericBase::determine_calculations
void determine_calculations()
Determine which values are to be calculated, for both the FE itself and for the FEMap.
Definition: fe_base.C:756
libMesh::FEGenericBase::get_dphidz
const std::vector< std::vector< OutputShape > > & get_dphidz() const
Definition: fe_base.h:254
libMesh::FEGenericBase::OutputNumberGradient
TensorTools::IncrementRank< OutputNumber >::type OutputNumberGradient
Definition: fe_base.h:123
libMesh::FETransformationBase
This class handles the computation of the shape functions in the physical domain.
Definition: fe_base.h:53
libMesh::FEGenericBase::d2phidy2
std::vector< std::vector< OutputShape > > d2phidy2
Shape function second derivatives in the y direction.
Definition: fe_base.h:600
libMesh::FEGenericBase::get_dphidzeta
const std::vector< std::vector< OutputShape > > & get_dphidzeta() const
Definition: fe_base.h:278
libMesh::FEGenericBase::get_phi
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:206
libMesh::FEGenericBase::get_d2phidx2
const std::vector< std::vector< OutputShape > > & get_d2phidx2() const
Definition: fe_base.h:296
libMesh::FEGenericBase::OutputNumberTensor
TensorTools::IncrementRank< OutputNumberGradient >::type OutputNumberTensor
Definition: fe_base.h:124
libMesh::FEGenericBase::get_d2phidetadzeta
const std::vector< std::vector< OutputShape > > & get_d2phidetadzeta() const
Definition: fe_base.h:376
libMesh::MeshTools::weight
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:236
libMesh::PointLocatorBase
This is the base class for point locators.
Definition: point_locator_base.h:62
libMesh::FEGenericBase::d2phidx2
std::vector< std::vector< OutputShape > > d2phidx2
Shape function second derivatives in the x direction.
Definition: fe_base.h:585
libMesh::FEGenericBase::init_base_shape_functions
virtual void init_base_shape_functions(const std::vector< Point > &qp, const Elem *e)=0
Initialize the data fields for the base of an an infinite element.
libMesh::FEAbstract::calculate_d2phi
bool calculate_d2phi
Should we calculate shape function hessians?
Definition: fe_abstract.h:557
libMesh::FEGenericBase::weight
std::vector< Real > weight
Used for certain infinite element families: the additional radial weight in local coordinates,...
Definition: fe_base.h:642
libMesh::FEGenericBase::OutputTensor
TensorTools::IncrementRank< OutputGradient >::type OutputTensor
Definition: fe_base.h:120
libMesh::DenseVector< Number >
libMesh::FEGenericBase::get_d2phidz2
const std::vector< std::vector< OutputShape > > & get_d2phidz2() const
Definition: fe_base.h:336
libMesh::TensorTools::MakeNumber::type
std::complex< T > type
Definition: tensor_tools.h:196
libMesh::FEGenericBase::get_Sobolev_weight
const std::vector< Real > & get_Sobolev_weight() const
Definition: fe_base.h:418