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
67 const unsigned int n_comp);
167 get_required_properties()
const;
186 adiabatic_heating_enabled()
const;
193 shear_heating_enabled()
const;
210 parse_parameters (ParameterHandler &prm)
override;
254 register_heating_model (
const std::string &name,
255 const std::string &description,
256 void (*declare_parameters_function) (ParameterHandler &),
268 const std::vector<std::string> &
269 get_active_heating_model_names ()
const;
279 const std::list<std::unique_ptr<Interface<dim>>> &
280 get_active_heating_models ()
const;
296 template <
typename HeatingModelType,
297 typename =
typename std::enable_if_t<std::is_base_of<Interface<dim>,HeatingModelType>::value>>
300 has_matching_heating_model ()
const;
318 template <
typename HeatingModelType,
319 typename =
typename std::enable_if_t<std::is_base_of<Interface<dim>,HeatingModelType>::value>>
321 const HeatingModelType &
322 get_matching_heating_model ()
const;
343 <<
"Could not find entry <" 345 <<
"> among the names of registered heating model objects.");
351 template <
typename HeatingModelType,
typename>
356 return this->
template has_matching_active_plugin<HeatingModelType>();
361 template <
typename HeatingModelType,
typename>
363 const HeatingModelType &
366 return this->
template get_matching_active_plugin<HeatingModelType>();
388 #define ASPECT_REGISTER_HEATING_MODEL(classname,name,description) \ 389 template class classname<2>; \ 390 template class classname<3>; \ 391 namespace ASPECT_REGISTER_HEATING_MODEL_ ## classname \ 393 aspect::internal::Plugins::RegisterHelper<aspect::HeatingModel::Interface<2>,classname<2>> \ 394 dummy_ ## classname ## _2d (&aspect::HeatingModel::Manager<2>::register_heating_model, \ 395 name, description); \ 396 aspect::internal::Plugins::RegisterHelper<aspect::HeatingModel::Interface<3>,classname<3>> \ 397 dummy_ ## classname ## _3d (&aspect::HeatingModel::Manager<3>::register_heating_model, \ 398 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()
HeatingModelOutputs(const unsigned int n_points, const unsigned int n_comp)
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.")