DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
source_damage_anisoBrittle.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_anisoBrittle.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_anisoBrittle.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 lattice
20  use results
21 
22  implicit none
23  private
24 
25  integer, dimension(:), allocatable :: &
26  source_damage_anisobrittle_offset, & !< which source is my current source mechanism?
28 
29  type :: tparameters
30  real(preal) :: &
31  sdot_0, &
32  n
33  real(preal), dimension(:), allocatable :: &
34  critdisp, &
35  critload
36  real(preal), dimension(:,:,:,:), allocatable :: &
37  cleavage_systems
38  integer :: &
39  sum_n_cl
40  character(len=pStringLen), allocatable, dimension(:) :: &
41  output
42  end type tparameters
43 
44  type(tparameters), dimension(:), allocatable :: param
45 
46 
47  public :: &
52 
53 contains
54 
55 
56 !--------------------------------------------------------------------------------------------------
59 !--------------------------------------------------------------------------------------------------
61 
62  integer :: ninstance,sourceoffset,nipcmyphase,p
63  integer, dimension(:), allocatable :: n_cl
64  character(len=pStringLen) :: extmsg = ''
65 
66  write(6,'(/,a)') ' <<<+- source_'//source_damage_anisobrittle_label//' init -+>>>'; flush(6)
67 
68  ninstance = count(phase_source == source_damage_anisobrittle_id)
70  write(6,'(a16,1x,i5,/)') '# instances:',ninstance
71 
72  allocate(source_damage_anisobrittle_offset(size(config_phase)), source=0)
73  allocate(source_damage_anisobrittle_instance(size(config_phase)), source=0)
74  allocate(param(ninstance))
75 
76  do p = 1, size(config_phase)
78  do sourceoffset = 1, phase_nsources(p)
79  if (phase_source(sourceoffset,p) == source_damage_anisobrittle_id) then
80  source_damage_anisobrittle_offset(p) = sourceoffset
81  exit
82  endif
83  enddo
84 
85  if (all(phase_source(:,p) /= source_damage_anisobrittle_id)) cycle
86  associate(prm => param(source_damage_anisobrittle_instance(p)), &
87  config => config_phase(p))
88 
89  prm%output = config%getStrings('(output)',defaultval=emptystringarray)
90 
91  n_cl = config%getInts('ncleavage',defaultval=emptyintarray)
92  prm%sum_N_cl = sum(abs(n_cl))
93 
94  prm%n = config%getFloat('anisobrittle_ratesensitivity')
95  prm%sdot_0 = config%getFloat('anisobrittle_sdot0')
96 
97  prm%critDisp = config%getFloats('anisobrittle_criticaldisplacement',requiredsize=size(n_cl))
98  prm%critLoad = config%getFloats('anisobrittle_criticalload', requiredsize=size(n_cl))
99 
100  prm%cleavage_systems = lattice_schmidmatrix_cleavage(n_cl,config%getString('lattice_structure'),&
101  config%getFloat('c/a',defaultval=0.0_preal))
102 
103  ! expand: family => system
104  prm%critDisp = math_expand(prm%critDisp,n_cl)
105  prm%critLoad = math_expand(prm%critLoad,n_cl)
106 
107  ! sanity checks
108  if (prm%n <= 0.0_preal) extmsg = trim(extmsg)//' anisobrittle_n'
109  if (prm%sdot_0 <= 0.0_preal) extmsg = trim(extmsg)//' anisobrittle_sdot0'
110  if (any(prm%critLoad < 0.0_preal)) extmsg = trim(extmsg)//' anisobrittle_critLoad'
111  if (any(prm%critDisp < 0.0_preal)) extmsg = trim(extmsg)//' anisobrittle_critDisp'
112 
113  nipcmyphase = count(material_phaseat==p) * discretization_nip
114  call material_allocatesourcestate(p,sourceoffset,nipcmyphase,1,1,0)
115  sourcestate(p)%p(sourceoffset)%atol = config%getFloat('anisobrittle_atol',defaultval=1.0e-3_preal)
116  if(any(sourcestate(p)%p(sourceoffset)%atol < 0.0_preal)) extmsg = trim(extmsg)//' anisobrittle_atol'
117 
118  end associate
119 
120 !--------------------------------------------------------------------------------------------------
121 ! exit if any parameter is out of range
122  if (extmsg /= '') call io_error(211,ext_msg=trim(extmsg)//'('//source_damage_anisobrittle_label//')')
123 
124 enddo
125 
126 end subroutine source_damage_anisobrittle_init
127 
128 
129 !--------------------------------------------------------------------------------------------------
131 !--------------------------------------------------------------------------------------------------
132 subroutine source_damage_anisobrittle_dotstate(S, ipc, ip, el)
133 
134  integer, intent(in) :: &
135  ipc, & !< component-ID of integration point
136  ip, & !< integration point
137  el
138  real(preal), intent(in), dimension(3,3) :: &
139  s
140 
141  integer :: &
142  phase, &
143  constituent, &
144  sourceoffset, &
145  damageoffset, &
146  homog, &
147  i
148  real(preal) :: &
149  traction_d, traction_t, traction_n, traction_crit
150 
151  phase = material_phaseat(ipc,el)
152  constituent = material_phasememberat(ipc,ip,el)
153  sourceoffset = source_damage_anisobrittle_offset(phase)
154  homog = material_homogenizationat(el)
155  damageoffset = damagemapping(homog)%p(ip,el)
156 
157  associate(prm => param(source_damage_anisobrittle_instance(phase)))
158  sourcestate(phase)%p(sourceoffset)%dotState(1,constituent) = 0.0_preal
159  do i = 1, prm%sum_N_cl
160  traction_d = math_tensordot(s,prm%cleavage_systems(1:3,1:3,1,i))
161  traction_t = math_tensordot(s,prm%cleavage_systems(1:3,1:3,2,i))
162  traction_n = math_tensordot(s,prm%cleavage_systems(1:3,1:3,3,i))
163 
164  traction_crit = prm%critLoad(i)*damage(homog)%p(damageoffset)**2.0_preal
165 
166  sourcestate(phase)%p(sourceoffset)%dotState(1,constituent) &
167  = sourcestate(phase)%p(sourceoffset)%dotState(1,constituent) &
168  + prm%sdot_0 / prm%critDisp(i) &
169  * ((max(0.0_preal, abs(traction_d) - traction_crit)/traction_crit)**prm%n + &
170  (max(0.0_preal, abs(traction_t) - traction_crit)/traction_crit)**prm%n + &
171  (max(0.0_preal, abs(traction_n) - traction_crit)/traction_crit)**prm%n)
172  enddo
173  end associate
174 
176 
177 
178 !--------------------------------------------------------------------------------------------------
180 !--------------------------------------------------------------------------------------------------
181 subroutine source_damage_anisobrittle_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
182 
183  integer, intent(in) :: &
184  phase, &
185  constituent
186  real(preal), intent(in) :: &
187  phi
188  real(preal), intent(out) :: &
189  localphidot, &
190  dlocalphidot_dphi
191 
192  integer :: &
193  sourceoffset
194 
195  sourceoffset = source_damage_anisobrittle_offset(phase)
196 
197  dlocalphidot_dphi = -sourcestate(phase)%p(sourceoffset)%state(1,constituent)
198 
199  localphidot = 1.0_preal &
200  + dlocalphidot_dphi*phi
201 
203 
204 
205 !--------------------------------------------------------------------------------------------------
207 !--------------------------------------------------------------------------------------------------
208 subroutine source_damage_anisobrittle_results(phase,group)
209 
210  integer, intent(in) :: phase
211  character(len=*), intent(in) :: group
212 
213  integer :: o
214 
215  associate(prm => param(source_damage_anisobrittle_instance(phase)), &
216  stt => sourcestate(phase)%p(source_damage_anisobrittle_offset(phase))%state)
217  outputsloop: do o = 1,size(prm%output)
218  select case(trim(prm%output(o)))
219  case ('anisobrittle_drivingforce')
220  call results_writedataset(group,stt,'tbd','driving force','tbd')
221  end select
222  enddo outputsloop
223  end associate
224 
226 
material::sourcestate
type(tsourcestate), dimension(:), allocatable, public sourcestate
Definition: material.f90:139
source_damage_anisobrittle::source_damage_anisobrittle_getrateanditstangent
subroutine, public source_damage_anisobrittle_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
returns local part of nonlocal damage driving force
Definition: source_damage_anisoBrittle.f90:182
material::material_phaseat
integer, dimension(:,:), allocatable, public, protected material_phaseat
phase ID of each element
Definition: material.f90:132
source_damage_anisobrittle::source_damage_anisobrittle_offset
integer, dimension(:), allocatable source_damage_anisobrittle_offset
which source is my current source mechanism?
Definition: source_damage_anisoBrittle.f90:25
math::math_expand
pure real(preal) function, dimension(sum(how)) math_expand(what, how)
vector expansion
Definition: math.f90:207
source_damage_anisobrittle::source_damage_anisobrittle_instance
integer, dimension(:), allocatable source_damage_anisobrittle_instance
instance of source mechanism
Definition: source_damage_anisoBrittle.f90:25
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
source_damage_anisobrittle::source_damage_anisobrittle_init
subroutine, public source_damage_anisobrittle_init
module initialization
Definition: source_damage_anisoBrittle.f90:61
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
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
material::source_damage_anisobrittle_id
@, public source_damage_anisobrittle_id
Definition: material.f90:87
source_damage_anisobrittle
material subroutine incorporating anisotropic brittle damage source mechanism
Definition: source_damage_anisoBrittle.f90:11
material::source_damage_anisobrittle_label
character(len= *), parameter, public source_damage_anisobrittle_label
Definition: material.f90:25
source_damage_anisobrittle::param
type(tparameters), dimension(:), allocatable param
containers of constitutive parameters (len Ninstance)
Definition: source_damage_anisoBrittle.f90:44
lattice
contains lattice structure definitions including Schmid matrices for slip, twin, trans,
Definition: lattice.f90:13
prec::emptyintarray
integer, dimension(0), parameter emptyintarray
Definition: prec.f90:84
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
results
Definition: results.f90:11
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::phase_source
integer(kind(source_undefined_id)), dimension(:,:), allocatable, public, protected phase_source
active sources mechanisms of each phase
Definition: material.f90:105
source_damage_anisobrittle::source_damage_anisobrittle_dotstate
subroutine, public source_damage_anisobrittle_dotstate(S, ipc, ip, el)
calculates derived quantities from state
Definition: source_damage_anisoBrittle.f90:133
material::phase_nsources
integer, dimension(:), allocatable, public, protected phase_nsources
number of source mechanisms active in each phase
Definition: material.f90:113
source_damage_anisobrittle::tparameters
container type for internal constitutive parameters
Definition: source_damage_anisoBrittle.f90:29
source_damage_anisobrittle::source_damage_anisobrittle_results
subroutine, public source_damage_anisobrittle_results(phase, group)
writes results to HDF5 output file
Definition: source_damage_anisoBrittle.f90:209