fe_base.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_BASE_H
21 #define LIBMESH_FE_BASE_H
22 
23 // Local includes
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/auto_ptr.h" // deprecated
26 #include "libmesh/compare_types.h"
27 #include "libmesh/fe_abstract.h"
29 #include "libmesh/point.h"
31 #include "libmesh/tensor_tools.h"
32 #include "libmesh/type_n_tensor.h"
33 #include "libmesh/vector_value.h"
34 
35 // C++ includes
36 #include <cstddef>
37 #include <vector>
38 #include <memory>
39 
40 namespace libMesh
41 {
42 
43 
44 // forward declarations
45 template <typename T> class DenseMatrix;
46 template <typename T> class DenseVector;
47 class BoundaryInfo;
48 class DofConstraints;
49 class DofMap;
50 class Elem;
51 class MeshBase;
52 template <typename T> class NumericVector;
53 class QBase;
54 template <typename T> class FETransformationBase;
55 class FEType;
56 
57 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
58 class NodeConstraints;
59 #endif
60 
61 #ifdef LIBMESH_ENABLE_PERIODIC
62 class PeriodicBoundaries;
63 class PointLocatorBase;
64 #endif
65 
66 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
67 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
68 class InfFE;
69 #endif
70 
84 template <typename OutputType>
85 class FEGenericBase : public FEAbstract
86 {
87 protected:
88 
94  FEGenericBase (const unsigned int dim,
95  const FEType & fet);
96 
97 public:
98 
102  virtual ~FEGenericBase();
103 
112  static std::unique_ptr<FEGenericBase> build (const unsigned int dim,
113  const FEType & type);
114 
119  typedef OutputType OutputShape;
127 
128 
129 
130 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
131 
140  static std::unique_ptr<FEGenericBase> build_InfFE (const unsigned int dim,
141  const FEType & type);
142 
143 #endif
144 
145 #ifdef LIBMESH_ENABLE_AMR
146 
153  static void compute_proj_constraints (DofConstraints & constraints,
154  DofMap & dof_map,
155  const unsigned int variable_number,
156  const Elem * elem);
157 
165  static void coarsened_dof_values(const NumericVector<Number> & global_vector,
166  const DofMap & dof_map,
167  const Elem * coarse_elem,
168  DenseVector<Number> & coarse_dofs,
169  const unsigned int var,
170  const bool use_old_dof_indices = false);
171 
178  static void coarsened_dof_values(const NumericVector<Number> & global_vector,
179  const DofMap & dof_map,
180  const Elem * coarse_elem,
181  DenseVector<Number> & coarse_dofs,
182  const bool use_old_dof_indices = false);
183 
184 #endif // #ifdef LIBMESH_ENABLE_AMR
185 
186 #ifdef LIBMESH_ENABLE_PERIODIC
187 
193  static void compute_periodic_constraints (DofConstraints & constraints,
194  DofMap & dof_map,
195  const PeriodicBoundaries & boundaries,
196  const MeshBase & mesh,
197  const PointLocatorBase * point_locator,
198  const unsigned int variable_number,
199  const Elem * elem);
200 
201 #endif // LIBMESH_ENABLE_PERIODIC
202 
207  const std::vector<std::vector<OutputShape>> & get_phi() const
208  { libmesh_assert(!calculations_started || calculate_phi);
209  calculate_phi = true; return phi; }
210 
215  const std::vector<std::vector<OutputGradient>> & get_dphi() const
216  { libmesh_assert(!calculations_started || calculate_dphi);
217  calculate_dphi = calculate_dphiref = true; return dphi; }
218 
223  const std::vector<std::vector<OutputShape>> & get_curl_phi() const
224  { libmesh_assert(!calculations_started || calculate_curl_phi);
225  calculate_curl_phi = calculate_dphiref = true; return curl_phi; }
226 
231  const std::vector<std::vector<OutputDivergence>> & get_div_phi() const
232  { libmesh_assert(!calculations_started || calculate_div_phi);
233  calculate_div_phi = calculate_dphiref = true; return div_phi; }
234 
239  const std::vector<std::vector<OutputShape>> & get_dphidx() const
240  { libmesh_assert(!calculations_started || calculate_dphi);
241  calculate_dphi = calculate_dphiref = true; return dphidx; }
242 
247  const std::vector<std::vector<OutputShape>> & get_dphidy() const
248  { libmesh_assert(!calculations_started || calculate_dphi);
249  calculate_dphi = calculate_dphiref = true; return dphidy; }
250 
255  const std::vector<std::vector<OutputShape>> & get_dphidz() const
256  { libmesh_assert(!calculations_started || calculate_dphi);
257  calculate_dphi = calculate_dphiref = true; return dphidz; }
258 
263  const std::vector<std::vector<OutputShape>> & get_dphidxi() const
264  { libmesh_assert(!calculations_started || calculate_dphiref);
265  calculate_dphiref = true; return dphidxi; }
266 
271  const std::vector<std::vector<OutputShape>> & get_dphideta() const
272  { libmesh_assert(!calculations_started || calculate_dphiref);
273  calculate_dphiref = true; return dphideta; }
274 
279  const std::vector<std::vector<OutputShape>> & get_dphidzeta() const
280  { libmesh_assert(!calculations_started || calculate_dphiref);
281  calculate_dphiref = true; return dphidzeta; }
282 
283 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
284 
289  const std::vector<std::vector<OutputTensor>> & get_d2phi() const
290  { libmesh_assert(!calculations_started || calculate_d2phi);
291  calculate_d2phi = calculate_dphiref = true; return d2phi; }
292 
297  const std::vector<std::vector<OutputShape>> & get_d2phidx2() const
298  { libmesh_assert(!calculations_started || calculate_d2phi);
299  calculate_d2phi = calculate_dphiref = true; return d2phidx2; }
300 
305  const std::vector<std::vector<OutputShape>> & get_d2phidxdy() const
306  { libmesh_assert(!calculations_started || calculate_d2phi);
307  calculate_d2phi = calculate_dphiref = true; return d2phidxdy; }
308 
313  const std::vector<std::vector<OutputShape>> & get_d2phidxdz() const
314  { libmesh_assert(!calculations_started || calculate_d2phi);
315  calculate_d2phi = calculate_dphiref = true; return d2phidxdz; }
316 
321  const std::vector<std::vector<OutputShape>> & get_d2phidy2() const
322  { libmesh_assert(!calculations_started || calculate_d2phi);
323  calculate_d2phi = calculate_dphiref = true; return d2phidy2; }
324 
329  const std::vector<std::vector<OutputShape>> & get_d2phidydz() const
330  { libmesh_assert(!calculations_started || calculate_d2phi);
331  calculate_d2phi = calculate_dphiref = true; return d2phidydz; }
332 
337  const std::vector<std::vector<OutputShape>> & get_d2phidz2() const
338  { libmesh_assert(!calculations_started || calculate_d2phi);
339  calculate_d2phi = calculate_dphiref = true; return d2phidz2; }
340 
345  const std::vector<std::vector<OutputShape>> & get_d2phidxi2() const
346  { libmesh_assert(!calculations_started || calculate_d2phi);
347  calculate_d2phi = calculate_dphiref = true; return d2phidxi2; }
348 
353  const std::vector<std::vector<OutputShape>> & get_d2phidxideta() const
354  { libmesh_assert(!calculations_started || calculate_d2phi);
356 
361  const std::vector<std::vector<OutputShape>> & get_d2phidxidzeta() const
362  { libmesh_assert(!calculations_started || calculate_d2phi);
364 
369  const std::vector<std::vector<OutputShape>> & get_d2phideta2() const
370  { libmesh_assert(!calculations_started || calculate_d2phi);
371  calculate_d2phi = calculate_dphiref = true; return d2phideta2; }
372 
377  const std::vector<std::vector<OutputShape>> & get_d2phidetadzeta() const
378  { libmesh_assert(!calculations_started || calculate_d2phi);
380 
385  const std::vector<std::vector<OutputShape>> & get_d2phidzeta2() const
386  { libmesh_assert(!calculations_started || calculate_d2phi);
387  calculate_d2phi = calculate_dphiref = true; return d2phidzeta2; }
388 
389 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
390 
391 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
392 
403  const std::vector<OutputGradient> & get_dphase() const
404  { return dphase; }
405 
406 
419  const std::vector<Real> & get_Sobolev_weight() const
420  { return weight; }
421 
427  const std::vector<RealGradient> & get_Sobolev_dweight() const
428  { return dweight; }
429 
430 #endif
431 
432 
436  void print_phi(std::ostream & os) const;
437 
442  void print_dphi(std::ostream & os) const;
443 
444 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
445 
450  void print_d2phi(std::ostream & os) const;
451 
452 #endif
453 
454 
455 protected:
456 
457 
458 
459 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
460 
466  virtual void init_base_shape_functions(const std::vector<Point> & qp,
467  const Elem * e) = 0;
468 
469 #endif
470 
475  void determine_calculations();
476 
487  virtual void compute_shape_functions(const Elem * elem, const std::vector<Point> & qp);
488 
493  std::unique_ptr<FETransformationBase<OutputType>> _fe_trans;
494 
498  std::vector<std::vector<OutputShape>> phi;
499 
503  std::vector<std::vector<OutputGradient>> dphi;
504 
508  std::vector<std::vector<OutputShape>> curl_phi;
509 
513  std::vector<std::vector<OutputDivergence>> div_phi;
514 
518  std::vector<std::vector<OutputShape>> dphidxi;
519 
523  std::vector<std::vector<OutputShape>> dphideta;
524 
528  std::vector<std::vector<OutputShape>> dphidzeta;
529 
533  std::vector<std::vector<OutputShape>> dphidx;
534 
538  std::vector<std::vector<OutputShape>> dphidy;
539 
543  std::vector<std::vector<OutputShape>> dphidz;
544 
545 
546 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
547 
551  std::vector<std::vector<OutputTensor>> d2phi;
552 
556  std::vector<std::vector<OutputShape>> d2phidxi2;
557 
561  std::vector<std::vector<OutputShape>> d2phidxideta;
562 
566  std::vector<std::vector<OutputShape>> d2phidxidzeta;
567 
571  std::vector<std::vector<OutputShape>> d2phideta2;
572 
576  std::vector<std::vector<OutputShape>> d2phidetadzeta;
577 
581  std::vector<std::vector<OutputShape>> d2phidzeta2;
582 
586  std::vector<std::vector<OutputShape>> d2phidx2;
587 
591  std::vector<std::vector<OutputShape>> d2phidxdy;
592 
596  std::vector<std::vector<OutputShape>> d2phidxdz;
597 
601  std::vector<std::vector<OutputShape>> d2phidy2;
602 
606  std::vector<std::vector<OutputShape>> d2phidydz;
607 
611  std::vector<std::vector<OutputShape>> d2phidz2;
612 
613 #endif
614 
615 
616 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
617 
618  //--------------------------------------------------------------
619  /* protected members for infinite elements, which are accessed
620  * from the outside through some inline functions
621  */
622 
623 
629  std::vector<OutputGradient> dphase;
630 
636  std::vector<RealGradient> dweight;
637 
643  std::vector<Real> weight;
644 
645 #endif
646 
647 private:
648 
649 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
650 
656  template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
657  friend class InfFE;
658 
659 #endif
660 
661 
662 };
663 
664 
665 // Typedefs for convenience and backwards compatibility
668 
669 
670 
671 
672 // ------------------------------------------------------------
673 // FEGenericBase class inline members
674 template <typename OutputType>
675 inline
677  const FEType & fet) :
678  FEAbstract(d,fet),
679  _fe_trans( FETransformationBase<OutputType>::build(fet) ),
680  phi(),
681  dphi(),
682  curl_phi(),
683  div_phi(),
684  dphidxi(),
685  dphideta(),
686  dphidzeta(),
687  dphidx(),
688  dphidy(),
689  dphidz()
690 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
691  ,d2phi(),
692  d2phidxi2(),
693  d2phidxideta(),
694  d2phidxidzeta(),
695  d2phideta2(),
696  d2phidetadzeta(),
697  d2phidzeta2(),
698  d2phidx2(),
699  d2phidxdy(),
700  d2phidxdz(),
701  d2phidy2(),
702  d2phidydz(),
703  d2phidz2()
704 #endif
705 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
706  ,dphase(),
707  dweight(),
708  weight()
709 #endif
710 {
711 }
712 
713 
714 
715 template <typename OutputType>
716 inline
718 {
719 }
720 
721 } // namespace libMesh
722 
723 #endif // LIBMESH_FE_BASE_H
FEGenericBase(const unsigned int dim, const FEType &fet)
Definition: fe_base.h:676
std::vector< std::vector< OutputTensor > > d2phi
Definition: fe_base.h:551
Manages the family, order, etc. parameters for a given FE.
Definition: fe_type.h:179
std::vector< std::vector< OutputShape > > dphidxi
Definition: fe_base.h:518
std::vector< std::vector< OutputShape > > d2phidxdz
Definition: fe_base.h:596
const std::vector< std::vector< OutputShape > > & get_d2phidxdy() const
Definition: fe_base.h:305
std::vector< std::vector< OutputShape > > dphidzeta
Definition: fe_base.h:528
std::vector< std::vector< OutputShape > > d2phidydz
Definition: fe_base.h:606
Base class for all the infinite geometric element types.
Definition: fe.h:40
std::unique_ptr< FETransformationBase< OutputType > > _fe_trans
Definition: fe_base.h:493
TensorTools::DecrementRank< OutputNumber >::type OutputNumberDivergence
Definition: fe_base.h:126
const std::vector< std::vector< OutputShape > > & get_d2phideta2() const
Definition: fe_base.h:369
virtual void init_base_shape_functions(const std::vector< Point > &qp, const Elem *e)=0
Maps between boundary ids and PeriodicBoundaryBase objects.
virtual void compute_shape_functions(const Elem *elem, const std::vector< Point > &qp)
Definition: fe_base.C:669
const std::vector< std::vector< OutputShape > > & get_d2phidxideta() const
Definition: fe_base.h:353
std::vector< std::vector< OutputShape > > d2phidxideta
Definition: fe_base.h:561
const std::vector< std::vector< OutputShape > > & get_d2phidydz() const
Definition: fe_base.h:329
The base class for all geometric element types.
Definition: elem.h:100
TensorTools::IncrementRank< OutputNumberGradient >::type OutputNumberTensor
Definition: fe_base.h:125
MeshBase & mesh
std::vector< Real > weight
Definition: fe_base.h:643
TensorTools::IncrementRank< OutputNumber >::type OutputNumberGradient
Definition: fe_base.h:124
TensorTools::IncrementRank< OutputShape >::type OutputGradient
Definition: fe_base.h:120
std::vector< std::vector< OutputShape > > d2phidx2
Definition: fe_base.h:586
std::vector< std::vector< OutputShape > > curl_phi
Definition: fe_base.h:508
const std::vector< std::vector< OutputShape > > & get_dphideta() const
Definition: fe_base.h:271
const std::vector< std::vector< OutputDivergence > > & get_div_phi() const
Definition: fe_base.h:231
std::vector< std::vector< OutputShape > > dphidy
Definition: fe_base.h:538
FEGenericBase< RealGradient > FEVectorBase
Definition: fe_base.h:667
TensorTools::DecrementRank< OutputShape >::type OutputDivergence
Definition: fe_base.h:122
std::vector< std::vector< OutputShape > > d2phidy2
Definition: fe_base.h:601
void print_phi(std::ostream &os) const
Definition: fe_base.C:707
Base class for Mesh.
Definition: mesh_base.h:77
const std::vector< std::vector< OutputShape > > & get_dphidzeta() const
Definition: fe_base.h:279
const std::vector< OutputGradient > & get_dphase() const
Definition: fe_base.h:403
std::vector< std::vector< OutputShape > > d2phidetadzeta
Definition: fe_base.h:576
const std::vector< std::vector< OutputShape > > & get_dphidx() const
Definition: fe_base.h:239
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)
Definition: fe_base.C:791
std::vector< std::vector< OutputShape > > d2phidxidzeta
Definition: fe_base.h:566
Manages the degrees of freedom (DOFs) in a simulation.
Definition: dof_map.h:176
std::vector< std::vector< OutputShape > > d2phidxdy
Definition: fe_base.h:591
std::vector< std::vector< OutputShape > > dphidx
Definition: fe_base.h:533
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:215
void print_d2phi(std::ostream &os) const
Definition: fe_base.C:776
std::vector< std::vector< OutputShape > > phi
Definition: fe_base.h:498
std::vector< std::vector< OutputShape > > d2phideta2
Definition: fe_base.h:571
const unsigned int dim
Definition: fe_abstract.h:531
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:233
TensorTools::IncrementRank< OutputGradient >::type OutputTensor
Definition: fe_base.h:121
const std::vector< std::vector< OutputShape > > & get_d2phidx2() const
Definition: fe_base.h:297
const std::vector< RealGradient > & get_Sobolev_dweight() const
Definition: fe_base.h:427
std::vector< OutputGradient > dphase
Definition: fe_base.h:629
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
const std::vector< std::vector< OutputShape > > & get_dphidy() const
Definition: fe_base.h:247
static void compute_proj_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
Definition: fe_base.C:1366
const std::vector< std::vector< OutputShape > > & get_d2phidy2() const
Definition: fe_base.h:321
const std::vector< std::vector< OutputTensor > > & get_d2phi() const
Definition: fe_base.h:289
OutputType OutputShape
Definition: fe_base.h:119
TensorTools::MakeNumber< OutputShape >::type OutputNumber
Definition: fe_base.h:123
std::vector< std::vector< OutputDivergence > > div_phi
Definition: fe_base.h:513
const std::vector< std::vector< OutputShape > > & get_d2phidxi2() const
Definition: fe_base.h:345
FEGenericBase< Real > FEBase
void determine_calculations()
Definition: fe_base.C:728
const std::vector< std::vector< OutputShape > > & get_d2phidxdz() const
Definition: fe_base.h:313
std::vector< std::vector< OutputGradient > > dphi
Definition: fe_base.h:503
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)
Definition: fe_base.C:1650
const std::vector< std::vector< OutputShape > > & get_d2phidxidzeta() const
Definition: fe_base.h:361
const std::vector< std::vector< OutputShape > > & get_d2phidetadzeta() const
Definition: fe_base.h:377
const std::vector< std::vector< OutputShape > > & get_curl_phi() const
Definition: fe_base.h:223
const std::vector< Real > & get_Sobolev_weight() const
Definition: fe_base.h:419
std::vector< std::vector< OutputShape > > d2phidz2
Definition: fe_base.h:611
const std::vector< std::vector< OutputShape > > & get_d2phidzeta2() const
Definition: fe_base.h:385
static std::unique_ptr< FEGenericBase > build_InfFE(const unsigned int dim, const FEType &type)
std::vector< std::vector< OutputShape > > d2phidxi2
Definition: fe_base.h:556
const std::vector< std::vector< OutputShape > > & get_dphidxi() const
Definition: fe_base.h:263
void print_dphi(std::ostream &os) const
Definition: fe_base.C:718
std::vector< std::vector< OutputShape > > d2phidzeta2
Definition: fe_base.h:581
std::vector< std::vector< OutputShape > > dphidz
Definition: fe_base.h:543
const std::vector< std::vector< OutputShape > > & get_d2phidz2() const
Definition: fe_base.h:337
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:207
std::vector< std::vector< OutputShape > > dphideta
Definition: fe_base.h:523
std::vector< RealGradient > dweight
Definition: fe_base.h:636
const std::vector< std::vector< OutputShape > > & get_dphidz() const
Definition: fe_base.h:255