40template <
int dim,
int spacedim,
int n_components>
46 ,
fe(
FE_Q<dim, spacedim>(
std::
max(
par.inclusion_degree, 1u)), n_components)
54template <
int dim,
int spacedim,
int n_components>
59 if (
par.inclusion_type ==
"hyper_ball")
61 if constexpr (dim == 1 && spacedim == 3)
70 else if constexpr (dim == spacedim - 1)
75 else if (par.inclusion_type ==
"hyper_cube")
80 ExcMessage(
"Unknown inclusion type. Please check the "
87template <
int dim,
int spacedim,
int n_components>
99 const auto n_global_q_points =
102 std::vector<Point<spacedim>> points;
103 std::vector<double> weights;
104 points.reserve(n_global_q_points);
105 weights.reserve(n_global_q_points);
113 for (
const auto &cell :
dof_handler.active_cell_iterators())
115 fe_values.reinit(cell);
116 points.insert(points.end(),
117 fe_values.get_quadrature_points().begin(),
118 fe_values.get_quadrature_points().end());
119 weights.insert(weights.end(),
120 fe_values.get_JxW_values().begin(),
121 fe_values.get_JxW_values().end());
122 for (
const auto &
w : fe_values.get_JxW_values())
129template <
int dim,
int spacedim,
int n_components>
137 const auto comp_i = i % n_components;
138 const auto poly_i = i / n_components;
148 std::vector<bool> is_zero(basis_functions.size(),
false);
149 for (
unsigned int i = 0; i < basis_functions.size(); ++i)
151 for (
unsigned int j = 0; j < i; ++j)
160 const auto coeff = mass_matrix.matrix_scalar_product(basis_functions[i],
172 for (
int i = is_zero.size() - 1; i >= 0; --i)
173 if (is_zero[i] ==
true)
174 basis_functions.erase(basis_functions.begin() + i);
177 if (par.selected_coefficients.empty())
179 par.selected_coefficients.resize(basis_functions.size());
180 std::iota(par.selected_coefficients.begin(),
181 par.selected_coefficients.end(),
186 for (
const auto &index : par.selected_coefficients)
193 selected_basis_functions.resize(par.selected_coefficients.size(),
197 for (
unsigned int i = 0; i < par.selected_coefficients.size(); ++i)
199 selected_basis_functions[i] =
200 basis_functions[par.selected_coefficients[i]];
205 basis_functions_on_qpoints.reinit(global_quadrature.size(),
206 selected_basis_functions.size() *
209 std::vector<Vector<double>>
values(quadrature_formula.size(),
217 for (
const auto &cell : dof_handler.active_cell_iterators())
219 fe_values.reinit(cell);
220 const unsigned int shift =
221 cell->global_active_cell_index() * quadrature_formula.size();
223 for (
unsigned int i = 0; i < selected_basis_functions.size(); ++i)
225 fe_values.get_function_values(selected_basis_functions[i], values);
226 for (
unsigned int q = 0; q < quadrature_formula.size(); ++q)
228 for (
unsigned int comp = 0; comp < n_components; ++comp)
230 basis_functions_on_qpoints(q +
shift,
231 i * n_components + comp) =
240template <
int dim,
int spacedim,
int n_components>
243 const double scale)
const
250template <
int dim,
int spacedim,
int n_components>
260template <
int dim,
int spacedim,
int n_components>
268template <
int dim,
int spacedim,
int n_components>
271 const unsigned int i,
272 const unsigned int q,
273 const unsigned int comp)
const
282 "The basis functions on quadrature points are empty."));
288template <
int dim,
int spacedim,
int n_components>
298template <
int dim,
int spacedim,
int n_components>
308template <
int dim,
int spacedim,
int n_components>
315template <
int dim,
int spacedim,
int n_components>
324 "This allows one to select a subset of the components of the "
325 "basis functions used for the representation of the data "
326 "(boundary data or forcing data). Notice that these indices are "
327 "w.r.t. to the total number of components of the problem, that "
328 "is, dimension of the space P^{inclusion_degree} number of Fourier coefficients x number of vector "
329 "components. In particular any entry of this list must be in "
330 "the set [0,inclusion_degree*n_components). ");
335template <
int dim,
int spacedim,
int n_components>
340 const double scale)
const
343 ExcMessage(
"The new vertical direction must be non-zero."));
347 vertical[spacedim - 1] = 1;
348 if constexpr (spacedim == 3)
351 const double axis_norm = axis.
norm();
352 if (axis_norm < 1
e-10)
355 for (
unsigned int i = 0; i < spacedim; ++i)
369 else if constexpr (spacedim == 2)
385 for (
unsigned int q = 0; q < original_points.size(); ++q)
387 transformed_points[q] =
388 rotation * (
scale * original_points[q]) + new_origin;
391 transformed_weights[q] = original_weights[q] * weight_scale_factor;
398template <
int dim,
int spacedim,
int n_components>
405template <
int dim,
int spacedim,
int n_components>
ParameterAcceptor(const std::string §ion_name="")
void add_parameter(const std::string &entry, ParameterType ¶meter, const std::string &documentation="", ParameterHandler &prm_=prm, const Patterns::PatternBase &pattern=*Patterns::Tools::Convert< ParameterType >::to_pattern())
numbers::NumberTraits< Number >::real_type norm() const
Handles the construction and management of a reference inclusion geometry and its basis.
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
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
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.
static ::ExceptionBase & ExcIndexRange(size_type arg1, size_type arg2, size_type arg3)
#define AssertIndexRange(index, range)
#define AssertThrow(cond, exc)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
Expression max(const Expression &a, const Expression &b)
void hyper_ball(Triangulation< dim, spacedim > &tria, const Point< spacedim > ¢er={}, const double radius=1., const bool attach_spherical_manifold_on_boundary_cells=false)
void hyper_sphere(Triangulation< spacedim - 1, spacedim > &tria, const Point< spacedim > ¢er=Point< spacedim >(), const double radius=1.)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
void flatten_triangulation(const Triangulation< dim, spacedim1 > &in_tria, Triangulation< dim, spacedim2 > &out_tria)
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Number signed_angle(const Tensor< 1, spacedim, Number > &a, const Tensor< 1, spacedim, Number > &b, const Tensor< 1, spacedim, Number > &axis)
Number angle(const Tensor< 1, spacedim, Number > &a, const Tensor< 1, spacedim, Number > &b)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
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.