15template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
32template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
38 par.tensor_product_space_parameters,
40 background_tria.get_communicator())
42background_tria.get_mpi_communicator()
47 background_tria.get_communicator()
49 background_tria.get_mpi_communicator()
53 , background_tria(&background_tria)
54 , immersed_partitioner(background_tria)
56 this->preprocess_serial_triangulation =
62 this->set_partitioner = [&](
auto &tria) {
63 tria.set_partitioner(immersed_partitioner,
68template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
76 auto locally_owned_dofs = this->get_dof_handler().locally_owned_dofs();
77 auto locally_relevant_dofs =
80 coupling_constraints.clear();
81 coupling_constraints.reinit(locally_owned_dofs, locally_relevant_dofs);
83 coupling_constraints);
84 coupling_constraints.close();
88 *this->background_tria, mapping);
91 const auto &qpoints = this->get_locally_owned_qpoints();
92 const auto &weights = this->get_locally_owned_weights();
93 auto q_index = this->insert_points(qpoints, weights);
94 this->update_local_dof_indices(q_index);
101 coupling_rhs = std::make_unique<FunctionParser<spacedim>>(
102 this->get_reference_cross_section().n_selected_basis());
112 this->get_dof_handler().get_fe().n_components());
117 <<
"Selected degree: "
119 <<
", selected basis functions: "
123 std::cout <<
"Reduced coupling initialized" << std::endl;
124 std::cout <<
"Reduced grid name: "
130template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
137 const auto &fe = dh.
get_fe();
138 std::vector<types::global_dof_index> background_dof_indices(
139 fe.n_dofs_per_cell());
141 auto particle = this->get_particles().begin();
142 while (particle != this->get_particles().
end())
144 const auto &cell = particle->get_surrounding_cell();
147 dh_cell->get_dof_indices(background_dof_indices);
149 const auto pic = this->get_particles().particles_in_cell(cell);
150 Assert(pic.begin() == particle, ExcInternalError());
153 for (
const auto &p : pic)
155 const auto [immersed_cell_id, immersed_q, section_q] =
156 this->particle_id_to_cell_and_qpoint_indices(p.get_id());
159 if (immersed_cell_id != previous_cell_id)
161 const auto &immersed_dof_indices =
162 this->get_dof_indices(immersed_cell_id);
165 coupling_constraints,
166 immersed_dof_indices,
169 previous_cell_id = immersed_cell_id;
172 particle = pic.end();
176template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
181 return coupling_constraints;
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 §ion_name="")
static ParameterHandler prm
void enter_subsection(const std::string &subsection)
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())
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)
VectorType::value_type * end(VectorType &V)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
constexpr unsigned int invalid_unsigned_int
unsigned int global_cell_index
Combines tensor product space and particle coupling for reduced Lagrange multipliers.
void assemble_coupling_sparsity(DynamicSparsityPattern &dsp, const DoFHandler< spacedim > &dh, const AffineConstraints< double > &constraints) const
Assemble the sparsity pattern for the coupling matrix.
ReducedCoupling(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.
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.
RefinementParameters refinement_parameters
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.
void adjust_grids(Triangulation< spacedim, spacedim > &space_triangulation, Triangulation< reduced_dim, spacedim > &embedded_triangulation, const RefinementParameters ¶meters=RefinementParameters())