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.cc
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 #include <deal.II/base/config.h>
17 
19 
21 
22 
24 
25 
26 namespace Particles
27 {
28  namespace Utilities
29  {
30  template <int dim, int spacedim, typename SparsityType, typename number>
31  void
33  const DoFHandler<dim, spacedim> & space_dh,
34  const Particles::ParticleHandler<dim, spacedim> &particle_handler,
35  SparsityType & sparsity,
36  const AffineConstraints<number> & constraints,
37  const ComponentMask & space_comps)
38  {
39  if (particle_handler.n_locally_owned_particles() == 0)
40  return; // nothing to do here
41 
42  const auto &tria = space_dh.get_triangulation();
43  const auto &fe = space_dh.get_fe();
44  const auto max_particles_per_cell =
45  particle_handler.n_global_max_particles_per_cell();
46 
47  // Take care of components
48  const ComponentMask comps =
49  (space_comps.size() == 0 ? ComponentMask(fe.n_components(), true) :
50  space_comps);
51  AssertDimension(comps.size(), fe.n_components());
52 
53  const auto n_comps = comps.n_selected_components();
54 
55  // Global to local indices
56  std::vector<unsigned int> space_gtl(fe.n_components(),
58  for (unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
59  if (comps[i])
60  space_gtl[i] = j++;
61 
62  // [TODO]: when the add_entries_local_to_global below will implement
63  // the version with the dof_mask, this should be uncommented.
64  // // Construct a dof_mask, used to distribute entries to the sparsity
65  // Table<2, bool> dof_mask(max_particles_per_cell * n_comps,
66  // fe.dofs_per_cell);
67  // dof_mask.fill(false);
68  // for (unsigned int i = 0; i < space_fe.dofs_per_cell; ++i)
69  // {
70  // const auto comp_i = space_fe.system_to_component_index(i).first;
71  // if (space_gtl[comp_i] != numbers::invalid_unsigned_int)
72  // for (unsigned int j = 0; j < max_particles_per_cell; ++j)
73  // dof_mask(i, j * n_comps + space_gtl[comp_i]) = true;
74  // }
75 
76  std::vector<types::global_dof_index> dof_indices(fe.dofs_per_cell);
77  std::vector<types::particle_index> particle_indices(
78  max_particles_per_cell * n_comps);
79 
80  auto particle = particle_handler.begin();
81  while (particle != particle_handler.end())
82  {
83  const auto &cell = particle->get_surrounding_cell(tria);
84  const auto dh_cell =
85  typename DoFHandler<dim, spacedim>::cell_iterator(*cell, &space_dh);
86  dh_cell->get_dof_indices(dof_indices);
87  const auto pic = particle_handler.particles_in_cell(cell);
88  const auto n_particles = particle_handler.n_particles_in_cell(cell);
89  particle_indices.resize(n_particles * n_comps);
90  Assert(pic.begin() == particle, ExcInternalError());
91  for (unsigned int i = 0; particle != pic.end(); ++particle, ++i)
92  {
93  const auto p_id = particle->get_id();
94  for (unsigned int j = 0; j < fe.dofs_per_cell; ++j)
95  {
96  const auto comp_j =
97  space_gtl[fe.system_to_component_index(j).first];
98  if (comp_j != numbers::invalid_unsigned_int)
99  constraints.add_entries_local_to_global(
100  {p_id * n_comps + comp_j}, {dof_indices[j]}, sparsity);
101  }
102  }
103  // [TODO]: when this works, use this:
104  // constraints.add_entries_local_to_global(particle_indices,
105  // dof_indices,
106  // sparsity,
107  // dof_mask);
108  }
109  }
110 
111 
112 
113  template <int dim, int spacedim, typename MatrixType>
114  void
116  const DoFHandler<dim, spacedim> & space_dh,
117  const Particles::ParticleHandler<dim, spacedim> &particle_handler,
118  MatrixType & matrix,
120  const ComponentMask & space_comps)
121  {
122  if (particle_handler.n_locally_owned_particles() == 0)
123  {
124  matrix.compress(VectorOperation::add);
125  return; // nothing else to do here
126  }
127 
128  AssertDimension(matrix.n(), space_dh.n_dofs());
129 
130  const auto &tria = space_dh.get_triangulation();
131  const auto &fe = space_dh.get_fe();
132  const auto max_particles_per_cell =
133  particle_handler.n_global_max_particles_per_cell();
134 
135  // Take care of components
136  const ComponentMask comps =
137  (space_comps.size() == 0 ? ComponentMask(fe.n_components(), true) :
138  space_comps);
139  AssertDimension(comps.size(), fe.n_components());
140  const auto n_comps = comps.n_selected_components();
141 
143  particle_handler.n_global_particles() * n_comps);
144 
145 
146  // Global to local indices
147  std::vector<unsigned int> space_gtl(fe.n_components(),
149  for (unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
150  if (comps[i])
151  space_gtl[i] = j++;
152 
153  // [TODO]: when the add_entries_local_to_global below will implement
154  // the version with the dof_mask, this should be uncommented.
155  // // Construct a dof_mask, used to distribute entries to the sparsity
156  // Table<2, bool> dof_mask(max_particles_per_cell * n_comps,
157  // fe.dofs_per_cell);
158  // dof_mask.fill(false);
159  // for (unsigned int i = 0; i < space_fe.dofs_per_cell; ++i)
160  // {
161  // const auto comp_i = space_fe.system_to_component_index(i).first;
162  // if (space_gtl[comp_i] != numbers::invalid_unsigned_int)
163  // for (unsigned int j = 0; j < max_particles_per_cell; ++j)
164  // dof_mask(i, j * n_comps + space_gtl[comp_i]) = true;
165  // }
166 
167  std::vector<types::global_dof_index> dof_indices(fe.dofs_per_cell);
168  std::vector<types::particle_index> particle_indices(
169  max_particles_per_cell * n_comps);
170 
172  max_particles_per_cell * n_comps, fe.dofs_per_cell);
173 
174  auto particle = particle_handler.begin();
175  while (particle != particle_handler.end())
176  {
177  const auto &cell = particle->get_surrounding_cell(tria);
178  const auto &dh_cell =
179  typename DoFHandler<dim, spacedim>::cell_iterator(*cell, &space_dh);
180  dh_cell->get_dof_indices(dof_indices);
181  const auto pic = particle_handler.particles_in_cell(cell);
182  const auto n_particles = particle_handler.n_particles_in_cell(cell);
183  particle_indices.resize(n_particles * n_comps);
184  local_matrix.reinit({n_particles * n_comps, fe.dofs_per_cell});
185  Assert(pic.begin() == particle, ExcInternalError());
186  for (unsigned int i = 0; particle != pic.end(); ++particle, ++i)
187  {
188  const auto &reference_location =
189  particle->get_reference_location();
190 
191  for (unsigned int d = 0; d < n_comps; ++d)
192  particle_indices[i * n_comps + d] =
193  particle->get_id() * n_comps + d;
194 
195  for (unsigned int j = 0; j < fe.dofs_per_cell; ++j)
196  {
197  const auto comp_j =
198  space_gtl[fe.system_to_component_index(j).first];
199  if (comp_j != numbers::invalid_unsigned_int)
200  local_matrix(i * n_comps + comp_j, j) =
201  fe.shape_value(j, reference_location);
202  }
203  }
204  constraints.distribute_local_to_global(local_matrix,
205  particle_indices,
206  dof_indices,
207  matrix);
208  }
209  matrix.compress(VectorOperation::add);
210  }
211 
212 #include "utilities.inst"
213 
214  } // namespace Utilities
215 } // namespace Particles
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
TableBase::reinit
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
VectorOperation::add
@ add
Definition: vector_operation.h:53
ComponentMask
Definition: component_mask.h:83
Particles
Definition: data_out.h:27
utilities.h
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
DoFHandler
Definition: dof_handler.h:205
ComponentMask::size
unsigned int size() const
AffineConstraints::distribute_local_to_global
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
Definition: affine_constraints.h:1859
AffineConstraints::add_entries_local_to_global
void add_entries_local_to_global(const std::vector< size_type > &local_dof_indices, SparsityPatternType &sparsity_pattern, const bool keep_constrained_entries=true, const Table< 2, bool > &dof_mask=Table< 2, bool >()) const
Definition: affine_constraints.h:2034
Particles::ParticleHandler::n_particles_in_cell
unsigned int n_particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
Definition: particle_handler.cc:691
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
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
Particles::ParticleHandler::n_global_particles
types::particle_index n_global_particles() const
Definition: particle_handler.cc:655
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
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
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
Particles::ParticleHandler
Definition: data_out.h:30
generic_linear_algebra.h
config.h
Particles::ParticleHandler::n_global_max_particles_per_cell
types::particle_index n_global_max_particles_per_cell() const
Definition: particle_handler.cc:664
FullMatrix< typename MatrixType::value_type >
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