Bembel
AnsatzSpace.hpp
1 // This file is part of Bembel, the higher order C++ boundary element library.
2 //
3 // Copyright (C) 2024 see <http://www.bembel.eu>
4 //
5 // It was written as part of a cooperation of J. Doelz, H. Harbrecht, S. Kurz,
6 // M. Multerer, S. Schoeps, and F. Wolf at Technische Universitaet Darmstadt,
7 // Universitaet Basel, and Universita della Svizzera italiana, Lugano. This
8 // source code is subject to the GNU General Public License version 3 and
9 // provided WITHOUT ANY WARRANTY, see <http://www.bembel.eu> for further
10 // information.
11 #ifndef BEMBEL_SRC_ANSATZSPACE_ANSATZSPACE_HPP_
12 #define BEMBEL_SRC_ANSATZSPACE_ANSATZSPACE_HPP_
13 
14 namespace Bembel {
24 template <typename Derived>
25 class AnsatzSpace {
26  public:
29  // constructors
31 
39  AnsatzSpace(const AnsatzSpace &other) {
40  super_space_ = other.super_space_;
41  knot_repetition_ = other.knot_repetition_;
42  transformation_matrix_ = other.transformation_matrix_;
43  }
49  super_space_ = other.super_space_;
50  knot_repetition_ = other.knot_repetition_;
51  transformation_matrix_ = other.transformation_matrix_;
52  }
63  super_space_ = other.super_space_;
64  knot_repetition_ = other.knot_repetition_;
65  transformation_matrix_ = other.transformation_matrix_;
66  return *this;
67  }
80  AnsatzSpace(const Geometry &geometry, int refinement_level,
81  int polynomial_degree, int knot_repetition = 1) {
82  init_AnsatzSpace(geometry, refinement_level, polynomial_degree,
83  knot_repetition);
84  return;
85  }
87  // init_Ansatzspace
89 
102  void init_AnsatzSpace(const Geometry &geometry, int refinement_level,
103  int polynomial_degree, int knot_repetition) {
104  knot_repetition_ = knot_repetition;
105  super_space_.init_SuperSpace(geometry, refinement_level, polynomial_degree);
106  Projector<Derived> proj(super_space_, knot_repetition_);
107  Glue<Derived> glue(super_space_, proj);
108  transformation_matrix_ =
109  proj.get_projection_matrix() * glue.get_glue_matrix();
110  return;
111  }
113  // getters
115 
121  const SuperSpace<Derived> &get_superspace() const { return super_space_; }
122 
128  int get_knot_repetition() const { return knot_repetition_; }
129 
135  int get_refinement_level() const {
136  return super_space_.get_refinement_level();
137  }
138 
144  int get_polynomial_degree() const {
145  return super_space_.get_polynomial_degree();
146  }
147 
155  return super_space_.get_number_of_elements();
156  }
157 
163  int get_number_of_patches() const {
164  return super_space_.get_number_of_patches();
165  }
166 
172  int get_number_of_dofs() const { return transformation_matrix_.cols(); }
173 
179  const PatchVector &get_geometry() const {
180  return super_space_.get_geometry();
181  }
182 
189  const Eigen::SparseMatrix<double> &get_transformation_matrix() const {
190  return transformation_matrix_;
191  }
193  // private member variables
195  private:
196  Eigen::SparseMatrix<double> transformation_matrix_;
197  SuperSpace<Derived> super_space_;
198  int knot_repetition_;
199 };
200 } // namespace Bembel
201 #endif // BEMBEL_SRC_ANSATZSPACE_ANSATZSPACE_HPP_
The AnsatzSpace is the class that handles the assembly of the discrete basis.
Definition: AnsatzSpace.hpp:25
int get_number_of_patches() const
Retrieves the number of patches of the underlying Geometry.
AnsatzSpace(const Geometry &geometry, int refinement_level, int polynomial_degree, int knot_repetition=1)
Constructor for AnsatzSpace.
Definition: AnsatzSpace.hpp:80
int get_polynomial_degree() const
Retrieves the polynomial degree of this AnsatzSpace.
void init_AnsatzSpace(const Geometry &geometry, int refinement_level, int polynomial_degree, int knot_repetition)
Initializes the AnsatzSpace object.
const Eigen::SparseMatrix< double > & get_transformation_matrix() const
Retrieves the transformation matrix associated with this AnsatzSpace.
const PatchVector & get_geometry() const
Retrieves the geometry associated with this AnsatzSpace.
int get_knot_repetition() const
Retrieves the knot repetition value of this AnsatzSpace.
AnsatzSpace()
Default constructor.
Definition: AnsatzSpace.hpp:34
int get_refinement_level() const
Retrieves the refinement level of this AnsatzSpace.
const SuperSpace< Derived > & get_superspace() const
Retrieves the reference to the SuperSpace associated with this AnsatzSpace.
AnsatzSpace(AnsatzSpace &&other)
Move constructor.
Definition: AnsatzSpace.hpp:48
AnsatzSpace(const AnsatzSpace &other)
Copy constructor.
Definition: AnsatzSpace.hpp:39
AnsatzSpace & operator=(AnsatzSpace other)
Copy assignment operator.
Definition: AnsatzSpace.hpp:62
int get_number_of_elements() const
Retrieves the number of elements of the underlying ElementTree in the SuperSpace of this AnsatzSpace.
int get_number_of_dofs() const
Retrieves the number of degrees of freedom of this AnsatzSpace.
this class wraps a GeometryVector and provides some basic functionality, like reading Geometry files
Definition: Geometry.hpp:20
This class takes care of identifying DOFs on different edges, which must be identified with one anoth...
Definition: Glue.hpp:76
Eigen::SparseMatrix< double > get_glue_matrix() const
Returns Glue matrix to assemble continuous B-splines.
Definition: Glue.hpp:120
The Projector provides routines to assemble smooth B-Splines on each patch.
Definition: Projector.hpp:40
const Eigen::SparseMatrix< double > & get_projection_matrix()
This function returns the transformation matrix to assemble smooth tensor product B-splines.
Definition: Projector.hpp:99
std::vector< Bembel::Patch > PatchVector
typedef for PatchVector
Definition: PatchVector.hpp:24
Routines for the evalutation of pointwise errors.
Definition: AnsatzSpace.hpp:14
struct containing specifications on the linear operator has to be specialized or derived for any part...
The superspace manages local polynomial bases on each element of the mesh and provides an interface t...
Definition: SuperSpace.hpp:22