ASPECT
interface.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2015 - 2022 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_interpolator_interface_h
22 #define _aspect_particle_interpolator_interface_h
23 
24 #include <aspect/plugins.h>
25 #include <aspect/global.h>
26 
27 #include <deal.II/particles/particle.h>
28 #include <deal.II/particles/particle_handler.h>
29 #include <deal.II/base/point.h>
30 #include <deal.II/base/parameter_handler.h>
31 #include <deal.II/distributed/tria.h>
32 #include <deal.II/fe/component_mask.h>
33 
34 namespace aspect
35 {
36  namespace Particle
37  {
38  namespace Interpolator
39  {
40  using namespace dealii;
41  using namespace ::Particles;
42 
49  template <int dim>
50  class Interface
51  {
52  public:
57  virtual ~Interface () = default;
58 
78  virtual
79  std::vector<std::vector<double>>
80  properties_at_points(const ParticleHandler<dim> &particle_handler,
81  const std::vector<Point<dim>> &positions,
82  const typename parallel::distributed::Triangulation<dim>::active_cell_iterator &cell = typename parallel::distributed::Triangulation<dim>::active_cell_iterator()) const;
83 
108  virtual
109  std::vector<std::vector<double>>
110  properties_at_points(const ParticleHandler<dim> &particle_handler,
111  const std::vector<Point<dim>> &positions,
112  const ComponentMask &selected_properties,
113  const typename parallel::distributed::Triangulation<dim>::active_cell_iterator &cell = typename parallel::distributed::Triangulation<dim>::active_cell_iterator()) const = 0;
114 
121  static
122  void
123  declare_parameters (ParameterHandler &prm);
124 
131  virtual
132  void
133  parse_parameters (ParameterHandler &prm);
134  };
135 
136 
141  std::string
143 
144 
160  template <int dim>
161  void
162  register_particle_interpolator (const std::string &name,
163  const std::string &description,
164  void (*declare_parameters_function) (ParameterHandler &),
165  std::unique_ptr<Interface<dim>> (*factory_function) ());
166 
177  template <int dim>
178  std::unique_ptr<Interface<dim>>
179  create_particle_interpolator (ParameterHandler &prm);
180 
181 
187  template <int dim>
188  void
189  declare_parameters (ParameterHandler &prm);
190 
191 
201  template <int dim>
202  void
203  write_plugin_graph (std::ostream &output_stream);
204 
212 #define ASPECT_REGISTER_PARTICLE_INTERPOLATOR(classname, name, description) \
213  template class classname<2>; \
214  template class classname<3>; \
215  namespace ASPECT_REGISTER_PARTICLE_INTERPOLATOR_ ## classname \
216  { \
217  aspect::internal::Plugins::RegisterHelper<aspect::Particle::Interpolator::Interface<2>,classname<2>> \
218  dummy_ ## classname ## _2d (&aspect::Particle::Interpolator::register_particle_interpolator<2>, \
219  name, description); \
220  aspect::internal::Plugins::RegisterHelper<aspect::Particle::Interpolator::Interface<3>,classname<3>> \
221  dummy_ ## classname ## _3d (&aspect::Particle::Interpolator::register_particle_interpolator<3>, \
222  name, description); \
223  }
224  }
225  }
226 }
227 
228 
229 #endif
std::unique_ptr< Interface< dim > > create_particle_interpolator(ParameterHandler &prm)
StructuredDataLookup< dim > DEAL_II_DEPRECATED
void write_plugin_graph(std::ostream &output_stream)
std::string interpolator_object_names()
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:42
void register_particle_interpolator(const std::string &name, const std::string &description, void(*declare_parameters_function)(ParameterHandler &), std::unique_ptr< Interface< dim >>(*factory_function)())