24#ifndef rdlm_reference_inclusion
25#define rdlm_reference_inclusion
60template <
int dim,
int spacedim = dim,
int n_components = 1>
91template <
int dim,
int spacedim = dim,
int n_components = 1>
104 const std::vector<Vector<double>> &
111 const unsigned int q,
112 const unsigned int comp)
const;
129 const double scale)
const;
ParameterAcceptor(const std::string §ion_name="")
void make_grid()
Builds the mesh for the reference inclusion domain.
Quadrature< spacedim > get_transformed_quadrature(const Point< spacedim > &new_origin, const Tensor< 1, spacedim > &new_vertical, const double scale) const
unsigned int max_n_basis() const
Returns the total number of available basis functions.
PolynomialsP< spacedim > polynomials
Polynomial space used for modal basis generation.
const SparseMatrix< double > & get_mass_matrix() const
Returns the mass matrix corresponding to selected basis functions.
std::vector< Vector< double > > selected_basis_functions
Subset of selected basis functions.
DoFHandler< dim, spacedim > dof_handler
Degree of freedom handler for the inclusion.
SparseMatrix< double > mass_matrix
Mass matrix associated with the selected basis functions.
void compute_basis()
Computes and stores all basis functions.
SparsityPattern sparsity_pattern
Sparsity pattern used for assembling the mass matrix.
std::vector< Vector< double > > basis_functions
List of all computed basis functions.
ReferenceCrossSection(const ReferenceCrossSectionParameters< dim, spacedim, n_components > &par)
Constructs the ReferenceCrossSection from parameters.
void initialize()
Initializes the reference inclusion domain.
MappingQ< dim, spacedim > mapping
const std::vector< Vector< double > > & get_basis_functions() const
Returns the list of selected basis functions.
double measure(const double scale=1.0) const
unsigned int n_quadrature_points() const
Triangulation< dim, spacedim > triangulation
Triangulation of the reference inclusion domain.
const Quadrature< spacedim > & get_global_quadrature() const
Returns the global quadrature object in the embedding space.
unsigned int n_selected_basis() const
Returns the number of selected basis functions.
void setup_dofs()
Initializes the degrees of freedom and finite element space.
FESystem< dim, spacedim > fe
Finite element space for the inclusion.
const ReferenceCrossSectionParameters< dim, spacedim, n_components > & par
Reference to parameter configuration.
FullMatrix< double > basis_functions_on_qpoints
The basis functions computed on the quadrature points.
Quadrature< spacedim > global_quadrature
Quadrature rule in the embedding space.
const double & shape_value(const unsigned int i, const unsigned int q, const unsigned int comp) const
QGauss< dim > quadrature_formula
Quadrature formula for integration on the reference domain.
double reference_measure
Reference measure for the inclusion.
Parameter configuration for a ReferenceCrossSection.
ReferenceCrossSectionParameters()
Constructor that registers parameters.
unsigned int refinement_level
Refinement level of the mesh.
std::vector< unsigned int > selected_coefficients
List of selected coefficient indices for reduced modeling.
std::string inclusion_type
Geometric type of inclusion ("hyper_ball", etc.).
unsigned int inclusion_degree
Degree of the polynomial basis for inclusion.