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

Classes

struct  Formulation
 

Public Member Functions

void initialize_phase_function (std::shared_ptr< MaterialUtilities::PhaseFunction< dim >> &phase_function)
 
void calculate_reaction_terms (const typename Interface< dim >::MaterialModelInputs &in, const std::vector< double > &adiabatic_pressure, const std::vector< unsigned int > &phase_indices, const std::function< double(const double, const double, const double, const SymmetricTensor< 2, dim > &, const unsigned int, const double, const double)> &dislocation_viscosity, const std::function< double(const double, const double, const double, const double, const double, const unsigned int)> &diffusion_viscosity, const double min_eta, const double max_eta, typename Interface< dim >::MaterialModelOutputs &out) const
 
void create_additional_named_outputs (MaterialModel::MaterialModelOutputs< dim > &out) const
 
void fill_additional_outputs (const typename MaterialModel::MaterialModelInputs< dim > &in, const typename MaterialModel::MaterialModelOutputs< dim > &out, const std::vector< unsigned int > &phase_indices, const std::vector< double > &dislocation_viscosities, std::vector< std::unique_ptr< MaterialModel::AdditionalMaterialOutputs< dim >>> &additional_outputs) 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 Member Functions

double compute_partitioning_fraction (const double temperature) const
 

Private Attributes

std::vector< double > grain_growth_activation_energy
 
std::vector< double > grain_growth_activation_volume
 
std::vector< double > grain_growth_rate_constant
 
std::vector< double > grain_growth_exponent
 
double minimum_grain_size
 
std::vector< double > reciprocal_required_strain
 
std::vector< double > recrystallized_grain_size
 
std::vector< double > grain_boundary_energy
 
std::vector< double > boundary_area_change_work_fraction
 
std::vector< double > geometric_constant
 
Formulation::Kind grain_size_evolution_formulation
 
double grain_size_reduction_work_fraction_exponent
 
double minimum_grain_size_reduction_work_fraction
 
double maximum_grain_size_reduction_work_fraction
 
double temperature_minimum_partitioning_power
 
double temperature_maximum_partitioning_power
 
double phase_distribution
 
double roughness_to_grain_size
 
std::shared_ptr< MaterialUtilities::PhaseFunction< dim > > phase_function
 
unsigned int n_phase_transitions
 

Detailed Description

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

A model to calculate the change of the average grain size according to different published grain size evolution models. We use the grain size evolution laws described in Behn et al., 2009. Implications of grain size evolution on the seismic structure of the oceanic upper mantle, Earth Planet. Sci. Letters, 282, 178–189. Other material parameters are either prescribed similar to the 'simple' material model, or read from data files that were generated by the Perplex or Hefesto software. The material model is described in more detail in Dannberg, J., Z. Eilon, U. Faul, R. Gassmöller, P. Moulik, and R. Myhill (2017), The importance of grain size to mantle dynamics and seismological observations, Geochem. Geophys. Geosyst., 18, 3034–3061, doi:10.1002/2017GC006944., which is the canonical reference for this reaction model.

Definition at line 52 of file grain_size_evolution.h.

Member Function Documentation

§ initialize_phase_function()

template<int dim>
void aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::initialize_phase_function ( std::shared_ptr< MaterialUtilities::PhaseFunction< dim >> &  phase_function)

Receive the phase function object that is used to determine if material crossed a phase transition and grain size needs to be reset.

§ calculate_reaction_terms()

template<int dim>
void aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::calculate_reaction_terms ( const typename Interface< dim >::MaterialModelInputs in,
const std::vector< double > &  adiabatic_pressure,
const std::vector< unsigned int > &  phase_indices,
const std::function< double(const double, const double, const double, const SymmetricTensor< 2, dim > &, const unsigned int, const double, const double)> &  dislocation_viscosity,
const std::function< double(const double, const double, const double, const double, const double, const unsigned int)> &  diffusion_viscosity,
const double  min_eta,
const double  max_eta,
typename Interface< dim >::MaterialModelOutputs out 
) const

Rate of grain size growth (Ostwald ripening) or reduction (due to dynamic recrystallization and phase transformations) depends on temperature, pressure, strain rate, mineral phase and creep regime. We use the grain size growth laws as described in Behn, M. D., Hirth, G., & Elsenbeck, J. R. (2009). Implications of grain size evolution on the seismic structure of the oceanic upper mantle. Earth and Planetary Science Letters, 282(1), 178-189.

For the rate of grain size reduction due to dynamic crystallization there is a choice between different models as described in the struct Formulation.

Note
The functions dislocation_creep and diffusion_creep currently follow the interface of the respective functions in the grain size material model. This may be changed in the future to match the interface of the DiffusionCreep and DislocationCreep rheology modules.

§ create_additional_named_outputs()

template<int dim>
void aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::create_additional_named_outputs ( MaterialModel::MaterialModelOutputs< dim > &  out) const

Create the additional material model outputs produced by this reaction model. In this case the function creates ShearHeatingOutputs, which reduce the amount of work done to produce shear heat to account for the surface energy used to reduce the grain size.

§ fill_additional_outputs()

template<int dim>
void aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::fill_additional_outputs ( const typename MaterialModel::MaterialModelInputs< dim > &  in,
const typename MaterialModel::MaterialModelOutputs< dim > &  out,
const std::vector< unsigned int > &  phase_indices,
const std::vector< double > &  dislocation_viscosities,
std::vector< std::unique_ptr< MaterialModel::AdditionalMaterialOutputs< dim >>> &  additional_outputs 
) const

Fill the additional output objects created by the function create_additional_named_outputs. Note that although the plugin gets write access to all AdditionalMaterialOutputs it should only modify the ones it created.

§ declare_parameters()

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

Declare the parameters this function takes through input files.

§ parse_parameters()

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

Read the parameters from the parameter file.

§ compute_partitioning_fraction()

template<int dim>
double aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::compute_partitioning_fraction ( const double  temperature) const
private

This function returns the fraction of shear heating energy partitioned into grain damage using the implementation by Mulyukova and Bercovici (2018) Collapse of passive margins by lithospheric damage and plunging grain size. Earth and Planetary Science Letters, 484, 341-352.

Member Data Documentation

§ grain_growth_activation_energy

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::grain_growth_activation_energy
private

Parameters controlling the grain size evolution.

Definition at line 132 of file grain_size_evolution.h.

§ grain_growth_activation_volume

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::grain_growth_activation_volume
private

Definition at line 133 of file grain_size_evolution.h.

§ grain_growth_rate_constant

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::grain_growth_rate_constant
private

Definition at line 134 of file grain_size_evolution.h.

§ grain_growth_exponent

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::grain_growth_exponent
private

Definition at line 135 of file grain_size_evolution.h.

§ minimum_grain_size

template<int dim>
double aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::minimum_grain_size
private

Definition at line 136 of file grain_size_evolution.h.

§ reciprocal_required_strain

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::reciprocal_required_strain
private

Definition at line 137 of file grain_size_evolution.h.

§ recrystallized_grain_size

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::recrystallized_grain_size
private

Definition at line 138 of file grain_size_evolution.h.

§ grain_boundary_energy

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::grain_boundary_energy
private

Parameters controlling the dynamic grain recrystallization.

Definition at line 143 of file grain_size_evolution.h.

§ boundary_area_change_work_fraction

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::boundary_area_change_work_fraction
private

Definition at line 144 of file grain_size_evolution.h.

§ geometric_constant

template<int dim>
std::vector<double> aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::geometric_constant
private

Definition at line 145 of file grain_size_evolution.h.

§ grain_size_evolution_formulation

template<int dim>
Formulation::Kind aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::grain_size_evolution_formulation
private

A variable that records the formulation of how to evolve grain size. See the type of this variable for a description of available options.

Definition at line 205 of file grain_size_evolution.h.

§ grain_size_reduction_work_fraction_exponent

template<int dim>
double aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::grain_size_reduction_work_fraction_exponent
private

Parameters controlling the partitioning of energy into grain damage in the pinned state.

Definition at line 219 of file grain_size_evolution.h.

§ minimum_grain_size_reduction_work_fraction

template<int dim>
double aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::minimum_grain_size_reduction_work_fraction
private

Definition at line 220 of file grain_size_evolution.h.

§ maximum_grain_size_reduction_work_fraction

template<int dim>
double aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::maximum_grain_size_reduction_work_fraction
private

Definition at line 221 of file grain_size_evolution.h.

§ temperature_minimum_partitioning_power

template<int dim>
double aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::temperature_minimum_partitioning_power
private

Definition at line 222 of file grain_size_evolution.h.

§ temperature_maximum_partitioning_power

template<int dim>
double aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::temperature_maximum_partitioning_power
private

Definition at line 223 of file grain_size_evolution.h.

§ phase_distribution

template<int dim>
double aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::phase_distribution
private

Functions and parameters controlling conversion from interface roughness to grain size, used in pinned state formulation of grain damage. This conversion depends on the proportion of the two mineral phases.

A detailed description of this approach can be found in Appendix H.1, in Equation (8) in the main manuscript, and in equation (F.28) of Bercovici and Ricard (2012). Mechanisms for the generation of plate tectonics by two-phase grain-damage and pinning. Physics of the Earth and Planetary Interiors 202 (2012): 27-55.

Definition at line 235 of file grain_size_evolution.h.

§ roughness_to_grain_size

template<int dim>
double aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::roughness_to_grain_size
private

The factor used to convert interface roughness into the equivalent mean grain size for a given volume fraction of a mineral in the two-phase damage model. "Interface roughness" (or radius of curvature) is a measure for how much interface area is present in a collection of mineral grains and therefore dependent on the average grain size. For a discussion see Bercovici and Ricard (2012), Appendix H.1.

Definition at line 245 of file grain_size_evolution.h.

§ phase_function

template<int dim>
std::shared_ptr<MaterialUtilities::PhaseFunction<dim> > aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::phase_function
private

Object that handles phase transitions. Allows it to compute the phase function for each individual phase transition in the model, given the temperature, pressure, depth, and density gradient. This variable is a shared pointer so that material model and reaction model can share the same phase function. You can initialize the phase function in the material model and pass it to the reaction model using the function initialize_phase_function.

Definition at line 256 of file grain_size_evolution.h.

§ n_phase_transitions

template<int dim>
unsigned int aspect::MaterialModel::ReactionModel::GrainSizeEvolution< dim >::n_phase_transitions
private

Cache the total number of phase transitions defined in phase_function.

Definition at line 261 of file grain_size_evolution.h.


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