Reduced Lagrange Multipliers
 
Loading...
Searching...
No Matches
reference_cross_section.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2024 by Luca Heltai
4//
5// This file is part of the reduced_lagrange_multipliers application, based on
6// the deal.II library.
7//
8// The reduced_lagrange_multipliers application is free software; you can use
9// it, redistribute it, and/or modify it under the terms of the Apache-2.0
10// License WITH LLVM-exception as published by the Free Software Foundation;
11// either version 3.0 of the License, or (at your option) any later version. The
12// full text of the license can be found in the file LICENSE.md at the top level
13// of the reduced_lagrange_multipliers distribution.
14//
15// ---------------------------------------------------------------------
16
23
24#ifndef rdlm_reference_inclusion
25#define rdlm_reference_inclusion
26
31
33
34#include <deal.II/fe/fe.h>
37
38#include <deal.II/grid/tria.h>
39
42#include <deal.II/lac/vector.h>
43
44#include <fstream>
45
46using namespace dealii;
47
60template <int dim, int spacedim = dim, int n_components = 1>
62{
65
67 unsigned int refinement_level = 1;
68
70 std::string inclusion_type = "hyper_ball";
71
73 unsigned int inclusion_degree = 0;
74
76 mutable std::vector<unsigned int> selected_coefficients;
77};
78
91template <int dim, int spacedim = dim, int n_components = 1>
93{
94public:
98
101 get_global_quadrature() const;
102
104 const std::vector<Vector<double>> &
105 get_basis_functions() const;
106
109 const double &
110 shape_value(const unsigned int i,
111 const unsigned int q,
112 const unsigned int comp) const;
113
116 get_mass_matrix() const;
117
119 unsigned int
120 n_selected_basis() const;
121
123 unsigned int
124 max_n_basis() const;
125
128 const Tensor<1, spacedim> &new_vertical,
129 const double scale) const;
130
131 unsigned int
132 n_quadrature_points() const;
133
135
136 void
137 initialize();
138
145 double
146 measure(const double scale = 1.0) const;
147
148private:
150 void
151 make_grid();
152
154 void
155 setup_dofs();
156
158 void
160
163
166
169
172
175
176 // Mapping space for the inclusion.
178
181
183 std::vector<Vector<double>> basis_functions;
184
186 std::vector<Vector<double>> selected_basis_functions;
187
190
193
196
199
202};
203
204#endif
ParameterAcceptor(const std::string &section_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.
void scale(const double scaling_factor, 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.