Reference documentation for deal.II version Git 082d75bebd 2019-10-16 19:44:02 +0200
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
particle_accessor.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2019 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_particle_accessor_h
17 #define dealii_particles_particle_accessor_h
18 
19 #include <deal.II/base/array_view.h>
20 
21 #include <deal.II/grid/tria.h>
22 
23 #include <deal.II/particles/particle.h>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 namespace Particles
28 {
29  // Forward declarations
30 #ifndef DOXYGEN
31  template <int, int>
32  class ParticleIterator;
33  template <int, int>
34  class ParticleHandler;
35 #endif
36 
40  template <int dim, int spacedim = dim>
42  {
43  public:
57  void
58  write_data(void *&data) const;
59 
66  void
67  set_location(const Point<spacedim> &new_location);
68 
74  const Point<spacedim> &
75  get_location() const;
76 
83  void
84  set_reference_location(const Point<dim> &new_reference_location);
85 
89  const Point<dim> &
90  get_reference_location() const;
91 
96  get_id() const;
97 
105  void
106  set_property_pool(PropertyPool &property_pool);
107 
112  bool
113  has_properties() const;
114 
121  void
122  set_properties(const std::vector<double> &new_properties);
123 
129  const ArrayView<double>
130  get_properties();
131 
138  get_properties() const;
139 
145  std::size_t
146  serialized_size_in_bytes() const;
147 
156  const Triangulation<dim, spacedim> &triangulation) const;
157 
161  template <class Archive>
162  void
163  serialize(Archive &ar, const unsigned int version);
164 
168  void
169  next();
170 
174  void
175  prev();
176 
180  bool
181  operator!=(const ParticleAccessor<dim, spacedim> &other) const;
182 
186  bool
187  operator==(const ParticleAccessor<dim, spacedim> &other) const;
188 
189  protected:
194 
201  const std::multimap<internal::LevelInd, Particle<dim, spacedim>> &map,
202  const typename std::multimap<internal::LevelInd,
203  Particle<dim, spacedim>>::iterator
204  &particle);
205 
206  private:
211  std::multimap<internal::LevelInd, Particle<dim, spacedim>> *map;
212 
217  typename std::multimap<internal::LevelInd,
219 
220  // Make ParticleIterator a friend to allow it constructing
221  // ParticleAccessors.
222  template <int, int>
223  friend class ParticleIterator;
224  template <int, int>
225  friend class ParticleHandler;
226  };
227 
228 
229 
230  template <int dim, int spacedim>
231  template <class Archive>
232  void
234  const unsigned int version)
235  {
236  return particle->second.serialize(ar, version);
237  }
238 
239 
240 } // namespace Particles
241 
242 DEAL_II_NAMESPACE_CLOSE
243 
244 #endif
std::multimap< internal::LevelInd, Particle< dim, spacedim > > * map
std::multimap< internal::LevelInd, Particle< dim, spacedim > >::iterator particle
const Point< spacedim > & get_location() const
void set_property_pool(PropertyPool &property_pool)
void set_reference_location(const Point< dim > &new_reference_location)
const ArrayView< double > get_properties()
Triangulation< dim, spacedim >::cell_iterator get_surrounding_cell(const Triangulation< dim, spacedim > &triangulation) const
void serialize(Archive &ar, const unsigned int version)
unsigned int particle_index
Definition: particle.h:64
const Point< dim > & get_reference_location() const
types::particle_index get_id() const
void write_data(void *&data) const
void set_location(const Point< spacedim > &new_location)
void set_properties(const std::vector< double > &new_properties)
std::size_t serialized_size_in_bytes() const
bool operator==(const ParticleAccessor< dim, spacedim > &other) const
bool operator!=(const ParticleAccessor< dim, spacedim > &other) const