void MisesMat :: give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) { this->giveLinearElasticMaterial()->give1dStressStiffMtrx(answer, mode, gp, tStep); MisesMatStatus *status = static_cast< MisesMatStatus * >( this->giveStatus(gp) ); double kappa = status->giveCumulativePlasticStrain(); // increment of cumulative plastic strain as an indicator of plastic loading double tempKappa = status->giveTempCumulativePlasticStrain(); double omega = status->giveTempDamage(); double E = answer.at(1, 1); if ( mode != TangentStiffness ) { return; } if ( tempKappa <= kappa ) { // elastic loading - elastic stiffness plays the role of tangent stiffness answer.times(1 - omega); return; } // === plastic loading === const FloatArray &stressVector = status->giveTempEffectiveStress(); double stress = stressVector.at(1); answer.resize(1, 1); answer.at(1, 1) = ( 1 - omega ) * E * H / ( E + H ) - computeDamageParamPrime(tempKappa) * E / ( E + H ) * stress * signum(stress); }
// returns the consistent (algorithmic) tangent stiffness matrix void MisesMat :: give3dSSMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *atTime) { // start from the elastic stiffness this->giveLinearElasticMaterial()->give3dMaterialStiffnessMatrix(answer, mode, gp, atTime); if ( mode != TangentStiffness ) { return; } MisesMatStatus *status = static_cast< MisesMatStatus * >( this->giveStatus(gp) ); double kappa = status->giveCumulativePlasticStrain(); double tempKappa = status->giveTempCumulativePlasticStrain(); // increment of cumulative plastic strain as an indicator of plastic loading double dKappa = tempKappa - kappa; if ( dKappa <= 0.0 ) { // elastic loading - elastic stiffness plays the role of tangent stiffness return; } // === plastic loading === // yield stress at the beginning of the step double sigmaY = sig0 + H * kappa; // trial deviatoric stress and its norm StressVector trialStressDev(_3dMat); double trialStressVol; status->giveTrialStressVol(trialStressVol); status->giveTrialStressDev(trialStressDev); double trialS = trialStressDev.computeStressNorm(); // one correction term FloatMatrix stiffnessCorrection(6, 6); stiffnessCorrection.beDyadicProductOf(trialStressDev, trialStressDev); double factor = -2. * sqrt(6.) * G * G / trialS; double factor1 = factor * sigmaY / ( ( H + 3. * G ) * trialS * trialS ); stiffnessCorrection.times(factor1); answer.add(stiffnessCorrection); // another correction term stiffnessCorrection.bePinvID(); double factor2 = factor * dKappa; stiffnessCorrection.times(factor2); answer.add(stiffnessCorrection); //influence of damage // double omega = computeDamageParam(tempKappa); double omega = status->giveTempDamage(); answer.times(1. - omega); FloatArray effStress; status->giveTempEffectiveStress(effStress); double omegaPrime = computeDamageParamPrime(tempKappa); double scalar = -omegaPrime *sqrt(6.) * G / ( 3. * G + H ) / trialS; stiffnessCorrection.beDyadicProductOf(effStress, trialStressDev); stiffnessCorrection.times(scalar); answer.add(stiffnessCorrection); }
int MisesMat :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) { MisesMatStatus *status = static_cast< MisesMatStatus * >( this->giveStatus(gp) ); if ( type == IST_PlasticStrainTensor ) { answer = status->givePlasDef(); return 1; } else if ( type == IST_MaxEquivalentStrainLevel ) { answer.resize(1); answer.at(1) = status->giveCumulativePlasticStrain(); return 1; } else if ( ( type == IST_DamageScalar ) || ( type == IST_DamageTensor ) ) { answer.resize(1); answer.at(1) = status->giveDamage(); return 1; } else { return StructuralMaterial :: giveIPValue(answer, gp, type, tStep); } }
int MisesMat :: giveIPValue(FloatArray &answer, GaussPoint *aGaussPoint, InternalStateType type, TimeStep *atTime) { MisesMatStatus *status = static_cast< MisesMatStatus * >( this->giveStatus(aGaussPoint) ); if ( type == IST_PlasticStrainTensor ) { const FloatArray &ep = status->givePlasDef(); ///@todo Fix this so that it doesn't just fill in zeros for plane stress: StructuralMaterial :: giveFullSymVectorForm(answer, ep, aGaussPoint->giveMaterialMode()); return 1; } else if ( type == IST_MaxEquivalentStrainLevel ) { answer.resize(1); answer.at(1) = status->giveCumulativePlasticStrain(); return 1; } else if ( ( type == IST_DamageScalar ) || ( type == IST_DamageTensor ) ) { answer.resize(1); answer.at(1) = status->giveDamage(); return 1; } else { return StructuralMaterial :: giveIPValue(answer, aGaussPoint, type, atTime); } }
void MisesMat :: give3dLSMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) { MisesMatStatus *status = static_cast< MisesMatStatus * >( this->giveStatus(gp) ); // start from the elastic stiffness FloatMatrix I(6, 6); I.at(1, 1) = I.at(2, 2) = I.at(3, 3) = 1; I.at(4, 4) = I.at(5, 5) = I.at(6, 6) = 0.5; FloatArray delta(6); delta.at(1) = delta.at(2) = delta.at(3) = 1; FloatMatrix F; F.beMatrixForm( status->giveTempFVector() ); double J = F.giveDeterminant(); const FloatArray &trialStressDev = status->giveTrialStressDev(); double trialStressVol = status->giveTrialStressVol(); double trialS = computeStressNorm(trialStressDev); FloatArray n = trialStressDev; if ( trialS == 0 ) { n.resize(6); } else { n.times(1 / trialS); } FloatMatrix C; FloatMatrix help; help.beDyadicProductOf(delta, delta); C = help; help.times(-1. / 3.); FloatMatrix C1 = I; C1.add(help); C1.times(2 * trialStressVol); FloatMatrix n1, n2; n1.beDyadicProductOf(n, delta); n2.beDyadicProductOf(delta, n); help = n1; help.add(n2); C1.add(-2. / 3. * trialS, help); C.times(K * J * J); C.add(-K * ( J * J - 1 ), I); C.add(C1); ////////////////////////////////////////////////////////////////////////////////////////////////////////// FloatMatrix invF; FloatMatrix T(6, 6); invF.beInverseOf(F); ////////////////////////////////////////////////// //first row of pull back transformation matrix T.at(1, 1) = invF.at(1, 1) * invF.at(1, 1); T.at(1, 2) = invF.at(1, 2) * invF.at(1, 2); T.at(1, 3) = invF.at(1, 3) * invF.at(1, 3); T.at(1, 4) = 2. * invF.at(1, 2) * invF.at(1, 3); T.at(1, 5) = 2. * invF.at(1, 1) * invF.at(1, 3); T.at(1, 6) = 2. * invF.at(1, 1) * invF.at(1, 2); //second row of pull back transformation matrix T.at(2, 1) = invF.at(2, 1) * invF.at(2, 1); T.at(2, 2) = invF.at(2, 2) * invF.at(2, 2); T.at(2, 3) = invF.at(2, 3) * invF.at(2, 3); T.at(2, 4) = 2. * invF.at(2, 2) * invF.at(2, 3); T.at(2, 5) = 2. * invF.at(2, 1) * invF.at(2, 3); T.at(2, 6) = 2. * invF.at(2, 1) * invF.at(2, 2); //third row of pull back transformation matrix T.at(3, 1) = invF.at(3, 1) * invF.at(3, 1); T.at(3, 2) = invF.at(3, 2) * invF.at(3, 2); T.at(3, 3) = invF.at(3, 3) * invF.at(3, 3); T.at(3, 4) = 2. * invF.at(3, 2) * invF.at(3, 3); T.at(3, 5) = 2. * invF.at(3, 1) * invF.at(3, 3); T.at(3, 6) = 2. * invF.at(3, 1) * invF.at(3, 2); //fourth row of pull back transformation matrix T.at(4, 1) = invF.at(2, 1) * invF.at(3, 1); T.at(4, 2) = invF.at(2, 2) * invF.at(3, 2); T.at(4, 3) = invF.at(2, 3) * invF.at(3, 3); T.at(4, 4) = ( invF.at(2, 2) * invF.at(3, 3) + invF.at(2, 3) * invF.at(3, 2) ); T.at(4, 5) = ( invF.at(2, 1) * invF.at(3, 3) + invF.at(2, 3) * invF.at(3, 1) ); T.at(4, 6) = ( invF.at(2, 1) * invF.at(3, 2) + invF.at(2, 2) * invF.at(3, 1) ); //fifth row of pull back transformation matrix T.at(5, 1) = invF.at(1, 1) * invF.at(3, 1); T.at(5, 2) = invF.at(1, 2) * invF.at(3, 2); T.at(5, 3) = invF.at(1, 3) * invF.at(3, 3); T.at(5, 4) = ( invF.at(1, 2) * invF.at(3, 3) + invF.at(1, 3) * invF.at(3, 2) ); T.at(5, 5) = ( invF.at(1, 1) * invF.at(3, 3) + invF.at(1, 3) * invF.at(3, 1) ); T.at(5, 6) = ( invF.at(1, 1) * invF.at(3, 2) + invF.at(1, 2) * invF.at(3, 1) ); //sixth row of pull back transformation matrix T.at(6, 1) = invF.at(1, 1) * invF.at(2, 1); T.at(6, 2) = invF.at(1, 2) * invF.at(2, 2); T.at(6, 3) = invF.at(1, 3) * invF.at(2, 3); T.at(6, 4) = ( invF.at(1, 2) * invF.at(2, 3) + invF.at(1, 3) * invF.at(2, 2) ); T.at(6, 5) = ( invF.at(1, 1) * invF.at(2, 3) + invF.at(1, 3) * invF.at(2, 1) ); T.at(6, 6) = ( invF.at(1, 1) * invF.at(2, 2) + invF.at(1, 2) * invF.at(2, 1) ); /////////////////////////////////////////// if ( mode != TangentStiffness ) { help.beProductTOf(C, T); answer.beProductOf(T, help); return; } double kappa = status->giveCumulativePlasticStrain(); // increment of cumulative plastic strain as an indicator of plastic loading double dKappa = sqrt(3. / 2.) * ( status->giveTempCumulativePlasticStrain() - kappa ); //double dKappa = ( status->giveTempCumulativePlasticStrain() - kappa); if ( dKappa <= 0.0 ) { // elastic loading - elastic stiffness plays the role of tangent stiffness help.beProductTOf(C, T); answer.beProductOf(T, help); return; } // === plastic loading === //dKappa = dKappa*sqrt(3./2.); // trial deviatoric stress and its norm double beta0, beta1, beta2, beta3, beta4; if ( trialS == 0 ) { beta1 = 0; } else { beta1 = 2 * trialStressVol * dKappa / trialS; } if ( trialStressVol == 0 ) { beta0 = 0; beta2 = 0; beta3 = beta1; beta4 = 0; } else { beta0 = 1 + H / 3 / trialStressVol; beta2 = ( 1 - 1 / beta0 ) * 2. / 3. * trialS * dKappa / trialStressVol; beta3 = 1 / beta0 - beta1 + beta2; beta4 = ( 1 / beta0 - beta1 ) * trialS / trialStressVol; } FloatMatrix N; N.beDyadicProductOf(n, n); N.times(-2 * trialStressVol * beta3); C1.times(-beta1); FloatMatrix mN(3, 3); mN.at(1, 1) = n.at(1); mN.at(1, 2) = n.at(6); mN.at(1, 3) = n.at(5); mN.at(2, 1) = n.at(6); mN.at(2, 2) = n.at(2); mN.at(2, 3) = n.at(4); mN.at(3, 1) = n.at(5); mN.at(3, 2) = n.at(4); mN.at(3, 3) = n.at(3); FloatMatrix mN2; mN2.beProductOf(mN, mN); double volN2 = 1. / 3. * ( mN2.at(1, 1) + mN2.at(2, 2) + mN2.at(3, 3) ); FloatArray devN2(6); devN2.at(1) = mN2.at(1, 1) - volN2; devN2.at(2) = mN2.at(2, 2) - volN2; devN2.at(3) = mN2.at(3, 3) - volN2; devN2.at(4) = mN2.at(2, 3); devN2.at(5) = mN2.at(1, 3); devN2.at(6) = mN2.at(1, 2); FloatMatrix nonSymPart; nonSymPart.beDyadicProductOf(n, devN2); //symP.beTranspositionOf(nonSymPart); //symP.add(nonSymPart); //symP.times(1./2.); nonSymPart.times(-2 * trialStressVol * beta4); C.add(C1); C.add(N); C.add(nonSymPart); help.beProductTOf(C, T); answer.beProductOf(T, help); }
void MisesMat :: performPlasticityReturn(GaussPoint *gp, const FloatArray &totalStrain) { MisesMatStatus *status = static_cast< MisesMatStatus * >( this->giveStatus(gp) ); double kappa; FloatArray plStrain; FloatArray fullStress; // get the initial plastic strain and initial kappa from the status plStrain = status->givePlasticStrain(); kappa = status->giveCumulativePlasticStrain(); // === radial return algorithm === if ( totalStrain.giveSize() == 1 ) { LinearElasticMaterial *lmat = this->giveLinearElasticMaterial(); double E = lmat->give('E', gp); /*trial stress*/ fullStress.resize(6); fullStress.at(1) = E * ( totalStrain.at(1) - plStrain.at(1) ); double trialS = fabs(fullStress.at(1)); /*yield function*/ double yieldValue = trialS - (sig0 + H * kappa); // === radial return algorithm === if ( yieldValue > 0 ) { double dKappa = yieldValue / ( H + E ); kappa += dKappa; plStrain.at(1) += dKappa * signum( fullStress.at(1) ); plStrain.at(2) -= 0.5 * dKappa * signum( fullStress.at(1) ); plStrain.at(3) -= 0.5 * dKappa * signum( fullStress.at(1) ); fullStress.at(1) -= dKappa * E * signum( fullStress.at(1) ); } } else { // elastic predictor FloatArray elStrain = totalStrain; elStrain.subtract(plStrain); FloatArray elStrainDev; double elStrainVol; elStrainVol = computeDeviatoricVolumetricSplit(elStrainDev, elStrain); FloatArray trialStressDev; applyDeviatoricElasticStiffness(trialStressDev, elStrainDev, G); /**************************************************************/ double trialStressVol = 3 * K * elStrainVol; /**************************************************************/ // store the deviatoric and trial stress (reused by algorithmic stiffness) status->letTrialStressDevBe(trialStressDev); status->setTrialStressVol(trialStressVol); // check the yield condition at the trial state double trialS = computeStressNorm(trialStressDev); double yieldValue = sqrt(3./2.) * trialS - (sig0 + H * kappa); if ( yieldValue > 0. ) { // increment of cumulative plastic strain double dKappa = yieldValue / ( H + 3. * G ); kappa += dKappa; FloatArray dPlStrain; // the following line is equivalent to multiplication by scaling matrix P applyDeviatoricElasticCompliance(dPlStrain, trialStressDev, 0.5); // increment of plastic strain plStrain.add(sqrt(3. / 2.) * dKappa / trialS, dPlStrain); // scaling of deviatoric trial stress trialStressDev.times(1. - sqrt(6.) * G * dKappa / trialS); } // assemble the stress from the elastically computed volumetric part // and scaled deviatoric part computeDeviatoricVolumetricSum(fullStress, trialStressDev, trialStressVol); } // store the effective stress in status status->letTempEffectiveStressBe(fullStress); // store the plastic strain and cumulative plastic strain status->letTempPlasticStrainBe(plStrain); status->setTempCumulativePlasticStrain(kappa); }
void MisesMat :: giveFirstPKStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &totalDefGradOOFEM, TimeStep *tStep) { MisesMatStatus *status = static_cast< MisesMatStatus * >( this->giveStatus(gp) ); double kappa, dKappa, yieldValue, mi; FloatMatrix F, oldF, invOldF; FloatArray s; F.beMatrixForm(totalDefGradOOFEM); //(method assumes full 3D) kappa = status->giveCumulativePlasticStrain(); oldF.beMatrixForm( status->giveFVector() ); invOldF.beInverseOf(oldF); //relative deformation radient FloatMatrix f; f.beProductOf(F, invOldF); //compute elastic predictor FloatMatrix trialLeftCauchyGreen, help; f.times( 1./cbrt(f.giveDeterminant()) ); help.beProductOf(f, status->giveTempLeftCauchyGreen()); trialLeftCauchyGreen.beProductTOf(help, f); FloatMatrix E; E.beTProductOf(F, F); E.at(1, 1) -= 1.0; E.at(2, 2) -= 1.0; E.at(3, 3) -= 1.0; E.times(0.5); FloatArray e; e.beSymVectorFormOfStrain(E); FloatArray leftCauchyGreen; FloatArray leftCauchyGreenDev; double leftCauchyGreenVol; leftCauchyGreen.beSymVectorFormOfStrain(trialLeftCauchyGreen); leftCauchyGreenVol = computeDeviatoricVolumetricSplit(leftCauchyGreenDev, leftCauchyGreen); FloatArray trialStressDev; applyDeviatoricElasticStiffness(trialStressDev, leftCauchyGreenDev, G / 2.); s = trialStressDev; //check for plastic loading double trialS = computeStressNorm(trialStressDev); double sigmaY = sig0 + H * kappa; //yieldValue = sqrt(3./2.)*trialS-sigmaY; yieldValue = trialS - sqrt(2. / 3.) * sigmaY; //store deviatoric trial stress(reused by algorithmic stiffness) status->letTrialStressDevBe(trialStressDev); //the return-mapping algorithm double J = F.giveDeterminant(); mi = leftCauchyGreenVol * G; if ( yieldValue > 0 ) { //dKappa =sqrt(3./2.)* yieldValue/(H + 3.*mi); //kappa = kappa + dKappa; //trialStressDev.times(1-sqrt(6.)*mi*dKappa/trialS); dKappa = ( yieldValue / ( 2 * mi ) ) / ( 1 + H / ( 3 * mi ) ); FloatArray n = trialStressDev; n.times(2 * mi * dKappa / trialS); ////return map s.beDifferenceOf(trialStressDev, n); kappa += sqrt(2. / 3.) * dKappa; //update of intermediate configuration trialLeftCauchyGreen.beMatrixForm(s); trialLeftCauchyGreen.times(1.0 / G); trialLeftCauchyGreen.at(1, 1) += leftCauchyGreenVol; trialLeftCauchyGreen.at(2, 2) += leftCauchyGreenVol; trialLeftCauchyGreen.at(2, 2) += leftCauchyGreenVol; trialLeftCauchyGreen.times(J * J); } //addition of the elastic mean stress FloatMatrix kirchhoffStress; kirchhoffStress.beMatrixForm(s); kirchhoffStress.at(1, 1) += 1. / 2. * K * ( J * J - 1 ); kirchhoffStress.at(2, 2) += 1. / 2. * K * ( J * J - 1 ); kirchhoffStress.at(3, 3) += 1. / 2. * K * ( J * J - 1 ); FloatMatrix iF, Ep(3, 3), S; FloatArray vF, vS, ep; //transform Kirchhoff stress into Second Piola - Kirchhoff stress iF.beInverseOf(F); help.beProductOf(iF, kirchhoffStress); S.beProductTOf(help, iF); this->computeGLPlasticStrain(F, Ep, trialLeftCauchyGreen, J); ep.beSymVectorFormOfStrain(Ep); vS.beSymVectorForm(S); vF.beVectorForm(F); answer.beVectorForm(kirchhoffStress); status->setTrialStressVol(mi); status->letTempLeftCauchyGreenBe(trialLeftCauchyGreen); status->setTempCumulativePlasticStrain(kappa); status->letTempStressVectorBe(answer); status->letTempStrainVectorBe(e); status->letTempPlasticStrainBe(ep); status->letTempPVectorBe(answer); status->letTempFVectorBe(vF); }
// returns the stress vector in 3d stress space // computed from the previous plastic strain and current total strain void MisesMat :: performPlasticityReturn(GaussPoint *gp, const FloatArray &totalStrain) { double kappa, yieldValue, dKappa; FloatArray reducedStress; FloatArray strain, plStrain; MaterialMode mode = gp->giveMaterialMode(); MisesMatStatus *status = static_cast< MisesMatStatus * >( this->giveStatus(gp) ); StressVector fullStress(mode); this->initTempStatus(gp); this->initGpForNewStep(gp); // get the initial plastic strain and initial kappa from the status status->givePlasticStrain(plStrain); kappa = status->giveCumulativePlasticStrain(); // === radial return algorithm === if ( mode == _1dMat ) { LinearElasticMaterial *lmat = this->giveLinearElasticMaterial(); double E = lmat->give('E', gp); /*trial stress*/ fullStress.at(1) = E * ( totalStrain.at(1) - plStrain.at(1) ); double sigmaY = sig0 + H * kappa; double trialS = fullStress.at(1); trialS = fabs(trialS); /*yield function*/ yieldValue = trialS - sigmaY; // === radial return algorithm === if ( yieldValue > 0 ) { dKappa = yieldValue / ( H + E ); kappa += dKappa; fullStress.at(1) = fullStress.at(1) - dKappa *E *signum( fullStress.at(1) ); plStrain.at(1) = plStrain.at(1) + dKappa *signum( fullStress.at(1) ); } } else if ( ( mode == _PlaneStrain ) || ( mode == _3dMat ) ) { // elastic predictor StrainVector elStrain(totalStrain, mode); elStrain.subtract(plStrain); StrainVector elStrainDev(mode); double elStrainVol; elStrain.computeDeviatoricVolumetricSplit(elStrainDev, elStrainVol); StressVector trialStressDev(mode); elStrainDev.applyDeviatoricElasticStiffness(trialStressDev, G); /**************************************************************/ double trialStressVol; trialStressVol = 3 * K * elStrainVol; /**************************************************************/ // store the deviatoric and trial stress (reused by algorithmic stiffness) status->letTrialStressDevBe(trialStressDev); status->setTrialStressVol(trialStressVol); // check the yield condition at the trial state double trialS = trialStressDev.computeStressNorm(); double sigmaY = sig0 + H * kappa; yieldValue = sqrt(3. / 2.) * trialS - sigmaY; if ( yieldValue > 0. ) { // increment of cumulative plastic strain dKappa = yieldValue / ( H + 3. * G ); kappa += dKappa; StrainVector dPlStrain(mode); // the following line is equivalent to multiplication by scaling matrix P trialStressDev.applyDeviatoricElasticCompliance(dPlStrain, 0.5); // increment of plastic strain dPlStrain.times(sqrt(3. / 2.) * dKappa / trialS); plStrain.add(dPlStrain); // scaling of deviatoric trial stress trialStressDev.times(1. - sqrt(6.) * G * dKappa / trialS); } // assemble the stress from the elastically computed volumetric part // and scaled deviatoric part double stressVol = 3. * K * elStrainVol; trialStressDev.computeDeviatoricVolumetricSum(fullStress, stressVol); } // store the effective stress in status status->letTempEffectiveStressBe(fullStress); // store the plastic strain and cumulative plastic strain status->letTempPlasticStrainBe(plStrain); status->setTempCumulativePlasticStrain(kappa); }