27   template<
size_t dimension>
 
   33                  bool use_threads = 
true, 
 
   34                  std::array<int,dimension> 
const * dqr_p = 
nullptr, 
 
   35                  std::ostream & output = std::cout 
 
   41       using FunctionType = std::function<Eigen::Vector<double,Dimension::ExtDim(k,dimension)>(
size_t map_id, 
const Eigen::Vector<double,dimension> &)>;
 
   47                                   std::array<int,dimension> 
const * dqr_p = 
nullptr ) 
const;
 
   51       const Eigen::MatrixXd & 
full_diff(
size_t k , 
size_t d ,
size_t i_cell ) 
const  
   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];
 
   58       Eigen::MatrixXd 
compose_diff(
size_t k , 
size_t d ,
size_t i_cell ) 
const; 
 
   61       const Eigen::MatrixXd & 
potential(
size_t k , 
size_t d ,
size_t i_cell ) 
const  
   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];
 
   87       struct DDR_Operators { 
 
   88         std::array<Eigen::MatrixXd,dimension+1> 
full_diff; 
 
   89         std::array<Eigen::MatrixXd,dimension+1> diff; 
 
   90         std::array<Eigen::MatrixXd,dimension+1> P; 
 
   93       const Mesh<dimension>* _mesh;
 
   96       std::unique_ptr<PEC<dimension>> _ddr, _ddr_po;
 
   97       std::array<DOFSpace<dimension>,dimension+1> _dofspace; 
 
   98       std::array<std::vector<DDR_Operators>,dimension+1> _ops; 
 
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:23
 
Discrete spaces for DDR-PEC.