2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
!##############################################################
|
2007-04-17 13:28:53 +05:30
|
|
|
MODULE CPFEM
|
2007-04-11 15:36:28 +05:30
|
|
|
!##############################################################
|
2007-04-17 13:28:53 +05:30
|
|
|
! *** CPFEM engine ***
|
|
|
|
!
|
|
|
|
use prec, only: pReal,pInt
|
|
|
|
implicit none
|
|
|
|
!
|
|
|
|
! ****************************************************************
|
|
|
|
! *** General variables for the material behaviour calculation ***
|
|
|
|
! ****************************************************************
|
2007-10-16 14:03:59 +05:30
|
|
|
real(pReal) CPFEM_Temperature
|
2007-04-17 13:28:53 +05:30
|
|
|
real(pReal), dimension (:,:,:), allocatable :: CPFEM_stress_all
|
|
|
|
real(pReal), dimension (:,:,:,:), allocatable :: CPFEM_jacobi_all
|
2007-03-26 15:57:34 +05:30
|
|
|
real(pReal), dimension (:,:,:,:), allocatable :: CPFEM_ffn_all
|
|
|
|
real(pReal), dimension (:,:,:,:), allocatable :: CPFEM_ffn1_all
|
2007-04-17 13:28:53 +05:30
|
|
|
real(pReal), dimension (:,:,:,:), allocatable :: CPFEM_results
|
|
|
|
real(pReal), dimension (:,:,:,:), allocatable :: CPFEM_ini_ori
|
|
|
|
real(pReal), dimension (:,:,:,:), allocatable :: CPFEM_sigma_old
|
|
|
|
real(pReal), dimension (:,:,:,:), allocatable :: CPFEM_sigma_new
|
|
|
|
real(pReal), dimension (:,:,:,:,:), allocatable :: CPFEM_Fp_old
|
|
|
|
real(pReal), dimension (:,:,:,:,:), allocatable :: CPFEM_Fp_new
|
|
|
|
real(pReal), dimension (:,:,:,:), allocatable :: CPFEM_jaco_old
|
|
|
|
integer(pInt) :: CPFEM_inc_old = 0_pInt
|
2007-03-22 17:39:37 +05:30
|
|
|
integer(pInt) :: CPFEM_subinc_old = 1_pInt
|
2007-03-26 15:57:34 +05:30
|
|
|
integer(pInt) :: CPFEM_Nresults = 3_pInt
|
2007-04-17 13:28:53 +05:30
|
|
|
logical :: CPFEM_first_call = .true.
|
2007-03-22 17:39:37 +05:30
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
CONTAINS
|
2007-03-22 17:39:37 +05:30
|
|
|
|
2007-04-11 15:36:28 +05:30
|
|
|
!*********************************************************
|
|
|
|
!*** allocate the arrays defined in module CPFEM ***
|
|
|
|
!*** and initialize them ***
|
|
|
|
!*********************************************************
|
|
|
|
SUBROUTINE CPFEM_init()
|
|
|
|
!
|
|
|
|
use prec, only: pReal,pInt
|
2007-04-11 20:58:46 +05:30
|
|
|
use math, only: math_EulertoR
|
2007-04-11 15:36:28 +05:30
|
|
|
use mesh
|
|
|
|
use constitutive
|
|
|
|
!
|
|
|
|
implicit none
|
2007-04-11 20:58:46 +05:30
|
|
|
|
|
|
|
integer(pInt) e,i,g
|
2007-04-11 15:36:28 +05:30
|
|
|
!
|
|
|
|
! *** mpie.marc parameters ***
|
2007-10-16 14:03:59 +05:30
|
|
|
CPFEM_Temperature = 0.0_pReal
|
2007-04-11 15:36:28 +05:30
|
|
|
allocate(CPFEM_ffn_all (3,3,mesh_maxNips,mesh_NcpElems)) ; CPFEM_ffn_all = 0.0_pReal
|
|
|
|
allocate(CPFEM_ffn1_all (3,3,mesh_maxNips,mesh_NcpElems)) ; CPFEM_ffn1_all = 0.0_pReal
|
|
|
|
allocate(CPFEM_stress_all( 6,mesh_maxNips,mesh_NcpElems)) ; CPFEM_stress_all = 0.0_pReal
|
|
|
|
allocate(CPFEM_jacobi_all(6,6,mesh_maxNips,mesh_NcpElems)) ; CPFEM_jacobi_all = 0.0_pReal
|
|
|
|
!
|
|
|
|
! *** User defined results !!! MISSING incorporate consti_Nresults ***
|
|
|
|
allocate(CPFEM_results(CPFEM_Nresults+constitutive_maxNresults,constitutive_maxNgrains,mesh_maxNips,mesh_NcpElems))
|
|
|
|
CPFEM_results = 0.0_pReal
|
|
|
|
!
|
|
|
|
! *** Second Piola-Kirchoff stress tensor at (t=t0) and (t=t1) ***
|
|
|
|
allocate(CPFEM_sigma_old(6,constitutive_maxNgrains,mesh_maxNips,mesh_NcpElems)) ; CPFEM_sigma_old = 0.0_pReal
|
|
|
|
allocate(CPFEM_sigma_new(6,constitutive_maxNgrains,mesh_maxNips,mesh_NcpElems)) ; CPFEM_sigma_new = 0.0_pReal
|
|
|
|
!
|
2007-04-17 13:28:53 +05:30
|
|
|
! *** Plastic deformation gradient at (t=t0) and (t=t1) ***
|
2007-04-11 20:58:46 +05:30
|
|
|
allocate(CPFEM_Fp_old(3,3,constitutive_maxNgrains,mesh_maxNips,mesh_NcpElems))
|
|
|
|
forall (e=1:mesh_NcpElems,i=1:mesh_maxNips,g=1:constitutive_maxNgrains) &
|
|
|
|
CPFEM_Fp_old(:,:,g,i,e) = math_EulerToR(constitutive_EulerAngles(:,g,i,e)) ! plastic def gradient reflects init orientation
|
2007-04-11 15:36:28 +05:30
|
|
|
allocate(CPFEM_Fp_new(3,3,constitutive_maxNgrains,mesh_maxNips,mesh_NcpElems)) ; CPFEM_Fp_new = 0.0_pReal
|
2007-04-17 13:28:53 +05:30
|
|
|
!
|
2007-04-11 15:36:28 +05:30
|
|
|
! *** Old jacobian (consistent tangent) ***
|
|
|
|
allocate(CPFEM_jaco_old(6,6,mesh_maxNips,mesh_NcpElems)) ; CPFEM_jaco_old = 0.0_pReal
|
|
|
|
!
|
|
|
|
! *** Output to MARC output file ***
|
|
|
|
write(6,*)
|
|
|
|
write(6,*) 'Arrays allocated:'
|
|
|
|
write(6,*) 'CPFEM_ffn_all: ', shape(CPFEM_ffn_all)
|
|
|
|
write(6,*) 'CPFEM_ffn1_all: ', shape(CPFEM_ffn1_all)
|
|
|
|
write(6,*) 'CPFEM_stress_all: ', shape(CPFEM_stress_all)
|
|
|
|
write(6,*) 'CPFEM_jacobi_all: ', shape(CPFEM_jacobi_all)
|
|
|
|
write(6,*) 'CPFEM_results: ', shape(CPFEM_results)
|
|
|
|
write(6,*) 'CPFEM_sigma_old: ', shape(CPFEM_sigma_old)
|
|
|
|
write(6,*) 'CPFEM_sigma_new: ', shape(CPFEM_sigma_new)
|
|
|
|
write(6,*) 'CPFEM_Fp_old: ', shape(CPFEM_Fp_old)
|
|
|
|
write(6,*) 'CPFEM_Fp_new: ', shape(CPFEM_Fp_new)
|
|
|
|
write(6,*) 'CPFEM_jaco_old: ', shape(CPFEM_jaco_old)
|
|
|
|
write(6,*)
|
|
|
|
call flush(6)
|
|
|
|
return
|
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
END SUBROUTINE
|
2007-04-11 15:36:28 +05:30
|
|
|
!
|
|
|
|
!
|
2007-03-22 17:39:37 +05:30
|
|
|
!***********************************************************************
|
2007-04-11 15:36:28 +05:30
|
|
|
!*** perform initialization at first call, update variables and ***
|
|
|
|
!*** call the actual material model ***
|
2007-03-22 17:39:37 +05:30
|
|
|
!***********************************************************************
|
2007-10-16 14:03:59 +05:30
|
|
|
SUBROUTINE CPFEM_general(ffn, ffn1, Temperature, CPFEM_inc, CPFEM_subinc, CPFEM_cn, CPFEM_dt, CPFEM_en, CPFEM_in)
|
2007-03-22 17:39:37 +05:30
|
|
|
!
|
|
|
|
use prec, only: pReal,pInt
|
2007-04-13 19:50:33 +05:30
|
|
|
use math, only: math_init
|
2007-04-23 18:53:03 +05:30
|
|
|
use mesh, only: mesh_init,mesh_FEasCP
|
2007-10-15 15:46:51 +05:30
|
|
|
use crystal, only: crystal_Init
|
2007-04-13 19:50:33 +05:30
|
|
|
use constitutive, only: constitutive_init,constitutive_state_old,constitutive_state_new
|
2007-03-22 17:39:37 +05:30
|
|
|
implicit none
|
|
|
|
!
|
2007-10-16 14:03:59 +05:30
|
|
|
real(pReal) ffn(3,3), ffn1(3,3), Temperature, CPFEM_dt
|
2007-04-23 18:53:03 +05:30
|
|
|
integer(pInt) CPFEM_inc, CPFEM_subinc, CPFEM_cn, CPFEM_en, CPFEM_in, cp_en
|
2007-03-22 17:39:37 +05:30
|
|
|
!
|
|
|
|
! initialization step
|
2007-03-22 20:19:42 +05:30
|
|
|
if (CPFEM_first_call) then
|
2007-03-22 17:39:37 +05:30
|
|
|
! three dimensional stress state ?
|
2007-04-11 15:36:28 +05:30
|
|
|
call math_init()
|
2007-03-22 20:19:42 +05:30
|
|
|
call mesh_init()
|
2007-10-15 15:46:51 +05:30
|
|
|
call crystal_Init()
|
2007-03-22 20:19:42 +05:30
|
|
|
call constitutive_init()
|
|
|
|
call CPFEM_init()
|
2007-03-26 15:57:34 +05:30
|
|
|
CPFEM_first_call = .false.
|
2007-03-22 17:39:37 +05:30
|
|
|
endif
|
2007-04-11 15:36:28 +05:30
|
|
|
if (CPFEM_inc==CPFEM_inc_old) then ! not a new increment
|
2007-03-22 17:39:37 +05:30
|
|
|
! case of a new subincrement:update starting with subinc 2
|
|
|
|
if (CPFEM_subinc > CPFEM_subinc_old) then
|
|
|
|
CPFEM_sigma_old = CPFEM_sigma_new
|
|
|
|
CPFEM_Fp_old = CPFEM_Fp_new
|
|
|
|
constitutive_state_old = constitutive_state_new
|
|
|
|
CPFEM_subinc_old = CPFEM_subinc
|
|
|
|
endif
|
2007-04-17 13:28:53 +05:30
|
|
|
else ! new increment
|
2007-03-22 17:39:37 +05:30
|
|
|
CPFEM_sigma_old = CPFEM_sigma_new
|
|
|
|
CPFEM_Fp_old = CPFEM_Fp_new
|
|
|
|
constitutive_state_old = constitutive_state_new
|
|
|
|
CPFEM_inc_old = CPFEM_inc
|
|
|
|
CPFEM_subinc_old = 1_pInt
|
|
|
|
endif
|
2007-04-23 18:53:03 +05:30
|
|
|
|
|
|
|
cp_en = mesh_FEasCP('elem',CPFEM_en)
|
2007-10-16 14:03:59 +05:30
|
|
|
CPFEM_Temperature = Temperature
|
2007-03-22 17:39:37 +05:30
|
|
|
CPFEM_ffn_all(:,:,CPFEM_in, cp_en) = ffn
|
|
|
|
CPFEM_ffn1_all(:,:,CPFEM_in, cp_en) = ffn1
|
2007-04-11 15:36:28 +05:30
|
|
|
call CPFEM_stressIP(CPFEM_cn, CPFEM_dt, cp_en, CPFEM_in)
|
2007-03-22 17:39:37 +05:30
|
|
|
return
|
2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
END SUBROUTINE
|
2007-03-22 17:39:37 +05:30
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
|
2007-04-11 15:36:28 +05:30
|
|
|
!**********************************************************
|
|
|
|
!*** calculate the material behaviour at IP level ***
|
|
|
|
!**********************************************************
|
2007-04-17 13:28:53 +05:30
|
|
|
SUBROUTINE CPFEM_stressIP(&
|
|
|
|
CPFEM_cn,& ! Cycle number
|
|
|
|
CPFEM_dt,& ! Time increment (dt)
|
|
|
|
cp_en,& ! Element number
|
|
|
|
CPFEM_in) ! Integration point number
|
2007-04-11 15:36:28 +05:30
|
|
|
|
2007-04-11 20:25:06 +05:30
|
|
|
use prec, only: pReal,pInt,ijaco,nCutback
|
2007-04-17 13:28:53 +05:30
|
|
|
use math, only: math_pDecomposition,math_RtoEuler
|
|
|
|
use IO, only: IO_error
|
|
|
|
use mesh, only: mesh_element
|
|
|
|
use constitutive
|
|
|
|
!
|
|
|
|
implicit none
|
2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
integer(pInt), parameter :: i_now = 1_pInt,i_then = 2_pInt
|
|
|
|
character(len=128) msg
|
|
|
|
integer(pInt) CPFEM_cn,cp_en,CPFEM_in,grain,i
|
2007-04-11 20:21:49 +05:30
|
|
|
logical updateJaco,error
|
2007-04-11 15:36:28 +05:30
|
|
|
real(pReal) CPFEM_dt,dt,t,volfrac
|
|
|
|
real(pReal), dimension(6) :: cs,Tstar_v
|
2007-04-11 20:21:49 +05:30
|
|
|
real(pReal), dimension(6,6) :: cd
|
|
|
|
real(pReal), dimension(3,3) :: Fe,U,R,deltaFg
|
2007-04-11 15:36:28 +05:30
|
|
|
real(pReal), dimension(3,3,2) :: Fg,Fp
|
2007-04-11 20:21:49 +05:30
|
|
|
real(pReal), dimension(constitutive_maxNstatevars,2) :: state
|
2007-04-11 15:36:28 +05:30
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
updateJaco = (mod(CPFEM_cn,ijaco)==0) ! update consistent tangent every ijaco'th iteration
|
2007-04-11 15:36:28 +05:30
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
CPFEM_stress_all(:,CPFEM_in,cp_en) = 0.0_pReal ! average Cauchy stress
|
|
|
|
if (updateJaco) CPFEM_jaco_old(:,:,CPFEM_in,cp_en) = 0.0_pReal ! average consistent tangent
|
2007-04-11 15:36:28 +05:30
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
! -------------- grain loop -----------------
|
2007-04-25 19:28:10 +05:30
|
|
|
do grain = 1,texture_Ngrains(mesh_element(4,cp_en))
|
2007-04-11 15:36:28 +05:30
|
|
|
! -------------------------------------------
|
2007-04-17 13:28:53 +05:30
|
|
|
|
2007-04-11 15:36:28 +05:30
|
|
|
i = 0_pInt ! cutback counter
|
|
|
|
state(:,i_now) = constitutive_state_old(:,grain,CPFEM_in,cp_en)
|
|
|
|
Fg(:,:,i_now) = CPFEM_ffn_all(:,:,CPFEM_in,cp_en)
|
|
|
|
Fp(:,:,i_now) = CPFEM_Fp_old(:,:,grain,CPFEM_in,cp_en)
|
2007-04-17 13:28:53 +05:30
|
|
|
|
2007-04-11 15:36:28 +05:30
|
|
|
deltaFg = CPFEM_ffn1_all(:,:,CPFEM_in,cp_en)-CPFEM_ffn_all(:,:,CPFEM_in,cp_en)
|
|
|
|
dt = CPFEM_dt
|
2007-04-17 13:28:53 +05:30
|
|
|
|
2007-06-06 20:08:06 +05:30
|
|
|
Tstar_v = CPFEM_sigma_old(:,grain,CPFEM_in,cp_en) ! use last result as initial guess
|
2007-04-11 15:36:28 +05:30
|
|
|
Fg(:,:,i_then) = Fg(:,:,i_now)
|
2007-04-17 13:28:53 +05:30
|
|
|
state(:,i_then) = 0.0_pReal ! state_old as initial guess
|
2007-04-11 15:36:28 +05:30
|
|
|
t = 0.0_pReal
|
2007-04-17 13:28:53 +05:30
|
|
|
|
2007-04-11 15:36:28 +05:30
|
|
|
! ------- crystallite integration -----------
|
|
|
|
do
|
|
|
|
! -------------------------------------------
|
|
|
|
if (t+dt < CPFEM_dt) then ! intermediate solution
|
|
|
|
t = t+dt ! next time inc
|
|
|
|
Fg(:,:,i_then) = Fg(:,:,i_then)+deltaFg ! corresponding Fg
|
|
|
|
else ! full step solution
|
|
|
|
t = CPFEM_dt ! final time
|
2007-05-14 17:43:36 +05:30
|
|
|
Fg(:,:,i_then) = CPFEM_ffn1_all(:,:,CPFEM_in,cp_en) ! final Fg
|
2007-04-11 15:36:28 +05:30
|
|
|
endif
|
2007-04-17 13:28:53 +05:30
|
|
|
|
2007-04-11 20:21:49 +05:30
|
|
|
call CPFEM_stressCrystallite(msg,cs,cd,Tstar_v,Fp(:,:,i_then),Fe,state(:,i_then),&
|
2007-04-11 15:36:28 +05:30
|
|
|
dt,cp_en,CPFEM_in,grain,updateJaco .and. t==CPFEM_dt,&
|
|
|
|
Fg(:,:,i_now),Fg(:,:,i_then),Fp(:,:,i_now),state(:,i_now))
|
|
|
|
if (msg == 'ok') then ! solution converged
|
|
|
|
if (t == CPFEM_dt) exit ! reached final "then"
|
|
|
|
else ! solution not found
|
|
|
|
i = i+1_pInt ! inc cutback counter
|
2007-08-07 13:32:31 +05:30
|
|
|
! write(6,*) 'ncut:', i
|
2007-04-11 15:36:28 +05:30
|
|
|
if (i > nCutback) then ! limit exceeded?
|
|
|
|
write(6,*) 'cutback limit --> '//msg
|
|
|
|
write(6,*) 'Grain: ',grain
|
|
|
|
write(6,*) 'Integration point: ',CPFEM_in
|
|
|
|
write(6,*) 'Element: ',mesh_element(1,cp_en)
|
|
|
|
call IO_error(600)
|
|
|
|
return ! byebye
|
|
|
|
else
|
|
|
|
t = t-dt ! rewind time
|
|
|
|
Fg(:,:,i_then) = Fg(:,:,i_then)-deltaFg ! rewind Fg
|
|
|
|
dt = 0.5_pReal*dt ! cut time-step in half
|
|
|
|
deltaFg = 0.5_pReal*deltaFg ! cut Fg-step in half
|
|
|
|
endif
|
|
|
|
endif
|
|
|
|
enddo ! crystallite integration (cutback loop)
|
2007-04-17 13:28:53 +05:30
|
|
|
|
|
|
|
! ---- update crystallite matrices at t = t1 ----
|
2007-04-11 15:36:28 +05:30
|
|
|
CPFEM_Fp_new(:,:,grain,CPFEM_in,cp_en) = Fp(:,:,i_then)
|
|
|
|
constitutive_state_new(:,grain,CPFEM_in,cp_en) = state(:,i_then)
|
|
|
|
CPFEM_sigma_new(:,grain,CPFEM_in,cp_en) = Tstar_v
|
2007-04-17 13:28:53 +05:30
|
|
|
! ---- update results plotted in MENTAT ----
|
2007-04-11 20:21:49 +05:30
|
|
|
call math_pDecomposition(Fe,U,R,error) ! polar decomposition
|
|
|
|
if (error) then
|
|
|
|
write(6,*) 'polar decomposition'
|
|
|
|
write(6,*) 'Grain: ',grain
|
|
|
|
write(6,*) 'Integration point: ',CPFEM_in
|
|
|
|
write(6,*) 'Element: ',mesh_element(1,cp_en)
|
|
|
|
call IO_error(600)
|
|
|
|
return
|
|
|
|
endif
|
|
|
|
CPFEM_results(1:3,grain,CPFEM_in,cp_en) = math_RtoEuler(transpose(R)) ! orientation
|
2007-04-11 15:36:28 +05:30
|
|
|
CPFEM_results(4:3+constitutive_Nresults(grain,CPFEM_in,cp_en),grain,CPFEM_in,cp_en) = &
|
2007-10-16 14:03:59 +05:30
|
|
|
constitutive_post_results(Tstar_v,state(:,i_then),CPFEM_dt,CPFEM_Temperature,grain,CPFEM_in,cp_en)
|
2007-04-11 15:36:28 +05:30
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
! ---- contribute to IP result ----
|
|
|
|
volfrac = constitutive_matVolFrac(grain,CPFEM_in,cp_en)*constitutive_texVolFrac(grain,CPFEM_in,cp_en)
|
2007-04-11 20:21:49 +05:30
|
|
|
CPFEM_stress_all(:,CPFEM_in,cp_en) = CPFEM_stress_all(:,CPFEM_in,cp_en)+volfrac*cs ! average Cauchy stress
|
|
|
|
if (updateJaco) CPFEM_jaco_old(:,:,CPFEM_in,cp_en) = CPFEM_jaco_old(:,:,CPFEM_in,cp_en)+volfrac*cd ! average consistent tangent
|
2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
enddo ! grain loop
|
2007-04-17 13:28:53 +05:30
|
|
|
|
|
|
|
return
|
|
|
|
END SUBROUTINE
|
2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
|
|
|
|
!********************************************************************
|
|
|
|
! Calculates the stress for a single component
|
|
|
|
! it is based on the paper by Kalidindi et al.:
|
|
|
|
! J. Mech. Phys, Solids Vol. 40, No. 3, pp. 537-569, 1992
|
|
|
|
! it is modified to use anisotropic elasticity matrix
|
|
|
|
!********************************************************************
|
2007-04-17 13:28:53 +05:30
|
|
|
subroutine CPFEM_stressCrystallite(&
|
2007-04-11 15:36:28 +05:30
|
|
|
msg,& ! return message
|
2007-03-22 17:39:37 +05:30
|
|
|
cs,& ! Cauchy stress vector
|
2007-04-11 15:36:28 +05:30
|
|
|
dcs_de,& ! consistent tangent
|
|
|
|
Tstar_v,& ! second Piola-Kirchoff stress tensor
|
|
|
|
Fp_new,& ! new plastic deformation gradient
|
2007-04-11 20:21:49 +05:30
|
|
|
Fe_new,& ! new "elastic" deformation gradient
|
2007-04-11 15:36:28 +05:30
|
|
|
state_new,& ! new state variable array
|
|
|
|
!
|
2007-04-17 13:28:53 +05:30
|
|
|
dt,& ! time increment
|
2007-04-11 15:36:28 +05:30
|
|
|
cp_en,& ! element number
|
|
|
|
CPFEM_in,& ! integration point number
|
|
|
|
grain,& ! grain number
|
|
|
|
updateJaco,& ! boolean to calculate Jacobi matrix
|
2007-04-17 13:28:53 +05:30
|
|
|
Fg_old,& ! old global deformation gradient
|
|
|
|
Fg_new,& ! new global deformation gradient
|
|
|
|
Fp_old,& ! old plastic deformation gradient
|
|
|
|
state_old) ! old state variable array
|
2007-04-11 15:36:28 +05:30
|
|
|
|
2007-03-29 12:45:12 +05:30
|
|
|
use prec, only: pReal,pInt,pert_e
|
2007-03-28 12:09:48 +05:30
|
|
|
use constitutive, only: constitutive_Nstatevars
|
2007-04-17 13:28:53 +05:30
|
|
|
use math, only: math_Mandel6to33,mapMandel
|
|
|
|
implicit none
|
2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
character(len=*) msg
|
2007-04-24 12:19:13 +05:30
|
|
|
logical updateJaco
|
2007-04-11 15:36:28 +05:30
|
|
|
integer(pInt) cp_en,CPFEM_in,grain,i
|
|
|
|
real(pReal) dt
|
2007-04-11 20:21:49 +05:30
|
|
|
real(pReal), dimension(3,3) :: Fg_old,Fg_new,Fg_pert,Fp_old,Fp_new,Fp_pert,Fe_new,Fe_pert,E_pert
|
|
|
|
real(pReal), dimension(6) :: cs,Tstar_v,Tstar_v_pert
|
2007-04-11 15:36:28 +05:30
|
|
|
real(pReal), dimension(6,6) :: dcs_de
|
2007-04-11 20:21:49 +05:30
|
|
|
real(pReal), dimension(constitutive_Nstatevars(grain,CPFEM_in,cp_en)) :: state_old,state_new,state_pert
|
|
|
|
|
|
|
|
call CPFEM_timeIntegration(msg,Fp_new,Fe_new,Tstar_v,state_new, & ! def gradients and PK2 at end of time step
|
2007-04-17 13:28:53 +05:30
|
|
|
dt,cp_en,CPFEM_in,grain,Fg_new,Fp_old,state_old)
|
2007-04-11 15:36:28 +05:30
|
|
|
if (msg /= 'ok') return
|
|
|
|
cs = CPFEM_CauchyStress(Tstar_v,Fe_new) ! Cauchy stress
|
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
if (updateJaco) then ! consistent tangent using numerical perturbation of Fg
|
|
|
|
do i = 1,6 ! Fg component
|
2007-04-11 15:36:28 +05:30
|
|
|
E_pert = 0.0_pReal
|
|
|
|
E_pert(mapMandel(1,i),mapMandel(2,i)) = E_pert(mapMandel(1,i),mapMandel(2,i)) + pert_e/2.0_pReal
|
|
|
|
E_pert(mapMandel(2,i),mapMandel(1,i)) = E_pert(mapMandel(2,i),mapMandel(1,i)) + pert_e/2.0_pReal
|
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
Fg_pert = Fg_new+matmul(E_pert,Fg_old) ! perturbated Fg
|
|
|
|
Tstar_v_pert = Tstar_v ! initial guess from end of time step
|
|
|
|
state_pert = state_new ! initial guess from end of time step
|
2007-04-11 15:36:28 +05:30
|
|
|
call CPFEM_timeIntegration(msg,Fp_pert,Fe_pert,Tstar_v_pert,state_pert, &
|
|
|
|
dt,cp_en,CPFEM_in,grain,Fg_pert,Fp_old,state_old)
|
|
|
|
if (msg /= 'ok') then
|
|
|
|
msg = 'consistent tangent --> '//msg
|
|
|
|
return
|
|
|
|
endif
|
2007-04-11 20:21:49 +05:30
|
|
|
! Remark: (perturbated) Cauchy stress is Mandel hence dcs_de(:,4:6) is too large by sqrt(2)
|
2007-04-17 13:28:53 +05:30
|
|
|
dcs_de(:,i) = (CPFEM_CauchyStress(Tstar_v_pert,Fe_pert)-cs)/pert_e
|
|
|
|
enddo
|
2007-04-11 15:36:28 +05:30
|
|
|
endif
|
2007-04-11 20:21:49 +05:30
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
return
|
2007-04-11 15:36:28 +05:30
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
END SUBROUTINE
|
2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
|
|
|
|
!***********************************************************************
|
|
|
|
!*** fully-implicit two-level time integration ***
|
|
|
|
!***********************************************************************
|
2007-04-17 13:28:53 +05:30
|
|
|
SUBROUTINE CPFEM_timeIntegration(&
|
2007-04-11 15:36:28 +05:30
|
|
|
msg,& ! return message
|
|
|
|
Fp_new,& ! new plastic deformation gradient
|
|
|
|
Fe_new,& ! new "elastic" deformation gradient
|
|
|
|
Tstar_v,& ! 2nd PK stress (taken as initial guess if /= 0)
|
|
|
|
state_new,& ! current microstructure at end of time inc (taken as guess if /= 0)
|
|
|
|
!
|
2007-04-17 13:28:53 +05:30
|
|
|
dt,& ! time increment
|
2007-04-11 15:36:28 +05:30
|
|
|
cp_en,& ! element number
|
|
|
|
CPFEM_in,& ! integration point number
|
|
|
|
grain,& ! grain number
|
2007-04-17 13:28:53 +05:30
|
|
|
Fg_new,& ! new total def gradient
|
|
|
|
Fp_old,& ! former plastic def gradient
|
2007-04-11 15:36:28 +05:30
|
|
|
state_old) ! former microstructure
|
2007-04-17 13:28:53 +05:30
|
|
|
|
2007-07-20 19:02:44 +05:30
|
|
|
use prec
|
2007-04-11 15:36:28 +05:30
|
|
|
use constitutive, only: constitutive_Nstatevars,&
|
2007-10-16 17:00:05 +05:30
|
|
|
constitutive_homogenizedC,constitutive_dotState,constitutive_LpAndItsTangent,&
|
2007-10-16 14:03:59 +05:30
|
|
|
constitutive_Microstructure
|
2007-04-17 13:28:53 +05:30
|
|
|
use math
|
2007-03-26 15:57:34 +05:30
|
|
|
implicit none
|
2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
character(len=*) msg
|
|
|
|
integer(pInt) cp_en, CPFEM_in, grain
|
|
|
|
integer(pInt) iState,iStress,dummy, i,j,k,l,m
|
2007-08-07 13:32:31 +05:30
|
|
|
real(pReal) dt,det, p_hydro
|
2007-10-15 15:46:51 +05:30
|
|
|
real(pReal), dimension(6) :: Tstar_v,dTstar_v,Rstress, T_elastic, Rstress_old
|
2007-04-24 12:19:13 +05:30
|
|
|
real(pReal), dimension(6,6) :: C_66,Jacobi,invJacobi
|
2007-04-11 15:36:28 +05:30
|
|
|
real(pReal), dimension(3,3) :: Fg_new,Fp_old,Fp_new,Fe_new,invFp_old,invFp_new,Lp,A,B,AB
|
|
|
|
real(pReal), dimension(3,3,3,3) :: dLp, LTL
|
|
|
|
real(pReal), dimension(constitutive_Nstatevars(grain, CPFEM_in, cp_en)) :: state_old,state_new,dstate,Rstate,RstateS
|
|
|
|
logical failed
|
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
msg = 'ok' ! error-free so far
|
|
|
|
|
|
|
|
call math_invert3x3(Fp_old,invFp_old,det,failed) ! inversion of Fp
|
|
|
|
if (failed) then
|
|
|
|
msg = 'inversion Fp_old'
|
|
|
|
return
|
|
|
|
endif
|
2007-04-11 15:36:28 +05:30
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
C_66 = constitutive_HomogenizedC(grain, CPFEM_in, cp_en)
|
2007-03-28 19:02:25 +05:30
|
|
|
A = matmul(Fg_new,invFp_old) ! actually Fe
|
2007-03-26 15:57:34 +05:30
|
|
|
A = matmul(transpose(A), A)
|
2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
! former state guessed, if none specified
|
|
|
|
if (all(state_new == 0.0_pReal)) state_new = state_old
|
|
|
|
RstateS = state_new
|
|
|
|
iState = 0_pInt
|
2007-08-07 13:32:31 +05:30
|
|
|
|
2007-04-11 15:36:28 +05:30
|
|
|
Rstress = Tstar_v
|
2007-08-07 13:32:31 +05:30
|
|
|
Rstress_old=Rstress
|
2007-03-28 19:02:25 +05:30
|
|
|
|
2007-04-11 15:36:28 +05:30
|
|
|
state: do ! outer iteration: state
|
|
|
|
iState = iState+1
|
|
|
|
if (iState > nState) then
|
|
|
|
msg = 'limit state iteration'
|
|
|
|
return
|
2007-04-17 13:28:53 +05:30
|
|
|
endif
|
2007-10-16 14:03:59 +05:30
|
|
|
call constitutive_Microstructure(state_new,CPFEM_Temperature,grain,CPFEM_in,cp_en)
|
2007-04-23 18:53:03 +05:30
|
|
|
iStress = 0_pInt
|
2007-04-17 13:28:53 +05:30
|
|
|
stress: do ! inner iteration: stress
|
2007-04-11 15:36:28 +05:30
|
|
|
iStress = iStress+1
|
|
|
|
if (iStress > nStress) then ! too many loops required
|
|
|
|
msg = 'limit stress iteration'
|
|
|
|
return
|
|
|
|
endif
|
2007-07-20 19:02:44 +05:30
|
|
|
p_hydro=(Tstar_v(1)+Tstar_v(2)+Tstar_v(3))/3.0_pReal
|
|
|
|
forall(i=1:3) Tstar_v(i)=Tstar_v(i)-p_hydro
|
2007-10-16 14:03:59 +05:30
|
|
|
call constitutive_LpAndItsTangent(Lp,dLp,Tstar_v,state_new,CPFEM_Temperature,grain,CPFEM_in,cp_en)
|
2007-04-11 15:36:28 +05:30
|
|
|
B = math_I3-dt*Lp
|
2007-08-07 13:32:31 +05:30
|
|
|
! B = B / math_det3x3(B)**(1.0_pReal/3.0_pReal)
|
2007-04-23 18:53:03 +05:30
|
|
|
AB = matmul(A,B)
|
2007-07-20 19:02:44 +05:30
|
|
|
T_elastic= 0.5_pReal*matmul(C_66,math_Mandel33to6(matmul(transpose(B),AB)-math_I3))
|
|
|
|
p_hydro=(T_elastic(1)+T_elastic(2)+T_elastic(3))/3.0_pReal
|
|
|
|
forall(i=1:3) T_elastic(i)=T_elastic(i)-p_hydro
|
|
|
|
Rstress = Tstar_v - T_elastic
|
2007-08-07 13:32:31 +05:30
|
|
|
! step size control: if residuum does not improve redo iteration with reduced step size
|
2007-08-13 12:42:24 +05:30
|
|
|
if(maxval(abs(Rstress)) > maxval(abs(Rstress_old)) .and. &
|
|
|
|
maxval(abs(Rstress)) > 1.0e-6 .and. iStress > 1) then
|
2007-08-07 13:32:31 +05:30
|
|
|
! write(6,*) 'Hallo', iStress
|
|
|
|
Tstar_v=Tstar_v+0.5*dTstar_v
|
|
|
|
dTstar_v=0.5*dTstar_v
|
|
|
|
cycle
|
|
|
|
endif
|
|
|
|
if (iStress > 1 .and. (maxval(abs(Tstar_v)) < 1.0e-3_pReal .or. maxval(abs(Rstress/maxval(abs(Tstar_v)))) < tol_Stress)) exit stress
|
2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
! update stress guess using inverse of dRes/dTstar (Newton--Raphson)
|
|
|
|
LTL = 0.0_pReal
|
|
|
|
do i=1,3
|
|
|
|
do j=1,3
|
|
|
|
do k=1,3
|
|
|
|
do l=1,3
|
|
|
|
do m=1,3
|
2007-08-07 13:32:31 +05:30
|
|
|
LTL(i,j,k,l) = LTL(i,j,k,l) + dLp(j,i,m,k)*AB(m,l) + AB(m,i)*dLp(m,j,k,l)
|
2007-04-11 15:36:28 +05:30
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
Jacobi = math_identity2nd(6) + 0.5_pReal*dt*matmul(C_66,math_Mandel3333to66(LTL))
|
2007-06-06 20:08:06 +05:30
|
|
|
j = 0_pInt
|
|
|
|
call math_invert6x6(Jacobi,invJacobi,dummy,failed)
|
2007-04-11 15:36:28 +05:30
|
|
|
do while (failed .and. j <= nReg)
|
2007-04-17 13:28:53 +05:30
|
|
|
forall (i=1:6) Jacobi(i,i) = 1.05_pReal*maxval(Jacobi(i,:)) ! regularization
|
2007-06-06 20:08:06 +05:30
|
|
|
call math_invert6x6(Jacobi,invJacobi,dummy,failed)
|
2007-04-17 13:28:53 +05:30
|
|
|
j = j+1
|
2007-04-11 15:36:28 +05:30
|
|
|
enddo
|
|
|
|
if (failed) then
|
|
|
|
msg = 'regularization Jacobi'
|
|
|
|
return
|
2007-04-17 13:28:53 +05:30
|
|
|
endif
|
|
|
|
dTstar_v = matmul(invJacobi,Rstress) ! correction to Tstar
|
2007-08-07 13:32:31 +05:30
|
|
|
Rstress_old=Rstress
|
2007-04-17 13:28:53 +05:30
|
|
|
Tstar_v = Tstar_v-dTstar_v
|
2007-08-07 13:32:31 +05:30
|
|
|
! write(999,*) Tstar_v, dTstar_v, Rstress
|
2007-04-11 15:36:28 +05:30
|
|
|
|
2007-04-17 13:28:53 +05:30
|
|
|
enddo stress
|
2007-06-06 20:08:06 +05:30
|
|
|
! write(6,*) 'istress', istress
|
2007-08-13 12:42:24 +05:30
|
|
|
Tstar_v = 0.5_pReal*matmul(C_66,math_Mandel33to6(matmul(transpose(B),AB)-math_I3))
|
2007-10-16 14:03:59 +05:30
|
|
|
dstate = dt*constitutive_dotState(Tstar_v,state_new,CPFEM_Temperature,grain,CPFEM_in,cp_en) ! evolution of microstructure
|
2007-04-17 13:28:53 +05:30
|
|
|
Rstate = state_new - (state_old+dstate)
|
2007-04-11 15:36:28 +05:30
|
|
|
RstateS = 0.0_pReal
|
|
|
|
forall (i=1:constitutive_Nstatevars(grain,CPFEM_in,cp_en), state_new(i)/=0.0_pReal) &
|
2007-04-13 17:00:49 +05:30
|
|
|
RstateS(i) = Rstate(i)/state_new(i)
|
2007-04-17 13:28:53 +05:30
|
|
|
state_new = state_old+dstate
|
2007-06-06 20:08:06 +05:30
|
|
|
if (maxval(abs(RstateS)) < tol_State) exit state
|
2007-04-11 15:36:28 +05:30
|
|
|
|
|
|
|
enddo state
|
2007-06-06 20:08:06 +05:30
|
|
|
! write(6,*) 'istate', istate
|
2007-08-07 13:32:31 +05:30
|
|
|
! write(999,*) 'Tstar_v raus', Tstar_v
|
|
|
|
! write(999,*)
|
2007-04-17 13:28:53 +05:30
|
|
|
|
|
|
|
invFp_new = matmul(invFp_old,B)
|
2007-04-11 15:36:28 +05:30
|
|
|
call math_invert3x3(invFp_new,Fp_new,det,failed)
|
|
|
|
if (failed) then
|
|
|
|
msg = 'inversion Fp_new'
|
2007-03-26 20:21:01 +05:30
|
|
|
return
|
|
|
|
endif
|
2007-04-17 13:28:53 +05:30
|
|
|
Fp_new = Fp_new*det**(1.0_pReal/3.0_pReal) ! det = det(InvFp_new) !!
|
|
|
|
Fe_new = matmul(Fg_new,invFp_new)
|
|
|
|
return
|
2007-04-11 15:36:28 +05:30
|
|
|
END SUBROUTINE
|
2007-04-17 13:28:53 +05:30
|
|
|
|
|
|
|
|
2007-04-11 15:36:28 +05:30
|
|
|
FUNCTION CPFEM_CauchyStress(PK_v,Fe)
|
2007-03-28 12:09:48 +05:30
|
|
|
!***********************************************************************
|
|
|
|
!*** Cauchy stress calculation ***
|
|
|
|
!***********************************************************************
|
|
|
|
use prec, only: pReal,pInt
|
2007-03-28 19:42:41 +05:30
|
|
|
use math, only: math_Mandel33to6,math_Mandel6to33,math_det3x3
|
2007-03-28 12:09:48 +05:30
|
|
|
implicit none
|
|
|
|
! *** Subroutine parameters ***
|
2007-04-11 15:36:28 +05:30
|
|
|
real(pReal) PK_v(6), Fe(3,3), CPFEM_CauchyStress(6)
|
2007-03-28 19:02:25 +05:30
|
|
|
|
2007-04-11 15:36:28 +05:30
|
|
|
CPFEM_CauchyStress = math_Mandel33to6(matmul(matmul(Fe,math_Mandel6to33(PK_v)),transpose(Fe))/math_det3x3(Fe))
|
|
|
|
return
|
|
|
|
END FUNCTION
|
2007-04-17 13:28:53 +05:30
|
|
|
|
|
|
|
|
2007-04-11 15:36:28 +05:30
|
|
|
END MODULE
|