Bembel
Bembel::AnsatzSpace< Derived > Class Template Reference

The AnsatzSpace is the class that handles the assembly of the discrete basis. More...

#include <AnsatzSpace.hpp>

Detailed Description

template<typename Derived>
class Bembel::AnsatzSpace< Derived >

The AnsatzSpace is the class that handles the assembly of the discrete basis.

It invokes a SuperSpace and uses the Glue and Projector class to assemble a transformation matrix, which relates the SuperSpace to the desired basis.

Definition at line 25 of file AnsatzSpace.hpp.

Public Member Functions

 AnsatzSpace ()
 Default constructor.
 
 AnsatzSpace (const AnsatzSpace &other)
 Copy constructor. More...
 
 AnsatzSpace (AnsatzSpace &&other)
 Move constructor. More...
 
AnsatzSpaceoperator= (AnsatzSpace other)
 Copy assignment operator. More...
 
 AnsatzSpace (const Geometry &geometry, int refinement_level, int polynomial_degree, int knot_repetition=1)
 Constructor for AnsatzSpace. More...
 
void init_AnsatzSpace (const Geometry &geometry, int refinement_level, int polynomial_degree, int knot_repetition)
 Initializes the AnsatzSpace object. More...
 
const SuperSpace< Derived > & get_superspace () const
 Retrieves the reference to the SuperSpace associated with this AnsatzSpace. More...
 
int get_knot_repetition () const
 Retrieves the knot repetition value of this AnsatzSpace. More...
 
int get_refinement_level () const
 Retrieves the refinement level of this AnsatzSpace. More...
 
int get_polynomial_degree () const
 Retrieves the polynomial degree of this AnsatzSpace. More...
 
int get_number_of_elements () const
 Retrieves the number of elements of the underlying ElementTree in the SuperSpace of this AnsatzSpace. More...
 
int get_number_of_patches () const
 Retrieves the number of patches of the underlying Geometry. More...
 
int get_number_of_dofs () const
 Retrieves the number of degrees of freedom of this AnsatzSpace. More...
 
const PatchVectorget_geometry () const
 Retrieves the geometry associated with this AnsatzSpace. More...
 
const Eigen::SparseMatrix< double > & get_transformation_matrix () const
 Retrieves the transformation matrix associated with this AnsatzSpace. More...
 

Public Types

enum  { Form = LinearOperatorTraits<Derived>::Form }
 

Constructor & Destructor Documentation

◆ AnsatzSpace() [1/3]

template<typename Derived >
Bembel::AnsatzSpace< Derived >::AnsatzSpace ( const AnsatzSpace< Derived > &  other)
inline

Copy constructor.

Parameters
otherThe object to copy from

Definition at line 39 of file AnsatzSpace.hpp.

◆ AnsatzSpace() [2/3]

template<typename Derived >
Bembel::AnsatzSpace< Derived >::AnsatzSpace ( AnsatzSpace< Derived > &&  other)
inline

Move constructor.

Parameters
otherThe object to move from

Definition at line 48 of file AnsatzSpace.hpp.

◆ AnsatzSpace() [3/3]

template<typename Derived >
Bembel::AnsatzSpace< Derived >::AnsatzSpace ( const Geometry geometry,
int  refinement_level,
int  polynomial_degree,
int  knot_repetition = 1 
)
inline

Constructor for AnsatzSpace.

This constructor initializes an AnsatzSpace object with the provided parameters.

Parameters
geometryThe geometry object defining the space.
refinement_levelThe refinement level of the space.
polynomial_degreeThe degree of polynomials used in the space.
knot_repetition(optional) The number of repetitions of knots in the space.

Definition at line 80 of file AnsatzSpace.hpp.

Member Function Documentation

◆ get_geometry()

template<typename Derived >
const PatchVector& Bembel::AnsatzSpace< Derived >::get_geometry ( ) const
inline

Retrieves the geometry associated with this AnsatzSpace.

Returns
A const reference to the geometry.

Definition at line 179 of file AnsatzSpace.hpp.

◆ get_knot_repetition()

template<typename Derived >
int Bembel::AnsatzSpace< Derived >::get_knot_repetition ( ) const
inline

Retrieves the knot repetition value of this AnsatzSpace.

Returns
The knot repetition value.

Definition at line 128 of file AnsatzSpace.hpp.

◆ get_number_of_dofs()

template<typename Derived >
int Bembel::AnsatzSpace< Derived >::get_number_of_dofs ( ) const
inline

Retrieves the number of degrees of freedom of this AnsatzSpace.

Returns
The number of degrees of freedom.

Definition at line 172 of file AnsatzSpace.hpp.

◆ get_number_of_elements()

template<typename Derived >
int Bembel::AnsatzSpace< Derived >::get_number_of_elements ( ) const
inline

Retrieves the number of elements of the underlying ElementTree in the SuperSpace of this AnsatzSpace.

Returns
The number of elements.

Definition at line 154 of file AnsatzSpace.hpp.

◆ get_number_of_patches()

template<typename Derived >
int Bembel::AnsatzSpace< Derived >::get_number_of_patches ( ) const
inline

Retrieves the number of patches of the underlying Geometry.

Returns
The number of patches.

Definition at line 163 of file AnsatzSpace.hpp.

◆ get_polynomial_degree()

template<typename Derived >
int Bembel::AnsatzSpace< Derived >::get_polynomial_degree ( ) const
inline

Retrieves the polynomial degree of this AnsatzSpace.

Returns
The polynomial degree.

Definition at line 144 of file AnsatzSpace.hpp.

◆ get_refinement_level()

template<typename Derived >
int Bembel::AnsatzSpace< Derived >::get_refinement_level ( ) const
inline

Retrieves the refinement level of this AnsatzSpace.

Returns
The refinement level.

Definition at line 135 of file AnsatzSpace.hpp.

◆ get_superspace()

template<typename Derived >
const SuperSpace<Derived>& Bembel::AnsatzSpace< Derived >::get_superspace ( ) const
inline

Retrieves the reference to the SuperSpace associated with this AnsatzSpace.

Returns
A const reference to the SuperSpace.

Definition at line 121 of file AnsatzSpace.hpp.

◆ get_transformation_matrix()

template<typename Derived >
const Eigen::SparseMatrix<double>& Bembel::AnsatzSpace< Derived >::get_transformation_matrix ( ) const
inline

Retrieves the transformation matrix associated with this AnsatzSpace.

Returns
A const reference to the transformation matrix.

Definition at line 189 of file AnsatzSpace.hpp.

◆ init_AnsatzSpace()

template<typename Derived >
void Bembel::AnsatzSpace< Derived >::init_AnsatzSpace ( const Geometry geometry,
int  refinement_level,
int  polynomial_degree,
int  knot_repetition 
)
inline

Initializes the AnsatzSpace object.

This function initializes the AnsatzSpace object with the provided parameters. It sets the knot repetition, initializes the super space, utilizes the Projector and Glue class to create a transformation matrix which assembles conforming B-Splines from local Bernstein polynomials.

Parameters
geometryThe geometry object.
refinement_levelThe refinement level of the space.
polynomial_degreeThe degree of polynomials used in the space.
knot_repetitionThe number of repetitions of knots in the space.

Definition at line 102 of file AnsatzSpace.hpp.

◆ operator=()

template<typename Derived >
AnsatzSpace& Bembel::AnsatzSpace< Derived >::operator= ( AnsatzSpace< Derived >  other)
inline

Copy assignment operator.

This operator assigns the contents of another AnsatzSpace object to this one.

Parameters
otherThe AnsatzSpace object to copy from.
Returns
A reference to the updated AnsatzSpace object.

Definition at line 62 of file AnsatzSpace.hpp.


The documentation for this class was generated from the following file: