Reduced Lagrange Multipliers
 
Loading...
Searching...
No Matches
ReducedCoupling< reduced_dim, dim, spacedim, n_components > Class Template Reference

Combines tensor product space and particle coupling for reduced Lagrange multipliers. More...

#include <reduced_coupling.h>

Inheritance diagram for ReducedCoupling< reduced_dim, dim, spacedim, n_components >:
TensorProductSpace< reduced_dim, dim, dim, 1 > ParticleCoupling< dim >

Public Member Functions

 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 initialize (const Mapping< spacedim > &mapping=StaticMappingQ1< spacedim >::mapping)
 Initialize the tensor product space and particle coupling.
 
void assemble_coupling_sparsity (DynamicSparsityPattern &dsp, const DoFHandler< spacedim > &dh, const AffineConstraints< double > &constraints) const
 Assemble the sparsity pattern for the coupling matrix.
 
template<typename MatrixType>
void assemble_coupling_matrix (MatrixType &coupling_matrix, const DoFHandler< spacedim > &dh, const AffineConstraints< double > &constraints) const
 Assemble the coupling matrix between background and reduced spaces.
 
template<typename MatrixType>
void assemble_coupling_mass_matrix (MatrixType &mass_matrix) const
 
template<typename VectorType>
void assemble_reduced_rhs (VectorType &reduced_rhs) const
 Assemble the right-hand side vector for the reduced space.
 
const AffineConstraints< double > & get_coupling_constraints () const
 Get the affine constraints associated with the coupling.
 
- Public Member Functions inherited from TensorProductSpace< reduced_dim, dim, dim, 1 >
 TensorProductSpace (const TensorProductSpaceParameters< reduced_dim, dim, spacedim, n_components > &par, MPI_Comm mpi_communicator=MPI_COMM_WORLD)
 
void initialize ()
 
const ReferenceCrossSection< dim - reduced_dim, spacedim, n_components > & get_reference_cross_section () const
 
void make_reduced_grid_and_properties ()
 
const DoFHandler< reduced_dim, spacedim > & get_dof_handler () const
 
const std::vector< Point< spacedim > > & get_locally_owned_qpoints () const
 
const std::vector< std::vector< double > > & get_locally_owned_weights () const
 
const std::vector< Point< spacedim > > & get_locally_owned_reduced_qpoints () const
 
const std::vector< std::vector< double > > & get_locally_owned_reduced_weights () const
 
const std::vector< std::vector< double > > & get_locally_owned_section_measure () const
 
void update_local_dof_indices (const std::map< unsigned int, IndexSet > &remote_q_point_indices)
 
const std::vector< types::global_dof_index > & get_dof_indices (const types::global_cell_index cell_index) const
 
std::tuple< unsigned int, unsigned int, unsigned int > particle_id_to_cell_and_qpoint_indices (const unsigned int qpoint_index) const
 
IndexSet locally_owned_qpoints () const
 
IndexSet locally_relevant_indices () const
 
auto get_quadrature () const -> const QGauss< reduced_dim > &
 
void compute_points_and_weights ()
 
const parallel::fullydistributed::Triangulation< reduced_dim, spacedim > & get_triangulation () const
 
double get_scaling (const unsigned int) const
 
const LinearAlgebra::distributed::Vector< double > & get_properties () const
 
const DoFHandler< reduced_dim, spacedim > & get_properties_dh () const
 
DoFHandler< reduced_dim, spacedim > & get_properties_dh ()
 
const std::vector< std::string > & get_properties_names () const
 
std::vector< std::string > & get_properties_names ()
 
- Public Member Functions inherited from ParticleCoupling< dim >
 ParticleCoupling (const ParticleCouplingParameters< dim > &par)
 Constructor.
 
void output_particles (const std::string &output_name) const
 Outputs the current state of the particles to a file.
 
void initialize_particle_handler (const parallel::TriangulationBase< dim > &tria_background, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping)
 
std::vector< std::vector< BoundingBox< dim > > > get_global_bounding_boxes () const
 
const Particles::ParticleHandler< dim > & get_particles () const
 
std::map< unsigned int, IndexSetinsert_points (const std::vector< Point< dim > > &points, const std::vector< std::vector< double > > &properties={})
 

Private Attributes

const MPI_Comm mpi_communicator
 The MPI communicator used for parallel operations.
 
const ReducedCouplingParameters< reduced_dim, dim, spacedim, n_components > & par
 Reference to the parameters used for this coupling.
 
SmartPointer< const parallel::TriangulationBase< spacedim > > background_tria
 The triangulation of the background domain.
 
AffineConstraints< double > coupling_constraints
 Affine constraints for the coupling.
 
std::unique_ptr< FunctionParser< spacedim > > coupling_rhs
 The right-hand side function for the coupling.
 
ImmersedRepartitioner< reduced_dim, spacedim > immersed_partitioner
 

Additional Inherited Members

- Public Attributes inherited from TensorProductSpace< reduced_dim, dim, dim, 1 >
std::function< void(parallel::fullydistributed::Triangulation< reduced_dim, spacedim > &)> set_partitioner
 
- Static Public Attributes inherited from TensorProductSpace< reduced_dim, dim, dim, 1 >
static constexpr int cross_section_dim
 
- Protected Member Functions inherited from TensorProductSpace< reduced_dim, dim, dim, 1 >
void setup_dofs ()
 
std::map< unsigned int, IndexSetlocal_q_point_indices_to_global_cell_indices (const std::map< unsigned int, IndexSet > &local_q_point_indices) const
 
- Protected Attributes inherited from TensorProductSpace< reduced_dim, dim, dim, 1 >
MPI_Comm mpi_communicator
 
const TensorProductSpaceParameters< reduced_dim, dim, spacedim, n_components > & par
 
ReferenceCrossSection< cross_section_dim, spacedim, n_components > reference_cross_section
 
parallel::fullydistributed::Triangulation< reduced_dim, spacedim > triangulation
 
FESystem< reduced_dim, spacedim > fe
 
QGauss< reduced_dim > quadrature_formula
 
DoFHandler< reduced_dim, spacedim > dof_handler
 
std::map< types::global_cell_index, std::vector< types::global_dof_index > > global_cell_to_dof_indices
 
std::vector< Point< spacedim > > all_qpoints
 
std::vector< std::vector< double > > all_weights
 
std::vector< Point< spacedim > > reduced_qpoints
 
std::vector< std::vector< double > > reduced_weights
 
LinearAlgebra::distributed::Vector< double > properties
 
DoFHandler< reduced_dim, spacedim > properties_dh
 
std::vector< std::string > properties_names
 
- Protected Attributes inherited from ParticleCoupling< dim >
const ParticleCouplingParameters< dim > & par
 Parameters for particle coupling.
 
MPI_Comm mpi_communicator
 Get the MPI communicator associated with the triangulation.
 
SmartPointer< const parallel::TriangulationBase< dim > > tria_background
 Smart pointer to the background triangulation.
 
SmartPointer< const Mapping< dim > > mapping
 Smart pointer to the mapping associated with the triangulation.
 
std::vector< std::vector< BoundingBox< dim > > > global_bounding_boxes
 
Particles::ParticleHandler< dim > particles
 Handler for managing particles in the simulation.
 

Detailed Description

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
class ReducedCoupling< reduced_dim, dim, spacedim, n_components >

Combines tensor product space and particle coupling for reduced Lagrange multipliers.

This class inherits from TensorProductSpace and ParticleCoupling, providing methods to initialize, assemble coupling matrices, and handle constraints for reduced coupling problems in the context of immersed or embedded finite element methods.

Template Parameters
reduced_dimThe reduced dimension of the problem.
dimThe dimension of the background domain.
spacedimThe space dimension (default: dim).
n_componentsNumber of components (default: 1).

Definition at line 119 of file reduced_coupling.h.

Constructor & Destructor Documentation

◆ ReducedCoupling()

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
ReducedCoupling< reduced_dim, dim, spacedim, n_components >::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.

Parameters
background_triaThe background domain triangulation.
parThe parameters for reduced coupling.

References background_tria, and par.

Member Function Documentation

◆ assemble_coupling_mass_matrix()

◆ assemble_coupling_matrix()

template<int reduced_dim, int dim, int spacedim, int n_components>
template<typename MatrixType>
void ReducedCoupling< reduced_dim, dim, spacedim, n_components >::assemble_coupling_matrix ( MatrixType & coupling_matrix,
const DoFHandler< spacedim > & dh,
const AffineConstraints< double > & constraints ) const
inline

◆ assemble_coupling_sparsity()

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
void ReducedCoupling< reduced_dim, dim, spacedim, n_components >::assemble_coupling_sparsity ( DynamicSparsityPattern & dsp,
const DoFHandler< spacedim > & dh,
const AffineConstraints< double > & constraints ) const

Assemble the sparsity pattern for the coupling matrix.

Parameters
dspThe dynamic sparsity pattern to fill.
dhThe DoFHandler for the background domain.
constraintsThe affine constraints to apply.

◆ assemble_reduced_rhs()

template<int reduced_dim, int dim, int spacedim, int n_components>
template<typename VectorType>
void ReducedCoupling< reduced_dim, dim, spacedim, n_components >::assemble_reduced_rhs ( VectorType & reduced_rhs) const
inline

Assemble the right-hand side vector for the reduced space.

Computes the right-hand side vector for the reduced space based on the following assumption: we would like to assemble the following $<Rg, w> = \sum_i \int_Gamma \phi_i\cdot g w_i \circ \Pi$ for a function $g$ defined on the full domain $\Gamma$.

Our assumption here is that the user does not provide $g$, but rather a function $\bar g$ defined on the reduced domain $\gamma$, such that $g = R^T \bar g$. With such a construction, we can specify, for example, a constant expression on the reduced domain by saying $g_0 = 1$, and we'd get $g = 1$ on the full domain. What we are actually assembling here is then $<R R^T \bar g, w> = <R^T \bar g, R^T w> = \sum_i \int_gamma \bar g_i \cdot w_i \int_D \phi_i^2 dD d\gamma$.

We know that $\int_D \phi_i^2 dD = \int_{\hat D} \hat{\phi}_i^2 J d\hat{D} = |\hat{D}| a^d = |D|$ is the scaling factor of the basis functions, where $a$ is the scaling of the cross-section. Notice that this means, in particular, that $||\phi_i||^2 = |D|$, and that this is the scaling that should be used in the mass matrix.

Template Parameters
VectorTypeThe vector type (e.g., Vector<double>).
Parameters
reduced_rhsThe right-hand side vector to assemble.

Definition at line 433 of file reduced_coupling.h.

References VectorOperation::add, coupling_constraints, coupling_rhs, FEValues< int dim, int spacedim >::dof_indices(), TensorProductSpace< reduced_dim, dim, dim, 1 >::get_dof_handler(), FEValues< int dim, int spacedim >::get_fe(), FEValues< int dim, int spacedim >::get_function_values(), FEValues< int dim, int spacedim >::get_JxW_values(), TensorProductSpace< reduced_dim, dim, dim, 1 >::get_quadrature(), FEValues< int dim, int spacedim >::get_quadrature_points(), TensorProductSpace< reduced_dim, dim, dim, 1 >::get_reference_cross_section(), numbers::invalid_unsigned_int, par, TensorProductSpace< reduced_dim, dim, dim, 1 >::properties, TensorProductSpace< reduced_dim, dim, dim, 1 >::properties_dh, TensorProductSpace< reduced_dim, dim, dim, 1 >::properties_names, FEValues< int dim, int spacedim >::quadrature_point_indices(), FEValues< int dim, int spacedim >::reinit(), FEValues< int dim, int spacedim >::shape_value(), and FiniteElement< int dim, int spacedim >::system_to_component_index().

◆ get_coupling_constraints()

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
const AffineConstraints< double > & ReducedCoupling< reduced_dim, dim, spacedim, n_components >::get_coupling_constraints ( ) const

Get the affine constraints associated with the coupling.

Returns
The affine constraints.

◆ initialize()

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
void ReducedCoupling< reduced_dim, dim, spacedim, n_components >::initialize ( const Mapping< spacedim > & mapping = StaticMappingQ1< spacedim >::mapping)

Initialize the tensor product space and particle coupling.

Parameters
mappingThe mapping to use (default: StaticMappingQ1).

References ParticleCoupling< dim >::mapping, and StaticMappingQ1< int dim, int spacedim >::mapping.

Member Data Documentation

◆ background_tria

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
SmartPointer<const parallel::TriangulationBase<spacedim> > ReducedCoupling< reduced_dim, dim, spacedim, n_components >::background_tria
private

The triangulation of the background domain.

Definition at line 223 of file reduced_coupling.h.

Referenced by ReducedCoupling().

◆ coupling_constraints

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
AffineConstraints<double> ReducedCoupling< reduced_dim, dim, spacedim, n_components >::coupling_constraints
private

Affine constraints for the coupling.

Definition at line 228 of file reduced_coupling.h.

Referenced by assemble_coupling_mass_matrix(), assemble_coupling_matrix(), and assemble_reduced_rhs().

◆ coupling_rhs

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
std::unique_ptr<FunctionParser<spacedim> > ReducedCoupling< reduced_dim, dim, spacedim, n_components >::coupling_rhs
private

The right-hand side function for the coupling.

Definition at line 233 of file reduced_coupling.h.

Referenced by assemble_reduced_rhs().

◆ immersed_partitioner

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
ImmersedRepartitioner<reduced_dim, spacedim> ReducedCoupling< reduced_dim, dim, spacedim, n_components >::immersed_partitioner
private

An ImmersedRepartitioner object that handles the repartitioning of the triangulation.

Definition at line 239 of file reduced_coupling.h.

◆ mpi_communicator

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
const MPI_Comm ReducedCoupling< reduced_dim, dim, spacedim, n_components >::mpi_communicator
private

The MPI communicator used for parallel operations.

Definition at line 212 of file reduced_coupling.h.

◆ par

template<int reduced_dim, int dim, int spacedim = dim, int n_components = 1>
const ReducedCouplingParameters<reduced_dim, dim, spacedim, n_components>& ReducedCoupling< reduced_dim, dim, spacedim, n_components >::par
private

Reference to the parameters used for this coupling.

Definition at line 218 of file reduced_coupling.h.

Referenced by assemble_coupling_mass_matrix(), assemble_reduced_rhs(), and ReducedCoupling().


The documentation for this class was generated from the following file: