HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
elementquad.hpp
Go to the documentation of this file.
1 // Class to create and store values of cell and edge basis functions on quadrature points
2 //
3 //
4 // Author: Jerome Droniou (jerome.droniou@monash.edu)
5 //
6 
7 /*
8  *
9  * This library was developed around HHO methods, although some parts of it have a more
10  * general purpose. If you use this code or part of it in a scientific publication,
11  * please mention the following book as a reference for the underlying principles
12  * of HHO schemes:
13  *
14  * The Hybrid High-Order Method for Polytopal Meshes: Design, Analysis, and Applications.
15  * D. A. Di Pietro and J. Droniou. Modeling, Simulation and Applications, vol. 19.
16  * Springer International Publishing, 2020, xxxi + 525p. doi: 10.1007/978-3-030-37203-3.
17  * url: https://hal.archives-ouvertes.fr/hal-02151813.
18  *
19  */
20 
21 #ifndef ELEMENTQUAD_HPP
22 #define ELEMENTQUAD_HPP
23 
24 #include <iostream>
25 #include <Eigen/Dense>
26 #include <hybridcore.hpp>
27 #include <quadraturerule.hpp>
28 
29 
35 namespace HArDCore2D {
36 
41  // ----------------------------------------------------------------------------
42  // Class definition
43  // ----------------------------------------------------------------------------
44 
49  class ElementQuad {
50 
51  public:
54  const HybridCore& hho,
55  const size_t iT,
56  const size_t doeT,
57  const size_t doeF
58  );
59 
61  inline const QuadratureRule & get_quadT() const
62  {
63  return m_quadT;
64  };
65 
67  inline const QuadratureRule & get_quadF(size_t ilF) const
68  {
69  return m_quadF[ilF];
70  };
71 
73  inline const boost::multi_array<double, 2> & get_phiT_quadT() const
74  {
75  return m_phiT_quadT;
76  };
77 
79  inline const boost::multi_array<VectorRd, 2> & get_dphiT_quadT() const
80  {
81  return m_dphiT_quadT;
82  };
83 
85  inline const boost::multi_array<double, 2> & get_phiT_quadF(size_t ilF) const
86  {
87  return m_phiT_quadF[ilF];
88  };
89 
91  inline const boost::multi_array<double, 2> get_phiF_quadF(size_t ilF) const
92  {
93  return m_phiF_quadF[ilF];
94  };
95 
97  inline const boost::multi_array<VectorRd, 2> & get_dphiT_quadF(size_t ilF) const
98  {
99  return m_dphiT_quadF[ilF];
100  };
101 
102 
104  boost::multi_array<VectorRd, 2> get_vec_phiT_quadT(
105  size_t degree
106  ) const;
107 
109  boost::multi_array<VectorRd, 2> get_vec_phiT_quadF(
110  size_t ilF,
111  size_t degree
112  ) const;
113 
115  boost::multi_array<VectorRd, 2> get_vec_phiF_quadF(
116  size_t ilF,
117  size_t degree
118  ) const;
119 
120 
121  private:
123  const HybridCore& m_hcore; // reference to the hybridcore instance
124  const size_t m_iT; // cell number
125  const size_t m_doeT; // degree of exactness of cell quadrature rules
126  const size_t m_doeF; // degree of exactness of edge quadrature rules
127 
129  const Cell& m_T;
130  QuadratureRule m_quadT;
131 
132  boost::multi_array<double, 2> m_phiT_quadT;
133  boost::multi_array<VectorRd, 2> m_dphiT_quadT;
134  boost::multi_array<VectorRd, 2> m_vec_phiT_quadT;
135 
137  std::vector<QuadratureRule> m_quadF;
138 
139  std::vector<boost::multi_array<double, 2>> m_phiT_quadF;
140  std::vector<boost::multi_array<double, 2>> m_phiF_quadF;
141  std::vector<boost::multi_array<VectorRd, 2>> m_dphiT_quadF;
142  std::vector<boost::multi_array<VectorRd, 2>> m_vec_phiT_quadF;
143  std::vector<boost::multi_array<VectorRd, 2>> m_vec_phiF_quadF;
144 
145  };
146 
147 
149 
150 } // end of namespace HArDCore2D
151 
152 #endif /* ELEMENTQUAD_HPP */
Definition: elementquad.hpp:49
Definition: hybridcore.hpp:179
const boost::multi_array< double, 2 > get_phiF_quadF(size_t ilF) const
Returns values of edge basis functions at edge quadrature nodes, for edge with local number ilF.
Definition: elementquad.hpp:91
const QuadratureRule & get_quadF(size_t ilF) const
Returns quadrature rules on edge with local number ilF.
Definition: elementquad.hpp:67
ElementQuad(const HybridCore &hho, const size_t iT, const size_t doeT, const size_t doeF)
Class constructor: loads the quadrature rules and values of basis functions/gradients at these points...
Definition: elementquad.cpp:17
const QuadratureRule & get_quadT() const
Returns quadrature rules in cell.
Definition: elementquad.hpp:61
boost::multi_array< VectorRd, 2 > get_vec_phiT_quadT(size_t degree) const
Builds on the fly the values of vector cell basis functions at cell quadrature nodes....
Definition: elementquad.cpp:62
const boost::multi_array< VectorRd, 2 > & get_dphiT_quadT() const
Returns values of gradients of cell basis functions at cell quadrature nodes.
Definition: elementquad.hpp:79
boost::multi_array< VectorRd, 2 > get_vec_phiF_quadF(size_t ilF, size_t degree) const
Builds on the fly the values of vector edge basis functions at edge quadrature nodes....
Definition: elementquad.cpp:110
const boost::multi_array< double, 2 > & get_phiT_quadT() const
Returns values of cell basis functions at cell quadrature nodes.
Definition: elementquad.hpp:73
const boost::multi_array< VectorRd, 2 > & get_dphiT_quadF(size_t ilF) const
Returns values of gradients of cell basis functions at edge quadrature nodes, for edge with local num...
Definition: elementquad.hpp:97
const boost::multi_array< double, 2 > & get_phiT_quadF(size_t ilF) const
Returns values of cell basis functions at edge quadrature nodes, for edge with local number ilF.
Definition: elementquad.hpp:85
boost::multi_array< VectorRd, 2 > get_vec_phiT_quadF(size_t ilF, size_t degree) const
Builds on the fly the values of vector cell basis functions at edge quadrature nodes....
Definition: elementquad.cpp:86
Polytope< DIMENSION > Cell
Definition: Polytope2D.hpp:151
std::vector< QuadratureNode > QuadratureRule
Definition: quadraturerule.hpp:55
Definition: ddr-klplate.hpp:27