diff --git a/code/constitutive_dislotwin.f90 b/code/constitutive_dislotwin.f90 index ebb2a266e..12d6e5747 100644 --- a/code/constitutive_dislotwin.f90 +++ b/code/constitutive_dislotwin.f90 @@ -2,237 +2,246 @@ !************************************ !* Module: CONSTITUTIVE * !************************************ -!* contains: * -!* - constitutive equations * -!* - parameters definition * -!* - orientations * -!************************************ - -! [TWIP steel FeMnC] - -! C11 175.0e9 # elastic constants in Pa -! C12 115.0e9 -! C44 135.0e9 -! lattice_structure fcc -! Nslip 12 -! Ntwin 12 -! constitution dislotwin -! (output) dislocationdensity -! (output) shearrate_slip -! (output) mfp_slip # mean free path -! (output) resolvedstress_slip -! (output) resistance_slip # passing stress -! (output) volumefraction -! (output) shearrate_twin -! (output) mfp_twin # mean free path -! (output) resolvedstress_twin -! (output) resistance_twin # "nucleation barrier" - -! ### dislocation density-based constitutive parameters ### -! burgers 2.56e-10 # Burgers vector [m] -! Qedge 5.5e-19 # Activation energy for dislocation glide [J/K] (0.5*G*b^3) -! grainsize 2.0e-5 # Average grain size [m] -! stacksize 5.0e-8 # Twin stack mean thickness [m] -! interaction_slipslip 1.0 2.2 3.0 1.6 3.8 4.5 # Dislocation interaction coefficients -! interaction_sliptwin 0.0 1.0 # Dislocation interaction coefficients -! interaction_twintwin 0.0 1.0 # Dislocation interaction coefficients -! # dislocation glide -! rho0 2.5e12 # Initial dislocation density [m/m³] -! Cmfpslip 1.0 # Adjustable parameter controlling dislocation mean free path -! Cactivolume 1.0 # Adjustable parameter controlling activation volume -! Cthresholdslip 0.1 # Adjustable parameter controlling threshold stress for dislocation motion -! Cstorage 0.02 # Adjustable parameter controlling dislocation storage -! Carecovery 15.0 # Adjustable parameter controlling athermal recovery -! # mechanical twinning -! Ndot0 0.0 # Number of potential twin source per volume per time [1/m³.s] -! fmax 1.0 # Maximum admissible twin volume fraction -! Cmfptwin 1.0 # Adjustable parameter controlling twin mean free path -! Cthresholdtwin 1.0 # Adjustable parameter controlling threshold stress for deformation twinning - MODULE constitutive_dislotwin -!*** Include other modules *** - use prec, only: pReal,pInt - implicit none +!* Include other modules +use prec, only: pReal,pInt +implicit none - character (len=*), parameter :: constitutive_dislotwin_label = 'dislotwin' - - integer(pInt), dimension(:), allocatable :: constitutive_dislotwin_sizeDotState, & - constitutive_dislotwin_sizeState, & - constitutive_dislotwin_sizePostResults - integer(pInt), dimension(:,:), allocatable,target :: constitutive_dislotwin_sizePostResult - character(len=64), dimension(:,:), allocatable,target :: constitutive_dislotwin_output - - character(len=32), dimension(:), allocatable :: constitutive_dislotwin_structureName - integer(pInt), dimension(:), allocatable :: constitutive_dislotwin_structure, & - constitutive_dislotwin_totalNslip, & - constitutive_dislotwin_totalNtwin - integer(pInt), dimension(:,:), allocatable :: constitutive_dislotwin_Nslip, & - constitutive_dislotwin_Ntwin, & - constitutive_dislotwin_slipFamily, & - constitutive_dislotwin_twinFamily - - real(pReal), dimension(:), allocatable :: constitutive_dislotwin_CoverA, & - constitutive_dislotwin_C11, & - constitutive_dislotwin_C12, & - constitutive_dislotwin_C13, & - constitutive_dislotwin_C33, & - constitutive_dislotwin_C44, & - constitutive_dislotwin_Gmod - real(pReal), dimension(:,:,:), allocatable :: constitutive_dislotwin_Cslip_66 - real(pReal), dimension(:,:,:,:), allocatable :: constitutive_dislotwin_Ctwin_66 - real(pReal), dimension(:,:,:,:,:), allocatable :: constitutive_dislotwin_Cslip_3333 - real(pReal), dimension(:,:,:,:,:,:), allocatable :: constitutive_dislotwin_Ctwin_3333 - real(pReal), dimension(:,:), allocatable :: constitutive_dislotwin_rho0, & - constitutive_dislotwin_Burgers, & - constitutive_dislotwin_Qedge, & - constitutive_dislotwin_stacksize, & - constitutive_dislotwin_Ndot0, & - - constitutive_dislotwin_interaction_slipslip, & - constitutive_dislotwin_interaction_sliptwin, & - constitutive_dislotwin_interaction_twinslip, & - constitutive_dislotwin_interaction_twintwin - real(pReal), dimension(:), allocatable :: constitutive_dislotwin_grainsize, & - constitutive_dislotwin_fmax, & - constitutive_dislotwin_Cmfpslip, & - constitutive_dislotwin_Cmfptwin, & - constitutive_dislotwin_Cthresholdslip, & - constitutive_dislotwin_Cthresholdtwin, & - constitutive_dislotwin_Cactivolume, & - constitutive_dislotwin_Carecovery, & - constitutive_dislotwin_Cstorage - - - real(pReal), dimension(:,:,:), allocatable :: constitutive_dislotwin_parall_interaction, & - constitutive_dislotwin_forest_interaction, & - constitutive_dislotwin_hardeningMatrix_sliptwin, & - constitutive_dislotwin_hardeningMatrix_twinslip, & - constitutive_dislotwin_hardeningMatrix_twintwin - -!************************************* -!* Definition of material properties * -!************************************* -!* Physical parameter, attack_frequency != Debye frequency -real(pReal), parameter :: attack_frequency = 1.0e10_pReal -!* Physical parameter, Boltzmann constant in J/Kelvin -real(pReal), parameter :: kB = 1.38e-23_pReal -!* Physical parameter, Avogadro number in 1/mol -real(pReal), parameter :: avogadro = 6.022e23_pReal -!* Physical parameter, Gas constant in J.mol/Kelvin -real(pReal), parameter :: Rgas = 8.314_pReal +!* Lists of states and physical parameters +character(len=*), parameter :: constitutive_dislotwin_label = 'dislotwin' +character(len=18), dimension(2), parameter:: constitutive_dislotwin_listBasicSlipStates = (/'rhoEdge ', & + 'rhoEdgeDip'/) +character(len=18), dimension(1), parameter:: constitutive_dislotwin_listBasicTwinStates = (/'twinFraction'/) +character(len=18), dimension(4), parameter:: constitutive_dislotwin_listDependentSlipStates =(/'invLambdaSlip ', & + 'invLambdaSlipTwin', & + 'meanFreePathSlip ', & + 'tauSlipThreshold '/) +character(len=18), dimension(4), parameter:: constitutive_dislotwin_listDependentTwinStates =(/'invLambdaTwin ', & + 'meanFreePathTwin', & + 'tauTwinThreshold', & + 'twinVolume '/) +real(pReal), parameter :: kB = 1.38e-23_pReal ! Boltzmann constant in J/Kelvin +!* Definition of global variables +integer(pInt), dimension(:), allocatable :: constitutive_dislotwin_sizeDotState, & ! number of dotStates + constitutive_dislotwin_sizeState, & ! total number of microstructural state variables + constitutive_dislotwin_sizePostResults ! cumulative size of post results +integer(pInt), dimension(:,:), allocatable, target :: constitutive_dislotwin_sizePostResult ! size of each post result output +character(len=64), dimension(:,:), allocatable, target :: constitutive_dislotwin_output ! name of each post result output +character(len=32), dimension(:), allocatable :: constitutive_dislotwin_structureName ! name of the lattice structure +integer(pInt), dimension(:), allocatable :: constitutive_dislotwin_structure, & ! number representing the kind of lattice structure + constitutive_dislotwin_totalNslip, & ! total number of active slip systems for each instance + constitutive_dislotwin_totalNtwin ! total number of active twin systems for each instance +integer(pInt), dimension(:,:), allocatable :: constitutive_dislotwin_Nslip, & ! number of active slip systems for each family and instance + constitutive_dislotwin_Ntwin, & ! number of active twin systems for each family and instance + constitutive_dislotwin_slipFamily, & ! lookup table relating active slip system to slip family for each instance + constitutive_dislotwin_twinFamily, & ! lookup table relating active twin system to twin family for each instance + constitutive_dislotwin_slipSystemLattice, & ! lookup table relating active slip system index to lattice slip system index for each instance + constitutive_dislotwin_twinSystemLattice ! lookup table relating active twin system index to lattice twin system index for each instance +real(pReal), dimension(:), allocatable :: constitutive_dislotwin_CoverA, & ! c/a ratio for hex type lattice + constitutive_dislotwin_C11, & ! C11 element in elasticity matrix + constitutive_dislotwin_C12, & ! C12 element in elasticity matrix + constitutive_dislotwin_C13, & ! C13 element in elasticity matrix + constitutive_dislotwin_C33, & ! C33 element in elasticity matrix + constitutive_dislotwin_C44, & ! C44 element in elasticity matrix + constitutive_dislotwin_Gmod, & ! shear modulus + constitutive_dislotwin_CAtomicVolume, & ! atomic volume in Bugers vector unit + constitutive_dislotwin_D0, & ! prefactor for self-diffusion coefficient + constitutive_dislotwin_Qsd, & ! activation energy for dislocation climb + constitutive_dislotwin_GrainSize, & ! grain size + constitutive_dislotwin_p, & ! p-exponent in glide velocity + constitutive_dislotwin_q, & ! q-exponent in glide velocity + constitutive_dislotwin_MaxTwinFraction, & ! maximum allowed total twin volume fraction + constitutive_dislotwin_r, & ! r-exponent in twin nucleation rate + constitutive_dislotwin_CEdgeDipMinDistance, & ! + constitutive_dislotwin_Cmfptwin, & ! + constitutive_dislotwin_Cthresholdtwin, & ! + constitutive_dislotwin_relevantRho ! dislocation density considered relevant +real(pReal), dimension(:,:,:), allocatable :: constitutive_dislotwin_Cslip_66 ! elasticity matrix in Mandel notation for each instance +real(pReal), dimension(:,:,:,:), allocatable :: constitutive_dislotwin_Ctwin_66 ! twin elasticity matrix in Mandel notation for each instance +real(pReal), dimension(:,:,:,:,:), allocatable :: constitutive_dislotwin_Cslip_3333 ! elasticity matrix for each instance +real(pReal), dimension(:,:,:,:,:,:), allocatable :: constitutive_dislotwin_Ctwin_3333 ! twin elasticity matrix for each instance +real(pReal), dimension(:,:), allocatable :: constitutive_dislotwin_rhoEdge0, & ! initial edge dislocation density per slip system for each family and instance + constitutive_dislotwin_rhoEdgeDip0, & ! initial edge dipole density per slip system for each family and instance + constitutive_dislotwin_burgersPerSlipFamily, & ! absolute length of burgers vector [m] for each slip family and instance + constitutive_dislotwin_burgersPerSlipSystem, & ! absolute length of burgers vector [m] for each slip system and instance + constitutive_dislotwin_burgersPerTwinFamily, & ! absolute length of burgers vector [m] for each twin family and instance + constitutive_dislotwin_burgersPerTwinSystem, & ! absolute length of burgers vector [m] for each twin system and instance + constitutive_dislotwin_QedgePerSlipFamily, & ! activation energy for glide [J] for each slip family and instance + constitutive_dislotwin_QedgePerSlipSystem, & ! activation energy for glide [J] for each slip system and instance + constitutive_dislotwin_v0PerSlipFamily, & ! dislocation velocity prefactor [m/s] for each family and instance + constitutive_dislotwin_v0PerSlipSystem, & ! dislocation velocity prefactor [m/s] for each slip system and instance + constitutive_dislotwin_Ndot0PerTwinFamily, & ! twin nucleation rate [1/m³s] for each twin family and instance + constitutive_dislotwin_Ndot0PerTwinSystem, & ! twin nucleation rate [1/m³s] for each twin system and instance + constitutive_dislotwin_twinsizePerTwinFamily, & ! twin thickness [m] for each twin family and instance + constitutive_dislotwin_twinsizePerTwinSystem, & ! twin thickness [m] for each twin system and instance + constitutive_dislotwin_CLambdaSlipPerSlipFamily, & ! Adj. parameter for distance between 2 forest dislocations for each slip family and instance + constitutive_dislotwin_CLambdaSlipPerSlipSystem, & ! Adj. parameter for distance between 2 forest dislocations for each slip system and instance + constitutive_dislotwin_interactionSlipSlip, & ! coefficients for slip-slip interaction for each interaction type and instance + constitutive_dislotwin_interactionSlipTwin, & ! coefficients for slip-twin interaction for each interaction type and instance + constitutive_dislotwin_interactionTwinSlip, & ! coefficients for twin-slip interaction for each interaction type and instance + constitutive_dislotwin_interactionTwinTwin ! coefficients for twin-twin interaction for each interaction type and instance +real(pReal), dimension(:,:,:), allocatable :: constitutive_dislotwin_interactionMatrixSlipSlip, & ! interaction matrix of the different slip systems for each instance + constitutive_dislotwin_interactionMatrixSlipTwin, & ! interaction matrix of slip systems with twin systems for each instance + constitutive_dislotwin_interactionMatrixTwinSlip, & ! interaction matrix of twin systems with slip systems for each instance + constitutive_dislotwin_interactionMatrixTwinTwin, & ! interaction matrix of the different twin systems for each instance + constitutive_dislotwin_forestProjectionEdge ! matrix of forest projections of edge dislocations for each instance CONTAINS !**************************************** -!* - constitutive_init -!* - constitutive_stateInit -!* - constitutive_homogenizedC -!* - constitutive_microstructure -!* - constitutive_LpAndItsTangent -!* - consistutive_dotState -!* - constitutive_dotTemperature -!* - consistutive_postResults +!* - constitutive_dislotwin_init +!* - constitutive_dislotwin_stateInit +!* - constitutive_dislotwin_relevantState +!* - constitutive_dislotwin_homogenizedC +!* - constitutive_dislotwin_microstructure +!* - constitutive_dislotwin_LpAndItsTangent +!* - consistutive_dislotwin_dotState +!* - constitutive_dislotwin_dotTemperature +!* - consistutive_dislotwin_postResults !**************************************** - subroutine constitutive_dislotwin_init(file) !************************************** !* Module initialization * !************************************** - use prec, only: pInt,pReal - use math, only: math_Mandel3333to66,math_Voigt66to3333,math_mul3x3 - use IO - use material - use lattice +use prec, only: pInt,pReal +use math, only: math_Mandel3333to66,math_Voigt66to3333,math_mul3x3 +use IO +use material +use lattice - integer(pInt), intent(in) :: file - integer(pInt), parameter :: maxNchunks = 21 - integer(pInt), dimension(1+2*maxNchunks) :: positions - integer(pInt) section,maxNinstance,i,j,k,l,m,n,o,p,q,r,s,output,mySize - character(len=64) tag - character(len=1024) line - real(pReal) x,y +!* Input variables +integer(pInt), intent(in) :: file +!* Local variables +integer(pInt), parameter :: maxNchunks = 21 +integer(pInt), dimension(1+2*maxNchunks) :: positions +integer(pInt) section,maxNinstance,f,i,j,k,l,m,n,o,p,q,r,s,s1,s2,t1,t2,ns,nt,output,mySize,myStructure,maxTotalNslip,maxTotalNtwin +character(len=64) tag +character(len=1024) line - write(6,*) - write(6,'(a20,a20,a12)') '<<<+- constitutive_',constitutive_dislotwin_label,' init -+>>>' - write(6,*) '$Id: constitutive_dislotwin.f90 412 2009-09-18 15:37:14Z MPIE\c.kords $' - write(6,*) +!write(6,*) +!write(6,'(a20,a20,a12)') '<<<+- constitutive_',constitutive_dislotwin_label,' init -+>>>' +!write(6,*) '$Id: constitutive_dislotwin.f90 412 2009-09-18 15:37:14Z MPIE\c.kords $' +!write(6,*) - maxNinstance = count(phase_constitution == constitutive_dislotwin_label) - if (maxNinstance == 0) return +maxNinstance = count(phase_constitution == constitutive_dislotwin_label) +if (maxNinstance == 0) return + +!* Space allocation for global variables +allocate(constitutive_dislotwin_sizeDotState(maxNinstance)) +allocate(constitutive_dislotwin_sizeState(maxNinstance)) +allocate(constitutive_dislotwin_sizePostResults(maxNinstance)) +allocate(constitutive_dislotwin_sizePostResult(maxval(phase_Noutput),maxNinstance)) +allocate(constitutive_dislotwin_output(maxval(phase_Noutput),maxNinstance)) +constitutive_dislotwin_sizeDotState = 0_pInt +constitutive_dislotwin_sizeState = 0_pInt +constitutive_dislotwin_sizePostResults = 0_pInt +constitutive_dislotwin_sizePostResult = 0_pInt +constitutive_dislotwin_output = '' - allocate(constitutive_dislotwin_sizeDotState(maxNinstance)) ; constitutive_dislotwin_sizeDotState = 0_pInt - allocate(constitutive_dislotwin_sizeState(maxNinstance)) ; constitutive_dislotwin_sizeState = 0_pInt - allocate(constitutive_dislotwin_sizePostResults(maxNinstance)) ; constitutive_dislotwin_sizePostResults = 0_pInt - allocate(constitutive_dislotwin_sizePostResult(maxval(phase_Noutput), & - maxNinstance)) ; constitutive_dislotwin_sizePostResult = 0_pInt - allocate(constitutive_dislotwin_output(maxval(phase_Noutput), & - maxNinstance)) ; constitutive_dislotwin_output = '' +allocate(constitutive_dislotwin_structureName(maxNinstance)) +allocate(constitutive_dislotwin_structure(maxNinstance)) +allocate(constitutive_dislotwin_Nslip(lattice_maxNslipFamily,maxNinstance)) +allocate(constitutive_dislotwin_Ntwin(lattice_maxNtwinFamily,maxNinstance)) +allocate(constitutive_dislotwin_slipFamily(lattice_maxNslip,maxNinstance)) +allocate(constitutive_dislotwin_twinFamily(lattice_maxNtwin,maxNinstance)) +allocate(constitutive_dislotwin_slipSystemLattice(lattice_maxNslip,maxNinstance)) +allocate(constitutive_dislotwin_twinSystemLattice(lattice_maxNtwin,maxNinstance)) +allocate(constitutive_dislotwin_totalNslip(maxNinstance)) +allocate(constitutive_dislotwin_totalNtwin(maxNinstance)) +constitutive_dislotwin_structureName = '' +constitutive_dislotwin_structure = 0_pInt +constitutive_dislotwin_Nslip = 0_pInt +constitutive_dislotwin_Ntwin = 0_pInt +constitutive_dislotwin_slipFamily = 0_pInt +constitutive_dislotwin_twinFamily = 0_pInt +constitutive_dislotwin_slipSystemLattice = 0.0_pReal +constitutive_dislotwin_twinSystemLattice = 0.0_pReal +constitutive_dislotwin_totalNslip = 0_pInt +constitutive_dislotwin_totalNtwin = 0_pInt - allocate(constitutive_dislotwin_structureName(maxNinstance)) ; constitutive_dislotwin_structureName = '' - allocate(constitutive_dislotwin_structure(maxNinstance)) ; constitutive_dislotwin_structure = 0_pInt - allocate(constitutive_dislotwin_Nslip(lattice_maxNslipFamily,& - maxNinstance)) ; constitutive_dislotwin_Nslip = 0_pInt - allocate(constitutive_dislotwin_Ntwin(lattice_maxNtwinFamily,& - maxNinstance)) ; constitutive_dislotwin_Ntwin = 0_pInt +allocate(constitutive_dislotwin_CoverA(maxNinstance)) +allocate(constitutive_dislotwin_C11(maxNinstance)) +allocate(constitutive_dislotwin_C12(maxNinstance)) +allocate(constitutive_dislotwin_C13(maxNinstance)) +allocate(constitutive_dislotwin_C33(maxNinstance)) +allocate(constitutive_dislotwin_C44(maxNinstance)) +allocate(constitutive_dislotwin_Gmod(maxNinstance)) +allocate(constitutive_dislotwin_CAtomicVolume(maxNinstance)) +allocate(constitutive_dislotwin_D0(maxNinstance)) +allocate(constitutive_dislotwin_Qsd(maxNinstance)) +allocate(constitutive_dislotwin_GrainSize(maxNinstance)) +allocate(constitutive_dislotwin_p(maxNinstance)) +allocate(constitutive_dislotwin_q(maxNinstance)) +allocate(constitutive_dislotwin_MaxTwinFraction(maxNinstance)) +allocate(constitutive_dislotwin_r(maxNinstance)) +allocate(constitutive_dislotwin_CEdgeDipMinDistance(maxNinstance)) +allocate(constitutive_dislotwin_Cmfptwin(maxNinstance)) +allocate(constitutive_dislotwin_Cthresholdtwin(maxNinstance)) +allocate(constitutive_dislotwin_relevantRho(maxNinstance)) +allocate(constitutive_dislotwin_Cslip_66(6,6,maxNinstance)) +allocate(constitutive_dislotwin_Cslip_3333(3,3,3,3,maxNinstance)) +constitutive_dislotwin_CoverA = 0.0_pReal +constitutive_dislotwin_C11 = 0.0_pReal +constitutive_dislotwin_C12 = 0.0_pReal +constitutive_dislotwin_C13 = 0.0_pReal +constitutive_dislotwin_C33 = 0.0_pReal +constitutive_dislotwin_C44 = 0.0_pReal +constitutive_dislotwin_Gmod = 0.0_pReal +constitutive_dislotwin_CAtomicVolume = 0.0_pReal +constitutive_dislotwin_D0 = 0.0_pReal +constitutive_dislotwin_Qsd = 0.0_pReal +constitutive_dislotwin_GrainSize = 0.0_pReal +constitutive_dislotwin_p = 0.0_pReal +constitutive_dislotwin_q = 0.0_pReal +constitutive_dislotwin_MaxTwinFraction = 0.0_pReal +constitutive_dislotwin_r = 0.0_pReal +constitutive_dislotwin_CEdgeDipMinDistance = 0.0_pReal +constitutive_dislotwin_Cmfptwin = 0.0_pReal +constitutive_dislotwin_Cthresholdtwin = 0.0_pReal +constitutive_dislotwin_relevantRho = 0.0_pReal +constitutive_dislotwin_Cslip_66 = 0.0_pReal +constitutive_dislotwin_Cslip_3333 = 0.0_pReal + +allocate(constitutive_dislotwin_rhoEdge0(lattice_maxNslipFamily,maxNinstance)) +allocate(constitutive_dislotwin_rhoEdgeDip0(lattice_maxNslipFamily,maxNinstance)) +allocate(constitutive_dislotwin_burgersPerSlipFamily(lattice_maxNslipFamily,maxNinstance)) +allocate(constitutive_dislotwin_burgersPerTwinFamily(lattice_maxNtwinFamily,maxNinstance)) +allocate(constitutive_dislotwin_QedgePerSlipFamily(lattice_maxNslipFamily,maxNinstance)) +allocate(constitutive_dislotwin_v0PerSlipFamily(lattice_maxNslipFamily,maxNinstance)) +allocate(constitutive_dislotwin_Ndot0PerTwinFamily(lattice_maxNtwinFamily,maxNinstance)) +allocate(constitutive_dislotwin_twinsizePerTwinFamily(lattice_maxNtwinFamily,maxNinstance)) +allocate(constitutive_dislotwin_CLambdaSlipPerSlipFamily(lattice_maxNslipFamily,maxNinstance)) +constitutive_dislotwin_rhoEdge0 = 0.0_pReal +constitutive_dislotwin_rhoEdgeDip0 = 0.0_pReal +constitutive_dislotwin_burgersPerSlipFamily = 0.0_pReal +constitutive_dislotwin_burgersPerTwinFamily = 0.0_pReal +constitutive_dislotwin_QedgePerSlipFamily = 0.0_pReal +constitutive_dislotwin_v0PerSlipFamily = 0.0_pReal +constitutive_dislotwin_Ndot0PerTwinFamily = 0.0_pReal +constitutive_dislotwin_twinsizePerTwinFamily = 0.0_pReal +constitutive_dislotwin_CLambdaSlipPerSlipFamily = 0.0_pReal + +allocate(constitutive_dislotwin_interactionSlipSlip(lattice_maxNinteraction,maxNinstance)) +allocate(constitutive_dislotwin_interactionSlipTwin(lattice_maxNinteraction,maxNinstance)) +allocate(constitutive_dislotwin_interactionTwinSlip(lattice_maxNinteraction,maxNinstance)) +allocate(constitutive_dislotwin_interactionTwinTwin(lattice_maxNinteraction,maxNinstance)) +constitutive_dislotwin_interactionSlipSlip = 0.0_pReal +constitutive_dislotwin_interactionSlipTwin = 0.0_pReal +constitutive_dislotwin_interactionTwinSlip = 0.0_pReal +constitutive_dislotwin_interactionTwinTwin = 0.0_pReal + +!* Readout data from material.config file +rewind(file) +line = '' +section = 0 - allocate(constitutive_dislotwin_slipFamily(lattice_maxNslip,& - maxNinstance)) ; constitutive_dislotwin_slipFamily = 0_pInt - allocate(constitutive_dislotwin_twinFamily(lattice_maxNtwin,& - maxNinstance)) ; constitutive_dislotwin_twinFamily = 0_pInt - - allocate(constitutive_dislotwin_totalNslip(maxNinstance)) ; constitutive_dislotwin_totalNslip = 0_pInt - allocate(constitutive_dislotwin_totalNtwin(maxNinstance)) ; constitutive_dislotwin_totalNtwin = 0_pInt - - allocate(constitutive_dislotwin_CoverA(maxNinstance)) ; constitutive_dislotwin_CoverA = 0.0_pReal - allocate(constitutive_dislotwin_C11(maxNinstance)) ; constitutive_dislotwin_C11 = 0.0_pReal - allocate(constitutive_dislotwin_C12(maxNinstance)) ; constitutive_dislotwin_C12 = 0.0_pReal - allocate(constitutive_dislotwin_C13(maxNinstance)) ; constitutive_dislotwin_C13 = 0.0_pReal - allocate(constitutive_dislotwin_C33(maxNinstance)) ; constitutive_dislotwin_C33 = 0.0_pReal - allocate(constitutive_dislotwin_C44(maxNinstance)) ; constitutive_dislotwin_C44 = 0.0_pReal - allocate(constitutive_dislotwin_Gmod(maxNinstance)) ; constitutive_dislotwin_Gmod = 0.0_pReal - allocate(constitutive_dislotwin_Cslip_66(6,6,maxNinstance)) ; constitutive_dislotwin_Cslip_66 = 0.0_pReal - allocate(constitutive_dislotwin_Cslip_3333(3,3,3,3,maxNinstance)) ; constitutive_dislotwin_Cslip_3333 = 0.0_pReal - - allocate(constitutive_dislotwin_rho0(lattice_maxNslipFamily, & - maxNinstance)) ; constitutive_dislotwin_rho0 = 0.0_pReal - allocate(constitutive_dislotwin_Burgers(lattice_maxNslipFamily, & - maxNinstance)) ; constitutive_dislotwin_Burgers = 0.0_pReal - allocate(constitutive_dislotwin_Qedge(lattice_maxNslipFamily, & - maxNinstance)) ; constitutive_dislotwin_Qedge = 0.0_pReal - allocate(constitutive_dislotwin_grainsize(maxNinstance)) ; constitutive_dislotwin_grainsize = 0.0_pReal - allocate(constitutive_dislotwin_stacksize(lattice_maxNtwinFamily, & - maxNinstance)) ; constitutive_dislotwin_stacksize = 0.0_pReal - allocate(constitutive_dislotwin_fmax(maxNinstance)) ; constitutive_dislotwin_fmax = 0.0_pReal - allocate(constitutive_dislotwin_Ndot0(lattice_maxNtwinFamily, & - maxNinstance)) ; constitutive_dislotwin_Ndot0 = 0.0_pReal - allocate(constitutive_dislotwin_Cmfpslip(maxNinstance)) ; constitutive_dislotwin_Cmfpslip = 0.0_pReal - allocate(constitutive_dislotwin_Cmfptwin(maxNinstance)) ; constitutive_dislotwin_Cmfptwin = 0.0_pReal - allocate(constitutive_dislotwin_Cthresholdslip(maxNinstance)) ; constitutive_dislotwin_Cthresholdslip = 0.0_pReal - allocate(constitutive_dislotwin_Cthresholdtwin(maxNinstance)) ; constitutive_dislotwin_Cthresholdtwin = 0.0_pReal - allocate(constitutive_dislotwin_Cactivolume(maxNinstance)) ; constitutive_dislotwin_Cactivolume = 0.0_pReal - allocate(constitutive_dislotwin_Carecovery(maxNinstance)) ; constitutive_dislotwin_Carecovery = 0.0_pReal - allocate(constitutive_dislotwin_Cstorage(maxNinstance)) ; constitutive_dislotwin_Cstorage = 0.0_pReal - - allocate(constitutive_dislotwin_interaction_slipslip(lattice_maxNinteraction,& - maxNinstance)) ; constitutive_dislotwin_interaction_slipslip = 0.0_pReal - allocate(constitutive_dislotwin_interaction_sliptwin(lattice_maxNinteraction,& - maxNinstance)) ; constitutive_dislotwin_interaction_sliptwin = 0.0_pReal - allocate(constitutive_dislotwin_interaction_twinslip(lattice_maxNinteraction,& - maxNinstance)) ; constitutive_dislotwin_interaction_twinslip = 0.0_pReal - allocate(constitutive_dislotwin_interaction_twintwin(lattice_maxNinteraction,& - maxNinstance)) ; constitutive_dislotwin_interaction_twintwin = 0.0_pReal - - rewind(file) - line = '' - section = 0 - - do while (IO_lc(IO_getTag(line,'<','>')) /= 'phase') ! wind forward to +do while (IO_lc(IO_getTag(line,'<','>')) /= 'phase') ! wind forward to read(file,'(a1024)',END=100) line - enddo +enddo - do ! read thru sections of phase part +do ! read thru sections of phase part read(file,'(a1024)',END=100) line if (IO_isBlank(line)) cycle ! skip empty lines if (IO_getTag(line,'<','>') /= '') exit ! stop at next part @@ -266,131 +275,197 @@ subroutine constitutive_dislotwin_init(file) forall (j = 1:lattice_maxNslipFamily) constitutive_dislotwin_Nslip(j,i) = IO_intValue(line,positions,1+j) case ('ntwin') forall (j = 1:lattice_maxNtwinFamily) constitutive_dislotwin_Ntwin(j,i) = IO_intValue(line,positions,1+j) - case ('rho0') - forall (j = 1:lattice_maxNslipFamily) constitutive_dislotwin_rho0(j,i) = IO_floatValue(line,positions,1+j) - case ('burgers') - forall (j = 1:lattice_maxNslipFamily) constitutive_dislotwin_Burgers(j,i) = IO_floatValue(line,positions,1+j) + case ('rhoedge0') + forall (j = 1:lattice_maxNslipFamily) constitutive_dislotwin_rhoEdge0(j,i) = IO_floatValue(line,positions,1+j) + case ('rhoedgedip0') + forall (j = 1:lattice_maxNslipFamily) constitutive_dislotwin_rhoEdgeDip0(j,i) = IO_floatValue(line,positions,1+j) + case ('slipburgers') + forall (j = 1:lattice_maxNslipFamily) constitutive_dislotwin_burgersPerSlipFamily(j,i) = IO_floatValue(line,positions,1+j) + case ('twinburgers') + forall (j = 1:lattice_maxNtwinFamily) constitutive_dislotwin_burgersPerTwinFamily(j,i) = IO_floatValue(line,positions,1+j) case ('qedge') - forall (j = 1:lattice_maxNslipFamily) constitutive_dislotwin_Qedge(j,i) = IO_floatValue(line,positions,1+j) - case ('grainsize') - constitutive_dislotwin_grainsize(i) = IO_floatValue(line,positions,2) - case ('stacksize') - forall (j = 1:lattice_maxNtwinFamily) constitutive_dislotwin_stacksize(j,i) = IO_floatValue(line,positions,1+j) - case ('fmax') - constitutive_dislotwin_fmax(i) = IO_floatValue(line,positions,2) - case ('ndot0') - forall (j = 1:lattice_maxNtwinFamily) constitutive_dislotwin_Ndot0(j,i) = IO_floatValue(line,positions,1+j) - case ('cmfpslip') - constitutive_dislotwin_Cmfpslip(i) = IO_floatValue(line,positions,2) + forall (j = 1:lattice_maxNslipFamily) constitutive_dislotwin_QedgePerSlipFamily(j,i) = IO_floatValue(line,positions,1+j) + case ('v0') + forall (j = 1:lattice_maxNslipFamily) constitutive_dislotwin_v0PerSlipFamily(j,i) = IO_floatValue(line,positions,1+j) + case ('ndot0') + forall (j = 1:lattice_maxNtwinFamily) constitutive_dislotwin_Ndot0PerTwinFamily(j,i) = IO_floatValue(line,positions,1+j) + case ('twinsize') + forall (j = 1:lattice_maxNtwinFamily) constitutive_dislotwin_twinsizePerTwinFamily(j,i) = IO_floatValue(line,positions,1+j) + case ('clambdaslip') + forall (j = 1:lattice_maxNslipFamily) constitutive_dislotwin_CLambdaSlipPerSlipFamily(j,i) = IO_floatValue(line,positions,1+j) + case ('grainsize') + constitutive_dislotwin_GrainSize(i) = IO_floatValue(line,positions,2) + case ('maxtwinfraction') + constitutive_dislotwin_MaxTwinFraction(i) = IO_floatValue(line,positions,2) + case ('pexponent') + constitutive_dislotwin_p(i) = IO_floatValue(line,positions,2) + case ('qexponent') + constitutive_dislotwin_q(i) = IO_floatValue(line,positions,2) + case ('rexponent') + constitutive_dislotwin_r(i) = IO_floatValue(line,positions,2) + case ('d0') + constitutive_dislotwin_D0(i) = IO_floatValue(line,positions,2) + case ('qsd') + constitutive_dislotwin_Qsd(i) = IO_floatValue(line,positions,2) + case ('relevantrho') + constitutive_dislotwin_relevantRho(i) = IO_floatValue(line,positions,2) case ('cmfptwin') constitutive_dislotwin_Cmfptwin(i) = IO_floatValue(line,positions,2) - case ('cthresholdslip') - constitutive_dislotwin_Cthresholdslip(i) = IO_floatValue(line,positions,2) case ('cthresholdtwin') constitutive_dislotwin_Cthresholdtwin(i) = IO_floatValue(line,positions,2) - case ('cactivolume') - constitutive_dislotwin_Cactivolume(i) = IO_floatValue(line,positions,2) - case ('carecovery') - constitutive_dislotwin_Carecovery(i) = IO_floatValue(line,positions,2) - case ('cstorage') - constitutive_dislotwin_Cstorage(i) = IO_floatValue(line,positions,2) - case ('interaction_slipslip') + case ('cedgedipmindistance') + constitutive_dislotwin_CEdgeDipMinDistance(i) = IO_floatValue(line,positions,2) + case ('catomicvolume') + constitutive_dislotwin_CAtomicVolume(i) = IO_floatValue(line,positions,2) + case ('interactionslipslip') forall (j = 1:lattice_maxNinteraction) & - constitutive_dislotwin_interaction_slipslip(j,i) = IO_floatValue(line,positions,1+j) - case ('interaction_sliptwin') + constitutive_dislotwin_interactionSlipSlip(j,i) = IO_floatValue(line,positions,1+j) + case ('interactionsliptwin') forall (j = 1:lattice_maxNinteraction) & - constitutive_dislotwin_interaction_sliptwin(j,i) = IO_floatValue(line,positions,1+j) - case ('interaction_twinslip') + constitutive_dislotwin_interactionSlipTwin(j,i) = IO_floatValue(line,positions,1+j) + case ('interactiontwinslip') forall (j = 1:lattice_maxNinteraction) & - constitutive_dislotwin_interaction_twinslip(j,i) = IO_floatValue(line,positions,1+j) - case ('interaction_twintwin') + constitutive_dislotwin_interactionTwinSlip(j,i) = IO_floatValue(line,positions,1+j) + case ('interactiontwintwin') forall (j = 1:lattice_maxNinteraction) & - constitutive_dislotwin_interaction_twintwin(j,i) = IO_floatValue(line,positions,1+j) + constitutive_dislotwin_interactionTwinTwin(j,i) = IO_floatValue(line,positions,1+j) end select endif - enddo +enddo 100 do i = 1,maxNinstance - constitutive_dislotwin_structure(i) = lattice_initializeStructure(constitutive_dislotwin_structureName(i), & - constitutive_dislotwin_CoverA(i)) - constitutive_dislotwin_Nslip(:,i) = min(lattice_NslipSystem(:,constitutive_dislotwin_structure(i)),& - constitutive_dislotwin_Nslip(:,i)) - constitutive_dislotwin_Ntwin(:,i) = min(lattice_NtwinSystem(:,constitutive_dislotwin_structure(i)),& - constitutive_dislotwin_Ntwin(:,i)) + constitutive_dislotwin_structure(i) = & + lattice_initializeStructure(constitutive_dislotwin_structureName(i),constitutive_dislotwin_CoverA(i)) + myStructure = constitutive_dislotwin_structure(i) + + !* Sanity checks + if (myStructure < 1 .or. myStructure > 3) call IO_error(205) + if (sum(constitutive_dislotwin_Nslip(:,i)) <= 0_pInt) call IO_error(225) + if (sum(constitutive_dislotwin_Ntwin(:,i)) < 0_pInt) call IO_error(225) !*** + do f = 1,lattice_maxNslipFamily + if (constitutive_dislotwin_Nslip(f,i) > 0_pInt) then + if (constitutive_dislotwin_rhoEdge0(f,i) < 0.0_pReal) call IO_error(220) + if (constitutive_dislotwin_rhoEdgeDip0(f,i) < 0.0_pReal) call IO_error(220) + if (constitutive_dislotwin_burgersPerSlipFamily(f,i) <= 0.0_pReal) call IO_error(221) + if (constitutive_dislotwin_v0PerSlipFamily(f,i) <= 0.0_pReal) call IO_error(226) + endif + enddo + do f = 1,lattice_maxNtwinFamily + if (constitutive_dislotwin_Nslip(f,i) > 0_pInt) then + if (constitutive_dislotwin_burgersPerTwinFamily(f,i) <= 0.0_pReal) call IO_error(221) !*** + if (constitutive_dislotwin_Ndot0PerTwinFamily(f,i) < 0.0_pReal) call IO_error(226) !*** + endif + enddo +! if (any(constitutive_dislotwin_interactionSlipSlip(1:maxval(lattice_interactionSlipSlip(:,:,myStructure)),i) < 1.0_pReal)) call IO_error(229) + if (constitutive_dislotwin_CAtomicVolume(i) <= 0.0_pReal) call IO_error(230) + if (constitutive_dislotwin_D0(i) <= 0.0_pReal) call IO_error(231) + if (constitutive_dislotwin_Qsd(i) <= 0.0_pReal) call IO_error(232) + if (constitutive_dislotwin_relevantRho(i) <= 0.0_pReal) call IO_error(233) + + !* Determine total number of active slip or twin systems + constitutive_dislotwin_Nslip(:,i) = min(lattice_NslipSystem(:,myStructure),constitutive_dislotwin_Nslip(:,i)) + constitutive_dislotwin_Ntwin(:,i) = min(lattice_NtwinSystem(:,myStructure),constitutive_dislotwin_Ntwin(:,i)) constitutive_dislotwin_totalNslip(i) = sum(constitutive_dislotwin_Nslip(:,i)) constitutive_dislotwin_totalNtwin(i) = sum(constitutive_dislotwin_Ntwin(:,i)) + +enddo + +!* Allocation of variables whose size depends on the total number of active slip systems +maxTotalNslip = maxval(constitutive_dislotwin_totalNslip) +maxTotalNtwin = maxval(constitutive_dislotwin_totalNtwin) -! sanity checks (still under construction) - if (constitutive_dislotwin_structure(i) < 1 .or. & ! sanity checks - constitutive_dislotwin_structure(i) > 3) call IO_error(205) - if (any(constitutive_dislotwin_rho0(:,i) < 0.0_pReal)) call IO_error(220) - if (any(constitutive_dislotwin_Burgers(:,i) <= 0.0_pReal .and. & - constitutive_dislotwin_Nslip(:,i) > 0)) call IO_error(221) - if (any(constitutive_dislotwin_Qedge(:,i) <= 0.0_pReal .and. & - constitutive_dislotwin_Nslip(:,i) > 0)) call IO_error(222) - enddo +allocate(constitutive_dislotwin_burgersPerSlipSystem(maxTotalNslip, maxNinstance)) +allocate(constitutive_dislotwin_burgersPerTwinSystem(maxTotalNtwin, maxNinstance)) +allocate(constitutive_dislotwin_QedgePerSlipSystem(maxTotalNslip, maxNinstance)) +allocate(constitutive_dislotwin_v0PerSlipSystem(maxTotalNslip, maxNinstance)) +allocate(constitutive_dislotwin_Ndot0PerTwinSystem(maxTotalNtwin, maxNinstance)) +allocate(constitutive_dislotwin_twinsizePerTwinSystem(maxTotalNtwin, maxNinstance)) +allocate(constitutive_dislotwin_CLambdaSlipPerSlipSystem(maxTotalNslip, maxNinstance)) +constitutive_dislotwin_burgersPerSlipSystem = 0.0_pReal +constitutive_dislotwin_burgersPerTwinSystem = 0.0_pReal +constitutive_dislotwin_QedgePerSlipSystem = 0.0_pReal +constitutive_dislotwin_v0PerSlipSystem = 0.0_pReal +constitutive_dislotwin_Ndot0PerTwinSystem = 0.0_pReal +constitutive_dislotwin_twinsizePerTwinSystem = 0.0_pReal +constitutive_dislotwin_CLambdaSlipPerSlipSystem = 0.0_pReal - allocate(constitutive_dislotwin_parall_interaction(maxval(constitutive_dislotwin_totalNslip),& - maxval(constitutive_dislotwin_totalNslip),& - maxNinstance)) - allocate(constitutive_dislotwin_forest_interaction(maxval(constitutive_dislotwin_totalNslip),& - maxval(constitutive_dislotwin_totalNslip),& - maxNinstance)) - allocate(constitutive_dislotwin_hardeningMatrix_sliptwin(maxval(constitutive_dislotwin_totalNslip),& - maxval(constitutive_dislotwin_totalNtwin),& - maxNinstance)) - allocate(constitutive_dislotwin_hardeningMatrix_twinslip(maxval(constitutive_dislotwin_totalNtwin),& - maxval(constitutive_dislotwin_totalNslip),& - maxNinstance)) - allocate(constitutive_dislotwin_hardeningMatrix_twintwin(maxval(constitutive_dislotwin_totalNtwin),& - maxval(constitutive_dislotwin_totalNtwin),& - maxNinstance)) - constitutive_dislotwin_parall_interaction = 0.0_pReal - constitutive_dislotwin_forest_interaction = 0.0_pReal - constitutive_dislotwin_hardeningMatrix_sliptwin = 0.0_pReal - constitutive_dislotwin_hardeningMatrix_twinslip = 0.0_pReal - constitutive_dislotwin_hardeningMatrix_twintwin = 0.0_pReal +allocate(constitutive_dislotwin_interactionMatrixSlipSlip(maxTotalNslip,maxTotalNslip,maxNinstance)) +allocate(constitutive_dislotwin_interactionMatrixSlipTwin(maxTotalNslip,maxTotalNtwin,maxNinstance)) +allocate(constitutive_dislotwin_interactionMatrixTwinSlip(maxTotalNtwin,maxTotalNslip,maxNinstance)) +allocate(constitutive_dislotwin_interactionMatrixTwinTwin(maxTotalNtwin,maxTotalNtwin,maxNinstance)) +allocate(constitutive_dislotwin_forestProjectionEdge(maxTotalNslip,maxTotalNslip,maxNinstance)) +constitutive_dislotwin_interactionMatrixSlipSlip = 0.0_pReal +constitutive_dislotwin_interactionMatrixSlipTwin = 0.0_pReal +constitutive_dislotwin_interactionMatrixTwinSlip = 0.0_pReal +constitutive_dislotwin_interactionMatrixTwinTwin = 0.0_pReal +constitutive_dislotwin_forestProjectionEdge = 0.0_pReal - allocate(constitutive_dislotwin_Ctwin_66(6,6,maxval(constitutive_dislotwin_totalNtwin),maxNinstance)) - constitutive_dislotwin_Ctwin_66 = 0.0_pReal +allocate(constitutive_dislotwin_Ctwin_66(6,6,maxTotalNtwin,maxNinstance)) +allocate(constitutive_dislotwin_Ctwin_3333(3,3,3,3,maxTotalNtwin,maxNinstance)) +constitutive_dislotwin_Ctwin_66 = 0.0_pReal +constitutive_dislotwin_Ctwin_3333 = 0.0_pReal - allocate(constitutive_dislotwin_Ctwin_3333(3,3,3,3,maxval(constitutive_dislotwin_totalNtwin),maxNinstance)) - constitutive_dislotwin_Ctwin_3333 = 0.0_pReal +do i = 1,maxNinstance + myStructure = constitutive_dislotwin_structure(i) - do i = 1,maxNinstance - do j = 1,maxval(phase_Noutput) - select case(constitutive_dislotwin_output(j,i)) - case('dislocationdensity', & - 'shearrate_slip', & + !* Inverse lookup of my slip system family + l = 0_pInt + do f = 1,lattice_maxNslipFamily + do k = 1,constitutive_dislotwin_Nslip(f,i) + l = l + 1 + constitutive_dislotwin_slipFamily(l,i) = f + constitutive_dislotwin_slipSystemLattice(l,i) = sum(lattice_NslipSystem(1:f-1,myStructure)) + k + enddo; enddo + + !* Inverse lookup of my twin system family + l = 0_pInt + do f = 1,lattice_maxNtwinFamily + do k = 1,constitutive_dislotwin_Ntwin(f,i) + l = l + 1 + constitutive_dislotwin_twinFamily(l,i) = f + constitutive_dislotwin_twinSystemLattice(l,i) = sum(lattice_NtwinSystem(1:f-1,myStructure)) + k + enddo; enddo + + !* Determine size of state array + ns = constitutive_dislotwin_totalNslip(i) + nt = constitutive_dislotwin_totalNtwin(i) + constitutive_dislotwin_sizeDotState(i) = & + size(constitutive_dislotwin_listBasicSlipStates)*ns+size(constitutive_dislotwin_listBasicTwinStates)*nt + constitutive_dislotwin_sizeState(i) = & + constitutive_dislotwin_sizeDotState(i)+ & + size(constitutive_dislotwin_listDependentSlipStates)*ns+size(constitutive_dislotwin_listDependentTwinStates)*nt + + !* Determine size of postResults array + do o = 1,maxval(phase_Noutput) + select case(constitutive_dislotwin_output(o,i)) + case('edge_density', & + 'dipole_density', & + 'shear_rate_slip', & 'mfp_slip', & - 'resolvedstress_slip', & - 'resistance_slip' & + 'resolved_stress_slip', & + 'threshold_stress_slip' & ) mySize = constitutive_dislotwin_totalNslip(i) - case('volumefraction', & - 'shearrate_twin', & + case('twin_fraction', & + 'shear_rate_twin', & 'mfp_twin', & - 'resolvedstress_twin', & - 'resistance_twin' & + 'resolved_stress_twin', & + 'threshold_stress_twin' & ) mySize = constitutive_dislotwin_totalNtwin(i) case default mySize = 0_pInt - end select + end select - if (mySize > 0_pInt) then ! any meaningful output found - constitutive_dislotwin_sizePostResult(j,i) = mySize - constitutive_dislotwin_sizePostResults(i) = constitutive_dislotwin_sizePostResults(i) + mySize - endif + if (mySize > 0_pInt) then ! any meaningful output found + constitutive_dislotwin_sizePostResult(o,i) = mySize + constitutive_dislotwin_sizePostResults(i) = constitutive_dislotwin_sizePostResults(i) + mySize + endif enddo - constitutive_dislotwin_sizeDotState(i) = constitutive_dislotwin_totalNslip(i) + constitutive_dislotwin_totalNtwin(i) - constitutive_dislotwin_sizeState(i) = 10*constitutive_dislotwin_totalNslip(i) + 5*constitutive_dislotwin_totalNtwin(i) - - constitutive_dislotwin_Gmod(i) = constitutive_dislotwin_C44(i) - - select case (constitutive_dislotwin_structure(i)) + !* Elasticity matrix and shear modulus according to material.config + select case (myStructure) case(1:2) ! cubic(s) forall(k=1:3) forall(j=1:3) & @@ -410,27 +485,12 @@ subroutine constitutive_dislotwin_init(file) constitutive_dislotwin_Cslip_66(3,2,i) = constitutive_dislotwin_C13(i) constitutive_dislotwin_Cslip_66(4,4,i) = constitutive_dislotwin_C44(i) constitutive_dislotwin_Cslip_66(5,5,i) = constitutive_dislotwin_C44(i) - constitutive_dislotwin_Cslip_66(6,6,i) = 0.5_pReal*(constitutive_dislotwin_C11(i)- & - constitutive_dislotwin_C12(i)) + constitutive_dislotwin_Cslip_66(6,6,i) = 0.5_pReal*(constitutive_dislotwin_C11(i)-constitutive_dislotwin_C12(i)) end select constitutive_dislotwin_Cslip_66(:,:,i) = math_Mandel3333to66(math_Voigt66to3333(constitutive_dislotwin_Cslip_66(:,:,i))) constitutive_dislotwin_Cslip_3333(:,:,:,:,i) = math_Voigt66to3333(constitutive_dislotwin_Cslip_66(:,:,i)) - - !* Inverse lookup of my slip system family - l = 0_pInt - do j = 1,lattice_maxNslipFamily - do k = 1,constitutive_dislotwin_Nslip(j,i) - l = l + 1 - constitutive_dislotwin_slipFamily(l,i) = j - enddo; enddo - - !* Inverse lookup of my twin system family - l = 0_pInt - do j = 1,lattice_maxNtwinFamily - do k = 1,constitutive_dislotwin_Ntwin(j,i) - l = l + 1 - constitutive_dislotwin_twinFamily(l,i) = j - enddo; enddo + constitutive_dislotwin_Gmod(i) = & + 0.2_pReal*(constitutive_dislotwin_C11(i)-constitutive_dislotwin_C12(i))+0.3_pReal*constitutive_dislotwin_C44(i) !* Construction of the twin elasticity matrices do j=1,lattice_maxNtwinFamily @@ -439,82 +499,76 @@ subroutine constitutive_dislotwin_init(file) constitutive_dislotwin_Ctwin_3333(l,m,n,o,sum(constitutive_dislotwin_Nslip(1:j-1,i))+k,i) = & constitutive_dislotwin_Ctwin_3333(l,m,n,o,sum(constitutive_dislotwin_Nslip(1:j-1,i))+k,i) + & constitutive_dislotwin_Cslip_3333(p,q,r,s,i)*& - lattice_Qtwin(l,p,sum(lattice_NslipSystem(1:j-1,constitutive_dislotwin_structure(i)))+k,constitutive_dislotwin_structure(i))* & - lattice_Qtwin(m,q,sum(lattice_NslipSystem(1:j-1,constitutive_dislotwin_structure(i)))+k,constitutive_dislotwin_structure(i))* & - lattice_Qtwin(n,r,sum(lattice_NslipSystem(1:j-1,constitutive_dislotwin_structure(i)))+k,constitutive_dislotwin_structure(i))* & - lattice_Qtwin(o,s,sum(lattice_NslipSystem(1:j-1,constitutive_dislotwin_structure(i)))+k,constitutive_dislotwin_structure(i)) + lattice_Qtwin(l,p,sum(lattice_NslipSystem(1:j-1,myStructure))+k,myStructure)* & + lattice_Qtwin(m,q,sum(lattice_NslipSystem(1:j-1,myStructure))+k,myStructure)* & + lattice_Qtwin(n,r,sum(lattice_NslipSystem(1:j-1,myStructure))+k,myStructure)* & + lattice_Qtwin(o,s,sum(lattice_NslipSystem(1:j-1,myStructure))+k,myStructure) enddo ; enddo ; enddo ; enddo ; enddo ; enddo ; enddo ; enddo constitutive_dislotwin_Ctwin_66(:,:,k,i) = math_Mandel3333to66(constitutive_dislotwin_Ctwin_3333(:,:,:,:,k,i)) enddo enddo - !* Construction of the hardening matrices - !* Iteration over the systems - do j=1,lattice_maxNslipFamily - do k=1,constitutive_dislotwin_Nslip(j,i) - do l=1,lattice_maxNslipFamily - do m=1,constitutive_dislotwin_Nslip(l,i) - !* Projection of the dislocation * - x = math_mul3x3(lattice_sn(:,sum(lattice_NslipSystem(1:j-1,constitutive_dislotwin_structure(i)))+k, & - constitutive_dislotwin_structure(i)), & - lattice_st(:,sum(lattice_NslipSystem(1:l-1,constitutive_dislotwin_structure(i)))+m, & - constitutive_dislotwin_structure(i))) - y = 1.0_pReal-x**(2.0_pReal) - !* Interaction matrix * - constitutive_dislotwin_forest_interaction(sum(constitutive_dislotwin_Nslip(1:j-1,i))+k, & - sum(constitutive_dislotwin_Nslip(1:l-1,i))+m,i) = & - abs(x)*constitutive_dislotwin_interaction_slipslip(lattice_interactionSlipSlip( & - sum(lattice_NslipSystem(1:j-1,constitutive_dislotwin_structure(i)))+k, & - sum(lattice_NslipSystem(1:l-1,constitutive_dislotwin_structure(i)))+m, & - constitutive_dislotwin_structure(i)),i) - if (y>0.0_pReal) & - constitutive_dislotwin_parall_interaction(sum(constitutive_dislotwin_Nslip(1:j-1,i))+k, & - sum(constitutive_dislotwin_Nslip(1:l-1,i))+m,i) = & - sqrt(y)*constitutive_dislotwin_interaction_slipslip(lattice_interactionSlipSlip( & - sum(lattice_NslipSystem(1:j-1,constitutive_dislotwin_structure(i)))+k, & - sum(lattice_NslipSystem(1:l-1,constitutive_dislotwin_structure(i)))+m, & - constitutive_dislotwin_structure(i)),i) - enddo; enddo; enddo; enddo + !* Burgers vector, dislocation velocity prefactor, mean free path prefactor and minimum dipole distance for each slip system + do s = 1,constitutive_dislotwin_totalNslip(i) + f = constitutive_dislotwin_slipFamily(s,i) + constitutive_dislotwin_burgersPerSlipSystem(s,i) = constitutive_dislotwin_burgersPerSlipFamily(f,i) + constitutive_dislotwin_QedgePerSlipSystem(s,i) = constitutive_dislotwin_QedgePerSlipFamily(f,i) + constitutive_dislotwin_v0PerSlipSystem(s,i) = constitutive_dislotwin_v0PerSlipFamily(f,i) + constitutive_dislotwin_CLambdaSlipPerSlipSystem(s,i) = constitutive_dislotwin_CLambdaSlipPerSlipFamily(f,i) + enddo + + !* Burgers vector, nucleation rate prefactor and twin size for each twin system + do s = 1,constitutive_dislotwin_totalNtwin(i) + f = constitutive_dislotwin_twinFamily(s,i) + constitutive_dislotwin_burgersPerTwinSystem(s,i) = constitutive_dislotwin_burgersPerTwinFamily(f,i) + constitutive_dislotwin_Ndot0PerTwinSystem(s,i) = constitutive_dislotwin_Ndot0PerTwinFamily(f,i) + constitutive_dislotwin_twinsizePerTwinSystem(s,i) = constitutive_dislotwin_twinsizePerTwinFamily(f,i) + enddo + + !* Construction of interaction matrices + do s1 = 1,constitutive_dislotwin_totalNslip(i) + do s2 = 1,constitutive_dislotwin_totalNslip(i) + constitutive_dislotwin_interactionMatrixSlipSlip(s1,s2,i) = & + constitutive_dislotwin_interactionSlipSlip(lattice_interactionSlipSlip(constitutive_dislotwin_slipSystemLattice(s1,i), & + constitutive_dislotwin_slipSystemLattice(s2,i), & + myStructure),i) + enddo; enddo + + do s1 = 1,constitutive_dislotwin_totalNslip(i) + do t2 = 1,constitutive_dislotwin_totalNtwin(i) + constitutive_dislotwin_interactionMatrixSlipTwin(s1,t2,i) = & + constitutive_dislotwin_interactionSlipTwin(lattice_interactionSlipTwin(constitutive_dislotwin_slipSystemLattice(s1,i), & + constitutive_dislotwin_twinSystemLattice(t2,i), & + myStructure),i) + enddo; enddo + + do t1 = 1,constitutive_dislotwin_totalNtwin(i) + do s2 = 1,constitutive_dislotwin_totalNslip(i) + constitutive_dislotwin_interactionMatrixTwinSlip(t1,s2,i) = & + constitutive_dislotwin_interactionTwinSlip(lattice_interactionTwinSlip(constitutive_dislotwin_twinSystemLattice(t1,i), & + constitutive_dislotwin_slipSystemLattice(s2,i), & + myStructure),i) + enddo; enddo - do j=1,lattice_maxNslipFamily - do k=1,constitutive_dislotwin_Nslip(j,i) - do l=1,lattice_maxNtwinFamily - do m=1,constitutive_dislotwin_Ntwin(l,i) - constitutive_dislotwin_hardeningMatrix_sliptwin(sum(constitutive_dislotwin_Nslip(1:j-1,i))+k,& - sum(constitutive_dislotwin_Ntwin(1:l-1,i))+m,i) = & - constitutive_dislotwin_interaction_sliptwin(lattice_interactionSlipTwin( & - sum(lattice_NslipSystem(1:j-1,constitutive_dislotwin_structure(i)))+k, & - sum(lattice_NtwinSystem(1:l-1,constitutive_dislotwin_structure(i)))+m, & - constitutive_dislotwin_structure(i)),i) - enddo; enddo; enddo; enddo + do t1 = 1,constitutive_dislotwin_totalNtwin(i) + do t2 = 1,constitutive_dislotwin_totalNtwin(i) + constitutive_dislotwin_interactionMatrixTwinTwin(t1,t2,i) = & + constitutive_dislotwin_interactionTwinTwin(lattice_interactionTwinTwin(constitutive_dislotwin_twinSystemLattice(t1,i), & + constitutive_dislotwin_twinSystemLattice(t2,i), & + myStructure),i) + enddo; enddo + + !* Calculation of forest projections for edge dislocations + do s1 = 1,constitutive_dislotwin_totalNslip(i) + do s2 = 1,constitutive_dislotwin_totalNslip(i) + constitutive_dislotwin_forestProjectionEdge(s1,s2,i) = & + abs(math_mul3x3(lattice_sn(:,constitutive_dislotwin_slipSystemLattice(s1,i),myStructure), & + lattice_st(:,constitutive_dislotwin_slipSystemLattice(s2,i),myStructure))) + enddo; enddo + +enddo - do j=1,lattice_maxNtwinFamily - do k=1,constitutive_dislotwin_Ntwin(j,i) - do l=1,lattice_maxNslipFamily - do m=1,constitutive_dislotwin_Nslip(l,i) - constitutive_dislotwin_hardeningMatrix_twinslip(sum(constitutive_dislotwin_Ntwin(1:j-1,i))+k,& - sum(constitutive_dislotwin_Nslip(1:l-1,i))+m,i) = & - constitutive_dislotwin_interaction_twinslip(lattice_interactionTwinSlip( & - sum(lattice_NtwinSystem(1:j-1,constitutive_dislotwin_structure(i)))+k, & - sum(lattice_NslipSystem(1:l-1,constitutive_dislotwin_structure(i)))+m, & - constitutive_dislotwin_structure(i)),i) - enddo; enddo; enddo; enddo - - do j=1,lattice_maxNtwinFamily - do k=1,constitutive_dislotwin_Ntwin(j,i) - do l=1,lattice_maxNtwinFamily - do m=1,constitutive_dislotwin_Ntwin(l,i) - constitutive_dislotwin_hardeningMatrix_twintwin(sum(constitutive_dislotwin_Ntwin(1:j-1,i))+k,& - sum(constitutive_dislotwin_Ntwin(1:l-1,i))+m,i) = & - constitutive_dislotwin_interaction_twintwin(lattice_interactionTwinTwin( & - sum(lattice_NtwinSystem(1:j-1,constitutive_dislotwin_structure(i)))+k, & - sum(lattice_NtwinSystem(1:l-1,constitutive_dislotwin_structure(i)))+m, & - constitutive_dislotwin_structure(i)), i ) - enddo; enddo; enddo; enddo - - enddo - - return +return end subroutine @@ -522,90 +576,135 @@ function constitutive_dislotwin_stateInit(myInstance) !********************************************************************* !* initial microstructural state * !********************************************************************* - use prec, only: pReal,pInt - use lattice, only: lattice_maxNslipFamily,lattice_maxNtwinFamily - implicit none +use prec, only: pReal,pInt +use math, only: pi +use lattice, only: lattice_maxNslipFamily,lattice_maxNtwinFamily +implicit none - !* Definition of variables - integer(pInt), intent(in) :: myInstance - integer(pInt) i - real(pReal), dimension(constitutive_dislotwin_sizeState(myInstance)) :: constitutive_dislotwin_stateInit - - constitutive_dislotwin_stateInit = 0.0_pReal - - do i = 1,lattice_maxNslipFamily - constitutive_dislotwin_stateInit(1+sum(constitutive_dislotwin_Nslip(1:i-1,myInstance)) : & - sum(constitutive_dislotwin_Nslip(1:i ,myInstance))) = & - constitutive_dislotwin_rho0(i,myInstance) - enddo - - return +!* Input-Output variables +integer(pInt) :: myInstance +real(pReal), dimension(constitutive_dislotwin_sizeState(myInstance)) :: constitutive_dislotwin_stateInit +!* Local variables +integer(pInt) s0,s1,s,t,f,ns,nt +real(pReal), dimension(constitutive_dislotwin_totalNslip(myInstance)) :: rhoEdge0,rhoEdgeDip0,invLambdaSlip0,MeanFreePathSlip0,tauSlipThreshold0 +real(pReal), dimension(constitutive_dislotwin_totalNtwin(myInstance)) :: MeanFreePathTwin0,TwinVolume0 + +ns = constitutive_dislotwin_totalNslip(myInstance) +nt = constitutive_dislotwin_totalNtwin(myInstance) +constitutive_dislotwin_stateInit = 0.0_pReal + +!* Initialize basic slip state variables +s1 = 0_pInt +do f = 1,lattice_maxNslipFamily + s0 = s1 + 1_pInt + s1 = s0 + constitutive_dislotwin_Nslip(f,myInstance) - 1_pInt + do s = s0,s1 + rhoEdge0(s) = constitutive_dislotwin_rhoEdge0(f,myInstance) + rhoEdgeDip0(s) = constitutive_dislotwin_rhoEdgeDip0(f,myInstance) + enddo +enddo +constitutive_dislotwin_stateInit(1:ns) = rhoEdge0 +constitutive_dislotwin_stateInit(ns+1:2*ns) = rhoEdgeDip0 + +!* Initialize dependent slip microstructural variables +forall (s = 1:ns) & +invLambdaSlip0(s) = sqrt(dot_product((rhoEdge0+rhoEdgeDip0),constitutive_dislotwin_forestProjectionEdge(1:ns,s,myInstance)))/ & +constitutive_dislotwin_CLambdaSlipPerSlipSystem(s,myInstance) +constitutive_dislotwin_stateInit(2*ns+nt+1:3*ns+nt) = invLambdaSlip0 + +forall (s = 1:ns) & +MeanFreePathSlip0(s) = & +constitutive_dislotwin_GrainSize(myInstance)/(1.0_pReal+invLambdaSlip0(s)*constitutive_dislotwin_GrainSize(myInstance)) +constitutive_dislotwin_stateInit(4*ns+2*nt+1:5*ns+2*nt) = MeanFreePathSlip0 + +forall (s = 1:ns) & +tauSlipThreshold0(s) = & +constitutive_dislotwin_Gmod(myInstance)*constitutive_dislotwin_burgersPerSlipSystem(s,myInstance)* & +sqrt(dot_product(rhoEdge0,constitutive_dislotwin_interactionMatrixSlipSlip(1:ns,s,myInstance))) +constitutive_dislotwin_stateInit(5*ns+3*nt+1:6*ns+3*nt) = tauSlipThreshold0 + +!* Initialize dependent twin microstructural variables +forall (t = 1:nt) & +MeanFreePathTwin0(t) = constitutive_dislotwin_GrainSize(myInstance) +constitutive_dislotwin_stateInit(5*ns+2*nt+1:5*ns+3*nt) = MeanFreePathTwin0 + +forall (t = 1:nt) & +TwinVolume0(t) = & +(pi/6.0_pReal)*constitutive_dislotwin_twinsizePerTwinSystem(t,myInstance)*MeanFreePathTwin0(t)**(2.0_pReal) +constitutive_dislotwin_stateInit(6*ns+4*nt+1:6*ns+5*nt) = TwinVolume0 + +!write(6,*) '#STATEINIT#' +!write(6,*) +!write(6,'(a,/,4(3(f30.20,x)/))') 'RhoEdge',rhoEdge0 +!write(6,'(a,/,4(3(f30.20,x)/))') 'RhoEdgedip',rhoEdgeDip0 +!write(6,'(a,/,4(3(f30.20,x)/))') 'invLambdaSlip',invLambdaSlip0 +!write(6,'(a,/,4(3(f30.20,x)/))') 'MeanFreePathSlip',MeanFreePathSlip0 +!write(6,'(a,/,4(3(f30.20,x)/))') 'tauSlipThreshold', tauSlipThreshold0 +!write(6,'(a,/,4(3(f30.20,x)/))') 'MeanFreePathTwin', MeanFreePathTwin0 +!write(6,'(a,/,4(3(f30.20,x)/))') 'TwinVolume', TwinVolume0 + +return end function +pure function constitutive_dislotwin_relevantState(myInstance) !********************************************************************* !* relevant microstructural state * !********************************************************************* -pure function constitutive_dislotwin_relevantState(myInstance) - -use prec, only: pReal, & - pInt +use prec, only: pReal, pInt implicit none -!*** input variables -integer(pInt), intent(in) :: myInstance ! number specifying the current instance of the constitution +!* Input-Output variables +integer(pInt), intent(in) :: myInstance +real(pReal), dimension(constitutive_dislotwin_sizeState(myInstance)) :: constitutive_dislotwin_relevantState -!*** output variables -real(pReal), dimension(constitutive_dislotwin_sizeState(myInstance)) :: & - constitutive_dislotwin_relevantState ! relevant state values for the current instance of this constitution - -!*** local variables - -constitutive_dislotwin_relevantState = 1.0e-200_pReal +constitutive_dislotwin_relevantState = constitutive_dislotwin_relevantRho(myInstance) +return endfunction -function constitutive_dislotwin_homogenizedC(state,ipc,ip,el) +pure function constitutive_dislotwin_homogenizedC(state,g,ip,el) !********************************************************************* !* calculates homogenized elacticity matrix * !* - state : microstructure quantities * -!* - ipc : component-ID of current integration point * +!* - g : component-ID of current integration point * !* - ip : current integration point * !* - el : current element * !********************************************************************* - use prec, only: pReal,pInt,p_vec - use mesh, only: mesh_NcpElems,mesh_maxNips - use material, only: homogenization_maxNgrains,material_phase,phase_constitutionInstance - implicit none +use prec, only: pReal,pInt,p_vec +use mesh, only: mesh_NcpElems,mesh_maxNips +use material, only: homogenization_maxNgrains,material_phase,phase_constitutionInstance +implicit none - !* Definition of variables - integer(pInt), intent(in) :: ipc,ip,el - integer(pInt) matID,ns,nt,i - real(pReal) sumf - real(pReal), dimension(6,6) :: constitutive_dislotwin_homogenizedC - type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems) :: state +!* Input-Output variables +integer(pInt), intent(in) :: g,ip,el +type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: state +real(pReal), dimension(6,6) :: constitutive_dislotwin_homogenizedC +!* Local variables +integer(pInt) myInstance,ns,nt,i +real(pReal) sumf - !* Shortened notation - matID = phase_constitutionInstance(material_phase(ipc,ip,el)) - ns = constitutive_dislotwin_totalNslip(matID) - nt = constitutive_dislotwin_totalNtwin(matID) +!* Shortened notation +myInstance = phase_constitutionInstance(material_phase(g,ip,el)) +ns = constitutive_dislotwin_totalNslip(myInstance) +nt = constitutive_dislotwin_totalNtwin(myInstance) - !* Total twin volume fraction - sumf = sum(state(ipc,ip,el)%p((ns+1):(ns+nt))) ! safe for nt == 0 +!* Total twin volume fraction +sumf = sum(state(g,ip,el)%p((2*ns+1):(2*ns+nt))) ! safe for nt == 0 - !* Homogenized elasticity matrix - constitutive_dislotwin_homogenizedC = (1.0_pReal-sumf)*constitutive_dislotwin_Cslip_66(:,:,matID) - do i=1,nt - constitutive_dislotwin_homogenizedC = constitutive_dislotwin_homogenizedC + & - state(ipc,ip,el)%p(ns+i)*constitutive_dislotwin_Ctwin_66(:,:,i,matID) - enddo +!* Homogenized elasticity matrix +constitutive_dislotwin_homogenizedC = (1.0_pReal-sumf)*constitutive_dislotwin_Cslip_66(:,:,myInstance) +do i=1,nt + constitutive_dislotwin_homogenizedC = & + constitutive_dislotwin_homogenizedC + state(g,ip,el)%p(2*ns+i)*constitutive_dislotwin_Ctwin_66(:,:,i,myInstance) +enddo - return +return end function -subroutine constitutive_dislotwin_microstructure(Temperature,state,ipc,ip,el) +subroutine constitutive_dislotwin_microstructure(Temperature,state,g,ip,el) !********************************************************************* !* calculates quantities characterizing the microstructure * !* - Temperature : temperature * @@ -614,132 +713,121 @@ subroutine constitutive_dislotwin_microstructure(Temperature,state,ipc,ip,el) !* - ip : current integration point * !* - el : current element * !********************************************************************* - use prec, only: pReal,pInt,p_vec - use math, only: pi - use mesh, only: mesh_NcpElems,mesh_maxNips - use material, only: homogenization_maxNgrains,material_phase,phase_constitutionInstance - use lattice, only: lattice_interactionSlipTwin,lattice_interactionTwinTwin - implicit none +use prec, only: pReal,pInt,p_vec +use math, only: pi +use mesh, only: mesh_NcpElems,mesh_maxNips +use material, only: homogenization_maxNgrains,material_phase,phase_constitutionInstance +use lattice, only: lattice_interactionSlipTwin,lattice_interactionTwinTwin +!use debug, only: debugger +implicit none - !* Definition of variables - integer(pInt), intent(in) :: ipc,ip,el - integer(pInt) matID,ns,nt,i - real(pReal) Temperature,sumf - real(pReal), dimension(constitutive_dislotwin_totalNtwin(phase_constitutionInstance(material_phase(ipc,ip,el)))) :: fOverStacksize - type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems) :: state - - !* Shortened notation - matID = phase_constitutionInstance(material_phase(ipc,ip,el)) - ns = constitutive_dislotwin_totalNslip(matID) - nt = constitutive_dislotwin_totalNtwin(matID) - !* State: 1 : ns rho_ssd - !* State: ns+1 : ns+nt f - !* State: ns+nt+1 : 2*ns+nt rho_forest - !* State: 2*ns+nt+1 : 3*ns+nt rho_parallel - !* State: 3*ns+nt+1 : 4*ns+nt 1/lambda_slip - !* State: 4*ns+nt+1 : 5*ns+nt 1/lambda_sliptwin - !* State: 5*ns+nt+1 : 5*ns+2*nt 1/lambda_twin - !* State: 5*ns+2*nt+1 : 6*ns+2*nt mfp_slip - !* State: 6*ns+2*nt+1 : 6*ns+3*nt mfp_twin - !* State: 6*ns+3*nt+1 : 7*ns+3*nt threshold_stress_slip - !* State: 7*ns+3*nt+1 : 7*ns+4*nt threshold_stress_twin - !* State: 7*ns+4*nt+1 : 8*ns+4*nt activation volume - !* State: 8*ns+4*nt+1 : 8*ns+5*nt twin volume - !* State: 8*ns+5*nt+1 : 9*ns+5*nt rho_mobile - !* State: 9*ns+5*nt+1 : 10*ns+5*nt initial shear rate - - !* Total twin volume fraction - sumf = sum(state(ipc,ip,el)%p((ns+1):(ns+nt))) ! safe for nt == 0 - - !* rescaled twin volume fraction for topology - forall (i = 1:nt) & - fOverStacksize(i) = state(ipc,ip,el)%p(ns+i)/constitutive_dislotwin_stacksize(constitutive_dislotwin_twinFamily(i,matID),matID) - - !* Forest and parallel dislocation densities - !$OMP CRITICAL (evilmatmul) - state(ipc,ip,el)%p((ns+nt+1):(2*ns+nt)) = & - matmul(constitutive_dislotwin_forest_interaction(1:ns,1:ns,matID),state(ipc,ip,el)%p(1:ns)) - state(ipc,ip,el)%p((2*ns+nt+1):(3*ns+nt)) = & - matmul(constitutive_dislotwin_parall_interaction(1:ns,1:ns,matID),state(ipc,ip,el)%p(1:ns)) - !$OMP END CRITICAL (evilmatmul) +!* Input-Output variables +integer(pInt), intent(in) :: g,ip,el +real(pReal), intent(in) :: Temperature +type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(inout) :: state +!* Local variables +integer(pInt) myInstance,myStructure,ns,nt,s,t,i +real(pReal) sumf,sfe +real(pReal), dimension(constitutive_dislotwin_totalNtwin(phase_constitutionInstance(material_phase(g,ip,el)))) :: fOverStacksize - !* 1/mean free distance between 2 forest dislocations seen by a moving dislocation - forall (i=1:ns) state(ipc,ip,el)%p(3*ns+nt+i) = sqrt(state(ipc,ip,el)%p(ns+nt+i)) +!* Shortened notation +myInstance = phase_constitutionInstance(material_phase(g,ip,el)) +myStructure = constitutive_dislotwin_structure(myInstance) +ns = constitutive_dislotwin_totalNslip(myInstance) +nt = constitutive_dislotwin_totalNtwin(myInstance) +!* State: 1 : ns rho_edge +!* State: ns+1 : 2*ns rho_dipole +!* State: 2*ns+1 : 2*ns+nt f +!* State: 2*ns+nt+1 : 3*ns+nt 1/lambda_slip +!* State: 3*ns+nt+1 : 4*ns+nt 1/lambda_sliptwin +!* State: 4*ns+nt+1 : 4*ns+2*nt 1/lambda_twin +!* State: 4*ns+2*nt+1 : 5*ns+2*nt mfp_slip +!* State: 5*ns+2*nt+1 : 5*ns+3*nt mfp_twin +!* State: 5*ns+3*nt+1 : 6*ns+3*nt threshold_stress_slip +!* State: 6*ns+3*nt+1 : 6*ns+4*nt threshold_stress_twin +!* State: 6*ns+4*nt+1 : 6*ns+5*nt twin volume - !* 1/mean free distance between 2 twin stacks from different systems seen by a moving dislocation - !$OMP CRITICAL (evilmatmul) - state(ipc,ip,el)%p((4*ns+nt+1):(5*ns+nt)) = 0.0_pReal - if (nt > 0_pInt) state(ipc,ip,el)%p((4*ns+nt+1):(5*ns+nt)) = & - matmul(constitutive_dislotwin_hardeningMatrix_sliptwin(1:ns,1:nt,matID),fOverStacksize(1:nt))/& - (2.0_pReal*(1.0_pReal-sumf)) - !$OMP END CRITICAL (evilmatmul) +!* Total twin volume fraction +sumf = sum(state(g,ip,el)%p((2*ns+1):(2*ns+nt))) ! safe for nt == 0 + +!* Stacking fault energy +sfe = 0.0002_pReal*Temperature-0.0396_pReal - !* 1/mean free distance between 2 twin stacks from different systems seen by a growing twin - !$OMP CRITICAL (evilmatmul) - if (nt > 0_pInt) state(ipc,ip,el)%p((5*ns+nt+1):(5*ns+2*nt)) = & - matmul(constitutive_dislotwin_hardeningMatrix_twintwin(1:nt,1:nt,matID),fOverStacksize(1:nt))/& - (2.0_pReal*(1.0_pReal-sumf)) - !$OMP END CRITICAL (evilmatmul) - - !* mean free path between 2 obstacles seen by a moving dislocation - do i=1,ns - if (nt > 0_pInt) then - state(ipc,ip,el)%p(5*ns+2*nt+i) = (constitutive_dislotwin_Cmfpslip(matID)*constitutive_dislotwin_grainsize(matID))/& - (1.0_pReal+constitutive_dislotwin_grainsize(matID)*& - (state(ipc,ip,el)%p(3*ns+nt+i)+state(ipc,ip,el)%p(4*ns+nt+i))) - else - state(ipc,ip,el)%p(5*ns+i) = (constitutive_dislotwin_Cmfpslip(matID)*constitutive_dislotwin_grainsize(matID))/& - (1.0_pReal+constitutive_dislotwin_grainsize(matID)*(state(ipc,ip,el)%p(3*ns+i))) - endif - enddo - - !* mean free path between 2 obstacles seen by a growing twin - forall (i=1:nt) & - state(ipc,ip,el)%p(6*ns+2*nt+i) = (constitutive_dislotwin_Cmfptwin(matID)*constitutive_dislotwin_grainsize(matID))/& - (1.0_pReal+constitutive_dislotwin_grainsize(matID)*state(ipc,ip,el)%p(5*ns+nt+i)) - - !* threshold stress for dislocation motion - forall (i=1:ns) & - state(ipc,ip,el)%p(6*ns+3*nt+i) = constitutive_dislotwin_Cthresholdslip(matID)*& - constitutive_dislotwin_Burgers(constitutive_dislotwin_slipFamily(i,matID),matID)*& - constitutive_dislotwin_Gmod(matID)*sqrt(state(ipc,ip,el)%p(2*ns+nt+i)) - - !* threshold stress for growing twin - forall (i=1:nt) & - state(ipc,ip,el)%p(7*ns+3*nt+i) = constitutive_dislotwin_Cthresholdtwin(matID)*(sqrt(3.0_pReal)/3.0_pReal)*(& - (0.0002_pReal*Temperature-0.0396_pReal)/constitutive_dislotwin_Burgers(constitutive_dislotwin_slipFamily(i,matID),matID)+& - (constitutive_dislotwin_Burgers(constitutive_dislotwin_slipFamily(i,matID),matID)*& - constitutive_dislotwin_Gmod(matID))/state(ipc,ip,el)%p(5*ns+2*nt+i)) - - !* activation volume for dislocation glide - forall (i=1:ns) & - state(ipc,ip,el)%p(7*ns+4*nt+i) = constitutive_dislotwin_Cactivolume(matID)*& - constitutive_dislotwin_Burgers(constitutive_dislotwin_slipFamily(i,matID),matID)**2*state(ipc,ip,el)%p(5*ns+2*nt+i) - - !* final twin volume after growth - forall (i=1:nt) & - state(ipc,ip,el)%p(8*ns+4*nt+i) = (pi/6.0_pReal)*& - constitutive_dislotwin_stacksize(constitutive_dislotwin_twinFamily(i,matID),matID)*& - state(ipc,ip,el)%p(6*ns+2*nt+i)*state(ipc,ip,el)%p(6*ns+2*nt+i) +!* rescaled twin volume fraction for topology +forall (t = 1:nt) & +fOverStacksize(t) = & +state(g,ip,el)%p(2*ns+t)/constitutive_dislotwin_twinsizePerTwinSystem(t,myInstance) - !* mobile dislocation densities - forall (i=1:ns) & - state(ipc,ip,el)%p(8*ns+5*nt+i) = (2.0_pReal*kB*Temperature*state(ipc,ip,el)%p(2*ns+nt+i))/& - (state(ipc,ip,el)%p(6*ns+3*nt+i)*state(ipc,ip,el)%p(7*ns+4*nt+i)) +!* 1/mean free distance between 2 forest dislocations seen by a moving dislocation +forall (s = 1:ns) & +state(g,ip,el)%p(2*ns+nt+s) = & +sqrt(dot_product((state(g,ip,el)%p(1:ns)+state(g,ip,el)%p(ns+1:2*ns)),constitutive_dislotwin_forestProjectionEdge(1:ns,s,myInstance)))/ & +constitutive_dislotwin_CLambdaSlipPerSlipSystem(s,myInstance) - !* initial shear rate for slip - forall (i=1:ns) & - state(ipc,ip,el)%p(9*ns+5*nt+i) = state(ipc,ip,el)%p(8*ns+5*nt+i)*& - constitutive_dislotwin_Burgers(constitutive_dislotwin_slipFamily(i,matID),matID)*& - attack_frequency*state(ipc,ip,el)%p(5*ns+2*nt+i)*& - exp(-constitutive_dislotwin_Qedge(constitutive_dislotwin_slipFamily(i,matID),matID)/& - ! -------------------- - (kB*Temperature)) - +!* 1/mean free distance between 2 twin stacks from different systems seen by a moving dislocation +!$OMP CRITICAL (evilmatmul) +state(g,ip,el)%p((3*ns+nt+1):(4*ns+nt)) = 0.0_pReal +if (nt > 0_pInt) state(g,ip,el)%p((3*ns+nt+1):(4*ns+nt)) = & +matmul(constitutive_dislotwin_interactionMatrixSlipTwin(1:ns,1:nt,myInstance),fOverStacksize(1:nt))/(1.0_pReal-sumf) +!$OMP END CRITICAL (evilmatmul) + +!* 1/mean free distance between 2 twin stacks from different systems seen by a growing twin +!$OMP CRITICAL (evilmatmul) +if (nt > 0_pInt) state(g,ip,el)%p((4*ns+nt+1):(4*ns+2*nt)) = & +matmul(constitutive_dislotwin_interactionMatrixTwinTwin(1:nt,1:nt,myInstance),fOverStacksize(1:nt))/(1.0_pReal-sumf) +!$OMP END CRITICAL (evilmatmul) + +!* mean free path between 2 obstacles seen by a moving dislocation +do s = 1,ns + if (nt > 0_pInt) then + state(g,ip,el)%p(4*ns+2*nt+s) = & + constitutive_dislotwin_GrainSize(myInstance)/(1.0_pReal+constitutive_dislotwin_GrainSize(myInstance)*& + (state(g,ip,el)%p(2*ns+nt+s)+state(g,ip,el)%p(3*ns+nt+s))) + else + state(g,ip,el)%p(4*ns+s) = & + constitutive_dislotwin_GrainSize(myInstance)/& + (1.0_pReal+constitutive_dislotwin_GrainSize(myInstance)*(state(g,ip,el)%p(2*ns+s))) + endif +enddo + +!* mean free path between 2 obstacles seen by a growing twin +forall (t = 1:nt) & +state(g,ip,el)%p(5*ns+2*nt+t) = & +(constitutive_dislotwin_Cmfptwin(myInstance)*constitutive_dislotwin_GrainSize(myInstance))/& +(1.0_pReal+constitutive_dislotwin_GrainSize(myInstance)*state(g,ip,el)%p(4*ns+nt+t)) + +!* threshold stress for dislocation motion +forall (s = 1:ns) & +state(g,ip,el)%p(5*ns+3*nt+s) = & +constitutive_dislotwin_Gmod(myInstance)*constitutive_dislotwin_burgersPerSlipSystem(s,myInstance)*& +sqrt(dot_product(state(g,ip,el)%p(1:ns),constitutive_dislotwin_interactionMatrixSlipSlip(1:ns,s,myInstance))) + +!* threshold stress for growing twin +forall (t = 1:nt) & +state(g,ip,el)%p(6*ns+3*nt+t) = & +constitutive_dislotwin_Cthresholdtwin(myInstance)*& +(sfe/(3.0_pReal*constitutive_dislotwin_burgersPerTwinSystem(t,myInstance))+& +(constitutive_dislotwin_burgersPerTwinSystem(t,myInstance)*constitutive_dislotwin_Gmod(myInstance))/& +state(g,ip,el)%p(5*ns+2*nt+t)) + +!* final twin volume after growth +forall (t = 1:nt) & +state(g,ip,el)%p(6*ns+4*nt+t) = & +(pi/6.0_pReal)*constitutive_dislotwin_twinsizePerTwinSystem(t,myInstance)*state(g,ip,el)%p(5*ns+2*nt+t)**(2.0_pReal) + +!if ((ip==1).and.(el==1)) then +! write(6,*) '#MICROSTRUCTURE#' +! write(6,*) +! write(6,'(a,/,4(3(f10.4,x)/))') 'rhoEdge',state(g,ip,el)%p(1:ns)/1e9 +! write(6,'(a,/,4(3(f10.4,x)/))') 'rhoEdgeDip',state(g,ip,el)%p(ns+1:2*ns)/1e9 +! write(6,'(a,/,4(3(f10.4,x)/))') 'Fraction',state(g,ip,el)%p(2*ns+1:2*ns+nt) +!endif + +return end subroutine -subroutine constitutive_dislotwin_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v,Temperature,state,ipc,ip,el) +subroutine constitutive_dislotwin_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v,Temperature,state,g,ip,el) !********************************************************************* !* calculates plastic velocity gradient and its tangent * !* INPUT: * @@ -753,107 +841,137 @@ subroutine constitutive_dislotwin_LpAndItsTangent(Lp,dLp_dTstar,Tstar_v,Temperat !* - Lp : plastic velocity gradient * !* - dLp_dTstar : derivative of Lp (4th-rank tensor) * !********************************************************************* - use prec, only: pReal,pInt,p_vec - use math, only: math_Plain3333to99 - use mesh, only: mesh_NcpElems,mesh_maxNips - use material, only: homogenization_maxNgrains,material_phase,phase_constitutionInstance - use lattice, only: lattice_Sslip,lattice_Sslip_v,lattice_Stwin,lattice_Stwin_v,lattice_maxNslipFamily,lattice_maxNtwinFamily, & - lattice_NslipSystem,lattice_NtwinSystem,lattice_shearTwin +use prec, only: pReal,pInt,p_vec +use math, only: math_Plain3333to99 +use mesh, only: mesh_NcpElems,mesh_maxNips +use material, only: homogenization_maxNgrains,material_phase,phase_constitutionInstance +use lattice, only: lattice_Sslip,lattice_Sslip_v,lattice_Stwin,lattice_Stwin_v,lattice_maxNslipFamily,lattice_maxNtwinFamily, & + lattice_NslipSystem,lattice_NtwinSystem,lattice_shearTwin +implicit none - implicit none - - !* Definition of variables - integer(pInt) ipc,ip,el - integer(pInt) matID,structID,ns,nt,f,i,j,k,l,m,n,index_myFamily - real(pReal) Temperature,sumf - type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems) :: state - real(pReal), dimension(6) :: Tstar_v - real(pReal), dimension(3,3) :: Lp - real(pReal), dimension(3,3,3,3) :: dLp_dTstar3333 - real(pReal), dimension(9,9) :: dLp_dTstar - real(pReal), dimension(constitutive_dislotwin_totalNslip(phase_constitutionInstance(material_phase(ipc,ip,el)))) :: & +!* Input-Output variables +integer(pInt), intent(in) :: g,ip,el +real(pReal), intent(in) :: Temperature +real(pReal), dimension(6), intent(in) :: Tstar_v +type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(inout) :: state +real(pReal), dimension(3,3), intent(out) :: Lp +real(pReal), dimension(9,9), intent(out) :: dLp_dTstar +!* Local variables +integer(pInt) myInstance,myStructure,ns,nt,f,i,j,k,l,m,n,index_myFamily +real(pReal) sumf,StressRatio_p,StressRatio_pminus1,StressRatio_r,BoltzmannRatio,DotGamma0 +real(pReal), dimension(3,3,3,3) :: dLp_dTstar3333 +real(pReal), dimension(constitutive_dislotwin_totalNslip(phase_constitutionInstance(material_phase(g,ip,el)))) :: & gdot_slip,dgdot_dtauslip,tau_slip - real(pReal), dimension(constitutive_dislotwin_totalNtwin(phase_constitutionInstance(material_phase(ipc,ip,el)))) :: & +real(pReal), dimension(constitutive_dislotwin_totalNtwin(phase_constitutionInstance(material_phase(g,ip,el)))) :: & gdot_twin,dgdot_dtautwin,tau_twin - !* Shortened notation - matID = phase_constitutionInstance(material_phase(ipc,ip,el)) - structID = constitutive_dislotwin_structure(matID) - ns = constitutive_dislotwin_totalNslip(matID) - nt = constitutive_dislotwin_totalNtwin(matID) +!* Shortened notation +myInstance = phase_constitutionInstance(material_phase(g,ip,el)) +myStructure = constitutive_dislotwin_structure(myInstance) +ns = constitutive_dislotwin_totalNslip(myInstance) +nt = constitutive_dislotwin_totalNtwin(myInstance) - !* Total twin volume fraction - sumf = sum(state(ipc,ip,el)%p((ns+1):(ns+nt))) ! safe for nt == 0 +!* Total twin volume fraction +sumf = sum(state(g,ip,el)%p((2*ns+1):(2*ns+nt))) ! safe for nt == 0 - Lp = 0.0_pReal - dLp_dTstar3333 = 0.0_pReal - dLp_dTstar = 0.0_pReal +Lp = 0.0_pReal +dLp_dTstar3333 = 0.0_pReal +dLp_dTstar = 0.0_pReal - !* Dislocation glide part - gdot_slip = 0.0_pReal - dgdot_dtauslip = 0.0_pReal - j = 0_pInt - do f = 1,lattice_maxNslipFamily ! loop over all slip families - index_myFamily = sum(lattice_NslipSystem(1:f-1,structID)) ! at which index starts my family - do i = 1,constitutive_dislotwin_Nslip(f,matID) ! process each (active) slip system in family - j = j+1_pInt +!* Dislocation glide part +gdot_slip = 0.0_pReal +dgdot_dtauslip = 0.0_pReal +j = 0_pInt +do f = 1,lattice_maxNslipFamily ! loop over all slip families + index_myFamily = sum(lattice_NslipSystem(1:f-1,myStructure)) ! at which index starts my family + do i = 1,constitutive_dislotwin_Nslip(f,myInstance) ! process each (active) slip system in family + j = j+1_pInt - !* Calculation of Lp - tau_slip(j) = dot_product(Tstar_v,lattice_Sslip_v(:,index_myFamily+i,structID)) - if ( abs(tau_slip(j)) > state(ipc,ip,el)%p(6*ns+3*nt+j) ) then - - gdot_slip(j) = state(ipc,ip,el)%p(9*ns+5*nt+j)*sign(1.0_pReal,tau_slip(j))*& - sinh(((abs(tau_slip(j))-state(ipc,ip,el)%p(6*ns+3*nt+j))*state(ipc,ip,el)%p(7*ns+4*nt+j))/(kB*Temperature)) - - dgdot_dtauslip(j) = (state(ipc,ip,el)%p(9*ns+5*nt+j)*state(ipc,ip,el)%p(7*ns+4*nt+j))/(kB*Temperature)*& - cosh(((abs(tau_slip(j))-state(ipc,ip,el)%p(6*ns+3*nt+j))*state(ipc,ip,el)%p(7*ns+4*nt+j))/(kB*Temperature)) - - endif - Lp = Lp + (1.0_pReal - sumf)*gdot_slip(j)*lattice_Sslip(:,:,index_myFamily+i,structID) + !* Calculation of Lp + + !* Resolved shear stress on slip system + tau_slip(j) = dot_product(Tstar_v,lattice_Sslip_v(:,index_myFamily+i,myStructure)) + !* Stress ratios + StressRatio_p = (abs(tau_slip(j))/state(g,ip,el)%p(5*ns+3*nt+j))**constitutive_dislotwin_p(myInstance) + StressRatio_pminus1 = (abs(tau_slip(j))/state(g,ip,el)%p(5*ns+3*nt+j))**(constitutive_dislotwin_p(myInstance)-1.0_pReal) + !* Boltzmann ratio + BoltzmannRatio = constitutive_dislotwin_QedgePerSlipSystem(f,myInstance)/(kB*Temperature) + !* Initial shear rates + DotGamma0 = & + state(g,ip,el)%p(j)*constitutive_dislotwin_burgersPerSlipSystem(f,myInstance)*constitutive_dislotwin_v0PerSlipSystem(f,myInstance) + + !* Shear rates due to slip + gdot_slip(j) = DotGamma0*exp(-BoltzmannRatio*(1-StressRatio_p)**constitutive_dislotwin_q(myInstance))*sign(1.0_pReal,tau_slip(j)) + + !* Derivatives of shear rates + dgdot_dtauslip(j) = & + ((gdot_slip(j)*BoltzmannRatio*& + constitutive_dislotwin_p(myInstance)*constitutive_dislotwin_q(myInstance))/state(g,ip,el)%p(5*ns+3*nt+j))*& + StressRatio_pminus1*(1-StressRatio_p)**(constitutive_dislotwin_q(myInstance)-1.0_pReal) + + !* Plastic velocity gradient for dislocation glide + Lp = Lp + (1.0_pReal - sumf)*gdot_slip(j)*lattice_Sslip(:,:,index_myFamily+i,myStructure) - !* Calculation of the tangent of Lp - forall (k=1:3,l=1:3,m=1:3,n=1:3) & - dLp_dTstar3333(k,l,m,n) = dLp_dTstar3333(k,l,m,n) + dgdot_dtauslip(j)*lattice_Sslip(k,l,index_myFamily+i,structID) & - *lattice_Sslip(m,n,index_myFamily+i,structID) - enddo - enddo + !* Calculation of the tangent of Lp + forall (k=1:3,l=1:3,m=1:3,n=1:3) & + dLp_dTstar3333(k,l,m,n) = & + dLp_dTstar3333(k,l,m,n) + dgdot_dtauslip(j)*& + lattice_Sslip(k,l,index_myFamily+i,myStructure)*& + lattice_Sslip(m,n,index_myFamily+i,myStructure) + enddo +enddo - !* Mechanical twinning part - gdot_twin = 0.0_pReal - dgdot_dtautwin = 0.0_pReal - j = 0_pInt - do f = 1,lattice_maxNtwinFamily ! loop over all slip families - index_myFamily = sum(lattice_NtwinSystem(1:f-1,structID)) ! at which index starts my family - do i = 1,constitutive_dislotwin_Ntwin(f,matID) ! process each (active) slip system in family - j = j+1_pInt +!* Mechanical twinning part +gdot_twin = 0.0_pReal +dgdot_dtautwin = 0.0_pReal +j = 0_pInt +do f = 1,lattice_maxNtwinFamily ! loop over all slip families + index_myFamily = sum(lattice_NtwinSystem(1:f-1,myStructure)) ! at which index starts my family + do i = 1,constitutive_dislotwin_Ntwin(f,myInstance) ! process each (active) slip system in family + j = j+1_pInt - !* Calculation of Lp - tau_twin(j) = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,structID)) - if ( tau_twin(j) > 0.0_pReal ) then - - gdot_twin(j) = (constitutive_dislotwin_fmax(matID) - sumf)*lattice_shearTwin(index_myFamily+i,structID)*& - state(ipc,ip,el)%p(8*ns+4*nt+j)*constitutive_dislotwin_Ndot0(f,matID)*& - exp(-(state(ipc,ip,el)%p(7*ns+3*nt+j)/tau_twin(j))**10.0_pReal) + !* Calculation of Lp + + !* Resolved shear stress on twin system + tau_twin(j) = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,myStructure)) + !* Stress ratios + StressRatio_r = (state(g,ip,el)%p(6*ns+3*nt+j)/tau_twin(j))**constitutive_dislotwin_r(myInstance) + + !* Shear rates and their derivatives due to twin + if ( tau_twin(j) > 0.0_pReal ) then + gdot_twin(j) = & + (constitutive_dislotwin_MaxTwinFraction(myInstance)-sumf)*lattice_shearTwin(index_myFamily+i,myStructure)*& + state(g,ip,el)%p(6*ns+4*nt+j)*constitutive_dislotwin_Ndot0PerTwinSystem(f,myInstance)*exp(-StressRatio_r) - dgdot_dtautwin(j) = (gdot_twin(j)*10.0_pReal*state(ipc,ip,el)%p(7*ns+3*nt+j)**10.0_pReal)/(tau_twin(j)**11.0_pReal) + dgdot_dtautwin(j) = ((gdot_twin(j)*constitutive_dislotwin_r(myInstance))/tau_twin(j))*StressRatio_r + endif + + !* Plastic velocity gradient for mechanical twinning + Lp = Lp + gdot_twin(j)*lattice_Stwin(:,:,index_myFamily+i,myStructure) - endif - Lp = Lp + gdot_twin(j)*lattice_Stwin(:,:,index_myFamily+i,structID) + !* Calculation of the tangent of Lp + forall (k=1:3,l=1:3,m=1:3,n=1:3) & + dLp_dTstar3333(k,l,m,n) = & + dLp_dTstar3333(k,l,m,n) + dgdot_dtautwin(j)*& + lattice_Stwin(k,l,index_myFamily+i,myStructure)*& + lattice_Stwin(m,n,index_myFamily+i,myStructure) + enddo +enddo - !* Calculation of the tangent of Lp - forall (k=1:3,l=1:3,m=1:3,n=1:3) & - dLp_dTstar3333(k,l,m,n) = dLp_dTstar3333(k,l,m,n) + dgdot_dtautwin(j)*lattice_Stwin(k,l,index_myFamily+i,structID) & - *lattice_Stwin(m,n,index_myFamily+i,structID) - enddo - enddo +dLp_dTstar = math_Plain3333to99(dLp_dTstar3333) + +!if ((ip==1).and.(el==1)) then +! write(6,*) '#MICROSTRUCTURE#' +! write(6,*) +! write(6,'(a,/,4(3(f10.4,x)/))') 'rhoEdge',state(g,ip,el)%p(1:12) +! write(6,'(a,/,4(3(f10.4,x)/))') 'rhoEdgeDip',state(g,ip,el)%p(13:24) +!endif - dLp_dTstar = math_Plain3333to99(dLp_dTstar3333) - - return +return end subroutine -function constitutive_dislotwin_dotState(Tstar_v,Temperature,state,ipc,ip,el) +function constitutive_dislotwin_dotState(Tstar_v,Temperature,state,g,ip,el) !********************************************************************* !* rate of change of microstructure * !* INPUT: * @@ -866,85 +984,152 @@ function constitutive_dislotwin_dotState(Tstar_v,Temperature,state,ipc,ip,el) !* OUTPUT: * !* - constitutive_dotState : evolution of state variable * !********************************************************************* - use prec, only: pReal,pInt,p_vec - use mesh, only: mesh_NcpElems,mesh_maxNips - use material, only: homogenization_maxNgrains,material_phase, phase_constitutionInstance - use lattice, only: lattice_Sslip,lattice_Sslip_v,lattice_Stwin,lattice_Stwin_v,lattice_maxNslipFamily,lattice_maxNtwinFamily, & +use prec, only: pReal,pInt,p_vec +use math, only: pi +use mesh, only: mesh_NcpElems,mesh_maxNips +use material, only: homogenization_maxNgrains,material_phase, phase_constitutionInstance +use lattice, only: lattice_Sslip,lattice_Sslip_v,lattice_Stwin,lattice_Stwin_v,lattice_maxNslipFamily,lattice_maxNtwinFamily, & lattice_NslipSystem,lattice_NtwinSystem,lattice_shearTwin - implicit none +implicit none -!* Definition of variables - integer(pInt) ipc,ip,el - integer(pInt) matID,structID,ns,nt,f,i,j,k,index_myFamily - real(pReal) Temperature,sumf - type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems) :: state - real(pReal), dimension(6) :: Tstar_v - real(pReal), dimension(constitutive_dislotwin_totalNslip(phase_constitutionInstance(material_phase(ipc,ip,el)))) :: & - gdot_slip,tau_slip,storage,arecovery - real(pReal), dimension(constitutive_dislotwin_totalNtwin(phase_constitutionInstance(material_phase(ipc,ip,el)))) :: & - gdot_twin,tau_twin - real(pReal), dimension(constitutive_dislotwin_sizeDotState(phase_constitutionInstance(material_phase(ipc,ip,el)))) :: & - constitutive_dislotwin_dotState +!* Input-Output variables +integer(pInt), intent(in) :: g,ip,el +real(pReal), intent(in) :: Temperature +real(pReal), dimension(6), intent(in) :: Tstar_v +type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: state +real(pReal), dimension(constitutive_dislotwin_sizeDotState(phase_constitutionInstance(material_phase(g,ip,el)))) :: & +constitutive_dislotwin_dotState +!* Local variables +integer(pInt) MyInstance,MyStructure,ns,nt,f,i,j,k,index_myFamily +real(pReal) sumf,StressRatio_p,StressRatio_pminus1,BoltzmannRatio,DotGamma0,EdgeDipMinDistance,AtomicVolume,VacancyDiffusion,StressRatio_r +real(pReal), dimension(constitutive_dislotwin_totalNslip(phase_constitutionInstance(material_phase(g,ip,el)))) :: & +gdot_slip,tau_slip,DotRhoMultiplication,EdgeDipDistance,DotRhoEdgeEdgeAnnihilation,DotRhoEdgeDipAnnihilation,& +ClimbVelocity,DotRhoEdgeDipClimb,DotRhoDipFormation +real(pReal), dimension(constitutive_dislotwin_totalNtwin(phase_constitutionInstance(material_phase(g,ip,el)))) :: gdot_twin,tau_twin - !* Shortened notation - matID = phase_constitutionInstance(material_phase(ipc,ip,el)) - structID = constitutive_dislotwin_structure(matID) - ns = constitutive_dislotwin_totalNslip(matID) - nt = constitutive_dislotwin_totalNtwin(matID) +!* Shortened notation +myInstance = phase_constitutionInstance(material_phase(g,ip,el)) +MyStructure = constitutive_dislotwin_structure(myInstance) +ns = constitutive_dislotwin_totalNslip(myInstance) +nt = constitutive_dislotwin_totalNtwin(myInstance) - !* Total twin volume fraction - sumf = sum(state(ipc,ip,el)%p((ns+1):(ns+nt))) ! safe for nt == 0 +!* Total twin volume fraction +sumf = sum(state(g,ip,el)%p((2*ns+1):(2*ns+nt))) ! safe for nt == 0 - constitutive_dislotwin_dotState = 0.0_pReal +constitutive_dislotwin_dotState = 0.0_pReal - !* Dislocation density evolution - gdot_slip = 0.0_pReal - j = 0_pInt - do f = 1,lattice_maxNslipFamily ! loop over all slip families - index_myFamily = sum(lattice_NslipSystem(1:f-1,structID)) ! at which index starts my family - do i = 1,constitutive_dislotwin_Nslip(f,matID) ! process each (active) slip system in family - j = j+1_pInt +!* Dislocation density evolution +gdot_slip = 0.0_pReal +j = 0_pInt +do f = 1,lattice_maxNslipFamily ! loop over all slip families + index_myFamily = sum(lattice_NslipSystem(1:f-1,MyStructure)) ! at which index starts my family + do i = 1,constitutive_dislotwin_Nslip(f,myInstance) ! process each (active) slip system in family + j = j+1_pInt + + !* Resolved shear stress on slip system + tau_slip(j) = dot_product(Tstar_v,lattice_Sslip_v(:,index_myFamily+i,myStructure)) + !* Stress ratios + StressRatio_p = (abs(tau_slip(j))/state(g,ip,el)%p(5*ns+3*nt+j))**constitutive_dislotwin_p(myInstance) + StressRatio_pminus1 = (abs(tau_slip(j))/state(g,ip,el)%p(5*ns+3*nt+j))**(constitutive_dislotwin_p(myInstance)-1.0_pReal) + !* Boltzmann ratio + BoltzmannRatio = constitutive_dislotwin_QedgePerSlipSystem(f,myInstance)/(kB*Temperature) + !* Initial shear rates + DotGamma0 = & + state(g,ip,el)%p(j)*constitutive_dislotwin_burgersPerSlipSystem(f,myInstance)*constitutive_dislotwin_v0PerSlipSystem(f,myInstance) + + !* Shear rates due to slip + gdot_slip(j) = DotGamma0*exp(-BoltzmannRatio*(1-StressRatio_p)**constitutive_dislotwin_q(myInstance))*sign(1.0_pReal,tau_slip(j)) + + !* Multiplication + DotRhoMultiplication(j) = abs(gdot_slip(j))/(constitutive_dislotwin_burgersPerSlipSystem(f,myInstance)*state(g,ip,el)%p(4*ns+2*nt+j)) + + !* Dipole formation + if (tau_slip(j) == 0.0_pReal) then + DotRhoDipFormation(j) = 0.0_pReal + else + EdgeDipDistance(j) = & + (3.0_pReal*constitutive_dislotwin_Gmod(myInstance)*constitutive_dislotwin_burgersPerSlipSystem(f,myInstance))/& + (16.0_pReal*pi*abs(tau_slip(j))) + DotRhoDipFormation(j) = & + ((2.0_pReal*EdgeDipDistance(j))/constitutive_dislotwin_burgersPerSlipSystem(f,myInstance))*& + state(g,ip,el)%p(j)*abs(gdot_slip(j)) + endif + + !* Spontaneous annihilation of 2 single edge dislocations + EdgeDipMinDistance = & + constitutive_dislotwin_CEdgeDipMinDistance(myInstance)*constitutive_dislotwin_burgersPerSlipSystem(f,myInstance) + DotRhoEdgeEdgeAnnihilation(j) = & + ((2.0_pReal*EdgeDipMinDistance)/constitutive_dislotwin_burgersPerSlipSystem(f,myInstance))*& + state(g,ip,el)%p(j)*abs(gdot_slip(j)) + + !* Spontaneous annihilation of a single edge dislocation with a dipole constituent + DotRhoEdgeDipAnnihilation(j) = & + ((2.0_pReal*EdgeDipMinDistance)/constitutive_dislotwin_burgersPerSlipSystem(f,myInstance))*& + state(g,ip,el)%p(ns+j)*abs(gdot_slip(j)) + + !* Dislocation dipole climb + AtomicVolume = & + constitutive_dislotwin_CAtomicVolume(myInstance)*constitutive_dislotwin_burgersPerSlipSystem(f,myInstance)**(3.0_pReal) + VacancyDiffusion = & + constitutive_dislotwin_D0(myInstance)*exp(-constitutive_dislotwin_Qsd(myInstance)/(kB*Temperature)) + ClimbVelocity(j) = & + ((3.0_pReal*constitutive_dislotwin_Gmod(myInstance)*VacancyDiffusion*AtomicVolume)/(2.0_pReal*pi*kB*Temperature))*& + (1/(EdgeDipDistance(j)+EdgeDipMinDistance)) + DotRhoEdgeDipClimb(j) = & + (4.0_pReal*ClimbVelocity(j)*state(g,ip,el)%p(ns+j))/(EdgeDipDistance(j)+EdgeDipMinDistance) + + !* Edge dislocation density rate of change + constitutive_dislotwin_dotState(j) = & + DotRhoMultiplication(j)-DotRhoDipFormation(j)-DotRhoEdgeEdgeAnnihilation(j) + + !* Edge dislocation dipole density rate of change + constitutive_dislotwin_dotState(ns+j) = & + DotRhoDipFormation(j)-DotRhoEdgeDipAnnihilation(j)-DotRhoEdgeDipClimb(j) - !* Calculation of Lp - tau_slip(j) = dot_product(Tstar_v,lattice_Sslip_v(:,index_myFamily+i,structID)) - if ( abs(tau_slip(j)) > state(ipc,ip,el)%p(6*ns+3*nt+j) ) then + enddo +enddo - gdot_slip(j) = state(ipc,ip,el)%p(9*ns+5*nt+j)*sign(1.0_pReal,tau_slip(j))* & - sinh(((abs(tau_slip(j))-state(ipc,ip,el)%p(6*ns+3*nt+j))*state(ipc,ip,el)%p(7*ns+4*nt+j))/(kB*Temperature)) - - storage(j) = (constitutive_dislotwin_Cstorage(matID)*abs(gdot_slip(j)))/& - (constitutive_dislotwin_Burgers(f,matID)*state(ipc,ip,el)%p(5*ns+2*nt+j)) - - arecovery(j) = constitutive_dislotwin_Carecovery(matID)*state(ipc,ip,el)%p(j)*abs(gdot_slip(j)) - - constitutive_dislotwin_dotState(j) = storage(j) - arecovery(j) - - endif - enddo - enddo - - !* Twin volume fraction evolution - gdot_twin = 0.0_pReal - j = 0_pInt - do f = 1,lattice_maxNtwinFamily ! loop over all twin families - index_myFamily = sum(lattice_NtwinSystem(1:f-1,structID)) ! at which index starts my family - do i = 1,constitutive_dislotwin_Ntwin(f,matID) ! process each (active) twin system in family - j = j+1_pInt - - !* Calculation of Lp - tau_twin(j) = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,structID)) - if ( tau_twin(j) > 0.0_pReal ) & - constitutive_dislotwin_dotState(ns+j) = (constitutive_dislotwin_fmax(matID) - sumf)* & - lattice_shearTwin(index_myFamily+i,structID)*state(ipc,ip,el)%p(8*ns+4*nt+j)*constitutive_dislotwin_Ndot0(f,matID)*& - exp(-(state(ipc,ip,el)%p(7*ns+3*nt+j)/tau_twin(j))**10.0_pReal) - enddo - enddo - - return +!* Twin volume fraction evolution +gdot_twin = 0.0_pReal +j = 0_pInt +do f = 1,lattice_maxNtwinFamily ! loop over all twin families + index_myFamily = sum(lattice_NtwinSystem(1:f-1,MyStructure)) ! at which index starts my family + do i = 1,constitutive_dislotwin_Ntwin(f,myInstance) ! process each (active) twin system in family + j = j+1_pInt + + !* Resolved shear stress on twin system + tau_twin(j) = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,myStructure)) + !* Stress ratios + StressRatio_r = (state(g,ip,el)%p(6*ns+3*nt+j)/tau_twin(j))**constitutive_dislotwin_r(myInstance) + + !* Shear rates and their derivatives due to twin + if ( tau_twin(j) > 0.0_pReal ) then + gdot_twin(j) = & + (constitutive_dislotwin_MaxTwinFraction(myInstance)-sumf)*& + state(g,ip,el)%p(6*ns+4*nt+j)*constitutive_dislotwin_Ndot0PerTwinSystem(f,myInstance)*exp(-StressRatio_r) + endif + enddo +enddo + +!if ((ip==1).and.(el==1)) then +! write(6,*) '#DOTSTATE#' +! write(6,*) +! write(6,'(a,/,4(3(f30.20,x)/))') 'tau slip',tau_slip +! write(6,'(a,/,4(3(f30.20,x)/))') 'gamma slip',gdot_slip +! write(6,'(a,/,4(3(f30.20,x)/))') 'RhoEdge',state(g,ip,el)%p(1:ns) +! write(6,'(a,/,4(3(f30.20,x)/))') 'Threshold Slip', state(g,ip,el)%p(5*ns+3*nt+1:6*ns+3*nt) +! write(6,'(a,/,4(3(f30.20,x)/))') 'Multiplication',DotRhoMultiplication +! write(6,'(a,/,4(3(f30.20,x)/))') 'DipFormation',DotRhoDipFormation +! write(6,'(a,/,4(3(f30.20,x)/))') 'SingleSingle',DotRhoEdgeEdgeAnnihilation +! write(6,'(a,/,4(3(f30.20,x)/))') 'SingleDipole',DotRhoEdgeDipAnnihilation +! write(6,'(a,/,4(3(f30.20,x)/))') 'DipClimb',DotRhoEdgeDipClimb +!endif + +return end function -function constitutive_dislotwin_dotTemperature(Tstar_v,Temperature,state,ipc,ip,el) +pure function constitutive_dislotwin_dotTemperature(Tstar_v,Temperature,state,g,ip,el) !********************************************************************* !* rate of change of microstructure * !* INPUT: * @@ -956,25 +1141,25 @@ function constitutive_dislotwin_dotTemperature(Tstar_v,Temperature,state,ipc,ip, !* OUTPUT: * !* - constitutive_dotTemperature : evolution of Temperature * !********************************************************************* - use prec, only: pReal,pInt,p_vec - use mesh, only: mesh_NcpElems,mesh_maxNips - use material, only: homogenization_maxNgrains - implicit none +use prec, only: pReal,pInt,p_vec +use mesh, only: mesh_NcpElems,mesh_maxNips +use material, only: homogenization_maxNgrains +implicit none -!* Definition of variables - integer(pInt) ipc,ip,el - real(pReal) Temperature - type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems) :: state - real(pReal), dimension(6) :: Tstar_v - real(pReal) constitutive_dislotwin_dotTemperature - - constitutive_dislotwin_dotTemperature = 0.0_pReal - - return +!* Input-Output variables +integer(pInt), intent(in) :: g,ip,el +real(pReal), intent(in) :: Temperature +real(pReal), dimension(6), intent(in) :: Tstar_v +type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: state +real(pReal) constitutive_dislotwin_dotTemperature + +constitutive_dislotwin_dotTemperature = 0.0_pReal + +return end function -pure function constitutive_dislotwin_postResults(Tstar_v,Temperature,dt,state,ipc,ip,el) +pure function constitutive_dislotwin_postResults(Tstar_v,Temperature,dt,state,g,ip,el) !********************************************************************* !* return array of constitutive results * !* INPUT: * @@ -985,124 +1170,140 @@ pure function constitutive_dislotwin_postResults(Tstar_v,Temperature,dt,state,ip !* - ip : current integration point * !* - el : current element * !********************************************************************* - use prec, only: pReal,pInt,p_vec - use mesh, only: mesh_NcpElems,mesh_maxNips - use material, only: homogenization_maxNgrains,material_phase,phase_constitutionInstance,phase_Noutput - use lattice, only: lattice_Sslip_v,lattice_Stwin_v,lattice_maxNslipFamily,lattice_maxNtwinFamily, & +use prec, only: pReal,pInt,p_vec +use mesh, only: mesh_NcpElems,mesh_maxNips +use material, only: homogenization_maxNgrains,material_phase,phase_constitutionInstance,phase_Noutput +use lattice, only: lattice_Sslip_v,lattice_Stwin_v,lattice_maxNslipFamily,lattice_maxNtwinFamily, & lattice_NslipSystem,lattice_NtwinSystem,lattice_shearTwin - implicit none +implicit none !* Definition of variables - integer(pInt), intent(in) :: ipc,ip,el - real(pReal), intent(in) :: dt,Temperature - real(pReal), dimension(6), intent(in) :: Tstar_v - type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: state - integer(pInt) matID,structID,ns,nt,f,o,i,c,j,index_myFamily - real(pReal) sumf,tau - real(pReal), dimension(constitutive_dislotwin_sizePostResults(phase_constitutionInstance(material_phase(ipc,ip,el)))) :: & - constitutive_dislotwin_postResults +integer(pInt), intent(in) :: g,ip,el +real(pReal), intent(in) :: dt,Temperature +real(pReal), dimension(6), intent(in) :: Tstar_v +type(p_vec), dimension(homogenization_maxNgrains,mesh_maxNips,mesh_NcpElems), intent(in) :: state +integer(pInt) myInstance,myStructure,ns,nt,f,o,i,c,j,index_myFamily +real(pReal) sumf,tau,StressRatio_p,StressRatio_pminus1,BoltzmannRatio,DotGamma0,StressRatio_r +real(pReal), dimension(constitutive_dislotwin_sizePostResults(phase_constitutionInstance(material_phase(g,ip,el)))) :: & +constitutive_dislotwin_postResults - !* Shortened notation - matID = phase_constitutionInstance(material_phase(ipc,ip,el)) - structID = constitutive_dislotwin_structure(matID) - ns = constitutive_dislotwin_totalNslip(matID) - nt = constitutive_dislotwin_totalNtwin(matID) +!* Shortened notation +myInstance = phase_constitutionInstance(material_phase(g,ip,el)) +myStructure = constitutive_dislotwin_structure(myInstance) +ns = constitutive_dislotwin_totalNslip(myInstance) +nt = constitutive_dislotwin_totalNtwin(myInstance) - !* Total twin volume fraction - sumf = sum(state(ipc,ip,el)%p((ns+1):(ns+nt))) ! safe for nt == 0 +!* Total twin volume fraction +sumf = sum(state(g,ip,el)%p((2*ns+1):(2*ns+nt))) ! safe for nt == 0 - !* Required output - c = 0_pInt - constitutive_dislotwin_postResults = 0.0_pReal +!* Required output +c = 0_pInt +constitutive_dislotwin_postResults = 0.0_pReal - do o = 1,phase_Noutput(material_phase(ipc,ip,el)) - select case(constitutive_dislotwin_output(o,matID)) +do o = 1,phase_Noutput(material_phase(g,ip,el)) + select case(constitutive_dislotwin_output(o,myInstance)) - case ('dislocationdensity') - constitutive_dislotwin_postResults(c+1:c+ns) = state(ipc,ip,el)%p(1:ns) + case ('edge_density') + constitutive_dislotwin_postResults(c+1:c+ns) = state(g,ip,el)%p(1:ns) c = c + ns + + case ('dipole_density') + constitutive_dislotwin_postResults(c+1:c+ns) = state(g,ip,el)%p(ns+1:2*ns) + c = c + ns - case ('shearrate_slip') + case ('shear_rate_slip') j = 0_pInt - do f = 1,lattice_maxNslipFamily ! loop over all slip families - index_myFamily = sum(lattice_NslipSystem(1:f-1,structID)) ! at which index starts my family - do i = 1,constitutive_dislotwin_Nslip(f,matID) ! process each (active) slip system in family + do f = 1,lattice_maxNslipFamily ! loop over all slip families + index_myFamily = sum(lattice_NslipSystem(1:f-1,myStructure)) ! at which index starts my family + do i = 1,constitutive_dislotwin_Nslip(f,myInstance) ! process each (active) slip system in family j = j + 1_pInt - tau = dot_product(Tstar_v,lattice_Sslip_v(:,index_myFamily+i,structID)) - if ( abs(tau) > state(ipc,ip,el)%p(6*ns+3*nt+j) ) then - constitutive_dislotwin_postResults(c+j) = state(ipc,ip,el)%p(9*ns+5*nt+j)*sign(1.0_pReal,tau)* & - sinh(((abs(tau)-state(ipc,ip,el)%p(6*ns+3*nt+j))*state(ipc,ip,el)%p(7*ns+4*nt+j))/(kB*Temperature)) - else - constitutive_dislotwin_postResults(c+j) = 0.0_pReal - endif - enddo ; enddo + + !* Resolved shear stress on slip system + tau = dot_product(Tstar_v,lattice_Sslip_v(:,index_myFamily+i,myStructure)) + !* Stress ratios + StressRatio_p = (abs(tau)/state(g,ip,el)%p(5*ns+3*nt+j))**constitutive_dislotwin_p(myInstance) + StressRatio_pminus1 = (abs(tau)/state(g,ip,el)%p(5*ns+3*nt+j))**(constitutive_dislotwin_p(myInstance)-1.0_pReal) + !* Boltzmann ratio + BoltzmannRatio = constitutive_dislotwin_QedgePerSlipSystem(f,myInstance)/(kB*Temperature) + !* Initial shear rates + DotGamma0 = & + state(g,ip,el)%p(j)*constitutive_dislotwin_burgersPerSlipSystem(f,myInstance)*constitutive_dislotwin_v0PerSlipSystem(f,myInstance) + + !* Shear rates due to slip + constitutive_dislotwin_postResults(c+j) = & + DotGamma0*exp(-BoltzmannRatio*(1-StressRatio_p)**constitutive_dislotwin_q(myInstance))*sign(1.0_pReal,tau) + enddo ; enddo c = c + ns case ('mfp_slip') - constitutive_dislotwin_postResults(c+1:c+ns) = state(ipc,ip,el)%p((5*ns+2*nt+1):(6*ns+2*nt)) + constitutive_dislotwin_postResults(c+1:c+ns) = state(g,ip,el)%p((4*ns+2*nt+1):(5*ns+2*nt)) c = c + ns - case ('resolvedstress_slip') + case ('resolved_stress_slip') j = 0_pInt - do f = 1,lattice_maxNslipFamily ! loop over all slip families - index_myFamily = sum(lattice_NslipSystem(1:f-1,structID)) ! at which index starts my family - do i = 1,constitutive_dislotwin_Nslip(f,matID) ! process each (active) slip system in family + do f = 1,lattice_maxNslipFamily ! loop over all slip families + index_myFamily = sum(lattice_NslipSystem(1:f-1,myStructure)) ! at which index starts my family + do i = 1,constitutive_dislotwin_Nslip(f,myInstance) ! process each (active) slip system in family j = j + 1_pInt - constitutive_dislotwin_postResults(c+j) = dot_product(Tstar_v,lattice_Sslip_v(:,index_myFamily+i,structID)) + constitutive_dislotwin_postResults(c+j) = dot_product(Tstar_v,lattice_Sslip_v(:,index_myFamily+i,myStructure)) enddo; enddo c = c + ns - case ('resistance_slip') - constitutive_dislotwin_postResults(c+1:c+ns) = state(ipc,ip,el)%p((6*ns+3*nt+1):(7*ns+3*nt)) + case ('threshold_stress_slip') + constitutive_dislotwin_postResults(c+1:c+ns) = state(g,ip,el)%p((5*ns+3*nt+1):(6*ns+3*nt)) c = c + ns - case ('volumefraction') - constitutive_dislotwin_postResults(c+1:c+nt) = state(ipc,ip,el)%p((ns+1):(ns+nt)) + case ('twin_fraction') + constitutive_dislotwin_postResults(c+1:c+nt) = state(g,ip,el)%p((2*ns+1):(2*ns+nt)) c = c + nt - case ('shearrate_twin') + case ('shear_rate_twin') if (nt > 0_pInt) then j = 0_pInt - do f = 1,lattice_maxNtwinFamily ! loop over all slip families - index_myFamily = sum(lattice_NtwinSystem(1:f-1,structID)) ! at which index starts my family - do i = 1,constitutive_dislotwin_Ntwin(f,matID) ! process each (active) slip system in family + do f = 1,lattice_maxNtwinFamily ! loop over all twin families + index_myFamily = sum(lattice_NtwinSystem(1:f-1,myStructure)) ! at which index starts my family + do i = 1,constitutive_dislotwin_Ntwin(f,myInstance) ! process each (active) twin system in family j = j + 1_pInt - tau = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,structID)) - if ( tau > 0.0_pReal ) then - constitutive_dislotwin_postResults(c+j) = (constitutive_dislotwin_fmax(matID) - sumf)* & - lattice_shearTwin(index_myFamily+i,structID)*state(ipc,ip,el)%p(8*ns+4*nt+j)* & - constitutive_dislotwin_Ndot0(f,matID)*exp(-(state(ipc,ip,el)%p(7*ns+3*nt+j)/tau)**10.0_pReal) - else - constitutive_dislotwin_postResults(c+j) = 0.0_pReal - endif + + !* Resolved shear stress on twin system + tau = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,myStructure)) + !* Stress ratios + StressRatio_r = (state(g,ip,el)%p(6*ns+3*nt+j)/tau)**constitutive_dislotwin_r(myInstance) + + !* Shear rates and their derivatives due to twin + if ( tau > 0.0_pReal ) then + constitutive_dislotwin_postResults(c+j) = & + (constitutive_dislotwin_MaxTwinFraction(myInstance)-sumf)*& + state(g,ip,el)%p(6*ns+4*nt+j)*constitutive_dislotwin_Ndot0PerTwinSystem(f,myInstance)*exp(-StressRatio_r) + endif enddo ; enddo endif c = c + nt case ('mfp_twin') - constitutive_dislotwin_postResults(c+1:c+nt) = state(ipc,ip,el)%p((6*ns+2*nt+1):(6*ns+3*nt)) + constitutive_dislotwin_postResults(c+1:c+nt) = state(g,ip,el)%p((5*ns+2*nt+1):(5*ns+3*nt)) c = c + nt - case ('resolvedstress_twin') + case ('resolved_stress_twin') if (nt > 0_pInt) then j = 0_pInt - do f = 1,lattice_maxNtwinFamily ! loop over all slip families - index_myFamily = sum(lattice_NtwinSystem(1:f-1,structID)) ! at which index starts my family - do i = 1,constitutive_dislotwin_Ntwin(f,matID) ! process each (active) slip system in family + do f = 1,lattice_maxNtwinFamily ! loop over all slip families + index_myFamily = sum(lattice_NtwinSystem(1:f-1,myStructure)) ! at which index starts my family + do i = 1,constitutive_dislotwin_Ntwin(f,myInstance) ! process each (active) slip system in family j = j + 1_pInt - constitutive_dislotwin_postResults(c+j) = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,structID)) + constitutive_dislotwin_postResults(c+j) = dot_product(Tstar_v,lattice_Stwin_v(:,index_myFamily+i,myStructure)) enddo; enddo endif c = c + nt - case ('resistance_twin') - constitutive_dislotwin_postResults(c+1:c+nt) = state(ipc,ip,el)%p((7*ns+3*nt+1):(7*ns+4*nt)) + case ('threshold_stress_twin') + constitutive_dislotwin_postResults(c+1:c+nt) = state(g,ip,el)%p((6*ns+3*nt+1):(6*ns+4*nt)) c = c + nt end select - enddo - - return +enddo + +return end function END MODULE \ No newline at end of file diff --git a/code/material.config b/code/material.config index 904b7ce4b..5f23f894f 100644 --- a/code/material.config +++ b/code/material.config @@ -209,6 +209,58 @@ interaction_twintwin 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 w0_slip 4.0 relevantResistance 1 + +[TWIP steel FeMnC] + +#(output) edge_density +#(output) dipole_density +#(output) shear_rate_slip +#(output) mfp_slip +#(output) resolved_stress_slip +#(output) threshold_stress_slip +#(output) twin_fraction +#(output) shear_rate_twin +#(output) mfp_twin +#(output) resolved_stress_twin +#(output) threshold_stress_twin + +### Material parameters ### +lattice_structure fcc +constitution dislotwin +C11 175.0e9 # From Music et al. Applied Physics Letters 91, 191904 (2007) +C12 115.0e9 +C44 135.0e9 +grainsize 2.0e-5 # Average grain size [m] + +### Dislocation glide adjusting parameters ### +Nslip 12 0 0 0 +slipburgers 2.56e-10 0 0 0 # Burgers vector of slip system [m] +rhoedge0 1.0e12 0 0 0 # Initial dislocation density [m/m³] +rhoedgedip0 1.0 0 0 0 # Initial dislocation density [m/m³] +Qedge 3.7e-19 0 0 0 # Activation energy for dislocation glide [J] +v0 1.0e-4 0 0 0 # Initial glide velocity [m/s] +CLambdaSlip 10.0 0 0 0 # Adj. parameter controlling dislocation mean free path +D0 4.0e-5 # Vacancy diffusion prefactor [m²/s] +Qsd 4.5e-19 # Activation energy for climb [J] +pexponent 1.0 # p-exponent in glide velocity +qexponent 1.0 # p-exponent in glide velocity +Catomicvolume 1.0 # Adj. parameter controlling the atomic volume [in b] +Cedgedipmindistance 1.0 # Adj. parameter controlling the minimum dipole distance [in b] +relevantRho 1.0e-200 +interactionSlipSlip 0.122 0.122 0.625 0.07 0.137 0.122 # Interaction coefficients (Kubin et al. 2008) + +### dislocation density-based constitutive parameters ### +Ntwin 12 +twinburgers 1.47e-10 # Burgers vector of twin system [m] +twinsize 5.0e-8 # Twin stack mean thickness [m] +fmax 1.0 # Maximum admissible twin volume fraction +Ndot0 0.0 # Number of potential source per volume per time [1/m³.s] +rexponent 10.0 # r-exponent in glide velocity +Cmfptwin 1.0 # Adj. parameter controlling twin mean free path +Cthresholdtwin 1.0 # Adj. parameter controlling slip threshold stress +interactionSlipTwin 0.0 1.0 # Dislocation interaction coefficients +interactionTwinTwin 0.0 1.0 # Dislocation interaction coefficients + ##################### #####################