DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
damage_nonlocal.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/damage_nonlocal.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/damage_nonlocal.f90"
5 !--------------------------------------------------------------------------------------------------
8 !--------------------------------------------------------------------------------------------------
10  use prec
11  use material
12  use config
13  use numerics
14  use crystallite
15  use lattice
20  use results
21 
22  implicit none
23  private
24 
25  type :: tparameters
26  character(len=pStringLen), allocatable, dimension(:) :: &
27  output
28  end type tparameters
29 
30  type(tparameters), dimension(:), allocatable :: &
31  param
32 
33  public :: &
40 
41 contains
42 
43 !--------------------------------------------------------------------------------------------------
46 !--------------------------------------------------------------------------------------------------
47 subroutine damage_nonlocal_init
48 
49  integer :: ninstance,nofmyhomog,h
50 
51  write(6,'(/,a)') ' <<<+- damage_'//damage_nonlocal_label//' init -+>>>'; flush(6)
52 
53  ninstance = count(damage_type == damage_nonlocal_id)
54  allocate(param(ninstance))
55 
56  do h = 1, size(config_homogenization)
57  if (damage_type(h) /= damage_nonlocal_id) cycle
58  associate(prm => param(damage_typeinstance(h)),config => config_homogenization(h))
59 
60  prm%output = config%getStrings('(output)',defaultval=emptystringarray)
61 
62  nofmyhomog = count(material_homogenizationat == h)
63  damagestate(h)%sizeState = 1
64  allocate(damagestate(h)%state0 (1,nofmyhomog), source=damage_initialphi(h))
65  allocate(damagestate(h)%subState0(1,nofmyhomog), source=damage_initialphi(h))
66  allocate(damagestate(h)%state (1,nofmyhomog), source=damage_initialphi(h))
67 
68  nullify(damagemapping(h)%p)
70  deallocate(damage(h)%p)
71  damage(h)%p => damagestate(h)%state(1,:)
72 
73  end associate
74  enddo
75 
76 end subroutine damage_nonlocal_init
77 
78 
79 !--------------------------------------------------------------------------------------------------
81 !--------------------------------------------------------------------------------------------------
82 subroutine damage_nonlocal_getsourceanditstangent(phiDot, dPhiDot_dPhi, phi, ip, el)
83 
84  integer, intent(in) :: &
85  ip, & !< integration point number
86  el
87  real(preal), intent(in) :: &
88  phi
89  integer :: &
90  phase, &
91  grain, &
92  source, &
93  constituent
94  real(preal) :: &
95  phidot, dphidot_dphi, localphidot, dlocalphidot_dphi
96 
97  phidot = 0.0_preal
98  dphidot_dphi = 0.0_preal
99  do grain = 1, homogenization_ngrains(material_homogenizationat(el))
100  phase = material_phaseat(grain,el)
101  constituent = material_phasememberat(grain,ip,el)
102  do source = 1, phase_nsources(phase)
103  select case(phase_source(source,phase))
104  case (source_damage_isobrittle_id)
105  call source_damage_isobrittle_getrateanditstangent (localphidot, dlocalphidot_dphi, phi, phase, constituent)
106 
107  case (source_damage_isoductile_id)
108  call source_damage_isoductile_getrateanditstangent (localphidot, dlocalphidot_dphi, phi, phase, constituent)
109 
110  case (source_damage_anisobrittle_id)
111  call source_damage_anisobrittle_getrateanditstangent(localphidot, dlocalphidot_dphi, phi, phase, constituent)
112 
113  case (source_damage_anisoductile_id)
114  call source_damage_anisoductile_getrateanditstangent(localphidot, dlocalphidot_dphi, phi, phase, constituent)
115 
116  case default
117  localphidot = 0.0_preal
118  dlocalphidot_dphi = 0.0_preal
119 
120  end select
121  phidot = phidot + localphidot
122  dphidot_dphi = dphidot_dphi + dlocalphidot_dphi
123  enddo
124  enddo
125 
126  phidot = phidot/real(homogenization_ngrains(material_homogenizationat(el)),preal)
127  dphidot_dphi = dphidot_dphi/real(homogenization_ngrains(material_homogenizationat(el)),preal)
128 
130 
131 
132 !--------------------------------------------------------------------------------------------------
134 !--------------------------------------------------------------------------------------------------
135 function damage_nonlocal_getdiffusion(ip,el)
136 
137  integer, intent(in) :: &
138  ip, & !< integration point number
139  el
140  real(preal), dimension(3,3) :: &
142  integer :: &
143  homog, &
144  grain
145 
146  homog = material_homogenizationat(el)
147  damage_nonlocal_getdiffusion = 0.0_preal
148  do grain = 1, homogenization_ngrains(homog)
150  crystallite_push33toref(grain,ip,el,lattice_damagediffusion(1:3,1:3,material_phaseat(grain,el)))
151  enddo
152 
154  charlength**2*damage_nonlocal_getdiffusion/real(homogenization_ngrains(homog),preal)
155 
156 end function damage_nonlocal_getdiffusion
157 
158 
159 !--------------------------------------------------------------------------------------------------
161 !--------------------------------------------------------------------------------------------------
162 real(pReal) function damage_nonlocal_getmobility(ip,el)
163 
164  integer, intent(in) :: &
165  ip, & !< integration point number
166  el
167  integer :: &
168  ipc
169 
170  damage_nonlocal_getmobility = 0.0_preal
171 
172  do ipc = 1, homogenization_ngrains(material_homogenizationat(el))
173  damage_nonlocal_getmobility = damage_nonlocal_getmobility + lattice_damagemobility(material_phaseat(ipc,el))
174  enddo
175 
177  real(homogenization_ngrains(material_homogenizationat(el)),preal)
178 
179 end function damage_nonlocal_getmobility
180 
181 
182 !--------------------------------------------------------------------------------------------------
184 !--------------------------------------------------------------------------------------------------
185 subroutine damage_nonlocal_putnonlocaldamage(phi,ip,el)
186 
187  integer, intent(in) :: &
188  ip, & !< integration point number
189  el
190  real(preal), intent(in) :: &
191  phi
192  integer :: &
193  homog, &
194  offset
195 
196  homog = material_homogenizationat(el)
197  offset = damagemapping(homog)%p(ip,el)
198  damage(homog)%p(offset) = phi
199 
201 
202 
203 !--------------------------------------------------------------------------------------------------
205 !--------------------------------------------------------------------------------------------------
206 subroutine damage_nonlocal_results(homog,group)
207 
208  integer, intent(in) :: homog
209  character(len=*), intent(in) :: group
210 
211  integer :: o
212 
213  associate(prm => param(damage_typeinstance(homog)))
214  outputsloop: do o = 1,size(prm%output)
215  select case(prm%output(o))
216  case ('damage')
217  call results_writedataset(group,damage(homog)%p,'phi',&
218  'damage indicator','-')
219  end select
220  enddo outputsloop
221  end associate
222 
223 end subroutine damage_nonlocal_results
224 
225 end module damage_nonlocal
source_damage_isobrittle
material subroutine incoprorating isotropic brittle damage source mechanism
Definition: source_damage_isoBrittle.f90:11
damage_nonlocal
material subroutine for non-locally evolving damage field
Definition: damage_nonlocal.f90:9
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
damage_nonlocal::damage_nonlocal_init
subroutine, public damage_nonlocal_init
module initialization
Definition: damage_nonlocal.f90:48
prec::emptystringarray
character(len=pstringlen), dimension(0), parameter emptystringarray
Definition: prec.f90:88
material::damage_nonlocal_id
@, public damage_nonlocal_id
Definition: material.f90:87
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
damage_nonlocal::damage_nonlocal_putnonlocaldamage
subroutine, public damage_nonlocal_putnonlocaldamage(phi, ip, el)
updated nonlocal damage field with solution from damage phase field PDE
Definition: damage_nonlocal.f90:186
material::damagestate
type(tstate), dimension(:), allocatable, public damagestate
Definition: material.f90:141
crystallite
crystallite state integration functions and reporting of results
Definition: crystallite.f90:15
prec
setting precision for real and int type
Definition: prec.f90:13
damage_nonlocal::damage_nonlocal_getmobility
real(preal) function, public damage_nonlocal_getmobility(ip, el)
Returns homogenized nonlocal damage mobility.
Definition: damage_nonlocal.f90:163
material::damage_typeinstance
integer, dimension(:), allocatable, public, protected damage_typeinstance
instance of particular type of each nonlocal damage
Definition: material.f90:113
source_damage_anisoductile
material subroutine incorporating anisotropic ductile damage source mechanism
Definition: source_damage_anisoDuctile.f90:11
damage_nonlocal::param
type(tparameters), dimension(:), allocatable param
Definition: damage_nonlocal.f90:30
source_damage_anisobrittle
material subroutine incorporating anisotropic brittle damage source mechanism
Definition: source_damage_anisoBrittle.f90:11
damage_nonlocal::tparameters
Definition: damage_nonlocal.f90:25
lattice
contains lattice structure definitions including Schmid matrices for slip, twin, trans,
Definition: lattice.f90:13
material::damage_nonlocal_label
character(len= *), parameter, public damage_nonlocal_label
Definition: material.f90:25
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_nonlocal::damage_nonlocal_results
subroutine, public damage_nonlocal_results(homog, group)
writes results to HDF5 output file
Definition: damage_nonlocal.f90:207
material::material_homogenizationmemberat
integer, dimension(:,:), allocatable, target, public material_homogenizationmemberat
position of the element within its homogenization instance
Definition: material.f90:130
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
numerics
Managing of parameters related to numerics.
Definition: numerics.f90:10
damage_nonlocal::damage_nonlocal_getdiffusion
real(preal) function, dimension(3, 3), public damage_nonlocal_getdiffusion(ip, el)
returns homogenized non local damage diffusion tensor in reference configuration
Definition: damage_nonlocal.f90:136
config::config_homogenization
type(tpartitionedstringlist), dimension(:), allocatable, public, protected config_homogenization
Definition: config.f90:23
damage_nonlocal::damage_nonlocal_getsourceanditstangent
subroutine, public damage_nonlocal_getsourceanditstangent(phiDot, dPhiDot_dPhi, phi, ip, el)
calculates homogenized damage driving forces
Definition: damage_nonlocal.f90:83