HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
Loading...
Searching...
No Matches
elementquad.hpp
Go to the documentation of this file.
1// Class to create and store values of cell and face 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 <cassert>
25//#include <cmath>
26
27//#include <functional>
28//#include <memory>
29#include <iostream>
30#include <Eigen/Dense>
31#include <hybridcore.hpp>
32#include <quadraturerule.hpp>
33#include <mesh.hpp>
34
35
41namespace HArDCore3D {
42
47 // ----------------------------------------------------------------------------
48 // Class definition
49 // ----------------------------------------------------------------------------
50
56
57 public:
60 const HybridCore& hho,
61 const size_t iT,
62 const size_t doeT,
63 const size_t doeF
64 );
65
67 inline const QuadratureRule & get_quadT() const
68 {
69 return m_quadT;
70 };
71
73 inline const QuadratureRule & get_quadF(size_t ilF) const
74 {
75 return m_quadF[ilF];
76 };
77
79 inline const boost::multi_array<double, 2> & get_phiT_quadT() const
80 {
81 return m_phiT_quadT;
82 };
83
85 inline const boost::multi_array<VectorRd, 2> & get_dphiT_quadT() const
86 {
87 return m_dphiT_quadT;
88 };
89
91 inline const boost::multi_array<double, 2> & get_phiT_quadF(size_t ilF) const
92 {
93 return m_phiT_quadF[ilF];
94 };
95
97 inline const boost::multi_array<double, 2> get_phiF_quadF(size_t ilF) const
98 {
99 return m_phiF_quadF[ilF];
100 };
101
103 inline const boost::multi_array<VectorRd, 2> & get_dphiT_quadF(size_t ilF) const
104 {
105 return m_dphiT_quadF[ilF];
106 };
107
108
110 boost::multi_array<VectorRd, 2> get_vec_phiT_quadT(
111 size_t degree
112 ) const;
113
115 boost::multi_array<VectorRd, 2> get_vec_phiT_quadF(
116 size_t ilF,
117 size_t degree
118 ) const;
119
121 boost::multi_array<VectorRd, 2> get_vec_phiF_quadF(
122 size_t ilF,
123 size_t degree
124 ) const;
125
126
127 private:
129 const HybridCore& m_hcore; // reference to the hybridcore instance
130 const size_t m_iT; // cell number
131 const size_t m_doeT; // degree of exactness of cell quadrature rules
132 const size_t m_doeF; // degree of exactness of face quadrature rules
133
135 const Cell& m_T;
136 QuadratureRule m_quadT;
137
138 boost::multi_array<double, 2> m_phiT_quadT;
139 boost::multi_array<VectorRd, 2> m_dphiT_quadT;
140 boost::multi_array<VectorRd, 2> m_vec_phiT_quadT;
141
143 std::vector<QuadratureRule> m_quadF;
144
145 std::vector<boost::multi_array<double, 2>> m_phiT_quadF;
146 std::vector<boost::multi_array<double, 2>> m_phiF_quadF;
147 std::vector<boost::multi_array<VectorRd, 2>> m_dphiT_quadF;
148 std::vector<boost::multi_array<VectorRd, 2>> m_vec_phiT_quadF;
149 std::vector<boost::multi_array<VectorRd, 2>> m_vec_phiF_quadF;
150
151 };
152
153
155
156} // end of namespace HArDCore3D
157
158#endif /* ELEMENTQUAD_HPP */
Definition elementquad.hpp:55
Definition hybridcore.hpp:174
@ Matrix
Definition basis.hpp:67
const QuadratureRule & get_quadF(size_t ilF) const
Returns quadrature rules on face with local number ilF.
Definition elementquad.hpp:73
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:61
const boost::multi_array< double, 2 > get_phiF_quadF(size_t ilF) const
Returns values of face basis functions at face quadrature nodes, for face with local number ilF.
Definition elementquad.hpp:97
boost::multi_array< VectorRd, 2 > get_vec_phiF_quadF(size_t ilF, size_t degree) const
Builds on the fly the values of vector face basis functions at face quadrature nodes....
Definition elementquad.cpp:109
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:85
const boost::multi_array< double, 2 > & get_phiT_quadT() const
Returns values of cell basis functions at cell quadrature nodes.
Definition elementquad.hpp:79
const boost::multi_array< VectorRd, 2 > & get_dphiT_quadF(size_t ilF) const
Returns values of gradients of cell basis functions at face quadrature nodes, for face with local num...
Definition elementquad.hpp:103
const QuadratureRule & get_quadT() const
Returns quadrature rules in cell.
Definition elementquad.hpp:67
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 face quadrature nodes....
Definition elementquad.cpp:85
const boost::multi_array< double, 2 > & get_phiT_quadF(size_t ilF) const
Returns values of cell basis functions at face quadrature nodes, for face with local number ilF.
Definition elementquad.hpp:91
std::vector< QuadratureNode > QuadratureRule
Definition quadraturerule.hpp:61
Definition ddr-magnetostatics.hpp:41