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_heating_model_interface_h
23 #define _aspect_heating_model_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 {
39  template <int dim> class SimulatorAccess;
46  namespace HeatingModel
47  {
48  using namespace dealii;
49 
57  {
68  HeatingModelOutputs (const unsigned int n_points,
69  const unsigned int n_comp);
70 
78  std::vector<double> heating_source_terms;
79 
99  std::vector<double> rates_of_temperature_change;
100 
105  std::vector<double> lhs_latent_heat_terms;
106 
112  void
113  reset ();
114  };
115 
121  template <int dim>
122  class Interface
123  {
124  public:
129  virtual ~Interface() = default;
130 
136  virtual
137  void
138  initialize ();
139 
151  virtual
152  void
153  update ();
154 
167  virtual
168  void
169  evaluate (const MaterialModel::MaterialModelInputs<dim> &material_model_inputs,
170  const MaterialModel::MaterialModelOutputs<dim> &material_model_outputs,
171  HeatingModel::HeatingModelOutputs &heating_model_outputs) const = 0;
172 
179  static
180  void
181  declare_parameters (ParameterHandler &prm);
182 
189  virtual
190  void
191  parse_parameters (ParameterHandler &prm);
192 
193 
200  virtual
201  void
202  create_additional_material_model_outputs(MaterialModel::MaterialModelOutputs<dim> &outputs) const;
203 
210  virtual
211  void
212  create_additional_material_model_inputs(MaterialModel::MaterialModelInputs<dim> &inputs) const;
213  };
214 
215 
222  template <int dim>
223  class Manager : public ::aspect::SimulatorAccess<dim>
224  {
225  public:
230  ~Manager () override;
231 
236  bool
237  adiabatic_heating_enabled() const;
238 
243  bool
244  shear_heating_enabled() const;
245 
250  static
251  void
252  declare_parameters (ParameterHandler &prm);
253 
254 
260  void
261  parse_parameters (ParameterHandler &prm);
262 
267  void
268  update ();
269 
270 
276  void
277  evaluate (const MaterialModel::MaterialModelInputs<dim> &material_model_inputs,
278  const MaterialModel::MaterialModelOutputs<dim> &material_model_outputs,
279  HeatingModel::HeatingModelOutputs &heating_model_outputs) const;
280 
286  virtual
287  void
288  create_additional_material_model_inputs_and_outputs(MaterialModel::MaterialModelInputs<dim> &material_model_inputs,
289  MaterialModel::MaterialModelOutputs<dim> &material_model_outputs) const;
290 
291 
311  static
312  void
313  register_heating_model (const std::string &name,
314  const std::string &description,
315  void (*declare_parameters_function) (ParameterHandler &),
316  std::unique_ptr<Interface<dim>> (*factory_function) ());
317 
318 
323  const std::vector<std::string> &
324  get_active_heating_model_names () const;
325 
330  const std::list<std::unique_ptr<Interface<dim>>> &
331  get_active_heating_models () const;
332 
342  template <typename HeatingModelType,
343  typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,HeatingModelType>::value>>
344  bool
345  has_matching_heating_model () const;
346 
358  template <typename HeatingModelType,
359  typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,HeatingModelType>::value>>
360  const HeatingModelType &
361  get_matching_heating_model () const;
362 
363 
373  static
374  void
375  write_plugin_graph (std::ostream &output_stream);
376 
380  DeclException1 (ExcHeatingModelNameNotFound,
381  std::string,
382  << "Could not find entry <"
383  << arg1
384  << "> among the names of registered heating model objects.");
385  private:
390  std::list<std::unique_ptr<Interface<dim>>> heating_model_objects;
391 
396  std::vector<std::string> model_names;
397  };
398 
399 
400 
401  template <int dim>
402  template <typename HeatingModelType, typename>
403  inline
404  bool
406  {
407  for (const auto &p : heating_model_objects)
408  if (Plugins::plugin_type_matches<HeatingModelType>(*p))
409  return true;
410  return false;
411  }
412 
413 
414  template <int dim>
415  template <typename HeatingModelType, typename>
416  inline
417  const HeatingModelType &
419  {
420  AssertThrow(has_matching_heating_model<HeatingModelType> (),
421  ExcMessage("You asked HeatingModel::Manager::get_heating_model() for a "
422  "heating model of type <" + boost::core::demangle(typeid(HeatingModelType).name()) + "> "
423  "that could not be found in the current model. Activate this "
424  "heating model in the input file."));
425 
426  typename std::list<std::unique_ptr<Interface<dim>>>::const_iterator heating_model;
427  for (typename std::list<std::unique_ptr<Interface<dim>>>::const_iterator
428  p = heating_model_objects.begin();
429  p != heating_model_objects.end(); ++p)
430  if (Plugins::plugin_type_matches<HeatingModelType>(*(*p)))
431  return Plugins::get_plugin_as_type<HeatingModelType>(*(*p));
432 
433  // We will never get here, because we had the Assert above. Just to avoid warnings.
434  return Plugins::get_plugin_as_type<HeatingModelType>(*(*heating_model));
435  }
436 
437 
444  template <int dim>
445  std::string
447 
448 
456 #define ASPECT_REGISTER_HEATING_MODEL(classname,name,description) \
457  template class classname<2>; \
458  template class classname<3>; \
459  namespace ASPECT_REGISTER_HEATING_MODEL_ ## classname \
460  { \
461  aspect::internal::Plugins::RegisterHelper<aspect::HeatingModel::Interface<2>,classname<2>> \
462  dummy_ ## classname ## _2d (&aspect::HeatingModel::Manager<2>::register_heating_model, \
463  name, description); \
464  aspect::internal::Plugins::RegisterHelper<aspect::HeatingModel::Interface<3>,classname<3>> \
465  dummy_ ## classname ## _3d (&aspect::HeatingModel::Manager<3>::register_heating_model, \
466  name, description); \
467  }
468  }
469 }
470 
471 
472 #endif
std::list< std::unique_ptr< Interface< dim > > > heating_model_objects
Definition: interface.h:390
std::vector< double > heating_source_terms
Definition: interface.h:78
std::vector< double > rates_of_temperature_change
Definition: interface.h:99
bool has_matching_heating_model() const
Definition: interface.h:405
void write_plugin_graph(std::ostream &output_stream)
const HeatingModelType & get_matching_heating_model() const
Definition: interface.h:418
std::string get_valid_model_names_pattern()
std::vector< double > lhs_latent_heat_terms
Definition: interface.h:105
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:42
std::vector< std::string > model_names
Definition: interface.h:396
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.")