Handles the construction and management of a reference inclusion geometry and its basis. More...
#include <reference_cross_section.h>
Public Member Functions | |
ReferenceCrossSection (const ReferenceCrossSectionParameters< dim, spacedim, n_components > &par) | |
Constructs the ReferenceCrossSection from parameters. | |
const Quadrature< spacedim > & | get_global_quadrature () const |
Returns the global quadrature object in the embedding space. | |
const std::vector< Vector< double > > & | get_basis_functions () const |
Returns the list of selected basis functions. | |
const double & | shape_value (const unsigned int i, const unsigned int q, const unsigned int comp) const |
const SparseMatrix< double > & | get_mass_matrix () const |
Returns the mass matrix corresponding to selected basis functions. | |
unsigned int | n_selected_basis () const |
Returns the number of selected basis functions. | |
unsigned int | max_n_basis () const |
Returns the total number of available basis functions. | |
Quadrature< spacedim > | get_transformed_quadrature (const Point< spacedim > &new_origin, const Tensor< 1, spacedim > &new_vertical, const double scale) const |
unsigned int | n_quadrature_points () const |
void | initialize () |
Initializes the reference inclusion domain. | |
double | measure (const double scale=1.0) const |
Private Member Functions | |
void | make_grid () |
Builds the mesh for the reference inclusion domain. | |
void | setup_dofs () |
Initializes the degrees of freedom and finite element space. | |
void | compute_basis () |
Computes and stores all basis functions. | |
Private Attributes | |
const ReferenceCrossSectionParameters< dim, spacedim, n_components > & | par |
Reference to parameter configuration. | |
PolynomialsP< spacedim > | polynomials |
Polynomial space used for modal basis generation. | |
Triangulation< dim, spacedim > | triangulation |
Triangulation of the reference inclusion domain. | |
QGauss< dim > | quadrature_formula |
Quadrature formula for integration on the reference domain. | |
FESystem< dim, spacedim > | fe |
Finite element space for the inclusion. | |
MappingQ< dim, spacedim > | mapping |
DoFHandler< dim, spacedim > | dof_handler |
Degree of freedom handler for the inclusion. | |
std::vector< Vector< double > > | basis_functions |
List of all computed basis functions. | |
std::vector< Vector< double > > | selected_basis_functions |
Subset of selected basis functions. | |
SparsityPattern | sparsity_pattern |
Sparsity pattern used for assembling the mass matrix. | |
SparseMatrix< double > | mass_matrix |
Mass matrix associated with the selected basis functions. | |
Quadrature< spacedim > | global_quadrature |
Quadrature rule in the embedding space. | |
double | reference_measure |
Reference measure for the inclusion. | |
FullMatrix< double > | basis_functions_on_qpoints |
The basis functions computed on the quadrature points. | |
Handles the construction and management of a reference inclusion geometry and its basis.
Used in reduced basis immersed boundary methods, this class initializes and stores a full finite element space, computes basis functions, and provides access to quadrature rules and mass matrices.
dim | The intrinsic dimension of the inclusion. |
spacedim | The embedding dimension. |
n_components | Number of components per field variable. |
Definition at line 92 of file reference_cross_section.h.
ReferenceCrossSection< dim, spacedim, n_components >::ReferenceCrossSection | ( | const ReferenceCrossSectionParameters< dim, spacedim, n_components > & | par | ) |
Constructs the ReferenceCrossSection from parameters.
Definition at line 41 of file reference_cross_section.cc.
References dof_handler, fe, initialize(), mapping, par, polynomials, quadrature_formula, and triangulation.
|
private |
Computes and stores all basis functions.
Definition at line 131 of file reference_cross_section.cc.
Referenced by initialize().
auto ReferenceCrossSection< dim, spacedim, n_components >::get_basis_functions | ( | ) | const |
Returns the list of selected basis functions.
Definition at line 262 of file reference_cross_section.cc.
References selected_basis_functions.
auto ReferenceCrossSection< dim, spacedim, n_components >::get_global_quadrature | ( | ) | const |
Returns the global quadrature object in the embedding space.
Definition at line 252 of file reference_cross_section.cc.
References global_quadrature.
auto ReferenceCrossSection< dim, spacedim, n_components >::get_mass_matrix | ( | ) | const |
Returns the mass matrix corresponding to selected basis functions.
Definition at line 290 of file reference_cross_section.cc.
References mass_matrix.
Quadrature< spacedim > ReferenceCrossSection< dim, spacedim, n_components >::get_transformed_quadrature | ( | const Point< spacedim > & | new_origin, |
const Tensor< 1, spacedim > & | new_vertical, | ||
const double | scale ) const |
Definition at line 337 of file reference_cross_section.cc.
References Physics::VectorRelations::angle(), global_quadrature, Tensor< int rank_, int dim, typename Number >::norm(), std::pow(), Physics::Transformations::Rotations::rotation_matrix_2d(), Physics::Transformations::Rotations::rotation_matrix_3d(), and Physics::VectorRelations::signed_angle().
void ReferenceCrossSection< dim, spacedim, n_components >::initialize | ( | ) |
Initializes the reference inclusion domain.
Definition at line 407 of file reference_cross_section.cc.
References compute_basis(), dof_handler, make_grid(), setup_dofs(), and triangulation.
Referenced by ReferenceCrossSection().
|
private |
Builds the mesh for the reference inclusion domain.
Definition at line 56 of file reference_cross_section.cc.
References par.
Referenced by initialize().
unsigned int ReferenceCrossSection< dim, spacedim, n_components >::max_n_basis | ( | ) | const |
Returns the total number of available basis functions.
Definition at line 310 of file reference_cross_section.cc.
References basis_functions.
double ReferenceCrossSection< dim, spacedim, n_components >::measure | ( | const double | scale = 1.0 | ) | const |
Compute the cross section measure for the inclusion. The returned value is exact w.r.t. to the used quadrature formula, i.e., it is the integral of one on the cross section. If a scale is provided, the measure is multiplied by the scale^dim.
scale | The scaling factor for the measure. |
Definition at line 242 of file reference_cross_section.cc.
References std::pow(), and reference_measure.
unsigned int ReferenceCrossSection< dim, spacedim, n_components >::n_quadrature_points | ( | ) | const |
Definition at line 400 of file reference_cross_section.cc.
References global_quadrature.
auto ReferenceCrossSection< dim, spacedim, n_components >::n_selected_basis | ( | ) | const |
Returns the number of selected basis functions.
Definition at line 300 of file reference_cross_section.cc.
References selected_basis_functions.
|
private |
Initializes the degrees of freedom and finite element space.
Definition at line 89 of file reference_cross_section.cc.
References dof_handler, fe, DoFTools::make_sparsity_pattern(), and sparsity_pattern.
Referenced by initialize().
const double & ReferenceCrossSection< dim, spacedim, n_components >::shape_value | ( | const unsigned int | i, |
const unsigned int | q, | ||
const unsigned int | comp ) const |
the component comp of the ith selected reference basis function, at the quadrature point index q,
Definition at line 270 of file reference_cross_section.cc.
References basis_functions_on_qpoints, global_quadrature, and selected_basis_functions.
|
private |
List of all computed basis functions.
Definition at line 183 of file reference_cross_section.h.
Referenced by max_n_basis().
|
private |
The basis functions computed on the quadrature points.
Definition at line 201 of file reference_cross_section.h.
Referenced by shape_value().
|
private |
Degree of freedom handler for the inclusion.
Definition at line 180 of file reference_cross_section.h.
Referenced by initialize(), ReferenceCrossSection(), and setup_dofs().
|
private |
Finite element space for the inclusion.
Definition at line 174 of file reference_cross_section.h.
Referenced by ReferenceCrossSection(), and setup_dofs().
|
private |
Quadrature rule in the embedding space.
Definition at line 195 of file reference_cross_section.h.
Referenced by get_global_quadrature(), get_transformed_quadrature(), n_quadrature_points(), and shape_value().
|
private |
Definition at line 177 of file reference_cross_section.h.
Referenced by ReferenceCrossSection().
|
private |
Mass matrix associated with the selected basis functions.
Definition at line 192 of file reference_cross_section.h.
Referenced by get_mass_matrix().
|
private |
Reference to parameter configuration.
Definition at line 162 of file reference_cross_section.h.
Referenced by make_grid(), and ReferenceCrossSection().
|
private |
Polynomial space used for modal basis generation.
Definition at line 165 of file reference_cross_section.h.
Referenced by ReferenceCrossSection().
|
private |
Quadrature formula for integration on the reference domain.
Definition at line 171 of file reference_cross_section.h.
Referenced by ReferenceCrossSection().
|
private |
Reference measure for the inclusion.
Definition at line 198 of file reference_cross_section.h.
Referenced by measure().
|
private |
Subset of selected basis functions.
Definition at line 186 of file reference_cross_section.h.
Referenced by get_basis_functions(), n_selected_basis(), and shape_value().
|
private |
Sparsity pattern used for assembling the mass matrix.
Definition at line 189 of file reference_cross_section.h.
Referenced by setup_dofs().
|
private |
Triangulation of the reference inclusion domain.
Definition at line 168 of file reference_cross_section.h.
Referenced by initialize(), and ReferenceCrossSection().