22 #ifndef _aspect_heating_model_interface_h 23 #define _aspect_heating_model_interface_h 30 #include <deal.II/base/point.h> 31 #include <deal.II/base/parameter_handler.h> 33 #include <boost/core/demangle.hpp> 39 template <
int dim>
class SimulatorAccess;
46 namespace HeatingModel
69 const unsigned int n_comp);
180 adiabatic_heating_enabled()
const;
187 shear_heating_enabled()
const;
204 parse_parameters (ParameterHandler &prm)
override;
248 register_heating_model (
const std::string &name,
249 const std::string &description,
250 void (*declare_parameters_function) (ParameterHandler &),
262 const std::vector<std::string> &
263 get_active_heating_model_names ()
const;
273 const std::list<std::unique_ptr<Interface<dim>>> &
274 get_active_heating_models ()
const;
290 template <
typename HeatingModelType,
291 typename =
typename std::enable_if_t<std::is_base_of<Interface<dim>,HeatingModelType>::value>>
294 has_matching_heating_model ()
const;
312 template <
typename HeatingModelType,
313 typename =
typename std::enable_if_t<std::is_base_of<Interface<dim>,HeatingModelType>::value>>
315 const HeatingModelType &
316 get_matching_heating_model ()
const;
337 <<
"Could not find entry <" 339 <<
"> among the names of registered heating model objects.");
345 template <
typename HeatingModelType,
typename>
350 return this->
template has_matching_active_plugin<HeatingModelType>();
355 template <
typename HeatingModelType,
typename>
357 const HeatingModelType &
360 return this->
template get_matching_active_plugin<HeatingModelType>();
382 #define ASPECT_REGISTER_HEATING_MODEL(classname,name,description) \ 383 template class classname<2>; \ 384 template class classname<3>; \ 385 namespace ASPECT_REGISTER_HEATING_MODEL_ ## classname \ 387 aspect::internal::Plugins::RegisterHelper<aspect::HeatingModel::Interface<2>,classname<2>> \ 388 dummy_ ## classname ## _2d (&aspect::HeatingModel::Manager<2>::register_heating_model, \ 389 name, description); \ 390 aspect::internal::Plugins::RegisterHelper<aspect::HeatingModel::Interface<3>,classname<3>> \ 391 dummy_ ## classname ## _3d (&aspect::HeatingModel::Manager<3>::register_heating_model, \ 392 name, description); \
std::vector< double > heating_source_terms
std::vector< double > rates_of_temperature_change
DEAL_II_DEPRECATED bool has_matching_heating_model() const
DEAL_II_DEPRECATED const HeatingModelType & get_matching_heating_model() const
void write_plugin_graph(std::ostream &output_stream)
std::string get_valid_model_names_pattern()
std::vector< double > lhs_latent_heat_terms
void declare_parameters(ParameterHandler &prm)
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.")