DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
source_damage_anisoDuctile.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_anisoDuctile.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_anisoDuctile.f90"
5 !--------------------------------------------------------------------------------------------------
10 !--------------------------------------------------------------------------------------------------
12  use prec
13  use debug
14  use io
15  use math
16  use discretization
17  use material
18  use config
19  use results
20 
21  implicit none
22  private
23 
24  integer, dimension(:), allocatable :: &
25  source_damage_anisoductile_offset, & !< which source is my current damage mechanism?
27 
28  type, private :: tparameters
29  real(preal) :: &
30  n
31  real(preal), dimension(:), allocatable :: &
32  critplasticstrain
33  character(len=pStringLen), allocatable, dimension(:) :: &
34  output
35  end type tparameters
36 
37  type(tparameters), dimension(:), allocatable, private :: param
38 
39 
40  public :: &
45 
46 contains
47 
48 
49 !--------------------------------------------------------------------------------------------------
52 !--------------------------------------------------------------------------------------------------
54 
55  integer :: ninstance,sourceoffset,nipcmyphase,p
56  integer, dimension(:), allocatable :: n_sl
57  character(len=pStringLen) :: extmsg = ''
58 
59  write(6,'(/,a)') ' <<<+- source_'//source_damage_anisoductile_label//' init -+>>>'; flush(6)
60 
61  ninstance = count(phase_source == source_damage_anisoductile_id)
63  write(6,'(a16,1x,i5,/)') '# instances:',ninstance
64 
65  allocate(source_damage_anisoductile_offset(size(config_phase)), source=0)
66  allocate(source_damage_anisoductile_instance(size(config_phase)), source=0)
67  allocate(param(ninstance))
68 
69  do p = 1, size(config_phase)
71  do sourceoffset = 1, phase_nsources(p)
72  if (phase_source(sourceoffset,p) == source_damage_anisoductile_id) then
73  source_damage_anisoductile_offset(p) = sourceoffset
74  exit
75  endif
76  enddo
77 
78  if (all(phase_source(:,p) /= source_damage_anisoductile_id)) cycle
79  associate(prm => param(source_damage_anisoductile_instance(p)), &
80  config => config_phase(p))
81 
82  prm%output = config%getStrings('(output)',defaultval=emptystringarray)
83 
84  n_sl = config%getInts('nslip',defaultval=emptyintarray)
85  prm%n = config%getFloat('anisoductile_ratesensitivity')
86  prm%critPlasticStrain = config%getFloats('anisoductile_criticalplasticstrain',requiredsize=size(n_sl))
87 
88  ! expand: family => system
89  prm%critPlasticStrain = math_expand(prm%critPlasticStrain,n_sl)
90 
91  ! sanity checks
92  if (prm%n <= 0.0_preal) extmsg = trim(extmsg)//' anisoductile_ratesensitivity'
93  if (any(prm%critPlasticStrain < 0.0_preal)) extmsg = trim(extmsg)//' anisoductile_criticalplasticstrain'
94 
95  nipcmyphase=count(material_phaseat==p) * discretization_nip
96  call material_allocatesourcestate(p,sourceoffset,nipcmyphase,1,1,0)
97  sourcestate(p)%p(sourceoffset)%atol = config%getFloat('anisoductile_atol',defaultval=1.0e-3_preal)
98  if(any(sourcestate(p)%p(sourceoffset)%atol < 0.0_preal)) extmsg = trim(extmsg)//' anisoductile_atol'
99 
100  end associate
101 
102 !--------------------------------------------------------------------------------------------------
103 ! exit if any parameter is out of range
104  if (extmsg /= '') call io_error(211,ext_msg=trim(extmsg)//'('//source_damage_anisoductile_label//')')
105 
106 enddo
107 
108 end subroutine source_damage_anisoductile_init
109 
110 
111 !--------------------------------------------------------------------------------------------------
113 !--------------------------------------------------------------------------------------------------
114 subroutine source_damage_anisoductile_dotstate(ipc, ip, el)
115 
116  integer, intent(in) :: &
117  ipc, & !< component-ID of integration point
118  ip, & !< integration point
119  el
120 
121  integer :: &
122  phase, &
123  constituent, &
124  sourceoffset, &
125  damageoffset, &
126  homog
127 
128  phase = material_phaseat(ipc,el)
129  constituent = material_phasememberat(ipc,ip,el)
130  sourceoffset = source_damage_anisoductile_offset(phase)
131  homog = material_homogenizationat(el)
132  damageoffset = damagemapping(homog)%p(ip,el)
133 
134  associate(prm => param(source_damage_anisoductile_instance(phase)))
135  sourcestate(phase)%p(sourceoffset)%dotState(1,constituent) &
136  = sum(plasticstate(phase)%slipRate(:,constituent)/(damage(homog)%p(damageoffset)**prm%n)/prm%critPlasticStrain)
137  end associate
138 
140 
141 
142 !--------------------------------------------------------------------------------------------------
144 !--------------------------------------------------------------------------------------------------
145 subroutine source_damage_anisoductile_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
146 
147  integer, intent(in) :: &
148  phase, &
149  constituent
150  real(preal), intent(in) :: &
151  phi
152  real(preal), intent(out) :: &
153  localphidot, &
154  dlocalphidot_dphi
155 
156  integer :: &
157  sourceoffset
158 
159  sourceoffset = source_damage_anisoductile_offset(phase)
160 
161  dlocalphidot_dphi = -sourcestate(phase)%p(sourceoffset)%state(1,constituent)
162 
163  localphidot = 1.0_preal &
164  + dlocalphidot_dphi*phi
165 
167 
168 
169 !--------------------------------------------------------------------------------------------------
171 !--------------------------------------------------------------------------------------------------
172 subroutine source_damage_anisoductile_results(phase,group)
173 
174  integer, intent(in) :: phase
175  character(len=*), intent(in) :: group
176 
177  integer :: o
178 
179  associate(prm => param(source_damage_anisoductile_instance(phase)), &
180  stt => sourcestate(phase)%p(source_damage_anisoductile_offset(phase))%state)
181  outputsloop: do o = 1,size(prm%output)
182  select case(trim(prm%output(o)))
183  case ('anisoductile_drivingforce')
184  call results_writedataset(group,stt,'tbd','driving force','tbd')
185  end select
186  enddo outputsloop
187  end associate
188 
190 
material::sourcestate
type(tsourcestate), dimension(:), allocatable, public sourcestate
Definition: material.f90:139
material::material_phaseat
integer, dimension(:,:), allocatable, public, protected material_phaseat
phase ID of each element
Definition: material.f90:132
math::math_expand
pure real(preal) function, dimension(sum(how)) math_expand(what, how)
vector expansion
Definition: math.f90:207
source_damage_anisoductile::source_damage_anisoductile_instance
integer, dimension(:), allocatable source_damage_anisoductile_instance
instance of damage source mechanism
Definition: source_damage_anisoDuctile.f90:24
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
config::config_phase
type(tpartitionedstringlist), dimension(:), allocatable, public, protected config_phase
Definition: config.f90:23
debug::debug_level
integer, dimension(debug_maxntype+2), public, protected debug_level
Definition: debug.f90:48
prec::emptystringarray
character(len=pstringlen), dimension(0), parameter emptystringarray
Definition: prec.f90:88
material
Parses material config file, either solverJobName.materialConfig or material.config.
Definition: material.f90:11
source_damage_anisoductile::source_damage_anisoductile_results
subroutine, public source_damage_anisoductile_results(phase, group)
writes results to HDF5 output file
Definition: source_damage_anisoDuctile.f90:173
config
Reads in the material configuration from file.
Definition: config.f90:13
prec
setting precision for real and int type
Definition: prec.f90:13
source_damage_anisoductile::source_damage_anisoductile_init
subroutine, public source_damage_anisoductile_init
module initialization
Definition: source_damage_anisoDuctile.f90:54
discretization
spatial discretization
Definition: discretization.f90:9
discretization::discretization_nip
integer, public, protected discretization_nip
Definition: discretization.f90:17
source_damage_anisoductile::param
type(tparameters), dimension(:), allocatable, private param
containers of constitutive parameters (len Ninstance)
Definition: source_damage_anisoDuctile.f90:37
io
input/output functions, partly depending on chosen solver
Definition: IO.f90:12
source_damage_anisoductile
material subroutine incorporating anisotropic ductile damage source mechanism
Definition: source_damage_anisoDuctile.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
debug
Reading in and interpretating the debugging settings for the various modules.
Definition: debug.f90:12
source_damage_anisoductile::source_damage_anisoductile_offset
integer, dimension(:), allocatable source_damage_anisoductile_offset
which source is my current damage mechanism?
Definition: source_damage_anisoDuctile.f90:24
material::material_allocatesourcestate
subroutine, public material_allocatesourcestate(phase, of, NipcMyPhase, sizeState, sizeDotState, sizeDeltaState)
allocates the source state of a phase
Definition: material.f90:750
prec::emptyintarray
integer, dimension(0), parameter emptyintarray
Definition: prec.f90:84
results
Definition: results.f90:11
source_damage_anisoductile::source_damage_anisoductile_getrateanditstangent
subroutine, public source_damage_anisoductile_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
returns local part of nonlocal damage driving force
Definition: source_damage_anisoDuctile.f90:146
math
Mathematical library, including random number generation and tensor representations.
Definition: math.f90:12
material::source_damage_anisoductile_label
character(len= *), parameter, public source_damage_anisoductile_label
Definition: material.f90:25
debug::debug_constitutive
integer, parameter, public debug_constitutive
stores debug level for constitutive part of DAMASK bitwise coded
Definition: debug.f90:32
material::phase_source
integer(kind(source_undefined_id)), dimension(:,:), allocatable, public, protected phase_source
active sources mechanisms of each phase
Definition: material.f90:105
material::phase_nsources
integer, dimension(:), allocatable, public, protected phase_nsources
number of source mechanisms active in each phase
Definition: material.f90:113
source_damage_anisoductile::tparameters
container type for internal constitutive parameters
Definition: source_damage_anisoDuctile.f90:28
source_damage_anisoductile::source_damage_anisoductile_dotstate
subroutine, public source_damage_anisoductile_dotstate(ipc, ip, el)
calculates derived quantities from state
Definition: source_damage_anisoDuctile.f90:115
material::source_damage_anisoductile_id
@, public source_damage_anisoductile_id
Definition: material.f90:87