HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
rtn-pk-pk-pk.hpp
Go to the documentation of this file.
1 #ifndef RTN_PK_PK_PK_HPP
2 #define RTN_PK_PK_PK_HPP
3 
4 #include <boost/fusion/include/map.hpp>
5 
6 #include "hypre.hpp"
7 
8 #include <globaldofspace.hpp>
9 
10 namespace HArDCore2D {
11 
12  namespace DSL {
13 
14  class RTNPkPkPk: public HYPRE {
15  public:
16  typedef boost::fusion::map<boost::fusion::pair<RolyCellType, std::string>,
17  boost::fusion::pair<RolyComplCellType, std::string> > CellDofsMapType;
20  const Mesh & mesh,
21  size_t degree,
22  const BoundaryConditions & bc,
23  const HYPREParameters & parameters
24  );
25 
27  Eigen::VectorXd interpolate(
28  const VelocityType & u,
29  const PressureType & p
30  ) const;
31 
32 
34  const CellDofsMapType & cellDofsMap() const
35  {
36  return m_cell_dofs_map;
37  }
38 
40  Eigen::VectorXd extendCellDofs(const Eigen::VectorXd & vh) const;
41 
43  const DiscreteSpace & extendedDofsSpace() const {
44  return *m_pkpo_pk;
45  }
46 
47  private:
49  template<typename F>
50  Eigen::VectorXd _interpolate_velocity(
51  const F & v,
52  const InterpolateParameters & parameters = {}
53  ) const;
54 
55  void _construct_local_operators(size_t iT);
56  void _construct_stabilization(size_t iT);
57  std::pair<Eigen::MatrixXd, Eigen::VectorXd> _velocity_pressure_coupling(
58  size_t iT,
60  const Eigen::VectorXd & uT_old,
61  const Eigen::VectorXd & pT_old
62  );
63  Eigen::MatrixXd _forcing_terms(
64  size_t iT,
66  );
67 
68 
69  CellDofsMapType m_cell_dofs_map;
70  GlobalDOFSpace m_rtn_dof_table;
71  std::unique_ptr<DiscreteSpace> m_pkpo_pk;
72  };
73 
74  //------------------------------------------------------------------------------
75 
76  template<typename F>
77  Eigen::VectorXd RTNPkPkPk::_interpolate_velocity(
78  const F & v,
79  const InterpolateParameters & parameters
80  ) const
81  {
82  Eigen::VectorXd vh = Eigen::VectorXd::Zero(m_velocity_space->dimension());
83 
84  // Degrees of quadrature rules
85  int doe_cell = (parameters.doe_cell > 0 ? doe_cell : 2 * m_degree + 1);
86  int doe_edge = (parameters.doe_edge > 0 ? doe_edge : 2 * m_degree + 3);
87 
88  // Interpolate at edges
89  std::function<void(size_t, size_t)> interpolate_edges
90  = [this, v, &vh, &doe_edge](size_t start, size_t end)->void
91  {
92  const auto & edge_dofs = m_velocity_space->get<PolynEdgeType>("EdgeDOFs");
93 
94  for (size_t iE = start; iE < end; iE++) {
95  const Edge & E = *this->m_velocity_space->mesh().edge(iE);
96  const auto & edge_dofs_E = edge_dofs[E.global_index()];
97  QuadratureRule quad = generate_quadrature_rule(E, doe_edge);
98  auto edge_dofs_quad = evaluate_quad<Function>::compute(*edge_dofs_E, quad);
99  vh.segment(this->m_velocity_space->globalOffset(E), this->m_velocity_space->numberOfLocalEdgeDofs())
100  = l2_projection(v, *edge_dofs_E, quad, edge_dofs_quad, GramMatrix(E, *edge_dofs_E));
101  } // for iE
102  };
103  parallel_for(m_velocity_space->mesh().n_edges(), interpolate_edges, parameters.use_threads);
104 
105  // Interpolate at cells
106  std::function<void(size_t, size_t)> interpolate_cells
107  = [this, &vh, v, &doe_cell, &doe_edge](size_t start, size_t end)->void
108  {
109  const auto & rtn_edge_dofs = m_velocity_space->get<PolyEdgeType>("RTNEdgeDOFs");
110  const auto & roly_cell_dofs = m_velocity_space->get<RolyCellType>("RolyCellDOFs");
111  const auto & roly_compl_cell_dofs = m_velocity_space->get<RolyComplCellType>("RolyComplCellDOFs");
112 
113  // Dimensions
114  size_t dim_rtn_cell_dofs = LocalSpaceDimensions<PolynCellType>::get(m_degree - 1);
115  size_t dim_rtn_edge_dofs = LocalSpaceDimensions<PolyEdgeType>::get(m_degree);
116  size_t dim_roly_cell_dofs = LocalSpaceDimensions<RolyCellType>::get(m_degree);
117  size_t dim_roly_compl_cell_dofs = LocalSpaceDimensions<RolyComplCellType>::get(m_degree + 1);
118  size_t dim_cell_dofs = dim_roly_cell_dofs + dim_roly_compl_cell_dofs;
119 
120  for (size_t iT = start; iT < end; iT++) {
121  const Cell & T = *this->m_velocity_space->mesh().cell(iT);
122 
123  // Interpolator
124  Eigen::MatrixXd MIT = Eigen::MatrixXd::Zero(dim_cell_dofs, dim_cell_dofs);
125  Eigen::VectorXd bIT = Eigen::VectorXd(dim_cell_dofs);
126 
127  MonomialCellIntegralsType int_mono_cell = IntegrateCellMonomials(T, 2 * (this->m_degree + 1));
128  QuadratureRule quad = generate_quadrature_rule(T, doe_cell);
129 
130  if (dim_rtn_cell_dofs > 0) {
131  const auto & rtn_cell_dofs = m_velocity_space->get<PolynCellType>("RTNCellDOFs")[iT];
132 
133  MIT.block(m_rtn_dof_table.localOffset(T), 0, dim_rtn_cell_dofs, dim_roly_cell_dofs)
134  = GramMatrix(T, *rtn_cell_dofs, *roly_cell_dofs[iT], int_mono_cell);
135  MIT.block(m_rtn_dof_table.localOffset(T), dim_roly_cell_dofs, dim_rtn_cell_dofs, dim_roly_compl_cell_dofs)
136  = GramMatrix(T, *rtn_cell_dofs, *roly_compl_cell_dofs[iT], int_mono_cell);
137 
138  auto rtn_cell_dofs_quad = evaluate_quad<Function>::compute(*rtn_cell_dofs, quad);
139  bIT.segment(m_rtn_dof_table.localOffset(T), dim_rtn_cell_dofs)
140  = l2_projection(v, *rtn_cell_dofs, quad, rtn_cell_dofs_quad, GramMatrix(T, *rtn_cell_dofs, int_mono_cell));
141  } // if
142 
143  for (size_t iE = 0; iE < T.n_edges(); iE++) {
144  const Edge & E = *T.edge(iE);
145  auto nTE = T.edge_normal(iE);
146 
147  const auto & rtn_edge_dofs_E = rtn_edge_dofs[E.global_index()];
148 
149  auto normal_component = [&nTE](const VectorRd & w)->double { return w.dot(nTE); };
150 
151  QuadratureRule quad_E = generate_quadrature_rule(E, doe_edge);
152  auto roly_cell_dofs_quad
153  = transform_values_quad<double>( evaluate_quad<Function>::compute(*roly_cell_dofs[iT], quad_E),
154  normal_component );
155  auto roly_compl_cell_dofs_quad
156  = transform_values_quad<double>( evaluate_quad<Function>::compute(*roly_compl_cell_dofs[iT], quad_E),
157  normal_component );
158 
159  auto rtn_edge_dofs_quad = evaluate_quad<Function>::compute(*rtn_edge_dofs[E.global_index()], quad_E);
160 
161  std::function<double(const VectorRd &)> evaluate_normal_component
162  = [v, &nTE](const VectorRd & x) { return v(x).dot(nTE); };
163 
164  MIT.block(m_rtn_dof_table.localOffset(T, E), 0, dim_rtn_edge_dofs, dim_roly_cell_dofs)
165  = compute_gram_matrix(rtn_edge_dofs_quad, roly_cell_dofs_quad, quad_E);
166  MIT.block(m_rtn_dof_table.localOffset(T, E), dim_roly_cell_dofs, dim_rtn_edge_dofs, dim_roly_compl_cell_dofs)
167  = compute_gram_matrix(rtn_edge_dofs_quad, roly_compl_cell_dofs_quad, quad_E);
168 
169  bIT.segment(m_rtn_dof_table.localOffset(T, E), dim_rtn_edge_dofs)
170  = l2_projection(evaluate_normal_component, *rtn_edge_dofs_E, quad_E, rtn_edge_dofs_quad, GramMatrix(E, *rtn_edge_dofs_E));
171  } // for iE
172 
173  vh.segment(this->m_velocity_space->globalOffset(T), this->m_velocity_space->numberOfLocalCellDofs())
174  = MIT.fullPivLu().solve(bIT);
175  } // for iT
176  };
177  parallel_for(m_velocity_space->mesh().n_cells(), interpolate_cells, m_use_threads);
178 
179  return vh;
180  }
181 
182  } // namespace DSL
183 
184 } // namespace HArDCore2D
185 
186 #endif
The BoundaryConditions class provides definition of boundary conditions.
Definition: BoundaryConditions.hpp:45
Definition: discrete-space.hpp:226
std::function< VectorRd(const VectorRd &)> VelocityType
Definition: hypre.hpp:35
size_t degree() const
Returns the degree.
Definition: hypre.hpp:62
const Mesh & mesh() const
Returns the mesh.
Definition: hypre.hpp:164
std::unique_ptr< DiscreteSpace > m_velocity_space
Definition: hypre.hpp:260
size_t m_degree
Definition: hypre.hpp:244
std::function< double(const VectorRd &)> PressureType
Definition: hypre.hpp:37
bool m_use_threads
Definition: hypre.hpp:247
Eigen::VectorXd extendCellDofs(const Eigen::VectorXd &vh) const
Extend cell DOFs from to .
const DiscreteSpace & extendedDofsSpace() const
Return the discrete space corresponding to the extended DOFs.
Definition: rtn-pk-pk-pk.hpp:43
RTNPkPkPk(const Mesh &mesh, size_t degree, const BoundaryConditions &bc, const HYPREParameters &parameters)
Constructor.
const CellDofsMapType & cellDofsMap() const
Return cell DOFs map.
Definition: rtn-pk-pk-pk.hpp:34
boost::fusion::map< boost::fusion::pair< RolyCellType, std::string >, boost::fusion::pair< RolyComplCellType, std::string > > CellDofsMapType
Definition: rtn-pk-pk-pk.hpp:17
Eigen::VectorXd interpolate(const VelocityType &u, const PressureType &p) const
Interpolate both velocity and pressure.
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition: globaldofspace.hpp:16
Definition: Mesh2D.hpp:26
function[maxeig, mineig, cond, mat]
Definition: compute_eigs.m:1
end
Definition: compute_eigs.m:12
Eigen::Vector2d VectorRd
Definition: basis.hpp:55
Eigen::MatrixXd compute_gram_matrix(const boost::multi_array< VectorRd, 2 > &B1, const boost::multi_array< double, 2 > &B2, const QuadratureRule &qr)
Compute the Gram-like matrix given a family of vector-valued and one of scalar-valued functions by te...
Definition: basis.cpp:225
static boost::multi_array< typename detail::basis_evaluation_traits< BasisType, BasisFunction >::ReturnValue, 2 > compute(const BasisType &basis, const QuadratureRule &quad)
Generic basis evaluation.
Definition: basis.hpp:2288
Eigen::VectorXd l2_projection(const std::function< typename BasisType::FunctionValue(const VectorRd &)> &f, const BasisType &basis, QuadratureRule &quad, const boost::multi_array< typename BasisType::FunctionValue, 2 > &basis_quad, const Eigen::MatrixXd &mass_basis=Eigen::MatrixXd::Zero(1, 1))
Compute the L2-projection of a function.
Definition: basis.hpp:2886
size_t localOffset(const Edge &E, const Vertex &V) const
Returns the local offset of the vertex V with respect to the edge E.
Definition: localdofspace.hpp:143
static void parallel_for(unsigned nb_elements, std::function< void(size_t start, size_t end)> functor, bool use_threads=true, unsigned nb_threads_max=1e9)
Generic function to execute threaded processes.
Definition: parallel_for.hpp:42
Polytope< DIMENSION > Cell
Definition: Polytope2D.hpp:151
Polytope< 1 > Edge
A Face is a Polytope with object_dim = DIMENSION - 1.
Definition: Polytope2D.hpp:147
std::vector< QuadratureNode > QuadratureRule
Definition: quadraturerule.hpp:55
std::unordered_map< VectorZd, double, VecHash > MonomialCellIntegralsType
Type for list of integrals of monomials.
Definition: GMpoly_cell.hpp:53
MonomialCellIntegralsType IntegrateCellMonomials(const Cell &T, const size_t maxdeg)
Compute all integrals on a cell of monomials up to a total degree, using vertex values.
Definition: GMpoly_cell.cpp:7
Eigen::MatrixXd GramMatrix(const Cell &T, const MonomialScalarBasisCell &basis1, const MonomialScalarBasisCell &basis2, MonomialCellIntegralsType mono_int_map={})
Computes the Gram Matrix of a pair of local scalar monomial bases.
Definition: GMpoly_cell.cpp:86
QuadratureRule generate_quadrature_rule(const Cell &T, const int doe, const bool force_split)
Generate quadrature rule on mesh element.
Definition: quadraturerule.cpp:10
if(strcmp(field, 'real')) % real valued entries T
Definition: mmread.m:93
const std::vector< std::unique_ptr< T > > & get(const DiscreteSpace &discrete_space, const std::string &name)
TensorizedVectorFamily< PolyCellType, dimspace > PolynCellType
Definition: discrete-space.hpp:22
TensorizedVectorFamily< Family< MonomialScalarBasisEdge >, dimspace > PolynEdgeType
Definition: discrete-space.hpp:23
Family< RolyComplBasisCell > RolyComplCellType
Definition: discrete-space.hpp:29
Family< CurlBasis< ShiftedBasis< MonomialScalarBasisCell > > > RolyCellType
Definition: discrete-space.hpp:28
Family< MonomialScalarBasisEdge > PolyEdgeType
Definition: discrete-space.hpp:21
Definition: ddr-klplate.hpp:27
static auto v
Definition: ddrcore-test.hpp:32
Definition: hypre.hpp:21
Definition: hho-interpolate.hpp:15
Evaluate a basis at quadrature nodes. 'BasisFunction' (=Function, Gradient, Curl, Divergence,...
Definition: basis.hpp:2284