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

Public Member Functions

void evaluate (const MaterialModel::MaterialModelInputs< dim > &in, MaterialModel::MaterialModelOutputs< dim > &out) const override
 
bool is_compressible () const override
 
void parse_parameters (ParameterHandler &prm) override
 
void create_additional_named_outputs (MaterialModel::MaterialModelOutputs< dim > &out) const override
 
double get_min_strain_rate () const
 
DEAL_II_DEPRECATED bool is_yielding (const double pressure, const double temperature, const std::vector< double > &composition, const SymmetricTensor< 2, dim > &strain_rate) const
 
bool is_yielding (const MaterialModelInputs< dim > &in) const
 
- Public Member Functions inherited from aspect::MaterialModel::Interface< dim >
virtual ~Interface ()=default
 
virtual void initialize ()
 
virtual void update ()
 
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::ModelDependenceget_model_dependence () 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
 
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::MaterialModel::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 Attributes

std::unique_ptr< Rheology::ViscoPlastic< dim > > rheology
 
std::vector< double > thermal_diffusivities
 
bool define_conductivities
 
std::vector< double > thermal_conductivities
 
std::vector< unsigned int > n_phase_transitions_for_each_chemical_composition
 
unsigned int n_phases
 
EquationOfState::MulticomponentIncompressible< dim > equation_of_state
 
MaterialUtilities::PhaseFunction< dim > phase_function
 

Additional Inherited Members

- Public Types inherited from aspect::MaterialModel::Interface< dim >
using MaterialModelInputs = MaterialModel::MaterialModelInputs< dim >
 
using MaterialModelOutputs = MaterialModel::MaterialModelOutputs< dim >
 
- Protected Attributes inherited from aspect::MaterialModel::Interface< dim >
NonlinearDependence::ModelDependence model_dependence
 

Detailed Description

template<int dim>
class aspect::MaterialModel::ViscoPlastic< dim >

A material model combining viscous and plastic deformation, with the option to also include viscoelastic deformation.

Viscous deformation is defined by a viscous flow law describing dislocation and diffusion creep: \( v = \frac{1}{2} A^{-\frac{1}{n}} d^{\frac{m}{n}} \dot{\varepsilon}_{ii}^{\frac{1-n}{n}} \exp\left(\frac{E + PV}{nRT}\right) \) where where \(A\) is the prefactor, \(n\) is the stress exponent, \(\dot{\varepsilon}_{ii}\) is the square root of the deviatoric strain rate tensor second invariant, \(d\) is grain size, \(m\) is the grain size exponent, \(E\) is activation energy, \(V\) is activation volume, \(P\) is pressure, \(R\) is the gas exponent and \(T\) is temperature.

One may select to use the diffusion ( \(v_{diff}\); \(n=1\), \(m!=0\)), dislocation ( \(v_{disl}\), \(n>1\), \(m=0\)) or composite \(\frac{v_{diff}v_{disl}}{v_{diff}+v_{disl}}\) equation form.

Viscous stress is limited by plastic deformation, which follows a Drucker Prager yield criterion: \(\sigma_y = C\cos(\phi) + P\sin(\phi)\) (2D) or in 3D \(\sigma_y = \frac{6C\cos(\phi) + 2P\sin(\phi)}{\sqrt{3}(3+\sin(\phi))}\) where \(\sigma_y\) is the yield stress, \(C\) is cohesion, \(phi\) is the angle of internal friction and \(P\) is pressure. If the viscous stress ( \(2v{\varepsilon}_{ii})\)) exceeds the yield stress ( \(\sigma_{y}\)), the viscosity is rescaled back to the yield surface: \(v_{y}=\sigma_{y}/(2{\varepsilon}_{ii})\)

When included, the viscoelastic rheology takes into account the elastic shear strength (e.g., shear modulus), while the tensile and volumetric strength (e.g., Young's and bulk modulus) are not considered. The model is incompressible and allows specifying an arbitrary number of compositional fields, where each field represents a different rock type or component of the viscoelastic stress tensor. The symmetric stress tensor in 2D and 3D, respectively, contains 3 or 6 components. The compositional fields representing these components must be named and listed in a very specific format, which is designed to minimize mislabeling stress tensor components as distinct 'compositional rock types' (or vice versa). For 2D models, the first three compositional fields must be labeled ve_stress_xx, ve_stress_yy and ve_stress_xy. In 3D, the first six compositional fields must be labeled ve_stress_xx, ve_stress_yy, ve_stress_zz, ve_stress_xy, ve_stress_xz, ve_stress_yz.

Combining this viscoelasticity implementation with non-linear viscous flow and plasticity produces a constitutive relationship commonly referred to as partial elastoviscoplastic (e.g., pEVP) in the geodynamics community. While extensively discussed and applied within the geodynamics literature, notable references include: Moresi et al. (2003), J. Comp. Phys., v. 184, p. 476-497. Gerya and Yuen (2007), Phys. Earth. Planet. Inter., v. 163, p. 83-105. Gerya (2010), Introduction to Numerical Geodynamic Modeling. Kaus (2010), Tectonophysics, v. 484, p. 36-47. Choi et al. (2013), J. Geophys. Res., v. 118, p. 2429-2444. Keller et al. (2013), Geophys. J. Int., v. 195, p. 1406-1442.

The overview below directly follows Moresi et al. (2003) eqns. 23-32. However, an important distinction between this material model and the studies above is the use of compositional fields, rather than particles, to track individual components of the viscoelastic stress tensor. The material model will be updated when an option to track and calculate viscoelastic stresses with particles is implemented. Moresi et al. (2003) begins (eqn. 23) by writing the deviatoric rate of deformation ( \(\hat{D}\)) as the sum of elastic ( \(\hat{D_{e}}\)) and viscous ( \(\hat{D_{v}}\)) components: \(\hat{D} = \hat{D_{e}} + \hat{D_{v}}\). These terms further decompose into \(\hat{D_{v}} = \frac{\tau}{2\eta}\) and \(\hat{D_{e}} = \frac{\overset{\triangledown}{\tau}}{2\mu}\), where \(\tau\) is the viscous deviatoric stress, \(\eta\) is the shear viscosity, \(\mu\) is the shear modulus and \(\overset{\triangledown}{\tau}\) is the Jaumann corotational stress rate. If plasticity is included the deviatoric rate of deformation may be written as: \(\hat{D} = \hat{D_{e}} + \hat{D_{v}} + \hat{D_{p}}\), where \(\hat{D_{p}}\) is the plastic component. As defined in the second paragraph, \(\hat{D_{p}}\) decomposes to \(\frac{\tau_{y}}{2\eta_{y}}\), where \(\tau_{y}\) is the yield stress and \(\eta_{y}\) is the viscosity rescaled to the yield surface.

Above, the Jaumann corotational stress rate (eqn. 24) from the elastic component contains the time derivative of the deviatoric stress ( \(\dot{\tau}\)) and terms that account for material spin (e.g., rotation) due to advection: \(\overset{\triangledown}{\tau} = \dot{\tau} + {\tau}W -W\tau\). Above, \(W\) is the material spin tensor (eqn. 25): \(W_{ij} = \frac{1}{2} \left (\frac{\partial V_{i}}{\partial x_{j}} - \frac{\partial V_{j}}{\partial x_{i}} \right )\).

The Jaumann stress-rate can also be approximated using terms from the time at the previous time step ( \(t\)) and current time step ( \(t + \Delta t^{e}\)): \(\smash[t]{\overset{\triangledown}{\tau}}^{t + \Delta t^{e}} \approx \frac{\tau^{t + \Delta t^{e} - \tau^{t}}}{\Delta t^{e}} - W^{t}\tau^{t} + \tau^{t}W^{t}\). In this material model, the size of the time step above ( \(\Delta t^{e}\)) can be specified as the numerical time step size or an independent fixed time step. If the latter case is selected, the user has an option to apply a stress averaging scheme to account for the differences between the numerical and fixed elastic time step (eqn. 32). If one selects to use a fixed elastic time step throughout the model run, an equal numerical and elastic time step can be achieved by using CFL and maximum time step values that restrict the numerical time step to the fixed elastic time step.

The formulation above allows rewriting the total rate of deformation (eqn. 29) as \(\tau^{t + \Delta t^{e}} = \eta_{eff} \left ( 2\hat{D}^{t + \triangle t^{e}} + \frac{\tau^{t}}{\mu \Delta t^{e}} + \frac{W^{t}\tau^{t} - \tau^{t}W^{t}}{\mu} \right )\).

The effective viscosity (eqn. 28) is a function of the viscosity ( \(\eta\)), elastic time step size ( \(\Delta t^{e}\)) and shear relaxation time ( \( \alpha = \frac{\eta}{\mu} \)): \(\eta_{eff} = \eta \frac{\Delta t^{e}}{\Delta t^{e} + \alpha}\) The magnitude of the shear modulus thus controls how much the effective viscosity is reduced relative to the initial viscosity.

Elastic effects are introduced into the governing Stokes equations through an elastic force term (eqn. 30) using stresses from the previous time step: \(F^{e,t} = -\frac{\eta_{eff}}{\mu \Delta t^{e}} \tau^{t}\). This force term is added onto the right-hand side force vector in the system of equations.

Several model parameters (reference densities, thermal expansivities thermal diffusivities, heat capacities and rheology parameters) can be defined per-compositional field. For each material parameter the user supplies a comma delimited list of length N+1, where N is the number of compositional fields. The additional field corresponds to the value for background material. They should be ordered ``background, composition1, composition2...''

If a list of values is given for the density, thermal expansivity, thermal diffusivity and heat capacity, the volume weighted sum of the values of each of the compositional fields is used in their place, for example \(\rho = \sum \left( \rho_i V_i \right)\)

The individual output viscosities for each compositional field are also averaged. The user can choose from a range of options for this viscosity averaging. If only one value is given for any of these parameters, all compositions are assigned the same value. The first value in the list is the value assigned to "background material" (regions where the sum of the compositional fields is < 1.0).

Definition at line 181 of file visco_plastic.h.

Member Function Documentation

§ evaluate()

template<int dim>
void aspect::MaterialModel::ViscoPlastic< dim >::evaluate ( const MaterialModel::MaterialModelInputs< dim > &  in,
MaterialModel::MaterialModelOutputs< dim > &  out 
) const
overridevirtual

Function to compute the material properties in out given the inputs in in.

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

§ is_compressible()

template<int dim>
bool aspect::MaterialModel::ViscoPlastic< dim >::is_compressible ( ) const
overridevirtual

Return whether the model is compressible or not. Incompressibility does not necessarily imply that the density is constant; rather, it may still depend on temperature or pressure. In the current context, compressibility means whether we should solve the continuity equation as \(\nabla \cdot (\rho \mathbf u)=0\) (compressible Stokes) or as \(\nabla \cdot \mathbf{u}=0\) (incompressible Stokes).

This material model is incompressible.

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

§ declare_parameters()

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

§ parse_parameters()

template<int dim>
void aspect::MaterialModel::ViscoPlastic< dim >::parse_parameters ( ParameterHandler &  prm)
overridevirtual

Read the parameters this class declares from the parameter file. The default implementation of this function does not read any parameters. Consequently, derived classes do not have to overload this function if they do not take any runtime parameters.

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

§ create_additional_named_outputs()

template<int dim>
void aspect::MaterialModel::ViscoPlastic< dim >::create_additional_named_outputs ( MaterialModel::MaterialModelOutputs< dim > &  outputs) const
overridevirtual

If this material model can produce additional named outputs that are derived from NamedAdditionalOutputs, create them in here. By default, this does nothing.

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

§ get_min_strain_rate()

template<int dim>
double aspect::MaterialModel::ViscoPlastic< dim >::get_min_strain_rate ( ) const

§ is_yielding() [1/2]

template<int dim>
DEAL_II_DEPRECATED bool aspect::MaterialModel::ViscoPlastic< dim >::is_yielding ( const double  pressure,
const double  temperature,
const std::vector< double > &  composition,
const SymmetricTensor< 2, dim > &  strain_rate 
) const

A function that returns whether the material is plastically yielding at the given pressure, temperature, composition, and strain rate.

Deprecated:
: Use the other function with this name instead, which allows to pass in more general input variables.

§ is_yielding() [2/2]

template<int dim>
bool aspect::MaterialModel::ViscoPlastic< dim >::is_yielding ( const MaterialModelInputs< dim > &  in) const

A function that returns whether the material is plastically yielding at the given input variables (pressure, temperature, composition, strain rate, and so on).

Member Data Documentation

§ rheology

template<int dim>
std::unique_ptr<Rheology::ViscoPlastic<dim> > aspect::MaterialModel::ViscoPlastic< dim >::rheology
private

Pointer to the object used to compute the rheological properties. In this case, the rheology in question is visco(elasto)plastic. The object contains functions for parameter declaration and parsing, and further functions that calculate viscosity and viscosity derivatives. It also contains functions that create and fill additional material model outputs, specifically plastic outputs. The rheology itself is a composite rheology, and so the object contains further objects and/or pointers to objects that provide functions and parameters for all subordinate rheologies.

Definition at line 247 of file visco_plastic.h.

§ thermal_diffusivities

template<int dim>
std::vector<double> aspect::MaterialModel::ViscoPlastic< dim >::thermal_diffusivities
private

Definition at line 249 of file visco_plastic.h.

§ define_conductivities

template<int dim>
bool aspect::MaterialModel::ViscoPlastic< dim >::define_conductivities
private

Whether to use user-defined thermal conductivities instead of thermal diffusivities.

Definition at line 254 of file visco_plastic.h.

§ thermal_conductivities

template<int dim>
std::vector<double> aspect::MaterialModel::ViscoPlastic< dim >::thermal_conductivities
private

Definition at line 256 of file visco_plastic.h.

§ n_phase_transitions_for_each_chemical_composition

template<int dim>
std::vector<unsigned int> aspect::MaterialModel::ViscoPlastic< dim >::n_phase_transitions_for_each_chemical_composition
private

Number of phase transitions for each chemical composition (including the background field).

Definition at line 261 of file visco_plastic.h.

§ n_phases

template<int dim>
unsigned int aspect::MaterialModel::ViscoPlastic< dim >::n_phases
private

Total number of phases.

Definition at line 266 of file visco_plastic.h.

§ equation_of_state

template<int dim>
EquationOfState::MulticomponentIncompressible<dim> aspect::MaterialModel::ViscoPlastic< dim >::equation_of_state
private

Object for computing the equation of state.

Definition at line 271 of file visco_plastic.h.

§ phase_function

template<int dim>
MaterialUtilities::PhaseFunction<dim> aspect::MaterialModel::ViscoPlastic< dim >::phase_function
private

Object that handles phase transitions.

Definition at line 276 of file visco_plastic.h.


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