DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
source_damage_isoBrittle.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_isoBrittle.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/source_damage_isoBrittle.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 results
20 
21  implicit none
22  private
23 
24  integer, dimension(:), allocatable :: &
27 
28  type, private :: tparameters
29  real(preal) :: &
30  critstrainenergy, &
31  n
32  character(len=pStringLen), allocatable, dimension(:) :: &
33  output
34  end type tparameters
35 
36  type(tparameters), dimension(:), allocatable :: param
37 
38 
39  public :: &
44 
45 contains
46 
47 
48 !--------------------------------------------------------------------------------------------------
51 !--------------------------------------------------------------------------------------------------
53 
54  integer :: ninstance,sourceoffset,nipcmyphase,p
55  character(len=pStringLen) :: extmsg = ''
56 
57  write(6,'(/,a)') ' <<<+- source_'//source_damage_isobrittle_label//' init -+>>>'; flush(6)
58 
59  ninstance = count(phase_source == source_damage_isobrittle_id)
61  write(6,'(a16,1x,i5,/)') '# instances:',ninstance
62 
63  allocate(source_damage_isobrittle_offset(size(config_phase)), source=0)
64  allocate(source_damage_isobrittle_instance(size(config_phase)), source=0)
65  allocate(param(ninstance))
66 
67  do p = 1, size(config_phase)
69  do sourceoffset = 1, phase_nsources(p)
70  if (phase_source(sourceoffset,p) == source_damage_isobrittle_id) then
71  source_damage_isobrittle_offset(p) = sourceoffset
72  exit
73  endif
74  enddo
75 
76  if (all(phase_source(:,p) /= source_damage_isobrittle_id)) cycle
77  associate(prm => param(source_damage_isobrittle_instance(p)), &
78  config => config_phase(p))
79 
80  prm%output = config%getStrings('(output)',defaultval=emptystringarray)
81 
82  prm%N = config%getFloat('isobrittle_n')
83  prm%critStrainEnergy = config%getFloat('isobrittle_criticalstrainenergy')
84 
85  ! sanity checks
86  if (prm%N <= 0.0_preal) extmsg = trim(extmsg)//' isobrittle_n'
87  if (prm%critStrainEnergy <= 0.0_preal) extmsg = trim(extmsg)//' isobrittle_criticalstrainenergy'
88 
89  nipcmyphase = count(material_phaseat==p) * discretization_nip
90  call material_allocatesourcestate(p,sourceoffset,nipcmyphase,1,1,1)
91  sourcestate(p)%p(sourceoffset)%atol = config%getFloat('isobrittle_atol',defaultval=1.0e-3_preal)
92  if(any(sourcestate(p)%p(sourceoffset)%atol < 0.0_preal)) extmsg = trim(extmsg)//' isobrittle_atol'
93 
94  end associate
95 
96 !--------------------------------------------------------------------------------------------------
97 ! exit if any parameter is out of range
98  if (extmsg /= '') call io_error(211,ext_msg=trim(extmsg)//'('//source_damage_isobrittle_label//')')
99 
100 enddo
101 
102 end subroutine source_damage_isobrittle_init
103 
104 
105 !--------------------------------------------------------------------------------------------------
107 !--------------------------------------------------------------------------------------------------
108 subroutine source_damage_isobrittle_deltastate(C, Fe, ipc, ip, el)
109 
110  integer, intent(in) :: &
111  ipc, & !< component-ID of integration point
112  ip, & !< integration point
113  el
114  real(preal), intent(in), dimension(3,3) :: &
115  fe
116  real(preal), intent(in), dimension(6,6) :: &
117  c
118 
119  integer :: &
120  phase, &
121  constituent, &
122  sourceoffset
123  real(preal), dimension(6) :: &
124  strain
125  real(preal) :: &
126  strainenergy
127 
128  phase = material_phaseat(ipc,el)
129  constituent = material_phasememberat(ipc,ip,el)
130  sourceoffset = source_damage_isobrittle_offset(phase)
131 
132  strain = 0.5_preal*math_sym33to6(matmul(transpose(fe),fe)-math_i3)
133 
134  associate(prm => param(source_damage_isobrittle_instance(phase)))
135  strainenergy = 2.0_preal*sum(strain*matmul(c,strain))/prm%critStrainEnergy
136  ! ToDo: check strainenergy = 2.0_pReal*dot_product(strain,matmul(C,strain))/param(instance)%critStrainEnergy
137 
138  if (strainenergy > sourcestate(phase)%p(sourceoffset)%subState0(1,constituent)) then
139  sourcestate(phase)%p(sourceoffset)%deltaState(1,constituent) = &
140  strainenergy - sourcestate(phase)%p(sourceoffset)%state(1,constituent)
141  else
142  sourcestate(phase)%p(sourceoffset)%deltaState(1,constituent) = &
143  sourcestate(phase)%p(sourceoffset)%subState0(1,constituent) - &
144  sourcestate(phase)%p(sourceoffset)%state(1,constituent)
145  endif
146  end associate
147 
149 
150 
151 !--------------------------------------------------------------------------------------------------
153 !--------------------------------------------------------------------------------------------------
154 subroutine source_damage_isobrittle_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
155 
156  integer, intent(in) :: &
157  phase, &
158  constituent
159  real(preal), intent(in) :: &
160  phi
161  real(preal), intent(out) :: &
162  localphidot, &
163  dlocalphidot_dphi
164 
165  integer :: &
166  sourceoffset
167 
168  sourceoffset = source_damage_isobrittle_offset(phase)
169 
170  associate(prm => param(source_damage_isobrittle_instance(phase)))
171  localphidot = (1.0_preal - phi)**(prm%n - 1.0_preal) &
172  - phi*sourcestate(phase)%p(sourceoffset)%state(1,constituent)
173  dlocalphidot_dphi = - (prm%n - 1.0_preal)* (1.0_preal - phi)**max(0.0_preal,prm%n - 2.0_preal) &
174  - sourcestate(phase)%p(sourceoffset)%state(1,constituent)
175  end associate
176 
178 
179 
180 !--------------------------------------------------------------------------------------------------
182 !--------------------------------------------------------------------------------------------------
183 subroutine source_damage_isobrittle_results(phase,group)
184 
185  integer, intent(in) :: phase
186  character(len=*), intent(in) :: group
187 
188  integer :: o
189 
190  associate(prm => param(source_damage_isobrittle_instance(phase)), &
191  stt => sourcestate(phase)%p(source_damage_isobrittle_offset(phase))%state)
192  outputsloop: do o = 1,size(prm%output)
193  select case(trim(prm%output(o)))
194  case ('isobrittle_drivingforce')
195  call results_writedataset(group,stt,'tbd','driving force','tbd')
196  end select
197  enddo outputsloop
198  end associate
199 
201 
202 end module source_damage_isobrittle
material::sourcestate
type(tsourcestate), dimension(:), allocatable, public sourcestate
Definition: material.f90:139
material::material_phaseat
integer, dimension(:,:), allocatable, public, protected material_phaseat
phase ID of each element
Definition: material.f90:132
source_damage_isobrittle
material subroutine incoprorating isotropic brittle damage source mechanism
Definition: source_damage_isoBrittle.f90:11
source_damage_isobrittle::source_damage_isobrittle_offset
integer, dimension(:), allocatable source_damage_isobrittle_offset
Definition: source_damage_isoBrittle.f90:24
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
material::source_damage_isobrittle_id
@, public source_damage_isobrittle_id
Definition: material.f90:87
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
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
source_damage_isobrittle::source_damage_isobrittle_instance
integer, dimension(:), allocatable source_damage_isobrittle_instance
Definition: source_damage_isoBrittle.f90:24
prec
setting precision for real and int type
Definition: prec.f90:13
material::source_damage_isobrittle_label
character(len= *), parameter, public source_damage_isobrittle_label
Definition: material.f90:25
discretization
spatial discretization
Definition: discretization.f90:9
discretization::discretization_nip
integer, public, protected discretization_nip
Definition: discretization.f90:17
source_damage_isobrittle::source_damage_isobrittle_deltastate
subroutine, public source_damage_isobrittle_deltastate(C, Fe, ipc, ip, el)
calculates derived quantities from state
Definition: source_damage_isoBrittle.f90:109
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_isobrittle::source_damage_isobrittle_results
subroutine, public source_damage_isobrittle_results(phase, group)
writes results to HDF5 output file
Definition: source_damage_isoBrittle.f90:184
source_damage_isobrittle::source_damage_isobrittle_getrateanditstangent
subroutine, public source_damage_isobrittle_getrateanditstangent(localphiDot, dLocalphiDot_dPhi, phi, phase, constituent)
returns local part of nonlocal damage driving force
Definition: source_damage_isoBrittle.f90:155
source_damage_isobrittle::source_damage_isobrittle_init
subroutine, public source_damage_isobrittle_init
module initialization
Definition: source_damage_isoBrittle.f90:53
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
material::phase_nsources
integer, dimension(:), allocatable, public, protected phase_nsources
number of source mechanisms active in each phase
Definition: material.f90:113
source_damage_isobrittle::tparameters
container type for internal constitutive parameters
Definition: source_damage_isoBrittle.f90:28
source_damage_isobrittle::param
type(tparameters), dimension(:), allocatable param
containers of constitutive parameters (len Ninstance)
Definition: source_damage_isoBrittle.f90:36