Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Functions
Particles::Utilities Namespace Reference

Functions

template<int dim, int spacedim, typename number = double>
void create_interpolation_sparsity_pattern (const DoFHandler< dim, spacedim > &space_dh, const Particles::ParticleHandler< dim, spacedim > &particle_handler, SparsityPatternBase &sparsity, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const ComponentMask &space_comps={})
 
template<int dim, int spacedim, typename MatrixType >
void create_interpolation_matrix (const DoFHandler< dim, spacedim > &space_dh, const Particles::ParticleHandler< dim, spacedim > &particle_handler, MatrixType &matrix, const AffineConstraints< typename MatrixType::value_type > &constraints=AffineConstraints< typename MatrixType::value_type >(), const ComponentMask &space_comps={})
 
template<int dim, int spacedim, typename InputVectorType , typename OutputVectorType >
void interpolate_field_on_particles (const DoFHandler< dim, spacedim > &field_dh, const Particles::ParticleHandler< dim, spacedim > &particle_handler, const InputVectorType &field_vector, OutputVectorType &interpolated_field, const ComponentMask &field_comps={})
 

Detailed Description

A namespace for functions offering tools to handle ParticleHandler objects and their coupling with DoFHandler objects.

Function Documentation

◆ create_interpolation_sparsity_pattern()

template<int dim, int spacedim, typename number = double>
void Particles::Utilities::create_interpolation_sparsity_pattern ( const DoFHandler< dim, spacedim > &  space_dh,
const Particles::ParticleHandler< dim, spacedim > &  particle_handler,
SparsityPatternBase sparsity,
const AffineConstraints< number > &  constraints = AffineConstraints<number>(),
const ComponentMask space_comps = {} 
)

Create an interpolation sparsity pattern for particles.

Given a triangulation representing the domain \(\Omega\), a particle handler of particles in \(\Omega\), and a scalar finite element space \(V(\Omega) = \text{span}\{v_j\}_{j=0}^n\), compute the sparsity pattern that would be necessary to assemble the matrix

\[ M_{i,j} \dealcoloneq v_j(x_i) , \]

where \(V(\Omega)\) is the finite element space associated with the space_dh, and the index i is given by the particle id whose position is x_i.

In the case of vector valued finite element spaces, the components on which interpolation must be performed can be selected using a component mask. Only primitive finite element spaces are supported.

When selecting more than one component, the resulting sparsity will have dimension equal to particle_handler.n_global_particles() * mask.n_selected_components() times space_dh.n_dofs(), and the corresponding matrix entries are given by

\[ M_{(i*n_comps+k),j} \dealcoloneq v_j(x_i) \cdot e_{comp_j}, \]

where comp_j is the only non zero component of the vector valued basis function v_j (equal to fe.system_to_component_index(j).first), k corresponds to its index within the selected components of the mask, and \(e_{comp_j}\) is the unit vector in the direction comp_j.

The sparsity is filled by locating the position of the particle with index i within the particle handler with respect to the embedding triangulation \(\Omega\), and coupling it with all the local degrees of freedom specified in the component mask space_comps, following the ordering in which they are selected in the mask space_comps.

If a particle does not fall within \(\Omega\), it is ignored, and the corresponding rows of the sparsity will be empty.

Constraints of the form supported by the AffineConstraints class may be supplied with the constraints argument. The method AffineConstraints::add_entries_local_to_global() is used to fill the final sparsity pattern.

Definition at line 31 of file utilities.cc.

◆ create_interpolation_matrix()

template<int dim, int spacedim, typename MatrixType >
void Particles::Utilities::create_interpolation_matrix ( const DoFHandler< dim, spacedim > &  space_dh,
const Particles::ParticleHandler< dim, spacedim > &  particle_handler,
MatrixType &  matrix,
const AffineConstraints< typename MatrixType::value_type > &  constraints = AffineConstraints<typename MatrixType::value_type>(),
const ComponentMask space_comps = {} 
)

Create an interpolation matrix for particles.

Given a triangulation representing the domains \(\Omega\), a particle handler of particles in \(\Omega\), and a scalar finite element space \(V(\Omega) = \text{span}\{v_j\}_{j=0}^n\), compute the matrix

\[ M_{ij} \dealcoloneq v_j(x_i) , \]

where \(V(\Omega)\) is the finite element space associated with the space_dh, and the index i is given by the particle id whose position is x_i.

In the case of vector valued finite element spaces, the components on which interpolation must be performed can be selected using a component mask. Only primitive finite element spaces are supported.

When selecting more than one component, the resulting sparsity will have dimension equal to particle_handler.n_global_particles() * mask.n_selected_components() times space_dh.n_dofs(), and the corresponding matrix entries are given by

\[ M_{(i*n_comps+k),j} \dealcoloneq v_j(x_i) \cdot e_{comp_j}, \]

where comp_j is the only non zero component of the vector valued basis function v_j (equal to fe.system_to_component_index(j).first), k corresponds to its index within the selected components of the mask, and \(e_{comp_j}\) is the unit vector in the direction comp_j.

The matrix is filled by locating the position of the particle with index i within the particle handler with respect to the embedding triangulation \(\Omega\), and coupling it with all the local degrees of freedom specified in the component mask space_comps, following the ordering in which they are selected in the mask space_comps.

If a particle does not fall within \(\Omega\), it is ignored, and the corresponding rows of the matrix will be zero.

Constraints of the form supported by the AffineConstraints class may be supplied with the constraints argument. The method AffineConstraints::distribute_local_to_global() is used to distribute the entries of the matrix to respect the given constraints.

Definition at line 113 of file utilities.cc.

◆ interpolate_field_on_particles()

template<int dim, int spacedim, typename InputVectorType , typename OutputVectorType >
void Particles::Utilities::interpolate_field_on_particles ( const DoFHandler< dim, spacedim > &  field_dh,
const Particles::ParticleHandler< dim, spacedim > &  particle_handler,
const InputVectorType &  field_vector,
OutputVectorType &  interpolated_field,
const ComponentMask field_comps = {} 
)

Given a DoFHandler and a ParticleHandler, interpolate a vector field at the position of the particles. The result is stored in an output vector whose size corresponds to the number of locally owned particles * number of active components

Parameters
[in]field_dhThe DOF Handler which was used to generate the field vector that is to be interpolated.
[in]particle_handlerThe particle handler whose particle serve as the interpolation points.
[in]field_vectorThe vector of the field to be interpolated. This vector must be coherent with the dof_handler provided
[in,out]interpolated_fieldThe interpolated value of the field at the position of the particles. The size of the vector must be n_locally_owned_particles times the n_components
[in]field_compsAn optional component mask that decides which subset of the vector fields are interpolated

Definition at line 183 of file utilities.h.