Manicore
Library to implement schemes on n-dimensionnal manifolds.
ddr_spaces.hpp
Go to the documentation of this file.
1 #ifndef DDR_SPACES_HPP
2 #define DDR_SPACES_HPP
3 
4 #include "dofspace.hpp"
5 #include "pec.hpp"
6 
7 #include <memory>
8 #include <iostream>
9 
13 namespace Manicore {
14 
17 
19 
27  template<size_t dimension>
28  class DDR_Spaces {
29  public:
32  int r,
33  bool use_threads = true,
34  std::array<int,dimension> const * dqr_p = nullptr,
35  std::ostream & output = std::cout
36  );
37 
39 
40  template<size_t k>
41  using FunctionType = std::function<Eigen::Vector<double,Dimension::ExtDim(k,dimension)>(size_t map_id, const Eigen::Vector<double,dimension> &)>;
42 
44 
45  template<size_t k >
46  Eigen::VectorXd interpolate(FunctionType<k> const & func,
47  std::array<int,dimension> const * dqr_p = nullptr ) const;
48 
50 
51  const Eigen::MatrixXd & full_diff(size_t k , size_t d ,size_t i_cell ) const // Return \star d in PL(r,d-k-1,d)
52  {
53  assert(d <= dimension && k < d && i_cell < _ops[d].size() && "Access of diff out of range");
54  return _ops[d][i_cell].full_diff[k];
55  }
57 
58  Eigen::MatrixXd compose_diff(size_t k , size_t d ,size_t i_cell ) const;
60 
61  const Eigen::MatrixXd & potential(size_t k , size_t d ,size_t i_cell ) const // Return \star P^k in PL(r,d-k,d)
62  {
63  assert(d <= dimension && k <= d && i_cell< _ops[d].size() && "Access of potential out of range");
64  return _ops[d][i_cell].P[k];
65  }
66 
68  DOFSpace<dimension> const & dofspace(size_t k) const {
69  return _dofspace[k];
70  }
72  const Mesh<dimension>* mesh() const {return _mesh;}
74  int degree() const {return _r;}
75 
77 
80  Eigen::MatrixXd computeL2Product(size_t k ,size_t i_cell ) const;
82 
84  double hmax() const;
85 
86  private:
87  struct DDR_Operators { // one for each form degree k
88  std::array<Eigen::MatrixXd,dimension+1> full_diff; // equal to \star d, PL(r,d-k-1,d) valued
89  std::array<Eigen::MatrixXd,dimension+1> diff; // equal to \star \ul{d}, PLtrimmed(r,d-k-1,d) valued
90  std::array<Eigen::MatrixXd,dimension+1> P; // equal to \star P, PL(r,d-k,d) valued
91  };
92 
93  const Mesh<dimension>* _mesh;
94  int _r;
95  bool _use_threads;
96  std::unique_ptr<PEC<dimension>> _ddr, _ddr_po;
97  std::array<DOFSpace<dimension>,dimension+1> _dofspace; // one for all form degree
98  std::array<std::vector<DDR_Operators>,dimension+1> _ops; // one for each cell dimension
99  };
101 }
102 #endif
103 
Implement the discrete operators of DDR-PEC.
Definition: ddr_spaces.hpp:28
Eigen::MatrixXd computeL2Product(size_t k, size_t i_cell) const
Compute the local product (including the cell boundary)
Definition: ddr_spaces.cpp:272
const Mesh< dimension > * mesh() const
Return the Mesh associated with the class.
Definition: ddr_spaces.hpp:72
DDR_Spaces(Mesh< dimension > const &mesh, int r, bool use_threads=true, std::array< int, dimension > const *dqr_p=nullptr, std::ostream &output=std::cout)
Definition: ddr_spaces.cpp:29
Eigen::VectorXd interpolate(FunctionType< k > const &func, std::array< int, dimension > const *dqr_p=nullptr) const
Interpolate the given function on the discrete spaces.
Definition: ddr_spaces.cpp:176
DOFSpace< dimension > const & dofspace(size_t k) const
Return DOFSpace associated to k forms.
Definition: ddr_spaces.hpp:68
std::function< Eigen::Vector< double, Dimension::ExtDim(k, dimension)>(size_t map_id, const Eigen::Vector< double, dimension > &)> FunctionType
Signature of the function to use in interpolate.
Definition: ddr_spaces.hpp:41
const Eigen::MatrixXd & potential(size_t k, size_t d, size_t i_cell) const
Compute the potential operator.
Definition: ddr_spaces.hpp:61
double hmax() const
Return the mesh size.
Definition: ddr_spaces.cpp:324
Eigen::MatrixXd compose_diff(size_t k, size_t d, size_t i_cell) const
Compute the projected differential operator.
Definition: ddr_spaces.cpp:248
int degree() const
Return the polynomial degree associated with the class.
Definition: ddr_spaces.hpp:74
const Eigen::MatrixXd & full_diff(size_t k, size_t d, size_t i_cell) const
Compute the full differential operator.
Definition: ddr_spaces.hpp:51
Convert between local and global data.
Definition: dofspace.hpp:19
Main data structure for the mesh.
Definition: mesh.hpp:54
Provides the Local-Global mapping of unknowns.
Definition: maxwell.hpp:21
Discrete spaces for DDR-PEC.