DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
thermal_adiabatic.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/thermal_adiabatic.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/thermal_adiabatic.f90"
5 !--------------------------------------------------------------------------------------------------
8 !--------------------------------------------------------------------------------------------------
10  use prec
11  use config
12  use numerics
13  use material
14  use results
17  use crystallite
18  use lattice
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 :: &
38 
39 contains
40 
41 
42 !--------------------------------------------------------------------------------------------------
45 !--------------------------------------------------------------------------------------------------
46 subroutine thermal_adiabatic_init
47 
48  integer :: maxninstance,h,nofmyhomog
49 
50  write(6,'(/,a)') ' <<<+- thermal_'//thermal_adiabatic_label//' init -+>>>'; flush(6)
51 
52  maxninstance = count(thermal_type == thermal_adiabatic_id)
53  if (maxninstance == 0) return
54 
55  allocate(param(maxninstance))
56 
57  do h = 1, size(thermal_type)
58  if (thermal_type(h) /= thermal_adiabatic_id) cycle
59  associate(prm => param(thermal_typeinstance(h)),config => config_homogenization(h))
60 
61  prm%output = config%getStrings('(output)',defaultval=emptystringarray)
62 
63  nofmyhomog=count(material_homogenizationat==h)
64  thermalstate(h)%sizeState = 1
65  allocate(thermalstate(h)%state0 (1,nofmyhomog), source=thermal_initialt(h))
66  allocate(thermalstate(h)%subState0(1,nofmyhomog), source=thermal_initialt(h))
67  allocate(thermalstate(h)%state (1,nofmyhomog), source=thermal_initialt(h))
68 
70  deallocate(temperature(h)%p)
71  temperature(h)%p => thermalstate(h)%state(1,:)
72  deallocate(temperaturerate(h)%p)
73  allocate (temperaturerate(h)%p(nofmyhomog), source=0.0_preal)
74 
75  end associate
76  enddo
77 
78 end subroutine thermal_adiabatic_init
79 
80 
81 !--------------------------------------------------------------------------------------------------
83 !--------------------------------------------------------------------------------------------------
84 function thermal_adiabatic_updatestate(subdt, ip, el)
85 
86  integer, intent(in) :: &
87  ip, & !< integration point number
88  el
89  real(preal), intent(in) :: &
90  subdt
91 
92  logical, dimension(2) :: &
94  integer :: &
95  homog, &
96  offset
97  real(preal) :: &
98  t, tdot, dtdot_dt
99 
100  homog = material_homogenizationat(el)
101  offset = material_homogenizationmemberat(ip,el)
102 
103  t = thermalstate(homog)%subState0(1,offset)
104  call thermal_adiabatic_getsourceanditstangent(tdot, dtdot_dt, t, ip, el)
106 
107  thermal_adiabatic_updatestate = [ abs(t - thermalstate(homog)%state(1,offset)) &
108  <= err_thermal_tolabs &
109  .or. abs(t - thermalstate(homog)%state(1,offset)) &
110  <= err_thermal_tolrel*abs(thermalstate(homog)%state(1,offset)), &
111  .true.]
112 
113  temperature(homog)%p(thermalmapping(homog)%p(ip,el)) = t
114  temperaturerate(homog)%p(thermalmapping(homog)%p(ip,el)) = &
115  (thermalstate(homog)%state(1,offset) - thermalstate(homog)%subState0(1,offset))/(subdt+tiny(0.0_preal))
116 
118 
119 
120 !--------------------------------------------------------------------------------------------------
122 !--------------------------------------------------------------------------------------------------
123 subroutine thermal_adiabatic_getsourceanditstangent(Tdot, dTdot_dT, T, ip, el)
124 
125  integer, intent(in) :: &
126  ip, & !< integration point number
127  el
128  real(preal), intent(in) :: &
129  t
130  real(preal), intent(out) :: &
131  tdot, dtdot_dt
132 
133  real(preal) :: &
134  my_tdot, my_dtdot_dt
135  integer :: &
136  phase, &
137  homog, &
138  instance, &
139  grain, &
140  source, &
141  constituent
142 
143  homog = material_homogenizationat(el)
144  instance = thermal_typeinstance(homog)
145 
146  tdot = 0.0_preal
147  dtdot_dt = 0.0_preal
148  do grain = 1, homogenization_ngrains(homog)
149  phase = material_phaseat(grain,el)
150  constituent = material_phasememberat(grain,ip,el)
151  do source = 1, phase_nsources(phase)
152  select case(phase_source(source,phase))
153  case (source_thermal_dissipation_id)
154  call source_thermal_dissipation_getrateanditstangent(my_tdot, my_dtdot_dt, &
155  crystallite_s(1:3,1:3,grain,ip,el), &
156  crystallite_lp(1:3,1:3,grain,ip,el), &
157  phase)
158 
159  case (source_thermal_externalheat_id)
160  call source_thermal_externalheat_getrateanditstangent(my_tdot, my_dtdot_dt, &
161  phase, constituent)
162 
163  case default
164  my_tdot = 0.0_preal
165  my_dtdot_dt = 0.0_preal
166  end select
167  tdot = tdot + my_tdot
168  dtdot_dt = dtdot_dt + my_dtdot_dt
169  enddo
170  enddo
171 
172  tdot = tdot/real(homogenization_ngrains(homog),preal)
173  dtdot_dt = dtdot_dt/real(homogenization_ngrains(homog),preal)
174 
176 
177 
178 !--------------------------------------------------------------------------------------------------
180 !--------------------------------------------------------------------------------------------------
181 function thermal_adiabatic_getspecificheat(ip,el)
182 
183  integer, intent(in) :: &
184  ip, & !< integration point number
185  el
186 
187  real(preal) :: &
189  integer :: &
190  grain
191 
193 
194  do grain = 1, homogenization_ngrains(material_homogenizationat(el))
196  + lattice_specificheat(material_phaseat(grain,el))
197  enddo
198 
200  / real(homogenization_ngrains(material_homogenizationat(el)),preal)
201 
203 
204 
205 !--------------------------------------------------------------------------------------------------
207 !--------------------------------------------------------------------------------------------------
208 function thermal_adiabatic_getmassdensity(ip,el)
209 
210  integer, intent(in) :: &
211  ip, & !< integration point number
212  el
213  real(preal) :: &
215  integer :: &
216  grain
217 
219 
220  do grain = 1, homogenization_ngrains(material_homogenizationat(el))
222  + lattice_massdensity(material_phaseat(grain,el))
223  enddo
224 
226  / real(homogenization_ngrains(material_homogenizationat(el)),preal)
227 
229 
230 
231 !--------------------------------------------------------------------------------------------------
233 !--------------------------------------------------------------------------------------------------
234 subroutine thermal_adiabatic_results(homog,group)
235 
236  integer, intent(in) :: homog
237  character(len=*), intent(in) :: group
238 
239  integer :: o
240 
241  associate(prm => param(damage_typeinstance(homog)))
242  outputsloop: do o = 1,size(prm%output)
243  select case(trim(prm%output(o)))
244  case('temperature') ! ToDo: should be 'T'
245  call results_writedataset(group,temperature(homog)%p,'T',&
246  'temperature','K')
247  end select
248  enddo outputsloop
249  end associate
250 
251 end subroutine thermal_adiabatic_results
252 
253 end module thermal_adiabatic
material::thermalmapping
type(thomogmapping), dimension(:), allocatable, public thermalmapping
mapping for thermal state/fields
Definition: material.f90:170
thermal_adiabatic::thermal_adiabatic_init
subroutine, public thermal_adiabatic_init
module initialization
Definition: thermal_adiabatic.f90:47
thermal_adiabatic::thermal_adiabatic_getsourceanditstangent
subroutine, public thermal_adiabatic_getsourceanditstangent(Tdot, dTdot_dT, T, ip, el)
returns heat generation rate
Definition: thermal_adiabatic.f90:124
material::thermal_adiabatic_id
@, public thermal_adiabatic_id
Definition: material.f90:87
source_thermal_dissipation
material subroutine for thermal source due to plastic dissipation
Definition: source_thermal_dissipation.f90:11
prec::emptystringarray
character(len=pstringlen), dimension(0), parameter emptystringarray
Definition: prec.f90:88
material::temperature
type(group_float), dimension(:), allocatable, public temperature
temperature field
Definition: material.f90:174
thermal_adiabatic
material subroutine for adiabatic temperature evolution
Definition: thermal_adiabatic.f90:9
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_thermal_dissipation::param
type(tparameters), dimension(:), allocatable param
containers of constitutive parameters (len Ninstance)
Definition: source_thermal_dissipation.f90:30
material::thermal_initialt
real(preal), dimension(:), allocatable, public, protected thermal_initialt
initial temperature per each homogenization
Definition: material.f90:124
crystallite
crystallite state integration functions and reporting of results
Definition: crystallite.f90:15
material::thermal_type
integer(kind(thermal_isothermal_id)), dimension(:), allocatable, public, protected thermal_type
thermal transport model
Definition: material.f90:94
prec
setting precision for real and int type
Definition: prec.f90:13
thermal_adiabatic::param
type(tparameters), dimension(:), allocatable param
Definition: thermal_adiabatic.f90:28
thermal_adiabatic::thermal_adiabatic_getspecificheat
real(preal) function, public thermal_adiabatic_getspecificheat(ip, el)
returns homogenized specific heat capacity
Definition: thermal_adiabatic.f90:182
thermal_adiabatic::thermal_adiabatic_getmassdensity
real(preal) function, public thermal_adiabatic_getmassdensity(ip, el)
returns homogenized mass density
Definition: thermal_adiabatic.f90:209
thermal_adiabatic::thermal_adiabatic_results
subroutine, public thermal_adiabatic_results(homog, group)
writes results to HDF5 output file
Definition: thermal_adiabatic.f90:235
material::thermalstate
type(tstate), dimension(:), allocatable, public thermalstate
Definition: material.f90:141
lattice
contains lattice structure definitions including Schmid matrices for slip, twin, trans,
Definition: lattice.f90:13
material::thermal_typeinstance
integer, dimension(:), allocatable, public, protected thermal_typeinstance
instance of particular type of each thermal transport
Definition: material.f90:113
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
source_thermal_externalheat
material subroutine for variable heat source
Definition: source_thermal_externalheat.f90:11
material::material_homogenizationmemberat
integer, dimension(:,:), allocatable, target, public material_homogenizationmemberat
position of the element within its homogenization instance
Definition: material.f90:130
thermal_adiabatic::tparameters
Definition: thermal_adiabatic.f90:23
material::thermal_adiabatic_label
character(len= *), parameter, public thermal_adiabatic_label
Definition: material.f90:25
numerics
Managing of parameters related to numerics.
Definition: numerics.f90:10
material::temperaturerate
type(group_float), dimension(:), allocatable, public temperaturerate
temperature change rate field
Definition: material.f90:174
thermal_adiabatic::thermal_adiabatic_updatestate
logical function, dimension(2), public thermal_adiabatic_updatestate(subdt, ip, el)
calculates adiabatic change in temperature based on local heat generation model
Definition: thermal_adiabatic.f90:85
config::config_homogenization
type(tpartitionedstringlist), dimension(:), allocatable, public, protected config_homogenization
Definition: config.f90:23