|
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/kinematics_slipplane_opening.f90"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/kinematics_slipplane_opening.f90"
31 real(
preal),
dimension(:),
allocatable :: &
33 real(
preal),
dimension(:,:,:),
allocatable :: &
54 integer :: ninstance,p,i
55 character(len=pStringLen) :: extmsg =
''
56 integer,
dimension(:),
allocatable :: n_sl
57 real(
preal),
dimension(:,:),
allocatable :: d,n,t
63 write(6,
'(a16,1x,i5,/)')
'# instances:',ninstance
66 allocate(
param(ninstance))
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))
80 config%getFloat(
'c/a',defaultval=0.0_preal))
82 config%getFloat(
'c/a',defaultval=0.0_preal))
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)))
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))
93 prm%critLoad =
config%getFloats(
'anisoductile_criticalload',requiredsize=
size(n_sl))
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'
118 integer,
intent(in) :: &
119 ipc, & !< grain number
120 ip, & !< integration point number
122 real(preal),
intent(in),
dimension(3,3) :: &
124 real(preal),
intent(out),
dimension(3,3) :: &
126 real(preal),
intent(out),
dimension(3,3,3,3) :: &
131 homog, damageoffset, &
134 traction_d, traction_t, traction_n, traction_crit, &
135 udotd, dudotd_dt, udott, dudott_dt, udotn, dudotn_dt
137 phase = material_phaseat(ipc,el)
139 homog = material_homogenizationat(el)
140 damageoffset = damagemapping(homog)%p(ip,el)
142 associate(prm =>
param(instance))
144 dld_dtstar = 0.0_preal
145 do i = 1, prm%sum_N_sl
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))
151 traction_crit = prm%critLoad(i)* damage(homog)%p(damageoffset)
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
160 if (dneq0(traction_d))
then
161 dudotd_dt = udotd*prm%n/traction_d
163 dudotd_dt = 0.0_preal
165 if (dneq0(traction_t))
then
166 dudott_dt = udott*prm%n/traction_t
168 dudott_dt = 0.0_preal
170 if (dneq0(traction_n))
then
171 dudotn_dt = udotn*prm%n/traction_n
173 dudotn_dt = 0.0_preal
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)
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)
pure real(preal) function, dimension(sum(how)) math_expand(what, how)
vector expansion
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
type(tpartitionedstringlist), dimension(:), allocatable, public, protected config_phase
integer, dimension(debug_maxntype+2), public, protected debug_level
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)
Parses material config file, either solverJobName.materialConfig or material.config.
Reads in the material configuration from file.
setting precision for real and int type
subroutine, public kinematics_slipplane_opening_init
module initialization
subroutine, public kinematics_slipplane_opening_lianditstangent(Ld, dLd_dTstar, S, ipc, ip, el)
contains the constitutive equation for calculating the velocity gradient
material subroutine incorporating kinematics resulting from opening of slip planes
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_direction(Nslip, structure, cOverA)
Slip direction of slip systems (|| b)
input/output functions, partly depending on chosen solver
integer, parameter, public debug_levelbasic
integer, parameter preal
number with 15 significant digits, up to 1e+-307 (typically 64 bit)
Reading in and interpretating the debugging settings for the various modules.
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_normal(Nslip, structure, cOverA)
Normal direction of slip systems (|| n)
@, public kinematics_slipplane_opening_id
integer(kind(source_undefined_id)), dimension(:,:), allocatable, public, protected phase_kinematics
active kinematic mechanisms of each phase
integer, dimension(:), allocatable kinematics_slipplane_opening_instance
type(tparameters), dimension(:), allocatable param
containers of constitutive parameters (len Ninstance)
real(preal) function, dimension(3, sum(nslip)), public lattice_slip_transverse(Nslip, structure, cOverA)
Transverse direction of slip systems ( || t = b x n)
contains lattice structure definitions including Schmid matrices for slip, twin, trans,
character(len= *), parameter, public kinematics_slipplane_opening_label
Mathematical library, including random number generation and tensor representations.
integer, parameter, public debug_constitutive
stores debug level for constitutive part of DAMASK bitwise coded
container type for internal constitutive parameters