Reduced Lagrange Multipliers
 
Loading...
Searching...
No Matches
reduced_coupling.cc
Go to the documentation of this file.
1#include "reduced_coupling.h"
2
4
7
8#include <fstream>
9
12#include "vtk_utils.h"
13
14template <int reduced_dim, int dim, int spacedim, int n_components>
17 : ParameterAcceptor("Reduced coupling/")
18{
19 this->enter_subsection("Representative domain");
20 this->add_parameter("Pre-refinement level", pre_refinement);
21 this->add_parameter("Reduced right hand side",
23 "",
24 this->prm,
26 1,
28 ";"));
29 this->leave_subsection();
30}
31
32template <int reduced_dim, int dim, int spacedim, int n_components>
34 const parallel::TriangulationBase<spacedim> &background_tria,
36 &par)
37 : TensorProductSpace<reduced_dim, dim, spacedim, n_components>(
38 par.tensor_product_space_parameters,
39#if DEAL_II_VERSION_GTE(9, 6, 0)
40 background_tria.get_communicator())
41#else
42background_tria.get_mpi_communicator()
43#endif
44 , ParticleCoupling<spacedim>(par.particle_coupling_parameters)
45 , mpi_communicator(
46#if DEAL_II_VERSION_GTE(9, 6, 0)
47 background_tria.get_communicator()
48#else
49 background_tria.get_mpi_communicator()
50#endif
51 )
52 , par(par)
53 , background_tria(&background_tria)
54 , immersed_partitioner(background_tria)
55{
56 this->set_partitioner = [&](auto &tria) {
57 tria.set_partitioner(immersed_partitioner,
59 };
60}
61
62template <int reduced_dim, int dim, int spacedim, int n_components>
63void
65 const Mapping<spacedim> &mapping)
66{
67 // Initialize the tensor product space
69
70 auto locally_owned_dofs = this->get_dof_handler().locally_owned_dofs();
71 auto locally_relevant_dofs =
73
74 coupling_constraints.clear();
75 coupling_constraints.reinit(locally_owned_dofs, locally_relevant_dofs);
76 DoFTools::make_hanging_node_constraints(this->get_dof_handler(),
77 coupling_constraints);
78 coupling_constraints.close();
79
80 // Initialize the particle coupling
82 *this->background_tria, mapping);
83
84 // Initialize the particles
85 const auto &qpoints = this->get_locally_owned_qpoints();
86 const auto &weights = this->get_locally_owned_weights();
87 auto q_index = this->insert_points(qpoints, weights);
88 this->update_local_dof_indices(q_index);
89
90 // Initialize the coupling rhs
91 typename FunctionParser<spacedim>::ConstMap constants;
92 constants["pi"] = numbers::PI;
93 constants["E"] = numbers::E;
94
95 coupling_rhs = std::make_unique<FunctionParser<spacedim>>(
96 this->get_reference_cross_section().n_selected_basis());
97
98 coupling_rhs->initialize(FunctionParser<spacedim>::default_variable_names() +
99 ",t",
101 constants,
102 true);
103
104 // This should be true. Let's double check
105 AssertDimension(coupling_rhs->n_components,
106 this->get_dof_handler().get_fe().n_components());
107
108 if (Utilities::MPI::this_mpi_process(mpi_communicator) == 0)
109 {
110 std::cout
111 << "Selected degree: "
112 << par.tensor_product_space_parameters.section.inclusion_degree
113 << ", selected basis functions: "
115 par.tensor_product_space_parameters.section.selected_coefficients)
116 << std::endl;
117 std::cout << "Reduced coupling initialized" << std::endl;
118 std::cout << "Reduced grid name: "
119 << par.tensor_product_space_parameters.reduced_grid_name
120 << std::endl;
121 }
122}
123
124template <int reduced_dim, int dim, int spacedim, int n_components>
125void
128 const DoFHandler<spacedim> &dh,
129 const AffineConstraints<double> &constraints) const
130{
131 const auto &fe = dh.get_fe();
132 std::vector<types::global_dof_index> background_dof_indices(
133 fe.n_dofs_per_cell());
134
135 auto particle = this->get_particles().begin();
136 while (particle != this->get_particles().end())
137 {
138 const auto &cell = particle->get_surrounding_cell();
139 const auto dh_cell =
140 typename DoFHandler<spacedim>::cell_iterator(*cell, &dh);
141 dh_cell->get_dof_indices(background_dof_indices);
142
143 const auto pic = this->get_particles().particles_in_cell(cell);
144 Assert(pic.begin() == particle, ExcInternalError());
145
147 for (const auto &p : pic)
148 {
149 const auto [immersed_cell_id, immersed_q, section_q] =
150 this->particle_id_to_cell_and_qpoint_indices(p.get_id());
151 // If cell id is the same, we can skip the rest of the loop. We
152 // already added these entries
153 if (immersed_cell_id != previous_cell_id)
154 {
155 const auto &immersed_dof_indices =
156 this->get_dof_indices(immersed_cell_id);
157
158 constraints.add_entries_local_to_global(background_dof_indices,
159 coupling_constraints,
160 immersed_dof_indices,
161 dsp);
162
163 previous_cell_id = immersed_cell_id;
164 }
165 }
166 particle = pic.end();
167 }
168}
169
170template <int reduced_dim, int dim, int spacedim, int n_components>
174{
175 return coupling_constraints;
176}
177
178// Explicit instantiations for ReducedCouplingParameters
183
188
189
190template struct ReducedCoupling<1, 2, 2, 1>;
191template struct ReducedCoupling<1, 2, 3, 1>;
192template struct ReducedCoupling<1, 3, 3, 1>;
193template struct ReducedCoupling<2, 3, 3, 1>;
194
195template struct ReducedCoupling<1, 2, 2, 2>;
196template struct ReducedCoupling<1, 2, 3, 3>;
197template struct ReducedCoupling<1, 3, 3, 3>;
198template struct ReducedCoupling<2, 3, 3, 3>;
void add_entries_local_to_global(const std::vector< size_type > &local_dof_indices, SparsityPatternBase &sparsity_pattern, const bool keep_constrained_entries=true, const Table< 2, bool > &dof_mask=Table< 2, bool >()) const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
std::map< std::string, double > ConstMap
static std::string default_variable_names()
ParameterAcceptor(const std::string &section_name="")
static ParameterHandler prm
void enter_subsection(const std::string &subsection)
void add_parameter(const std::string &entry, ParameterType &parameter, const std::string &documentation="", ParameterHandler &prm_=prm, const Patterns::PatternBase &pattern=*Patterns::Tools::Convert< ParameterType >::to_pattern())
static const unsigned int max_int_value
Manages the coupling of particles with a finite element background mesh.
void initialize_particle_handler(const parallel::TriangulationBase< dim > &tria_background, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping)
#define DEAL_II_VERSION_GTE(major, minor, subminor)
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
std::string to_string(const T &t)
VectorType::value_type * end(VectorType &V)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
constexpr double E
constexpr double PI
constexpr unsigned int invalid_unsigned_int
unsigned int global_cell_index
Combines tensor product space and particle coupling for reduced Lagrange multipliers.
ReducedCoupling(const parallel::TriangulationBase< spacedim > &background_tria, const ReducedCouplingParameters< reduced_dim, dim, spacedim, n_components > &par)
Constructor that initializes the ReducedCoupling object with background triangulation and parameters.
void assemble_coupling_sparsity(DynamicSparsityPattern &dsp, const DoFHandler< spacedim > &dh, const AffineConstraints< double > &constraints) const
Assemble the sparsity pattern for the coupling matrix.
const AffineConstraints< double > & get_coupling_constraints() const
Get the affine constraints associated with the coupling.
Parameter structure for configuring ReducedCoupling objects.
ReducedCouplingParameters()
Constructor that registers parameters with the ParameterAcceptor.
TensorProductSpaceParameters< reduced_dim, dim, spacedim, n_components > tensor_product_space_parameters
Parameters for the tensor product space.
std::vector< std::string > coupling_rhs_expressions
Right hand side expressions for the reduced coupling.
unsigned int pre_refinement
Number of pre-refinements to apply to the grid before distribution.