DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
kinematics_cleavage_opening.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/kinematics_cleavage_opening.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/kinematics_cleavage_opening.f90"
5 !--------------------------------------------------------------------------------------------------
10 !--------------------------------------------------------------------------------------------------
12  use prec
13  use io
14  use config
15  use debug
16  use math
17  use lattice
18  use material
19 
20  implicit none
21  private
22 
23  integer, dimension(:), allocatable :: kinematics_cleavage_opening_instance
24 
25  type :: tparameters
26  integer :: &
27  sum_n_cl
28  real(preal) :: &
29  sdot0, &
30  n
31  real(preal), dimension(:), allocatable :: &
32  critload
33  real(preal), dimension(:,:,:,:), allocatable :: &
34  cleavage_systems
35  end type tparameters
36 
37  type(tparameters), dimension(:), allocatable :: param
38 
39  public :: &
42 
43 contains
44 
45 
46 !--------------------------------------------------------------------------------------------------
49 !--------------------------------------------------------------------------------------------------
51 
52  integer :: ninstance,p
53  integer, dimension(:), allocatable :: n_cl
54  character(len=pStringLen) :: extmsg = ''
55 
56  write(6,'(/,a)') ' <<<+- kinematics_'//kinematics_cleavage_opening_label//' init -+>>>'; flush(6)
57 
60  write(6,'(a16,1x,i5,/)') '# instances:',ninstance
61 
62  allocate(kinematics_cleavage_opening_instance(size(config_phase)), source=0)
63  allocate(param(ninstance))
64 
65  do p = 1, size(config_phase)
67  if (all(phase_kinematics(:,p) /= kinematics_cleavage_opening_id)) cycle
68 
69  associate(prm => param(kinematics_cleavage_opening_instance(p)), &
70  config => config_phase(p))
71 
72  n_cl = config%getInts('ncleavage')
73  prm%sum_N_cl = sum(abs(n_cl))
74 
75  prm%n = config%getFloat('anisobrittle_ratesensitivity')
76  prm%sdot0 = config%getFloat('anisobrittle_sdot0')
77 
78  prm%critLoad = config%getFloats('anisobrittle_criticalload',requiredsize=size(n_cl))
79 
80  prm%cleavage_systems = lattice_schmidmatrix_cleavage(n_cl,config%getString('lattice_structure'),&
81  config%getFloat('c/a',defaultval=0.0_preal))
82 
83  ! expand: family => system
84  prm%critLoad = math_expand(prm%critLoad,n_cl)
85 
86  ! sanity checks
87  if (prm%n <= 0.0_preal) extmsg = trim(extmsg)//' anisobrittle_n'
88  if (prm%sdot0 <= 0.0_preal) extmsg = trim(extmsg)//' anisobrittle_sdot0'
89  if (any(prm%critLoad < 0.0_preal)) extmsg = trim(extmsg)//' anisobrittle_critLoad'
90 
91 !--------------------------------------------------------------------------------------------------
92 ! exit if any parameter is out of range
93  if (extmsg /= '') call io_error(211,ext_msg=trim(extmsg)//'('//kinematics_cleavage_opening_label//')')
94 
95  end associate
96  enddo
97 
99 
100 
101 !--------------------------------------------------------------------------------------------------
103 !--------------------------------------------------------------------------------------------------
104 subroutine kinematics_cleavage_opening_lianditstangent(Ld, dLd_dTstar, S, ipc, ip, el)
105 
106  integer, intent(in) :: &
107  ipc, & !< grain number
108  ip, & !< integration point number
109  el
110  real(preal), intent(in), dimension(3,3) :: &
111  s
112  real(preal), intent(out), dimension(3,3) :: &
113  ld
114  real(preal), intent(out), dimension(3,3,3,3) :: &
115  dld_dtstar
116 
117  integer :: &
118  homog, damageoffset, &
119  i, k, l, m, n
120  real(preal) :: &
121  traction_d, traction_t, traction_n, traction_crit, &
122  udotd, dudotd_dt, udott, dudott_dt, udotn, dudotn_dt
123 
124  homog = material_homogenizationat(el)
125  damageoffset = damagemapping(homog)%p(ip,el)
126 
127  ld = 0.0_preal
128  dld_dtstar = 0.0_preal
129  associate(prm => param(kinematics_cleavage_opening_instance(material_phaseat(ipc,el))))
130  do i = 1,prm%sum_N_cl
131  traction_crit = prm%critLoad(i)* damage(homog)%p(damageoffset)**2.0_preal
132 
133  traction_d = math_tensordot(s,prm%cleavage_systems(1:3,1:3,1,i))
134  if (abs(traction_d) > traction_crit + tol_math_check) then
135  udotd = sign(1.0_preal,traction_d)* prm%sdot0 * ((abs(traction_d) - traction_crit)/traction_crit)**prm%n
136  ld = ld + udotd*prm%cleavage_systems(1:3,1:3,1,i)
137  dudotd_dt = sign(1.0_preal,traction_d)*udotd*prm%n / (abs(traction_d) - traction_crit)
138  forall (k=1:3,l=1:3,m=1:3,n=1:3) &
139  dld_dtstar(k,l,m,n) = dld_dtstar(k,l,m,n) &
140  + dudotd_dt*prm%cleavage_systems(k,l,1,i) * prm%cleavage_systems(m,n,1,i)
141  endif
142 
143  traction_t = math_tensordot(s,prm%cleavage_systems(1:3,1:3,2,i))
144  if (abs(traction_t) > traction_crit + tol_math_check) then
145  udott = sign(1.0_preal,traction_t)* prm%sdot0 * ((abs(traction_t) - traction_crit)/traction_crit)**prm%n
146  ld = ld + udott*prm%cleavage_systems(1:3,1:3,2,i)
147  dudott_dt = sign(1.0_preal,traction_t)*udott*prm%n / (abs(traction_t) - traction_crit)
148  forall (k=1:3,l=1:3,m=1:3,n=1:3) &
149  dld_dtstar(k,l,m,n) = dld_dtstar(k,l,m,n) &
150  + dudott_dt*prm%cleavage_systems(k,l,2,i) * prm%cleavage_systems(m,n,2,i)
151  endif
152 
153  traction_n = math_tensordot(s,prm%cleavage_systems(1:3,1:3,3,i))
154  if (abs(traction_n) > traction_crit + tol_math_check) then
155  udotn = sign(1.0_preal,traction_n)* prm%sdot0 * ((abs(traction_n) - traction_crit)/traction_crit)**prm%n
156  ld = ld + udotn*prm%cleavage_systems(1:3,1:3,3,i)
157  dudotn_dt = sign(1.0_preal,traction_n)*udotn*prm%n / (abs(traction_n) - traction_crit)
158  forall (k=1:3,l=1:3,m=1:3,n=1:3) &
159  dld_dtstar(k,l,m,n) = dld_dtstar(k,l,m,n) &
160  + dudotn_dt*prm%cleavage_systems(k,l,3,i) * prm%cleavage_systems(m,n,3,i)
161  endif
162  enddo
163  end associate
164 
166 
kinematics_cleavage_opening::kinematics_cleavage_opening_init
subroutine, public kinematics_cleavage_opening_init
module initialization
Definition: kinematics_cleavage_opening.f90:51
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
kinematics_cleavage_opening
material subroutine incorporating kinematics resulting from opening of cleavage planes
Definition: kinematics_cleavage_opening.f90:11
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
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
material::kinematics_cleavage_opening_label
character(len= *), parameter, public kinematics_cleavage_opening_label
Definition: material.f90:25
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
kinematics_cleavage_opening::tparameters
container type for internal constitutive parameters
Definition: kinematics_cleavage_opening.f90:25
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_cleavage_opening::kinematics_cleavage_opening_instance
integer, dimension(:), allocatable kinematics_cleavage_opening_instance
Definition: kinematics_cleavage_opening.f90:23
lattice
contains lattice structure definitions including Schmid matrices for slip, twin, trans,
Definition: lattice.f90:13
lattice::lattice_schmidmatrix_cleavage
real(preal) function, dimension(3, 3, 3, sum(ncleavage)), public lattice_schmidmatrix_cleavage(Ncleavage, structure, cOverA)
Schmid matrix for cleavage details only active cleavage systems are considered.
Definition: lattice.f90:1546
kinematics_cleavage_opening::kinematics_cleavage_opening_lianditstangent
subroutine, public kinematics_cleavage_opening_lianditstangent(Ld, dLd_dTstar, S, ipc, ip, el)
contains the constitutive equation for calculating the velocity gradient
Definition: kinematics_cleavage_opening.f90:105
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
material::kinematics_cleavage_opening_id
@, public kinematics_cleavage_opening_id
Definition: material.f90:87
kinematics_cleavage_opening::param
type(tparameters), dimension(:), allocatable param
containers of constitutive parameters (len Ninstance)
Definition: kinematics_cleavage_opening.f90:37