11 #ifndef BEMBEL_SRC_ANSATZSPACE_ANSATZSPACE_HPP_
12 #define BEMBEL_SRC_ANSATZSPACE_ANSATZSPACE_HPP_
24 template <
typename Derived>
40 super_space_ = other.super_space_;
41 knot_repetition_ = other.knot_repetition_;
42 transformation_matrix_ = other.transformation_matrix_;
49 super_space_ = other.super_space_;
50 knot_repetition_ = other.knot_repetition_;
51 transformation_matrix_ = other.transformation_matrix_;
63 super_space_ = other.super_space_;
64 knot_repetition_ = other.knot_repetition_;
65 transformation_matrix_ = other.transformation_matrix_;
81 int polynomial_degree,
int knot_repetition = 1) {
103 int polynomial_degree,
int knot_repetition) {
104 knot_repetition_ = knot_repetition;
105 super_space_.init_SuperSpace(geometry, refinement_level, polynomial_degree);
108 transformation_matrix_ =
136 return super_space_.get_refinement_level();
145 return super_space_.get_polynomial_degree();
155 return super_space_.get_number_of_elements();
164 return super_space_.get_number_of_patches();
180 return super_space_.get_geometry();
190 return transformation_matrix_;
196 Eigen::SparseMatrix<double> transformation_matrix_;
198 int knot_repetition_;
The AnsatzSpace is the class that handles the assembly of the discrete basis.
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.
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.
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.
AnsatzSpace(const AnsatzSpace &other)
Copy constructor.
AnsatzSpace & operator=(AnsatzSpace other)
Copy assignment operator.
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
This class takes care of identifying DOFs on different edges, which must be identified with one anoth...
Eigen::SparseMatrix< double > get_glue_matrix() const
Returns Glue matrix to assemble continuous B-splines.
The Projector provides routines to assemble smooth B-Splines on each patch.
const Eigen::SparseMatrix< double > & get_projection_matrix()
This function returns the transformation matrix to assemble smooth tensor product B-splines.
std::vector< Bembel::Patch > PatchVector
typedef for PatchVector
Routines for the evalutation of pointwise errors.
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...