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

Public Member Functions

double melt_fraction (const MaterialModel::MaterialModelInputs< dim > &in, const unsigned int porosity_idx, unsigned int q) const
 
std::vector< double > tian_equilibrium_bound_water_content (const MaterialModel::MaterialModelInputs< dim > &in, unsigned int q) const
 
void parse_parameters (ParameterHandler &prm)
 
- 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 tian_max_peridotite_water
 
double tian_max_gabbro_water
 
double tian_max_MORB_water
 
double tian_max_sediment_water
 
std::vector< double > LR_peridotite_poly_coeffs {-19.0609, 168.983, -630.032, 1281.84, -1543.14, 1111.88, -459.142, 95.4143, 1.97246}
 
std::vector< double > csat_peridotite_poly_coeffs {0.00115628, 2.42179}
 
std::vector< double > Td_peridotite_poly_coeffs {-15.4627, 94.9716, 636.603}
 
std::vector< double > LR_gabbro_poly_coeffs {-1.81745, 7.67198, -10.8507, 5.09329, 8.14519}
 
std::vector< double > csat_gabbro_poly_coeffs {-0.0176673, 0.0893044, 1.52732}
 
std::vector< double > Td_gabbro_poly_coeffs {-1.72277, 20.5898, 637.517}
 
std::vector< double > LR_MORB_poly_coeffs {-1.78177, 7.50871, -10.4840, 5.19725, 7.96365}
 
std::vector< double > csat_MORB_poly_coeffs {0.0102725, -0.115390, 0.324452, 1.41588}
 
std::vector< double > Td_MORB_poly_coeffs {-3.81280, 22.7809, 638.049}
 
std::vector< double > LR_sediment_poly_coeffs {-2.03283, 10.8186, -21.2119, 18.3351, -6.48711, 8.32459}
 
std::vector< double > csat_sediment_poly_coeffs {-0.150662, 0.301807, 1.01867}
 
std::vector< double > Td_sediment_poly_coeffs {2.83277, -24.7593, 85.9090, 524.898}
 
const std::array< double, 4 > pressure_cutoffs {{10, 26, 16, 50}}
 
std::vector< std::vector< double > > devolatilization_enthalpy_changes
 
std::vector< std::vector< double > > water_mass_fractions
 
std::vector< std::vector< double > > devolatilization_onset_temperatures
 

Detailed Description

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

A melt model that calculates the solubility of water according to parameterized phase diagrams for four lithologies: 1) sediment 2) mid-ocean ridge basalt (MORB) 3) gabbro 4) peridotite from Tian, 2019 https://doi.org/10.1029/2019GC008488.

These functions can be used in the calculation of reactive fluid transport of water.

Definition at line 54 of file tian2019_solubility.h.

Member Function Documentation

§ melt_fraction()

template<int dim>
double aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::melt_fraction ( const MaterialModel::MaterialModelInputs< dim > &  in,
const unsigned int  porosity_idx,
unsigned int  q 
) const

Compute the free fluid fraction that is present in the material based on the fluid content of the material and the fluid solubility for the given input conditions. in and melt_fraction need to have the same size.

Parameters
inObject that contains the current conditions.
porosity_idxthe index of the "porosity" composition
qthe quadrature point index

§ tian_equilibrium_bound_water_content()

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::tian_equilibrium_bound_water_content ( const MaterialModel::MaterialModelInputs< dim > &  in,
unsigned int  q 
) const

Compute the maximum allowed bound water content at the input pressure and temperature conditions. This is used to determine how free water interacts with the solid phase.

Parameters
inObject that contains the current conditions.
qthe quadrature point index

§ declare_parameters()

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

Declare the parameters this function takes through input files.

§ parse_parameters()

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

Read the parameters from the parameter file.

Member Data Documentation

§ tian_max_peridotite_water

template<int dim>
double aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::tian_max_peridotite_water
private

The maximum water content for each of the 4 rock types in the tian approximation method. These are important for keeping the polynomial bounded within reasonable values.

Definition at line 102 of file tian2019_solubility.h.

§ tian_max_gabbro_water

template<int dim>
double aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::tian_max_gabbro_water
private

Definition at line 103 of file tian2019_solubility.h.

§ tian_max_MORB_water

template<int dim>
double aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::tian_max_MORB_water
private

Definition at line 104 of file tian2019_solubility.h.

§ tian_max_sediment_water

template<int dim>
double aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::tian_max_sediment_water
private

Definition at line 105 of file tian2019_solubility.h.

§ LR_peridotite_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::LR_peridotite_poly_coeffs {-19.0609, 168.983, -630.032, 1281.84, -1543.14, 1111.88, -459.142, 95.4143, 1.97246}
private

The following coefficients are taken from a publication from Tian et al., 2019, and can be found in Table 3 (Gabbro), Table B1 (MORB), Table B2 (Sediments) and Table B3 (peridotite). LR refers to the effective enthalpy change for devolatilization reactions, csat is the saturated mass fraction of water in the solid, and Td is the onset temperature of devolatilization for water.

Definition at line 115 of file tian2019_solubility.h.

§ csat_peridotite_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::csat_peridotite_poly_coeffs {0.00115628, 2.42179}
private

Definition at line 116 of file tian2019_solubility.h.

§ Td_peridotite_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::Td_peridotite_poly_coeffs {-15.4627, 94.9716, 636.603}
private

Definition at line 117 of file tian2019_solubility.h.

§ LR_gabbro_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::LR_gabbro_poly_coeffs {-1.81745, 7.67198, -10.8507, 5.09329, 8.14519}
private

Definition at line 119 of file tian2019_solubility.h.

§ csat_gabbro_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::csat_gabbro_poly_coeffs {-0.0176673, 0.0893044, 1.52732}
private

Definition at line 120 of file tian2019_solubility.h.

§ Td_gabbro_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::Td_gabbro_poly_coeffs {-1.72277, 20.5898, 637.517}
private

Definition at line 121 of file tian2019_solubility.h.

§ LR_MORB_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::LR_MORB_poly_coeffs {-1.78177, 7.50871, -10.4840, 5.19725, 7.96365}
private

Definition at line 123 of file tian2019_solubility.h.

§ csat_MORB_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::csat_MORB_poly_coeffs {0.0102725, -0.115390, 0.324452, 1.41588}
private

Definition at line 124 of file tian2019_solubility.h.

§ Td_MORB_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::Td_MORB_poly_coeffs {-3.81280, 22.7809, 638.049}
private

Definition at line 125 of file tian2019_solubility.h.

§ LR_sediment_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::LR_sediment_poly_coeffs {-2.03283, 10.8186, -21.2119, 18.3351, -6.48711, 8.32459}
private

Definition at line 127 of file tian2019_solubility.h.

§ csat_sediment_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::csat_sediment_poly_coeffs {-0.150662, 0.301807, 1.01867}
private

Definition at line 128 of file tian2019_solubility.h.

§ Td_sediment_poly_coeffs

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::Td_sediment_poly_coeffs {2.83277, -24.7593, 85.9090, 524.898}
private

Definition at line 129 of file tian2019_solubility.h.

§ pressure_cutoffs

template<int dim>
const std::array<double,4 > aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::pressure_cutoffs {{10, 26, 16, 50}}
private

The polynomials breakdown above certain pressures, 10 GPa for peridotite, 26 GPa for gabbro, 16 GPa for MORB, and 50 GPa for sediment. These cutoff pressures were determined by extending the pressure range in Tian et al. (2019) and observing where the maximum allowed water contents jump towards infinite values.

Definition at line 136 of file tian2019_solubility.h.

§ devolatilization_enthalpy_changes

template<int dim>
std::vector<std::vector<double> > aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::devolatilization_enthalpy_changes
private

§ water_mass_fractions

template<int dim>
std::vector<std::vector<double> > aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::water_mass_fractions
private

§ devolatilization_onset_temperatures

template<int dim>
std::vector<std::vector<double> > aspect::MaterialModel::ReactionModel::Tian2019Solubility< dim >::devolatilization_onset_temperatures
private

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