1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
17 integer,
dimension(:),
allocatable :: &
22 real(pReal),
dimension(:),
allocatable :: &
27 character(len=pStringLen),
allocatable,
dimension(:) :: &
32 real(pReal),
pointer,
dimension(:) :: &
35 real(pReal),
pointer,
dimension(:,:) :: &
39 type :: trgcdependentstate
40 real(pReal),
allocatable,
dimension(:) :: &
44 real(pReal),
allocatable,
dimension(:,:) :: &
46 real(pReal),
allocatable,
dimension(:,:,:) :: &
48 end type trgcdependentstate
52 atol, & !< absolute tolerance of RGC residuum
53 rtol, & !< relative tolerance of RGC residuum
54 absMax, & !< absolute maximum of RGC residuum
55 relMax, & !< relative maximum of RGC residuum
56 pPert, & !< perturbation for computing RGC penalty tangent
57 xSmoo, & !< RGC penalty smoothing parameter (hyperbolic tangent)
58 viscPower, & !< power (sensitivity rate) of numerical viscosity in RGC scheme, Default 1.0e0: Newton viscosity (linear model)
59 viscModus, & !< stress modulus of RGC numerical viscosity, Default 0.0e0: No viscosity is applied
60 refRelaxRate, & !< reference relaxation rate in RGC viscosity
61 maxdRelax, & !< threshold of maximum relaxation vector increment (if exceed this then cutback)
62 maxVolDiscr, & !< threshold of maximum volume discrepancy allowed
63 volDiscrMod, & !< stiffness of RGC volume discrepancy (zero = without volume discrepancy constraint)
65 end type tnumerics_rgc
67 type(tparameters),
dimension(:),
allocatable :: &
69 type(tRGCstate),
dimension(:),
allocatable :: &
72 type(tRGCdependentState),
dimension(:),
allocatable :: &
74 type(tNumerics_RGC) :: &
82 module subroutine mech_rgc_init
88 sizeState, nIntFaceTot
90 write(6,
'(/,a)')
' <<<+- homogenization_'//homogenization_rgc_label//
' init -+>>>';
flush(6)
92 write(6,
'(/,a)')
' Tjahjanto et al., International Journal of Material Forming 2(1):939–942, 2009'
93 write(6,
'(a)')
' https://doi.org/10.1007/s12289-009-0619-1'
95 write(6,
'(/,a)')
' Tjahjanto et al., Modelling and Simulation in Materials Science and Engineering 18:015006, 2010'
96 write(6,
'(a)')
' https://doi.org/10.1088/0965-0393/18/1/015006'
98 ninstance = count(homogenization_type == homogenization_rgc_id)
99 if (iand(debug_level(debug_homogenization),debug_levelbasic) /= 0) &
100 write(6,
'(a16,1x,i5,/)')
'# instances:',ninstance
102 allocate(param(ninstance))
103 allocate(state(ninstance))
104 allocate(state0(ninstance))
105 allocate(dependentstate(ninstance))
107 num%atol = config_numerics%getFloat(
'atol_rgc', defaultval=1.0e+4_preal)
108 num%rtol = config_numerics%getFloat(
'rtol_rgc', defaultval=1.0e-3_preal)
109 num%absMax = config_numerics%getFloat(
'amax_rgc', defaultval=1.0e+10_preal)
110 num%relMax = config_numerics%getFloat(
'rmax_rgc', defaultval=1.0e+2_preal)
111 num%pPert = config_numerics%getFloat(
'perturbpenalty_rgc', defaultval=1.0e-7_preal)
112 num%xSmoo = config_numerics%getFloat(
'relvantmismatch_rgc', defaultval=1.0e-5_preal)
113 num%viscPower = config_numerics%getFloat(
'viscositypower_rgc', defaultval=1.0e+0_preal)
114 num%viscModus = config_numerics%getFloat(
'viscositymodulus_rgc', defaultval=0.0e+0_preal)
115 num%refRelaxRate = config_numerics%getFloat(
'refrelaxationrate_rgc',defaultval=1.0e-3_preal)
116 num%maxdRelax = config_numerics%getFloat(
'maxrelaxationrate_rgc',defaultval=1.0e+0_preal)
117 num%maxVolDiscr = config_numerics%getFloat(
'maxvoldiscrepancy_rgc',defaultval=1.0e-5_preal)
118 num%volDiscrMod = config_numerics%getFloat(
'voldiscrepancymod_rgc',defaultval=1.0e+12_preal)
119 num%volDiscrPow = config_numerics%getFloat(
'dicrepancypower_rgc', defaultval=5.0_preal)
121 if (num%atol <= 0.0_preal)
call io_error(301,ext_msg=
'absTol_RGC')
122 if (num%rtol <= 0.0_preal)
call io_error(301,ext_msg=
'relTol_RGC')
123 if (num%absMax <= 0.0_preal)
call io_error(301,ext_msg=
'absMax_RGC')
124 if (num%relMax <= 0.0_preal)
call io_error(301,ext_msg=
'relMax_RGC')
125 if (num%pPert <= 0.0_preal)
call io_error(301,ext_msg=
'pPert_RGC')
126 if (num%xSmoo <= 0.0_preal)
call io_error(301,ext_msg=
'xSmoo_RGC')
127 if (num%viscPower < 0.0_preal)
call io_error(301,ext_msg=
'viscPower_RGC')
128 if (num%viscModus < 0.0_preal)
call io_error(301,ext_msg=
'viscModus_RGC')
129 if (num%refRelaxRate <= 0.0_preal)
call io_error(301,ext_msg=
'refRelaxRate_RGC')
130 if (num%maxdRelax <= 0.0_preal)
call io_error(301,ext_msg=
'maxdRelax_RGC')
131 if (num%maxVolDiscr <= 0.0_preal)
call io_error(301,ext_msg=
'maxVolDiscr_RGC')
132 if (num%volDiscrMod < 0.0_preal)
call io_error(301,ext_msg=
'volDiscrMod_RGC')
133 if (num%volDiscrPow <= 0.0_preal)
call io_error(301,ext_msg=
'volDiscrPw_RGC')
135 do h = 1,
size(homogenization_type)
136 if (homogenization_type(h) /= homogenization_rgc_id) cycle
137 associate(prm => param(homogenization_typeinstance(h)), &
138 stt => state(homogenization_typeinstance(h)), &
139 st0 => state0(homogenization_typeinstance(h)), &
140 dst => dependentstate(homogenization_typeinstance(h)), &
141 config => config_homogenization(h))
149 prm%output =
config%getStrings(
'(output)',defaultval=emptystringarray)
151 prm%Nconstituents =
config%getInts(
'clustersize',requiredsize=3)
152 if (homogenization_ngrains(h) /= product(prm%Nconstituents)) &
153 call io_error(211,ext_msg=
'clustersize ('//homogenization_rgc_label//
')')
155 prm%xiAlpha =
config%getFloat(
'scalingparameter')
156 prm%ciAlpha =
config%getFloat(
'overproportionality')
158 prm%dAlpha =
config%getFloats(
'grainsize', requiredsize=3)
159 prm%angles =
config%getFloats(
'clusterorientation',requiredsize=3)
161 nofmyhomog = count(material_homogenizationat == h)
162 nintfacetot = 3*( (prm%Nconstituents(1)-1)*prm%Nconstituents(2)*prm%Nconstituents(3) &
163 + prm%Nconstituents(1)*(prm%Nconstituents(2)-1)*prm%Nconstituents(3) &
164 + prm%Nconstituents(1)*prm%Nconstituents(2)*(prm%Nconstituents(3)-1))
165 sizestate = nintfacetot &
166 +
size([
'avg constitutive work ',
'average penalty energy'])
168 homogstate(h)%sizeState = sizestate
169 allocate(homogstate(h)%state0 (sizestate,nofmyhomog), source=0.0_preal)
170 allocate(homogstate(h)%subState0(sizestate,nofmyhomog), source=0.0_preal)
171 allocate(homogstate(h)%state (sizestate,nofmyhomog), source=0.0_preal)
173 stt%relaxationVector => homogstate(h)%state(1:nintfacetot,:)
174 st0%relaxationVector => homogstate(h)%state0(1:nintfacetot,:)
175 stt%work => homogstate(h)%state(nintfacetot+1,:)
176 stt%penaltyEnergy => homogstate(h)%state(nintfacetot+2,:)
178 allocate(dst%volumeDiscrepancy( nofmyhomog))
179 allocate(dst%relaxationRate_avg( nofmyhomog))
180 allocate(dst%relaxationRate_max( nofmyhomog))
181 allocate(dst%mismatch( 3,nofmyhomog))
185 dependentstate(homogenization_typeinstance(h))%orientation = spread(
eu2om(prm%angles*inrad),3,nofmyhomog)
192 end subroutine mech_rgc_init
198 module subroutine mech_rgc_partitiondeformation(f,avgf,instance,of)
200 real(pReal),
dimension (:,:,:),
intent(out) :: F
202 real(pReal),
dimension (3,3),
intent(in) :: avgF
203 integer,
intent(in) :: &
207 real(pReal),
dimension(3) :: aVect,nVect
208 integer,
dimension(4) :: intFace
209 integer,
dimension(3) :: iGrain3
210 integer :: iGrain,iFace,i,j
212 associate(prm => param(instance))
217 do igrain = 1,product(prm%Nconstituents)
218 igrain3 = grain1to3(igrain,prm%Nconstituents)
220 intface = getinterface(iface,igrain3)
221 avect = relaxationvector(intface,instance,of)
222 nvect = interfacenormal(intface,instance,of)
223 forall (i=1:3,j=1:3) &
224 f(i,j,igrain) = f(i,j,igrain) + avect(i)*nvect(j)
226 f(1:3,1:3,igrain) = f(1:3,1:3,igrain) + avgf
228 # 234 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
233 end subroutine mech_rgc_partitiondeformation
240 module procedure mech_rgc_updatestate
242 integer,
dimension(4) :: intFaceN,intFaceP,faceID
243 integer,
dimension(3) :: nGDim,iGr3N,iGr3P
244 integer :: instance,iNum,i,j,nIntFaceTot,iGrN,iGrP,iMun,iFace,k,l,ipert,iGrain,nGrain, of
245 real(pReal),
dimension(3,3,size(P,3)) :: R,pF,pR,D,pD
246 real(pReal),
dimension(3,size(P,3)) :: NN,devNull
247 real(pReal),
dimension(3) :: normP,normN,mornP,mornN
248 real(pReal) :: residMax,stresMax
250 real(pReal),
dimension(:,:),
allocatable :: tract,jmatrix,jnverse,smatrix,pmatrix,rmatrix
251 real(pReal),
dimension(:),
allocatable :: resid,relax,p_relax,p_resid,drelax
257 zerotimestep:
if(deq0(dt))
then
258 mech_rgc_updatestate = .true.
262 instance = homogenization_typeinstance(material_homogenizationat(el))
263 of = material_homogenizationmemberat(ip,el)
265 associate(stt => state(instance), st0 => state0(instance), dst => dependentstate(instance), prm => param(instance))
269 ngdim = prm%Nconstituents
270 ngrain = product(ngdim)
271 nintfacetot = (ngdim(1)-1)*ngdim(2)*ngdim(3) &
272 + ngdim(1)*(ngdim(2)-1)*ngdim(3) &
273 + ngdim(1)*ngdim(2)*(ngdim(3)-1)
277 allocate(resid(3*nintfacetot), source=0.0_preal)
278 allocate(tract(nintfacetot,3), source=0.0_preal)
279 relax = stt%relaxationVector(:,of)
280 drelax = stt%relaxationVector(:,of) - st0%relaxationVector(:,of)
282 # 296 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
286 call stresspenalty(r,nn,avgf,f,ip,el,instance,of)
290 call volumepenalty(d,dst%volumeDiscrepancy(of),avgf,f,ngrain,instance,of)
292 # 320 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
296 do inum = 1,nintfacetot
297 faceid = interface1to4(inum,param(instance)%Nconstituents)
302 igrn = grain3to1(igr3n,param(instance)%Nconstituents)
303 intfacen = getinterface(2*faceid(1),igr3n)
304 normn = interfacenormal(intfacen,instance,of)
309 igr3p(faceid(1)) = igr3n(faceid(1))+1
310 igrp = grain3to1(igr3p,param(instance)%Nconstituents)
311 intfacep = getinterface(2*faceid(1)-1,igr3p)
312 normp = interfacenormal(intfacep,instance,of)
317 tract(inum,i) = sign(num%viscModus*(abs(drelax(i+3*(inum-1)))/(num%refRelaxRate*dt))**num%viscPower, &
318 drelax(i+3*(inum-1)))
320 tract(inum,i) = tract(inum,i) + (p(i,j,igrp) + r(i,j,igrp) + d(i,j,igrp))*normp(j) &
321 + (p(i,j,igrn) + r(i,j,igrn) + d(i,j,igrn))*normn(j)
322 resid(i+3*(inum-1)) = tract(inum,i)
337 stresmax = maxval(abs(p))
338 residmax = maxval(abs(tract))
340 # 380 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
342 mech_rgc_updatestate = .false.
346 if (residmax < num%rtol*stresmax .or. residmax < num%atol)
then
347 mech_rgc_updatestate = .true.
355 do igrain = 1,product(prm%Nconstituents)
356 do i = 1,3;
do j = 1,3
357 stt%work(of) = stt%work(of) &
358 + p(i,j,igrain)*(f(i,j,igrain) - f0(i,j,igrain))/real(ngrain,preal)
359 stt%penaltyEnergy(of) = stt%penaltyEnergy(of) &
360 + r(i,j,igrain)*(f(i,j,igrain) - f0(i,j,igrain))/real(ngrain,preal)
364 dst%mismatch(1:3,of) = sum(nn,2)/real(ngrain,preal)
365 dst%relaxationRate_avg(of) = sum(abs(drelax))/dt/real(3*nintfacetot,preal)
366 dst%relaxationRate_max(of) = maxval(abs(drelax))/dt
368 # 420 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
374 elseif (residmax > num%relMax*stresmax .or. residmax > num%absMax)
then
375 mech_rgc_updatestate = [.true.,.false.]
398 allocate(smatrix(3*nintfacetot,3*nintfacetot), source=0.0_preal)
399 do inum = 1,nintfacetot
400 faceid = interface1to4(inum,param(instance)%Nconstituents)
405 igrn = grain3to1(igr3n,param(instance)%Nconstituents)
406 intfacen = getinterface(2*faceid(1),igr3n)
407 normn = interfacenormal(intfacen,instance,of)
409 intfacen = getinterface(iface,igr3n)
410 mornn = interfacenormal(intfacen,instance,of)
411 imun = interface4to1(intfacen,param(instance)%Nconstituents)
413 do i=1,3;
do j=1,3;
do k=1,3;
do l=1,3
414 smatrix(3*(inum-1)+i,3*(imun-1)+j) = smatrix(3*(inum-1)+i,3*(imun-1)+j) &
415 + dpdf(i,k,j,l,igrn)*normn(k)*mornn(l)
416 enddo;enddo;enddo;
enddo
425 igr3p(faceid(1)) = igr3n(faceid(1))+1
426 igrp = grain3to1(igr3p,param(instance)%Nconstituents)
427 intfacep = getinterface(2*faceid(1)-1,igr3p)
428 normp = interfacenormal(intfacep,instance,of)
430 intfacep = getinterface(iface,igr3p)
431 mornp = interfacenormal(intfacep,instance,of)
432 imun = interface4to1(intfacep,param(instance)%Nconstituents)
434 do i=1,3;
do j=1,3;
do k=1,3;
do l=1,3
435 smatrix(3*(inum-1)+i,3*(imun-1)+j) = smatrix(3*(inum-1)+i,3*(imun-1)+j) &
436 + dpdf(i,k,j,l,igrp)*normp(k)*mornp(l)
437 enddo;enddo;enddo;
enddo
442 # 503 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
447 allocate(pmatrix(3*nintfacetot,3*nintfacetot), source=0.0_preal)
448 allocate(p_relax(3*nintfacetot), source=0.0_preal)
449 allocate(p_resid(3*nintfacetot), source=0.0_preal)
451 do ipert = 1,3*nintfacetot
453 p_relax(ipert) = relax(ipert) + num%pPert
454 stt%relaxationVector(:,of) = p_relax
455 call graindeformation(pf,avgf,instance,of)
456 call stresspenalty(pr,devnull, avgf,pf,ip,el,instance,of)
457 call volumepenalty(pd,devnull(1,1), avgf,pf,ngrain,instance,of)
462 do inum = 1,nintfacetot
463 faceid = interface1to4(inum,param(instance)%Nconstituents)
468 igrn = grain3to1(igr3n,param(instance)%Nconstituents)
469 intfacen = getinterface(2*faceid(1),igr3n)
470 normn = interfacenormal(intfacen,instance,of)
475 igr3p(faceid(1)) = igr3n(faceid(1))+1
476 igrp = grain3to1(igr3p,param(instance)%Nconstituents)
477 intfacep = getinterface(2*faceid(1)-1,igr3p)
478 normp = interfacenormal(intfacep,instance,of)
483 do i = 1,3;
do j = 1,3
484 p_resid(i+3*(inum-1)) = p_resid(i+3*(inum-1)) + (pr(i,j,igrp) - r(i,j,igrp))*normp(j) &
485 + (pr(i,j,igrn) - r(i,j,igrn))*normn(j) &
486 + (pd(i,j,igrp) - d(i,j,igrp))*normp(j) &
487 + (pd(i,j,igrn) - d(i,j,igrn))*normn(j)
490 pmatrix(:,ipert) = p_resid/num%pPert
493 # 563 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
497 allocate(rmatrix(3*nintfacetot,3*nintfacetot),source=0.0_preal)
499 rmatrix(i,i) = num%viscModus*num%viscPower/(num%refRelaxRate*dt)* &
500 (abs(drelax(i))/(num%refRelaxRate*dt))**(num%viscPower - 1.0_preal)
503 # 582 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
507 allocate(jmatrix(3*nintfacetot,3*nintfacetot)); jmatrix = smatrix + pmatrix + rmatrix
509 # 597 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
513 allocate(jnverse(3*nintfacetot,3*nintfacetot),source=0.0_preal)
514 call math_invert(jnverse,error,jmatrix)
516 # 613 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
521 do i = 1,3*nintfacetot;
do j = 1,3*nintfacetot
522 drelax(i) = drelax(i) - jnverse(i,j)*resid(j)
524 stt%relaxationVector(:,of) = relax + drelax
525 if (any(abs(drelax) > num%maxdRelax))
then
526 mech_rgc_updatestate = [.true.,.false.]
528 write(6,
'(1x,a,1x,i3,1x,a,1x,i3,1x,a)')
'RGC_updateState: ip',ip,
'| el',el,
'enforces cutback'
529 write(6,
'(1x,a,1x,e15.8)')
'due to large relaxation change =',maxval(abs(drelax))
534 # 640 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
542 subroutine stresspenalty(rPen,nMis,avgF,fDef,ip,el,instance,of)
544 real(pReal),
dimension (:,:,:),
intent(out) :: rPen
545 real(pReal),
dimension (:,:),
intent(out) :: nMis
547 real(pReal),
dimension (:,:,:),
intent(in) :: fDef
548 real(pReal),
dimension (3,3),
intent(in) :: avgF
549 integer,
intent(in) :: ip,el,instance,of
551 integer,
dimension (4) :: intFace
552 integer,
dimension (3) :: iGrain3,iGNghb3,nGDim
553 real(pReal),
dimension (3,3) :: gDef,nDef
554 real(pReal),
dimension (3) :: nVect,surfCorr
555 real(pReal),
dimension (2) :: Gmoduli
556 integer :: iGrain,iGNghb,iFace,i,j,k,l
557 real(pReal) :: muGrain,muGNghb,nDefNorm,bgGrain,bgGNghb
558 real(pReal),
parameter :: nDefToler = 1.0e-10_preal
563 ngdim = param(instance)%Nconstituents
571 surfcorr = surfacecorrection(avgf,instance,of)
573 associate(prm => param(instance))
575 # 688 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization_mech_RGC.f90"
579 grainloop:
do igrain = 1,product(prm%Nconstituents)
580 gmoduli = equivalentmoduli(igrain,ip,el)
583 igrain3 = grain1to3(igrain,prm%Nconstituents)
585 interfaceloop:
do iface = 1,6
586 intface = getinterface(iface,igrain3)
587 nvect = interfacenormal(intface,instance,of)
589 ignghb3(abs(intface(1))) = ignghb3(abs(intface(1))) &
590 + int(real(intface(1),preal)/real(abs(intface(1)),preal))
591 where(ignghb3 < 1) ignghb3 = ngdim
592 where(ignghb3 >ngdim) ignghb3 = 1
593 ignghb = grain3to1(ignghb3,prm%Nconstituents)
594 gmoduli = equivalentmoduli(ignghb,ip,el)
597 gdef = 0.5_preal*(fdef(1:3,1:3,ignghb) - fdef(1:3,1:3,igrain))
603 do i = 1,3;
do j = 1,3
604 do k = 1,3;
do l = 1,3
605 ndef(i,j) = ndef(i,j) - nvect(k)*gdef(i,l)*math_levicivita(j,k,l)
607 ndefnorm = ndefnorm + ndef(i,j)**2.0_preal
609 ndefnorm = max(ndeftoler,sqrt(ndefnorm))
610 nmis(abs(intface(1)),igrain) = nmis(abs(intface(1)),igrain) + ndefnorm
621 do i = 1,3;
do j = 1,3;
do k = 1,3;
do l = 1,3
622 rpen(i,j,igrain) = rpen(i,j,igrain) + 0.5_preal*(mugrain*bggrain + mugnghb*bggnghb)*prm%xiAlpha &
623 *surfcorr(abs(intface(1)))/prm%dAlpha(abs(intface(1))) &
624 *cosh(prm%ciAlpha*ndefnorm) &
625 *0.5_preal*nvect(l)*ndef(i,k)/ndefnorm*math_levicivita(k,l,j) &
626 *tanh(ndefnorm/num%xSmoo)
627 enddo; enddo;enddo;
enddo
640 end subroutine stresspenalty
646 subroutine volumepenalty(vPen,vDiscrep,fAvg,fDef,nGrain,instance,of)
648 real(pReal),
dimension (:,:,:),
intent(out) :: vPen
649 real(pReal),
intent(out) :: vDiscrep
651 real(pReal),
dimension (:,:,:),
intent(in) :: fDef
652 real(pReal),
dimension (3,3),
intent(in) :: fAvg
653 integer,
intent(in) :: &
658 real(pReal),
dimension(size(vPen,3)) :: gVol
663 vdiscrep = math_det33(favg)
665 gvol(i) = math_det33(fdef(1:3,1:3,i))
666 vdiscrep = vdiscrep - gvol(i)/real(ngrain,preal)
674 vpen(:,:,i) = -1.0_preal/real(ngrain,preal)*num%volDiscrMod*num%volDiscrPow/num%maxVolDiscr* &
675 sign((abs(vdiscrep)/num%maxVolDiscr)**(num%volDiscrPow - 1.0),vdiscrep)* &
676 gvol(i)*transpose(math_inv33(fdef(:,:,i)))
687 end subroutine volumepenalty
694 function surfacecorrection(avgF,instance,of)
696 real(pReal),
dimension(3) :: surfaceCorrection
698 real(pReal),
dimension(3,3),
intent(in) :: avgF
699 integer,
intent(in) :: &
702 real(pReal),
dimension(3,3) :: invC
703 real(pReal),
dimension(3) :: nVect
708 call math_invert33(invc,detf,error,matmul(transpose(avgf),avgf))
710 surfacecorrection = 0.0_preal
712 nvect = interfacenormal([ibase,1,1,1],instance,of)
713 do i = 1,3;
do j = 1,3
714 surfacecorrection(ibase) = surfacecorrection(ibase) + invc(i,j)*nvect(i)*nvect(j)
716 surfacecorrection(ibase) = sqrt(surfacecorrection(ibase))*detf
719 end function surfacecorrection
725 function equivalentmoduli(grainID,ip,el)
727 real(pReal),
dimension(2) :: equivalentModuli
729 integer,
intent(in) :: &
731 ip, & !< integration point number
733 real(pReal),
dimension(6,6) :: elasTens
739 elastens = constitutive_homogenizedc(grainid,ip,el)
743 cequiv_11 = (elastens(1,1) + elastens(2,2) + elastens(3,3))/3.0_preal
744 cequiv_12 = (elastens(1,2) + elastens(2,3) + elastens(3,1) + &
745 elastens(1,3) + elastens(2,1) + elastens(3,2))/6.0_preal
746 cequiv_44 = (elastens(4,4) + elastens(5,5) + elastens(6,6))/3.0_preal
747 equivalentmoduli(1) = 0.2_preal*(cequiv_11 - cequiv_12) + 0.6_preal*cequiv_44
751 equivalentmoduli(2) = 2.5e-10_preal
753 end function equivalentmoduli
760 subroutine graindeformation(F, avgF, instance, of)
762 real(pReal),
dimension(:,:,:),
intent(out) :: F
764 real(pReal),
dimension(:,:),
intent(in) :: avgF
765 integer,
intent(in) :: &
769 real(pReal),
dimension(3) :: aVect,nVect
770 integer,
dimension(4) :: intFace
771 integer,
dimension(3) :: iGrain3
772 integer :: iGrain,iFace,i,j
777 associate(prm => param(instance))
780 do igrain = 1,product(prm%Nconstituents)
781 igrain3 = grain1to3(igrain,prm%Nconstituents)
783 intface = getinterface(iface,igrain3)
784 avect = relaxationvector(intface,instance,of)
785 nvect = interfacenormal(intface,instance,of)
786 forall (i=1:3,j=1:3) &
787 f(i,j,igrain) = f(i,j,igrain) + avect(i)*nvect(j)
789 f(1:3,1:3,igrain) = f(1:3,1:3,igrain) + avgf
794 end subroutine graindeformation
796 end procedure mech_rgc_updatestate
802 module subroutine mech_rgc_averagestressanditstangent(avgp,davgpdavgf,p,dpdf,instance)
804 real(pReal),
dimension (3,3),
intent(out) :: avgP
805 real(pReal),
dimension (3,3,3,3),
intent(out) :: dAvgPdAvgF
807 real(pReal),
dimension (:,:,:),
intent(in) :: P
808 real(pReal),
dimension (:,:,:,:,:),
intent(in) :: dPdF
809 integer,
intent(in) :: instance
811 avgp = sum(p,3) /real(product(param(instance)%Nconstituents),preal)
812 davgpdavgf = sum(dpdf,5)/real(product(param(instance)%Nconstituents),preal)
814 end subroutine mech_rgc_averagestressanditstangent
820 module subroutine mech_rgc_results(instance,group)
822 integer,
intent(in) :: instance
823 character(len=*),
intent(in) :: group
827 associate(stt => state(instance), dst => dependentstate(instance), prm => param(instance))
828 outputsloop:
do o = 1,
size(prm%output)
829 select case(trim(prm%output(o)))
830 case(
'constitutivework')
831 call results_writedataset(group,stt%work,
'W',&
832 'work density',
'J/m³')
833 case(
'magnitudemismatch')
834 call results_writedataset(group,dst%mismatch,
'N',&
835 'average mismatch tensor',
'1')
836 case(
'penaltyenergy')
837 call results_writedataset(group,stt%penaltyEnergy,
'R',&
838 'mismatch penalty density',
'J/m³')
839 case(
'volumediscrepancy')
840 call results_writedataset(group,dst%volumeDiscrepancy,
'Delta_V',&
841 'volume discrepancy',
'm³')
842 case(
'maximumrelaxrate')
843 call results_writedataset(group,dst%relaxationrate_max,
'max_alpha_dot',&
844 'maximum relaxation rate',
'm/s')
845 case(
'averagerelaxrate')
846 call results_writedataset(group,dst%relaxationrate_avg,
'avg_alpha_dot',&
847 'average relaxation rate',
'm/s')
852 end subroutine mech_rgc_results
858 pure function relaxationvector(intFace,instance,of)
860 real(pReal),
dimension (3) :: relaxationVector
862 integer,
intent(in) :: instance,of
863 integer,
dimension(4),
intent(in) :: intFace
870 inum = interface4to1(intface,param(instance)%Nconstituents)
872 relaxationvector = state(instance)%relaxationVector((3*inum-2):(3*inum),of)
874 relaxationvector = 0.0_preal
877 end function relaxationvector
883 pure function interfacenormal(intFace,instance,of)
885 real(pReal),
dimension(3) :: interfaceNormal
887 integer,
dimension(4),
intent(in) :: intFace
888 integer,
intent(in) :: &
896 interfacenormal = 0.0_preal
897 npos = abs(intface(1))
898 interfacenormal(npos) = real(intface(1)/abs(intface(1)),preal)
900 interfacenormal = matmul(dependentstate(instance)%orientation(1:3,1:3,of),interfacenormal)
902 end function interfacenormal
908 pure function getinterface(iFace,iGrain3)
910 integer,
dimension(4) :: getInterface
912 integer,
dimension(3),
intent(in) :: iGrain3
913 integer,
intent(in) :: iFace
917 idir = (int(real(iface-1,preal)/2.0_preal)+1)*(-1)**iface
918 getinterface(1) = idir
922 getinterface(2:4) = igrain3
923 if (idir < 0) getinterface(1-idir) = getinterface(1-idir)-1
925 end function getinterface
931 pure function grain1to3(grain1,nGDim)
933 integer,
dimension(3) :: grain1to3
935 integer,
intent(in) :: grain1
936 integer,
dimension(3),
intent(in) :: nGDim
938 grain1to3 = 1 + [mod((grain1-1), ngdim(1)), &
939 mod((grain1-1)/ ngdim(1),ngdim(2)), &
940 (grain1-1)/(ngdim(1)*ngdim(2))]
942 end function grain1to3
948 integer pure function grain3to1(grain3,nGDim)
950 integer,
dimension(3),
intent(in) :: grain3
951 integer,
dimension(3),
intent(in) :: nGDim
953 grain3to1 = grain3(1) &
954 + ngdim(1)*(grain3(2)-1) &
955 + ngdim(1)*ngdim(2)*(grain3(3)-1)
957 end function grain3to1
963 integer pure function interface4to1(iFace4D, nGDim)
965 integer,
dimension(4),
intent(in) :: iFace4D
966 integer,
dimension(3),
intent(in) :: nGDim
969 select case(abs(iface4d(1)))
972 if ((iface4d(2) == 0) .or. (iface4d(2) == ngdim(1)))
then
975 interface4to1 = iface4d(3) + ngdim(2)*(iface4d(4)-1) &
976 + ngdim(2)*ngdim(3)*(iface4d(2)-1)
980 if ((iface4d(3) == 0) .or. (iface4d(3) == ngdim(2)))
then
983 interface4to1 = iface4d(4) + ngdim(3)*(iface4d(2)-1) &
984 + ngdim(3)*ngdim(1)*(iface4d(3)-1) &
985 + (ngdim(1)-1)*ngdim(2)*ngdim(3)
989 if ((iface4d(4) == 0) .or. (iface4d(4) == ngdim(3)))
then
992 interface4to1 = iface4d(2) + ngdim(1)*(iface4d(3)-1) &
993 + ngdim(1)*ngdim(2)*(iface4d(4)-1) &
994 + (ngdim(1)-1)*ngdim(2)*ngdim(3) &
995 + ngdim(1)*(ngdim(2)-1)*ngdim(3)
1003 end function interface4to1
1009 pure function interface1to4(iFace1D, nGDim)
1011 integer,
dimension(4) :: interface1to4
1013 integer,
intent(in) :: iFace1D
1014 integer,
dimension(3),
intent(in) :: nGDim
1015 integer,
dimension(3) :: nIntFace
1019 nintface = [(ngdim(1)-1)*ngdim(2)*ngdim(3), &
1020 ngdim(1)*(ngdim(2)-1)*ngdim(3), &
1021 ngdim(1)*ngdim(2)*(ngdim(3)-1)]
1025 if (iface1d > 0 .and. iface1d <= nintface(1))
then
1026 interface1to4(1) = 1
1027 interface1to4(3) = mod((iface1d-1),ngdim(2))+1
1028 interface1to4(4) = mod(int(real(iface1d-1,preal)/real(ngdim(2),preal)),ngdim(3))+1
1029 interface1to4(2) = int(real(iface1d-1,preal)/real(ngdim(2),preal)/real(ngdim(3),preal))+1
1030 elseif (iface1d > nintface(1) .and. iface1d <= (nintface(2) + nintface(1)))
then
1031 interface1to4(1) = 2
1032 interface1to4(4) = mod((iface1d-nintface(1)-1),ngdim(3))+1
1033 interface1to4(2) = mod(int(real(iface1d-nintface(1)-1,preal)/real(ngdim(3),preal)),ngdim(1))+1
1034 interface1to4(3) = int(real(iface1d-nintface(1)-1,preal)/real(ngdim(3),preal)/real(ngdim(1),preal))+1
1035 elseif (iface1d > nintface(2) + nintface(1) .and. iface1d <= (nintface(3) + nintface(2) + nintface(1)))
then
1036 interface1to4(1) = 3
1037 interface1to4(2) = mod((iface1d-nintface(2)-nintface(1)-1),ngdim(1))+1
1038 interface1to4(3) = mod(int(real(iface1d-nintface(2)-nintface(1)-1,preal)/real(ngdim(1),preal)),ngdim(2))+1
1039 interface1to4(4) = int(real(iface1d-nintface(2)-nintface(1)-1,preal)/real(ngdim(1),preal)/real(ngdim(2),preal))+1
1042 end function interface1to4
1045 end submodule homogenization_mech_rgc