ASPECT
interface.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 2023 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 
22 #ifndef _aspect_boundary_velocity_interface_h
23 #define _aspect_boundary_velocity_interface_h
24 
25 #include <aspect/plugins.h>
28 #include <aspect/utilities.h>
29 
30 #include <deal.II/base/point.h>
31 #include <deal.II/base/parameter_handler.h>
32 
33 #include <boost/core/demangle.hpp>
34 #include <typeinfo>
35 
36 
37 namespace aspect
38 {
45  namespace BoundaryVelocity
46  {
47  using namespace dealii;
48 
54  template <int dim>
55  class Interface
56  {
57  public:
62  virtual ~Interface() = default;
63 
69  virtual
70  void
71  initialize ();
72 
84  virtual
85  void
86  update ();
87 
100  virtual
101  Tensor<1,dim>
102  boundary_velocity (const types::boundary_id boundary_indicator,
103  const Point<dim> &position) const = 0;
104 
111  static
112  void
113  declare_parameters (ParameterHandler &prm);
114 
121  virtual
122  void
123  parse_parameters (ParameterHandler &prm);
124  };
125 
131  template <int dim>
132  class Manager : public ::aspect::SimulatorAccess<dim>
133  {
134  public:
139  ~Manager () override;
140 
150  virtual
151  void
152  update ();
153 
159  Tensor<1,dim>
160  boundary_velocity (const types::boundary_id boundary_indicator,
161  const Point<dim> &position) const;
162 
180  static
181  void
182  register_boundary_velocity (const std::string &name,
183  const std::string &description,
184  void (*declare_parameters_function) (ParameterHandler &),
185  std::unique_ptr<Interface<dim>> (*factory_function) ());
186 
187 
200  const std::map<types::boundary_id, std::pair<std::string,std::vector<std::string>>> &
201  get_active_boundary_velocity_names () const;
202 
212  const std::map<types::boundary_id,std::vector<std::unique_ptr<BoundaryVelocity::Interface<dim>>>> &
213  get_active_boundary_velocity_conditions () const;
214 
219  const std::set<types::boundary_id> &
220  get_zero_boundary_velocity_indicators () const;
221 
226  const std::set<types::boundary_id> &
227  get_tangential_boundary_velocity_indicators () const;
228 
233  static
234  void
235  declare_parameters (ParameterHandler &prm);
236 
242  void
243  parse_parameters (ParameterHandler &prm);
244 
254  template <typename BoundaryVelocityType,
255  typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,BoundaryVelocityType>::value>>
256  bool
257  has_matching_boundary_velocity_model () const;
258 
270  template <typename BoundaryVelocityType,
271  typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,BoundaryVelocityType>::value>>
272  const BoundaryVelocityType &
273  get_matching_boundary_velocity_model () const;
274 
284  static
285  void
286  write_plugin_graph (std::ostream &output_stream);
287 
288 
292  DeclException1 (ExcBoundaryVelocityNameNotFound,
293  std::string,
294  << "Could not find entry <"
295  << arg1
296  << "> among the names of registered boundary velocity objects.");
297  private:
302  std::map<types::boundary_id,std::vector<std::unique_ptr<BoundaryVelocity::Interface<dim>>>> boundary_velocity_objects;
303 
312  std::map<types::boundary_id, std::pair<std::string,std::vector<std::string>>> boundary_velocity_indicators;
313 
318  std::set<types::boundary_id> zero_velocity_boundary_indicators;
319 
324  std::set<types::boundary_id> tangential_velocity_boundary_indicators;
325  };
326 
327 
328 
329  template <int dim>
330  template <typename BoundaryVelocityType, typename>
331  inline
332  bool
334  {
335  for (const auto &boundary : boundary_velocity_objects)
336  for (const auto &p : boundary.second)
337  if (Plugins::plugin_type_matches<BoundaryVelocityType>(*p))
338  return true;
339  return false;
340  }
341 
342 
343  template <int dim>
344  template <typename BoundaryVelocityType, typename>
345  inline
346  const BoundaryVelocityType &
348  {
349  AssertThrow(has_matching_boundary_velocity_model<BoundaryVelocityType> (),
350  ExcMessage("You asked BoundaryVelocity::Manager::get_boundary_velocity_model() for a "
351  "boundary velocity model of type <" + boost::core::demangle(typeid(BoundaryVelocityType).name()) + "> "
352  "that could not be found in the current model. Activate this "
353  "boundary velocity model in the input file."));
354 
355  typename std::map<types::boundary_id,std::vector<std::unique_ptr<BoundaryVelocity::Interface<dim>>>>::const_iterator boundary_velocity_model;
356  for (const auto &boundary : boundary_velocity_objects)
357  for (const auto &p : boundary)
358  if (Plugins::plugin_type_matches<BoundaryVelocityType>(*p))
359  return Plugins::get_plugin_as_type<BoundaryVelocityType>(*p);
360 
361  // We will never get here, because we had the Assert above. Just to avoid warnings.
362  return Plugins::get_plugin_as_type<BoundaryVelocityType>(*(*boundary_velocity_model));
363  }
364 
365 
372  template <int dim>
373  std::string
375 
376 
384 #define ASPECT_REGISTER_BOUNDARY_VELOCITY_MODEL(classname, name, description) \
385  template class classname<2>; \
386  template class classname<3>; \
387  namespace ASPECT_REGISTER_BOUNDARY_VELOCITY_MODEL_ ## classname \
388  { \
389  aspect::internal::Plugins::RegisterHelper<aspect::BoundaryVelocity::Interface<2>,classname<2>> \
390  dummy_ ## classname ## _2d (&aspect::BoundaryVelocity::Manager<2>::register_boundary_velocity, \
391  name, description); \
392  aspect::internal::Plugins::RegisterHelper<aspect::BoundaryVelocity::Interface<3>,classname<3>> \
393  dummy_ ## classname ## _3d (&aspect::BoundaryVelocity::Manager<3>::register_boundary_velocity, \
394  name, description); \
395  }
396  }
397 }
398 
399 
400 #endif
void write_plugin_graph(std::ostream &output_stream)
std::set< types::boundary_id > tangential_velocity_boundary_indicators
Definition: interface.h:324
std::map< types::boundary_id, std::pair< std::string, std::vector< std::string > > > boundary_velocity_indicators
Definition: interface.h:312
std::string get_valid_model_names_pattern()
std::map< types::boundary_id, std::vector< std::unique_ptr< BoundaryVelocity::Interface< dim > > > > boundary_velocity_objects
Definition: interface.h:302
bool has_matching_boundary_velocity_model() const
Definition: interface.h:333
const BoundaryVelocityType & get_matching_boundary_velocity_model() const
Definition: interface.h:347
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:42
std::set< types::boundary_id > zero_velocity_boundary_indicators
Definition: interface.h:318
DeclException1(ProbabilityFunctionNegative, Point< dim >,<< "Your probability density function in the particle generator " "returned a negative probability density for the following position: "<< arg1<< ". Please check your function expression.")