30 template <
int dim,
int spacedim,
typename number>
42 const auto &fe = space_dh.
get_fe();
43 const auto max_particles_per_cell =
55 std::vector<unsigned int> space_gtl(fe.n_components(),
57 for (
unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
75 std::vector<types::global_dof_index> dof_indices(fe.n_dofs_per_cell());
76 std::vector<types::particle_index> particle_indices(
77 max_particles_per_cell * n_comps);
79 auto particle = particle_handler.
begin();
80 while (particle != particle_handler.
end())
82 const auto &cell = particle->get_surrounding_cell();
85 dh_cell->get_dof_indices(dof_indices);
88 particle_indices.resize(n_particles * n_comps);
90 for (
unsigned int i = 0; particle != pic.end(); ++particle, ++i)
92 const auto p_id = particle->get_id();
93 for (
unsigned int j = 0; j < fe.n_dofs_per_cell(); ++j)
96 space_gtl[fe.system_to_component_index(j).first];
99 {p_id * n_comps + comp_j}, {dof_indices[j]}, sparsity);
112 template <
int dim,
int spacedim,
typename MatrixType>
129 const auto &fe = space_dh.
get_fe();
130 const auto max_particles_per_cell =
145 std::vector<unsigned int> space_gtl(fe.n_components(),
147 for (
unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
165 std::vector<types::global_dof_index> dof_indices(fe.n_dofs_per_cell());
166 std::vector<types::particle_index> particle_indices(
167 max_particles_per_cell * n_comps);
170 max_particles_per_cell * n_comps, fe.n_dofs_per_cell());
172 auto particle = particle_handler.
begin();
173 while (particle != particle_handler.
end())
175 const auto &cell = particle->get_surrounding_cell();
176 const auto &dh_cell =
178 dh_cell->get_dof_indices(dof_indices);
181 particle_indices.resize(n_particles * n_comps);
182 local_matrix.reinit({n_particles * n_comps, fe.n_dofs_per_cell()});
184 for (
unsigned int i = 0; particle != pic.end(); ++particle, ++i)
186 const auto &reference_location =
187 particle->get_reference_location();
189 for (
unsigned int d = 0; d < n_comps; ++d)
190 particle_indices[i * n_comps + d] =
191 particle->get_id() * n_comps + d;
193 for (
unsigned int j = 0; j < fe.n_dofs_per_cell(); ++j)
196 space_gtl[fe.system_to_component_index(j).first];
198 local_matrix(i * n_comps + comp_j, j) =
199 fe.shape_value(j, reference_location);
210#include "utilities.inst"
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
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
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 types::fe_index index=0) const
types::global_dof_index n_dofs() const
types::particle_index n_global_particles() const
particle_iterator begin() const
particle_iterator end() const
types::particle_index n_locally_owned_particles() const
types::particle_index n_particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
particle_iterator_range particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
types::particle_index n_global_max_particles_per_cell() 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 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())
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=ComponentMask())
static const unsigned int invalid_unsigned_int