inf_fe.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_INF_FE_H
21 #define LIBMESH_INF_FE_H
22 
23 #include "libmesh/libmesh_config.h"
24 
25 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
26 
27 // Local includes
28 #include "libmesh/fe_base.h"
29 
30 // C++ includes
31 #include <cstddef>
32 
33 namespace libMesh
34 {
35 
36 
37 // forward declarations
38 class Elem;
39 class FEComputeData;
40 
41 
42 
74 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
75 class InfFE : public FEBase
76 {
77 
78  /*
79  * Protect the nested class
80  */
81 protected:
82 
93  class Radial
94  {
95  private:
96 
100  Radial() {}
101 
102  public:
103 
108  static Real decay (const Real v);
109 
114  static Real decay_deriv (const Real) { return -.5; }
115 
120  static Real D (const Real v) { return (1.-v)*(1.-v)/4.; }
121 
125  static Real D_deriv (const Real v) { return (v-1.)/2.; }
126 
131  static Order mapping_order() { return FIRST; }
132 
147  static unsigned int n_dofs (const Order o_radial)
148  { return static_cast<unsigned int>(o_radial)+1; }
149 
159  static unsigned int n_dofs_at_node (const Order o_radial,
160  const unsigned int n_onion);
161 
170  static unsigned int n_dofs_per_elem (const Order o_radial)
171  { return static_cast<unsigned int>(o_radial)+1; }
172 
173  };
174 
175 
176 
185  class Base
186  {
187  private:
188 
192  Base() {}
193 
194  public:
195 
201  static Elem * build_elem (const Elem * inf_elem);
202 
208  static ElemType get_elem_type (const ElemType type);
209 
215  static unsigned int n_base_mapping_sf (const ElemType base_elem_type,
216  const Order base_mapping_order);
217 
218  };
219 
220 
221 
222 
223 
224 
225 
226 public:
227 
228  // InfFE continued
229 
243  explicit
244  InfFE(const FEType & fet);
245  ~InfFE() {}
246 
247  // The static public members for access from FEInterface etc
248 
261  static Real shape(const FEType & fet,
262  const ElemType t,
263  const unsigned int i,
264  const Point & p);
265 
278  static Real shape(const FEType & fet,
279  const Elem * elem,
280  const unsigned int i,
281  const Point & p);
282 
293  static void compute_data(const FEType & fe_t,
294  const Elem * inf_elem,
295  FEComputeData & data);
296 
301  static unsigned int n_shape_functions (const FEType & fet,
302  const ElemType t)
303  { return n_dofs(fet, t); }
304 
311  static unsigned int n_dofs(const FEType & fet,
312  const ElemType inf_elem_type);
313 
318  static unsigned int n_dofs_at_node(const FEType & fet,
319  const ElemType inf_elem_type,
320  const unsigned int n);
321 
326  static unsigned int n_dofs_per_elem(const FEType & fet,
327  const ElemType inf_elem_type);
328 
332  virtual FEContinuity get_continuity() const override
333  { return C_ZERO; } // FIXME - is this true??
334 
339  virtual bool is_hierarchic() const override
340  { return false; } // FIXME - Inf FEs don't handle p elevation yet
341 
349  static void nodal_soln(const FEType & fet,
350  const Elem * elem,
351  const std::vector<Number> & elem_soln,
352  std::vector<Number> & nodal_soln);
353 
358  static Point map (const Elem * inf_elem,
359  const Point & reference_point);
360 
376  static Point inverse_map (const Elem * elem,
377  const Point & p,
378  const Real tolerance = TOLERANCE,
379  const bool secure = true);
380 
381 
388  static void inverse_map (const Elem * elem,
389  const std::vector<Point> & physical_points,
390  std::vector<Point> & reference_points,
391  const Real tolerance = TOLERANCE,
392  const bool secure = true);
393 
394 
395  // The workhorses of InfFE. These are often used during matrix assembly.
396 
404  virtual void reinit (const Elem * elem,
405  const std::vector<Point> * const pts = nullptr,
406  const std::vector<Real> * const weights = nullptr) override;
407 
413  virtual void reinit (const Elem * elem,
414  const unsigned int side,
415  const Real tolerance = TOLERANCE,
416  const std::vector<Point> * const pts = nullptr,
417  const std::vector<Real> * const weights = nullptr) override;
418 
424  virtual void edge_reinit (const Elem * elem,
425  const unsigned int edge,
426  const Real tolerance = TOLERANCE,
427  const std::vector<Point> * const pts = nullptr,
428  const std::vector<Real> * const weights = nullptr) override;
429 
434  virtual void side_map (const Elem * /* elem */,
435  const Elem * /* side */,
436  const unsigned int /* s */,
437  const std::vector<Point> & /* reference_side_points */,
438  std::vector<Point> & /* reference_points */) override
439  {
440  libmesh_not_implemented();
441  }
442 
454  virtual void attach_quadrature_rule (QBase * q) override;
455 
460  virtual unsigned int n_shape_functions () const override
461  { return _n_total_approx_sf; }
462 
468  virtual unsigned int n_quadrature_points () const override
469  { libmesh_assert(radial_qrule); return _n_total_qp; }
470 
471 
472 protected:
473 
474  // static members used by the workhorses
475 
492  static Real eval(Real v,
493  Order o_radial,
494  unsigned int i);
495 
501  static Real eval_deriv(Real v,
502  Order o_radial,
503  unsigned int i);
504 
505 
506 
507  // Non-static members used by the workhorses
508 
513  void update_base_elem (const Elem * inf_elem);
514 
518  virtual void init_base_shape_functions(const std::vector<Point> &,
519  const Elem *) override
520  { libmesh_not_implemented(); }
521 
527  void init_radial_shape_functions(const Elem * inf_elem,
528  const std::vector<Point> * radial_pts = nullptr);
529 
536  void init_shape_functions(const std::vector<Point> & radial_qp,
537  const std::vector<Point> & base_qp,
538  const Elem * inf_elem);
539 
544  void init_face_shape_functions (const std::vector<Point> & qp,
545  const Elem * side);
546 
555  void combine_base_radial(const Elem * inf_elem);
556 
567  virtual void compute_shape_functions(const Elem *, const std::vector<Point> &) override;
568 
569 
570 
571  // Miscellaneous static members
572 
579  static void compute_node_indices (const ElemType inf_elem_type,
580  const unsigned int outer_node_index,
581  unsigned int & base_node,
582  unsigned int & radial_node);
583 
592  static void compute_node_indices_fast (const ElemType inf_elem_type,
593  const unsigned int outer_node_index,
594  unsigned int & base_node,
595  unsigned int & radial_node);
596 
603  static void compute_shape_indices (const FEType & fet,
604  const ElemType inf_elem_type,
605  const unsigned int i,
606  unsigned int & base_shape,
607  unsigned int & radial_shape);
608 
612  std::vector<Real> dist;
613 
620  std::vector<Real> dweightdv;
621 
629  std::vector<Real> som;
634  std::vector<Real> dsomdv;
635 
640  std::vector<std::vector<Real>> mode;
641 
646  std::vector<std::vector<Real>> dmodedv;
647 
651  std::vector<std::vector<Real>> radial_map;
652 
653 
657  std::vector<std::vector<Real>> dradialdv_map;
658 
664  std::vector<Real> dphasedxi;
665 
671  std::vector<Real> dphasedeta;
672 
678  std::vector<Real> dphasedzeta;
679 
680 
681 
682 
683  // numbering scheme maps
684 
693  std::vector<unsigned int> _radial_node_index;
694 
703  std::vector<unsigned int> _base_node_index;
704 
713  std::vector<unsigned int> _radial_shape_index;
714 
723  std::vector<unsigned int> _base_shape_index;
724 
725 
726 
727 
728  // some more protected members
729 
734  unsigned int _n_total_approx_sf;
735 
740  unsigned int _n_total_qp;
741 
746  std::vector<Real> _total_qrule_weights;
747 
752  std::unique_ptr<QBase> base_qrule;
753 
758  std::unique_ptr<QBase> radial_qrule;
759 
764  std::unique_ptr<Elem> base_elem;
765 
772  std::unique_ptr<FEBase> base_fe;
773 
783 
784 
785 private:
786 
790  virtual bool shapes_need_reinit() const override;
791 
800 
801 
802 #ifdef DEBUG
803 
808  static bool _warned_for_shape;
809 
810 #endif
811 
817  template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
818  friend class InfFE;
819 
820 };
821 
822 
823 
824 // InfFE::Radial class inline members
825 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
826 inline
828 {
829  switch (Dim)
830  //TODO:[DD] What decay do i have in 2D and 1D?
831  {
832  case 3:
833  return (1.-v)/2.;
834 
835  case 2:
836  return 0.;
837 
838  case 1:
839  return 0.;
840 
841  default:
842  libmesh_error_msg("Invalid Dim = " << Dim);
843  }
844 }
845 
846 } // namespace libMesh
847 
848 
849 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
850 
851 
852 #endif // LIBMESH_INF_FE_H
void init_face_shape_functions(const std::vector< Point > &qp, const Elem *side)
static void compute_node_indices_fast(const ElemType inf_elem_type, const unsigned int outer_node_index, unsigned int &base_node, unsigned int &radial_node)
Manages the family, order, etc. parameters for a given FE.
Definition: fe_type.h:179
virtual void init_base_shape_functions(const std::vector< Point > &, const Elem *) override
Definition: inf_fe.h:518
static unsigned int n_dofs_per_elem(const FEType &fet, const ElemType inf_elem_type)
static Real decay_deriv(const Real)
Definition: inf_fe.h:114
std::vector< Real > dsomdv
Definition: inf_fe.h:634
static bool _warned_for_shape
Definition: inf_fe.h:808
std::vector< std::vector< Real > > dradialdv_map
Definition: inf_fe.h:657
std::vector< std::vector< Real > > radial_map
Definition: inf_fe.h:651
static Real D(const Real v)
Definition: inf_fe.h:120
Base class for all the infinite geometric element types.
Definition: fe.h:40
std::vector< Real > dphasedxi
Definition: inf_fe.h:664
static unsigned int n_dofs_at_node(const FEType &fet, const ElemType inf_elem_type, const unsigned int n)
Definition: inf_fe_static.C:72
void combine_base_radial(const Elem *inf_elem)
Definition: inf_fe.C:733
static Real eval_deriv(Real v, Order o_radial, unsigned int i)
std::vector< Real > dweightdv
Definition: inf_fe.h:620
Helper class used with FEInterface::compute_data().
std::unique_ptr< QBase > radial_qrule
Definition: inf_fe.h:758
std::vector< unsigned int > _base_node_index
Definition: inf_fe.h:703
void init_shape_functions(const std::vector< Point > &radial_qp, const std::vector< Point > &base_qp, const Elem *inf_elem)
Definition: inf_fe.C:404
unsigned short int side
Definition: xdr_io.C:50
unsigned int _n_total_qp
Definition: inf_fe.h:740
static Elem * build_elem(const Elem *inf_elem)
The base class for all geometric element types.
Definition: elem.h:100
static Real shape(const FEType &fet, const ElemType t, const unsigned int i, const Point &p)
static unsigned int n_dofs_at_node(const Order o_radial, const unsigned int n_onion)
static const Real TOLERANCE
static unsigned int n_shape_functions(const FEType &fet, const ElemType t)
Definition: inf_fe.h:301
static ElemType get_elem_type(const ElemType type)
static void compute_node_indices(const ElemType inf_elem_type, const unsigned int outer_node_index, unsigned int &base_node, unsigned int &radial_node)
virtual FEContinuity get_continuity() const override
Definition: inf_fe.h:332
std::vector< std::vector< Real > > mode
Definition: inf_fe.h:640
virtual void compute_shape_functions(const Elem *, const std::vector< Point > &) override
Definition: inf_fe.C:878
friend class InfFE
Definition: inf_fe.h:818
static unsigned int n_dofs_per_elem(const Order o_radial)
Definition: inf_fe.h:170
static unsigned int n_dofs(const FEType &fet, const ElemType inf_elem_type)
Definition: inf_fe_static.C:55
virtual unsigned int n_quadrature_points() const override
Definition: inf_fe.h:468
std::unique_ptr< Elem > base_elem
Definition: inf_fe.h:764
static unsigned int n_base_mapping_sf(const ElemType base_elem_type, const Order base_mapping_order)
std::unique_ptr< FEBase > base_fe
Definition: inf_fe.h:772
static Real decay(const Real v)
Definition: inf_fe.h:827
std::vector< Real > dphasedzeta
Definition: inf_fe.h:678
static Real eval(Real v, Order o_radial, unsigned int i)
FEType current_fe_type
Definition: inf_fe.h:782
std::vector< Real > dphasedeta
Definition: inf_fe.h:671
static Real D_deriv(const Real v)
Definition: inf_fe.h:125
static unsigned int n_dofs(const Order o_radial)
Definition: inf_fe.h:147
FEGenericBase< Real > FEBase
virtual bool shapes_need_reinit() const override
Definition: inf_fe.C:980
unsigned int _n_total_approx_sf
Definition: inf_fe.h:734
void init_radial_shape_functions(const Elem *inf_elem, const std::vector< Point > *radial_pts=nullptr)
Definition: inf_fe.C:320
std::vector< unsigned int > _radial_node_index
Definition: inf_fe.h:693
static Point inverse_map(const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: inf_fe_map.C:90
std::vector< Real > som
Definition: inf_fe.h:629
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void edge_reinit(const Elem *elem, const unsigned int edge, const Real tolerance=TOLERANCE, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
std::unique_ptr< QBase > base_qrule
Definition: inf_fe.h:752
static void compute_shape_indices(const FEType &fet, const ElemType inf_elem_type, const unsigned int i, unsigned int &base_shape, unsigned int &radial_shape)
static Point map(const Elem *inf_elem, const Point &reference_point)
Definition: inf_fe_map.C:40
std::vector< unsigned int > _radial_shape_index
Definition: inf_fe.h:713
static bool _warned_for_nodal_soln
Definition: inf_fe.h:807
static void compute_data(const FEType &fe_t, const Elem *inf_elem, FEComputeData &data)
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
Definition: inf_fe.C:109
static void nodal_soln(const FEType &fet, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
IterBase * data
std::vector< unsigned int > _base_shape_index
Definition: inf_fe.h:723
virtual void attach_quadrature_rule(QBase *q) override
Definition: inf_fe.C:69
std::vector< Real > dist
Definition: inf_fe.h:612
virtual void side_map(const Elem *, const Elem *, const unsigned int, const std::vector< Point > &, std::vector< Point > &) override
Definition: inf_fe.h:434
static Order mapping_order()
Definition: inf_fe.h:131
static ElemType _compute_node_indices_fast_current_elem_type
Definition: inf_fe.h:799
std::vector< std::vector< Real > > dmodedv
Definition: inf_fe.h:646
A geometric point in (x,y,z) space.
Definition: point.h:38
virtual unsigned int n_shape_functions() const override
Definition: inf_fe.h:460
virtual bool is_hierarchic() const override
Definition: inf_fe.h:339
Base class for all quadrature families and orders.
Definition: quadrature.h:62
std::vector< Real > _total_qrule_weights
Definition: inf_fe.h:746
void update_base_elem(const Elem *inf_elem)
Definition: inf_fe.C:98