![]() |
ASPECT
|
Static Public Member Functions | |
static void | declare_parameters (ParameterHandler &prm) |
![]() | |
static void | declare_parameters (ParameterHandler &prm) |
![]() | |
static void | get_composition_values_at_q_point (const std::vector< std::vector< double >> &composition_values, const unsigned int q, std::vector< double > &composition_values_at_q_point) |
Private Types | |
enum | StrainRateDistribution { constant, latitudinal_variation } |
Private Attributes | |
double | tidal_frequency |
double | elastic_shear_modulus |
double | constant_tidal_strain_rate |
enum aspect::HeatingModel::TidalHeating::StrainRateDistribution | strain_rate_distribution |
double | maximum_tidal_strain_rate |
double | minimum_tidal_strain_rate |
A class that implements tidal heating.
Definition at line 39 of file tidal_heating.h.
|
private |
Specify the distribution of time-averaged tidal strain rate.
Enumerator | |
---|---|
constant | |
latitudinal_variation |
Definition at line 100 of file tidal_heating.h.
|
overridevirtual |
Initialize function, which sets the start time and start timestep of tidal heating.
Reimplemented from aspect::Plugins::InterfaceBase.
|
overridevirtual |
Return the tidal heating terms.
Implements aspect::HeatingModel::Interface< dim >.
|
overridevirtual |
Specify which material model outputs the heating model requires for computing the heating terms.
Reimplemented from aspect::HeatingModel::Interface< dim >.
|
static |
Declare the parameters this class takes through input files.
|
overridevirtual |
Read the parameters this class declares from the parameter file.
Reimplemented from aspect::Plugins::InterfaceBase.
|
private |
Parameters used for tidal heating (H), which is defined using the following Equation is from Tobie et al. (2003) (https://doi.org/10.1029/2003JE002099) H = 2*(viscosity)*(time-averaged tidal strain rate)^2/(1+((viscosity)*(tidal frequency)/(elastic shear modulus))^2)) viscosity (Pa s) = viscosity calculated by the selected material in ASPECT time-averaged strain rate = constant_tidal_strain_rate tidal frequency = tidal_frequency elastic shear modulus = elastic_shear_modulus
strain_rate_distribution lets selection for distribution of tidal strain rate. If 'constant', the tidal strain rate is fixed to 'Constant tidal strain rate'. If 'latitudinal variation', 'Maximum tidal strain rate' and 'Minimum tidal strain rate' are used. Then, the equation follows as (maximum_tidal_strain_rate - minimum_tidal_strain_rate)*cos(theta/2)+(maximum_tidal_strain_rate+minimum_tidal_strain_rate)/2. This is shown in Fig.3 of Nimmo et al. (2007) (https://doi.org/10.1016/j.icarus.2007.04.021).
Definition at line 93 of file tidal_heating.h.
|
private |
Definition at line 94 of file tidal_heating.h.
|
private |
Definition at line 95 of file tidal_heating.h.
|
private |
|
private |
Definition at line 106 of file tidal_heating.h.
|
private |
Definition at line 107 of file tidal_heating.h.