DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
kinematics_slipplane_opening.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/kinematics_slipplane_opening.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/kinematics_slipplane_opening.f90"
5 !--------------------------------------------------------------------------------------------------
10 !--------------------------------------------------------------------------------------------------
12  use prec
13  use config
14  use io
15  use debug
16  use math
17  use lattice
18  use material
19 
20  implicit none
21  private
22 
23  integer, dimension(:), allocatable :: kinematics_slipplane_opening_instance
24 
25  type :: tparameters
26  integer :: &
27  sum_n_sl
28  real(preal) :: &
29  sdot0, &
30  n
31  real(preal), dimension(:), allocatable :: &
32  critload
33  real(preal), dimension(:,:,:), allocatable :: &
34  p_d, &
35  p_t, &
36  p_n
37  end type tparameters
38 
39  type(tparameters), dimension(:), allocatable :: param
40 
41  public :: &
44 
45 contains
46 
47 
48 !--------------------------------------------------------------------------------------------------
51 !--------------------------------------------------------------------------------------------------
53 
54  integer :: ninstance,p,i
55  character(len=pStringLen) :: extmsg = ''
56  integer, dimension(:), allocatable :: n_sl
57  real(preal), dimension(:,:), allocatable :: d,n,t
58 
59  write(6,'(/,a)') ' <<<+- kinematics_'//kinematics_slipplane_opening_label//' init -+>>>'; flush(6)
60 
63  write(6,'(a16,1x,i5,/)') '# instances:',ninstance
64 
65  allocate(kinematics_slipplane_opening_instance(size(config_phase)), source=0)
66  allocate(param(ninstance))
67 
68  do p = 1, size(config_phase)
71  associate(prm => param(kinematics_slipplane_opening_instance(p)), &
72  config => config_phase(p))
73 
74  prm%sdot0 = config%getFloat('anisoductile_sdot0')
75  prm%n = config%getFloat('anisoductile_ratesensitivity')
76  n_sl = config%getInts('nslip')
77  prm%sum_N_sl = sum(abs(n_sl))
78 
79  d = lattice_slip_direction(n_sl,config%getString('lattice_structure'),&
80  config%getFloat('c/a',defaultval=0.0_preal))
81  t = lattice_slip_transverse(n_sl,config%getString('lattice_structure'),&
82  config%getFloat('c/a',defaultval=0.0_preal))
83  n = lattice_slip_normal(n_sl,config%getString('lattice_structure'),&
84  config%getFloat('c/a',defaultval=0.0_preal))
85  allocate(prm%P_d(3,3,size(d,2)),prm%P_t(3,3,size(t,2)),prm%P_n(3,3,size(n,2)))
86 
87  do i=1, size(n,2)
88  prm%P_d(1:3,1:3,i) = math_outer(d(1:3,i), n(1:3,i))
89  prm%P_t(1:3,1:3,i) = math_outer(t(1:3,i), n(1:3,i))
90  prm%P_n(1:3,1:3,i) = math_outer(n(1:3,i), n(1:3,i))
91  enddo
92 
93  prm%critLoad = config%getFloats('anisoductile_criticalload',requiredsize=size(n_sl))
94 
95  ! expand: family => system
96  prm%critLoad = math_expand(prm%critLoad,n_sl)
97 
98  ! sanity checks
99  if (prm%n <= 0.0_preal) extmsg = trim(extmsg)//' anisoDuctile_n'
100  if (prm%sdot0 <= 0.0_preal) extmsg = trim(extmsg)//' anisoDuctile_sdot0'
101  if (any(prm%critLoad < 0.0_preal)) extmsg = trim(extmsg)//' anisoDuctile_critLoad'
102 
103 !--------------------------------------------------------------------------------------------------
104 ! exit if any parameter is out of range
105  if (extmsg /= '') call io_error(211,ext_msg=trim(extmsg)//'('//kinematics_slipplane_opening_label//')')
106 
107  end associate
108  enddo
109 
111 
112 
113 !--------------------------------------------------------------------------------------------------
115 !--------------------------------------------------------------------------------------------------
116 subroutine kinematics_slipplane_opening_lianditstangent(Ld, dLd_dTstar, S, ipc, ip, el)
117 
118  integer, intent(in) :: &
119  ipc, & !< grain number
120  ip, & !< integration point number
121  el
122  real(preal), intent(in), dimension(3,3) :: &
123  s
124  real(preal), intent(out), dimension(3,3) :: &
125  ld
126  real(preal), intent(out), dimension(3,3,3,3) :: &
127  dld_dtstar
128 
129  integer :: &
130  instance, phase, &
131  homog, damageoffset, &
132  i, k, l, m, n
133  real(preal) :: &
134  traction_d, traction_t, traction_n, traction_crit, &
135  udotd, dudotd_dt, udott, dudott_dt, udotn, dudotn_dt
136 
137  phase = material_phaseat(ipc,el)
138  instance = kinematics_slipplane_opening_instance(phase)
139  homog = material_homogenizationat(el)
140  damageoffset = damagemapping(homog)%p(ip,el)
141 
142  associate(prm => param(instance))
143  ld = 0.0_preal
144  dld_dtstar = 0.0_preal
145  do i = 1, prm%sum_N_sl
146 
147  traction_d = math_tensordot(s,prm%P_d(1:3,1:3,i))
148  traction_t = math_tensordot(s,prm%P_t(1:3,1:3,i))
149  traction_n = math_tensordot(s,prm%P_n(1:3,1:3,i))
150 
151  traction_crit = prm%critLoad(i)* damage(homog)%p(damageoffset) ! degrading critical load carrying capacity by damage
152 
153  udotd = sign(1.0_preal,traction_d)* prm%sdot0* ( abs(traction_d)/traction_crit &
154  - abs(traction_d)/prm%critLoad(i))**prm%n
155  udott = sign(1.0_preal,traction_t)* prm%sdot0* ( abs(traction_t)/traction_crit &
156  - abs(traction_t)/prm%critLoad(i))**prm%n
157  udotn = prm%sdot0* ( max(0.0_preal,traction_n)/traction_crit &
158  - max(0.0_preal,traction_n)/prm%critLoad(i))**prm%n
159 
160  if (dneq0(traction_d)) then
161  dudotd_dt = udotd*prm%n/traction_d
162  else
163  dudotd_dt = 0.0_preal
164  endif
165  if (dneq0(traction_t)) then
166  dudott_dt = udott*prm%n/traction_t
167  else
168  dudott_dt = 0.0_preal
169  endif
170  if (dneq0(traction_n)) then
171  dudotn_dt = udotn*prm%n/traction_n
172  else
173  dudotn_dt = 0.0_preal
174  endif
175 
176  forall (k=1:3,l=1:3,m=1:3,n=1:3) &
177  dld_dtstar(k,l,m,n) = dld_dtstar(k,l,m,n) &
178  + dudotd_dt*prm%P_d(k,l,i)*prm%P_d(m,n,i) &
179  + dudott_dt*prm%P_t(k,l,i)*prm%P_t(m,n,i) &
180  + dudotn_dt*prm%P_n(k,l,i)*prm%P_n(m,n,i)
181 
182  ld = ld &
183  + udotd*prm%P_d(1:3,1:3,i) &
184  + udott*prm%P_t(1:3,1:3,i) &
185  + udotn*prm%P_n(1:3,1:3,i)
186  enddo
187 
188  end associate
189 
191 
math::math_expand
pure real(preal) function, dimension(sum(how)) math_expand(what, how)
vector expansion
Definition: math.f90:207
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
math::math_outer
pure real(preal) function, dimension(size(a, 1), size(b, 1)) math_outer(A, B)
outer product of arbitrary sized vectors (A ⊗ B / i,j)
Definition: math.f90:327
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
kinematics_slipplane_opening::kinematics_slipplane_opening_init
subroutine, public kinematics_slipplane_opening_init
module initialization
Definition: kinematics_slipplane_opening.f90:53
kinematics_slipplane_opening::kinematics_slipplane_opening_lianditstangent
subroutine, public kinematics_slipplane_opening_lianditstangent(Ld, dLd_dTstar, S, ipc, ip, el)
contains the constitutive equation for calculating the velocity gradient
Definition: kinematics_slipplane_opening.f90:117
kinematics_slipplane_opening
material subroutine incorporating kinematics resulting from opening of slip planes
Definition: kinematics_slipplane_opening.f90:11
lattice::lattice_slip_direction
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_direction(Nslip, structure, cOverA)
Slip direction of slip systems (|| b)
Definition: lattice.f90:1594
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
lattice::lattice_slip_normal
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_normal(Nslip, structure, cOverA)
Normal direction of slip systems (|| n)
Definition: lattice.f90:1612
material::kinematics_slipplane_opening_id
@, public kinematics_slipplane_opening_id
Definition: material.f90:87
material::phase_kinematics
integer(kind(source_undefined_id)), dimension(:,:), allocatable, public, protected phase_kinematics
active kinematic mechanisms of each phase
Definition: material.f90:105
kinematics_slipplane_opening::kinematics_slipplane_opening_instance
integer, dimension(:), allocatable kinematics_slipplane_opening_instance
Definition: kinematics_slipplane_opening.f90:23
kinematics_slipplane_opening::param
type(tparameters), dimension(:), allocatable param
containers of constitutive parameters (len Ninstance)
Definition: kinematics_slipplane_opening.f90:39
lattice::lattice_slip_transverse
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_transverse(Nslip, structure, cOverA)
Transverse direction of slip systems ( || t = b x n)
Definition: lattice.f90:1630
lattice
contains lattice structure definitions including Schmid matrices for slip, twin, trans,
Definition: lattice.f90:13
material::kinematics_slipplane_opening_label
character(len= *), parameter, public kinematics_slipplane_opening_label
Definition: material.f90:25
math
Mathematical library, including random number generation and tensor representations.
Definition: math.f90:12
debug::debug_constitutive
integer, parameter, public debug_constitutive
stores debug level for constitutive part of DAMASK bitwise coded
Definition: debug.f90:32
kinematics_slipplane_opening::tparameters
container type for internal constitutive parameters
Definition: kinematics_slipplane_opening.f90:25