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

Public Member Functions

 Katz2003MantleMelting ()
 
void parse_parameters (ParameterHandler &prm)
 
double melt_fraction (const double temperature, const double pressure) const
 
double entropy_change (const double temperature, const double pressure, const double maximum_melt_fraction, const NonlinearDependence::Dependence dependence) const
 
void calculate_reaction_rate_outputs (const typename Interface< dim >::MaterialModelInputs &in, typename Interface< dim >::MaterialModelOutputs &out) const
 
void calculate_fluid_outputs (const typename Interface< dim >::MaterialModelInputs &in, typename Interface< dim >::MaterialModelOutputs &out, const double reference_T) const
 
double reference_darcy_coefficient () const
 
- 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
 
unsigned int n_particle_managers () const
 
const Particle::Manager< dim > & get_particle_manager (const unsigned int particle_manager_index) const
 
Particle::Manager< dim > & get_particle_manager (const unsigned int particle_manager_index)
 
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::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 Attributes

double reference_rho_fluid
 
double xi_0
 
double viscosity_fluid
 
double thermal_bulk_viscosity_exponent
 
double alpha_phi
 
double extraction_depth
 
double melt_compressibility
 
double melt_bulk_modulus_derivative
 
double depletion_solidus_change
 
bool fractional_melting
 
double freezing_rate
 
double melting_time_scale
 
double reference_permeability
 
double A1
 
double A2
 
double A3
 
double B1
 
double B2
 
double B3
 
double C1
 
double C2
 
double C3
 
double r1
 
double r2
 
double M_cpx
 
double beta
 
double peridotite_melting_entropy_change
 

Detailed Description

template<int dim>
class aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >

A melt model that calculates melt fraction and entropy change according to the melting model for dry peridotite of Katz, 2003. This also includes a computation of the latent heat of melting (if the latent heat heating model is active).

These functions can be used in the calculation of melting and melt transport in the melt_simple material model and can be extended to other material models

Definition at line 50 of file katz2003_mantle_melting.h.

Constructor & Destructor Documentation

§ Katz2003MantleMelting()

Member Function Documentation

§ declare_parameters()

template<int dim>
static void aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare the parameters this function takes through input files.

§ parse_parameters()

template<int dim>
void aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::parse_parameters ( ParameterHandler &  prm)

Read the parameters from the parameter file.

§ melt_fraction()

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::melt_fraction ( const double  temperature,
const double  pressure 
) const

Percentage of material that is molten for a given temperature and pressure (assuming equilibrium conditions). Melting model after Katz, 2003, for dry peridotite.

§ entropy_change()

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::entropy_change ( const double  temperature,
const double  pressure,
const double  maximum_melt_fraction,
const NonlinearDependence::Dependence  dependence 
) const

Compute the change in entropy due to melting for a given temperature and pressure, and under the assumption that a fraction maximum_melt_fraction of the material has already been molten previously. The entropy change is computed with respect to temperature or pressure, depending on dependence. This is needed to calculate the latent heat of melt.

§ calculate_reaction_rate_outputs()

template<int dim>
void aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::calculate_reaction_rate_outputs ( const typename Interface< dim >::MaterialModelInputs in,
typename Interface< dim >::MaterialModelOutputs out 
) const

Compute all the reaction rate variables needed for a reactive transport model based on the Katz 2003 formulation. Takes the material model inputs in to compute the material model outputs out. This function mainly fills the reaction_rate_out object but populates out.reaction_terms, out.entropy_derivative_pressure and entropy_derivative_temperature

§ calculate_fluid_outputs()

template<int dim>
void aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::calculate_fluid_outputs ( const typename Interface< dim >::MaterialModelInputs in,
typename Interface< dim >::MaterialModelOutputs out,
const double  reference_T 
) const

Compute all the fluid variables needed for a reactive transport model based on the Katz 2003 formulation. This function fills melt outputs, the out object should already contain outputs for the solid and this function uses the inputs in and the solid outputs out to fill MeltOutputs. Solid outputs such as out.Thermal_expansion_coefficients are expected to have already been computed when this function is called. Solid viscosities are also modified in the out object here because the presence of melt weakens the material.

§ reference_darcy_coefficient()

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::reference_darcy_coefficient ( ) const

Member Data Documentation

§ reference_rho_fluid

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::reference_rho_fluid
private

Parameters for anhydrous melting of peridotite after Katz, 2003

Definition at line 123 of file katz2003_mantle_melting.h.

§ xi_0

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::xi_0
private

Definition at line 124 of file katz2003_mantle_melting.h.

§ viscosity_fluid

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::viscosity_fluid
private

Definition at line 125 of file katz2003_mantle_melting.h.

§ thermal_bulk_viscosity_exponent

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::thermal_bulk_viscosity_exponent
private

Definition at line 126 of file katz2003_mantle_melting.h.

§ alpha_phi

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::alpha_phi
private

Definition at line 127 of file katz2003_mantle_melting.h.

§ extraction_depth

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::extraction_depth
private

Definition at line 128 of file katz2003_mantle_melting.h.

§ melt_compressibility

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::melt_compressibility
private

Definition at line 129 of file katz2003_mantle_melting.h.

§ melt_bulk_modulus_derivative

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::melt_bulk_modulus_derivative
private

Definition at line 130 of file katz2003_mantle_melting.h.

§ depletion_solidus_change

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::depletion_solidus_change
private

Definition at line 131 of file katz2003_mantle_melting.h.

§ fractional_melting

template<int dim>
bool aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::fractional_melting
private

Definition at line 132 of file katz2003_mantle_melting.h.

§ freezing_rate

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::freezing_rate
private

Definition at line 133 of file katz2003_mantle_melting.h.

§ melting_time_scale

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::melting_time_scale
private

Definition at line 134 of file katz2003_mantle_melting.h.

§ reference_permeability

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::reference_permeability
private

Definition at line 135 of file katz2003_mantle_melting.h.

§ A1

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::A1
private

Definition at line 138 of file katz2003_mantle_melting.h.

§ A2

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::A2
private

Definition at line 139 of file katz2003_mantle_melting.h.

§ A3

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::A3
private

Definition at line 140 of file katz2003_mantle_melting.h.

§ B1

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::B1
private

Definition at line 143 of file katz2003_mantle_melting.h.

§ B2

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::B2
private

Definition at line 144 of file katz2003_mantle_melting.h.

§ B3

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::B3
private

Definition at line 145 of file katz2003_mantle_melting.h.

§ C1

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::C1
private

Definition at line 148 of file katz2003_mantle_melting.h.

§ C2

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::C2
private

Definition at line 149 of file katz2003_mantle_melting.h.

§ C3

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::C3
private

Definition at line 150 of file katz2003_mantle_melting.h.

§ r1

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::r1
private

Definition at line 153 of file katz2003_mantle_melting.h.

§ r2

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::r2
private

Definition at line 154 of file katz2003_mantle_melting.h.

§ M_cpx

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::M_cpx
private

Definition at line 155 of file katz2003_mantle_melting.h.

§ beta

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::beta
private

Definition at line 158 of file katz2003_mantle_melting.h.

§ peridotite_melting_entropy_change

template<int dim>
double aspect::MaterialModel::ReactionModel::Katz2003MantleMelting< dim >::peridotite_melting_entropy_change
private

Definition at line 161 of file katz2003_mantle_melting.h.


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