![]() |
ASPECT
|
Public Member Functions | |
PrescribedPlasticDilation (const unsigned int n_points) | |
std::vector< double > | get_nth_output (const unsigned int idx) const override |
![]() | |
NamedAdditionalMaterialOutputs (const std::vector< std::string > &output_names) | |
NamedAdditionalMaterialOutputs (const std::vector< std::string > &output_names, const unsigned int n_points) | |
~NamedAdditionalMaterialOutputs () override | |
const std::vector< std::string > & | get_names () const |
void | average (const MaterialAveraging::AveragingOperation, const FullMatrix< double > &, const FullMatrix< double > &) override |
![]() | |
virtual | ~AdditionalMaterialOutputs ()=default |
Public Attributes | |
std::vector< double > | dilation_lhs_term |
std::vector< double > | dilation_rhs_term |
![]() | |
std::vector< std::vector< double > > | output_values |
An AdditionalOutput that allows prescribing a dilation applied to the Stokes solution.
This is typically used in a MaterialModel to add dilation when plastic failure occurs, but can be used for other forms of dilation as well. When plastic dilation is included, a term \(\bar\alpha\gamma\) should be added to the right-hand side of the mass conservation equation, where \(\bar\alpha\) is the negative derivative of plastic potential with respect to the pressure ( \(\sin\psi\) in 2D case), and \(\gamma\) is the plastic multiplier. The plastic multiplier is given by \(\gamma = (\tau_{II} - \alpha p - k) / \eta^{ve}\), where \(\tau_{II}\) is the second invariant of the deviatoric stress, \(\alpha\) is the negative derivative of yield function with respect to the pressure ( \(\sin\phi\) in 2D case), \(k\) is cohesion, and \(\eta^{ve}\), is the pre-yielding viscosity. When the Picard method or Defect Correction Method is applied, the term \(\bar\alpha\gamma\) should be split into two terms: \(\bar\alpha\gamma = \bar\alpha\alpha p / \eta^{ve} + \bar\alpha(\tau_{II} - k) / \eta^{ve}\), the former of which should be moved to the left-hand side in order to guarantee the stability of the nonlinear solver. Therefore, this output provides two terms: dilation_lhs_term corresponds to \(\bar\alpha\alpha / \eta^{ve}\) (p is replaced by the shape function), and dilation_rhs_term corresponds to \((\tau_{II} - k) / \eta^{ve}\).
Definition at line 1230 of file interface.h.
|
explicit |
Constructor
|
overridevirtual |
Function for NamedAdditionalMaterialOutputs interface
Reimplemented from aspect::MaterialModel::NamedAdditionalMaterialOutputs< dim >.
std::vector<double> aspect::MaterialModel::PrescribedPlasticDilation< dim >::dilation_lhs_term |
A scalar value per evaluation point corresponding to the LHS term due to plastic dilation.
Definition at line 1247 of file interface.h.
std::vector<double> aspect::MaterialModel::PrescribedPlasticDilation< dim >::dilation_rhs_term |
A scalar value per evaluation point corresponding to the RHS term due to plastic dilation.
Definition at line 1253 of file interface.h.