ASPECT
Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
aspect::BoundaryTemperature::DynamicCore< dim > Class Template Reference
Inheritance diagram for aspect::BoundaryTemperature::DynamicCore< dim >:
Inheritance graph
[legend]

Classes

struct  str_data_OES
 

Public Member Functions

 DynamicCore ()
 
double boundary_temperature (const types::boundary_id boundary_indicator, const Point< dim > &location) const override
 
double minimal_temperature (const std::set< types::boundary_id > &fixed_boundary_ids) const override
 
double maximal_temperature (const std::set< types::boundary_id > &fixed_boundary_ids) const override
 
void parse_parameters (ParameterHandler &prm) override
 
void update () override
 
const internal::CoreDataget_core_data () const
 
bool is_OES_used () const
 
- Public Member Functions inherited from aspect::BoundaryTemperature::Interface< dim >
virtual ~Interface ()=default
 
virtual void initialize ()
 
- Public Member Functions inherited from aspect::SimulatorAccess< dim >
 SimulatorAccess ()
 
 SimulatorAccess (const Simulator< dim > &simulator_object)
 
virtual ~SimulatorAccess ()=default
 
virtual void initialize_simulator (const Simulator< dim > &simulator_object)
 
const Introspection< dim > & introspection () const
 
const Simulator< dim > & get_simulator () const
 
const Parameters< dim > & get_parameters () const
 
SimulatorSignals< dim > & get_signals () const
 
MPI_Comm get_mpi_communicator () const
 
TimerOutput & get_computing_timer () const
 
const ConditionalOStream & get_pcout () const
 
double get_time () const
 
double get_timestep () const
 
double get_old_timestep () const
 
unsigned int get_timestep_number () const
 
const TimeStepping::Manager< dim > & get_timestepping_manager () const
 
unsigned int get_nonlinear_iteration () const
 
const parallel::distributed::Triangulation< dim > & get_triangulation () const
 
double get_volume () const
 
const Mapping< dim > & get_mapping () const
 
std::string get_output_directory () const
 
bool include_adiabatic_heating () const
 
bool include_latent_heat () const
 
bool include_melt_transport () const
 
int get_stokes_velocity_degree () const
 
double get_adiabatic_surface_temperature () const
 
double get_surface_pressure () const
 
bool convert_output_to_years () const
 
unsigned int get_pre_refinement_step () const
 
unsigned int n_compositional_fields () const
 
double get_end_time () const
 
void get_refinement_criteria (Vector< float > &estimated_error_per_cell) const
 
void get_artificial_viscosity (Vector< float > &viscosity_per_cell, const bool skip_interior_cells=false) const
 
void get_artificial_viscosity_composition (Vector< float > &viscosity_per_cell, const unsigned int compositional_variable) const
 
const LinearAlgebra::BlockVectorget_current_linearization_point () const
 
const LinearAlgebra::BlockVectorget_solution () const
 
const LinearAlgebra::BlockVectorget_old_solution () const
 
const LinearAlgebra::BlockVectorget_old_old_solution () const
 
const LinearAlgebra::BlockVectorget_reaction_vector () const
 
const LinearAlgebra::BlockVectorget_mesh_velocity () const
 
const DoFHandler< dim > & get_dof_handler () const
 
const FiniteElement< dim > & get_fe () const
 
const LinearAlgebra::BlockSparseMatrixget_system_matrix () const
 
const LinearAlgebra::BlockSparseMatrixget_system_preconditioner_matrix () const
 
const MaterialModel::Interface< dim > & get_material_model () const
 
const GravityModel::Interface< dim > & get_gravity_model () const
 
const InitialTopographyModel::Interface< dim > & get_initial_topography_model () const
 
const GeometryModel::Interface< dim > & get_geometry_model () const
 
const AdiabaticConditions::Interface< dim > & get_adiabatic_conditions () const
 
bool has_boundary_temperature () const
 
const BoundaryTemperature::Manager< dim > & get_boundary_temperature_manager () const
 
const BoundaryHeatFlux::Interface< dim > & get_boundary_heat_flux () const
 
bool has_boundary_composition () const
 
const BoundaryComposition::Manager< dim > & get_boundary_composition_manager () const
 
const BoundaryTraction::Manager< dim > & get_boundary_traction_manager () const
 
std::shared_ptr< const InitialTemperature::Manager< dim > > get_initial_temperature_manager_pointer () const
 
const InitialTemperature::Manager< dim > & get_initial_temperature_manager () const
 
std::shared_ptr< const InitialComposition::Manager< dim > > get_initial_composition_manager_pointer () const
 
const InitialComposition::Manager< dim > & get_initial_composition_manager () const
 
const std::set< types::boundary_id > & get_fixed_temperature_boundary_indicators () const
 
const std::set< types::boundary_id > & get_fixed_heat_flux_boundary_indicators () const
 
const std::set< types::boundary_id > & get_fixed_composition_boundary_indicators () const
 
const std::set< types::boundary_id > & get_mesh_deformation_boundary_indicators () const
 
const BoundaryVelocity::Manager< dim > & get_boundary_velocity_manager () const
 
const HeatingModel::Manager< dim > & get_heating_model_manager () const
 
const MeshRefinement::Manager< dim > & get_mesh_refinement_manager () const
 
const MeltHandler< dim > & get_melt_handler () const
 
const VolumeOfFluidHandler< dim > & get_volume_of_fluid_handler () const
 
const NewtonHandler< dim > & get_newton_handler () const
 
const MeshDeformation::MeshDeformationHandler< dim > & get_mesh_deformation_handler () const
 
const LateralAveraging< dim > & get_lateral_averaging () const
 
const AffineConstraints< double > & get_current_constraints () const
 
bool simulator_is_past_initialization () const
 
double get_pressure_scaling () const
 
bool pressure_rhs_needs_compatibility_modification () const
 
bool model_has_prescribed_stokes_solution () const
 
TableHandler & get_statistics_object () const
 
const Postprocess::Manager< dim > & get_postprocess_manager () const
 
const Particle::World< dim > & get_particle_world () const
 
Particle::World< dim > & get_particle_world ()
 
bool is_stokes_matrix_free ()
 
const StokesMatrixFreeHandler< dim > & get_stokes_matrix_free () const
 
RotationProperties< dim > compute_net_angular_momentum (const bool use_constant_density, const LinearAlgebra::BlockVector &solution, const bool limit_to_top_faces=false) const
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 
- Static Public Member Functions inherited from aspect::BoundaryTemperature::Interface< dim >
static void declare_parameters (ParameterHandler &prm)
 
- Static Public Member Functions inherited from aspect::SimulatorAccess< dim >
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 Member Functions

void read_data_OES ()
 
double get_OES (double t) const
 
bool solve_time_step (double &X, double &T, double &R)
 
double get_dT (double r) const
 
double get_Tc (double r) const
 
double get_Ts (double r) const
 
double get_solidus (double X, double p) const
 
double get_initial_Ri (double T)
 
double get_X (double r) const
 
double get_Mass (double r) const
 
double fun_Sn (double B, double R, double n) const
 
double get_Rho (double r) const
 
double get_g (double r) const
 
double get_T (double Tc, double r) const
 
double get_Pressure (double r) const
 
double get_gravity_potential (double r) const
 
void get_specific_heating (double Tc, double &Qs, double &Es)
 
void get_radio_heating (double Tc, double &Qr, double &Er)
 
void get_gravity_heating (double Tc, double r, double X, double &Qg, double &Eg)
 
void get_adiabatic_heating (double Tc, double &Ek, double &Qk)
 
void get_latent_heating (double Tc, double r, double &El, double &Ql)
 
void get_heat_solution (double Tc, double r, double X, double &Eh)
 
double get_radioheating_rate () const
 
void update_core_data ()
 

Private Attributes

internal::CoreData core_data
 
double inner_temperature
 
double outer_temperature
 
double init_dT_dt
 
double init_dR_dt
 
double init_dX_dt
 
bool is_first_call
 
double Rc
 
double CpRho
 
double X_init
 
double Delta
 
double g
 
double P_CMB
 
double P_Core
 
double Tm0
 
double Tm1
 
double Tm2
 
double Theta
 
bool composition_dependency
 
bool use_bw11
 
double K0
 
double Alpha
 
double Rho_0
 
double Rho_cen
 
double Lh
 
double Rh
 
double Beta_c
 
double k_c
 
double Cp
 
unsigned int n_radioheating_elements
 
std::vector< double > heating_rate
 
std::vector< double > half_life
 
std::vector< double > initial_concentration
 
double L
 
double D
 
double Mc
 
int max_steps
 
double dTa
 
std::string name_OES
 
std::vector< struct str_data_OESdata_OES
 

Detailed Description

template<int dim>
class aspect::BoundaryTemperature::DynamicCore< dim >

A class that implements a temperature boundary condition for a spherical shell geometry in which the temperature at the outer surfaces are constant and the core-mantle boundaries (CMB) temperature is calculated by core energy balance. The formulation of the core energy balance is from [1] .

Definition at line 98 of file dynamic_core.h.

Constructor & Destructor Documentation

§ DynamicCore()

Constructor

Member Function Documentation

§ boundary_temperature()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::boundary_temperature ( const types::boundary_id  boundary_indicator,
const Point< dim > &  location 
) const
overridevirtual

Return the temperature that is to hold at a particular location on the boundary of the domain. This function returns the temperatures at the inner and outer boundaries.

Parameters
boundary_indicatorThe boundary indicator of the part of the boundary of the domain on which the point is located at which we are requesting the temperature.
locationThe location of the point at which we ask for the temperature.

Implements aspect::BoundaryTemperature::Interface< dim >.

§ minimal_temperature()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::minimal_temperature ( const std::set< types::boundary_id > &  fixed_boundary_ids) const
overridevirtual

Return the minimal temperature on that part of the boundary on which Dirichlet conditions are posed.

This value is used in computing dimensionless numbers such as the Nusselt number indicating heat flux.

Implements aspect::BoundaryTemperature::Interface< dim >.

§ maximal_temperature()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::maximal_temperature ( const std::set< types::boundary_id > &  fixed_boundary_ids) const
overridevirtual

Return the maximal temperature on that part of the boundary on which Dirichlet conditions are posed.

This value is used in computing dimensionless numbers such as the Nusselt number indicating heat flux.

Implements aspect::BoundaryTemperature::Interface< dim >.

§ declare_parameters()

template<int dim>
static void aspect::BoundaryTemperature::DynamicCore< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare the parameters this class takes through input files. This class declares the inner and outer boundary temperatures.

§ parse_parameters()

template<int dim>
void aspect::BoundaryTemperature::DynamicCore< dim >::parse_parameters ( ParameterHandler &  prm)
overridevirtual

Read the parameters this class declares from the parameter file.

Reimplemented from aspect::BoundaryTemperature::Interface< dim >.

§ update()

template<int dim>
void aspect::BoundaryTemperature::DynamicCore< dim >::update ( )
overridevirtual

This function update the core-mantle boundary (CMB) temperature by the core energy balance solver using the core-mantle boundary heat flux.

Reimplemented from aspect::BoundaryTemperature::Interface< dim >.

§ get_core_data()

template<int dim>
const internal::CoreData& aspect::BoundaryTemperature::DynamicCore< dim >::get_core_data ( ) const

Pass core data to other modules

§ is_OES_used()

template<int dim>
bool aspect::BoundaryTemperature::DynamicCore< dim >::is_OES_used ( ) const

Check if other energy source in the core is in use. The 'other energy source' is used for external core energy source. For example if someone want to test the early lunar core powered by precession (Dwyer, C. A., et al. (2011). "A long-lived lunar dynamo driven by continuous mechanical stirring." Nature 479(7372): 212-214.)

§ read_data_OES()

template<int dim>
void aspect::BoundaryTemperature::DynamicCore< dim >::read_data_OES ( )
private

§ get_OES()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_OES ( double  t) const
private

§ solve_time_step()

template<int dim>
bool aspect::BoundaryTemperature::DynamicCore< dim >::solve_time_step ( double &  X,
double &  T,
double &  R 
)
private

Solve core energy balance for each time step. When solving the change in core-mantle boundary temperature T, inner core radius R, and light component (e.g. S, O, Si) composition X, the following relations has to be respected:

  1. At the inner core boundary the adiabatic temperature should be equal to solidus temperature
  2. The following energy production rate should be balanced in core: Heat flux at core-mantle boundary Q Specific heat Qs*dT/dt Radioactive heating Qr Gravitational contribution Qg*dR/dt Latent heat Ql*dR/dt So that Q+Qs*dT/dt+Qr+Qg*dR/dt*Ql*dR/dt=0
  3. The light component composition X depends on inner core radius (See function get_X() ), and core solidus may dependent on X as well. This becomes a small nonlinear problem. Directly iterate through the above three equations doesn't converge well. Alternatively we solve the inner core radius by bisection method. A single solution between fully liquid and fully solid core is expected. Otherwise this function will throw exception and terminate.

    At Earth core condition, a inner core is forming at the center of the Earth and surrounded by a liquid outer core. However, the core solidus is influenced by light components (e.g. S) and its slope is very closed to core adiabatic. So there is an alternative scenario that the crystallization happens first at the core mantle boundary instead of at the center, which is called a 'snowing core' (Stewart, A. J., et al. (2007). "Mars: a new core-crystallization regime." Science 316(5829): 1323-1325.). This also provides a valid solution for the solver. So the returning bool is set to true for normal core, and false for 'snowing core'. TODO: The current code is only able to treat normal core scenario, treating 'snowing core' scenario may be possible and could be added.

§ get_dT()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_dT ( double  r) const
private

Compute the difference between solidus and adiabatic at inner core boundary for a given inner core radius.

§ get_Tc()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_Tc ( double  r) const
private

Use energy balance to calculate core mantle boundary temperature with a given inner core radius.

§ get_Ts()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_Ts ( double  r) const
private

Get the solidus temperature at inner core boundary with a given inner core radius.

§ get_solidus()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_solidus ( double  X,
double  p 
) const
private

Compute the core solidus at certain pressure

§ get_initial_Ri()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_initial_Ri ( double  T)
private

Get initial inner core radius with given initial core mantle temperature.

§ get_X()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_X ( double  r) const
private

Get the light composition concentration in the outer core from given inner core radius r

§ get_Mass()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_Mass ( double  r) const
private

Compute the mass inside certain radius within the core_data.

§ fun_Sn()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::fun_Sn ( double  B,
double  R,
double  n 
) const
private

Calculate Sn(B,R), referring to [1] .

§ get_Rho()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_Rho ( double  r) const
private

Calculate density at given r

§ get_g()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_g ( double  r) const
private

Calculate gravitational acceleration at given r

§ get_T()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_T ( double  Tc,
double  r 
) const
private

Calculate the core temperature at given r Tc is the temperature at CMB

§ get_Pressure()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_Pressure ( double  r) const
private

Calculate pressure at given r

§ get_gravity_potential()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_gravity_potential ( double  r) const
private

Calculate the gravitational potential at given r

§ get_specific_heating()

template<int dim>
void aspect::BoundaryTemperature::DynamicCore< dim >::get_specific_heating ( double  Tc,
double &  Qs,
double &  Es 
)
private

Calculate energy and entropy change rate factor (regarding the core cooling rated Tc/dt) Qs and Es with given core-mantle boundary (CMB) temperature Tc

§ get_radio_heating()

template<int dim>
void aspect::BoundaryTemperature::DynamicCore< dim >::get_radio_heating ( double  Tc,
double &  Qr,
double &  Er 
)
private

Calculate energy and entropy change rate factor (regarding the radioactive heating rate H) Qr and Er with given CMB temperature Tc

§ get_gravity_heating()

template<int dim>
void aspect::BoundaryTemperature::DynamicCore< dim >::get_gravity_heating ( double  Tc,
double  r,
double  X,
double &  Qg,
double &  Eg 
)
private

Calculate energy and entropy change rate factor (regarding the inner core growth rate dR/dt) Qg and Eg with given Tc(CMB temperature), r(inner core radius), X(light composition concentration)

§ get_adiabatic_heating()

template<int dim>
void aspect::BoundaryTemperature::DynamicCore< dim >::get_adiabatic_heating ( double  Tc,
double &  Ek,
double &  Qk 
)
private

Calculate energy and entropy change rate factor (regarding the core cooling rate Tc/dt) Qk and Ek with given Tc(CMB temperature)

§ get_latent_heating()

template<int dim>
void aspect::BoundaryTemperature::DynamicCore< dim >::get_latent_heating ( double  Tc,
double  r,
double &  El,
double &  Ql 
)
private

Calculate energy and entropy change rate factor (regarding the inner core growth rate dR/dt) Ql and El with given Tc(CMB temperature), r(inner core radius)

§ get_heat_solution()

template<int dim>
void aspect::BoundaryTemperature::DynamicCore< dim >::get_heat_solution ( double  Tc,
double  r,
double  X,
double &  Eh 
)
private

Calculate entropy of heat of solution Eh

§ get_radioheating_rate()

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::get_radioheating_rate ( ) const
private

return radio heating rate at certain time

§ update_core_data()

template<int dim>
void aspect::BoundaryTemperature::DynamicCore< dim >::update_core_data ( )
private

Update the data for core dynamic simulation, the data will be used in the next timestep and for postprocess.

Member Data Documentation

§ core_data

template<int dim>
internal::CoreData aspect::BoundaryTemperature::DynamicCore< dim >::core_data
private

Data for core energy balance it get updated each time step.

Definition at line 179 of file dynamic_core.h.

§ inner_temperature

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::inner_temperature
private

Temperature at the inner boundary.

Definition at line 184 of file dynamic_core.h.

§ outer_temperature

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::outer_temperature
private

Temperatures at the outer boundaries.

Definition at line 189 of file dynamic_core.h.

§ init_dT_dt

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::init_dT_dt
private

Initial CMB temperature changing rate

Definition at line 194 of file dynamic_core.h.

§ init_dR_dt

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::init_dR_dt
private

Initial inner core radius changing rate

Definition at line 199 of file dynamic_core.h.

§ init_dX_dt

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::init_dX_dt
private

Initial light composition changing rate

Definition at line 204 of file dynamic_core.h.

§ is_first_call

template<int dim>
bool aspect::BoundaryTemperature::DynamicCore< dim >::is_first_call
private

Flag for determining the initial call for update().

Definition at line 209 of file dynamic_core.h.

§ Rc

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Rc
private

Core radius

Definition at line 214 of file dynamic_core.h.

§ CpRho

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::CpRho
private

(Heat capacity) * density

Definition at line 219 of file dynamic_core.h.

§ X_init

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::X_init
private

Initial light composition concentration

Definition at line 224 of file dynamic_core.h.

§ Delta

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Delta
private

Partition coefficient of the light element

Definition at line 229 of file dynamic_core.h.

§ g

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::g
private

Gravitational acceleration

Definition at line 234 of file dynamic_core.h.

§ P_CMB

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::P_CMB
private

Pressure at the core mantle boundary

Definition at line 239 of file dynamic_core.h.

§ P_Core

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::P_Core
private

Pressure at the center of the core

Definition at line 244 of file dynamic_core.h.

§ Tm0

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Tm0
private

Parameters for core solidus following: if not dependent on composition Tm(p)= Tm0*(1-Theta)*(1+Tm1*p+Tm2*p^2) if depend on composition X Tm(p)= Tm0*(1-Theta*X)*(1+Tm1*p+Tm2*p^2)

Definition at line 253 of file dynamic_core.h.

§ Tm1

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Tm1
private

Definition at line 254 of file dynamic_core.h.

§ Tm2

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Tm2
private

Definition at line 255 of file dynamic_core.h.

§ Theta

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Theta
private

Definition at line 256 of file dynamic_core.h.

§ composition_dependency

template<int dim>
bool aspect::BoundaryTemperature::DynamicCore< dim >::composition_dependency
private

Definition at line 257 of file dynamic_core.h.

§ use_bw11

template<int dim>
bool aspect::BoundaryTemperature::DynamicCore< dim >::use_bw11
private

If using the Fe-FeS system solidus from Buono & Walker (2011) instead.

Definition at line 262 of file dynamic_core.h.

§ K0

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::K0
private

Compressibility at zero pressure

Definition at line 268 of file dynamic_core.h.

§ Alpha

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Alpha
private

Thermal expansivity

Definition at line 273 of file dynamic_core.h.

§ Rho_0

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Rho_0
private

Density at zero pressure

Definition at line 278 of file dynamic_core.h.

§ Rho_cen

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Rho_cen
private

Density at the center of the planet

Definition at line 283 of file dynamic_core.h.

§ Lh

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Lh
private

Latent heat of fusion

Definition at line 288 of file dynamic_core.h.

§ Rh

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Rh
private

Heat of reaction

Definition at line 293 of file dynamic_core.h.

§ Beta_c

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Beta_c
private

Compositional expansion coefficient

Definition at line 298 of file dynamic_core.h.

§ k_c

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::k_c
private

Heat conductivity of the core

Definition at line 303 of file dynamic_core.h.

§ Cp

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Cp
private

Heat capacity

Definition at line 308 of file dynamic_core.h.

§ n_radioheating_elements

template<int dim>
unsigned int aspect::BoundaryTemperature::DynamicCore< dim >::n_radioheating_elements
private

Number of radioheating element in core

Definition at line 313 of file dynamic_core.h.

§ heating_rate

template<int dim>
std::vector<double> aspect::BoundaryTemperature::DynamicCore< dim >::heating_rate
private

Heating rates of different elements

Definition at line 318 of file dynamic_core.h.

§ half_life

template<int dim>
std::vector<double> aspect::BoundaryTemperature::DynamicCore< dim >::half_life
private

Half life of different elements

Definition at line 323 of file dynamic_core.h.

§ initial_concentration

template<int dim>
std::vector<double> aspect::BoundaryTemperature::DynamicCore< dim >::initial_concentration
private

Initial concentration of different elements

Definition at line 328 of file dynamic_core.h.

§ L

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::L
private

Two length scales in [1] .

Definition at line 333 of file dynamic_core.h.

§ D

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::D
private

Definition at line 334 of file dynamic_core.h.

§ Mc

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::Mc
private

Mass of the core

Definition at line 339 of file dynamic_core.h.

§ max_steps

template<int dim>
int aspect::BoundaryTemperature::DynamicCore< dim >::max_steps
private

Max iterations for the core energy balance solver.

Definition at line 344 of file dynamic_core.h.

§ dTa

template<int dim>
double aspect::BoundaryTemperature::DynamicCore< dim >::dTa
private

Temperature correction value for adiabatic

Definition at line 349 of file dynamic_core.h.

§ name_OES

template<int dim>
std::string aspect::BoundaryTemperature::DynamicCore< dim >::name_OES
private

Other energy source into the core, e.g. the mechanical stirring of the moon.

Definition at line 354 of file dynamic_core.h.

§ data_OES

template<int dim>
std::vector<struct str_data_OES> aspect::BoundaryTemperature::DynamicCore< dim >::data_OES
private

Definition at line 360 of file dynamic_core.h.


The documentation for this class was generated from the following file: