DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
homogenization.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization.f90"
5 !--------------------------------------------------------------------------------------------------
10 !--------------------------------------------------------------------------------------------------
12  use prec
13  use io
14  use config
15  use debug
16  use math
17  use material
18  use numerics
19  use constitutive
20  use crystallite
21  use fesolving
22  use discretization
26  use damage_none
27  use damage_local
28  use damage_nonlocal
29  use results
30 
31  implicit none
32  private
33 
34 !--------------------------------------------------------------------------------------------------
35 ! General variables for the homogenization at a material point
36  real(preal), dimension(:,:,:,:), allocatable, public :: &
37  materialpoint_f0, & !< def grad of IP at start of FE increment
38  materialpoint_f, & !< def grad of IP to be reached at end of FE increment
40  real(preal), dimension(:,:,:,:,:,:), allocatable, public :: &
42 
43  real(preal), dimension(:,:,:,:), allocatable :: &
44  materialpoint_subf0, & !< def grad of IP at beginning of homogenization increment
46  real(preal), dimension(:,:), allocatable :: &
50  logical, dimension(:,:), allocatable :: &
53  logical, dimension(:,:,:), allocatable :: &
55 
56  type :: tnumerics
57  integer :: &
58  nmpstate
59  real(preal) :: &
60  substepminhomog, & !< minimum (relative) size of sub-step allowed during cutback in homogenization
61  substepsizehomog, & !< size of first substep when cutback in homogenization
62  stepincreasehomog
63  end type tnumerics
64 
65  type(tnumerics) :: num
66 
67  interface
68 
69  module subroutine mech_none_init
70  end subroutine mech_none_init
71 
72  module subroutine mech_isostrain_init
73  end subroutine mech_isostrain_init
74 
75  module subroutine mech_rgc_init
76  end subroutine mech_rgc_init
77 
78 
79  module subroutine mech_isostrain_partitiondeformation(f,avgf)
80  real(preal), dimension (:,:,:), intent(out) :: f
81  real(preal), dimension (3,3), intent(in) :: avgf
82  end subroutine mech_isostrain_partitiondeformation
83 
84  module subroutine mech_rgc_partitiondeformation(f,avgf,instance,of)
85  real(preal), dimension (:,:,:), intent(out) :: f
86  real(preal), dimension (3,3), intent(in) :: avgf
87  integer, intent(in) :: &
88  instance, &
89  of
90  end subroutine mech_rgc_partitiondeformation
91 
92 
93  module subroutine mech_isostrain_averagestressanditstangent(avgp,davgpdavgf,p,dpdf,instance)
94  real(preal), dimension (3,3), intent(out) :: avgp
95  real(preal), dimension (3,3,3,3), intent(out) :: davgpdavgf
96 
97  real(preal), dimension (:,:,:), intent(in) :: p
98  real(preal), dimension (:,:,:,:,:), intent(in) :: dpdf
99  integer, intent(in) :: instance
100  end subroutine mech_isostrain_averagestressanditstangent
101 
102  module subroutine mech_rgc_averagestressanditstangent(avgp,davgpdavgf,p,dpdf,instance)
103  real(preal), dimension (3,3), intent(out) :: avgp
104  real(preal), dimension (3,3,3,3), intent(out) :: davgpdavgf
105 
106  real(preal), dimension (:,:,:), intent(in) :: p
107  real(preal), dimension (:,:,:,:,:), intent(in) :: dpdf
108  integer, intent(in) :: instance
109  end subroutine mech_rgc_averagestressanditstangent
110 
111 
112  module function mech_rgc_updatestate(p,f,f0,avgf,dt,dpdf,ip,el)
113  logical, dimension(2) :: mech_rgc_updatestate
114  real(preal), dimension(:,:,:), intent(in) :: &
115  p,& !< partitioned stresses
116  f,& !< partitioned deformation gradients
117  f0
118  real(preal), dimension(:,:,:,:,:), intent(in) :: dpdf
119  real(preal), dimension(3,3), intent(in) :: avgf
120  real(preal), intent(in) :: dt
121  integer, intent(in) :: &
122  ip, & !< integration point number
123  el
124  end function mech_rgc_updatestate
125 
126 
127  module subroutine mech_rgc_results(instance,group)
128  integer, intent(in) :: instance
129  character(len=*), intent(in) :: group
130  end subroutine mech_rgc_results
131 
132  end interface
133 
134  public :: &
135  homogenization_init, &
136  materialpoint_stressanditstangent, &
137  homogenization_results
138 
139 contains
140 
141 
142 !--------------------------------------------------------------------------------------------------
144 !--------------------------------------------------------------------------------------------------
145 subroutine homogenization_init
146 
147  if (any(homogenization_type == homogenization_none_id)) call mech_none_init
148  if (any(homogenization_type == homogenization_isostrain_id)) call mech_isostrain_init
149  if (any(homogenization_type == homogenization_rgc_id)) call mech_rgc_init
150 
154 
158 
159  call config_deallocate('material.config/homogenization')
160 
161 !--------------------------------------------------------------------------------------------------
162 ! allocate and initialize global variables
163  allocate(materialpoint_dpdf(3,3,3,3,discretization_nip,discretization_nelem), source=0.0_preal)
164  materialpoint_f0 = spread(spread(math_i3,3,discretization_nip),4,discretization_nelem) ! initialize to identity
165  materialpoint_f = materialpoint_f0 ! initialize to identity
166  allocate(materialpoint_subf0(3,3,discretization_nip,discretization_nelem), source=0.0_preal)
167  allocate(materialpoint_subf(3,3,discretization_nip,discretization_nelem), source=0.0_preal)
168  allocate(materialpoint_p(3,3,discretization_nip,discretization_nelem), source=0.0_preal)
169  allocate(materialpoint_subfrac(discretization_nip,discretization_nelem), source=0.0_preal)
170  allocate(materialpoint_substep(discretization_nip,discretization_nelem), source=0.0_preal)
171  allocate(materialpoint_subdt(discretization_nip,discretization_nelem), source=0.0_preal)
172  allocate(materialpoint_requested(discretization_nip,discretization_nelem), source=.false.)
173  allocate(materialpoint_converged(discretization_nip,discretization_nelem), source=.true.)
174  allocate(materialpoint_doneandhappy(2,discretization_nip,discretization_nelem), source=.true.)
175 
176  write(6,'(/,a)') ' <<<+- homogenization init -+>>>'; flush(6)
177 
178  if (iand(debug_level(debug_homogenization), debug_levelbasic) /= 0) then
179  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_dPdF: ', shape(materialpoint_dpdf)
180  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_F0: ', shape(materialpoint_f0)
181  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_F: ', shape(materialpoint_f)
182  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_subF0: ', shape(materialpoint_subf0)
183  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_subF: ', shape(materialpoint_subf)
184  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_P: ', shape(materialpoint_p)
185  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_subFrac: ', shape(materialpoint_subfrac)
186  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_subStep: ', shape(materialpoint_substep)
187  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_subdt: ', shape(materialpoint_subdt)
188  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_requested: ', shape(materialpoint_requested)
189  write(6,'(a32,1x,7(i8,1x))') 'materialpoint_converged: ', shape(materialpoint_converged)
190  write(6,'(a32,1x,7(i8,1x),/)') 'materialpoint_doneAndHappy: ', shape(materialpoint_doneandhappy)
191  endif
192  flush(6)
193 
195  call io_error(602,ext_msg='constituent', el=debug_e, g=debug_g)
196 
197  num%nMPstate = config_numerics%getInt( 'nmpstate', defaultval=10)
198  num%subStepMinHomog = config_numerics%getFloat('substepminhomog', defaultval=1.0e-3_preal)
199  num%subStepSizeHomog = config_numerics%getFloat('substepsizehomog', defaultval=0.25_preal)
200  num%stepIncreaseHomog = config_numerics%getFloat('stepincreasehomog', defaultval=1.5_preal)
201  if (num%nMPstate < 1) call io_error(301,ext_msg='nMPstate')
202  if (num%subStepMinHomog <= 0.0_preal) call io_error(301,ext_msg='subStepMinHomog')
203  if (num%subStepSizeHomog <= 0.0_preal) call io_error(301,ext_msg='subStepSizeHomog')
204  if (num%stepIncreaseHomog <= 0.0_preal) call io_error(301,ext_msg='stepIncreaseHomog')
205 
206 end subroutine homogenization_init
207 
208 
209 !--------------------------------------------------------------------------------------------------
211 !--------------------------------------------------------------------------------------------------
212 subroutine materialpoint_stressanditstangent(updateJaco,dt)
213 
214  real(preal), intent(in) :: dt
215  logical, intent(in) :: updatejaco
216  integer :: &
217  niterationhomog, &
218  niterationmpstate, &
219  g, & !< grain number
220  i, & !< integration point number
221  e, & !< element number
222  mysource, &
223  myngrains
224 
225 # 231 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization.f90"
226 
227 !--------------------------------------------------------------------------------------------------
228 ! initialize restoration points of ...
231  do i = fesolving_execip(1),fesolving_execip(2);
232  do g = 1,myngrains
233 
234  plasticstate(material_phaseat(g,e))%partionedState0(:,material_phasememberat(g,i,e)) = &
236  do mysource = 1, phase_nsources(material_phaseat(g,e))
237  sourcestate(material_phaseat(g,e))%p(mysource)%partionedState0(:,material_phasememberat(g,i,e)) = &
238  sourcestate(material_phaseat(g,e))%p(mysource)%state0( :,material_phasememberat(g,i,e))
239  enddo
240 
241  crystallite_partionedfp0(1:3,1:3,g,i,e) = crystallite_fp0(1:3,1:3,g,i,e)
242  crystallite_partionedlp0(1:3,1:3,g,i,e) = crystallite_lp0(1:3,1:3,g,i,e)
243  crystallite_partionedfi0(1:3,1:3,g,i,e) = crystallite_fi0(1:3,1:3,g,i,e)
244  crystallite_partionedli0(1:3,1:3,g,i,e) = crystallite_li0(1:3,1:3,g,i,e)
245  crystallite_partionedf0(1:3,1:3,g,i,e) = crystallite_f0(1:3,1:3,g,i,e)
246  crystallite_partioneds0(1:3,1:3,g,i,e) = crystallite_s0(1:3,1:3,g,i,e)
247 
248  enddo
249 
250 
251  materialpoint_subf0(1:3,1:3,i,e) = materialpoint_f0(1:3,1:3,i,e)
252  materialpoint_subfrac(i,e) = 0.0_preal
253  materialpoint_substep(i,e) = 1.0_preal/num%subStepSizeHomog ! <<added to adopt flexibility in cutback size>>
254  materialpoint_converged(i,e) = .false. ! pretend failed step of twice the required size
255  materialpoint_requested(i,e) = .true. ! everybody requires calculation
256 
257  if (homogstate(material_homogenizationat(e))%sizeState > 0) &
259  homogstate(material_homogenizationat(e))%State0( :,material_homogenizationmemberat(i,e)) ! ...internal homogenization state
260 
261  if (thermalstate(material_homogenizationat(e))%sizeState > 0) &
263  thermalstate(material_homogenizationat(e))%State0( :,material_homogenizationmemberat(i,e)) ! ...internal thermal state
264 
265  if (damagestate(material_homogenizationat(e))%sizeState > 0) &
267  damagestate(material_homogenizationat(e))%State0( :,material_homogenizationmemberat(i,e)) ! ...internal damage state
268  enddo
269  enddo
270 
271  niterationhomog = 0
272 
273  cutbacklooping: do while (.not. terminallyill .and. &
274  any(materialpoint_substep(:,fesolving_execelem(1):fesolving_execelem(2)) > num%subStepMinHomog))
275 
276  !$OMP PARALLEL DO PRIVATE(myNgrains)
277  elementlooping1: do e = fesolving_execelem(1),fesolving_execelem(2)
279  iplooping1: do i = fesolving_execip(1),fesolving_execip(2)
280 
281  converged: if (materialpoint_converged(i,e)) then
282 # 296 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization.f90"
283 
284 !---------------------------------------------------------------------------------------------------
285 ! calculate new subStep and new subFrac
286  materialpoint_subfrac(i,e) = materialpoint_subfrac(i,e) + materialpoint_substep(i,e)
287  materialpoint_substep(i,e) = min(1.0_preal-materialpoint_subfrac(i,e), &
288  num%stepIncreaseHomog*materialpoint_substep(i,e)) ! introduce flexibility for step increase/acceleration
289 
290  steppingneeded: if (materialpoint_substep(i,e) > num%subStepMinHomog) then
291 
292  ! wind forward grain starting point of...
293  crystallite_partionedf0(1:3,1:3,1:myngrains,i,e) = &
294  crystallite_partionedf(1:3,1:3,1:myngrains,i,e)
295 
296  crystallite_partionedfp0(1:3,1:3,1:myngrains,i,e) = &
297  crystallite_fp(1:3,1:3,1:myngrains,i,e)
298 
299  crystallite_partionedlp0(1:3,1:3,1:myngrains,i,e) = &
300  crystallite_lp(1:3,1:3,1:myngrains,i,e)
301 
302  crystallite_partionedfi0(1:3,1:3,1:myngrains,i,e) = &
303  crystallite_fi(1:3,1:3,1:myngrains,i,e)
304 
305  crystallite_partionedli0(1:3,1:3,1:myngrains,i,e) = &
306  crystallite_li(1:3,1:3,1:myngrains,i,e)
307 
308  crystallite_partioneds0(1:3,1:3,1:myngrains,i,e) = &
309  crystallite_s(1:3,1:3,1:myngrains,i,e)
310 
311  do g = 1,myngrains
312  plasticstate(material_phaseat(g,e))%partionedState0(:,material_phasememberat(g,i,e)) = &
314  do mysource = 1, phase_nsources(material_phaseat(g,e))
315  sourcestate(material_phaseat(g,e))%p(mysource)%partionedState0(:,material_phasememberat(g,i,e)) = &
316  sourcestate(material_phaseat(g,e))%p(mysource)%state (:,material_phasememberat(g,i,e))
317  enddo
318  enddo
319 
320  if(homogstate(material_homogenizationat(e))%sizeState > 0) &
323  if(thermalstate(material_homogenizationat(e))%sizeState > 0) &
326  if(damagestate(material_homogenizationat(e))%sizeState > 0) &
329 
330  materialpoint_subf0(1:3,1:3,i,e) = materialpoint_subf(1:3,1:3,i,e)
331 
332  endif steppingneeded
333 
334  else converged
335  if ( (myngrains == 1 .and. materialpoint_substep(i,e) <= 1.0 ) .or. & ! single grain already tried internal subStepping in crystallite
336  num%subStepSizeHomog * materialpoint_substep(i,e) <= num%subStepMinHomog ) then ! would require too small subStep
337  ! cutback makes no sense
338  !$OMP FLUSH(terminallyIll)
339  if (.not. terminallyill) then ! so first signals terminally ill...
340  !$OMP CRITICAL (write2out)
341  write(6,*) 'Integration point ', i,' at element ', e, ' terminally ill'
342  !$OMP END CRITICAL (write2out)
343  endif
344  terminallyill = .true. ! ...and kills all others
345  else ! cutback makes sense
346  materialpoint_substep(i,e) = num%subStepSizeHomog * materialpoint_substep(i,e) ! crystallite had severe trouble, so do a significant cutback
347 
348 # 370 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization.f90"
349 
350 !--------------------------------------------------------------------------------------------------
351 ! restore...
352  if (materialpoint_substep(i,e) < 1.0_preal) then ! protect against fake cutback from \Delta t = 2 to 1. Maybe that "trick" is not necessary anymore at all? I.e. start with \Delta t = 1
353  crystallite_lp(1:3,1:3,1:myngrains,i,e) = &
354  crystallite_partionedlp0(1:3,1:3,1:myngrains,i,e)
355  crystallite_li(1:3,1:3,1:myngrains,i,e) = &
356  crystallite_partionedli0(1:3,1:3,1:myngrains,i,e)
357  endif ! maybe protecting everything from overwriting (not only L) makes even more sense
358  crystallite_fp(1:3,1:3,1:myngrains,i,e) = &
359  crystallite_partionedfp0(1:3,1:3,1:myngrains,i,e)
360  crystallite_fi(1:3,1:3,1:myngrains,i,e) = &
361  crystallite_partionedfi0(1:3,1:3,1:myngrains,i,e)
362  crystallite_s(1:3,1:3,1:myngrains,i,e) = &
363  crystallite_partioneds0(1:3,1:3,1:myngrains,i,e)
364  do g = 1, myngrains
365  plasticstate(material_phaseat(g,e))%state( :,material_phasememberat(g,i,e)) = &
366  plasticstate(material_phaseat(g,e))%partionedState0(:,material_phasememberat(g,i,e))
367  do mysource = 1, phase_nsources(material_phaseat(g,e))
368  sourcestate(material_phaseat(g,e))%p(mysource)%state( :,material_phasememberat(g,i,e)) = &
369  sourcestate(material_phaseat(g,e))%p(mysource)%partionedState0(:,material_phasememberat(g,i,e))
370  enddo
371  enddo
372  if(homogstate(material_homogenizationat(e))%sizeState > 0) &
375  if(thermalstate(material_homogenizationat(e))%sizeState > 0) &
378  if(damagestate(material_homogenizationat(e))%sizeState > 0) &
381  endif
382  endif converged
383 
384  if (materialpoint_substep(i,e) > num%subStepMinHomog) then
385  materialpoint_requested(i,e) = .true.
386  materialpoint_subf(1:3,1:3,i,e) = materialpoint_subf0(1:3,1:3,i,e) &
387  + materialpoint_substep(i,e) * (materialpoint_f(1:3,1:3,i,e) &
388  - materialpoint_f0(1:3,1:3,i,e))
389  materialpoint_subdt(i,e) = materialpoint_substep(i,e) * dt
390  materialpoint_doneandhappy(1:2,i,e) = [.false.,.true.]
391  endif
392  enddo iplooping1
393  enddo elementlooping1
394  !$OMP END PARALLEL DO
395 
396  niterationmpstate = 0
397 
398  convergencelooping: do while (.not. terminallyill .and. &
399  any( materialpoint_requested(:,fesolving_execelem(1):fesolving_execelem(2)) &
400  .and. .not. materialpoint_doneandhappy(1,:,fesolving_execelem(1):fesolving_execelem(2)) &
401  ) .and. &
402  niterationmpstate < num%nMPstate)
403  niterationmpstate = niterationmpstate + 1
404 
405 !--------------------------------------------------------------------------------------------------
406 ! deformation partitioning
407 ! based on materialpoint_subF0,.._subF,crystallite_partionedF0, and homogenization_state,
408 ! results in crystallite_partionedF
409  !$OMP PARALLEL DO PRIVATE(myNgrains)
410  elementlooping2: do e = fesolving_execelem(1),fesolving_execelem(2)
412  iplooping2: do i = fesolving_execip(1),fesolving_execip(2)
413  if ( materialpoint_requested(i,e) .and. & ! process requested but...
414  .not. materialpoint_doneandhappy(1,i,e)) then ! ...not yet done material points
415  call partitiondeformation(i,e) ! partition deformation onto constituents
416  crystallite_dt(1:myngrains,i,e) = materialpoint_subdt(i,e) ! propagate materialpoint dt to grains
417  crystallite_requested(1:myngrains,i,e) = .true. ! request calculation for constituents
418  else
419  crystallite_requested(1:myngrains,i,e) = .false. ! calculation for constituents not required anymore
420  endif
421  enddo iplooping2
422  enddo elementlooping2
423  !$OMP END PARALLEL DO
424 
425 !--------------------------------------------------------------------------------------------------
426 ! crystallite integration
427 ! based on crystallite_partionedF0,.._partionedF
428 ! incrementing by crystallite_dt
429 
430  materialpoint_converged = crystallite_stress() !ToDo: MD not sure if that is the best logic
431 
432 !--------------------------------------------------------------------------------------------------
433 ! state update
434  !$OMP PARALLEL DO
435  elementlooping3: do e = fesolving_execelem(1),fesolving_execelem(2)
436  iplooping3: do i = fesolving_execip(1),fesolving_execip(2)
437  if ( materialpoint_requested(i,e) .and. &
438  .not. materialpoint_doneandhappy(1,i,e)) then
439  if (.not. materialpoint_converged(i,e)) then
440  materialpoint_doneandhappy(1:2,i,e) = [.true.,.false.]
441  else
442  materialpoint_doneandhappy(1:2,i,e) = updatestate(i,e)
443  materialpoint_converged(i,e) = all(materialpoint_doneandhappy(1:2,i,e)) ! converged if done and happy
444  endif
445  endif
446  enddo iplooping3
447  enddo elementlooping3
448  !$OMP END PARALLEL DO
449 
450  enddo convergencelooping
451 
452  niterationhomog = niterationhomog + 1
453 
454  enddo cutbacklooping
455 
456  if(updatejaco) call crystallite_stresstangent
457 
458  if (.not. terminallyill ) then
459  call crystallite_orientations() ! calculate crystal orientations
460  !$OMP PARALLEL DO
461  elementlooping4: do e = fesolving_execelem(1),fesolving_execelem(2)
462  iplooping4: do i = fesolving_execip(1),fesolving_execip(2)
463  call averagestressanditstangent(i,e)
464  enddo iplooping4
465  enddo elementlooping4
466  !$OMP END PARALLEL DO
467  else
468  write(6,'(/,a,/)') '<< HOMOG >> Material Point terminally ill'
469  endif
470 
471 end subroutine materialpoint_stressanditstangent
472 
473 
474 !--------------------------------------------------------------------------------------------------
476 !--------------------------------------------------------------------------------------------------
477 subroutine partitiondeformation(ip,el)
478 
479  integer, intent(in) :: &
480  ip, & !< integration point
481  el
482 
483  chosenhomogenization: select case(homogenization_type(material_homogenizationat(el)))
484 
485  case (homogenization_none_id) chosenhomogenization
486  crystallite_partionedf(1:3,1:3,1,ip,el) = materialpoint_subf(1:3,1:3,ip,el)
487 
488  case (homogenization_isostrain_id) chosenhomogenization
489  call mech_isostrain_partitiondeformation(&
491  materialpoint_subf(1:3,1:3,ip,el))
492 
493  case (homogenization_rgc_id) chosenhomogenization
494  call mech_rgc_partitiondeformation(&
496  materialpoint_subf(1:3,1:3,ip,el),&
497  ip, &
498  el)
499  end select chosenhomogenization
500 
501 end subroutine partitiondeformation
502 
503 
504 !--------------------------------------------------------------------------------------------------
507 !--------------------------------------------------------------------------------------------------
508 function updatestate(ip,el)
509 
510  integer, intent(in) :: &
511  ip, & !< integration point
512  el
513  logical, dimension(2) :: updatestate
514 
515  updatestate = .true.
516  chosenhomogenization: select case(homogenization_type(material_homogenizationat(el)))
517  case (homogenization_rgc_id) chosenhomogenization
518  updatestate = &
519  updatestate .and. &
520  mech_rgc_updatestate(crystallite_p(1:3,1:3,1:homogenization_ngrains(material_homogenizationat(el)),ip,el), &
523  materialpoint_subf(1:3,1:3,ip,el),&
524  materialpoint_subdt(ip,el), &
525  crystallite_dpdf(1:3,1:3,1:3,1:3,1:homogenization_ngrains(material_homogenizationat(el)),ip,el), &
526  ip, &
527  el)
528  end select chosenhomogenization
529 
530  chosenthermal: select case (thermal_type(material_homogenizationat(el)))
531  case (thermal_adiabatic_id) chosenthermal
532  updatestate = &
533  updatestate .and. &
534  thermal_adiabatic_updatestate(materialpoint_subdt(ip,el), &
535  ip, &
536  el)
537  end select chosenthermal
538 
539  chosendamage: select case (damage_type(material_homogenizationat(el)))
540  case (damage_local_id) chosendamage
541  updatestate = &
542  updatestate .and. &
543  damage_local_updatestate(materialpoint_subdt(ip,el), &
544  ip, &
545  el)
546  end select chosendamage
547 
548 end function updatestate
549 
550 
551 !--------------------------------------------------------------------------------------------------
553 !--------------------------------------------------------------------------------------------------
554 subroutine averagestressanditstangent(ip,el)
555 
556  integer, intent(in) :: &
557  ip, & !< integration point
558  el
559 
560  chosenhomogenization: select case(homogenization_type(material_homogenizationat(el)))
561  case (homogenization_none_id) chosenhomogenization
562  materialpoint_p(1:3,1:3,ip,el) = crystallite_p(1:3,1:3,1,ip,el)
563  materialpoint_dpdf(1:3,1:3,1:3,1:3,ip,el) = crystallite_dpdf(1:3,1:3,1:3,1:3,1,ip,el)
564 
565  case (homogenization_isostrain_id) chosenhomogenization
566  call mech_isostrain_averagestressanditstangent(&
567  materialpoint_p(1:3,1:3,ip,el), &
568  materialpoint_dpdf(1:3,1:3,1:3,1:3,ip,el),&
570  crystallite_dpdf(1:3,1:3,1:3,1:3,1:homogenization_ngrains(material_homogenizationat(el)),ip,el), &
572 
573  case (homogenization_rgc_id) chosenhomogenization
574  call mech_rgc_averagestressanditstangent(&
575  materialpoint_p(1:3,1:3,ip,el), &
576  materialpoint_dpdf(1:3,1:3,1:3,1:3,ip,el),&
578  crystallite_dpdf(1:3,1:3,1:3,1:3,1:homogenization_ngrains(material_homogenizationat(el)),ip,el), &
580  end select chosenhomogenization
581 
582 end subroutine averagestressanditstangent
583 
584 
585 !--------------------------------------------------------------------------------------------------
587 !--------------------------------------------------------------------------------------------------
588 subroutine homogenization_results
589  use material, only: &
590  material_homogenization_type => homogenization_type
591 
592  integer :: p
593  character(len=pStringLen) :: group_base,group
594 
595  !real(pReal), dimension(:,:,:), allocatable :: temp
596 
597  do p=1,size(config_name_homogenization)
598  group_base = 'current/materialpoint/'//trim(config_name_homogenization(p))
599  call results_closegroup(results_addgroup(group_base))
600 
601  group = trim(group_base)//'/generic'
603  !temp = reshape(materialpoint_F,[3,3,discretization_nIP*discretization_nElem])
604  !call results_writeDataset(group,temp,'F',&
605  ! 'deformation gradient','1')
606  !temp = reshape(materialpoint_P,[3,3,discretization_nIP*discretization_nElem])
607  !call results_writeDataset(group,temp,'P',&
608  ! '1st Piola-Kirchoff stress','Pa')
609 
610  group = trim(group_base)//'/mech'
612  select case(material_homogenization_type(p))
613  case(homogenization_rgc_id)
614  call mech_rgc_results(homogenization_typeinstance(p),group)
615  end select
616 
617  group = trim(group_base)//'/damage'
619  select case(damage_type(p))
620  case(damage_local_id)
621  call damage_local_results(p,group)
622  case(damage_nonlocal_id)
623  call damage_nonlocal_results(p,group)
624  end select
625 
626  group = trim(group_base)//'/thermal'
628  select case(thermal_type(p))
629  case(thermal_adiabatic_id)
630  call thermal_adiabatic_results(p,group)
631  case(thermal_conduction_id)
632  call thermal_conduction_results(p,group)
633  end select
634 
635  enddo
636 
637 end subroutine homogenization_results
638 
639 end module homogenization
material::sourcestate
type(tsourcestate), dimension(:), allocatable, public sourcestate
Definition: material.f90:139
fesolving::fesolving_execip
integer, dimension(2) fesolving_execip
for ping-pong scheme always range to max IP, otherwise one specific IP
Definition: FEsolving.f90:18
crystallite::crystallite_s0
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_s0
2nd Piola-Kirchhoff stress vector at start of FE inc
Definition: crystallite.f90:43
thermal_adiabatic::thermal_adiabatic_init
subroutine, public thermal_adiabatic_init
module initialization
Definition: thermal_adiabatic.f90:47
material::material_phaseat
integer, dimension(:,:), allocatable, public, protected material_phaseat
phase ID of each element
Definition: material.f90:132
crystallite::crystallite_dt
real(preal), dimension(:,:,:), allocatable, public crystallite_dt
requested time increment of each grain
Definition: crystallite.f90:35
homogenization::materialpoint_subdt
real(preal), dimension(:,:), allocatable materialpoint_subdt
Definition: homogenization.f90:46
damage_nonlocal
material subroutine for non-locally evolving damage field
Definition: damage_nonlocal.f90:9
homogenization::materialpoint_dpdf
real(preal), dimension(:,:,:,:,:,:), allocatable, public materialpoint_dpdf
tangent of first P–K stress at IP
Definition: homogenization.f90:40
crystallite::crystallite_stress
logical function, dimension(discretization_nip, discretization_nelem), public crystallite_stress(dummyArgumentToPreventInternalCompilerErrorWithGCC)
calculate stress (P)
Definition: crystallite.f90:284
homogenization::materialpoint_p
real(preal), dimension(:,:,:,:), allocatable, public materialpoint_p
first P–K stress of IP
Definition: homogenization.f90:36
homogenization::materialpoint_converged
logical, dimension(:,:), allocatable materialpoint_converged
Definition: homogenization.f90:50
thermal_isothermal::thermal_isothermal_init
subroutine thermal_isothermal_init
allocates all neccessary fields, reads information from material configuration file
Definition: thermal_isothermal.f90:22
homogenization::materialpoint_substep
real(preal), dimension(:,:), allocatable materialpoint_substep
Definition: homogenization.f90:46
damage_nonlocal::damage_nonlocal_init
subroutine, public damage_nonlocal_init
module initialization
Definition: damage_nonlocal.f90:48
material::material_phasememberat
integer, dimension(:,:,:), allocatable, public, protected material_phasememberat
position of the element within its phase instance
Definition: material.f90:134
io::io_error
subroutine, public io_error(error_ID, el, ip, g, instance, ext_msg)
write error statements to standard out and terminate the Marc/spectral run with exit #9xxx
Definition: IO.f90:305
crystallite::crystallite_partionedfp0
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedfp0
plastic def grad at start of homog inc
Definition: crystallite.f90:52
homogenization::materialpoint_f0
real(preal), dimension(:,:,:,:), allocatable, public materialpoint_f0
def grad of IP at start of FE increment
Definition: homogenization.f90:36
debug::debug_level
integer, dimension(debug_maxntype+2), public, protected debug_level
Definition: debug.f90:48
crystallite::crystallite_p
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_p
1st Piola-Kirchhoff stress per grain
Definition: crystallite.f90:43
material::thermal_adiabatic_id
@, public thermal_adiabatic_id
Definition: material.f90:87
homogenization::materialpoint_subfrac
real(preal), dimension(:,:), allocatable materialpoint_subfrac
Definition: homogenization.f90:46
material::damage_nonlocal_id
@, public damage_nonlocal_id
Definition: material.f90:87
damage_local
material subroutine for locally evolving damage field
Definition: damage_local.f90:9
thermal_adiabatic
material subroutine for adiabatic temperature evolution
Definition: thermal_adiabatic.f90:9
homogenization::materialpoint_f
real(preal), dimension(:,:,:,:), allocatable, public materialpoint_f
def grad of IP to be reached at end of FE increment
Definition: homogenization.f90:36
crystallite::crystallite_fi0
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_fi0
intermediate def grad at start of FE inc
Definition: crystallite.f90:43
crystallite::crystallite_fi
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_fi
current intermediate def grad (end of converged time step)
Definition: crystallite.f90:52
crystallite::crystallite_li
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_li
current intermediate velocitiy grad (end of converged time step)
Definition: crystallite.f90:52
crystallite::crystallite_s
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_s
current 2nd Piola-Kirchhoff stress vector (end of converged time step)
Definition: crystallite.f90:52
material
Parses material config file, either solverJobName.materialConfig or material.config.
Definition: material.f90:11
thermal_conduction::thermal_conduction_results
subroutine, public thermal_conduction_results(homog, group)
writes results to HDF5 output file
Definition: thermal_conduction.f90:244
crystallite::crystallite_fp0
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_fp0
plastic def grad at start of FE inc
Definition: crystallite.f90:43
config
Reads in the material configuration from file.
Definition: config.f90:13
crystallite::crystallite_stresstangent
subroutine, public crystallite_stresstangent
calculate tangent (dPdF)
Definition: crystallite.f90:437
material::damagestate
type(tstate), dimension(:), allocatable, public damagestate
Definition: material.f90:141
homogenization::materialpoint_subf
real(preal), dimension(:,:,:,:), allocatable materialpoint_subf
def grad of IP to be reached at end of homog inc
Definition: homogenization.f90:43
config::config_deallocate
subroutine, public config_deallocate(what)
deallocates the linked lists that store the content of the configuration files
Definition: config.f90:290
material::homogenization_ngrains
integer, dimension(:), allocatable, public, protected homogenization_ngrains
number of grains in each homogenization
Definition: material.f90:113
crystallite::num
type(tnumerics) num
Definition: crystallite.f90:103
crystallite::crystallite_lp0
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_lp0
plastic velocitiy grad at start of FE inc
Definition: crystallite.f90:43
crystallite::crystallite_dpdf
real(preal), dimension(:,:,:,:,:,:,:), allocatable, public, protected crystallite_dpdf
current individual dPdF per grain (end of converged time step)
Definition: crystallite.f90:72
crystallite::crystallite_partionedlp0
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedlp0
plastic velocity grad at start of homog inc
Definition: crystallite.f90:52
material::thermal_isothermal_id
@, public thermal_isothermal_id
Definition: material.f90:87
damage_local::damage_local_init
subroutine, public damage_local_init
module initialization
Definition: damage_local.f90:43
material::damage_local_id
@, public damage_local_id
Definition: material.f90:87
crystallite
crystallite state integration functions and reporting of results
Definition: crystallite.f90:15
material::thermal_type
integer(kind(thermal_isothermal_id)), dimension(:), allocatable, public, protected thermal_type
thermal transport model
Definition: material.f90:94
prec
setting precision for real and int type
Definition: prec.f90:13
material::homogenization_typeinstance
integer, dimension(:), allocatable, public, protected homogenization_typeinstance
instance of particular type of each homogenization
Definition: material.f90:113
material::homogstate
type(tstate), dimension(:), allocatable, public homogstate
Definition: material.f90:141
crystallite::crystallite_li0
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_li0
intermediate velocitiy grad at start of FE inc
Definition: crystallite.f90:43
discretization
spatial discretization
Definition: discretization.f90:9
discretization::discretization_nip
integer, public, protected discretization_nip
Definition: discretization.f90:17
material::homogenization_isostrain_id
@, public homogenization_isostrain_id
Definition: material.f90:87
io
input/output functions, partly depending on chosen solver
Definition: IO.f90:12
debug::debug_homogenization
integer, parameter, public debug_homogenization
Definition: debug.f90:32
homogenization
homogenization manager, organizing deformation partitioning and stress homogenization
Definition: homogenization.f90:11
debug::debug_levelbasic
integer, parameter, public debug_levelbasic
Definition: debug.f90:19
prec::preal
integer, parameter preal
number with 15 significant digits, up to 1e+-307 (typically 64 bit)
Definition: prec.f90:20
thermal_adiabatic::thermal_adiabatic_results
subroutine, public thermal_adiabatic_results(homog, group)
writes results to HDF5 output file
Definition: thermal_adiabatic.f90:235
damage_local::damage_local_results
subroutine, public damage_local_results(homog, group)
writes results to HDF5 output file
Definition: damage_local.f90:166
debug
Reading in and interpretating the debugging settings for the various modules.
Definition: debug.f90:12
damage_none::damage_none_init
subroutine damage_none_init
allocates all neccessary fields, reads information from material configuration file
Definition: damage_none.f90:22
crystallite::crystallite_partioneds0
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partioneds0
2nd Piola-Kirchhoff stress vector at start of homog inc
Definition: crystallite.f90:52
thermal_isothermal
material subroutine for isothermal temperature field
Definition: thermal_isothermal.f90:9
damage_none
material subroutine for constant damage field
Definition: damage_none.f90:9
crystallite::crystallite_partionedfi0
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedfi0
intermediate def grad at start of homog inc
Definition: crystallite.f90:52
results::results_addgroup
integer(hid_t) function, public results_addgroup(groupName)
adds a new group to the results file
Definition: results.f90:154
config::config_name_homogenization
character(len=pstringlen), dimension(:), allocatable, public, protected config_name_homogenization
name of each homogenization
Definition: config.f90:34
material::thermalstate
type(tstate), dimension(:), allocatable, public thermalstate
Definition: material.f90:141
material::damage_none_id
@, public damage_none_id
Definition: material.f90:87
thermal_conduction
material subroutine for temperature evolution from heat conduction
Definition: thermal_conduction.f90:9
constitutive
elasticity, plasticity, internal microstructure state
Definition: constitutive.f90:10
homogenization::materialpoint_subf0
real(preal), dimension(:,:,:,:), allocatable materialpoint_subf0
def grad of IP at beginning of homogenization increment
Definition: homogenization.f90:43
crystallite::crystallite_partionedf0
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedf0
def grad at start of homog inc
Definition: crystallite.f90:52
homogenization::tnumerics
Definition: homogenization.f90:56
crystallite::crystallite_fp
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_fp
current plastic def grad (end of converged time step)
Definition: crystallite.f90:52
results
Definition: results.f90:11
config::config_numerics
type(tpartitionedstringlist), public, protected config_numerics
Definition: config.f90:30
material::material_homogenizationat
integer, dimension(:), allocatable, public, protected material_homogenizationat
homogenization ID of each element (copy of discretization_homogenizationAt)
Definition: material.f90:128
material::plasticstate
type(tplasticstate), dimension(:), allocatable, public plasticstate
Definition: material.f90:137
math
Mathematical library, including random number generation and tensor representations.
Definition: math.f90:12
fesolving::terminallyill
logical terminallyill
at least one material point is terminally ill
Definition: FEsolving.f90:15
crystallite::crystallite_f0
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_f0
def grad at start of FE inc
Definition: crystallite.f90:43
crystallite::converged
logical pure function converged(residuum, state, atol)
determines whether a point is converged
Definition: crystallite.f90:1686
crystallite::crystallite_partionedli0
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedli0
intermediate velocity grad at start of homog inc
Definition: crystallite.f90:52
discretization::discretization_nelem
integer, public, protected discretization_nelem
Definition: discretization.f90:17
damage_nonlocal::damage_nonlocal_results
subroutine, public damage_nonlocal_results(homog, group)
writes results to HDF5 output file
Definition: damage_nonlocal.f90:207
debug::debug_e
integer, public, protected debug_e
Definition: debug.f90:51
homogenization::materialpoint_doneandhappy
logical, dimension(:,:,:), allocatable materialpoint_doneandhappy
Definition: homogenization.f90:53
crystallite::crystallite_requested
logical, dimension(:,:,:), allocatable, public crystallite_requested
used by upper level (homogenization) to request crystallite calculation
Definition: crystallite.f90:74
results::results_closegroup
subroutine, public results_closegroup(group_id)
close a group
Definition: results.f90:166
material::material_homogenizationmemberat
integer, dimension(:,:), allocatable, target, public material_homogenizationmemberat
position of the element within its homogenization instance
Definition: material.f90:130
crystallite::crystallite_partionedf
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedf
def grad to be reached at end of homog inc
Definition: crystallite.f90:52
crystallite::crystallite_orientations
subroutine, public crystallite_orientations
calculates orientations
Definition: crystallite.f90:573
fesolving::fesolving_execelem
integer, dimension(2) fesolving_execelem
for ping-pong scheme always whole range, otherwise one specific element
Definition: FEsolving.f90:18
thermal_conduction::thermal_conduction_init
subroutine, public thermal_conduction_init
module initialization
Definition: thermal_conduction.f90:47
material::homogenization_type
integer(kind(homogenization_undefined_id)), dimension(:), allocatable, public, protected homogenization_type
type of each homogenization
Definition: material.f90:98
material::phase_nsources
integer, dimension(:), allocatable, public, protected phase_nsources
number of source mechanisms active in each phase
Definition: material.f90:113
material::damage_type
integer(kind(damage_none_id)), dimension(:), allocatable, public, protected damage_type
nonlocal damage model
Definition: material.f90:96
crystallite::crystallite_lp
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_lp
current plastic velocitiy grad (end of converged time step)
Definition: crystallite.f90:52
debug::debug_g
integer, public, protected debug_g
Definition: debug.f90:51
damage_local::damage_local_updatestate
logical function, dimension(2), public damage_local_updatestate(subdt, ip, el)
calculates local change in damage field
Definition: damage_local.f90:78
numerics
Managing of parameters related to numerics.
Definition: numerics.f90:10
homogenization::materialpoint_requested
logical, dimension(:,:), allocatable materialpoint_requested
Definition: homogenization.f90:50
fesolving
global variables for flow control
Definition: FEsolving.f90:10
thermal_adiabatic::thermal_adiabatic_updatestate
logical function, dimension(2), public thermal_adiabatic_updatestate(subdt, ip, el)
calculates adiabatic change in temperature based on local heat generation model
Definition: thermal_adiabatic.f90:85
material::homogenization_rgc_id
@, public homogenization_rgc_id
Definition: material.f90:87
material::homogenization_none_id
@, public homogenization_none_id
Definition: material.f90:87
material::thermal_conduction_id
@, public thermal_conduction_id
Definition: material.f90:87
math::math_i3
real(preal), dimension(3, 3), parameter math_i3
3x3 Identity
Definition: math.f90:32