Reference documentation for deal.II version 9.2.0
\(\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\}}\)
utilities.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_particles_utilities
17 #define dealii_particles_utilities
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/index_set.h>
22 #include <deal.II/base/point.h>
24 
26 
28 #include <deal.II/fe/fe.h>
29 #include <deal.II/fe/mapping_q1.h>
30 
32 
34 
36 
37 
39 
40 
41 
42 namespace Particles
43 {
48  namespace Utilities
49  {
96  template <int dim,
97  int spacedim,
98  typename SparsityType,
99  typename number = double>
100  void
102  const DoFHandler<dim, spacedim> & space_dh,
103  const Particles::ParticleHandler<dim, spacedim> &particle_handler,
104  SparsityType & sparsity,
105  const AffineConstraints<number> & constraints =
107  const ComponentMask &space_comps = ComponentMask());
108 
154  template <int dim, int spacedim, typename MatrixType>
155  void
157  const DoFHandler<dim, spacedim> & space_dh,
158  const Particles::ParticleHandler<dim, spacedim> &particle_handler,
159  MatrixType & matrix,
162  const ComponentMask &space_comps = ComponentMask());
163 
188  template <int dim,
189  int spacedim,
190  typename InputVectorType,
191  typename OutputVectorType>
192  void
194  const DoFHandler<dim, spacedim> & field_dh,
195  const Particles::ParticleHandler<dim, spacedim> &particle_handler,
196  const InputVectorType & field_vector,
197  OutputVectorType & interpolated_field,
198  const ComponentMask &field_comps = ComponentMask())
199  {
200  if (particle_handler.n_locally_owned_particles() == 0)
201  {
202  interpolated_field.compress(VectorOperation::add);
203  return; // nothing else to do here
204  }
205 
206  const auto &tria = field_dh.get_triangulation();
207  const auto &fe = field_dh.get_fe();
208  auto particle = particle_handler.begin();
209 
210  // Take care of components
211  const ComponentMask comps =
212  (field_comps.size() == 0 ? ComponentMask(fe.n_components(), true) :
213  field_comps);
214  AssertDimension(comps.size(), fe.n_components());
215  const auto n_comps = comps.n_selected_components();
216 
217  AssertDimension(field_vector.size(), field_dh.n_dofs());
218  AssertDimension(interpolated_field.size(),
219  particle_handler.get_next_free_particle_index() *
220  n_comps);
221 
222  // Global to local indices
223  std::vector<unsigned int> space_gtl(fe.n_components(),
225  for (unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
226  if (comps[i])
227  space_gtl[i] = j++;
228 
229  std::vector<types::global_dof_index> dof_indices(fe.dofs_per_cell);
230 
231  while (particle != particle_handler.end())
232  {
233  const auto &cell = particle->get_surrounding_cell(tria);
234  const auto &dh_cell =
235  typename DoFHandler<dim, spacedim>::cell_iterator(*cell, &field_dh);
236  dh_cell->get_dof_indices(dof_indices);
237  const auto pic = particle_handler.particles_in_cell(cell);
238 
239  Assert(pic.begin() == particle, ExcInternalError());
240  for (unsigned int i = 0; particle != pic.end(); ++particle, ++i)
241  {
242  const auto &reference_location =
243  particle->get_reference_location();
244 
245  const auto id = particle->get_id();
246 
247  for (unsigned int j = 0; j < fe.dofs_per_cell; ++j)
248  {
249  const auto comp_j =
250  space_gtl[fe.system_to_component_index(j).first];
251  if (comp_j != numbers::invalid_unsigned_int)
252  interpolated_field[id * n_comps + comp_j] +=
253  fe.shape_value(j, reference_location) *
254  field_vector(dof_indices[j]);
255  }
256  }
257  }
258  interpolated_field.compress(VectorOperation::add);
259  }
260 
261  } // namespace Utilities
262 } // namespace Particles
264 
265 #endif
Particles::ParticleHandler::particles_in_cell
particle_iterator_range particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
Definition: particle_handler.cc:313
DoFHandler::cell_iterator
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:357
ComponentMask::n_selected_components
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
particle_handler.h
Particles::ParticleHandler::get_next_free_particle_index
types::particle_index get_next_free_particle_index() const
Definition: particle_handler.cc:712
mapping_q1.h
VectorOperation::add
@ add
Definition: vector_operation.h:53
ComponentMask
Definition: component_mask.h:83
Particles
Definition: data_out.h:27
DoFHandler
Definition: dof_handler.h:205
ComponentMask::size
unsigned int size() const
Particles::Utilities::create_interpolation_sparsity_pattern
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())
Definition: utilities.cc:32
DoFHandler::get_triangulation
const Triangulation< dim, spacedim > & get_triangulation() const
index_set.h
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
LAPACKSupport::matrix
@ matrix
Contents is actually a matrix.
Definition: lapack_support.h:60
fe.h
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
component_mask.h
Particles::ParticleHandler::end
particle_iterator end() const
Definition: particle_handler.cc:247
Particles::ParticleHandler::n_locally_owned_particles
types::particle_index n_locally_owned_particles() const
Definition: particle_handler.cc:673
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
AffineConstraints
Definition: affine_constraints.h:180
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
Particles::Utilities::create_interpolation_matrix
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())
Definition: utilities.cc:115
dof_handler.h
affine_constraints.h
grid_tools_cache.h
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
Particles::ParticleHandler
Definition: data_out.h:30
Particles::Utilities::interpolate_field_on_particles
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())
Definition: utilities.h:193
quadrature.h
config.h
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Particles::ParticleHandler::begin
particle_iterator begin() const
Definition: particle_handler.cc:229
DoFHandler::get_fe
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
Utilities
Definition: cuda.h:31
DoFHandler::n_dofs
types::global_dof_index n_dofs() const
point.h