|
DAMASK with grid solvers
Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
|
|
Go to the documentation of this file. 1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization.f90"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization.f90"
36 real(
preal),
dimension(:,:,:,:),
allocatable,
public :: &
38 materialpoint_f, & !< def grad of IP to be reached at end of FE increment
40 real(
preal),
dimension(:,:,:,:,:,:),
allocatable,
public :: &
43 real(
preal),
dimension(:,:,:,:),
allocatable :: &
46 real(
preal),
dimension(:,:),
allocatable :: &
50 logical,
dimension(:,:),
allocatable :: &
53 logical,
dimension(:,:,:),
allocatable :: &
60 substepminhomog, & !< minimum (relative) size of sub-step allowed during cutback in homogenization
61 substepsizehomog, & !< size of first substep when cutback in homogenization
69 module subroutine mech_none_init
70 end subroutine mech_none_init
72 module subroutine mech_isostrain_init
73 end subroutine mech_isostrain_init
75 module subroutine mech_rgc_init
76 end subroutine mech_rgc_init
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
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) :: &
90 end subroutine mech_rgc_partitiondeformation
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
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
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
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
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
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
124 end function mech_rgc_updatestate
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
135 homogenization_init, &
136 materialpoint_stressanditstangent, &
137 homogenization_results
145 subroutine homogenization_init
165 materialpoint_f = materialpoint_f0
176 write(6,
'(/,a)')
' <<<+- homogenization init -+>>>';
flush(6)
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)
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')
206 end subroutine homogenization_init
212 subroutine materialpoint_stressanditstangent(updateJaco,dt)
214 real(
preal),
intent(in) :: dt
215 logical,
intent(in) :: updatejaco
220 i, & !< integration point number
221 e, & !< element number
225 # 231 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization.f90"
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
254 materialpoint_converged(i,e) = .false.
255 materialpoint_requested(i,e) = .true.
281 converged:
if (materialpoint_converged(i,e))
then
282 # 296 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization.f90"
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))
290 steppingneeded:
if (materialpoint_substep(i,e) >
num%subStepMinHomog)
then
330 materialpoint_subf0(1:3,1:3,i,e) = materialpoint_subf(1:3,1:3,i,e)
335 if ( (myngrains == 1 .and. materialpoint_substep(i,e) <= 1.0 ) .or. &
336 num%subStepSizeHomog * materialpoint_substep(i,e) <=
num%subStepMinHomog )
then
341 write(6,*)
'Integration point ', i,
' at element ', e,
' terminally ill'
346 materialpoint_substep(i,e) =
num%subStepSizeHomog * materialpoint_substep(i,e)
348 # 370 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/homogenization.f90"
352 if (materialpoint_substep(i,e) < 1.0_preal)
then
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.]
393 enddo elementlooping1
396 niterationmpstate = 0
402 niterationmpstate <
num%nMPstate)
403 niterationmpstate = niterationmpstate + 1
413 if ( materialpoint_requested(i,e) .and. &
414 .not. materialpoint_doneandhappy(1,i,e))
then
415 call partitiondeformation(i,e)
422 enddo elementlooping2
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.]
442 materialpoint_doneandhappy(1:2,i,e) = updatestate(i,e)
443 materialpoint_converged(i,e) = all(materialpoint_doneandhappy(1:2,i,e))
447 enddo elementlooping3
450 enddo convergencelooping
452 niterationhomog = niterationhomog + 1
463 call averagestressanditstangent(i,e)
465 enddo elementlooping4
468 write(6,
'(/,a,/)')
'<< HOMOG >> Material Point terminally ill'
471 end subroutine materialpoint_stressanditstangent
477 subroutine partitiondeformation(ip,el)
479 integer,
intent(in) :: &
480 ip, & !< integration point
489 call mech_isostrain_partitiondeformation(&
491 materialpoint_subf(1:3,1:3,ip,el))
494 call mech_rgc_partitiondeformation(&
496 materialpoint_subf(1:3,1:3,ip,el),&
499 end select chosenhomogenization
501 end subroutine partitiondeformation
508 function updatestate(ip,el)
510 integer,
intent(in) :: &
511 ip, & !< integration point
513 logical,
dimension(2) :: updatestate
523 materialpoint_subf(1:3,1:3,ip,el),&
524 materialpoint_subdt(ip,el), &
528 end select chosenhomogenization
537 end select chosenthermal
546 end select chosendamage
548 end function updatestate
554 subroutine averagestressanditstangent(ip,el)
556 integer,
intent(in) :: &
557 ip, & !< integration point
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)
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),&
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),&
580 end select chosenhomogenization
582 end subroutine averagestressanditstangent
588 subroutine homogenization_results
593 character(len=pStringLen) :: group_base,group
601 group = trim(group_base)//
'/generic'
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)
617 group = trim(group_base)//
'/damage'
619 select case(damage_type(p))
620 case(damage_local_id)
622 case(damage_nonlocal_id)
626 group = trim(group_base)//
'/thermal'
628 select case(thermal_type(p))
629 case(thermal_adiabatic_id)
631 case(thermal_conduction_id)
637 end subroutine homogenization_results
type(tsourcestate), dimension(:), allocatable, public sourcestate
integer, dimension(2) fesolving_execip
for ping-pong scheme always range to max IP, otherwise one specific IP
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_s0
2nd Piola-Kirchhoff stress vector at start of FE inc
subroutine, public thermal_adiabatic_init
module initialization
integer, dimension(:,:), allocatable, public, protected material_phaseat
phase ID of each element
real(preal), dimension(:,:,:), allocatable, public crystallite_dt
requested time increment of each grain
real(preal), dimension(:,:), allocatable materialpoint_subdt
material subroutine for non-locally evolving damage field
real(preal), dimension(:,:,:,:,:,:), allocatable, public materialpoint_dpdf
tangent of first P–K stress at IP
logical function, dimension(discretization_nip, discretization_nelem), public crystallite_stress(dummyArgumentToPreventInternalCompilerErrorWithGCC)
calculate stress (P)
real(preal), dimension(:,:,:,:), allocatable, public materialpoint_p
first P–K stress of IP
logical, dimension(:,:), allocatable materialpoint_converged
subroutine thermal_isothermal_init
allocates all neccessary fields, reads information from material configuration file
real(preal), dimension(:,:), allocatable materialpoint_substep
subroutine, public damage_nonlocal_init
module initialization
integer, dimension(:,:,:), allocatable, public, protected material_phasememberat
position of the element within its phase instance
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
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedfp0
plastic def grad at start of homog inc
real(preal), dimension(:,:,:,:), allocatable, public materialpoint_f0
def grad of IP at start of FE increment
integer, dimension(debug_maxntype+2), public, protected debug_level
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_p
1st Piola-Kirchhoff stress per grain
@, public thermal_adiabatic_id
real(preal), dimension(:,:), allocatable materialpoint_subfrac
@, public damage_nonlocal_id
material subroutine for locally evolving damage field
material subroutine for adiabatic temperature evolution
real(preal), dimension(:,:,:,:), allocatable, public materialpoint_f
def grad of IP to be reached at end of FE increment
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_fi0
intermediate def grad at start of FE inc
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_fi
current intermediate def grad (end of converged time step)
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_li
current intermediate velocitiy grad (end of converged time step)
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_s
current 2nd Piola-Kirchhoff stress vector (end of converged time step)
Parses material config file, either solverJobName.materialConfig or material.config.
subroutine, public thermal_conduction_results(homog, group)
writes results to HDF5 output file
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_fp0
plastic def grad at start of FE inc
Reads in the material configuration from file.
subroutine, public crystallite_stresstangent
calculate tangent (dPdF)
type(tstate), dimension(:), allocatable, public damagestate
real(preal), dimension(:,:,:,:), allocatable materialpoint_subf
def grad of IP to be reached at end of homog inc
subroutine, public config_deallocate(what)
deallocates the linked lists that store the content of the configuration files
integer, dimension(:), allocatable, public, protected homogenization_ngrains
number of grains in each homogenization
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_lp0
plastic velocitiy grad at start of FE inc
real(preal), dimension(:,:,:,:,:,:,:), allocatable, public, protected crystallite_dpdf
current individual dPdF per grain (end of converged time step)
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedlp0
plastic velocity grad at start of homog inc
@, public thermal_isothermal_id
subroutine, public damage_local_init
module initialization
@, public damage_local_id
crystallite state integration functions and reporting of results
integer(kind(thermal_isothermal_id)), dimension(:), allocatable, public, protected thermal_type
thermal transport model
setting precision for real and int type
integer, dimension(:), allocatable, public, protected homogenization_typeinstance
instance of particular type of each homogenization
type(tstate), dimension(:), allocatable, public homogstate
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_li0
intermediate velocitiy grad at start of FE inc
integer, public, protected discretization_nip
@, public homogenization_isostrain_id
input/output functions, partly depending on chosen solver
integer, parameter, public debug_homogenization
homogenization manager, organizing deformation partitioning and stress homogenization
integer, parameter, public debug_levelbasic
integer, parameter preal
number with 15 significant digits, up to 1e+-307 (typically 64 bit)
subroutine, public thermal_adiabatic_results(homog, group)
writes results to HDF5 output file
subroutine, public damage_local_results(homog, group)
writes results to HDF5 output file
Reading in and interpretating the debugging settings for the various modules.
subroutine damage_none_init
allocates all neccessary fields, reads information from material configuration file
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partioneds0
2nd Piola-Kirchhoff stress vector at start of homog inc
material subroutine for isothermal temperature field
material subroutine for constant damage field
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedfi0
intermediate def grad at start of homog inc
integer(hid_t) function, public results_addgroup(groupName)
adds a new group to the results file
character(len=pstringlen), dimension(:), allocatable, public, protected config_name_homogenization
name of each homogenization
type(tstate), dimension(:), allocatable, public thermalstate
material subroutine for temperature evolution from heat conduction
elasticity, plasticity, internal microstructure state
real(preal), dimension(:,:,:,:), allocatable materialpoint_subf0
def grad of IP at beginning of homogenization increment
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedf0
def grad at start of homog inc
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_fp
current plastic def grad (end of converged time step)
type(tpartitionedstringlist), public, protected config_numerics
integer, dimension(:), allocatable, public, protected material_homogenizationat
homogenization ID of each element (copy of discretization_homogenizationAt)
type(tplasticstate), dimension(:), allocatable, public plasticstate
Mathematical library, including random number generation and tensor representations.
logical terminallyill
at least one material point is terminally ill
real(preal), dimension(:,:,:,:,:), allocatable, public, protected crystallite_f0
def grad at start of FE inc
logical pure function converged(residuum, state, atol)
determines whether a point is converged
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedli0
intermediate velocity grad at start of homog inc
integer, public, protected discretization_nelem
subroutine, public damage_nonlocal_results(homog, group)
writes results to HDF5 output file
integer, public, protected debug_e
logical, dimension(:,:,:), allocatable materialpoint_doneandhappy
logical, dimension(:,:,:), allocatable, public crystallite_requested
used by upper level (homogenization) to request crystallite calculation
subroutine, public results_closegroup(group_id)
close a group
integer, dimension(:,:), allocatable, target, public material_homogenizationmemberat
position of the element within its homogenization instance
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_partionedf
def grad to be reached at end of homog inc
subroutine, public crystallite_orientations
calculates orientations
integer, dimension(2) fesolving_execelem
for ping-pong scheme always whole range, otherwise one specific element
subroutine, public thermal_conduction_init
module initialization
integer(kind(homogenization_undefined_id)), dimension(:), allocatable, public, protected homogenization_type
type of each homogenization
integer, dimension(:), allocatable, public, protected phase_nsources
number of source mechanisms active in each phase
integer(kind(damage_none_id)), dimension(:), allocatable, public, protected damage_type
nonlocal damage model
real(preal), dimension(:,:,:,:,:), allocatable, public crystallite_lp
current plastic velocitiy grad (end of converged time step)
integer, public, protected debug_g
logical function, dimension(2), public damage_local_updatestate(subdt, ip, el)
calculates local change in damage field
Managing of parameters related to numerics.
logical, dimension(:,:), allocatable materialpoint_requested
global variables for flow control
logical function, dimension(2), public thermal_adiabatic_updatestate(subdt, ip, el)
calculates adiabatic change in temperature based on local heat generation model
@, public homogenization_rgc_id
@, public homogenization_none_id
@, public thermal_conduction_id
real(preal), dimension(3, 3), parameter math_i3
3x3 Identity