DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
source_damage_isoDuctile.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_isoDuctile.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_isoDuctile.f90"
5 !--------------------------------------------------------------------------------------------------
10 !--------------------------------------------------------------------------------------------------
12  use prec
13  use debug
14  use io
15  use discretization
16  use material
17  use config
18  use results
19 
20  implicit none
21  private
22 
23  integer, dimension(:), allocatable :: &
24  source_damage_isoductile_offset, & !< which source is my current damage mechanism?
26 
27  type, private :: tparameters
28  real(preal) :: &
29  critplasticstrain, &
30  n
31  character(len=pStringLen), allocatable, dimension(:) :: &
32  output
33  end type tparameters
34 
35  type(tparameters), dimension(:), allocatable, private :: param
36 
37 
38  public :: &
43 
44 contains
45 
46 
47 !--------------------------------------------------------------------------------------------------
50 !--------------------------------------------------------------------------------------------------
52 
53  integer :: ninstance,sourceoffset,nipcmyphase,p
54  character(len=pStringLen) :: extmsg = ''
55 
56  write(6,'(/,a)') ' <<<+- source_'//source_damage_isoductile_label//' init -+>>>'; flush(6)
57 
58  ninstance = count(phase_source == source_damage_isoductile_id)
60  write(6,'(a16,1x,i5,/)') '# instances:',ninstance
61 
62  allocate(source_damage_isoductile_offset(size(config_phase)), source=0)
63  allocate(source_damage_isoductile_instance(size(config_phase)), source=0)
64  allocate(param(ninstance))
65 
66  do p = 1, size(config_phase)
68  do sourceoffset = 1, phase_nsources(p)
69  if (phase_source(sourceoffset,p) == source_damage_isoductile_id) then
70  source_damage_isoductile_offset(p) = sourceoffset
71  exit
72  endif
73  enddo
74 
75  if (all(phase_source(:,p) /= source_damage_isoductile_id)) cycle
76  associate(prm => param(source_damage_isoductile_instance(p)), &
77  config => config_phase(p))
78 
79  prm%output = config%getStrings('(output)',defaultval=emptystringarray)
80 
81  prm%N = config%getFloat('isoductile_ratesensitivity')
82  prm%critPlasticStrain = config%getFloat('isoductile_criticalplasticstrain')
83 
84  ! sanity checks
85  if (prm%N <= 0.0_preal) extmsg = trim(extmsg)//' isoductile_ratesensitivity'
86  if (prm%critPlasticStrain <= 0.0_preal) extmsg = trim(extmsg)//' isoductile_criticalplasticstrain'
87 
88  nipcmyphase=count(material_phaseat==p) * discretization_nip
89  call material_allocatesourcestate(p,sourceoffset,nipcmyphase,1,1,0)
90  sourcestate(p)%p(sourceoffset)%atol = config%getFloat('isoductile_atol',defaultval=1.0e-3_preal)
91  if(any(sourcestate(p)%p(sourceoffset)%atol < 0.0_preal)) extmsg = trim(extmsg)//' isoductile_atol'
92 
93  end associate
94 
95 !--------------------------------------------------------------------------------------------------
96 ! exit if any parameter is out of range
97  if (extmsg /= '') call io_error(211,ext_msg=trim(extmsg)//'('//source_damage_isoductile_label//')')
98 
99 enddo
100 
101 end subroutine source_damage_isoductile_init
102 
103 
104 !--------------------------------------------------------------------------------------------------
106 !--------------------------------------------------------------------------------------------------
107 subroutine source_damage_isoductile_dotstate(ipc, ip, el)
108 
109  integer, intent(in) :: &
110  ipc, & !< component-ID of integration point
111  ip, & !< integration point
112  el
113 
114  integer :: &
115  phase, &
116  constituent, &
117  sourceoffset, &
118  damageoffset, &
119  homog
120 
121  phase = material_phaseat(ipc,el)
122  constituent = material_phasememberat(ipc,ip,el)
123  sourceoffset = source_damage_isoductile_offset(phase)
124  homog = material_homogenizationat(el)
125  damageoffset = damagemapping(homog)%p(ip,el)
126 
127  associate(prm => param(source_damage_isoductile_instance(phase)))
128  sourcestate(phase)%p(sourceoffset)%dotState(1,constituent) = &
129  sum(plasticstate(phase)%slipRate(:,constituent))/(damage(homog)%p(damageoffset)**prm%N)/prm%critPlasticStrain
130  end associate
131 
133 
134 
135 !--------------------------------------------------------------------------------------------------
137 !--------------------------------------------------------------------------------------------------
138 subroutine source_damage_isoductile_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
139 
140  integer, intent(in) :: &
141  phase, &
142  constituent
143  real(preal), intent(in) :: &
144  phi
145  real(preal), intent(out) :: &
146  localphidot, &
147  dlocalphidot_dphi
148 
149  integer :: &
150  sourceoffset
151 
152  sourceoffset = source_damage_isoductile_offset(phase)
153 
154  dlocalphidot_dphi = -sourcestate(phase)%p(sourceoffset)%state(1,constituent)
155 
156  localphidot = 1.0_preal &
157  + dlocalphidot_dphi*phi
158 
160 
161 
162 !--------------------------------------------------------------------------------------------------
164 !--------------------------------------------------------------------------------------------------
165 subroutine source_damage_isoductile_results(phase,group)
166 
167  integer, intent(in) :: phase
168  character(len=*), intent(in) :: group
169 
170  integer :: o
171 
172  associate(prm => param(source_damage_isoductile_instance(phase)), &
173  stt => sourcestate(phase)%p(source_damage_isoductile_offset(phase))%state)
174  outputsloop: do o = 1,size(prm%output)
175  select case(trim(prm%output(o)))
176  case ('isoductile_drivingforce')
177  call results_writedataset(group,stt,'tbd','driving force','tbd')
178  end select
179  enddo outputsloop
180  end associate
181 
183 
184 end module source_damage_isoductile
material::sourcestate
type(tsourcestate), dimension(:), allocatable, public sourcestate
Definition: material.f90:139
source_damage_isoductile::source_damage_isoductile_dotstate
subroutine, public source_damage_isoductile_dotstate(ipc, ip, el)
calculates derived quantities from state
Definition: source_damage_isoDuctile.f90:108
material::material_phaseat
integer, dimension(:,:), allocatable, public, protected material_phaseat
phase ID of each element
Definition: material.f90:132
material::source_damage_isoductile_id
@, public source_damage_isoductile_id
Definition: material.f90:87
source_damage_isoductile
material subroutine incoprorating isotropic ductile damage source mechanism
Definition: source_damage_isoDuctile.f90:11
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
source_damage_isoductile::source_damage_isoductile_results
subroutine, public source_damage_isoductile_results(phase, group)
writes results to HDF5 output file
Definition: source_damage_isoDuctile.f90:166
prec::emptystringarray
character(len=pstringlen), dimension(0), parameter emptystringarray
Definition: prec.f90:88
source_damage_isoductile::source_damage_isoductile_init
subroutine, public source_damage_isoductile_init
module initialization
Definition: source_damage_isoDuctile.f90:52
material
Parses material config file, either solverJobName.materialConfig or material.config.
Definition: material.f90:11
config
Reads in the material configuration from file.
Definition: config.f90:13
prec
setting precision for real and int type
Definition: prec.f90:13
discretization
spatial discretization
Definition: discretization.f90:9
source_damage_isoductile::source_damage_isoductile_instance
integer, dimension(:), allocatable source_damage_isoductile_instance
instance of damage source mechanism
Definition: source_damage_isoDuctile.f90:23
discretization::discretization_nip
integer, public, protected discretization_nip
Definition: discretization.f90:17
io
input/output functions, partly depending on chosen solver
Definition: IO.f90:12
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
material::material_allocatesourcestate
subroutine, public material_allocatesourcestate(phase, of, NipcMyPhase, sizeState, sizeDotState, sizeDeltaState)
allocates the source state of a phase
Definition: material.f90:750
source_damage_isoductile::param
type(tparameters), dimension(:), allocatable, private param
containers of constitutive parameters (len Ninstance)
Definition: source_damage_isoDuctile.f90:35
source_damage_isoductile::source_damage_isoductile_getrateanditstangent
subroutine, public source_damage_isoductile_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
returns local part of nonlocal damage driving force
Definition: source_damage_isoDuctile.f90:139
source_damage_isoductile::tparameters
container type for internal constitutive parameters
Definition: source_damage_isoDuctile.f90:27
results
Definition: results.f90:11
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
material::source_damage_isoductile_label
character(len= *), parameter, public source_damage_isoductile_label
Definition: material.f90:25
source_damage_isoductile::source_damage_isoductile_offset
integer, dimension(:), allocatable source_damage_isoductile_offset
which source is my current damage mechanism?
Definition: source_damage_isoDuctile.f90:23