Extraction crystal structure information from constitutive.f90. Creation of crystal.f90

This commit is contained in:
Luc Hantcherli 2007-10-11 11:36:09 +00:00
parent 7bc27e02ca
commit a6fef91c6b
3 changed files with 209 additions and 174 deletions

View File

@ -3,8 +3,8 @@
!* Module: CONSTITUTIVE * !* Module: CONSTITUTIVE *
!************************************ !************************************
!* contains: * !* contains: *
!* - parameters definition *
!* - constitutive equations * !* - constitutive equations *
!* - Schmid matrices calculation *
!* - Hardening matrices definition * !* - Hardening matrices definition *
!* - Parameters definition * !* - Parameters definition *
!* - orientations * !* - orientations *
@ -14,6 +14,7 @@ MODULE constitutive
!*** Include other modules *** !*** Include other modules ***
use prec, only: pReal,pInt use prec, only: pReal,pInt
use crystal, only: crystal_MaxMaxNslipOfStructure,crystal_MaxCrystalStructure
implicit none implicit none
! MISSING consistency check after reading 'mattex.mpie' ! MISSING consistency check after reading 'mattex.mpie'
@ -88,145 +89,19 @@ integer(pInt) constitutive_maxNresults
integer(pInt), dimension(:,:,:), allocatable :: constitutive_Nresults integer(pInt), dimension(:,:,:), allocatable :: constitutive_Nresults
real(pReal), dimension(:,:,:,:), allocatable :: constitutive_results real(pReal), dimension(:,:,:,:), allocatable :: constitutive_results
!************************************
!* Crystal structures *
!************************************
!* Number of crystal structures (1-FCC,2-BCC,3-HCP)
integer(pInt), parameter :: constitutive_MaxCrystalStructure = 3
!* Total number of slip systems per crystal structure
!* (as to be changed according the definition of slip systems)
integer(pInt), dimension(constitutive_MaxCrystalStructure), parameter :: constitutive_MaxNslipOfStructure = &
reshape((/12,48,12/),(/constitutive_MaxCrystalStructure/))
!* Maximum number of slip systems over crystal structures
integer(pInt), parameter :: constitutive_MaxMaxNslipOfStructure = 48
!* Slip direction, slip normales and Schmid matrices
real(pReal), dimension(3,3,constitutive_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: constitutive_Sslip
real(pReal), dimension(6,constitutive_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: constitutive_Sslip_v
real(pReal), dimension(3,constitutive_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: constitutive_sn
real(pReal), dimension(3,constitutive_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: constitutive_sd
!*** Slip systems for FCC structures (1) ***
!* System {111}<110> Sort according Eisenlohr&Hantcherli
data constitutive_sd(:, 1,1)/ 0, 1,-1/ ; data constitutive_sn(:, 1,1)/ 1, 1, 1/
data constitutive_sd(:, 2,1)/-1, 0, 1/ ; data constitutive_sn(:, 2,1)/ 1, 1, 1/
data constitutive_sd(:, 3,1)/ 1,-1, 0/ ; data constitutive_sn(:, 3,1)/ 1, 1, 1/
data constitutive_sd(:, 4,1)/ 0,-1,-1/ ; data constitutive_sn(:, 4,1)/-1,-1, 1/
data constitutive_sd(:, 5,1)/ 1, 0, 1/ ; data constitutive_sn(:, 5,1)/-1,-1, 1/
data constitutive_sd(:, 6,1)/-1, 1, 0/ ; data constitutive_sn(:, 6,1)/-1,-1, 1/
data constitutive_sd(:, 7,1)/ 0,-1, 1/ ; data constitutive_sn(:, 7,1)/ 1,-1,-1/
data constitutive_sd(:, 8,1)/-1, 0,-1/ ; data constitutive_sn(:, 8,1)/ 1,-1,-1/
data constitutive_sd(:, 9,1)/ 1, 1, 0/ ; data constitutive_sn(:, 9,1)/ 1,-1,-1/
data constitutive_sd(:,10,1)/ 0, 1, 1/ ; data constitutive_sn(:,10,1)/-1, 1,-1/
data constitutive_sd(:,11,1)/ 1, 0,-1/ ; data constitutive_sn(:,11,1)/-1, 1,-1/
data constitutive_sd(:,12,1)/-1,-1, 0/ ; data constitutive_sn(:,12,1)/-1, 1,-1/
!*** Slip systems for BCC structures (2) ***
!* System {110}<111>
!* Sort?
data constitutive_sd(:, 1,2)/ 1,-1, 1/ ; data constitutive_sn(:, 1,2)/ 0, 1, 1/
data constitutive_sd(:, 2,2)/-1,-1, 1/ ; data constitutive_sn(:, 2,2)/ 0, 1, 1/
data constitutive_sd(:, 3,2)/ 1, 1, 1/ ; data constitutive_sn(:, 3,2)/ 0,-1, 1/
data constitutive_sd(:, 4,2)/-1, 1, 1/ ; data constitutive_sn(:, 4,2)/ 0,-1, 1/
data constitutive_sd(:, 5,2)/-1, 1, 1/ ; data constitutive_sn(:, 5,2)/ 1, 0, 1/
data constitutive_sd(:, 6,2)/-1,-1, 1/ ; data constitutive_sn(:, 6,2)/ 1, 0, 1/
data constitutive_sd(:, 7,2)/ 1, 1, 1/ ; data constitutive_sn(:, 7,2)/-1, 0, 1/
data constitutive_sd(:, 8,2)/ 1,-1, 1/ ; data constitutive_sn(:, 8,2)/-1, 0, 1/
data constitutive_sd(:, 9,2)/-1, 1, 1/ ; data constitutive_sn(:, 9,2)/ 1, 1, 0/
data constitutive_sd(:,10,2)/-1, 1,-1/ ; data constitutive_sn(:,10,2)/ 1, 1, 0/
data constitutive_sd(:,11,2)/ 1, 1, 1/ ; data constitutive_sn(:,11,2)/-1, 1, 0/
data constitutive_sd(:,12,2)/ 1, 1,-1/ ; data constitutive_sn(:,12,2)/-1, 1, 0/
!* System {112}<111>
!* Sort?
data constitutive_sd(:,13,2)/-1, 1, 1/ ; data constitutive_sn(:,13,2)/ 2, 1, 1/
data constitutive_sd(:,14,2)/ 1, 1, 1/ ; data constitutive_sn(:,14,2)/-2, 1, 1/
data constitutive_sd(:,15,2)/ 1, 1,-1/ ; data constitutive_sn(:,15,2)/ 2,-1, 1/
data constitutive_sd(:,16,2)/ 1,-1, 1/ ; data constitutive_sn(:,16,2)/ 2, 1,-1/
data constitutive_sd(:,17,2)/ 1,-1, 1/ ; data constitutive_sn(:,17,2)/ 1, 2, 1/
data constitutive_sd(:,18,2)/ 1, 1,-1/ ; data constitutive_sn(:,18,2)/-1, 2, 1/
data constitutive_sd(:,19,2)/ 1, 1, 1/ ; data constitutive_sn(:,19,2)/ 1,-2, 1/
data constitutive_sd(:,20,2)/-1, 1, 1/ ; data constitutive_sn(:,20,2)/ 1, 2,-1/
data constitutive_sd(:,21,2)/ 1, 1,-1/ ; data constitutive_sn(:,21,2)/ 1, 1, 2/
data constitutive_sd(:,22,2)/ 1,-1, 1/ ; data constitutive_sn(:,22,2)/-1, 1, 2/
data constitutive_sd(:,23,2)/-1, 1, 1/ ; data constitutive_sn(:,23,2)/ 1,-1, 2/
data constitutive_sd(:,24,2)/ 1, 1, 1/ ; data constitutive_sn(:,24,2)/ 1, 1,-2/
!* System {123}<111>
!* Sort?
data constitutive_sd(:,25,2)/ 1, 1,-1/ ; data constitutive_sn(:,25,2)/ 1, 2, 3/
data constitutive_sd(:,26,2)/ 1,-1, 1/ ; data constitutive_sn(:,26,2)/-1, 2, 3/
data constitutive_sd(:,27,2)/-1, 1, 1/ ; data constitutive_sn(:,27,2)/ 1,-2, 3/
data constitutive_sd(:,28,2)/ 1, 1, 1/ ; data constitutive_sn(:,28,2)/ 1, 2,-3/
data constitutive_sd(:,29,2)/ 1,-1, 1/ ; data constitutive_sn(:,29,2)/ 1, 3, 2/
data constitutive_sd(:,30,2)/ 1, 1,-1/ ; data constitutive_sn(:,30,2)/-1, 3, 2/
data constitutive_sd(:,31,2)/ 1, 1, 1/ ; data constitutive_sn(:,31,2)/ 1,-3, 2/
data constitutive_sd(:,32,2)/-1, 1, 1/ ; data constitutive_sn(:,32,2)/ 1, 3,-2/
data constitutive_sd(:,33,2)/ 1, 1,-1/ ; data constitutive_sn(:,33,2)/ 2, 1, 3/
data constitutive_sd(:,34,2)/ 1,-1, 1/ ; data constitutive_sn(:,34,2)/-2, 1, 3/
data constitutive_sd(:,35,2)/-1, 1, 1/ ; data constitutive_sn(:,35,2)/ 2,-1, 3/
data constitutive_sd(:,36,2)/ 1, 1, 1/ ; data constitutive_sn(:,36,2)/ 2, 1,-3/
data constitutive_sd(:,37,2)/ 1,-1, 1/ ; data constitutive_sn(:,37,2)/ 2, 3, 1/
data constitutive_sd(:,38,2)/ 1, 1,-1/ ; data constitutive_sn(:,38,2)/-2, 3, 1/
data constitutive_sd(:,39,2)/ 1, 1, 1/ ; data constitutive_sn(:,39,2)/ 2,-3, 1/
data constitutive_sd(:,40,2)/-1, 1, 1/ ; data constitutive_sn(:,40,2)/ 2, 3,-1/
data constitutive_sd(:,41,2)/-1, 1, 1/ ; data constitutive_sn(:,41,2)/ 3, 1, 2/
data constitutive_sd(:,42,2)/ 1, 1, 1/ ; data constitutive_sn(:,42,2)/-3, 1, 2/
data constitutive_sd(:,43,2)/ 1, 1,-1/ ; data constitutive_sn(:,43,2)/ 3,-1, 2/
data constitutive_sd(:,44,2)/ 1,-1, 1/ ; data constitutive_sn(:,44,2)/ 3, 1,-2/
data constitutive_sd(:,45,2)/-1, 1, 1/ ; data constitutive_sn(:,45,2)/ 3, 2, 1/
data constitutive_sd(:,46,2)/ 1, 1, 1/ ; data constitutive_sn(:,46,2)/-3, 2, 1/
data constitutive_sd(:,47,2)/ 1, 1,-1/ ; data constitutive_sn(:,47,2)/ 3,-2, 1/
data constitutive_sd(:,48,2)/ 1,-1, 1/ ; data constitutive_sn(:,48,2)/ 3, 2,-1/
!*** Slip systems for HCP structures (3) ***
!* Basal systems {0001}<1120> (independent of c/a-ratio)
!* 1- (0 0 0 1)[-2 1 1 0]
!* 2- (0 0 0 1)[ 1 -2 1 0]
!* 3- (0 0 0 1)[ 1 1 -2 0]
!* Plane (hkil)->(hkl)
!* Direction [uvtw]->[(u-t) (v-t) w]
!* Automatical transformation from Bravais to Miller
!* not done for the moment
!* Sort?
data constitutive_sd(:, 1,3)/-1, 0, 0/ ; data constitutive_sn(:, 1,3)/ 0, 0, 1/
data constitutive_sd(:, 2,3)/ 0,-1, 0/ ; data constitutive_sn(:, 2,3)/ 0, 0, 1/
data constitutive_sd(:, 3,3)/ 1, 1, 0/ ; data constitutive_sn(:, 3,3)/ 0, 0, 1/
!* 1st type prismatic systems {1010}<1120> (independent of c/a-ratio)
!* 1- ( 0 1 -1 0)[-2 1 1 0]
!* 2- ( 1 0 -1 0)[ 1 -2 1 0]
!* 3- (-1 1 0 0)[ 1 1 -2 0]
!* Sort?
data constitutive_sd(:, 4,3)/-1, 0, 0/ ; data constitutive_sn(:, 4,3)/ 0, 1, 0/
data constitutive_sd(:, 5,3)/ 0,-1, 0/ ; data constitutive_sn(:, 5,3)/ 1, 0, 0/
data constitutive_sd(:, 6,3)/ 1, 1, 0/ ; data constitutive_sn(:, 6,3)/-1, 1, 0/
!* 1st type 1st order pyramidal systems {1011}<1120>
!* plane normales depend on the c/a-ratio
!* 1- ( 0 -1 1 1)[-2 1 1 0]
!* 2- ( 0 1 -1 1)[-2 1 1 0]
!* 3- (-1 0 1 1)[ 1 -2 1 0]
!* 4- ( 1 0 -1 1)[ 1 -2 1 0]
!* 5- (-1 1 0 1)[ 1 1 -2 0]
!* 6- ( 1 -1 0 1)[ 1 1 -2 0]
!* Sort?
data constitutive_sd(:, 7,3)/-1, 0, 0/ ; data constitutive_sn(:, 7,3)/ 0,-1, 1/
data constitutive_sd(:, 8,3)/ 0,-1, 0/ ; data constitutive_sn(:, 8,3)/ 0, 1, 1/
data constitutive_sd(:, 9,3)/ 1, 1, 0/ ; data constitutive_sn(:, 9,3)/-1, 0, 1/
data constitutive_sd(:,10,3)/-1, 0, 0/ ; data constitutive_sn(:,10,3)/ 1, 0, 1/
data constitutive_sd(:,11,3)/ 0,-1, 0/ ; data constitutive_sn(:,11,3)/-1, 1, 1/
data constitutive_sd(:,12,3)/ 1, 1, 0/ ; data constitutive_sn(:,12,3)/ 1,-1, 1/
!*********************************************** !***********************************************
!* slip-slip interaction * !* slip-slip interaction *
!*********************************************** !***********************************************
!* (defined for the moment as crystal structure property and not as material property) !* (defined for the moment as crystal structure property and not as material property)
!* (may be changed in the future) !* (may be changed in the future)
real(pReal), dimension(constitutive_MaxMaxNslipOfStructure,constitutive_MaxMaxNslipOfStructure,& real(pReal), dimension(crystal_MaxMaxNslipOfStructure,crystal_MaxMaxNslipOfStructure,&
constitutive_MaxCrystalStructure) :: constitutive_HardeningMatrix crystal_MaxCrystalStructure) :: constitutive_HardeningMatrix
real(pReal), parameter :: constitutive_LatentHardening=1.4_pReal real(pReal), parameter :: constitutive_LatentHardening=1.4_pReal
CONTAINS CONTAINS
!**************************************** !****************************************
!* - constitutive_Init !* - constitutive_Init
!* - constitutive_SchmidMatrices
!* - constitutive_HardeningMatrices !* - constitutive_HardeningMatrices
!* - constitutive_CountSections !* - constitutive_CountSections
!* - constitutive_Parse_UnknownPart !* - constitutive_Parse_UnknownPart
@ -244,50 +119,12 @@ subroutine constitutive_Init()
!************************************** !**************************************
!* Module initialization * !* Module initialization *
!************************************** !**************************************
call constitutive_SchmidMatrices()
call constitutive_HardeningMatrices() call constitutive_HardeningMatrices()
call constitutive_Parse_MatTexDat(mattexFile) call constitutive_Parse_MatTexDat(mattexFile)
call constitutive_Assignment() call constitutive_Assignment()
end subroutine end subroutine
subroutine constitutive_SchmidMatrices()
!**************************************
!* Calculation of Schmid matrices *
!**************************************
use prec, only: pReal,pInt
implicit none
!* Definition of variables
integer(pInt) i,j,k,l
real(pReal) invNorm
!* Iteration over the crystal structures
do l=1,3
!* Iteration over the systems
do k=1,constitutive_MaxNslipOfStructure(l)
!* Defintion of Schmid matrix
forall (i=1:3,j=1:3)
constitutive_Sslip(i,j,k,l)=constitutive_sd(i,k,l)*constitutive_sn(j,k,l)
endforall
!* Normalization of Schmid matrix
invNorm=dsqrt(1.0_pReal/((constitutive_sn(1,k,l)**2+constitutive_sn(2,k,l)**2+constitutive_sn(3,k,l)**2)*&
(constitutive_sd(1,k,l)**2+constitutive_sd(2,k,l)**2+constitutive_sd(3,k,l)**2)))
constitutive_Sslip(:,:,k,l)=constitutive_Sslip(:,:,k,l)*invNorm
!* Vectorization of normalized Schmid matrix
constitutive_Sslip_v(1,k,l)=constitutive_Sslip(1,1,k,l)
constitutive_Sslip_v(2,k,l)=constitutive_Sslip(2,2,k,l)
constitutive_Sslip_v(3,k,l)=constitutive_Sslip(3,3,k,l)
!* be compatible with Mandel notation of Tstar
constitutive_Sslip_v(4,k,l)=(constitutive_Sslip(1,2,k,l)+constitutive_Sslip(2,1,k,l))/dsqrt(2.0_pReal)
constitutive_Sslip_v(5,k,l)=(constitutive_Sslip(2,3,k,l)+constitutive_Sslip(3,2,k,l))/dsqrt(2.0_pReal)
constitutive_Sslip_v(6,k,l)=(constitutive_Sslip(1,3,k,l)+constitutive_Sslip(3,1,k,l))/dsqrt(2.0_pReal)
enddo
enddo
end subroutine
subroutine constitutive_HardeningMatrices() subroutine constitutive_HardeningMatrices()
!**************************************** !****************************************
!* Hardening matrix (see Kalidindi) * !* Hardening matrix (see Kalidindi) *
@ -915,6 +752,7 @@ subroutine constitutive_LpAndItsTangent(Lp,dLp_dTstar, Tstar_v,state,ipc,ip,el)
!* - dLp_dTstar : derivative of Lp (4th-order tensor) * !* - dLp_dTstar : derivative of Lp (4th-order tensor) *
!********************************************************************* !*********************************************************************
use prec, only: pReal,pInt use prec, only: pReal,pInt
use crystal, only: crystal_Sslip,crystal_Sslip_v
implicit none implicit none
!* Definition of variables !* Definition of variables
@ -931,10 +769,10 @@ matID = constitutive_matID(ipc,ip,el)
!* Calculation of Lp !* Calculation of Lp
Lp = 0.0_pReal Lp = 0.0_pReal
do i=1,material_Nslip(matID) do i=1,material_Nslip(matID)
tau_slip(i)=dot_product(Tstar_v,constitutive_Sslip_v(:,i,material_CrystalStructure(matID))) tau_slip(i)=dot_product(Tstar_v,crystal_Sslip_v(:,i,material_CrystalStructure(matID)))
gdot_slip(i)=material_gdot0_slip(matID)*(abs(tau_slip(i))/state(i))**& gdot_slip(i)=material_gdot0_slip(matID)*(abs(tau_slip(i))/state(i))**&
material_n_slip(matID)*sign(1.0_pReal,tau_slip(i)) material_n_slip(matID)*sign(1.0_pReal,tau_slip(i))
Lp=Lp+gdot_slip(i)*constitutive_Sslip(:,:,i,material_CrystalStructure(matID)) Lp=Lp+gdot_slip(i)*crystal_Sslip(:,:,i,material_CrystalStructure(matID))
enddo enddo
!* Calculation of the tangent of Lp !* Calculation of the tangent of Lp
@ -944,9 +782,9 @@ do i=1,material_Nslip(matID)
(material_n_slip(matID)-1.0_pReal)*material_n_slip(matID)/state(i) (material_n_slip(matID)-1.0_pReal)*material_n_slip(matID)/state(i)
forall (k=1:3,l=1:3,m=1:3,n=1:3) forall (k=1:3,l=1:3,m=1:3,n=1:3)
dLp_dTstar(k,l,m,n) = dLp_dTstar(k,l,m,n)+ & dLp_dTstar(k,l,m,n) = dLp_dTstar(k,l,m,n)+ &
dgdot_dtauslip(i)*constitutive_Sslip(k,l,i,material_CrystalStructure(matID))* & dgdot_dtauslip(i)*crystal_Sslip(k,l,i,material_CrystalStructure(matID))* &
(constitutive_Sslip(m,n,i,material_CrystalStructure(matID))+ & (crystal_Sslip(m,n,i,material_CrystalStructure(matID))+ &
constitutive_Sslip(n,m,i,material_CrystalStructure(matID)))/2.0_pReal ! force m,n symmetry crystal_Sslip(n,m,i,material_CrystalStructure(matID)))/2.0_pReal ! force m,n symmetry
endforall endforall
enddo enddo
@ -968,6 +806,7 @@ function constitutive_dotState(Tstar_v,state,ipc,ip,el)
!* - constitutive_DotState : evolution of state variable * !* - constitutive_DotState : evolution of state variable *
!********************************************************************* !*********************************************************************
use prec, only: pReal,pInt use prec, only: pReal,pInt
use crystal, only: crystal_Sslip_v
implicit none implicit none
!* Definition of variables !* Definition of variables
@ -982,7 +821,7 @@ matID = constitutive_matID(ipc,ip,el)
!* Self-Hardening of each system !* Self-Hardening of each system
do i=1,constitutive_Nstatevars(ipc,ip,el) do i=1,constitutive_Nstatevars(ipc,ip,el)
tau_slip = dot_product(Tstar_v,constitutive_Sslip_v(:,i,material_CrystalStructure(matID))) tau_slip = dot_product(Tstar_v,crystal_Sslip_v(:,i,material_CrystalStructure(matID)))
gdot_slip = material_gdot0_slip(matID)*(abs(tau_slip)/state(i))**& gdot_slip = material_gdot0_slip(matID)*(abs(tau_slip)/state(i))**&
material_n_slip(matID)*sign(1.0_pReal,tau_slip) material_n_slip(matID)*sign(1.0_pReal,tau_slip)
self_hardening(i)=material_h0(matID)*(1.0_pReal-state(i)/& self_hardening(i)=material_h0(matID)*(1.0_pReal-state(i)/&
@ -1009,6 +848,7 @@ function constitutive_post_results(Tstar_v,state,dt,ipc,ip,el)
!* - el : current element * !* - el : current element *
!********************************************************************* !*********************************************************************
use prec, only: pReal,pInt use prec, only: pReal,pInt
use crystal, only: crystal_Sslip_v
implicit none implicit none
!* Definition of variables !* Definition of variables
@ -1027,7 +867,7 @@ if(constitutive_Nresults(ipc,ip,el)==0) return
do i=1,material_Nslip(matID) do i=1,material_Nslip(matID)
constitutive_post_results(i) = state(i) constitutive_post_results(i) = state(i)
tau_slip=dot_product(Tstar_v,constitutive_Sslip_v(:,i,material_CrystalStructure(matID))) tau_slip=dot_product(Tstar_v,crystal_Sslip_v(:,i,material_CrystalStructure(matID)))
constitutive_post_results(i+material_Nslip(matID)) = & constitutive_post_results(i+material_Nslip(matID)) = &
dt*material_gdot0_slip(matID)*(abs(tau_slip)/state(i))**material_n_slip(matID)*sign(1.0_pReal,tau_slip) dt*material_gdot0_slip(matID)*(abs(tau_slip)/state(i))**material_n_slip(matID)*sign(1.0_pReal,tau_slip)
enddo enddo

194
trunk/crystal.f90 Normal file
View File

@ -0,0 +1,194 @@
!************************************
!* Module: CRYSTAL *
!************************************
!* contains: *
!* - Crystal structure definition *
!* - Slip system definition *
!* - Schmid matrices calculation *
!************************************
MODULE crystal
!*** Include other modules ***
use prec, only: pReal,pInt
implicit none
!************************************
!* Crystal structures *
!************************************
!* Number of crystal structures (1-FCC,2-BCC,3-HCP)
integer(pInt), parameter :: crystal_MaxCrystalStructure = 3
!* Total number of slip systems per crystal structure
!* (as to be changed according the definition of slip systems)
integer(pInt), dimension(crystal_MaxCrystalStructure), parameter :: crystal_MaxNslipOfStructure = &
reshape((/12,48,12/),(/crystal_MaxCrystalStructure/))
!* Maximum number of slip systems over crystal structures
integer(pInt), parameter :: crystal_MaxMaxNslipOfStructure = 48
!* Slip direction, slip normales and Schmid matrices
real(pReal), dimension(3,3,crystal_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: crystal_Sslip
real(pReal), dimension(6,crystal_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: crystal_Sslip_v
real(pReal), dimension(3,crystal_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: crystal_sn
real(pReal), dimension(3,crystal_MaxMaxNslipOfStructure,constitutive_MaxCrystalStructure) :: crystal_sd
!*** Slip systems for FCC structures (1) ***
!* System {111}<110> Sort according Eisenlohr&Hantcherli
data crystal_sd(:, 1,1)/ 0, 1,-1/ ; data crystal_sn(:, 1,1)/ 1, 1, 1/
data crystal_sd(:, 2,1)/-1, 0, 1/ ; data crystal_sn(:, 2,1)/ 1, 1, 1/
data crystal_sd(:, 3,1)/ 1,-1, 0/ ; data crystal_sn(:, 3,1)/ 1, 1, 1/
data crystal_sd(:, 4,1)/ 0,-1,-1/ ; data crystal_sn(:, 4,1)/-1,-1, 1/
data crystal_sd(:, 5,1)/ 1, 0, 1/ ; data crystal_sn(:, 5,1)/-1,-1, 1/
data crystal_sd(:, 6,1)/-1, 1, 0/ ; data crystal_sn(:, 6,1)/-1,-1, 1/
data crystal_sd(:, 7,1)/ 0,-1, 1/ ; data crystal_sn(:, 7,1)/ 1,-1,-1/
data crystal_sd(:, 8,1)/-1, 0,-1/ ; data crystal_sn(:, 8,1)/ 1,-1,-1/
data crystal_sd(:, 9,1)/ 1, 1, 0/ ; data crystal_sn(:, 9,1)/ 1,-1,-1/
data crystal_sd(:,10,1)/ 0, 1, 1/ ; data crystal_sn(:,10,1)/-1, 1,-1/
data crystal_sd(:,11,1)/ 1, 0,-1/ ; data crystal_sn(:,11,1)/-1, 1,-1/
data crystal_sd(:,12,1)/-1,-1, 0/ ; data crystal_sn(:,12,1)/-1, 1,-1/
!*** Slip systems for BCC structures (2) ***
!* System {110}<111>
!* Sort?
data crystal_sd(:, 1,2)/ 1,-1, 1/ ; data crystal_sn(:, 1,2)/ 0, 1, 1/
data crystal_sd(:, 2,2)/-1,-1, 1/ ; data crystal_sn(:, 2,2)/ 0, 1, 1/
data crystal_sd(:, 3,2)/ 1, 1, 1/ ; data crystal_sn(:, 3,2)/ 0,-1, 1/
data crystal_sd(:, 4,2)/-1, 1, 1/ ; data crystal_sn(:, 4,2)/ 0,-1, 1/
data crystal_sd(:, 5,2)/-1, 1, 1/ ; data crystal_sn(:, 5,2)/ 1, 0, 1/
data crystal_sd(:, 6,2)/-1,-1, 1/ ; data crystal_sn(:, 6,2)/ 1, 0, 1/
data crystal_sd(:, 7,2)/ 1, 1, 1/ ; data crystal_sn(:, 7,2)/-1, 0, 1/
data crystal_sd(:, 8,2)/ 1,-1, 1/ ; data crystal_sn(:, 8,2)/-1, 0, 1/
data crystal_sd(:, 9,2)/-1, 1, 1/ ; data crystal_sn(:, 9,2)/ 1, 1, 0/
data crystal_sd(:,10,2)/-1, 1,-1/ ; data crystal_sn(:,10,2)/ 1, 1, 0/
data crystal_sd(:,11,2)/ 1, 1, 1/ ; data crystal_sn(:,11,2)/-1, 1, 0/
data crystal_sd(:,12,2)/ 1, 1,-1/ ; data crystal_sn(:,12,2)/-1, 1, 0/
!* System {112}<111>
!* Sort?
data crystal_sd(:,13,2)/-1, 1, 1/ ; data crystal_sn(:,13,2)/ 2, 1, 1/
data crystal_sd(:,14,2)/ 1, 1, 1/ ; data crystal_sn(:,14,2)/-2, 1, 1/
data crystal_sd(:,15,2)/ 1, 1,-1/ ; data crystal_sn(:,15,2)/ 2,-1, 1/
data crystal_sd(:,16,2)/ 1,-1, 1/ ; data crystal_sn(:,16,2)/ 2, 1,-1/
data crystal_sd(:,17,2)/ 1,-1, 1/ ; data crystal_sn(:,17,2)/ 1, 2, 1/
data crystal_sd(:,18,2)/ 1, 1,-1/ ; data crystal_sn(:,18,2)/-1, 2, 1/
data crystal_sd(:,19,2)/ 1, 1, 1/ ; data crystal_sn(:,19,2)/ 1,-2, 1/
data crystal_sd(:,20,2)/-1, 1, 1/ ; data crystal_sn(:,20,2)/ 1, 2,-1/
data crystal_sd(:,21,2)/ 1, 1,-1/ ; data crystal_sn(:,21,2)/ 1, 1, 2/
data crystal_sd(:,22,2)/ 1,-1, 1/ ; data crystal_sn(:,22,2)/-1, 1, 2/
data crystal_sd(:,23,2)/-1, 1, 1/ ; data crystal_sn(:,23,2)/ 1,-1, 2/
data crystal_sd(:,24,2)/ 1, 1, 1/ ; data crystal_sn(:,24,2)/ 1, 1,-2/
!* System {123}<111>
!* Sort?
data crystal_sd(:,25,2)/ 1, 1,-1/ ; data crystal_sn(:,25,2)/ 1, 2, 3/
data crystal_sd(:,26,2)/ 1,-1, 1/ ; data crystal_sn(:,26,2)/-1, 2, 3/
data crystal_sd(:,27,2)/-1, 1, 1/ ; data crystal_sn(:,27,2)/ 1,-2, 3/
data crystal_sd(:,28,2)/ 1, 1, 1/ ; data crystal_sn(:,28,2)/ 1, 2,-3/
data crystal_sd(:,29,2)/ 1,-1, 1/ ; data crystal_sn(:,29,2)/ 1, 3, 2/
data crystal_sd(:,30,2)/ 1, 1,-1/ ; data crystal_sn(:,30,2)/-1, 3, 2/
data crystal_sd(:,31,2)/ 1, 1, 1/ ; data crystal_sn(:,31,2)/ 1,-3, 2/
data crystal_sd(:,32,2)/-1, 1, 1/ ; data crystal_sn(:,32,2)/ 1, 3,-2/
data crystal_sd(:,33,2)/ 1, 1,-1/ ; data crystal_sn(:,33,2)/ 2, 1, 3/
data crystal_sd(:,34,2)/ 1,-1, 1/ ; data crystal_sn(:,34,2)/-2, 1, 3/
data crystal_sd(:,35,2)/-1, 1, 1/ ; data crystal_sn(:,35,2)/ 2,-1, 3/
data crystal_sd(:,36,2)/ 1, 1, 1/ ; data crystal_sn(:,36,2)/ 2, 1,-3/
data crystal_sd(:,37,2)/ 1,-1, 1/ ; data crystal_sn(:,37,2)/ 2, 3, 1/
data crystal_sd(:,38,2)/ 1, 1,-1/ ; data crystal_sn(:,38,2)/-2, 3, 1/
data crystal_sd(:,39,2)/ 1, 1, 1/ ; data crystal_sn(:,39,2)/ 2,-3, 1/
data crystal_sd(:,40,2)/-1, 1, 1/ ; data crystal_sn(:,40,2)/ 2, 3,-1/
data crystal_sd(:,41,2)/-1, 1, 1/ ; data crystal_sn(:,41,2)/ 3, 1, 2/
data crystal_sd(:,42,2)/ 1, 1, 1/ ; data crystal_sn(:,42,2)/-3, 1, 2/
data crystal_sd(:,43,2)/ 1, 1,-1/ ; data crystal_sn(:,43,2)/ 3,-1, 2/
data crystal_sd(:,44,2)/ 1,-1, 1/ ; data crystal_sn(:,44,2)/ 3, 1,-2/
data crystal_sd(:,45,2)/-1, 1, 1/ ; data crystal_sn(:,45,2)/ 3, 2, 1/
data crystal_sd(:,46,2)/ 1, 1, 1/ ; data crystal_sn(:,46,2)/-3, 2, 1/
data crystal_sd(:,47,2)/ 1, 1,-1/ ; data crystal_sn(:,47,2)/ 3,-2, 1/
data crystal_sd(:,48,2)/ 1,-1, 1/ ; data crystal_sn(:,48,2)/ 3, 2,-1/
!*** Slip systems for HCP structures (3) ***
!* Basal systems {0001}<1120> (independent of c/a-ratio)
!* 1- (0 0 0 1)[-2 1 1 0]
!* 2- (0 0 0 1)[ 1 -2 1 0]
!* 3- (0 0 0 1)[ 1 1 -2 0]
!* Plane (hkil)->(hkl)
!* Direction [uvtw]->[(u-t) (v-t) w]
!* Automatical transformation from Bravais to Miller
!* not done for the moment
!* Sort?
data crystal_sd(:, 1,3)/-1, 0, 0/ ; data crystal_sn(:, 1,3)/ 0, 0, 1/
data crystal_sd(:, 2,3)/ 0,-1, 0/ ; data crystal_sn(:, 2,3)/ 0, 0, 1/
data crystal_sd(:, 3,3)/ 1, 1, 0/ ; data crystal_sn(:, 3,3)/ 0, 0, 1/
!* 1st type prismatic systems {1010}<1120> (independent of c/a-ratio)
!* 1- ( 0 1 -1 0)[-2 1 1 0]
!* 2- ( 1 0 -1 0)[ 1 -2 1 0]
!* 3- (-1 1 0 0)[ 1 1 -2 0]
!* Sort?
data crystal_sd(:, 4,3)/-1, 0, 0/ ; data crystal_sn(:, 4,3)/ 0, 1, 0/
data crystal_sd(:, 5,3)/ 0,-1, 0/ ; data crystal_sn(:, 5,3)/ 1, 0, 0/
data crystal_sd(:, 6,3)/ 1, 1, 0/ ; data crystal_sn(:, 6,3)/-1, 1, 0/
!* 1st type 1st order pyramidal systems {1011}<1120>
!* plane normales depend on the c/a-ratio
!* 1- ( 0 -1 1 1)[-2 1 1 0]
!* 2- ( 0 1 -1 1)[-2 1 1 0]
!* 3- (-1 0 1 1)[ 1 -2 1 0]
!* 4- ( 1 0 -1 1)[ 1 -2 1 0]
!* 5- (-1 1 0 1)[ 1 1 -2 0]
!* 6- ( 1 -1 0 1)[ 1 1 -2 0]
!* Sort?
data crystal_sd(:, 7,3)/-1, 0, 0/ ; data crystal_sn(:, 7,3)/ 0,-1, 1/
data crystal_sd(:, 8,3)/ 0,-1, 0/ ; data crystal_sn(:, 8,3)/ 0, 1, 1/
data crystal_sd(:, 9,3)/ 1, 1, 0/ ; data crystal_sn(:, 9,3)/-1, 0, 1/
data crystal_sd(:,10,3)/-1, 0, 0/ ; data crystal_sn(:,10,3)/ 1, 0, 1/
data crystal_sd(:,11,3)/ 0,-1, 0/ ; data crystal_sn(:,11,3)/-1, 1, 1/
data crystal_sd(:,12,3)/ 1, 1, 0/ ; data crystal_sn(:,12,3)/ 1,-1, 1/
CONTAINS
!****************************************
!* - crystal_Init
!* - crystal_SchmidMatrices
!****************************************
subroutine crystal_Init()
!**************************************
!* Module initialization *
!**************************************
call crystal_SchmidMatrices()
end subroutine
subroutine crystal_SchmidMatrices()
!**************************************
!* Calculation of Schmid matrices *
!**************************************
use prec, only: pReal,pInt
implicit none
!* Definition of variables
integer(pInt) i,j,k,l
real(pReal) invNorm
!* Iteration over the crystal structures
do l=1,3
!* Iteration over the systems
do k=1,crystal_MaxNslipOfStructure(l)
!* Defintion of Schmid matrix
forall (i=1:3,j=1:3)
crystal_Sslip(i,j,k,l)=crystal_sd(i,k,l)*crystal_sn(j,k,l)
endforall
!* Normalization of Schmid matrix
invNorm=dsqrt(1.0_pReal/((crystal_sn(1,k,l)**2+crystal_sn(2,k,l)**2+crystal_sn(3,k,l)**2)*&
(crystal_sd(1,k,l)**2+crystal_sd(2,k,l)**2+crystal_sd(3,k,l)**2)))
crystal_Sslip(:,:,k,l)=crystal_Sslip(:,:,k,l)*invNorm
!* Vectorization of normalized Schmid matrix
crystal_Sslip_v(1,k,l)=crystal_Sslip(1,1,k,l)
crystal_Sslip_v(2,k,l)=crystal_Sslip(2,2,k,l)
crystal_Sslip_v(3,k,l)=crystal_Sslip(3,3,k,l)
!* be compatible with Mandel notation of Tstar
crystal_Sslip_v(4,k,l)=(crystal_Sslip(1,2,k,l)+crystal_Sslip(2,1,k,l))/dsqrt(2.0_pReal)
crystal_Sslip_v(5,k,l)=(crystal_Sslip(2,3,k,l)+crystal_Sslip(3,2,k,l))/dsqrt(2.0_pReal)
crystal_Sslip_v(6,k,l)=(crystal_Sslip(1,3,k,l)+crystal_Sslip(3,1,k,l))/dsqrt(2.0_pReal)
enddo
enddo
end subroutine
END MODULE

View File

@ -31,6 +31,7 @@
include "math.f90" include "math.f90"
include "IO.f90" include "IO.f90"
include "mesh.f90" include "mesh.f90"
include "crystal.f90"
include "constitutive.f90" include "constitutive.f90"
include "CPFEM.f90" include "CPFEM.f90"
! !