16#ifndef dealii_particles_utilities
17#define dealii_particles_utilities
96 typename SparsityType,
97 typename number =
double>
102 SparsityType & sparsity,
150 template <
int dim,
int spacedim,
typename MatrixType>
184 typename InputVectorType,
185 typename OutputVectorType>
190 const InputVectorType & field_vector,
191 OutputVectorType & interpolated_field,
200 const auto &fe = field_dh.
get_fe();
201 auto particle = particle_handler.
begin();
205 (field_comps.size() == 0 ?
ComponentMask(fe.n_components(),
true) :
216 std::vector<unsigned int> space_gtl(fe.n_components(),
218 for (
unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
222 std::vector<types::global_dof_index> dof_indices(fe.n_dofs_per_cell());
224 while (particle != particle_handler.
end())
226 const auto &cell = particle->get_surrounding_cell();
227 const auto &dh_cell =
229 dh_cell->get_dof_indices(dof_indices);
233 for (
unsigned int i = 0; particle != pic.end(); ++particle, ++i)
236 particle->get_reference_location();
238 const auto id = particle->get_id();
240 for (
unsigned int j = 0; j < fe.n_dofs_per_cell(); ++j)
243 space_gtl[fe.system_to_component_index(j).first];
245 interpolated_field[
id * n_comps + comp_j] +=
246 fe.shape_value(j, reference_location) *
247 field_vector(dof_indices[j]);
unsigned int size() const
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
types::global_dof_index n_dofs() const
particle_iterator begin() const
particle_iterator end() const
types::particle_index n_locally_owned_particles() const
particle_iterator_range particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
types::particle_index get_next_free_particle_index() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
typename ActiveSelector::cell_iterator cell_iterator
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=ComponentMask())
void create_interpolation_sparsity_pattern(const DoFHandler< dim, spacedim > &space_dh, const Particles::ParticleHandler< dim, spacedim > &particle_handler, SparsityType &sparsity, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const ComponentMask &space_comps=ComponentMask())
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=ComponentMask())
static const unsigned int invalid_unsigned_int