DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
damage_local.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/damage_local.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/damage_local.f90"
5 !--------------------------------------------------------------------------------------------------
8 !--------------------------------------------------------------------------------------------------
10  use prec
11  use material
12  use config
13  use numerics
18  use results
19 
20  implicit none
21  private
22 
23  type :: tparameters
24  character(len=pStringLen), allocatable, dimension(:) :: &
25  output
26  end type tparameters
27 
28  type(tparameters), dimension(:), allocatable :: &
29  param
30 
31  public :: &
35 
36 contains
37 
38 !--------------------------------------------------------------------------------------------------
41 !--------------------------------------------------------------------------------------------------
42 subroutine damage_local_init
43 
44  integer :: ninstance,nofmyhomog,h
45 
46  write(6,'(/,a)') ' <<<+- damage_'//damage_local_label//' init -+>>>'; flush(6)
47 
48  ninstance = count(damage_type == damage_local_id)
49  allocate(param(ninstance))
50 
51  do h = 1, size(config_homogenization)
52  if (damage_type(h) /= damage_local_id) cycle
53  associate(prm => param(damage_typeinstance(h)),config => config_homogenization(h))
54 
55  prm%output = config%getStrings('(output)',defaultval=emptystringarray)
56 
57  nofmyhomog = count(material_homogenizationat == h)
58  damagestate(h)%sizeState = 1
59  allocate(damagestate(h)%state0 (1,nofmyhomog), source=damage_initialphi(h))
60  allocate(damagestate(h)%subState0(1,nofmyhomog), source=damage_initialphi(h))
61  allocate(damagestate(h)%state (1,nofmyhomog), source=damage_initialphi(h))
62 
63  nullify(damagemapping(h)%p)
65  deallocate(damage(h)%p)
66  damage(h)%p => damagestate(h)%state(1,:)
67 
68  end associate
69  enddo
70 
71 end subroutine damage_local_init
72 
73 
74 !--------------------------------------------------------------------------------------------------
76 !--------------------------------------------------------------------------------------------------
77 function damage_local_updatestate(subdt, ip, el)
78 
79  integer, intent(in) :: &
80  ip, & !< integration point number
81  el
82  real(preal), intent(in) :: &
83  subdt
84  logical, dimension(2) :: &
86  integer :: &
87  homog, &
88  offset
89  real(preal) :: &
90  phi, phidot, dphidot_dphi
91 
92  homog = material_homogenizationat(el)
93  offset = material_homogenizationmemberat(ip,el)
94  phi = damagestate(homog)%subState0(1,offset)
95  call damage_local_getsourceanditstangent(phidot, dphidot_dphi, phi, ip, el)
96  phi = max(residualstiffness,min(1.0_preal,phi + subdt*phidot))
97 
98  damage_local_updatestate = [ abs(phi - damagestate(homog)%state(1,offset)) &
99  <= err_damage_tolabs &
100  .or. abs(phi - damagestate(homog)%state(1,offset)) &
101  <= err_damage_tolrel*abs(damagestate(homog)%state(1,offset)), &
102  .true.]
103 
104  damagestate(homog)%state(1,offset) = phi
105 
106 end function damage_local_updatestate
107 
108 
109 !--------------------------------------------------------------------------------------------------
111 !--------------------------------------------------------------------------------------------------
112 subroutine damage_local_getsourceanditstangent(phiDot, dPhiDot_dPhi, phi, ip, el)
113 
114  integer, intent(in) :: &
115  ip, & !< integration point number
116  el
117  real(pReal), intent(in) :: &
118  phi
119  integer :: &
120  phase, &
121  grain, &
122  source, &
123  constituent
124  real(pReal) :: &
125  phiDot, dPhiDot_dPhi, localphiDot, dLocalphiDot_dPhi
126 
127  phidot = 0.0_preal
128  dphidot_dphi = 0.0_preal
129  do grain = 1, homogenization_ngrains(material_homogenizationat(el))
130  phase = material_phaseat(grain,el)
131  constituent = material_phasememberat(grain,ip,el)
132  do source = 1, phase_nsources(phase)
133  select case(phase_source(source,phase))
134  case (source_damage_isobrittle_id)
135  call source_damage_isobrittle_getrateanditstangent (localphidot, dlocalphidot_dphi, phi, phase, constituent)
136 
137  case (source_damage_isoductile_id)
138  call source_damage_isoductile_getrateanditstangent (localphidot, dlocalphidot_dphi, phi, phase, constituent)
139 
140  case (source_damage_anisobrittle_id)
141  call source_damage_anisobrittle_getrateanditstangent(localphidot, dlocalphidot_dphi, phi, phase, constituent)
142 
143  case (source_damage_anisoductile_id)
144  call source_damage_anisoductile_getrateanditstangent(localphidot, dlocalphidot_dphi, phi, phase, constituent)
145 
146  case default
147  localphidot = 0.0_preal
148  dlocalphidot_dphi = 0.0_preal
149 
150  end select
151  phidot = phidot + localphidot
152  dphidot_dphi = dphidot_dphi + dlocalphidot_dphi
153  enddo
154  enddo
155 
156  phidot = phidot/real(homogenization_ngrains(material_homogenizationat(el)),preal)
157  dphidot_dphi = dphidot_dphi/real(homogenization_ngrains(material_homogenizationat(el)),preal)
158 
160 
161 
162 !--------------------------------------------------------------------------------------------------
164 !--------------------------------------------------------------------------------------------------
165 subroutine damage_local_results(homog,group)
166 
167  integer, intent(in) :: homog
168  character(len=*), intent(in) :: group
169 
170  integer :: o
171 
172  associate(prm => param(damage_typeinstance(homog)))
173  outputsloop: do o = 1,size(prm%output)
174  select case(prm%output(o))
175  case ('damage')
176  call results_writedataset(group,damage(homog)%p,'phi',&
177  'damage indicator','-')
178  end select
179  enddo outputsloop
180  end associate
181 
182 end subroutine damage_local_results
183 
184 
185 end module damage_local
material::damage_local_label
character(len= *), parameter, public damage_local_label
Definition: material.f90:25
source_damage_isobrittle
material subroutine incoprorating isotropic brittle damage source mechanism
Definition: source_damage_isoBrittle.f90:11
source_damage_isoductile
material subroutine incoprorating isotropic ductile damage source mechanism
Definition: source_damage_isoDuctile.f90:11
material::damagemapping
type(thomogmapping), dimension(:), allocatable, public damagemapping
mapping for damage state/fields
Definition: material.f90:170
prec::emptystringarray
character(len=pstringlen), dimension(0), parameter emptystringarray
Definition: prec.f90:88
damage_local
material subroutine for locally evolving damage field
Definition: damage_local.f90:9
material::damage_initialphi
real(preal), dimension(:), allocatable, public, protected damage_initialphi
initial damage per each homogenization
Definition: material.f90:124
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
material::damagestate
type(tstate), dimension(:), allocatable, public damagestate
Definition: material.f90:141
damage_local::damage_local_init
subroutine, public damage_local_init
module initialization
Definition: damage_local.f90:43
material::damage_local_id
@, public damage_local_id
Definition: material.f90:87
prec
setting precision for real and int type
Definition: prec.f90:13
material::damage_typeinstance
integer, dimension(:), allocatable, public, protected damage_typeinstance
instance of particular type of each nonlocal damage
Definition: material.f90:113
damage_local::tparameters
Definition: damage_local.f90:23
source_damage_anisoductile
material subroutine incorporating anisotropic ductile damage source mechanism
Definition: source_damage_anisoDuctile.f90:11
damage_local::damage_local_results
subroutine, public damage_local_results(homog, group)
writes results to HDF5 output file
Definition: damage_local.f90:166
source_damage_anisobrittle
material subroutine incorporating anisotropic brittle damage source mechanism
Definition: source_damage_anisoBrittle.f90:11
results
Definition: results.f90:11
material::material_homogenizationat
integer, dimension(:), allocatable, public, protected material_homogenizationat
homogenization ID of each element (copy of discretization_homogenizationAt)
Definition: material.f90:128
damage_local::damage_local_getsourceanditstangent
subroutine damage_local_getsourceanditstangent(phiDot, dPhiDot_dPhi, phi, ip, el)
calculates homogenized local damage driving forces
Definition: damage_local.f90:113
material::material_homogenizationmemberat
integer, dimension(:,:), allocatable, target, public material_homogenizationmemberat
position of the element within its homogenization instance
Definition: material.f90:130
damage_local::param
type(tparameters), dimension(:), allocatable param
Definition: damage_local.f90:28
material::damage
type(group_float), dimension(:), allocatable, public damage
damage field
Definition: material.f90:174
material::damage_type
integer(kind(damage_none_id)), dimension(:), allocatable, public, protected damage_type
nonlocal damage model
Definition: material.f90:96
damage_local::damage_local_updatestate
logical function, dimension(2), public damage_local_updatestate(subdt, ip, el)
calculates local change in damage field
Definition: damage_local.f90:78
numerics
Managing of parameters related to numerics.
Definition: numerics.f90:10
config::config_homogenization
type(tpartitionedstringlist), dimension(:), allocatable, public, protected config_homogenization
Definition: config.f90:23