![]() |
ASPECT
|
Public Member Functions | |
CrystalPreferredOrientation () | |
virtual void | initialize () override |
virtual void | initialize_one_particle_property (const Point< dim > &position, std::vector< double > &particle_properties) const override |
virtual void | update_one_particle_property (const unsigned int data_position, const Point< dim > &position, const Vector< double > &solution, const std::vector< Tensor< 1, dim >> &gradients, const ArrayView< double > &particle_properties) const override |
UpdateTimeFlags | need_update () const override |
InitializationModeForLateParticles | late_initialization_mode () const override |
virtual UpdateFlags | get_needed_update_flags () const override |
virtual std::vector< std::pair< std::string, unsigned int > > | get_property_information () const override |
std::pair< std::vector< double >, std::vector< Tensor< 2, 3 > > > | compute_derivatives (const unsigned int cpo_index, const ArrayView< double > &data, const unsigned int mineral_i, const SymmetricTensor< 2, 3 > &strain_rate_3d, const Tensor< 2, 3 > &velocity_gradient_tensor, const Point< dim > &position, const double temperature, const double pressure, const Tensor< 1, dim > &velocity, const std::vector< double > &compositions, const SymmetricTensor< 2, dim > &strain_rate, const SymmetricTensor< 2, dim > &deviatoric_strain_rate, const double water_content) const |
Computes the volume fraction and grain orientation derivatives of all the grains of a mineral. More... | |
std::pair< std::vector< double >, std::vector< Tensor< 2, 3 > > > | compute_derivatives_drex_2004 (const unsigned int cpo_index, const ArrayView< double > &data, const unsigned int mineral_i, const SymmetricTensor< 2, 3 > &strain_rate_3d, const Tensor< 2, 3 > &velocity_gradient_tensor, const std::array< double, 4 > ref_resolved_shear_stress, const bool prevent_nondimensionalization=false) const |
Computes the CPO derivatives with the D-Rex 2004 algorithm. More... | |
void | parse_parameters (ParameterHandler &prm) override |
unsigned int | get_number_of_grains () const |
unsigned int | get_number_of_minerals () const |
DeformationType | determine_deformation_type (const DeformationTypeSelector deformation_type_selector, const Point< dim > &position, const double temperature, const double pressure, const Tensor< 1, dim > &velocity, const std::vector< double > &compositions, const SymmetricTensor< 2, dim > &strain_rate, const SymmetricTensor< 2, dim > &deviatoric_strain_rate, const double water_content) const |
Determines the deformation type from the deformation type selector. If the provided deformation_type_selector is a specific deformation type, the function will simply return the corresponding deformation type. However, if the deformation_type_selector is an algorithm to determine the current deformation type (e.g. based on measured lab data or analytical models), then the function computes the appropriate deformation type at the given conditions and returns the compute deformation type. More... | |
DeformationType | determine_deformation_type_karato_2008 (const double stress, const double water_content) const |
Computes the deformation type given the stress and water content according to the table in Karato 2008. More... | |
std::array< double, 4 > | reference_resolved_shear_stress_from_deformation_type (DeformationType deformation_type, double max_value=1e60) const |
Computes the reference resolved shear stress (RRSS) based on the selected deformation type. More... | |
double | get_deformation_type (const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i) const |
Returns the value in the data array representing the deformation type. More... | |
void | set_deformation_type (const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const double deformation_type) const |
Sets the value in the data array representing the deformation type. More... | |
double | get_volume_fraction_mineral (const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i) const |
Returns the value in the data array representing the volume fraction of a mineral. More... | |
void | set_volume_fraction_mineral (const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const double volume_fraction_mineral) const |
Sets the value in the data array representing the volume fraction of a mineral. More... | |
double | get_volume_fractions_grains (const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const unsigned int grain_i) const |
Returns the value in the data array representing the volume fraction of a grain. More... | |
void | set_volume_fractions_grains (const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const unsigned int grain_i, const double volume_fractions_grains) const |
Sets the value in the data array representing the volume fraction of a grain. More... | |
Tensor< 2, 3 > | get_rotation_matrix_grains (const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const unsigned int grain_i) const |
Gets the rotation matrix for a grain in a mineral. More... | |
void | set_rotation_matrix_grains (const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const unsigned int grain_i, const Tensor< 2, 3 > &rotation_matrix) const |
Sets the rotation matrix for a grain in a mineral. More... | |
![]() | |
virtual | ~Interface ()=default |
virtual void | update_particle_property (const unsigned int data_position, const Vector< double > &solution, const std::vector< Tensor< 1, dim >> &gradients, typename ParticleHandler< dim >::particle_iterator &particle) 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 |
Static Public Member Functions | |
static void | declare_parameters (ParameterHandler &prm) |
![]() | |
static void | declare_parameters (ParameterHandler &prm) |
![]() | |
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 | |
void | compute_random_rotation_matrix (Tensor< 2, 3 > &rotation_matrix) const |
double | advect_forward_euler (const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const double dt, const std::pair< std::vector< double >, std::vector< Tensor< 2, 3 >>> &derivatives) const |
Updates the volume fractions and rotation matrices with a Forward Euler scheme. More... | |
double | advect_backward_euler (const unsigned int cpo_data_position, const ArrayView< double > &data, const unsigned int mineral_i, const double dt, const std::pair< std::vector< double >, std::vector< Tensor< 2, 3 >>> &derivatives) const |
Updates the volume fractions and rotation matrices with a Backward Euler scheme. More... | |
std::pair< std::vector< double >, std::vector< Tensor< 2, 3 > > > | compute_derivatives_spin_tensor (const Tensor< 2, 3 > &velocity_gradient_tensor) const |
Private Attributes | |
boost::mt19937 | random_number_generator |
unsigned int | random_number_seed |
unsigned int | n_grains |
unsigned int | n_minerals |
unsigned int | water_index |
std::vector< DeformationTypeSelector > | deformation_type_selector |
std::vector< double > | volume_fractions_minerals |
AdvectionMethod | advection_method |
CPODerivativeAlgorithm | cpo_derivative_algorithm |
double | property_advection_tolerance |
unsigned int | property_advection_max_iterations |
Tensor< 3, 3 > | permutation_operator_3d |
D-Rex variables | |
double | stress_exponent |
double | nucleation_efficiency |
double | exponent_p |
double | threshold_GBS |
double | mobility |
The plugin manages and computes the evolution of Lattice/Crystal Preferred Orientations (LPO/CPO) on particles. Each ASPECT particle represents many grains. Each grain is assigned a size and a orientation matrix. This allows tracking the LPO evolution with kinematic polycrystal CPO evolution models such as D-Rex (Kaminski and Ribe, 2001; ́Kaminski et al., 2004).
This plugin stores M minerals and for each mineral it stores N grains. The total amount of memory stored is 2 doubles per mineral plus 10 doubles per grain, resulting in a total memory of M * (2 + N * 10) . The layout of the data stored is the following (note that for this plugin the following dims are always 3):
Last used data entry is n_minerals * (n_grains * 10 + 2).
We store the same number of grains for all minerals (e.g. olivine and enstatite grains), although their volume fractions may not be the same. This is because we need a minimum number of grains per tracer to perform reliable statistics on it. This minimum should be the same for all minerals.
Definition at line 119 of file crystal_preferred_orientation.h.
aspect::Particle::Property::CrystalPreferredOrientation< dim >::CrystalPreferredOrientation | ( | ) |
Constructor
|
overridevirtual |
Initialization function. This function is called once at the beginning of the program after parse_parameters is run.
Reimplemented from aspect::Particle::Property::Interface< dim >.
|
overridevirtual |
Initialization function. This function is called once at the creation of every particle for every property to initialize its value.
[in] | position | The current particle position. |
[in,out] | particle_properties | The properties of the particle that is initialized within the call of this function. The purpose of this function should be to extend this vector by a number of properties. |
Reimplemented from aspect::Particle::Property::Interface< dim >.
|
overridevirtual |
Update function. This function is called every time an update is request by need_update() for every particle for every property.
[in] | data_position | An unsigned integer that denotes which component of the particle property vector is associated with the current property. For properties that own several components it denotes the first component of this property, all other components fill consecutive entries in the particle_properties vector. |
[in] | position | The current particle position. |
[in] | solution | The values of the solution variables at the current particle position. |
[in] | gradients | The gradients of the solution variables at the current particle position. |
[in,out] | particle_properties | The properties of the particle that is updated within the call of this function. |
Reimplemented from aspect::Particle::Property::Interface< dim >.
|
overridevirtual |
This implementation tells the particle manager that we need to update particle properties every time step.
Reimplemented from aspect::Particle::Property::Interface< dim >.
|
overridevirtual |
The CPO of late particles is initialized by interpolating from existing particles.
Reimplemented from aspect::Particle::Property::Interface< dim >.
|
overridevirtual |
Return which data has to be provided to update the property. The integrated strains needs the gradients of the velocity.
Reimplemented from aspect::Particle::Property::Interface< dim >.
|
overridevirtual |
Set up the information about the names and number of components this property requires.
Implements aspect::Particle::Property::Interface< dim >.
std::pair<std::vector<double>, std::vector<Tensor<2,3> > > aspect::Particle::Property::CrystalPreferredOrientation< dim >::compute_derivatives | ( | const unsigned int | cpo_index, |
const ArrayView< double > & | data, | ||
const unsigned int | mineral_i, | ||
const SymmetricTensor< 2, 3 > & | strain_rate_3d, | ||
const Tensor< 2, 3 > & | velocity_gradient_tensor, | ||
const Point< dim > & | position, | ||
const double | temperature, | ||
const double | pressure, | ||
const Tensor< 1, dim > & | velocity, | ||
const std::vector< double > & | compositions, | ||
const SymmetricTensor< 2, dim > & | strain_rate, | ||
const SymmetricTensor< 2, dim > & | deviatoric_strain_rate, | ||
const double | water_content | ||
) | const |
Computes the volume fraction and grain orientation derivatives of all the grains of a mineral.
cpo_index | The location where the CPO data starts in the data array. |
data | The data array containing the CPO data. |
mineral_i | The mineral for which to compute the derivatives for. |
strain_rate_3d | The 3D strain rate at the location where the derivative is requested. |
velocity_gradient_tensor | The velocity gradient tensor at the location where the derivative is requested. |
position | the location for which the derivative is requested. |
temperature | The temperature at the location where the derivative is requested. |
pressure | The pressure at the location where the derivative is requested. |
velocity | The veloicty at the location where the derivative is requested. |
compositions | The compositios at the location where the derivative is requested. |
strain_rate | The strain-rate at the location where the derivative is requested. |
deviatoric_strain_rate | The deviatoric strain-rate at the location where the derivative is requested. |
water_content | The water content at the location where the derivative is requested. |
std::pair<std::vector<double>, std::vector<Tensor<2,3> > > aspect::Particle::Property::CrystalPreferredOrientation< dim >::compute_derivatives_drex_2004 | ( | const unsigned int | cpo_index, |
const ArrayView< double > & | data, | ||
const unsigned int | mineral_i, | ||
const SymmetricTensor< 2, 3 > & | strain_rate_3d, | ||
const Tensor< 2, 3 > & | velocity_gradient_tensor, | ||
const std::array< double, 4 > | ref_resolved_shear_stress, | ||
const bool | prevent_nondimensionalization = false |
||
) | const |
Computes the CPO derivatives with the D-Rex 2004 algorithm.
cpo_index | The location where the CPO data starts in the data array. |
data | The data array containing the CPO data. |
mineral_i | The mineral for which to compute the derivatives for. |
strain_rate_3d | The 3D strain rate |
velocity_gradient_tensor | The velocity gradient tensor |
ref_resolved_shear_stress | Represent one value per slip plane. The planes are ordered from weakest to strongest with relative values, where the inactive plane is infinity strong. So it is a measure of strength on each slip plane. |
prevent_nondimensionalization | Prevent nondimensializing values internally. Only for unit testing purposes. |
|
static |
Declare the parameters this class takes through input files.
|
overridevirtual |
Read the parameters this class declares from the parameter file.
Reimplemented from aspect::Particle::Property::Interface< dim >.
unsigned int aspect::Particle::Property::CrystalPreferredOrientation< dim >::get_number_of_grains | ( | ) | const |
Return the number of grains per particle
unsigned int aspect::Particle::Property::CrystalPreferredOrientation< dim >::get_number_of_minerals | ( | ) | const |
Return the number of minerals per particle
DeformationType aspect::Particle::Property::CrystalPreferredOrientation< dim >::determine_deformation_type | ( | const DeformationTypeSelector | deformation_type_selector, |
const Point< dim > & | position, | ||
const double | temperature, | ||
const double | pressure, | ||
const Tensor< 1, dim > & | velocity, | ||
const std::vector< double > & | compositions, | ||
const SymmetricTensor< 2, dim > & | strain_rate, | ||
const SymmetricTensor< 2, dim > & | deviatoric_strain_rate, | ||
const double | water_content | ||
) | const |
Determines the deformation type from the deformation type selector. If the provided deformation_type_selector
is a specific deformation type, the function will simply return the corresponding deformation type. However, if the deformation_type_selector
is an algorithm to determine the current deformation type (e.g. based on measured lab data or analytical models), then the function computes the appropriate deformation type at the given conditions and returns the compute deformation type.
DeformationType aspect::Particle::Property::CrystalPreferredOrientation< dim >::determine_deformation_type_karato_2008 | ( | const double | stress, |
const double | water_content | ||
) | const |
Computes the deformation type given the stress and water content according to the table in Karato 2008.
std::array<double,4> aspect::Particle::Property::CrystalPreferredOrientation< dim >::reference_resolved_shear_stress_from_deformation_type | ( | DeformationType | deformation_type, |
double | max_value = 1e60 |
||
) | const |
Computes the reference resolved shear stress (RRSS) based on the selected deformation type.
The inactive plane should theoretically be infinitely strong, but this is nummerically not desirable, so an optional max_value can be set to indicate an inactive plane.
It is currently designed to return the relative strength of the slip planes for olivine, which are are 4, but this could be generalized.
|
inline |
Returns the value in the data array representing the deformation type.
cpo_data_position | The starting index/position of the cpo data in the particle data vector. |
data | The particle data vector. |
mineral_i | The mineral to get the value of the deformation type for. |
Definition at line 343 of file crystal_preferred_orientation.h.
|
inline |
Sets the value in the data array representing the deformation type.
cpo_data_position | The starting index/position of the cpo data in the particle data vector. |
data | The particle data vector. |
mineral_i | The mineral to set the value deformation type for. |
deformation_type | The value of the of the deformation type to set. |
Definition at line 359 of file crystal_preferred_orientation.h.
|
inline |
Returns the value in the data array representing the volume fraction of a mineral.
cpo_data_position | The starting index/position of the cpo data in the particle data vector. |
data | The particle data vector. |
mineral_i | The mineral to get the value of the volume fraction of a mineral for. |
Definition at line 375 of file crystal_preferred_orientation.h.
|
inline |
Sets the value in the data array representing the volume fraction of a mineral.
cpo_data_position | The starting index/position of the cpo data in the particle data vector. |
data | The particle data vector. |
mineral_i | The mineral to set the value of the volume fraction of a mineral for. |
volume_fraction_mineral | The value of the of the volume fraction of a mineral to set. |
Definition at line 391 of file crystal_preferred_orientation.h.
|
inline |
Returns the value in the data array representing the volume fraction of a grain.
cpo_data_position | The starting index/position of the cpo data in the particle data vector. |
data | The particle data vector. |
mineral_i | The mineral to get the value of the volume fraction of a grain for. |
grain_i | The grain to get the value of the volume fraction of. |
Definition at line 408 of file crystal_preferred_orientation.h.
|
inline |
Sets the value in the data array representing the volume fraction of a grain.
cpo_data_position | The starting index/position of the cpo data in the particle data vector. |
data | The particle data vector. |
mineral_i | The mineral to set the value of the volume fraction of a grain for. |
grain_i | The grain to set the value of the volume fraction of. |
volume_fractions_grains | The value of the of the volume fraction of a grain to set. |
Definition at line 426 of file crystal_preferred_orientation.h.
|
inline |
Gets the rotation matrix for a grain in a mineral.
cpo_data_position | The starting index/position of the cpo data in the particle data vector. |
data | The particle data vector. |
mineral_i | The mineral to get the value of the rotation matrix of a grain for. |
grain_i | The grain to get the value of the rotation matrix of. |
Definition at line 445 of file crystal_preferred_orientation.h.
|
inline |
Sets the rotation matrix for a grain in a mineral.
cpo_data_position | The starting index/position of the cpo data in the particle data vector. |
data | The particle data vector. |
mineral_i | The mineral to set the value of the rotation matrix of a grain for. |
grain_i | The grain to get the value of the rotation matrix of. |
rotation_matrix | The rotation matrix to set for the grain in the mineral. |
Definition at line 469 of file crystal_preferred_orientation.h.
|
private |
Computes a random rotation matrix.
|
private |
Updates the volume fractions and rotation matrices with a Forward Euler scheme.
Updates the volume fractions and rotation matrices with a Forward Euler scheme: \(x_t = x_{t-1} + dt * x_{t-1} * \frac{dx_t}{dt}\). The function returns the sum of the new volume fractions.
cpo_data_position | The starting index/position of the cpo data in the particle data vector. |
data | The particle data vector. |
mineral_i | Which mineral to advect for. |
dt | The time step used for the advection step |
derivatives | A pair containing the derivatives for the volume fractions and orientations respectively. |
|
private |
Updates the volume fractions and rotation matrices with a Backward Euler scheme.
Updates the volume fractions and rotation matrices with a Backward Euler scheme: \(x_t = x_{t-1} + dt * x_{t-1} * \frac{dx_t}{dt}\). The function returns the sum of the new volume fractions.
cpo_data_position | The starting index/position of the cpo data in the particle data vector. |
data | The particle data vector. |
mineral_i | Which mineral to advect for. |
dt | The time step used for the advection step |
derivatives | A pair containing the derivatives for the volume fractions and orientations respectively. |
|
private |
Computes and returns the volume fraction and grain orientation derivatives such that the grains stay the same size and the orientations rotating passively with the particle.
velocity_gradient_tensor | is the velocity gradient tensor at the location of the particle. |
|
mutableprivate |
Random number generator used for initialization of particles
Definition at line 547 of file crystal_preferred_orientation.h.
|
private |
Definition at line 548 of file crystal_preferred_orientation.h.
|
private |
Definition at line 550 of file crystal_preferred_orientation.h.
|
private |
Definition at line 552 of file crystal_preferred_orientation.h.
|
private |
The index of the water composition.
Definition at line 557 of file crystal_preferred_orientation.h.
|
private |
A vector containing the deformation type selectors provided by the user. Should be one of the following: "Olivine: Karato 2008", "Olivine: A-fabric", "Olivine: B-fabric", "Olivine: C-fabric", "Olivine: D-fabric", "Olivine: E-fabric", "Enstatite" or "Passive".
Definition at line 565 of file crystal_preferred_orientation.h.
|
private |
Store the volume fraction for each mineral.
Definition at line 570 of file crystal_preferred_orientation.h.
|
private |
Advection method for particle properties
Definition at line 575 of file crystal_preferred_orientation.h.
|
private |
What algorithm to use to compute the derivatives
Definition at line 580 of file crystal_preferred_orientation.h.
|
private |
This value determines the tolerance used for the Backward Euler and Crank-Nicolson iterations.
Definition at line 586 of file crystal_preferred_orientation.h.
|
private |
This value determines the the maximum number of iterations used for the Backward Euler and Crank-Nicolson iterations.
Definition at line 592 of file crystal_preferred_orientation.h.
|
private |
The tensor representation of the permutation symbol.
Definition at line 597 of file crystal_preferred_orientation.h.
|
private |
Stress exponent
Definition at line 606 of file crystal_preferred_orientation.h.
|
private |
efficiency of nucleation parameter. lambda_m in equation 8 of Kaminski et al. (2004, Geophys. J. Int)
Definition at line 612 of file crystal_preferred_orientation.h.
|
private |
An exponent described in equation 10 of Kaminski and Ribe (2001, EPSL)
Definition at line 617 of file crystal_preferred_orientation.h.
|
private |
The Dimensionless Grain Boundary Sliding (GBS) threshold. This is a grain size threshold below which grain deform by GBS and become strain-free grains.
Definition at line 624 of file crystal_preferred_orientation.h.
|
private |
Dimensionless grain boundary mobility as described by equation 14 in Kaminski and Ribe (2001, EPSL).
Definition at line 630 of file crystal_preferred_orientation.h.