DAMASK_EICMD/src/homogenization_mech_RGC.f90

1313 lines
59 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @author Denny Tjahjanto, Max-Planck-Institut für Eisenforschung GmbH
!> @author Franz Roters, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Relaxed grain cluster (RGC) homogenization scheme
!> Nconstituents is defined as p x q x r (cluster)
!--------------------------------------------------------------------------------------------------
2019-04-06 00:28:56 +05:30
module homogenization_mech_RGC
use prec, only: &
2019-04-06 00:18:20 +05:30
pReal
implicit none
private
2019-04-06 00:15:56 +05:30
integer, dimension(:,:), allocatable,target, public :: &
homogenization_RGC_sizePostResult
2019-04-06 00:15:56 +05:30
character(len=64), dimension(:,:), allocatable,target, public :: &
homogenization_RGC_output ! name of each post result output
2018-08-28 16:27:22 +05:30
2018-08-31 04:54:28 +05:30
enum, bind(c)
2019-01-13 02:34:03 +05:30
enumerator :: &
undefined_ID, &
constitutivework_ID, &
penaltyenergy_ID, &
volumediscrepancy_ID, &
averagerelaxrate_ID,&
maximumrelaxrate_ID,&
magnitudemismatch_ID
2018-08-31 04:54:28 +05:30
end enum
2019-01-13 02:34:03 +05:30
type, private :: tParameters
2019-04-06 00:15:56 +05:30
integer, dimension(:), allocatable :: &
2018-08-28 16:27:22 +05:30
Nconstituents
real(pReal) :: &
xiAlpha, &
ciAlpha
real(pReal), dimension(:), allocatable :: &
dAlpha, &
angles
2019-04-06 00:15:56 +05:30
integer :: &
of_debug = 0
2018-08-31 04:54:28 +05:30
integer(kind(undefined_ID)), dimension(:), allocatable :: &
2019-01-13 02:34:03 +05:30
outputID
2019-03-08 03:49:39 +05:30
end type tParameters
2018-08-28 16:27:22 +05:30
2018-11-04 02:06:49 +05:30
type, private :: tRGCstate
2018-11-03 21:20:43 +05:30
real(pReal), pointer, dimension(:) :: &
work, &
2019-01-13 03:52:13 +05:30
penaltyEnergy
2018-11-03 21:20:43 +05:30
real(pReal), pointer, dimension(:,:) :: &
2019-01-13 03:52:13 +05:30
relaxationVector
2018-11-04 02:06:49 +05:30
end type tRGCstate
type, private :: tRGCdependentState
2019-01-13 03:52:13 +05:30
real(pReal), allocatable, dimension(:) :: &
volumeDiscrepancy, &
relaxationRate_avg, &
relaxationRate_max
real(pReal), allocatable, dimension(:,:) :: &
mismatch
2018-11-04 02:06:49 +05:30
real(pReal), allocatable, dimension(:,:,:) :: &
orientation
end type tRGCdependentState
2019-01-13 03:52:13 +05:30
type(tparameters), dimension(:), allocatable, private :: &
param
type(tRGCstate), dimension(:), allocatable, private :: &
state, &
state0
2019-01-13 03:52:13 +05:30
type(tRGCdependentState), dimension(:), allocatable, private :: &
dependentState
2018-11-04 02:06:49 +05:30
public :: &
homogenization_RGC_init, &
homogenization_RGC_partitionDeformation, &
homogenization_RGC_averageStressAndItsTangent, &
homogenization_RGC_updateState, &
homogenization_RGC_postResults
private :: &
relaxationVector, &
interfaceNormal, &
getInterface, &
grain1to3, &
grain3to1, &
interface4to1, &
interface1to4
contains
!--------------------------------------------------------------------------------------------------
2018-04-17 18:12:04 +05:30
!> @brief allocates all necessary fields, reads information from material configuration file
!--------------------------------------------------------------------------------------------------
2018-11-03 21:20:43 +05:30
subroutine homogenization_RGC_init()
use debug, only: &
#ifdef DEBUG
debug_i, &
debug_e, &
#endif
debug_level, &
debug_homogenization, &
debug_levelBasic
use math, only: &
math_EulerToR, &
INRAD
2019-01-13 02:34:03 +05:30
use IO, only: &
2019-03-08 03:49:39 +05:30
IO_error
use material, only: &
2019-01-13 14:23:37 +05:30
#ifdef DEBUG
mappingHomogenization, &
#endif
homogenization_type, &
2019-03-10 15:32:32 +05:30
material_homogenizationAt, &
homogState, &
HOMOGENIZATION_RGC_ID, &
HOMOGENIZATION_RGC_LABEL, &
homogenization_typeInstance, &
homogenization_Noutput, &
homogenization_Ngrains
2019-01-13 02:34:03 +05:30
use config, only: &
config_homogenization
implicit none
2019-04-06 00:15:56 +05:30
integer :: &
2019-01-13 02:34:03 +05:30
Ninstance, &
h, i, &
2019-01-13 02:34:03 +05:30
NofMyHomog, outputSize, &
sizeState, nIntFaceTot
2018-08-31 04:54:28 +05:30
character(len=65536), dimension(0), parameter :: emptyStringArray = [character(len=65536)::]
2019-01-13 02:34:03 +05:30
integer(kind(undefined_ID)) :: &
outputID
character(len=65536), dimension(:), allocatable :: &
outputs
write(6,'(/,a)') ' <<<+- homogenization_'//HOMOGENIZATION_RGC_label//' init -+>>>'
2019-03-09 15:32:12 +05:30
write(6,'(/,a)') ' Tjahjanto et al., International Journal of Material Forming 2(1):939942, 2009'
2018-08-31 04:54:28 +05:30
write(6,'(a)') ' https://doi.org/10.1007/s12289-009-0619-1'
2019-03-09 15:32:12 +05:30
write(6,'(/,a)') ' Tjahjanto et al., Modelling and Simulation in Materials Science and Engineering 18:015006, 2010'
2018-08-31 04:54:28 +05:30
write(6,'(a)') ' https://doi.org/10.1088/0965-0393/18/1/015006'
2019-03-09 15:32:12 +05:30
Ninstance = count(homogenization_type == HOMOGENIZATION_RGC_ID)
if (iand(debug_level(debug_HOMOGENIZATION),debug_levelBasic) /= 0) &
2019-01-13 02:34:03 +05:30
write(6,'(a16,1x,i5,/)') '# instances:',Ninstance
2019-01-13 02:34:03 +05:30
allocate(param(Ninstance))
allocate(state(Ninstance))
allocate(state0(Ninstance))
2019-01-13 02:34:03 +05:30
allocate(dependentState(Ninstance))
2019-04-06 00:15:56 +05:30
allocate(homogenization_RGC_sizePostResult(maxval(homogenization_Noutput),Ninstance),source=0)
2019-01-13 02:34:03 +05:30
allocate(homogenization_RGC_output(maxval(homogenization_Noutput),Ninstance))
homogenization_RGC_output=''
2019-04-06 00:15:56 +05:30
do h = 1, size(homogenization_type)
if (homogenization_type(h) /= HOMOGENIZATION_RGC_ID) cycle
2019-01-13 02:34:03 +05:30
associate(prm => param(homogenization_typeInstance(h)), &
stt => state(homogenization_typeInstance(h)), &
st0 => state0(homogenization_typeInstance(h)), &
2019-01-13 02:34:03 +05:30
dst => dependentState(homogenization_typeInstance(h)), &
config => config_homogenization(h))
#ifdef DEBUG
2019-01-13 14:23:37 +05:30
if (h==material_homogenizationAt(debug_e)) then
prm%of_debug = mappingHomogenization(1,debug_i,debug_e)
endif
2019-01-13 02:34:03 +05:30
#endif
2018-08-31 04:54:28 +05:30
prm%Nconstituents = config%getInts('clustersize',requiredSize=3)
if (homogenization_Ngrains(h) /= product(prm%Nconstituents)) &
2019-04-06 00:15:56 +05:30
call IO_error(211,ext_msg='clustersize ('//HOMOGENIZATION_RGC_label//')')
2019-01-13 02:34:03 +05:30
prm%xiAlpha = config%getFloat('scalingparameter')
prm%ciAlpha = config%getFloat('overproportionality')
prm%dAlpha = config%getFloats('grainsize', requiredSize=3)
prm%angles = config%getFloats('clusterorientation',requiredSize=3)
2018-08-29 17:06:46 +05:30
2019-01-13 02:34:03 +05:30
outputs = config%getStrings('(output)',defaultVal=emptyStringArray)
2018-08-31 04:54:28 +05:30
allocate(prm%outputID(0))
2019-04-06 00:15:56 +05:30
do i=1, size(outputs)
2018-08-31 04:54:28 +05:30
outputID = undefined_ID
select case(outputs(i))
2019-01-13 02:34:03 +05:30
2018-08-31 04:54:28 +05:30
case('constitutivework')
outputID = constitutivework_ID
2019-04-06 00:15:56 +05:30
outputSize = 1
2018-08-31 04:54:28 +05:30
case('penaltyenergy')
outputID = penaltyenergy_ID
2019-04-06 00:15:56 +05:30
outputSize = 1
2018-08-31 04:54:28 +05:30
case('volumediscrepancy')
outputID = volumediscrepancy_ID
2019-04-06 00:15:56 +05:30
outputSize = 1
2018-08-31 04:54:28 +05:30
case('averagerelaxrate')
outputID = averagerelaxrate_ID
2019-04-06 00:15:56 +05:30
outputSize = 1
2018-08-31 04:54:28 +05:30
case('maximumrelaxrate')
outputID = maximumrelaxrate_ID
2019-04-06 00:15:56 +05:30
outputSize = 1
2018-08-31 04:54:28 +05:30
case('magnitudemismatch')
outputID = magnitudemismatch_ID
2019-04-06 00:15:56 +05:30
outputSize = 3
2019-01-13 02:34:03 +05:30
2018-08-31 04:54:28 +05:30
end select
2019-01-13 02:34:03 +05:30
if (outputID /= undefined_ID) then
2019-01-13 02:34:03 +05:30
homogenization_RGC_output(i,homogenization_typeInstance(h)) = outputs(i)
homogenization_RGC_sizePostResult(i,homogenization_typeInstance(h)) = outputSize
prm%outputID = [prm%outputID , outputID]
endif
2019-01-13 02:34:03 +05:30
2018-08-31 04:54:28 +05:30
enddo
2018-08-29 17:06:46 +05:30
2019-03-10 15:32:32 +05:30
NofMyHomog = count(material_homogenizationAt == h)
2019-04-06 00:15:56 +05:30
nIntFaceTot = 3*( (prm%Nconstituents(1)-1)*prm%Nconstituents(2)*prm%Nconstituents(3) &
+ prm%Nconstituents(1)*(prm%Nconstituents(2)-1)*prm%Nconstituents(3) &
+ prm%Nconstituents(1)*prm%Nconstituents(2)*(prm%Nconstituents(3)-1))
2019-01-13 02:34:03 +05:30
sizeState = nIntFaceTot &
2019-01-13 03:52:13 +05:30
+ size(['avg constitutive work ','average penalty energy'])
2019-01-13 02:34:03 +05:30
homogState(h)%sizeState = sizeState
homogState(h)%sizePostResults = sum(homogenization_RGC_sizePostResult(:,homogenization_typeInstance(h)))
allocate(homogState(h)%state0 (sizeState,NofMyHomog), source=0.0_pReal)
allocate(homogState(h)%subState0(sizeState,NofMyHomog), source=0.0_pReal)
allocate(homogState(h)%state (sizeState,NofMyHomog), source=0.0_pReal)
stt%relaxationVector => homogState(h)%state(1:nIntFaceTot,:)
st0%relaxationVector => homogState(h)%state0(1:nIntFaceTot,:)
2019-01-13 02:34:03 +05:30
stt%work => homogState(h)%state(nIntFaceTot+1,:)
2019-01-13 03:52:13 +05:30
stt%penaltyEnergy => homogState(h)%state(nIntFaceTot+2,:)
2019-01-13 02:34:03 +05:30
allocate(dst%volumeDiscrepancy( NofMyHomog))
allocate(dst%relaxationRate_avg( NofMyHomog))
allocate(dst%relaxationRate_max( NofMyHomog))
allocate(dst%mismatch( 3,NofMyHomog))
2018-11-04 02:06:49 +05:30
!--------------------------------------------------------------------------------------------------
2019-01-13 02:34:03 +05:30
! assigning cluster orientations
dependentState(homogenization_typeInstance(h))%orientation = spread(math_EulerToR(prm%angles*inRad),3,NofMyHomog)
!dst%orientation = spread(math_EulerToR(prm%angles*inRad),3,NofMyHomog) ifort version 18.0.1 crashes (for whatever reason)
2019-01-13 02:34:03 +05:30
end associate
2019-01-13 02:34:03 +05:30
enddo
end subroutine homogenization_RGC_init
!--------------------------------------------------------------------------------------------------
!> @brief partitions the deformation gradient onto the constituents
!--------------------------------------------------------------------------------------------------
2018-11-04 03:53:52 +05:30
subroutine homogenization_RGC_partitionDeformation(F,avgF,instance,of)
2019-01-13 02:34:03 +05:30
#ifdef DEBUG
use debug, only: &
debug_level, &
debug_homogenization, &
debug_levelExtensive
2019-01-13 02:34:03 +05:30
#endif
implicit none
real(pReal), dimension (:,:,:), intent(out) :: F !< partioned F per grain
2019-01-13 02:34:03 +05:30
real(pReal), dimension (:,:), intent(in) :: avgF !< averaged F
2019-04-06 00:15:56 +05:30
integer, intent(in) :: &
2018-11-04 03:53:52 +05:30
instance, &
2019-01-13 02:34:03 +05:30
of
2019-04-06 00:15:56 +05:30
real(pReal), dimension(3) :: aVect,nVect
integer, dimension(4) :: intFace
integer, dimension(3) :: iGrain3
integer :: iGrain,iFace,i,j
2018-11-04 03:43:20 +05:30
associate(prm => param(instance))
2019-01-13 02:34:03 +05:30
!--------------------------------------------------------------------------------------------------
! compute the deformation gradient of individual grains due to relaxations
F = 0.0_pReal
2019-04-06 00:15:56 +05:30
do iGrain = 1,product(prm%Nconstituents)
2018-11-04 12:41:35 +05:30
iGrain3 = grain1to3(iGrain,prm%Nconstituents)
2019-04-06 00:15:56 +05:30
do iFace = 1,6
2018-11-04 03:43:20 +05:30
intFace = getInterface(iFace,iGrain3) ! identifying 6 interfaces of each grain
aVect = relaxationVector(intFace,instance,of) ! get the relaxation vectors for each interface from global relaxation vector array
nVect = interfaceNormal(intFace,instance,of)
2019-04-06 00:15:56 +05:30
forall (i=1:3,j=1:3) &
F(i,j,iGrain) = F(i,j,iGrain) + aVect(i)*nVect(j) ! calculating deformation relaxations due to interface relaxation
enddo
F(1:3,1:3,iGrain) = F(1:3,1:3,iGrain) + avgF ! resulting relaxed deformation gradient
2019-01-13 02:34:03 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0) then
write(6,'(1x,a32,1x,i3)')'Deformation gradient of grain: ',iGrain
2019-04-06 00:15:56 +05:30
do i = 1,3
write(6,'(1x,3(e15.8,1x))')(F(i,j,iGrain), j = 1,3)
enddo
write(6,*)' '
flush(6)
endif
2019-01-13 02:34:03 +05:30
#endif
enddo
2019-01-13 02:34:03 +05:30
2018-11-04 03:43:20 +05:30
end associate
end subroutine homogenization_RGC_partitionDeformation
!--------------------------------------------------------------------------------------------------
!> @brief update the internal state of the homogenization scheme and tell whether "done" and
! "happy" with result
!--------------------------------------------------------------------------------------------------
function homogenization_RGC_updateState(P,F,F0,avgF,dt,dPdF,ip,el)
2016-05-29 14:15:03 +05:30
use prec, only: &
dEq0
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
use debug, only: &
debug_level, &
debug_homogenization,&
2019-01-13 14:23:37 +05:30
debug_levelExtensive
2019-01-13 03:37:35 +05:30
#endif
use math, only: &
math_invert2
use material, only: &
2018-11-04 13:30:35 +05:30
material_homogenizationAt, &
homogenization_typeInstance, &
2019-01-13 14:03:47 +05:30
mappingHomogenization
use numerics, only: &
absTol_RGC, &
relTol_RGC, &
absMax_RGC, &
relMax_RGC, &
pPert_RGC, &
maxdRelax_RGC, &
viscPower_RGC, &
viscModus_RGC, &
refRelaxRate_RGC
implicit none
2019-04-06 00:15:56 +05:30
real(pReal), dimension(:,:,:), intent(in) :: &
P,& !< array of P
F,& !< array of F
F0 !< array of initial F
2019-04-06 00:15:56 +05:30
real(pReal), dimension(:,:,:,:,:), intent(in) :: dPdF !< array of current grain stiffness
real(pReal), dimension(3,3), intent(in) :: avgF !< average F
real(pReal), intent(in) :: dt !< time increment
integer, intent(in) :: &
ip, & !< integration point number
el !< element number
2018-11-03 21:20:43 +05:30
2019-04-06 00:15:56 +05:30
logical, dimension(2) :: homogenization_RGC_updateState
2019-01-13 14:03:47 +05:30
2019-04-06 00:15:56 +05:30
integer, dimension(4) :: intFaceN,intFaceP,faceID
integer, dimension(3) :: nGDim,iGr3N,iGr3P
integer :: instance,iNum,i,j,nIntFaceTot,iGrN,iGrP,iMun,iFace,k,l,ipert,iGrain,nGrain, of
real(pReal), dimension(3,3,size(P,3)) :: R,pF,pR,D,pD
real(pReal), dimension(3,size(P,3)) :: NN,devNull
real(pReal), dimension(3) :: normP,normN,mornP,mornN
2019-01-13 14:03:47 +05:30
real(pReal) :: residMax,stresMax
logical :: error
real(pReal), dimension(:,:), allocatable :: tract,jmatrix,jnverse,smatrix,pmatrix,rmatrix
real(pReal), dimension(:), allocatable :: resid,relax,p_relax,p_resid,drelax
2019-01-13 14:03:47 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
integer, dimension(3) :: stresLoc
integer, dimension(2) :: residLoc
2019-01-13 14:03:47 +05:30
#endif
zeroTimeStep: if(dEq0(dt)) then
2016-05-29 14:15:03 +05:30
homogenization_RGC_updateState = .true. ! pretend everything is fine and return
return
2016-05-29 14:15:03 +05:30
endif zeroTimeStep
2018-11-04 13:30:35 +05:30
instance = homogenization_typeInstance(material_homogenizationAt(el))
2018-11-04 00:30:02 +05:30
of = mappingHomogenization(1,ip,el)
2019-01-13 03:37:35 +05:30
associate(stt => state(instance), st0 => state0(instance), dst => dependentState(instance), prm => param(instance))
2019-01-13 03:37:35 +05:30
!--------------------------------------------------------------------------------------------------
! get the dimension of the cluster (grains and interfaces)
nGDim = prm%Nconstituents
nGrain = product(nGDim)
2019-04-06 00:15:56 +05:30
nIntFaceTot = (nGDim(1)-1)*nGDim(2)*nGDim(3) &
+ nGDim(1)*(nGDim(2)-1)*nGDim(3) &
+ nGDim(1)*nGDim(2)*(nGDim(3)-1)
!--------------------------------------------------------------------------------------------------
! allocate the size of the global relaxation arrays/jacobian matrices depending on the size of the cluster
2019-04-06 00:15:56 +05:30
allocate(resid(3*nIntFaceTot), source=0.0_pReal)
allocate(tract(nIntFaceTot,3), source=0.0_pReal)
2019-01-13 03:37:35 +05:30
relax = stt%relaxationVector(:,of)
drelax = stt%relaxationVector(:,of) - st0%relaxationVector(:,of)
2019-01-13 02:34:03 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0) then
write(6,'(1x,a30)')'Obtained state: '
2019-04-06 00:15:56 +05:30
do i = 1,size(stt%relaxationVector(:,of))
2019-01-13 03:37:35 +05:30
write(6,'(1x,2(e15.8,1x))') stt%relaxationVector(i,of)
enddo
write(6,*)' '
endif
2019-01-13 02:34:03 +05:30
#endif
!--------------------------------------------------------------------------------------------------
! computing interface mismatch and stress penalty tensor for all interfaces of all grains
2019-01-13 03:37:35 +05:30
call stressPenalty(R,NN,avgF,F,ip,el,instance,of)
!--------------------------------------------------------------------------------------------------
! calculating volume discrepancy and stress penalty related to overall volume discrepancy
2019-01-13 03:52:13 +05:30
call volumePenalty(D,dst%volumeDiscrepancy(of),avgF,F,nGrain,instance,of)
2019-01-13 02:34:03 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0) then
do iGrain = 1,nGrain
write(6,'(1x,a30,1x,i3,1x,a4,3(1x,e15.8))')'Mismatch magnitude of grain(',iGrain,') :',&
NN(1,iGrain),NN(2,iGrain),NN(3,iGrain)
write(6,'(/,1x,a30,1x,i3)')'Stress and penalties of grain: ',iGrain
2019-04-06 00:15:56 +05:30
do i = 1,3
write(6,'(1x,3(e15.8,1x),1x,3(e15.8,1x),1x,3(e15.8,1x))')(P(i,j,iGrain), j = 1,3), &
(R(i,j,iGrain), j = 1,3), &
(D(i,j,iGrain), j = 1,3)
enddo
write(6,*)' '
enddo
endif
2019-01-13 02:34:03 +05:30
#endif
!------------------------------------------------------------------------------------------------
! computing the residual stress from the balance of traction at all (interior) interfaces
2019-04-06 00:15:56 +05:30
do iNum = 1,nIntFaceTot
faceID = interface1to4(iNum,param(instance)%Nconstituents) ! identifying the interface ID in local coordinate system (4-dimensional index)
!--------------------------------------------------------------------------------------------------
! identify the left/bottom/back grain (-|N)
iGr3N = faceID(2:4) ! identifying the grain ID in local coordinate system (3-dimensional index)
2019-04-06 00:15:56 +05:30
iGrN = grain3to1(iGr3N,param(instance)%Nconstituents) ! translate the local grain ID into global coordinate system (1-dimensional index)
intFaceN = getInterface(2*faceID(1),iGr3N)
2018-11-04 03:43:20 +05:30
normN = interfaceNormal(intFaceN,instance,of)
!--------------------------------------------------------------------------------------------------
! identify the right/up/front grain (+|P)
iGr3P = iGr3N
2019-04-06 00:15:56 +05:30
iGr3P(faceID(1)) = iGr3N(faceID(1))+1 ! identifying the grain ID in local coordinate system (3-dimensional index)
iGrP = grain3to1(iGr3P,param(instance)%Nconstituents) ! translate the local grain ID into global coordinate system (1-dimensional index)
intFaceP = getInterface(2*faceID(1)-1,iGr3P)
2018-11-04 03:43:20 +05:30
normP = interfaceNormal(intFaceP,instance,of)
!--------------------------------------------------------------------------------------------------
! compute the residual of traction at the interface (in local system, 4-dimensional index)
2019-04-06 00:15:56 +05:30
do i = 1,3
tract(iNum,i) = sign(viscModus_RGC*(abs(drelax(i+3*(iNum-1)))/(refRelaxRate_RGC*dt))**viscPower_RGC, &
drelax(i+3*(iNum-1))) ! contribution from the relaxation viscosity
do j = 1,3
tract(iNum,i) = tract(iNum,i) + (P(i,j,iGrP) + R(i,j,iGrP) + D(i,j,iGrP))*normP(j) & ! contribution from material stress P, mismatch penalty R, and volume penalty D projected into the interface
+ (P(i,j,iGrN) + R(i,j,iGrN) + D(i,j,iGrN))*normN(j)
2019-04-06 00:15:56 +05:30
resid(i+3*(iNum-1)) = tract(iNum,i) ! translate the local residual into global 1-dimensional residual array
enddo
enddo
2019-01-13 02:34:03 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0) then
write(6,'(1x,a30,1x,i3)')'Traction at interface: ',iNum
2019-04-06 00:15:56 +05:30
write(6,'(1x,3(e15.8,1x))')(tract(iNum,j), j = 1,3)
write(6,*)' '
endif
2019-01-13 02:34:03 +05:30
#endif
enddo
!--------------------------------------------------------------------------------------------------
! convergence check for stress residual
stresMax = maxval(abs(P)) ! get the maximum of first Piola-Kirchhoff (material) stress
residMax = maxval(abs(tract)) ! get the maximum of the residual
2019-01-13 02:34:03 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0 &
2019-01-13 14:23:37 +05:30
.and. prm%of_debug == of) then
2019-04-06 00:15:56 +05:30
stresLoc = maxloc(abs(P))
residLoc = maxloc(abs(tract))
write(6,'(1x,a)')' '
write(6,'(1x,a,1x,i2,1x,i4)')'RGC residual check ...',ip,el
write(6,'(1x,a15,1x,e15.8,1x,a7,i3,1x,a12,i2,i2)')'Max stress: ',stresMax, &
'@ grain',stresLoc(3),'in component',stresLoc(1),stresLoc(2)
write(6,'(1x,a15,1x,e15.8,1x,a7,i3,1x,a12,i2)')'Max residual: ',residMax, &
'@ iface',residLoc(1),'in direction',residLoc(2)
flush(6)
endif
2019-01-13 02:34:03 +05:30
#endif
homogenization_RGC_updateState = .false.
!--------------------------------------------------------------------------------------------------
! If convergence reached => done and happy
if (residMax < relTol_RGC*stresMax .or. residMax < absTol_RGC) then
homogenization_RGC_updateState = .true.
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0 &
2019-01-13 14:23:37 +05:30
.and. prm%of_debug == of) write(6,'(1x,a55,/)')'... done and happy'
flush(6)
2019-01-13 03:37:35 +05:30
#endif
!--------------------------------------------------------------------------------------------------
! compute/update the state for postResult, i.e., all energy densities computed by time-integration
2019-04-06 00:15:56 +05:30
do iGrain = 1,product(prm%Nconstituents)
do i = 1,3;do j = 1,3
2019-01-13 03:37:35 +05:30
stt%work(of) = stt%work(of) &
+ P(i,j,iGrain)*(F(i,j,iGrain) - F0(i,j,iGrain))/real(nGrain,pReal)
stt%penaltyEnergy(of) = stt%penaltyEnergy(of) &
+ R(i,j,iGrain)*(F(i,j,iGrain) - F0(i,j,iGrain))/real(nGrain,pReal)
2018-11-04 03:53:52 +05:30
enddo; enddo
enddo
2018-11-04 02:06:49 +05:30
2019-01-13 03:52:13 +05:30
dst%mismatch(1:3,of) = sum(NN,2)/real(nGrain,pReal)
2019-04-06 00:15:56 +05:30
dst%relaxationRate_avg(of) = sum(abs(drelax))/dt/real(3*nIntFaceTot,pReal)
2019-01-13 03:52:13 +05:30
dst%relaxationRate_max(of) = maxval(abs(drelax))/dt
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0 &
2019-01-13 14:23:37 +05:30
.and. prm%of_debug == of) then
2019-01-13 03:37:35 +05:30
write(6,'(1x,a30,1x,e15.8)') 'Constitutive work: ',stt%work(of)
2019-01-13 03:52:13 +05:30
write(6,'(1x,a30,3(1x,e15.8))')'Magnitude mismatch: ',dst%mismatch(1,of), &
dst%mismatch(2,of), &
dst%mismatch(3,of)
2019-01-13 03:37:35 +05:30
write(6,'(1x,a30,1x,e15.8)') 'Penalty energy: ', stt%penaltyEnergy(of)
2019-01-13 03:52:13 +05:30
write(6,'(1x,a30,1x,e15.8,/)') 'Volume discrepancy: ', dst%volumeDiscrepancy(of)
write(6,'(1x,a30,1x,e15.8)') 'Maximum relaxation rate: ', dst%relaxationRate_max(of)
write(6,'(1x,a30,1x,e15.8,/)') 'Average relaxation rate: ', dst%relaxationRate_avg(of)
flush(6)
endif
2019-01-13 03:37:35 +05:30
#endif
return
!--------------------------------------------------------------------------------------------------
! if residual blows-up => done but unhappy
elseif (residMax > relMax_RGC*stresMax .or. residMax > absMax_RGC) then ! try to restart when residual blows up exceeding maximum bound
homogenization_RGC_updateState = [.true.,.false.] ! with direct cut-back
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0 &
2019-01-13 14:23:37 +05:30
.and. prm%of_debug == of) write(6,'(1x,a,/)') '... broken'
flush(6)
2019-01-13 03:37:35 +05:30
#endif
return
else ! proceed with computing the Jacobian and state update
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0 &
2019-01-13 14:23:37 +05:30
.and. prm%of_debug == of) write(6,'(1x,a,/)') '... not yet done'
flush(6)
2019-01-13 03:37:35 +05:30
#endif
endif
!---------------------------------------------------------------------------------------------------
! construct the global Jacobian matrix for updating the global relaxation vector array when
! convergence is not yet reached ...
!--------------------------------------------------------------------------------------------------
! ... of the constitutive stress tangent, assembled from dPdF or material constitutive model "smatrix"
allocate(smatrix(3*nIntFaceTot,3*nIntFaceTot), source=0.0_pReal)
2019-04-06 00:15:56 +05:30
do iNum = 1,nIntFaceTot
2018-12-18 22:47:06 +05:30
faceID = interface1to4(iNum,param(instance)%Nconstituents) ! assembling of local dPdF into global Jacobian matrix
!--------------------------------------------------------------------------------------------------
! identify the left/bottom/back grain (-|N)
iGr3N = faceID(2:4) ! identifying the grain ID in local coordinate sytem
2018-12-18 22:47:06 +05:30
iGrN = grain3to1(iGr3N,param(instance)%Nconstituents) ! translate into global grain ID
2019-04-06 00:15:56 +05:30
intFaceN = getInterface(2*faceID(1),iGr3N) ! identifying the connecting interface in local coordinate system
2018-11-04 03:43:20 +05:30
normN = interfaceNormal(intFaceN,instance,of)
2019-04-06 00:15:56 +05:30
do iFace = 1,6
2018-12-18 22:47:06 +05:30
intFaceN = getInterface(iFace,iGr3N) ! identifying all interfaces that influence relaxation of the above interface
2018-11-04 03:43:20 +05:30
mornN = interfaceNormal(intFaceN,instance,of)
2018-11-04 12:26:27 +05:30
iMun = interface4to1(intFaceN,param(instance)%Nconstituents) ! translate the interfaces ID into local 4-dimensional index
if (iMun > 0) then ! get the corresponding tangent
2019-04-06 00:15:56 +05:30
do i=1,3; do j=1,3; do k=1,3; do l=1,3
2018-11-04 03:53:52 +05:30
smatrix(3*(iNum-1)+i,3*(iMun-1)+j) = smatrix(3*(iNum-1)+i,3*(iMun-1)+j) &
+ dPdF(i,k,j,l,iGrN)*normN(k)*mornN(l)
enddo;enddo;enddo;enddo
! projecting the material tangent dPdF into the interface
! to obtain the Jacobian matrix contribution of dPdF
endif
enddo
!--------------------------------------------------------------------------------------------------
! identify the right/up/front grain (+|P)
iGr3P = iGr3N
2019-04-06 00:15:56 +05:30
iGr3P(faceID(1)) = iGr3N(faceID(1))+1 ! identifying the grain ID in local coordinate sytem
2018-12-18 22:47:06 +05:30
iGrP = grain3to1(iGr3P,param(instance)%Nconstituents) ! translate into global grain ID
2019-04-06 00:15:56 +05:30
intFaceP = getInterface(2*faceID(1)-1,iGr3P) ! identifying the connecting interface in local coordinate system
2018-11-04 03:43:20 +05:30
normP = interfaceNormal(intFaceP,instance,of)
2019-04-06 00:15:56 +05:30
do iFace = 1,6
2018-12-18 22:47:06 +05:30
intFaceP = getInterface(iFace,iGr3P) ! identifying all interfaces that influence relaxation of the above interface
2018-11-04 03:43:20 +05:30
mornP = interfaceNormal(intFaceP,instance,of)
2018-12-18 22:47:06 +05:30
iMun = interface4to1(intFaceP,param(instance)%Nconstituents) ! translate the interfaces ID into local 4-dimensional index
2019-04-06 00:15:56 +05:30
if (iMun > 0) then ! get the corresponding tangent
do i=1,3; do j=1,3; do k=1,3; do l=1,3
2018-11-04 13:19:40 +05:30
smatrix(3*(iNum-1)+i,3*(iMun-1)+j) = smatrix(3*(iNum-1)+i,3*(iMun-1)+j) &
+ dPdF(i,k,j,l,iGrP)*normP(k)*mornP(l)
enddo;enddo;enddo;enddo
endif
enddo
enddo
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0) then
write(6,'(1x,a30)')'Jacobian matrix of stress'
2019-04-06 00:15:56 +05:30
do i = 1,3*nIntFaceTot
write(6,'(1x,100(e11.4,1x))')(smatrix(i,j), j = 1,3*nIntFaceTot)
enddo
write(6,*)' '
flush(6)
endif
2019-01-13 03:37:35 +05:30
#endif
!--------------------------------------------------------------------------------------------------
! ... of the stress penalty tangent (mismatch penalty and volume penalty, computed using numerical
! perturbation method) "pmatrix"
allocate(pmatrix(3*nIntFaceTot,3*nIntFaceTot), source=0.0_pReal)
allocate(p_relax(3*nIntFaceTot), source=0.0_pReal)
allocate(p_resid(3*nIntFaceTot), source=0.0_pReal)
2018-11-04 13:19:40 +05:30
2019-04-06 00:15:56 +05:30
do ipert = 1,3*nIntFaceTot
p_relax = relax
p_relax(ipert) = relax(ipert) + pPert_RGC ! perturb the relaxation vector
2019-01-13 03:37:35 +05:30
stt%relaxationVector(:,of) = p_relax
2018-12-18 22:47:06 +05:30
call grainDeformation(pF,avgF,instance,of) ! rain deformation from perturbed state
2019-01-13 14:03:47 +05:30
call stressPenalty(pR,DevNull, avgF,pF,ip,el,instance,of) ! stress penalty due to interface mismatch from perturbed state
call volumePenalty(pD,devNull(1,1), avgF,pF,nGrain,instance,of) ! stress penalty due to volume discrepancy from perturbed state
!--------------------------------------------------------------------------------------------------
! computing the global stress residual array from the perturbed state
p_resid = 0.0_pReal
2019-04-06 00:15:56 +05:30
do iNum = 1,nIntFaceTot
2018-11-04 12:41:35 +05:30
faceID = interface1to4(iNum,param(instance)%Nconstituents) ! identifying the interface ID in local coordinate system (4-dimensional index)
!--------------------------------------------------------------------------------------------------
! identify the left/bottom/back grain (-|N)
2018-12-18 22:47:06 +05:30
iGr3N = faceID(2:4) ! identify the grain ID in local coordinate system (3-dimensional index)
iGrN = grain3to1(iGr3N,param(instance)%Nconstituents) ! translate the local grain ID into global coordinate system (1-dimensional index)
2019-04-06 00:15:56 +05:30
intFaceN = getInterface(2*faceID(1),iGr3N) ! identify the interface ID of the grain
2018-11-04 03:43:20 +05:30
normN = interfaceNormal(intFaceN,instance,of)
!--------------------------------------------------------------------------------------------------
! identify the right/up/front grain (+|P)
iGr3P = iGr3N
2019-04-06 00:15:56 +05:30
iGr3P(faceID(1)) = iGr3N(faceID(1))+1 ! identify the grain ID in local coordinate system (3-dimensional index)
2018-12-18 22:47:06 +05:30
iGrP = grain3to1(iGr3P,param(instance)%Nconstituents) ! translate the local grain ID into global coordinate system (1-dimensional index)
2019-04-06 00:15:56 +05:30
intFaceP = getInterface(2*faceID(1)-1,iGr3P) ! identify the interface ID of the grain
2018-11-04 03:43:20 +05:30
normP = interfaceNormal(intFaceP,instance,of)
!--------------------------------------------------------------------------------------------------
! compute the residual stress (contribution of mismatch and volume penalties) from perturbed state
! at all interfaces
2019-04-06 00:15:56 +05:30
do i = 1,3; do j = 1,3
p_resid(i+3*(iNum-1)) = p_resid(i+3*(iNum-1)) + (pR(i,j,iGrP) - R(i,j,iGrP))*normP(j) &
+ (pR(i,j,iGrN) - R(i,j,iGrN))*normN(j) &
+ (pD(i,j,iGrP) - D(i,j,iGrP))*normP(j) &
+ (pD(i,j,iGrN) - D(i,j,iGrN))*normN(j)
enddo; enddo
enddo
pmatrix(:,ipert) = p_resid/pPert_RGC
enddo
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization), debug_levelExtensive) /= 0) then
write(6,'(1x,a30)')'Jacobian matrix of penalty'
2019-04-06 00:15:56 +05:30
do i = 1,3*nIntFaceTot
write(6,'(1x,100(e11.4,1x))')(pmatrix(i,j), j = 1,3*nIntFaceTot)
enddo
write(6,*)' '
flush(6)
endif
2019-01-13 03:37:35 +05:30
#endif
!--------------------------------------------------------------------------------------------------
! ... of the numerical viscosity traction "rmatrix"
allocate(rmatrix(3*nIntFaceTot,3*nIntFaceTot),source=0.0_pReal)
2019-04-06 00:15:56 +05:30
forall (i=1:3*nIntFaceTot) &
rmatrix(i,i) = viscModus_RGC*viscPower_RGC/(refRelaxRate_RGC*dt)* & ! tangent due to numerical viscosity traction appears
(abs(drelax(i))/(refRelaxRate_RGC*dt))**(viscPower_RGC - 1.0_pReal) ! only in the main diagonal term
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization), debug_levelExtensive) /= 0) then
write(6,'(1x,a30)')'Jacobian matrix of penalty'
2019-04-06 00:15:56 +05:30
do i = 1,3*nIntFaceTot
write(6,'(1x,100(e11.4,1x))')(rmatrix(i,j), j = 1,3*nIntFaceTot)
enddo
write(6,*)' '
flush(6)
endif
2019-01-13 03:37:35 +05:30
#endif
!--------------------------------------------------------------------------------------------------
! The overall Jacobian matrix summarizing contributions of smatrix, pmatrix, rmatrix
allocate(jmatrix(3*nIntFaceTot,3*nIntFaceTot)); jmatrix = smatrix + pmatrix + rmatrix
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization), debug_levelExtensive) /= 0) then
write(6,'(1x,a30)')'Jacobian matrix (total)'
2019-04-06 00:15:56 +05:30
do i = 1,3*nIntFaceTot
write(6,'(1x,100(e11.4,1x))')(jmatrix(i,j), j = 1,3*nIntFaceTot)
enddo
write(6,*)' '
flush(6)
endif
2019-01-13 03:37:35 +05:30
#endif
!--------------------------------------------------------------------------------------------------
! computing the update of the state variable (relaxation vectors) using the Jacobian matrix
2019-04-06 00:15:56 +05:30
allocate(jnverse(3*nIntFaceTot,3*nIntFaceTot),source=0.0_pReal)
call math_invert2(jnverse,error,jmatrix)
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization), debug_levelExtensive) /= 0) then
write(6,'(1x,a30)')'Jacobian inverse'
2019-04-06 00:15:56 +05:30
do i = 1,3*nIntFaceTot
write(6,'(1x,100(e11.4,1x))')(jnverse(i,j), j = 1,3*nIntFaceTot)
enddo
write(6,*)' '
flush(6)
endif
2019-01-13 03:37:35 +05:30
#endif
!--------------------------------------------------------------------------------------------------
! calculate the state update (global relaxation vectors) for the next Newton-Raphson iteration
drelax = 0.0_pReal
2019-04-06 00:15:56 +05:30
do i = 1,3*nIntFaceTot;do j = 1,3*nIntFaceTot
2018-11-04 13:19:40 +05:30
drelax(i) = drelax(i) - jnverse(i,j)*resid(j) ! Calculate the correction for the state variable
enddo; enddo
stt%relaxationVector(:,of) = relax + drelax ! Updateing the state variable for the next iteration
if (any(abs(drelax) > maxdRelax_RGC)) then ! Forcing cutback when the incremental change of relaxation vector becomes too large
homogenization_RGC_updateState = [.true.,.false.]
!$OMP CRITICAL (write2out)
write(6,'(1x,a,1x,i3,1x,a,1x,i3,1x,a)')'RGC_updateState: ip',ip,'| el',el,'enforces cutback'
write(6,'(1x,a,1x,e15.8)')'due to large relaxation change =',maxval(abs(drelax))
flush(6)
!$OMP END CRITICAL (write2out)
endif
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_homogenization, debug_levelExtensive) > 0) then
write(6,'(1x,a30)')'Returned state: '
2019-04-06 00:15:56 +05:30
do i = 1,size(stt%relaxationVector(:,of))
write(6,'(1x,2(e15.8,1x))') stt%relaxationVector(i,of)
enddo
write(6,*)' '
flush(6)
endif
2019-01-13 03:37:35 +05:30
#endif
end associate
2018-11-04 01:41:43 +05:30
contains
!--------------------------------------------------------------------------------------------------
!> @brief calculate stress-like penalty due to deformation mismatch
!--------------------------------------------------------------------------------------------------
2019-01-13 03:37:35 +05:30
subroutine stressPenalty(rPen,nMis,avgF,fDef,ip,el,instance,of)
2018-11-04 01:41:43 +05:30
use math, only: &
math_civita
use numerics, only: &
xSmoo_RGC
implicit none
2019-01-13 14:03:47 +05:30
real(pReal), dimension (:,:,:), intent(out) :: rPen !< stress-like penalty
real(pReal), dimension (:,:), intent(out) :: nMis !< total amount of mismatch
2019-01-13 03:37:35 +05:30
2019-01-13 14:03:47 +05:30
real(pReal), dimension (:,:,:), intent(in) :: fDef !< deformation gradients
real(pReal), dimension (3,3), intent(in) :: avgF !< initial effective stretch tensor
2019-04-06 00:15:56 +05:30
integer, intent(in) :: ip,el,instance,of
2019-01-13 03:37:35 +05:30
2019-04-06 00:15:56 +05:30
integer, dimension (4) :: intFace
integer, dimension (3) :: iGrain3,iGNghb3,nGDim
2018-11-04 01:41:43 +05:30
real(pReal), dimension (3,3) :: gDef,nDef
real(pReal), dimension (3) :: nVect,surfCorr
real(pReal), dimension (2) :: Gmoduli
2019-04-06 00:15:56 +05:30
integer :: iGrain,iGNghb,iFace,i,j,k,l
2018-11-04 01:41:43 +05:30
real(pReal) :: muGrain,muGNghb,nDefNorm,bgGrain,bgGNghb
real(pReal), parameter :: nDefToler = 1.0e-10_pReal
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2018-11-04 03:53:52 +05:30
logical :: debugActive
2019-01-13 03:37:35 +05:30
#endif
2018-11-04 03:53:52 +05:30
2018-11-04 01:41:43 +05:30
nGDim = param(instance)%Nconstituents
rPen = 0.0_pReal
nMis = 0.0_pReal
!--------------------------------------------------------------------------------------------------
! get the correction factor the modulus of penalty stress representing the evolution of area of
! the interfaces due to deformations
2018-11-04 03:43:20 +05:30
surfCorr = surfaceCorrection(avgF,instance,of)
2019-01-13 03:37:35 +05:30
2018-11-04 01:41:43 +05:30
associate(prm => param(instance))
2018-11-04 03:53:52 +05:30
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
debugActive = iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0 &
.and. prm%of_debug == of
2019-01-13 03:37:35 +05:30
2018-11-04 03:53:52 +05:30
if (debugActive) then
write(6,'(1x,a20,2(1x,i3))')'Correction factor: ',ip,el
write(6,*) surfCorr
2018-11-04 01:41:43 +05:30
endif
2019-01-13 03:37:35 +05:30
#endif
2018-11-04 01:41:43 +05:30
!--------------------------------------------------------------------------------------------------
! computing the mismatch and penalty stress tensor of all grains
2019-04-06 00:15:56 +05:30
grainLoop: do iGrain = 1,product(prm%Nconstituents)
2018-11-04 01:41:43 +05:30
Gmoduli = equivalentModuli(iGrain,ip,el)
2018-12-18 22:47:06 +05:30
muGrain = Gmoduli(1) ! collecting the equivalent shear modulus of grain
bgGrain = Gmoduli(2) ! and the lengthh of Burgers vector
iGrain3 = grain1to3(iGrain,prm%Nconstituents) ! get the grain ID in local 3-dimensional index (x,y,z)-position
2018-11-04 01:41:43 +05:30
2019-04-06 00:15:56 +05:30
interfaceLoop: do iFace = 1,6
2018-12-18 22:47:06 +05:30
intFace = getInterface(iFace,iGrain3) ! get the 4-dimensional index of the interface in local numbering system of the grain
2018-11-04 03:43:20 +05:30
nVect = interfaceNormal(intFace,instance,of)
2018-12-18 22:47:06 +05:30
iGNghb3 = iGrain3 ! identify the neighboring grain across the interface
2018-11-04 02:06:49 +05:30
iGNghb3(abs(intFace(1))) = iGNghb3(abs(intFace(1))) &
2019-04-06 00:18:20 +05:30
+ int(real(intFace(1),pReal)/real(abs(intFace(1)),pReal))
2018-11-04 02:06:49 +05:30
where(iGNghb3 < 1) iGNghb3 = nGDim
2019-04-06 00:15:56 +05:30
where(iGNghb3 >nGDim) iGNghb3 = 1
2018-12-18 22:47:06 +05:30
iGNghb = grain3to1(iGNghb3,prm%Nconstituents) ! get the ID of the neighboring grain
Gmoduli = equivalentModuli(iGNghb,ip,el) ! collect the shear modulus and Burgers vector of the neighbor
2018-11-04 01:41:43 +05:30
muGNghb = Gmoduli(1)
bgGNghb = Gmoduli(2)
2018-12-18 22:47:06 +05:30
gDef = 0.5_pReal*(fDef(1:3,1:3,iGNghb) - fDef(1:3,1:3,iGrain)) ! difference/jump in deformation gradeint across the neighbor
2018-11-04 01:41:43 +05:30
!--------------------------------------------------------------------------------------------------
! compute the mismatch tensor of all interfaces
nDefNorm = 0.0_pReal
nDef = 0.0_pReal
2019-04-06 00:15:56 +05:30
do i = 1,3; do j = 1,3
do k = 1,3; do l = 1,3
2018-12-18 22:47:06 +05:30
nDef(i,j) = nDef(i,j) - nVect(k)*gDef(i,l)*math_civita(j,k,l) ! compute the interface mismatch tensor from the jump of deformation gradient
2018-11-04 01:41:43 +05:30
enddo; enddo
2018-12-18 22:47:06 +05:30
nDefNorm = nDefNorm + nDef(i,j)**2.0_pReal ! compute the norm of the mismatch tensor
2018-11-04 01:41:43 +05:30
enddo; enddo
2018-12-18 22:47:06 +05:30
nDefNorm = max(nDefToler,sqrt(nDefNorm)) ! approximation to zero mismatch if mismatch is zero (singularity)
nMis(abs(intFace(1)),iGrain) = nMis(abs(intFace(1)),iGrain) + nDefNorm ! total amount of mismatch experienced by the grain (at all six interfaces)
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2018-11-04 13:19:40 +05:30
if (debugActive) then
write(6,'(1x,a20,i2,1x,a20,1x,i3)')'Mismatch to face: ',intFace(1),'neighbor grain: ',iGNghb
write(6,*) transpose(nDef)
write(6,'(1x,a20,e11.4)')'with magnitude: ',nDefNorm
endif
2019-01-13 03:37:35 +05:30
#endif
2018-11-04 01:41:43 +05:30
!--------------------------------------------------------------------------------------------------
! compute the stress penalty of all interfaces
2019-04-06 00:15:56 +05:30
do i = 1,3; do j = 1,3; do k = 1,3; do l = 1,3
2018-11-04 13:19:40 +05:30
rPen(i,j,iGrain) = rPen(i,j,iGrain) + 0.5_pReal*(muGrain*bgGrain + muGNghb*bgGNghb)*prm%xiAlpha &
*surfCorr(abs(intFace(1)))/prm%dAlpha(abs(intFace(1))) &
*cosh(prm%ciAlpha*nDefNorm) &
*0.5_pReal*nVect(l)*nDef(i,k)/nDefNorm*math_civita(k,l,j) &
*tanh(nDefNorm/xSmoo_RGC)
enddo; enddo;enddo; enddo
2019-01-13 03:37:35 +05:30
enddo interfaceLoop
#ifdef DEBUG
2018-11-04 03:53:52 +05:30
if (debugActive) then
2019-01-13 03:37:35 +05:30
write(6,'(1x,a20,i2)')'Penalty of grain: ',iGrain
write(6,*) transpose(rPen(1:3,1:3,iGrain))
2018-11-04 01:41:43 +05:30
endif
2019-01-13 03:37:35 +05:30
#endif
2018-11-04 01:41:43 +05:30
2019-01-13 03:37:35 +05:30
enddo grainLoop
2018-11-04 01:41:43 +05:30
end associate
end subroutine stressPenalty
!--------------------------------------------------------------------------------------------------
!> @brief calculate stress-like penalty due to volume discrepancy
!--------------------------------------------------------------------------------------------------
2019-01-13 03:37:35 +05:30
subroutine volumePenalty(vPen,vDiscrep,fAvg,fDef,nGrain,instance,of)
2018-11-04 01:41:43 +05:30
use math, only: &
math_det33, &
math_inv33
use numerics, only: &
maxVolDiscr_RGC,&
volDiscrMod_RGC,&
volDiscrPow_RGC
implicit none
2019-01-13 14:03:47 +05:30
real(pReal), dimension (:,:,:), intent(out) :: vPen ! stress-like penalty due to volume
real(pReal), intent(out) :: vDiscrep ! total volume discrepancy
2019-01-13 03:37:35 +05:30
2019-01-13 14:03:47 +05:30
real(pReal), dimension (:,:,:), intent(in) :: fDef ! deformation gradients
real(pReal), dimension (3,3), intent(in) :: fAvg ! overall deformation gradient
2019-04-06 00:15:56 +05:30
integer, intent(in) :: &
2019-01-13 03:37:35 +05:30
Ngrain, &
instance, &
of
2019-01-13 14:03:47 +05:30
real(pReal), dimension(size(vPen,3)) :: gVol
2019-04-06 00:15:56 +05:30
integer :: i
2018-11-04 03:53:52 +05:30
2018-11-04 01:41:43 +05:30
!--------------------------------------------------------------------------------------------------
! compute the volumes of grains and of cluster
2019-01-13 03:37:35 +05:30
vDiscrep = math_det33(fAvg) ! compute the volume of the cluster
2019-04-06 00:15:56 +05:30
do i = 1,nGrain
2019-01-13 03:37:35 +05:30
gVol(i) = math_det33(fDef(1:3,1:3,i)) ! compute the volume of individual grains
vDiscrep = vDiscrep - gVol(i)/real(nGrain,pReal) ! calculate the difference/dicrepancy between
! the volume of the cluster and the the total volume of grains
2018-11-04 01:41:43 +05:30
enddo
!--------------------------------------------------------------------------------------------------
! calculate the stress and penalty due to volume discrepancy
vPen = 0.0_pReal
2019-04-06 00:15:56 +05:30
do i = 1,nGrain
2019-01-13 03:37:35 +05:30
vPen(:,:,i) = -1.0_pReal/real(nGrain,pReal)*volDiscrMod_RGC*volDiscrPow_RGC/maxVolDiscr_RGC* &
2018-11-04 01:41:43 +05:30
sign((abs(vDiscrep)/maxVolDiscr_RGC)**(volDiscrPow_RGC - 1.0),vDiscrep)* &
2019-01-13 03:37:35 +05:30
gVol(i)*transpose(math_inv33(fDef(:,:,i)))
2018-11-04 01:41:43 +05:30
2019-01-13 03:37:35 +05:30
#ifdef DEBUG
2019-04-06 00:15:56 +05:30
if (iand(debug_level(debug_homogenization),debug_levelExtensive) /= 0 &
2019-01-13 03:37:35 +05:30
.and. param(instance)%of_debug == of) then
write(6,'(1x,a30,i2)')'Volume penalty of grain: ',i
write(6,*) transpose(vPen(:,:,i))
2018-11-04 01:41:43 +05:30
endif
2019-01-13 03:37:35 +05:30
#endif
2018-11-04 01:41:43 +05:30
enddo
2018-11-04 02:06:49 +05:30
2018-11-04 01:41:43 +05:30
end subroutine volumePenalty
2018-11-04 02:06:49 +05:30
2018-11-04 01:41:43 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief compute the correction factor accouted for surface evolution (area change) due to
! deformation
!--------------------------------------------------------------------------------------------------
2018-11-04 03:43:20 +05:30
function surfaceCorrection(avgF,instance,of)
2018-11-04 01:41:43 +05:30
use math, only: &
math_invert33
2018-11-04 01:41:43 +05:30
implicit none
real(pReal), dimension(3) :: surfaceCorrection
2019-04-06 00:15:56 +05:30
2018-11-04 01:41:43 +05:30
real(pReal), dimension(3,3), intent(in) :: avgF !< average F
2019-04-06 00:15:56 +05:30
integer, intent(in) :: &
2018-11-04 03:43:20 +05:30
instance, &
of
2018-11-04 01:41:43 +05:30
real(pReal), dimension(3,3) :: invC
real(pReal), dimension(3) :: nVect
real(pReal) :: detF
2019-04-06 00:15:56 +05:30
integer :: i,j,iBase
2018-11-04 01:41:43 +05:30
logical :: error
call math_invert33(matmul(transpose(avgF),avgF),invC,detF,error)
2018-11-04 01:41:43 +05:30
surfaceCorrection = 0.0_pReal
2019-04-06 00:15:56 +05:30
do iBase = 1,3
nVect = interfaceNormal([iBase,1,1,1],instance,of)
do i = 1,3; do j = 1,3
surfaceCorrection(iBase) = surfaceCorrection(iBase) + invC(i,j)*nVect(i)*nVect(j) ! compute the component of (the inverse of) the stretch in the direction of the normal
2018-11-04 01:41:43 +05:30
enddo; enddo
surfaceCorrection(iBase) = sqrt(surfaceCorrection(iBase))*detF ! get the surface correction factor (area contraction/enlargement)
2018-11-04 01:41:43 +05:30
enddo
end function surfaceCorrection
2018-11-04 02:06:49 +05:30
2018-11-04 01:41:43 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief compute the equivalent shear and bulk moduli from the elasticity tensor
!--------------------------------------------------------------------------------------------------
function equivalentModuli(grainID,ip,el)
use constitutive, only: &
constitutive_homogenizedC
implicit none
2018-11-04 12:41:35 +05:30
real(pReal), dimension(2) :: equivalentModuli
2019-04-06 00:15:56 +05:30
integer, intent(in) :: &
2018-11-04 01:41:43 +05:30
grainID,&
ip, & !< integration point number
el !< element number
2018-11-04 12:41:35 +05:30
real(pReal), dimension(6,6) :: elasTens
2018-11-04 01:41:43 +05:30
real(pReal) :: &
cEquiv_11, &
cEquiv_12, &
cEquiv_44
elasTens = constitutive_homogenizedC(grainID,ip,el)
!--------------------------------------------------------------------------------------------------
! compute the equivalent shear modulus after Turterltaub and Suiker, JMPS (2005)
cEquiv_11 = (elasTens(1,1) + elasTens(2,2) + elasTens(3,3))/3.0_pReal
cEquiv_12 = (elasTens(1,2) + elasTens(2,3) + elasTens(3,1) + &
elasTens(1,3) + elasTens(2,1) + elasTens(3,2))/6.0_pReal
cEquiv_44 = (elasTens(4,4) + elasTens(5,5) + elasTens(6,6))/3.0_pReal
equivalentModuli(1) = 0.2_pReal*(cEquiv_11 - cEquiv_12) + 0.6_pReal*cEquiv_44
!--------------------------------------------------------------------------------------------------
! obtain the length of Burgers vector (could be model dependend)
equivalentModuli(2) = 2.5e-10_pReal
end function equivalentModuli
!--------------------------------------------------------------------------------------------------
!> @brief calculating the grain deformation gradient (the same with
! homogenization_RGC_partitionDeformation, but used only for perturbation scheme)
2018-11-04 01:41:43 +05:30
!--------------------------------------------------------------------------------------------------
subroutine grainDeformation(F, avgF, instance, of)
2018-11-04 01:41:43 +05:30
implicit none
2019-04-06 00:15:56 +05:30
real(pReal), dimension(:,:,:), intent(out) :: F !< partioned F per grain
2019-04-06 00:15:56 +05:30
real(pReal), dimension(:,:), intent(in) :: avgF !< averaged F
integer, intent(in) :: &
instance, &
of
2019-01-13 03:37:35 +05:30
2019-04-06 00:15:56 +05:30
real(pReal), dimension(3) :: aVect,nVect
integer, dimension(4) :: intFace
integer, dimension(3) :: iGrain3
integer :: iGrain,iFace,i,j
2018-11-04 01:41:43 +05:30
!-------------------------------------------------------------------------------------------------
2018-11-04 03:43:20 +05:30
! compute the deformation gradient of individual grains due to relaxations
2019-01-13 03:37:35 +05:30
associate(prm => param(instance))
F = 0.0_pReal
2019-04-06 00:15:56 +05:30
do iGrain = 1,product(prm%Nconstituents)
2019-01-13 03:37:35 +05:30
iGrain3 = grain1to3(iGrain,prm%Nconstituents)
2019-04-06 00:15:56 +05:30
do iFace = 1,6
intFace = getInterface(iFace,iGrain3)
aVect = relaxationVector(intFace,instance,of)
nVect = interfaceNormal(intFace,instance,of)
2019-04-06 00:15:56 +05:30
forall (i=1:3,j=1:3) &
F(i,j,iGrain) = F(i,j,iGrain) + aVect(i)*nVect(j) ! effective relaxations
enddo
F(1:3,1:3,iGrain) = F(1:3,1:3,iGrain) + avgF ! relaxed deformation gradient
enddo
2018-11-04 03:43:20 +05:30
2019-01-13 03:37:35 +05:30
end associate
2018-11-04 01:41:43 +05:30
end subroutine grainDeformation
end function homogenization_RGC_updateState
!--------------------------------------------------------------------------------------------------
!> @brief derive average stress and stiffness from constituent quantities
!--------------------------------------------------------------------------------------------------
subroutine homogenization_RGC_averageStressAndItsTangent(avgP,dAvgPdAvgF,P,dPdF,instance)
implicit none
2019-04-06 00:15:56 +05:30
real(pReal), dimension (3,3), intent(out) :: avgP !< average stress at material point
real(pReal), dimension (3,3,3,3), intent(out) :: dAvgPdAvgF !< average stiffness at material point
2019-04-06 00:15:56 +05:30
real(pReal), dimension (:,:,:), intent(in) :: P !< partitioned stresses
real(pReal), dimension (:,:,:,:,:), intent(in) :: dPdF !< partitioned stiffnesses
integer, intent(in) :: instance
avgP = sum(P,3) /real(product(param(instance)%Nconstituents),pReal)
dAvgPdAvgF = sum(dPdF,5)/real(product(param(instance)%Nconstituents),pReal)
end subroutine homogenization_RGC_averageStressAndItsTangent
!--------------------------------------------------------------------------------------------------
!> @brief return array of homogenization results for post file inclusion
!--------------------------------------------------------------------------------------------------
2019-01-13 03:37:35 +05:30
pure function homogenization_RGC_postResults(instance,of) result(postResults)
implicit none
2019-04-06 00:15:56 +05:30
integer, intent(in) :: &
2019-01-13 03:37:35 +05:30
instance, &
of
2019-04-06 00:15:56 +05:30
integer :: &
2019-01-13 03:37:35 +05:30
o,c
real(pReal), dimension(sum(homogenization_RGC_sizePostResult(:,instance))) :: &
postResults
2019-01-13 03:52:13 +05:30
associate(stt => state(instance), dst => dependentState(instance), prm => param(instance))
2019-04-06 00:15:56 +05:30
c = 0
2019-01-13 03:37:35 +05:30
2019-04-06 00:15:56 +05:30
outputsLoop: do o = 1,size(prm%outputID)
select case(prm%outputID(o))
2019-01-13 03:37:35 +05:30
case (constitutivework_ID)
2019-01-13 03:37:35 +05:30
postResults(c+1) = stt%work(of)
2019-04-06 00:15:56 +05:30
c = c + 1
case (magnitudemismatch_ID)
2019-01-13 03:52:13 +05:30
postResults(c+1:c+3) = dst%mismatch(1:3,of)
2019-04-06 00:15:56 +05:30
c = c + 3
case (penaltyenergy_ID)
2019-01-13 03:37:35 +05:30
postResults(c+1) = stt%penaltyEnergy(of)
2019-04-06 00:15:56 +05:30
c = c + 1
case (volumediscrepancy_ID)
2019-01-13 03:52:13 +05:30
postResults(c+1) = dst%volumeDiscrepancy(of)
2019-04-06 00:15:56 +05:30
c = c + 1
case (averagerelaxrate_ID)
2019-01-13 03:52:13 +05:30
postResults(c+1) = dst%relaxationrate_avg(of)
2019-04-06 00:15:56 +05:30
c = c + 1
case (maximumrelaxrate_ID)
2019-01-13 03:52:13 +05:30
postResults(c+1) = dst%relaxationrate_max(of)
2019-04-06 00:15:56 +05:30
c = c + 1
end select
2019-01-13 03:37:35 +05:30
enddo outputsLoop
2019-01-13 03:37:35 +05:30
2018-08-30 14:16:30 +05:30
end associate
2019-01-13 03:37:35 +05:30
end function homogenization_RGC_postResults
!--------------------------------------------------------------------------------------------------
!> @brief collect relaxation vectors of an interface
!--------------------------------------------------------------------------------------------------
2018-11-04 12:03:57 +05:30
pure function relaxationVector(intFace,instance,of)
implicit none
2019-04-06 00:15:56 +05:30
real(pReal), dimension (3) :: relaxationVector
integer, intent(in) :: instance,of
integer, dimension(4), intent(in) :: intFace !< set of interface ID in 4D array (normal and position)
integer :: iNum
2019-01-13 03:37:35 +05:30
2019-04-06 00:15:56 +05:30
!--------------------------------------------------------------------------------------------------
! collect the interface relaxation vector from the global state array
2019-01-13 03:37:35 +05:30
2018-11-04 12:26:27 +05:30
iNum = interface4to1(intFace,param(instance)%Nconstituents) ! identify the position of the interface in global state array
2019-04-06 00:15:56 +05:30
if (iNum > 0) then
2019-01-13 03:37:35 +05:30
relaxationVector = state(instance)%relaxationVector((3*iNum-2):(3*iNum),of)
else
relaxationVector = 0.0_pReal
endif
end function relaxationVector
!--------------------------------------------------------------------------------------------------
!> @brief identify the normal of an interface
!--------------------------------------------------------------------------------------------------
2018-11-04 12:03:57 +05:30
pure function interfaceNormal(intFace,instance,of)
implicit none
2019-04-06 00:15:56 +05:30
real(pReal), dimension(3) :: interfaceNormal
integer, dimension(4), intent(in) :: intFace !< interface ID in 4D array (normal and position)
integer, intent(in) :: &
2018-11-04 03:43:20 +05:30
instance, &
of
2019-04-06 00:15:56 +05:30
integer :: nPos
!--------------------------------------------------------------------------------------------------
! get the normal of the interface, identified from the value of intFace(1)
interfaceNormal = 0.0_pReal
nPos = abs(intFace(1)) ! identify the position of the interface in global state array
2018-08-30 21:21:31 +05:30
interfaceNormal(nPos) = real(intFace(1)/abs(intFace(1)),pReal) ! get the normal vector w.r.t. cluster axis
2019-04-06 00:15:56 +05:30
interfaceNormal = matmul(dependentState(instance)%orientation(1:3,1:3,of),interfaceNormal) ! map the normal vector into sample coordinate system (basis)
end function interfaceNormal
!--------------------------------------------------------------------------------------------------
!> @brief collect six faces of a grain in 4D (normal and position)
!--------------------------------------------------------------------------------------------------
2018-11-04 12:03:57 +05:30
pure function getInterface(iFace,iGrain3)
implicit none
2019-04-06 00:15:56 +05:30
integer, dimension(4) :: getInterface
integer, dimension(3), intent(in) :: iGrain3 !< grain ID in 3D array
integer, intent(in) :: iFace !< face index (1..6) mapped like (-e1,-e2,-e3,+e1,+e2,+e3) or iDir = (-1,-2,-3,1,2,3)
integer :: iDir
!* Direction of interface normal
2019-04-06 00:15:56 +05:30
iDir = (int(real(iFace-1,pReal)/2.0_pReal)+1)*(-1)**iFace
getInterface(1) = iDir
!--------------------------------------------------------------------------------------------------
! identify the interface position by the direction of its normal
getInterface(2:4) = iGrain3
2019-04-06 00:15:56 +05:30
if (iDir < 0) getInterface(1-iDir) = getInterface(1-iDir)-1 ! to have a correlation with coordinate/position in real space
end function getInterface
2018-11-04 12:41:35 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief map grain ID from in 1D (global array) to in 3D (local position)
!--------------------------------------------------------------------------------------------------
2018-11-04 12:41:35 +05:30
pure function grain1to3(grain1,nGDim)
implicit none
2019-04-06 00:15:56 +05:30
integer, dimension(3) :: grain1to3
integer, intent(in) :: grain1 !< grain ID in 1D array
integer, dimension(3), intent(in) :: nGDim
2019-04-06 00:15:56 +05:30
grain1to3 = 1 + [mod((grain1-1),nGDim(1)), &
mod((grain1-1)/nGDim(1),nGDim(2)), &
(grain1-1)/(nGDim(1)*nGDim(2))]
end function grain1to3
!--------------------------------------------------------------------------------------------------
!> @brief map grain ID from in 3D (local position) to in 1D (global array)
!--------------------------------------------------------------------------------------------------
2019-04-06 00:15:56 +05:30
integer pure function grain3to1(grain3,nGDim)
implicit none
2019-04-06 00:15:56 +05:30
integer, dimension(3), intent(in) :: grain3 !< grain ID in 3D array (pos.x,pos.y,pos.z)
integer, dimension(3), intent(in) :: nGDim
2018-11-04 12:41:35 +05:30
grain3to1 = grain3(1) &
2019-04-06 00:15:56 +05:30
+ nGDim(1)*(grain3(2)-1) &
+ nGDim(1)*nGDim(2)*(grain3(3)-1)
end function grain3to1
!--------------------------------------------------------------------------------------------------
!> @brief maps interface ID from 4D (normal and local position) into 1D (global array)
!--------------------------------------------------------------------------------------------------
2019-04-06 00:15:56 +05:30
integer pure function interface4to1(iFace4D, nGDim)
implicit none
2019-04-06 00:15:56 +05:30
integer, dimension(4), intent(in) :: iFace4D !< interface ID in 4D array (n.dir,pos.x,pos.y,pos.z)
integer, dimension(3), intent(in) :: nGDim
2018-11-04 13:49:24 +05:30
select case(abs(iFace4D(1)))
2019-04-06 00:15:56 +05:30
case(1)
if ((iFace4D(2) == 0) .or. (iFace4D(2) == nGDim(1))) then
interface4to1 = 0
2018-11-04 13:49:24 +05:30
else
2019-04-06 00:15:56 +05:30
interface4to1 = iFace4D(3) + nGDim(2)*(iFace4D(4)-1) &
+ nGDim(2)*nGDim(3)*(iFace4D(2)-1)
2018-11-04 13:49:24 +05:30
endif
2019-04-06 00:15:56 +05:30
case(2)
if ((iFace4D(3) == 0) .or. (iFace4D(3) == nGDim(2))) then
interface4to1 = 0
2018-11-04 13:49:24 +05:30
else
2019-04-06 00:15:56 +05:30
interface4to1 = iFace4D(4) + nGDim(3)*(iFace4D(2)-1) &
+ nGDim(3)*nGDim(1)*(iFace4D(3)-1) &
+ (nGDim(1)-1)*nGDim(2)*nGDim(3) ! total number of interfaces normal //e1
2018-11-04 13:49:24 +05:30
endif
2019-04-06 00:15:56 +05:30
case(3)
if ((iFace4D(4) == 0) .or. (iFace4D(4) == nGDim(3))) then
interface4to1 = 0
2018-11-04 13:49:24 +05:30
else
2019-04-06 00:15:56 +05:30
interface4to1 = iFace4D(2) + nGDim(1)*(iFace4D(3)-1) &
+ nGDim(1)*nGDim(2)*(iFace4D(4)-1) &
+ (nGDim(1)-1)*nGDim(2)*nGDim(3) & ! total number of interfaces normal //e1
+ nGDim(1)*(nGDim(2)-1)*nGDim(3) ! total number of interfaces normal //e2
2018-11-04 13:49:24 +05:30
endif
case default
2019-04-06 00:15:56 +05:30
interface4to1 = -1
2018-11-04 13:49:24 +05:30
end select
end function interface4to1
!--------------------------------------------------------------------------------------------------
!> @brief maps interface ID from 1D (global array) into 4D (normal and local position)
!--------------------------------------------------------------------------------------------------
2018-11-04 12:41:35 +05:30
pure function interface1to4(iFace1D, nGDim)
implicit none
2019-04-06 00:15:56 +05:30
integer, dimension(4) :: interface1to4
integer, intent(in) :: iFace1D !< interface ID in 1D array
integer, dimension(3), intent(in) :: nGDim
integer, dimension(3) :: nIntFace
!--------------------------------------------------------------------------------------------------
! compute the total number of interfaces, which ...
2019-04-06 00:15:56 +05:30
nIntFace = [(nGDim(1)-1)*nGDim(2)*nGDim(3), & ! ... normal //e1
nGDim(1)*(nGDim(2)-1)*nGDim(3), & ! ... normal //e2
nGDim(1)*nGDim(2)*(nGDim(3)-1)] ! ... normal //e3
!--------------------------------------------------------------------------------------------------
! get the corresponding interface ID in 4D (normal and local position)
if (iFace1D > 0 .and. iFace1D <= nIntFace(1)) then ! interface with normal //e1
2019-04-06 00:15:56 +05:30
interface1to4(1) = 1
interface1to4(3) = mod((iFace1D-1),nGDim(2))+1
interface1to4(4) = mod(int(real(iFace1D-1,pReal)/real(nGDim(2),pReal)),nGDim(3))+1
interface1to4(2) = int(real(iFace1D-1,pReal)/real(nGDim(2),pReal)/real(nGDim(3),pReal))+1
elseif (iFace1D > nIntFace(1) .and. iFace1D <= (nIntFace(2) + nIntFace(1))) then ! interface with normal //e2
2019-04-06 00:15:56 +05:30
interface1to4(1) = 2
interface1to4(4) = mod((iFace1D-nIntFace(1)-1),nGDim(3))+1
interface1to4(2) = mod(int(real(iFace1D-nIntFace(1)-1,pReal)/real(nGDim(3),pReal)),nGDim(1))+1
interface1to4(3) = int(real(iFace1D-nIntFace(1)-1,pReal)/real(nGDim(3),pReal)/real(nGDim(1),pReal))+1
elseif (iFace1D > nIntFace(2) + nIntFace(1) .and. iFace1D <= (nIntFace(3) + nIntFace(2) + nIntFace(1))) then ! interface with normal //e3
2019-04-06 00:15:56 +05:30
interface1to4(1) = 3
interface1to4(2) = mod((iFace1D-nIntFace(2)-nIntFace(1)-1),nGDim(1))+1
interface1to4(3) = mod(int(real(iFace1D-nIntFace(2)-nIntFace(1)-1,pReal)/real(nGDim(1),pReal)),nGDim(2))+1
interface1to4(4) = int(real(iFace1D-nIntFace(2)-nIntFace(1)-1,pReal)/real(nGDim(1),pReal)/real(nGDim(2),pReal))+1
endif
end function interface1to4
2019-04-06 00:28:56 +05:30
end module homogenization_mech_RGC