ASPECT
interface.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2013 - 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_composition_interface_h
23 #define _aspect_boundary_composition_interface_h
24 
25 #include <aspect/plugins.h>
27 #include <aspect/utilities.h>
29 
30 #include <deal.II/base/parameter_handler.h>
31 #include <deal.II/distributed/tria.h>
32 
33 #include <boost/core/demangle.hpp>
34 #include <typeinfo>
35 
36 
37 namespace aspect
38 {
45  namespace BoundaryComposition
46  {
47  using namespace dealii;
48 
54  template <int dim>
55  class Interface
56  {
57  public:
62  virtual ~Interface() = default;
63 
69  virtual void initialize ();
70 
82  virtual
83  void
84  update ();
85 
100  virtual
101  double
102  boundary_composition (const types::boundary_id boundary_indicator,
103  const Point<dim> &position,
104  const unsigned int compositional_field) const = 0;
105 
112  static
113  void
114  declare_parameters (ParameterHandler &prm);
115 
122  virtual
123  void
124  parse_parameters (ParameterHandler &prm);
125  };
126 
132  template <int dim>
133  class Manager : public ::aspect::SimulatorAccess<dim>
134  {
135  public:
140  ~Manager () override;
141 
151  virtual
152  void
153  update ();
154 
159  static
160  void
161  declare_parameters (ParameterHandler &prm);
162 
168  void
169  parse_parameters (ParameterHandler &prm);
170 
176  double
177  boundary_composition (const types::boundary_id boundary_indicator,
178  const Point<dim> &position,
179  const unsigned int compositional_field) const;
180 
198  static
199  void
200  register_boundary_composition (const std::string &name,
201  const std::string &description,
202  void (*declare_parameters_function) (ParameterHandler &),
203  std::unique_ptr<Interface<dim>> (*factory_function) ());
204 
205 
210  const std::vector<std::string> &
211  get_active_boundary_composition_names () const;
212 
217  const std::vector<std::unique_ptr<Interface<dim>>> &
218  get_active_boundary_composition_conditions () const;
219 
229  template <typename BoundaryCompositionType,
230  typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,BoundaryCompositionType>::value>>
231  bool
232  has_matching_boundary_composition_model () const;
233 
245  template <typename BoundaryCompositionType,
246  typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,BoundaryCompositionType>::value>>
247  const BoundaryCompositionType &
248  get_matching_boundary_composition_model () const;
249 
250  /*
251  * Return a set of boundary indicators for which boundary
252  * compositions are prescribed.
253  */
254  const std::set<types::boundary_id> &
255  get_fixed_composition_boundary_indicators() const;
256 
257  /*
258  * Return whether Dirichlet boundary conditions will be applied
259  * on parts of the boundaries where material flows out.
260  */
261  bool
262  allows_fixed_composition_on_outflow_boundaries() const;
263 
273  static
274  void
275  write_plugin_graph (std::ostream &output_stream);
276 
277 
281  DeclException1 (ExcBoundaryCompositionNameNotFound,
282  std::string,
283  << "Could not find entry <"
284  << arg1
285  << "> among the names of registered boundary composition objects.");
286  private:
291  std::vector<std::unique_ptr<Interface<dim>>> boundary_composition_objects;
292 
297  std::vector<std::string> model_names;
298 
305  std::vector<aspect::Utilities::Operator> model_operators;
306 
311  std::set<types::boundary_id> fixed_composition_boundary_indicators;
312 
318  };
319 
320 
321 
322  template <int dim>
323  template <typename BoundaryCompositionType, typename>
324  inline
325  bool
327  {
328  for (const auto &p : boundary_composition_objects)
329  if (Plugins::plugin_type_matches<BoundaryCompositionType>(*p))
330  return true;
331 
332  return false;
333  }
334 
335 
336 
337  template <int dim>
338  template <typename BoundaryCompositionType, typename>
339  inline
340  const BoundaryCompositionType &
342  {
343  AssertThrow(has_matching_boundary_composition_model<BoundaryCompositionType> (),
344  ExcMessage("You asked BoundaryComposition::Manager::get_boundary_composition_model() for a "
345  "boundary composition model of type <" + boost::core::demangle(typeid(BoundaryCompositionType).name()) + "> "
346  "that could not be found in the current model. Activate this "
347  "boundary composition model in the input file."));
348 
349  typename std::vector<std::unique_ptr<Interface<dim>>>::const_iterator boundary_composition_model;
350  for (typename std::vector<std::unique_ptr<Interface<dim>>>::const_iterator
351  p = boundary_composition_objects.begin();
352  p != boundary_composition_objects.end(); ++p)
353  if (Plugins::plugin_type_matches<BoundaryCompositionType>(*(*p)))
354  return Plugins::get_plugin_as_type<BoundaryCompositionType>(*(*p));
355 
356  // We will never get here, because we had the Assert above. Just to avoid warnings.
357  return Plugins::get_plugin_as_type<BoundaryCompositionType>(*(*boundary_composition_model));
358  }
359 
360 
361 
362 
369  template <int dim>
370  std::string
372 
373 
381 #define ASPECT_REGISTER_BOUNDARY_COMPOSITION_MODEL(classname, name, description) \
382  template class classname<2>; \
383  template class classname<3>; \
384  namespace ASPECT_REGISTER_BOUNDARY_COMPOSITION_MODEL_ ## classname \
385  { \
386  aspect::internal::Plugins::RegisterHelper<aspect::BoundaryComposition::Interface<2>,classname<2>> \
387  dummy_ ## classname ## _2d (&aspect::BoundaryComposition::Manager<2>::register_boundary_composition, \
388  name, description); \
389  aspect::internal::Plugins::RegisterHelper<aspect::BoundaryComposition::Interface<3>,classname<3>> \
390  dummy_ ## classname ## _3d (&aspect::BoundaryComposition::Manager<3>::register_boundary_composition, \
391  name, description); \
392  }
393  }
394 }
395 
396 
397 #endif
std::vector< aspect::Utilities::Operator > model_operators
Definition: interface.h:305
bool has_matching_boundary_composition_model() const
Definition: interface.h:326
std::vector< std::unique_ptr< Interface< dim > > > boundary_composition_objects
Definition: interface.h:291
void write_plugin_graph(std::ostream &output_stream)
const BoundaryCompositionType & get_matching_boundary_composition_model() const
Definition: interface.h:341
std::vector< std::string > model_names
Definition: interface.h:297
void declare_parameters(ParameterHandler &prm)
std::string get_valid_model_names_pattern()
Definition: compat.h:42
std::set< types::boundary_id > fixed_composition_boundary_indicators
Definition: interface.h:311
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.")