![]() |
OGS
|
|
Classes | |
| class | AT_DegradationDerivative |
| class | COHESIVE_DegradationDerivative |
| class | DegradationDerivative |
Enumerations | |
| enum class | PhaseFieldModel { AT1 , AT2 , COHESIVE } |
| enum class | SofteningCurve { Linear , Exponential } |
| enum class | EnergySplitModel { Isotropic , VolDev , Spectral , EffectiveStress , OrthoVolDev , OrthoMasonry } |
Functions | |
| double | evaluateHTensSpectral (int const i, int const j, Eigen::Matrix< double, 3, 1 > const &principal_strain) |
| H terms in the Spectral decomposition: | |
| double | evaluateHCompSpectral (int const i, int const j, Eigen::Matrix< double, 3, 1 > const &principal_strain) |
| template<> | |
| MathLib::KelvinVector::KelvinMatrixType< 3 > | aOdotB< 3 > (MathLib::KelvinVector::KelvinVectorType< 3 > const &A, MathLib::KelvinVector::KelvinVectorType< 3 > const &B) |
| template<> | |
| MathLib::KelvinVector::KelvinMatrixType< 2 > | aOdotB< 2 > (MathLib::KelvinVector::KelvinVectorType< 2 > const &A, MathLib::KelvinVector::KelvinVectorType< 2 > const &B) |
| template<int DisplacementDim> | |
| MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > | aOdotB (MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &A, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &B) |
| double | heaviside (double const v) |
| double | macaulayTensile (double const v) |
| double | macaulayCompressive (double v) |
| template<int DisplacementDim> | |
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > | calculateVolDevDegradedStress (double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps) |
| template<int DisplacementDim> | |
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > | calculateSpectralDegradedStress (double const degradation, double const lambda, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps) |
| template<int DisplacementDim> | |
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > | calculateIsotropicDegradedStress (double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps) |
| template<int DisplacementDim> | |
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > | calculateIsotropicDegradedStressWithRankineEnergy (double const degradation, double const bulk_modulus, double const mu, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps) |
| template<int DisplacementDim> | |
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > | calculateOrthoVolDevDegradedStress (double const degradation, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const &C_ortho) |
| template<int DisplacementDim> | |
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > | calculateOrthoMasonryDegradedStress (double const degradation, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const &C_ortho) |
| template<int DisplacementDim> | |
| PhaseFieldModel | convertStringToPhaseFieldModel (std::string const &phasefield_model) |
| template<int DisplacementDim> | |
| SofteningCurve | convertStringToSofteningCurve (std::optional< std::string > const &softening_curve) |
| template<int DisplacementDim> | |
| EnergySplitModel | convertStringToEnergySplitModel (std::string const &energy_split_model) |
| template<int DisplacementDim> | |
| std::unique_ptr< DegradationDerivative > | creatDegradationDerivative (PhaseFieldModel const &phasefield_model, double const lch, SofteningCurve const &softening_curve) |
| template<typename T_DNDX, typename T_N, typename T_W, typename T_D, typename T_LOCAL_JAC, typename T_LOCAL_RHS> | |
| void | calculateCrackLocalJacobianAndResidual (T_DNDX &dNdx, T_N &N, T_W &w, T_D &d, T_LOCAL_JAC &local_Jac, T_LOCAL_RHS local_rhs, double const gc, double const ls, PhaseFieldModel &phasefield_model) |
| template<typename T_VECTOR, typename T_MATRIX, int DisplacementDim> | |
| void | calculateStress (T_VECTOR &sigma, T_VECTOR &sigma_tensile, T_VECTOR &sigma_compressive, T_VECTOR &eps_tensile, T_MATRIX &D, T_MATRIX &C_tensile, T_MATRIX &C_compressive, double &strain_energy_tensile, double &elastic_energy, double const degradation, T_VECTOR const &eps, EnergySplitModel const &energy_split_model, double const t, ParameterLib::SpatialPosition const &x, MaterialLib::Solids::MechanicsBase< DisplacementDim > const &solid_material) |
|
strong |
| Enumerator | |
|---|---|
| Isotropic | |
| VolDev | |
| Spectral | |
| EffectiveStress | |
| OrthoVolDev | |
| OrthoMasonry | |
Definition at line 85 of file PhaseFieldBase.h.
|
strong |
| Enumerator | |
|---|---|
| AT1 | |
| AT2 | |
| COHESIVE | |
Definition at line 28 of file PhaseFieldBase.h.
|
strong |
| Enumerator | |
|---|---|
| Linear | |
| Exponential | |
Definition at line 57 of file PhaseFieldBase.h.
| MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > MaterialLib::Solids::Phasefield::aOdotB | ( | MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & | A, |
| MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & | B ) |
Decompose the stiffness into tensile and compressive part. Judging by the physical observations, compression perpendicular to a crack does not cause crack propagation. Thus, the phase-field parameter is only involved into the tensile part to degrade the elastic strain energy.
Referenced by calculateSpectralDegradedStress().
| MathLib::KelvinVector::KelvinMatrixType< 2 > MaterialLib::Solids::Phasefield::aOdotB< 2 > | ( | MathLib::KelvinVector::KelvinVectorType< 2 > const & | A, |
| MathLib::KelvinVector::KelvinVectorType< 2 > const & | B ) |
Definition at line 33 of file LinearElasticIsotropicPhaseField.cpp.
References heaviside(), and macaulayCompressive().
| MathLib::KelvinVector::KelvinMatrixType< 3 > MaterialLib::Solids::Phasefield::aOdotB< 3 > | ( | MathLib::KelvinVector::KelvinVectorType< 3 > const & | A, |
| MathLib::KelvinVector::KelvinVectorType< 3 > const & | B ) |
Double–minor symmetrized tensor product
aOdotB computes the symmetric 4th-order tensor product
\[(A \odot B)_{ijkl} = \tfrac{1}{4} (A_{ik}B_{jl} + A_{il}B_{jk} + A_{jk}B_{il} + A_{jl}B_{ik}), \]
which is symmetric in both index pairs \((i,j)\) and \((k,l)\) (“double-minor symmetry”).
Reference Mehrabadi, M.M., & Cowin, S.C. (1990). Eigentensors of linear anisotropic elastic materials. Q. J. Mech. Appl. Math., 43(1), 15–41.
Definition at line 33 of file LinearElasticIsotropicPhaseField.cpp.
| void MaterialLib::Solids::Phasefield::calculateCrackLocalJacobianAndResidual | ( | T_DNDX & | dNdx, |
| T_N & | N, | ||
| T_W & | w, | ||
| T_D & | d, | ||
| T_LOCAL_JAC & | local_Jac, | ||
| T_LOCAL_RHS | local_rhs, | ||
| double const | gc, | ||
| double const | ls, | ||
| PhaseFieldModel & | phasefield_model ) |
Definition at line 297 of file PhaseFieldBase.h.
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateIsotropicDegradedStress | ( | double const | degradation, |
| double const | bulk_modulus, | ||
| double const | mu, | ||
| MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & | eps ) |
Definition at line 269 of file LinearElasticIsotropicPhaseField.h.
References MathLib::KelvinVector::kelvin_vector_dimensions().
Referenced by calculateStress().
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateIsotropicDegradedStressWithRankineEnergy | ( | double const | degradation, |
| double const | bulk_modulus, | ||
| double const | mu, | ||
| MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & | eps ) |
Definition at line 322 of file LinearElasticIsotropicPhaseField.h.
References heaviside(), MathLib::KelvinVector::kelvin_vector_dimensions(), and MathLib::KelvinVector::kelvinVectorToTensor().
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateOrthoMasonryDegradedStress | ( | double const | degradation, |
| MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & | eps, | ||
| MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const & | C_ortho ) |
Definition at line 127 of file LinearElasticOrthotropicPhaseField.h.
References heaviside(), MathLib::KelvinVector::kelvinVectorToTensor(), macaulayCompressive(), macaulayTensile(), and MathLib::KelvinVector::tensorToKelvin().
Referenced by calculateStress().
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateOrthoVolDevDegradedStress | ( | double const | degradation, |
| MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & | eps, | ||
| MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > const & | C_ortho ) |
Definition at line 38 of file LinearElasticOrthotropicPhaseField.h.
References MathLib::KelvinVector::kelvin_vector_dimensions().
Referenced by calculateStress().
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateSpectralDegradedStress | ( | double const | degradation, |
| double const | lambda, | ||
| double const | mu, | ||
| MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & | eps ) |
Definition at line 152 of file LinearElasticIsotropicPhaseField.h.
References aOdotB(), evaluateHCompSpectral(), evaluateHTensSpectral(), heaviside(), MathLib::KelvinVector::kelvin_vector_dimensions(), MathLib::KelvinVector::kelvinVectorToTensor(), macaulayCompressive(), macaulayTensile(), and MathLib::KelvinVector::tensorToKelvin().
Referenced by calculateStress().
| void MaterialLib::Solids::Phasefield::calculateStress | ( | T_VECTOR & | sigma, |
| T_VECTOR & | sigma_tensile, | ||
| T_VECTOR & | sigma_compressive, | ||
| T_VECTOR & | eps_tensile, | ||
| T_MATRIX & | D, | ||
| T_MATRIX & | C_tensile, | ||
| T_MATRIX & | C_compressive, | ||
| double & | strain_energy_tensile, | ||
| double & | elastic_energy, | ||
| double const | degradation, | ||
| T_VECTOR const & | eps, | ||
| EnergySplitModel const & | energy_split_model, | ||
| double const | t, | ||
| ParameterLib::SpatialPosition const & | x, | ||
| MaterialLib::Solids::MechanicsBase< DisplacementDim > const & | solid_material ) |
Definition at line 347 of file PhaseFieldBase.h.
References calculateIsotropicDegradedStress(), calculateOrthoMasonryDegradedStress(), calculateOrthoVolDevDegradedStress(), calculateSpectralDegradedStress(), calculateVolDevDegradedStress(), EffectiveStress, Isotropic, OGS_FATAL, OrthoMasonry, OrthoVolDev, Spectral, and VolDev.
Referenced by ProcessLib::HMPhaseField::IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >::updateConstitutiveRelation(), and ProcessLib::PhaseField::IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >::updateConstitutiveRelation().
| std::tuple< MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinVectorType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, double, double, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim >, MathLib::KelvinVector::KelvinMatrixType< DisplacementDim > > MaterialLib::Solids::Phasefield::calculateVolDevDegradedStress | ( | double const | degradation, |
| double const | bulk_modulus, | ||
| double const | mu, | ||
| MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const & | eps ) |
Definition at line 75 of file LinearElasticIsotropicPhaseField.h.
References heaviside(), MathLib::KelvinVector::kelvin_vector_dimensions(), macaulayCompressive(), and macaulayTensile().
Referenced by calculateStress().
| EnergySplitModel MaterialLib::Solids::Phasefield::convertStringToEnergySplitModel | ( | std::string const & | energy_split_model | ) |
Definition at line 96 of file PhaseFieldBase.h.
References EffectiveStress, Isotropic, OGS_FATAL, OrthoMasonry, OrthoVolDev, Spectral, and VolDev.
Referenced by ProcessLib::HMPhaseField::createHMPhaseFieldProcess(), and ProcessLib::PhaseField::createPhaseFieldProcess().
| PhaseFieldModel MaterialLib::Solids::Phasefield::convertStringToPhaseFieldModel | ( | std::string const & | phasefield_model | ) |
Definition at line 36 of file PhaseFieldBase.h.
References AT1, AT2, COHESIVE, and OGS_FATAL.
Referenced by ProcessLib::HMPhaseField::createHMPhaseFieldProcess(), and ProcessLib::PhaseField::createPhaseFieldProcess().
| SofteningCurve MaterialLib::Solids::Phasefield::convertStringToSofteningCurve | ( | std::optional< std::string > const & | softening_curve | ) |
Definition at line 64 of file PhaseFieldBase.h.
References Exponential, Linear, and OGS_FATAL.
Referenced by ProcessLib::HMPhaseField::createHMPhaseFieldProcess(), and ProcessLib::PhaseField::createPhaseFieldProcess().
| std::unique_ptr< DegradationDerivative > MaterialLib::Solids::Phasefield::creatDegradationDerivative | ( | PhaseFieldModel const & | phasefield_model, |
| double const | lch, | ||
| SofteningCurve const & | softening_curve ) |
Definition at line 281 of file PhaseFieldBase.h.
References COHESIVE.
| double MaterialLib::Solids::Phasefield::evaluateHCompSpectral | ( | int const | i, |
| int const | j, | ||
| Eigen::Matrix< double, 3, 1 > const & | principal_strain ) |
Definition at line 33 of file LinearElasticIsotropicPhaseField.cpp.
Referenced by calculateSpectralDegradedStress().
| double MaterialLib::Solids::Phasefield::evaluateHTensSpectral | ( | int const | i, |
| int const | j, | ||
| Eigen::Matrix< double, 3, 1 > const & | principal_strain ) |
H terms in the Spectral decomposition:
Definition at line 14 of file LinearElasticIsotropicPhaseField.cpp.
References heaviside(), and macaulayTensile().
Referenced by calculateSpectralDegradedStress().
|
inline |
heaviside function returns 1.0 if the argument is positive and 0.0 if negative
Definition at line 37 of file LinearElasticIsotropicPhaseField.h.
Referenced by aOdotB< 2 >(), calculateIsotropicDegradedStressWithRankineEnergy(), calculateOrthoMasonryDegradedStress(), calculateSpectralDegradedStress(), calculateVolDevDegradedStress(), evaluateHTensSpectral(), macaulayCompressive(), and macaulayTensile().
|
inline |
Definition at line 48 of file LinearElasticIsotropicPhaseField.h.
References heaviside().
Referenced by aOdotB< 2 >(), calculateOrthoMasonryDegradedStress(), calculateSpectralDegradedStress(), and calculateVolDevDegradedStress().
|
inline |
Macaulay brackets: positive strain is tensile and negative strain for compressive
Definition at line 44 of file LinearElasticIsotropicPhaseField.h.
References heaviside().
Referenced by calculateOrthoMasonryDegradedStress(), calculateSpectralDegradedStress(), calculateVolDevDegradedStress(), and evaluateHTensSpectral().