void SimpleCrossSection :: giveGeneralizedStress_Beam2d(FloatArray &answer, GaussPoint *gp, const FloatArray &strain, TimeStep *tStep) { /**Note: (by bp): This assumes that the behaviour is elastic there exist a nuumber of nonlinear integral material models for beams defined directly in terms of integral forces and moments and corresponding deformations and curvatures. This would require to implement support at material model level. Mikael: That would not be a continuum material model, but it would highly depend on the cross-section shape, thus, it should be a special cross-section model instead. This cross-section assumes you can split the response into inertia moments and pure material response. This is only possible for a constant, elastic response (i.e. elastic). Therefore, this cross-section may only be allowed to give the elastic response. */ StructuralMaterial *mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) ); FloatArray elasticStrain, et, e0; FloatMatrix tangent; elasticStrain = strain; this->giveTemperatureVector(et, gp, tStep); if ( et.giveSize() > 0 ) { mat->giveThermalDilatationVector(e0, gp, tStep); double thick = this->give(CS_Thickness, gp); elasticStrain.at(1) -= e0.at(1) * ( et.at(1) - mat->giveReferenceTemperature() ); if ( et.giveSize() > 1 ) { elasticStrain.at(2) -= e0.at(1) * et.at(2) / thick; // kappa_x } } this->give2dBeamStiffMtrx(tangent, ElasticStiffness, gp, tStep); answer.beProductOf(tangent, elasticStrain); StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) ); status->letTempStrainVectorBe(strain); status->letTempStressVectorBe(answer); }
void SimpleVitrificationMaterial :: giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep) { FloatArray strainVector; FloatMatrix d; FloatArray deltaStrain; StructuralMaterialStatus *status = dynamic_cast< StructuralMaterialStatus * >( this->giveStatus(gp) ); this->giveStressDependentPartOfStrainVector(strainVector, gp, reducedStrain, tStep, VM_Total); deltaStrain.beDifferenceOf( strainVector, status->giveStrainVector() ); this->give3dMaterialStiffnessMatrix(d, TangentStiffness, gp, tStep); FloatArray deltaStress; deltaStress.beProductOf(d, deltaStrain); answer = status->giveStressVector(); answer.add(deltaStress); // update gp status->letTempStrainVectorBe(reducedStrain); status->letTempStressVectorBe(answer); }
void PhaseFieldElement :: computeStiffnessMatrix_ud(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep) { FloatMatrix B, DB, N_d, DN, S(3,1); FloatArray stress; NLStructuralElement *el = this->giveElement(); StructuralCrossSection *cs = dynamic_cast<StructuralCrossSection* > (el->giveCrossSection() ); answer.clear(); for ( auto &gp: el->giveIntegrationRule(0) ) { StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() ); double dV = el->computeVolumeAround(gp); // compute int_V ( B^t * D_B * B )dV this->computeNd_matrixAt(*gp->giveNaturalCoordinates(), N_d); // stress FloatArray reducedStrain, a_u; FloatMatrix B_u; el->computeBmatrixAt( gp, B_u ); stress = matStat->giveTempStressVector(); stress.times( this->computeGPrim( gp, VM_Total, tStep ) ); S.setColumn(stress,1); DN.beProductOf(S, N_d); answer.plusProductUnsym(B_u, DN, dV); } }
void HyperElasticMaterial :: giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &totalStrain, TimeStep *tStep) { double J2; FloatMatrix C(3, 3); FloatMatrix invC; FloatArray strainVector; StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( this->giveStatus(gp) ); this->giveStressDependentPartOfStrainVector_3d(strainVector, gp, totalStrain, tStep, VM_Total); C.at(1, 1) = 1. + 2. * strainVector.at(1); C.at(2, 2) = 1. + 2. * strainVector.at(2); C.at(3, 3) = 1. + 2. * strainVector.at(3); C.at(1, 2) = C.at(2, 1) = strainVector.at(6); C.at(1, 3) = C.at(3, 1) = strainVector.at(5); C.at(2, 3) = C.at(3, 2) = strainVector.at(4); invC.beInverseOf(C); J2 = C.giveDeterminant(); answer.resize(6); double aux = ( K - 2. / 3. * G ) * ( J2 - 1. ) / 2. - G; answer.at(1) = aux * invC.at(1, 1) + G; answer.at(2) = aux * invC.at(2, 2) + G; answer.at(3) = aux * invC.at(3, 3) + G; answer.at(4) = aux * invC.at(2, 3); answer.at(5) = aux * invC.at(1, 3); answer.at(6) = aux * invC.at(1, 2); // update gp status->letTempStrainVectorBe(totalStrain); status->letTempStressVectorBe(answer); }
void SimpleCrossSection :: giveGeneralizedStress_Shell(FloatArray &answer, GaussPoint *gp, const FloatArray &strain, TimeStep *tStep) { /**Note: (by bp): This assumes that the behaviour is elastic there exist a nuumber of nonlinear integral material models for beams/plates/shells defined directly in terms of integral forces and moments and corresponding deformations and curvatures. This would require to implement support at material model level. Mikael: See earlier response to comment */ StructuralMaterial *mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) ); FloatArray elasticStrain, et, e0; FloatMatrix tangent; elasticStrain = strain; this->giveTemperatureVector(et, gp, tStep); if ( et.giveSize() ) { double thick = this->give(CS_Thickness, gp); mat->giveThermalDilatationVector(e0, gp, tStep); elasticStrain.at(1) -= e0.at(1) * ( et.at(1) - mat->giveReferenceTemperature() ); elasticStrain.at(2) -= e0.at(2) * ( et.at(1) - mat->giveReferenceTemperature() ); if ( et.giveSize() > 1 ) { elasticStrain.at(4) -= e0.at(1) * et.at(2) / thick; // kappa_x elasticStrain.at(5) -= e0.at(2) * et.at(2) / thick; // kappa_y } } this->give3dShellStiffMtrx(tangent, ElasticStiffness, gp, tStep); answer.beProductOf(tangent, elasticStrain); StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) ); status->letTempStrainVectorBe(strain); status->letTempStressVectorBe(answer); }
double PhaseFieldElement :: computeFreeEnergy(GaussPoint *gp, TimeStep *tStep) { StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() ); FloatArray strain, stress; stress = matStat->giveTempStressVector(); strain = matStat->giveTempStrainVector(); return 0.5 * stress.dotProduct( strain ); }
void CCTPlate3d :: giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep) // returns characteristic tensor of the receiver at given gp and tStep // strain vector = (Kappa_x, Kappa_y, Kappa_xy, Gamma_zx, Gamma_zy) { FloatArray charVect; Material *mat = this->giveMaterial(); StructuralMaterialStatus *ms = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) ); answer.resize(3, 3); answer.zero(); if ( ( type == LocalForceTensor ) || ( type == GlobalForceTensor ) ) { //this->computeStressVector(charVect, gp, tStep); charVect = ms->giveStressVector(); answer.at(1, 3) = charVect.at(4); answer.at(3, 1) = charVect.at(4); answer.at(2, 3) = charVect.at(5); answer.at(3, 2) = charVect.at(5); } else if ( ( type == LocalMomentumTensor ) || ( type == GlobalMomentumTensor ) ) { //this->computeStressVector(charVect, gp, tStep); charVect = ms->giveStressVector(); answer.at(1, 1) = charVect.at(1); answer.at(2, 2) = charVect.at(2); answer.at(1, 2) = charVect.at(3); answer.at(2, 1) = charVect.at(3); } else if ( ( type == LocalStrainTensor ) || ( type == GlobalStrainTensor ) ) { //this->computeStrainVector(charVect, gp, tStep); charVect = ms->giveStrainVector(); answer.at(1, 3) = charVect.at(4) / 2.; answer.at(3, 1) = charVect.at(4) / 2.; answer.at(2, 3) = charVect.at(5) / 2.; answer.at(3, 2) = charVect.at(5) / 2.; } else if ( ( type == LocalCurvatureTensor ) || ( type == GlobalCurvatureTensor ) ) { //this->computeStrainVector(charVect, gp, tStep); charVect = ms->giveStrainVector(); answer.at(1, 1) = charVect.at(1); answer.at(2, 2) = charVect.at(2); answer.at(1, 2) = charVect.at(3) / 2.; answer.at(2, 1) = charVect.at(3) / 2.; } else { _error("GiveCharacteristicTensor: unsupported tensor mode"); exit(1); } if ( ( type == GlobalForceTensor ) || ( type == GlobalMomentumTensor ) || ( type == GlobalStrainTensor ) || ( type == GlobalCurvatureTensor ) ) { this->computeGtoLRotationMatrix(); answer.rotatedWith(* GtoLRotationMatrix); } }
void WinklerPasternakMaterial::giveRealStressVector_2dPlateSubSoil(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) { FloatMatrix tangent; this->give2dPlateSubSoilStiffMtrx(tangent, ElasticStiffness, gp, tStep); answer.beProductOf(tangent, reducedE); StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( this->giveStatus(gp) ); status->letTempStrainVectorBe(reducedE); status->letTempStressVectorBe(answer); }
void HyperElasticMaterial :: give3dMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) { double J2, c11, c22, c33, c12, c13, c23, A, B; FloatMatrix C(3, 3); FloatMatrix invC; StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( this->giveStatus(gp) ); C.at(1, 1) = 1. + 2. * status->giveTempStrainVector().at(1); C.at(2, 2) = 1. + 2. * status->giveTempStrainVector().at(2); C.at(3, 3) = 1. + 2. * status->giveTempStrainVector().at(3); C.at(2, 3) = C.at(3, 2) = status->giveTempStrainVector().at(4); C.at(1, 3) = C.at(3, 1) = status->giveTempStrainVector().at(5); C.at(1, 2) = C.at(2, 1) = status->giveTempStrainVector().at(6); invC.beInverseOf(C); J2 = C.giveDeterminant(); c11 = invC.at(1, 1); c22 = invC.at(2, 2); c33 = invC.at(3, 3); c12 = invC.at(1, 2); c13 = invC.at(1, 3); c23 = invC.at(2, 3); A = ( K - 2. / 3. * G ) * J2; B = -( K - 2. / 3. * G ) * ( J2 - 1. ) + 2. * G; answer.resize(6, 6); answer.at(1, 1) = ( A + B ) * c11 * c11; answer.at(2, 2) = ( A + B ) * c22 * c22; answer.at(3, 3) = ( A + B ) * c33 * c33; answer.at(4, 4) = A * c23 * c23 + B / 2. * ( c22 * c33 + c23 * c23 ); answer.at(5, 5) = A * c13 * c13 + B / 2. * ( c11 * c33 + c13 * c13 ); answer.at(6, 6) = A * c12 * c12 + B / 2. * ( c11 * c22 + c12 * c12 ); answer.at(1, 2) = answer.at(2, 1) = A * c11 * c22 + B * c12 * c12; answer.at(1, 3) = answer.at(3, 1) = A * c11 * c33 + B * c13 * c13; answer.at(1, 4) = answer.at(4, 1) = A * c11 * c23 + B * c12 * c13; answer.at(1, 5) = answer.at(5, 1) = A * c11 * c13 + B * c11 * c13; answer.at(1, 6) = answer.at(6, 1) = A * c11 * c12 + B * c11 * c12; answer.at(2, 3) = answer.at(3, 2) = A * c22 * c33 + B * c23 * c23; answer.at(2, 4) = answer.at(4, 2) = A * c22 * c23 + B * c22 * c23; answer.at(2, 5) = answer.at(5, 2) = A * c22 * c13 + B * c12 * c23; answer.at(2, 6) = answer.at(6, 2) = A * c22 * c12 + B * c22 * c12; answer.at(3, 4) = answer.at(4, 3) = A * c33 * c23 + B * c33 * c23; answer.at(3, 5) = answer.at(5, 3) = A * c33 * c13 + B * c33 * c13; answer.at(3, 6) = answer.at(6, 3) = A * c33 * c12 + B * c13 * c23; answer.at(4, 5) = answer.at(5, 4) = A * c23 * c13 + B / 2. * ( c12 * c33 + c13 * c23 ); answer.at(4, 6) = answer.at(6, 4) = A * c23 * c12 + B / 2. * ( c12 * c23 + c22 * c13 ); answer.at(5, 6) = answer.at(6, 5) = A * c13 * c12 + B / 2. * ( c11 * c23 + c12 * c13 ); }
void FiberedCrossSection :: giveGeneralizedStress_Beam3d(FloatArray &answer, GaussPoint *gp, const FloatArray &strain, TimeStep *tStep) { double fiberThick, fiberWidth, fiberZCoord, fiberYCoord; FloatArray fiberStrain, reducedFiberStress; StructuralElement *element = static_cast< StructuralElement * >( gp->giveElement() ); FiberedCrossSectionInterface *interface; if ( ( interface = static_cast< FiberedCrossSectionInterface * >( element->giveInterface(FiberedCrossSectionInterfaceType) ) ) == NULL ) { OOFEM_ERROR("element with no fiber support encountered"); } answer.resize(6); answer.zero(); for ( int i = 1; i <= numberOfFibers; i++ ) { GaussPoint *fiberGp = this->giveSlaveGaussPoint(gp, i - 1); StructuralMaterial *fiberMat = static_cast< StructuralMaterial * >( domain->giveMaterial( fiberMaterials.at(i) ) ); // the question is whether this function should exist ? // if yes the element details will be hidden. // good idea also should be existence of element::GiveBmatrixOfLayer // and computing strains here - but first idea looks better // but treating of geometric non-linearities may become more complicated // another approach - use several functions with assumed kinematic constraints // resolve current layer z-coordinate fiberThick = this->fiberThicks.at(i); fiberWidth = this->fiberWidths.at(i); fiberYCoord = fiberGp->giveNaturalCoordinate(1); fiberZCoord = fiberGp->giveNaturalCoordinate(2); interface->FiberedCrossSectionInterface_computeStrainVectorInFiber(fiberStrain, strain, fiberGp, tStep); fiberMat->giveRealStressVector_Fiber(reducedFiberStress, fiberGp, fiberStrain, tStep); // perform integration // 1) membrane terms N, Qz, Qy answer.at(1) += reducedFiberStress.at(1) * fiberWidth * fiberThick; answer.at(2) += reducedFiberStress.at(2) * fiberWidth * fiberThick; answer.at(3) += reducedFiberStress.at(3) * fiberWidth * fiberThick; // 2) bending terms mx, my, mxy answer.at(4) += ( reducedFiberStress.at(2) * fiberWidth * fiberThick * fiberYCoord - reducedFiberStress.at(3) * fiberWidth * fiberThick * fiberZCoord ); answer.at(5) += reducedFiberStress.at(1) * fiberWidth * fiberThick * fiberZCoord; answer.at(6) -= reducedFiberStress.at(1) * fiberWidth * fiberThick * fiberYCoord; } // now we must update master gp ///@ todo simply chosen the first fiber material as master material /JB StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * > ( domain->giveMaterial( fiberMaterials.at(1) )->giveStatus(gp) ); status->letTempStrainVectorBe(strain); status->letTempStressVectorBe(answer); }
void StructuralMaterialSettable :: giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &totalStrain, TimeStep *atTime) { StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( this->giveStatus(gp) ); const FloatArray& stressVector = status->giveStressVector(); status->letTempStrainVectorBe(totalStrain); status->letTempStressVectorBe(stressVector); answer = stressVector; }
int FiberedCrossSection :: giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) { Material *mat = this->giveDomain()->giveMaterial( fiberMaterials.at(1) ); ///@todo For now, create material status according to the first fiber material StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) ); if ( type == IST_BeamForceMomentTensor ) { answer = status->giveStressVector(); return 1; } else if ( type == IST_BeamStrainCurvatureTensor ) { answer = status->giveStrainVector(); return 1; } return CrossSection :: giveIPValue(answer, gp, type, tStep); }
void SimpleCrossSection :: giveGeneralizedStress_MembraneRot(FloatArray &answer, GaussPoint *gp, const FloatArray &strain, TimeStep *tStep) { FloatMatrix tangent; this->giveMembraneRotStiffMtrx(tangent, ElasticStiffness, gp, tStep); answer.beProductOf(tangent, strain); StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( this->giveMaterial(gp)->giveStatus(gp) ); status->letTempStrainVectorBe(strain); status->letTempStressVectorBe(answer); ///@todo We should support nonlinear behavior for the membrane part. In fact, should be even bundle the rotation part with the membrane? /// We gain nothing from this design anyway as the rotation field is always separate. Separate manual integration by the element would be an option. }
void LIBeam3dNL :: computeTempCurv(FloatArray &answer, TimeStep *tStep) { IntegrationRule *iRule = this->giveDefaultIntegrationRulePtr(); GaussPoint *gp = iRule->getIntegrationPoint(0); FloatArray ui(3), ac(3), PrevEpsilon; FloatMatrix sc(3, 3), tmid(3, 3); // update curvature at midpoint // first, compute Tmid // ask increments this->computeVectorOf(VM_Incremental, tStep, ui); ac.at(1) = 0.5 * ( ui.at(10) - ui.at(4) ); ac.at(2) = 0.5 * ( ui.at(11) - ui.at(5) ); ac.at(3) = 0.5 * ( ui.at(12) - ui.at(6) ); this->computeSMtrx(sc, ac); sc.times(1. / 2.); // compute I+sc sc.at(1, 1) += 1.0; sc.at(2, 2) += 1.0; sc.at(3, 3) += 1.0; tmid.beProductOf(sc, this->tc); // update curvature at centre ac.at(1) = ( ui.at(10) - ui.at(4) ); ac.at(2) = ( ui.at(11) - ui.at(5) ); ac.at(3) = ( ui.at(12) - ui.at(6) ); answer.beTProductOf(tmid, ac); answer.times(1 / this->l0); // ask for previous kappa StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() ); if ( matStat ) { PrevEpsilon = matStat->giveStrainVector(); } if ( PrevEpsilon.giveSize() ) { answer.at(1) += PrevEpsilon.at(4); answer.at(2) += PrevEpsilon.at(5); answer.at(3) += PrevEpsilon.at(6); } }
void LinQuad3DPlaneStress :: giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep) // returns characteristic tensor of the receiver at given gp and tStep // strain vector = (Eps_X, Eps_y, Gamma_xy, Kappa_z) { FloatArray charVect; StructuralMaterialStatus *ms = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() ); answer.resize(3, 3); answer.zero(); if ( ( type == LocalForceTensor ) || ( type == GlobalForceTensor ) ) { //this->computeStressVector(charVect, gp, tStep); charVect = ms->giveStressVector(); answer.at(1, 1) = charVect.at(1); answer.at(2, 2) = charVect.at(2); answer.at(1, 2) = charVect.at(3); answer.at(2, 1) = charVect.at(3); } else if ( ( type == LocalMomentTensor ) || ( type == GlobalMomentTensor ) ) { } else if ( ( type == LocalStrainTensor ) || ( type == GlobalStrainTensor ) ) { //this->computeStrainVector(charVect, gp, tStep); charVect = ms->giveStrainVector(); answer.at(1, 1) = charVect.at(1); answer.at(2, 2) = charVect.at(2); answer.at(1, 2) = charVect.at(3) / 2.; answer.at(2, 1) = charVect.at(3) / 2.; } else if ( ( type == LocalCurvatureTensor ) || ( type == GlobalCurvatureTensor ) ) { } else { OOFEM_ERROR("unsupported tensor mode"); exit(1); } if ( ( type == GlobalForceTensor ) || ( type == GlobalMomentTensor ) || ( type == GlobalStrainTensor ) || ( type == GlobalCurvatureTensor ) ) { this->computeGtoLRotationMatrix(); answer.rotatedWith(* GtoLRotationMatrix); } }
void StructuralMaterialEvaluator :: solveYourself() { Domain *d = this->giveDomain(1); MaterialMode mode = _3dMat; FloatArray initialStrain(6); gps.clear(); gps.reserve(d->giveNumberOfMaterialModels()); for ( int i = 1; i <= d->giveNumberOfMaterialModels(); i++ ) { std :: unique_ptr< GaussPoint > gp = std::make_unique<GaussPoint>(nullptr, i, FloatArray(0), 1, mode); gps.emplace_back( std :: move(gp) ); // Initialize the strain vector; StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( d->giveMaterial(i)->giveStatus( gps[i-1].get() ) ); status->letStrainVectorBe(initialStrain); } std :: string outname = this->giveOutputBaseFileName() + ".matdata"; this->outfile.open( outname.c_str() ); this->timer.startTimer(EngngModelTimer :: EMTT_AnalysisTimer); TimeStep *tStep = giveNextStep(); // Note, strain == strain-rate (kept as strain for brevity) int maxiter = 100; // User input? FloatArray stressC, deltaStrain, strain, stress, res; stressC.resize( sControl.giveSize() ); res.resize( sControl.giveSize() ); FloatMatrix tangent, reducedTangent; for ( int istep = 1; istep <= this->numberOfSteps; ++istep ) { this->timer.startTimer(EngngModelTimer :: EMTT_SolutionStepTimer); for ( int imat = 1; imat <= d->giveNumberOfMaterialModels(); ++imat ) { GaussPoint *gp = gps[imat-1].get(); StructuralMaterial *mat = static_cast< StructuralMaterial * >( d->giveMaterial(imat) ); StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) ); strain = status->giveStrainVector(); // Update the controlled parts for ( int j = 1; j <= eControl.giveSize(); ++j ) { int p = eControl.at(j); strain.at(p) = d->giveFunction( cmpntFunctions.at(p) )->evaluateAtTime( tStep->giveIntrinsicTime() ); } for ( int j = 1; j <= sControl.giveSize(); ++j ) { int p = sControl.at(j); stressC.at(j) = d->giveFunction( cmpntFunctions.at(p) )->evaluateAtTime( tStep->giveIntrinsicTime() ); } //strain.add(-100, {6.27e-06, 6.27e-06, 6.27e-06, 0, 0, 0}); for ( int iter = 1; iter < maxiter; iter++ ) { #if 0 // Debugging: mat->give3dMaterialStiffnessMatrix(tangent, TangentStiffness, gp, tStep); tangent.printYourself("# tangent"); strain.zero(); mat->giveRealStressVector_3d(stress, gp, strain, tStep); FloatArray strain2; tangent.solveForRhs(stress, strain2); strain2.printYourself("# thermal expansion"); break; #endif strain.printYourself("Macro strain guess"); mat->giveRealStressVector_3d(stress, gp, strain, tStep); for ( int j = 1; j <= sControl.giveSize(); ++j ) { res.at(j) = stressC.at(j) - stress.at( sControl.at(j) ); } OOFEM_LOG_INFO("*** Time step: %d (t = %.2e), Material %d, Iteration: %d, Residual = %e (tolerance %.2e)\n", istep, tStep->giveIntrinsicTime(), imat, iter, res.computeNorm(), tolerance); if ( res.computeNorm() <= tolerance ) { break; } else { if ( tangent.giveNumberOfRows() == 0 || !keepTangent ) { mat->give3dMaterialStiffnessMatrix(tangent, TangentStiffness, gp, tStep); } // Pick out the stress-controlled part; reducedTangent.beSubMatrixOf(tangent, sControl, sControl); // Update stress-controlled part of the strain reducedTangent.solveForRhs(res, deltaStrain); //deltaStrain.printYourself("deltaStrain"); for ( int j = 1; j <= sControl.giveSize(); ++j ) { strain.at( sControl.at(j) ) += deltaStrain.at(j); } } } if ( res.computeNorm() > tolerance ) { OOFEM_WARNING("Residual did not converge!"); } // This material model has converged, so we update it and go on to the next. gp->updateYourself(tStep); } this->timer.stopTimer(EngngModelTimer :: EMTT_SolutionStepTimer); this->doStepOutput(tStep); tStep = giveNextStep(); } this->timer.stopTimer(EngngModelTimer :: EMTT_AnalysisTimer); this->outfile.close(); }
std::vector<std::unique_ptr<EnrichmentItem>> NCPrincipalStress::nucleateEnrichmentItems() { SpatialLocalizer *octree = this->mpDomain->giveSpatialLocalizer(); XfemManager *xMan = mpDomain->giveXfemManager(); std::vector<std::unique_ptr<EnrichmentItem>> eiList; // Center coordinates of newly inserted cracks std::vector<FloatArray> center_coord_inserted_cracks; // Loop over all elements and all bulk GP. for(auto &el : mpDomain->giveElements() ) { int numIR = el->giveNumberOfIntegrationRules(); int csNum = el->giveCrossSection()->giveNumber(); if(csNum == mCrossSectionInd || true) { for(int irInd = 0; irInd < numIR; irInd++) { IntegrationRule *ir = el->giveIntegrationRule(irInd); int numGP = ir->giveNumberOfIntegrationPoints(); for(int gpInd = 0; gpInd < numGP; gpInd++) { GaussPoint *gp = ir->getIntegrationPoint(gpInd); // int csNum = gp->giveCrossSection()->giveNumber(); // printf("csNum: %d\n", csNum); StructuralMaterialStatus *ms = dynamic_cast<StructuralMaterialStatus*>(gp->giveMaterialStatus()); if(ms != NULL) { const FloatArray &stress = ms->giveTempStressVector(); FloatArray principalVals; FloatMatrix principalDirs; StructuralMaterial::computePrincipalValDir(principalVals, principalDirs, stress, principal_stress); if(principalVals[0] > mStressThreshold) { // printf("\nFound GP with stress above threshold.\n"); // printf("principalVals: "); principalVals.printYourself(); FloatArray crackNormal; crackNormal.beColumnOf(principalDirs, 1); // printf("crackNormal: "); crackNormal.printYourself(); FloatArray crackTangent = {-crackNormal(1), crackNormal(0)}; crackTangent.normalize(); // printf("crackTangent: "); crackTangent.printYourself(); // Create geometry FloatArray pc = {gp->giveGlobalCoordinates()(0), gp->giveGlobalCoordinates()(1)}; // printf("Global coord: "); pc.printYourself(); FloatArray ps = pc; ps.add(-0.5*mInitialCrackLength, crackTangent); FloatArray pe = pc; pe.add(0.5*mInitialCrackLength, crackTangent); if(mCutOneEl) { // If desired, ensure that the crack cuts exactly one element. Line line(ps, pe); std::vector<FloatArray> intersecPoints; // line.computeIntersectionPoints(el.get(), intersecPoints); for ( int i = 1; i <= el->giveNumberOfDofManagers(); i++ ) { // int n1 = i; // int n2 = 0; // if ( i < el->giveNumberOfDofManagers() ) { // n2 = i + 1; // } else { // n2 = 1; // } // const FloatArray &p1 = *(el->giveDofManager(n1)->giveCoordinates()); // const FloatArray &p2 = *(el->giveDofManager(n2)->giveCoordinates()); } // printf("intersecPoints.size(): %lu\n", intersecPoints.size()); if(intersecPoints.size() == 2) { ps = std::move(intersecPoints[0]); pe = std::move(intersecPoints[1]); } else { OOFEM_ERROR("intersecPoints.size() != 2") } } FloatArray points = {ps(0), ps(1), pc(0), pc(1), pe(0), pe(1)}; // double diffX = 0.5*(ps(0) + pe(0)) - pc(0); // printf("diffX: %e\n", diffX); // double diffY = 0.5*(ps(1) + pe(1)) - pc(1); // printf("diffY: %e\n", diffY); // TODO: Check if nucleation is allowed, by checking for already existing cracks close to the GP. // Idea: Nucleation is not allowed if we are within an enriched element. In this way, branching is not // completely prohibited, but we avoid initiating multiple similar cracks. bool insertionAllowed = true; Element *el_s = octree->giveElementContainingPoint(ps); if(el_s) { if( xMan->isElementEnriched(el_s) ) { insertionAllowed = false; } } Element *el_c = octree->giveElementContainingPoint(pc); if(el_c) { if( xMan->isElementEnriched(el_c) ) { insertionAllowed = false; } } Element *el_e = octree->giveElementContainingPoint(pe); if(el_e) { if( xMan->isElementEnriched(el_e) ) { insertionAllowed = false; } } for(const auto &x: center_coord_inserted_cracks) { if( x.distance(pc) < 2.0*mInitialCrackLength) { insertionAllowed = false; break; printf("Preventing insertion.\n"); } } if(insertionAllowed) { int n = xMan->giveNumberOfEnrichmentItems() + 1; std::unique_ptr<Crack> crack = std::make_unique<Crack>(n, xMan, mpDomain); // Geometry std::unique_ptr<BasicGeometry> geom = std::make_unique<PolygonLine>(); geom->insertVertexBack(ps); geom->insertVertexBack(pc); geom->insertVertexBack(pe); crack->setGeometry(std::move(geom)); // Enrichment function EnrichmentFunction *ef = new HeavisideFunction(1, mpDomain); crack->setEnrichmentFunction(ef); // Enrichment fronts // EnrichmentFront *efStart = new EnrFrontLinearBranchFuncOneEl(); EnrichmentFront *efStart = new EnrFrontCohesiveBranchFuncOneEl(); crack->setEnrichmentFrontStart(efStart); // EnrichmentFront *efEnd = new EnrFrontLinearBranchFuncOneEl(); EnrichmentFront *efEnd = new EnrFrontCohesiveBranchFuncOneEl(); crack->setEnrichmentFrontEnd(efEnd); /////////////////////////////////////// // Propagation law // Options // double radius = 0.5*mInitialCrackLength, angleInc = 10.0, incrementLength = 0.5*mInitialCrackLength, hoopStressThreshold = 0.0; // bool useRadialBasisFunc = true; // PLHoopStressCirc *pl = new PLHoopStressCirc(); // pl->setRadius(radius); // pl->setAngleInc(angleInc); // pl->setIncrementLength(incrementLength); // pl->setHoopStressThreshold(hoopStressThreshold); // pl->setUseRadialBasisFunc(useRadialBasisFunc); // PLDoNothing *pl = new PLDoNothing(); PLMaterialForce *pl = new PLMaterialForce(); pl->setRadius(mMatForceRadius); pl->setIncrementLength(mIncrementLength); // pl->setIncrementLength(0.25); // pl->setCrackPropThreshold(0.25); pl->setCrackPropThreshold(mCrackPropThreshold); crack->setPropagationLaw(pl); crack->updateDofIdPool(); center_coord_inserted_cracks.push_back(pc); eiList.push_back( std::unique_ptr<EnrichmentItem>(std::move(crack)) ); // printf("Nucleating a crack in NCPrincipalStress::nucleateEnrichmentItems.\n"); // printf("el->giveGlobalNumber(): %d\n", el->giveGlobalNumber() ); // We only introduce one crack per element in a single time step. break; } } } } } } // If correct csNum }
void NLStructuralElement :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord) { FloatMatrix B; FloatArray vStress, vStrain, u; // This function can be quite costly to do inside the loops when one has many slave dofs. this->computeVectorOf(VM_Total, tStep, u); // subtract initial displacements, if defined if ( initialDisplacements ) { u.subtract(* initialDisplacements); } // zero answer will resize accordingly when adding first contribution answer.clear(); for ( auto &gp: *this->giveDefaultIntegrationRulePtr() ) { StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() ); // Engineering (small strain) stress if ( nlGeometry == 0 ) { this->computeBmatrixAt(gp, B); if ( useUpdatedGpRecord == 1 ) { vStress = matStat->giveStressVector(); } else { ///@todo Is this really what we should do for inactive elements? if ( !this->isActivated(tStep) ) { vStrain.resize( StructuralMaterial :: giveSizeOfVoigtSymVector( gp->giveMaterialMode() ) ); vStrain.zero(); } vStrain.beProductOf(B, u); this->computeStressVector(vStress, vStrain, gp, tStep); } } else if ( nlGeometry == 1 ) { // First Piola-Kirchhoff stress if ( this->domain->giveEngngModel()->giveFormulation() == AL ) { // Cauchy stress if ( useUpdatedGpRecord == 1 ) { vStress = matStat->giveCVector(); } else { this->computeCauchyStressVector(vStress, gp, tStep); } this->computeBmatrixAt(gp, B); } else { // First Piola-Kirchhoff stress if ( useUpdatedGpRecord == 1 ) { vStress = matStat->givePVector(); } else { this->computeFirstPKStressVector(vStress, gp, tStep); ///@todo This is actaully inefficient since it constructs B and twice and collects the nodal unknowns over and over. } this->computeBHmatrixAt(gp, B); } } if ( vStress.giveSize() == 0 ) { /// @todo is this check really necessary? break; } // Compute nodal internal forces at nodes as f = B^T*Stress dV double dV = this->computeVolumeAround(gp); if ( nlGeometry == 1 ) { // First Piola-Kirchhoff stress if ( vStress.giveSize() == 9 ) { FloatArray stressTemp; StructuralMaterial :: giveReducedVectorForm( stressTemp, vStress, gp->giveMaterialMode() ); answer.plusProduct(B, stressTemp, dV); } else { answer.plusProduct(B, vStress, dV); } } else { if ( vStress.giveSize() == 6 ) { // It may happen that e.g. plane strain is computed // using the default 3D implementation. If so, // the stress needs to be reduced. // (Note that no reduction will take place if // the simulation is actually 3D.) FloatArray stressTemp; StructuralMaterial :: giveReducedSymVectorForm( stressTemp, vStress, gp->giveMaterialMode() ); answer.plusProduct(B, stressTemp, dV); } else { answer.plusProduct(B, vStress, dV); } } } // If inactive: update fields but do not give any contribution to the internal forces if ( !this->isActivated(tStep) ) { answer.zero(); return; } }
void NLStructuralElement :: giveInternalForcesVector_withIRulesAsSubcells(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord) { /* * Returns nodal representation of real internal forces computed from first Piola-Kirchoff stress * if useGpRecord == 1 then stresses stored in the gp are used, otherwise stresses are computed * this must be done if you want internal forces after element->updateYourself() has been called * for the same time step. * The integration procedure uses an integrationRulesArray for numerical integration. * Each integration rule is considered to represent a separate sub-cell/element. Typically this would be used when * integration of the element domain needs special treatment, e.g. when using the XFEM. */ FloatMatrix B; FloatArray vStress, vStrain; IntArray irlocnum; FloatArray *m = & answer, temp; if ( this->giveInterpolation() && this->giveInterpolation()->hasSubPatchFormulation() ) { m = & temp; } // zero answer will resize accordingly when adding first contribution answer.clear(); // loop over individual integration rules for ( auto &iRule: integrationRulesArray ) { for ( GaussPoint *gp: *iRule ) { StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() ); if ( nlGeometry == 0 ) { this->computeBmatrixAt(gp, B); if ( useUpdatedGpRecord == 1 ) { vStress = matStat->giveStressVector(); } else { this->computeStrainVector(vStrain, gp, tStep); this->computeStressVector(vStress, vStrain, gp, tStep); } } else if ( nlGeometry == 1 ) { if ( this->domain->giveEngngModel()->giveFormulation() == AL ) { // Cauchy stress if ( useUpdatedGpRecord == 1 ) { vStress = matStat->giveCVector(); } else { this->computeCauchyStressVector(vStress, gp, tStep); } this->computeBmatrixAt(gp, B); } else { // First Piola-Kirchhoff stress if ( useUpdatedGpRecord == 1 ) { vStress = matStat->givePVector(); } else { this->computeFirstPKStressVector(vStress, gp, tStep); } this->computeBHmatrixAt(gp, B); } } if ( vStress.giveSize() == 0 ) { //@todo is this really necessary? break; } // compute nodal representation of internal forces at nodes as f = B^T*stress dV double dV = this->computeVolumeAround(gp); m->plusProduct(B, vStress, dV); // localize irule contribution into element matrix if ( this->giveIntegrationRuleLocalCodeNumbers(irlocnum, *iRule) ) { answer.assemble(* m, irlocnum); m->clear(); } } } // if inactive: update fields but do not give any contribution to the structure if ( !this->isActivated(tStep) ) { answer.zero(); return; } }
void MicroplaneMaterial_Bazant :: giveRealStressVector(FloatArray &answer, GaussPoint *gp, const FloatArray &totalStrain, TimeStep *tStep) { int i, mPlaneIndex, mPlaneIndex1; double SvDash, SvSum = 0.; double SD; FloatArray mPlaneNormalStress(numberOfMicroplanes), mPlaneShear_L_Stress(numberOfMicroplanes), mPlaneShear_M_Stress(numberOfMicroplanes); double mPlaneIntegrationWeight; Microplane *mPlane; FloatArray mPlaneStressCmpns, mPlaneStrainCmpns; FloatArray stressIncrement; answer.resize(6); answer.zero(); StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( this->giveStatus(gp) ); this->initTempStatus(gp); for ( mPlaneIndex = 0; mPlaneIndex < numberOfMicroplanes; mPlaneIndex++ ) { mPlane = this->giveMicroplane(mPlaneIndex, gp); mPlaneIndex1 = mPlaneIndex + 1; // compute strain projections on mPlaneIndex-th microplane computeStrainVectorComponents(mPlaneStrainCmpns, mPlane, totalStrain); // compute real stresses on this microplane giveRealMicroplaneStressVector(mPlaneStressCmpns, mPlane, mPlaneStrainCmpns, tStep); mPlaneNormalStress.at(mPlaneIndex1) = mPlaneStressCmpns.at(2); mPlaneShear_L_Stress.at(mPlaneIndex1) = mPlaneStressCmpns.at(3); mPlaneShear_M_Stress.at(mPlaneIndex1) = mPlaneStressCmpns.at(4); mPlaneIntegrationWeight = this->giveMicroplaneIntegrationWeight(mPlane); SvSum += mPlaneNormalStress.at(mPlaneIndex1) * mPlaneIntegrationWeight; SD = mPlaneNormalStress.at(mPlaneIndex1) - mPlaneStressCmpns.at(1); //SDSum += SD* mPlaneIntegrationWeight; // perform homogenization // mPlaneStressCmpns.at(1) je SVdash // mPlaneStressCmpns.at(2) je SN // mPlaneStressCmpns.at(3) je SL // mPlaneStressCmpns.at(4) je SM // answer (1 az 6) for ( i = 0; i < 6; i++ ) { answer.at(i + 1) += ( ( N [ mPlaneIndex ] [ i ] - Kronecker [ i ] / 3. ) * SD + L [ mPlaneIndex ] [ i ] * mPlaneShear_L_Stress.at(mPlaneIndex1) + M [ mPlaneIndex ] [ i ] * mPlaneShear_M_Stress.at(mPlaneIndex1) ) * mPlaneIntegrationWeight; } } SvSum = SvSum * 6.; //nakonec answer take *6 SvDash = mPlaneStressCmpns.at(1); //volumetric stress is the same for all mplanes //and does not need to be homogenized . //Only updating accordinging to mean normal stress must be done. //Use updateVolumetricStressTo() if necessary // sv=min(integr(sn)/2PI,SvDash) if ( SvDash > SvSum / 3. ) { SvDash = SvSum / 3.; answer.zero(); for ( mPlaneIndex = 0; mPlaneIndex < numberOfMicroplanes; mPlaneIndex++ ) { mPlane = this->giveMicroplane(mPlaneIndex, gp); mPlaneIndex1 = mPlaneIndex + 1; updateVolumetricStressTo(mPlane, SvDash); SD = mPlaneNormalStress.at(mPlaneIndex1) - SvDash; mPlaneIntegrationWeight = this->giveMicroplaneIntegrationWeight(mPlane); for ( i = 0; i < 6; i++ ) { answer.at(i + 1) += ( ( N [ mPlaneIndex ] [ i ] - Kronecker [ i ] / 3. ) * SD + L [ mPlaneIndex ] [ i ] * mPlaneShear_L_Stress.at(mPlaneIndex1) + M [ mPlaneIndex ] [ i ] * mPlaneShear_M_Stress.at(mPlaneIndex1) ) * mPlaneIntegrationWeight; } } } answer.times(6.0); //2nd constraint, addition of volumetric part answer.at(1) += SvDash; answer.at(2) += SvDash; answer.at(3) += SvDash; // uncomment this //status -> letStrainIncrementVectorBe (reducedStrainIncrement); status->letTempStrainVectorBe(totalStrain); // uncomment this // stressIncrement = answer; // crossSection->giveReducedCharacteristicVector(stressIncrement, gp, answer); // stressIncrement.subtract (status -> giveStressVector()); // status -> letStressIncrementVectorBe (stressIncrement); status->letTempStressVectorBe(answer); }
void PLHoopStressCirc :: propagateInterfaces(Domain &iDomain, EnrichmentDomain &ioEnrDom) { // Fetch crack tip data TipInfo tipInfoStart, tipInfoEnd; ioEnrDom.giveTipInfos(tipInfoStart, tipInfoEnd); std :: vector< TipInfo >tipInfo = {tipInfoStart, tipInfoEnd}; SpatialLocalizer *localizer = iDomain.giveSpatialLocalizer(); for ( size_t tipIndex = 0; tipIndex < tipInfo.size(); tipIndex++ ) { // Construct circle points on an arc from -90 to 90 degrees double angle = -90.0 + mAngleInc; std :: vector< double >angles; while ( angle <= ( 90.0 - mAngleInc ) ) { angles.push_back(angle * M_PI / 180.0); angle += mAngleInc; } const FloatArray &xT = tipInfo [ tipIndex ].mGlobalCoord; const FloatArray &t = tipInfo [ tipIndex ].mTangDir; const FloatArray &n = tipInfo [ tipIndex ].mNormalDir; // It is meaningless to propagate a tip that is not inside any element Element *el = localizer->giveElementContainingPoint(tipInfo [ tipIndex ].mGlobalCoord); if ( el != NULL ) { std :: vector< FloatArray >circPoints; for ( size_t i = 0; i < angles.size(); i++ ) { FloatArray tangent(2); tangent.zero(); tangent.add(cos(angles [ i ]), t); tangent.add(sin(angles [ i ]), n); tangent.normalize(); FloatArray x(xT); x.add(mRadius, tangent); circPoints.push_back(x); } std :: vector< double >sigTTArray, sigRTArray; // Loop over circle points for ( size_t pointIndex = 0; pointIndex < circPoints.size(); pointIndex++ ) { FloatArray stressVec; if ( mUseRadialBasisFunc ) { // Interpolate stress with radial basis functions // Choose a cut-off length l: // take the distance between two nodes in the element containing the // crack tip multiplied by a constant factor. // ( This choice implies that we hope that the element has reasonable // aspect ratio.) const FloatArray &x1 = * ( el->giveDofManager(1)->giveCoordinates() ); const FloatArray &x2 = * ( el->giveDofManager(2)->giveCoordinates() ); const double l = 1.0 * x1.distance(x2); // Use the octree to get all elements that have // at least one Gauss point in a certain region around the tip. const double searchRadius = 3.0 * l; std :: set< int >elIndices; localizer->giveAllElementsWithIpWithinBox(elIndices, circPoints [ pointIndex ], searchRadius); // Loop over the elements and Gauss points obtained. // Evaluate the interpolation. FloatArray sumQiWiVi; double sumWiVi = 0.0; for ( int elIndex: elIndices ) { Element *gpEl = iDomain.giveElement(elIndex); IntegrationRule *iRule = gpEl->giveDefaultIntegrationRulePtr(); for ( GaussPoint *gp_i: *iRule ) { //////////////////////////////////////// // Compute global gp coordinates FloatArray N; FEInterpolation *interp = gpEl->giveInterpolation(); interp->evalN( N, * ( gp_i->giveCoordinates() ), FEIElementGeometryWrapper(gpEl) ); // Compute global coordinates of Gauss point FloatArray globalCoord(2); globalCoord.zero(); for ( int i = 1; i <= gpEl->giveNumberOfDofManagers(); i++ ) { DofManager *dMan = gpEl->giveDofManager(i); globalCoord.at(1) += N.at(i) * dMan->giveCoordinate(1); globalCoord.at(2) += N.at(i) * dMan->giveCoordinate(2); } //////////////////////////////////////// // Compute weight of kernel function FloatArray tipToGP; tipToGP.beDifferenceOf(globalCoord, xT); bool inFrontOfCrack = true; if ( tipToGP.dotProduct(t) < 0.0 ) { inFrontOfCrack = false; } double r = circPoints [ pointIndex ].distance(globalCoord); if ( r < l && inFrontOfCrack ) { double w = ( ( l - r ) / ( pow(2.0 * M_PI, 1.5) * pow(l, 3) ) ) * exp( -0.5 * pow(r, 2) / pow(l, 2) ); // Compute gp volume double V = gpEl->computeVolumeAround(gp_i); // Get stress StructuralMaterialStatus *ms = dynamic_cast< StructuralMaterialStatus * >( gp_i->giveMaterialStatus() ); if ( ms == NULL ) { OOFEM_ERROR("failed to fetch MaterialStatus."); } FloatArray stressVecGP = ms->giveStressVector(); if ( sumQiWiVi.giveSize() != stressVecGP.giveSize() ) { sumQiWiVi.resize( stressVecGP.giveSize() ); sumQiWiVi.zero(); } // Add to numerator sumQiWiVi.add(w * V, stressVecGP); // Add to denominator sumWiVi += w * V; } } } if ( fabs(sumWiVi) > 1.0e-12 ) { stressVec.beScaled(1.0 / sumWiVi, sumQiWiVi); } else { // Take stress from closest Gauss point int region = 1; bool useCZGP = false; GaussPoint &gp = * ( localizer->giveClosestIP(circPoints [ pointIndex ], region, useCZGP) ); // Compute stresses StructuralMaterialStatus *ms = dynamic_cast< StructuralMaterialStatus * >( gp.giveMaterialStatus() ); if ( ms == NULL ) { OOFEM_ERROR("failed to fetch MaterialStatus."); } stressVec = ms->giveStressVector(); } } else { // Take stress from closest Gauss point int region = 1; bool useCZGP = false; GaussPoint &gp = * ( localizer->giveClosestIP(circPoints [ pointIndex ], region, useCZGP) ); // Compute stresses StructuralMaterialStatus *ms = dynamic_cast< StructuralMaterialStatus * >( gp.giveMaterialStatus() ); if ( ms == NULL ) { OOFEM_ERROR("failed to fetch MaterialStatus."); } stressVec = ms->giveStressVector(); } FloatMatrix stress(2, 2); int shearPos = stressVec.giveSize(); stress.at(1, 1) = stressVec.at(1); stress.at(1, 2) = stressVec.at(shearPos); stress.at(2, 1) = stressVec.at(shearPos); stress.at(2, 2) = stressVec.at(2); // Rotation matrix FloatMatrix rot(2, 2); rot.at(1, 1) = cos(angles [ pointIndex ]); rot.at(1, 2) = -sin(angles [ pointIndex ]); rot.at(2, 1) = sin(angles [ pointIndex ]); rot.at(2, 2) = cos(angles [ pointIndex ]); FloatArray tRot, nRot; tRot.beProductOf(rot, t); nRot.beProductOf(rot, n); FloatMatrix rotTot(2, 2); rotTot.setColumn(tRot, 1); rotTot.setColumn(nRot, 2); FloatMatrix tmp, stressRot; tmp.beTProductOf(rotTot, stress); stressRot.beProductOf(tmp, rotTot); const double sigThetaTheta = stressRot.at(2, 2); sigTTArray.push_back(sigThetaTheta); const double sigRTheta = stressRot.at(1, 2); sigRTArray.push_back(sigRTheta); } ////////////////////////////// // Compute propagation angle // Find angles that fulfill sigRT = 0 const double stressTol = 1.0e-9; double maxSigTT = 0.0, maxAngle = 0.0; bool foundZeroLevel = false; for ( size_t segIndex = 0; segIndex < ( circPoints.size() - 1 ); segIndex++ ) { // If the shear stress sigRT changes sign over the segment if ( sigRTArray [ segIndex ] * sigRTArray [ segIndex + 1 ] < stressTol ) { // Compute location of zero level double xi = EnrichmentItem :: calcXiZeroLevel(sigRTArray [ segIndex ], sigRTArray [ segIndex + 1 ]); double theta = 0.5 * ( 1.0 - xi ) * angles [ segIndex ] + 0.5 * ( 1.0 + xi ) * angles [ segIndex + 1 ]; double sigThetaTheta = 0.5 * ( 1.0 - xi ) * sigTTArray [ segIndex ] + 0.5 * ( 1.0 + xi ) * sigTTArray [ segIndex + 1 ]; // printf("Found candidate: theta: %e sigThetaTheta: %e\n", theta, sigThetaTheta); if ( sigThetaTheta > maxSigTT ) { foundZeroLevel = true; maxSigTT = sigThetaTheta; maxAngle = theta; } } } if ( !foundZeroLevel ) { printf("No zero level was found.\n"); } if ( iDomain.giveXfemManager()->giveVtkDebug() ) { XFEMDebugTools :: WriteArrayToMatlab("sigTTvsAngle.m", angles, sigTTArray); XFEMDebugTools :: WriteArrayToMatlab("sigRTvsAngle.m", angles, sigRTArray); XFEMDebugTools :: WriteArrayToGnuplot("sigTTvsAngle.dat", angles, sigTTArray); XFEMDebugTools :: WriteArrayToGnuplot("sigRTvsAngle.dat", angles, sigRTArray); } // Compare with threshold if ( maxSigTT > mHoopStressThreshold && foundZeroLevel ) { // Rotation matrix FloatMatrix rot(2, 2); rot.at(1, 1) = cos(maxAngle); rot.at(1, 2) = -sin(maxAngle); rot.at(2, 1) = sin(maxAngle); rot.at(2, 2) = cos(maxAngle); FloatArray dir; dir.beProductOf(rot, tipInfo [ tipIndex ].mTangDir); // Fill up struct std :: vector< TipPropagation >tipPropagations; TipPropagation tipProp; tipProp.mTipIndex = tipIndex; tipProp.mPropagationDir = dir; tipProp.mPropagationLength = mIncrementLength; tipPropagations.push_back(tipProp); // Propagate ioEnrDom.propagateTips(tipPropagations); } } } }