ASPECT
world.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2012 - 2021 by the authors of the ASPECT code.
3 
4  This file is part of ASPECT.
5 
6  ASPECT is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2, or (at your option)
9  any later version.
10 
11  ASPECT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with ASPECT; see the file LICENSE. If not see
18  <http://www.gnu.org/licenses/>.
19  */
20 
21 #ifndef _aspect_particle_world_h
22 #define _aspect_particle_world_h
23 
24 #include <aspect/global.h>
25 
26 #include <deal.II/particles/particle.h>
27 #include <deal.II/particles/particle_accessor.h>
28 #include <deal.II/particles/particle_iterator.h>
29 #include <deal.II/particles/particle_handler.h>
30 #include <deal.II/particles/property_pool.h>
31 
32 #include <deal.II/matrix_free/fe_point_evaluation.h>
33 
38 
41 
42 #include <deal.II/base/timer.h>
43 #include <deal.II/base/array_view.h>
44 
45 #include <boost/serialization/unique_ptr.hpp>
46 
47 namespace aspect
48 {
49  template<int dim>
50  struct SimulatorSignals;
51 
52  namespace Particle
53  {
54  using namespace dealii;
55  using namespace ::Particles;
56 
57  namespace Generator
58  {
59  template <int dim>
60  class Interface;
61  }
62 
63 
64  namespace Property
65  {
66  template <int dim>
67  class Manager;
68  }
69 
70  namespace internal
71  {
77  template <int dim>
79  {
80  public:
84  virtual ~SolutionEvaluators() = default;
85 
92  virtual void
93  reinit(const typename DoFHandler<dim>::active_cell_iterator &cell,
94  const ArrayView<Point<dim>> &positions,
95  const ArrayView<double> &solution_values,
96  const UpdateFlags update_flags) = 0;
97 
104  virtual void get_solution(const unsigned int evaluation_point,
105  Vector<double> &solution) = 0;
106 
113  virtual void get_gradients(const unsigned int evaluation_point,
114  std::vector<Tensor<1, dim>> &gradients) = 0;
115 
120  virtual FEPointEvaluation<dim, dim> &
121  get_velocity_or_fluid_velocity_evaluator(const bool use_fluid_velocity) = 0;
122 
126  virtual NonMatching::MappingInfo<dim> &
127  get_mapping_info() = 0;
128  };
129 
133  template <int dim>
134  std::unique_ptr<internal::SolutionEvaluators<dim>>
135  construct_solution_evaluators(const SimulatorAccess<dim> &simulator_access,
136  const UpdateFlags update_flags);
137  }
138 
150  template <int dim>
151  class World : public SimulatorAccess<dim>
152  {
153  public:
157  World();
158 
162  ~World() override;
163 
167  void initialize();
168 
174  const Property::Manager<dim> &
175  get_property_manager() const;
176 
182  const Particles::ParticleHandler<dim> &
183  get_particle_handler() const;
184 
196  Particles::ParticleHandler<dim> &
197  get_particle_handler();
198 
217  void copy_particle_handler (const Particles::ParticleHandler<dim> &from_particle_handler,
218  Particles::ParticleHandler<dim> &to_particle_handler) const;
219 
225  void backup_particles ();
226 
233  void restore_particles ();
234 
235 
239  void setup_initial_state ();
240 
247  get_interpolator() const;
248 
252  void generate_particles();
256  void initialize_particles();
257 
265  void advance_timestep();
266 
278  types::particle_index n_global_particles() const;
279 
290  void
291  connect_to_signals(aspect::SimulatorSignals<dim> &signals);
292 
298 #if DEAL_II_VERSION_GTE(9,6,0)
299  unsigned int
300  cell_weight(const typename parallel::distributed::Triangulation<dim>::cell_iterator &cell,
301  const CellStatus status);
302 #else
303  unsigned int
304  cell_weight(const typename parallel::distributed::Triangulation<dim>::cell_iterator &cell,
305  const typename parallel::distributed::Triangulation<dim>::CellStatus status);
306 #endif
307 
311  void update_particles();
312 
316  template <class Archive>
317  void serialize (Archive &ar, const unsigned int version);
318 
322  virtual
323  void
324  save (std::ostringstream &os) const;
325 
329  virtual
330  void
331  load (std::istringstream &is);
332 
336  static
337  void
338  declare_parameters (ParameterHandler &prm);
339 
343  virtual
344  void
345  parse_parameters (ParameterHandler &prm);
346 
347  private:
349  {
350  enum Kind
351  {
352  no_balancing = 0x0,
353  remove_particles = 0x1,
354  add_particles = 0x2,
355  repartition = 0x4,
356  remove_and_add_particles = remove_particles | add_particles
357  };
358  };
359 
363  std::unique_ptr<Generator::Interface<dim>> generator;
364 
368  std::unique_ptr<Integrator::Interface<dim>> integrator;
369 
373  std::unique_ptr<Interpolator::Interface<dim>> interpolator;
374 
379  std::unique_ptr<Particles::ParticleHandler<dim>> particle_handler;
380 
387  Particles::ParticleHandler<dim> particle_handler_backup;
388 
394  std::unique_ptr<Property::Manager<dim>> property_manager;
395 
400 
414 
427 
436  unsigned int particle_weight;
437 
446 
453  std::map<types::subdomain_id, unsigned int>
454  get_subdomain_id_to_neighbor_map() const;
455 
461  void
462  apply_particle_per_cell_bounds();
463 
468  void advect_particles();
469 
473  void
474  local_initialize_particles(const typename ParticleHandler<dim>::particle_iterator &begin_particle,
475  const typename ParticleHandler<dim>::particle_iterator &end_particle);
476 
480  void
481  local_update_particles(const typename DoFHandler<dim>::active_cell_iterator &cell,
482  const typename ParticleHandler<dim>::particle_iterator &begin_particle,
483  const typename ParticleHandler<dim>::particle_iterator &end_particle);
484 
493  void
494  local_update_particles(const typename DoFHandler<dim>::active_cell_iterator &cell,
495  const typename ParticleHandler<dim>::particle_iterator &begin_particle,
496  const typename ParticleHandler<dim>::particle_iterator &end_particle,
498 
507  void
508  local_advect_particles(const typename DoFHandler<dim>::active_cell_iterator &cell,
509  const typename ParticleHandler<dim>::particle_iterator &begin_particle,
510  const typename ParticleHandler<dim>::particle_iterator &end_particle);
511 
525  void
526  local_advect_particles(const typename DoFHandler<dim>::active_cell_iterator &cell,
527  const typename ParticleHandler<dim>::particle_iterator &begin_particle,
528  const typename ParticleHandler<dim>::particle_iterator &end_particle,
530 
535  void
536  connect_particle_handler_signals(aspect::SimulatorSignals<dim> &signals,
537  ParticleHandler<dim> &particle_handler,
538  const bool connect_to_checkpoint_signals = true) const;
539  };
540 
541  /* -------------------------- inline and template functions ---------------------- */
542 
543  template <int dim>
544  template <class Archive>
545  void World<dim>::serialize (Archive &ar, const unsigned int)
546  {
547  // Note that although Boost claims to handle serialization of pointers
548  // correctly, at least for the case of unique_ptr it seems to not work.
549  // It works correctly when archiving the content of the pointer instead.
550  ar
551  &(*particle_handler)
552  ;
553  }
554  }
555 }
556 
557 #endif
void serialize(Archive &ar, const unsigned int version)
Definition: world.h:545
std::unique_ptr< internal::SolutionEvaluators< dim > > construct_solution_evaluators(const SimulatorAccess< dim > &simulator_access, const UpdateFlags update_flags)
unsigned int min_particles_per_cell
Definition: world.h:413
std::unique_ptr< Particles::ParticleHandler< dim > > particle_handler
Definition: world.h:379
unsigned int max_particles_per_cell
Definition: world.h:426
std::unique_ptr< Interpolator::Interface< dim > > interpolator
Definition: world.h:373
unsigned int particle_weight
Definition: world.h:436
std::unique_ptr< Generator::Interface< dim > > generator
Definition: world.h:363
void declare_parameters(ParameterHandler &prm)
std::unique_ptr< Integrator::Interface< dim > > integrator
Definition: world.h:368
std::unique_ptr< Property::Manager< dim > > property_manager
Definition: world.h:394
ParticleLoadBalancing::Kind particle_load_balancing
Definition: world.h:399
bool update_ghost_particles
Definition: world.h:445
Particles::ParticleHandler< dim > particle_handler_backup
Definition: world.h:387
Definition: compat.h:42