|
void | initialize () override |
|
bool | is_compressible () const override |
|
void | evaluate (const typename Interface< dim >::MaterialModelInputs &in, typename Interface< dim >::MaterialModelOutputs &out) const override |
|
double | enthalpy (const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position) const |
|
std::array< std::pair< double, unsigned int >, 2 > | enthalpy_derivative (const typename Interface< dim >::MaterialModelInputs &in) const |
|
virtual | ~Interface ()=default |
|
virtual void | update () |
|
virtual void | evaluate (const MaterialModel::MaterialModelInputs< dim > &in, MaterialModel::MaterialModelOutputs< dim > &out) const =0 |
|
virtual void | fill_additional_material_model_inputs (MaterialModel::MaterialModelInputs< dim > &input, const LinearAlgebra::BlockVector &solution, const FEValuesBase< dim > &fe_values, const Introspection< dim > &introspection) const |
|
const NonlinearDependence::ModelDependence & | get_model_dependence () const |
|
| 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 |
|
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 |
|
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::BlockVector & | get_current_linearization_point () const |
|
const LinearAlgebra::BlockVector & | get_solution () const |
|
const LinearAlgebra::BlockVector & | get_old_solution () const |
|
const LinearAlgebra::BlockVector & | get_old_old_solution () const |
|
const LinearAlgebra::BlockVector & | get_reaction_vector () const |
|
const LinearAlgebra::BlockVector & | get_mesh_velocity () const |
|
const DoFHandler< dim > & | get_dof_handler () const |
|
const FiniteElement< dim > & | get_fe () const |
|
const LinearAlgebra::BlockSparseMatrix & | get_system_matrix () const |
|
const LinearAlgebra::BlockSparseMatrix & | get_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 std::map< types::boundary_id, std::unique_ptr< BoundaryTraction::Interface< dim > > > & | get_boundary_traction () 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 |
|
|
double | diffusion_viscosity (const double temperature, const double adiabatic_temperature, const double adiabatic_pressure, const double grain_size, const double second_strain_rate_invariant, const Point< dim > &position) const |
|
double | dislocation_viscosity (const double temperature, const double adiabatic_temperature, const double adiabatic_pressure, const SymmetricTensor< 2, dim > &strain_rate, const Point< dim > &position, const double diffusion_viscosity, const double viscosity_guess=0) const |
|
double | density (const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position) const |
|
double | compressibility (const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position) const |
|
double | specific_heat (const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position) const |
|
double | thermal_expansion_coefficient (const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position) const |
|
double | seismic_Vp (const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position) const |
|
double | seismic_Vs (const double temperature, const double pressure, const std::vector< double > &compositional_fields, const Point< dim > &position) const |
|
double | grain_size_change (const double temperature, const double pressure, const std::vector< double > &compositional_fields, const SymmetricTensor< 2, dim > &strain_rate, const Tensor< 1, dim > &velocity, const Point< dim > &position, const unsigned int grain_size_index, const int crossed_transition) const |
|
double | phase_function (const Point< dim > &position, const double temperature, const double pressure, const unsigned int phase) const |
|
unsigned int | get_phase_index (const Point< dim > &position, const double temperature, const double pressure) const |
|
void | convert_log_grain_size (std::vector< double > &compositional_fields) const |
|
template<int dim>
class aspect::MaterialModel::GrainSize< dim >
A material model that relies on compositional fields that stand for average grain sizes of a mineral phase and source terms for them that determine the grain size evolution in dependence of the strain rate, temperature, phase transitions, and the creep regime. This material model only works if a compositional field named 'grain_size' is present. In the diffusion creep regime, the viscosity depends on this grain size. 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 material model.
Definition at line 90 of file grain_size.h.
template<int dim>
double aspect::MaterialModel::GrainSize< dim >::grain_size_change |
( |
const double |
temperature, |
|
|
const double |
pressure, |
|
|
const std::vector< double > & |
compositional_fields, |
|
|
const SymmetricTensor< 2, dim > & |
strain_rate, |
|
|
const Tensor< 1, dim > & |
velocity, |
|
|
const Point< dim > & |
position, |
|
|
const unsigned int |
grain_size_index, |
|
|
const int |
crossed_transition |
|
) |
| const |
|
protected |
Rate of grain size growth (Ostwald ripening) or reduction (due to dynamic recrystallization and phase transformations) in dependence on temperature, pressure, strain rate, mineral phase and creep regime. We use the grain size growth laws as for example 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 the choice between the paleowattmeter (Austins and Evans, 2007) and the paleopiezometer (Hall and Parmentier, 2003) as described in the parameter use_paleowattmeter.