DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
thermal_conduction.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/thermal_conduction.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/thermal_conduction.f90"
5 !--------------------------------------------------------------------------------------------------
8 !--------------------------------------------------------------------------------------------------
10  use prec
11  use material
12  use config
13  use lattice
14  use results
15  use crystallite
18 
19  implicit none
20  private
21 
22  type :: tparameters
23  character(len=pStringLen), allocatable, dimension(:) :: &
24  output
25  end type tparameters
26 
27  type(tparameters), dimension(:), allocatable :: &
28  param
29 
30  public :: &
38 
39 contains
40 
41 
42 !--------------------------------------------------------------------------------------------------
45 !--------------------------------------------------------------------------------------------------
46 subroutine thermal_conduction_init
47 
48  integer :: ninstance,nofmyhomog,h
49 
50  write(6,'(/,a)') ' <<<+- thermal_'//thermal_conduction_label//' init -+>>>'; flush(6)
51 
52  ninstance = count(thermal_type == thermal_conduction_id)
53  allocate(param(ninstance))
54 
55  do h = 1, size(config_homogenization)
56  if (thermal_type(h) /= thermal_conduction_id) cycle
57  associate(prm => param(thermal_typeinstance(h)),config => config_homogenization(h))
58 
59  prm%output = config%getStrings('(output)',defaultval=emptystringarray)
60 
61  nofmyhomog=count(material_homogenizationat==h)
62  thermalstate(h)%sizeState = 0
63  allocate(thermalstate(h)%state0 (0,nofmyhomog))
64  allocate(thermalstate(h)%subState0(0,nofmyhomog))
65  allocate(thermalstate(h)%state (0,nofmyhomog))
66 
68  deallocate(temperature(h)%p)
69  allocate (temperature(h)%p(nofmyhomog), source=thermal_initialt(h))
70  deallocate(temperaturerate(h)%p)
71  allocate (temperaturerate(h)%p(nofmyhomog), source=0.0_preal)
72 
73  end associate
74  enddo
75 
76 end subroutine thermal_conduction_init
77 
78 
79 !--------------------------------------------------------------------------------------------------
81 !--------------------------------------------------------------------------------------------------
82 subroutine thermal_conduction_getsourceanditstangent(Tdot, dTdot_dT, T, ip, el)
83 
84  integer, intent(in) :: &
85  ip, & !< integration point number
86  el
87  real(preal), intent(in) :: &
88  t
89  real(preal), intent(out) :: &
90  tdot, dtdot_dt
91  real(preal) :: &
92  my_tdot, my_dtdot_dt
93  integer :: &
94  phase, &
95  homog, &
96  offset, &
97  instance, &
98  grain, &
99  source, &
100  constituent
101 
102  homog = material_homogenizationat(el)
103  offset = material_homogenizationmemberat(ip,el)
104  instance = thermal_typeinstance(homog)
105 
106  tdot = 0.0_preal
107  dtdot_dt = 0.0_preal
108  do grain = 1, homogenization_ngrains(homog)
109  phase = material_phaseat(grain,el)
110  constituent = material_phasememberat(grain,ip,el)
111  do source = 1, phase_nsources(phase)
112  select case(phase_source(source,phase))
113  case (source_thermal_dissipation_id)
114  call source_thermal_dissipation_getrateanditstangent(my_tdot, my_dtdot_dt, &
115  crystallite_s(1:3,1:3,grain,ip,el), &
116  crystallite_lp(1:3,1:3,grain,ip,el), &
117  phase)
118 
119  case (source_thermal_externalheat_id)
120  call source_thermal_externalheat_getrateanditstangent(my_tdot, my_dtdot_dt, &
121  phase, constituent)
122  case default
123  my_tdot = 0.0_preal
124  my_dtdot_dt = 0.0_preal
125 
126  end select
127  tdot = tdot + my_tdot
128  dtdot_dt = dtdot_dt + my_dtdot_dt
129  enddo
130  enddo
131 
132  tdot = tdot/real(homogenization_ngrains(homog),preal)
133  dtdot_dt = dtdot_dt/real(homogenization_ngrains(homog),preal)
134 
136 
137 
138 !--------------------------------------------------------------------------------------------------
140 !--------------------------------------------------------------------------------------------------
142 
143  integer, intent(in) :: &
144  ip, & !< integration point number
145  el
146  real(preal), dimension(3,3) :: &
148  integer :: &
149  grain
150 
151 
153  do grain = 1, homogenization_ngrains(material_homogenizationat(el))
155  crystallite_push33toref(grain,ip,el,lattice_thermalconductivity(:,:,material_phaseat(grain,el)))
156  enddo
157 
159  / real(homogenization_ngrains(material_homogenizationat(el)),preal)
160 
162 
163 
164 !--------------------------------------------------------------------------------------------------
166 !--------------------------------------------------------------------------------------------------
168 
169  integer, intent(in) :: &
170  ip, & !< integration point number
171  el
172  real(preal) :: &
174  integer :: &
175  grain
176 
178 
179  do grain = 1, homogenization_ngrains(material_homogenizationat(el))
181  + lattice_specificheat(material_phaseat(grain,el))
182  enddo
183 
185  / real(homogenization_ngrains(material_homogenizationat(el)),preal)
186 
188 
189 
190 !--------------------------------------------------------------------------------------------------
192 !--------------------------------------------------------------------------------------------------
193 function thermal_conduction_getmassdensity(ip,el)
194 
195  integer, intent(in) :: &
196  ip, & !< integration point number
197  el
198  real(preal) :: &
200  integer :: &
201  grain
202 
204 
205 
206  do grain = 1, homogenization_ngrains(material_homogenizationat(el))
208  + lattice_massdensity(material_phaseat(grain,el))
209  enddo
210 
212  / real(homogenization_ngrains(material_homogenizationat(el)),preal)
213 
215 
216 
217 !--------------------------------------------------------------------------------------------------
219 !--------------------------------------------------------------------------------------------------
220 subroutine thermal_conduction_puttemperatureanditsrate(T,Tdot,ip,el)
221 
222  integer, intent(in) :: &
223  ip, & !< integration point number
224  el
225  real(preal), intent(in) :: &
226  t, &
227  tdot
228  integer :: &
229  homog, &
230  offset
231 
232  homog = material_homogenizationat(el)
233  offset = thermalmapping(homog)%p(ip,el)
234  temperature(homog)%p(offset) = t
235  temperaturerate(homog)%p(offset) = tdot
236 
238 
239 
240 !--------------------------------------------------------------------------------------------------
242 !--------------------------------------------------------------------------------------------------
243 subroutine thermal_conduction_results(homog,group)
244 
245  integer, intent(in) :: homog
246  character(len=*), intent(in) :: group
247 
248  integer :: o
249 
250  associate(prm => param(damage_typeinstance(homog)))
251  outputsloop: do o = 1,size(prm%output)
252  select case(trim(prm%output(o)))
253  case('temperature') ! ToDo: should be 'T'
254  call results_writedataset(group,temperature(homog)%p,'T',&
255  'temperature','K')
256  end select
257  enddo outputsloop
258  end associate
259 
260 end subroutine thermal_conduction_results
261 
262 end module thermal_conduction
material::thermalmapping
type(thomogmapping), dimension(:), allocatable, public thermalmapping
mapping for thermal state/fields
Definition: material.f90:170
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_conduction::thermal_conduction_getmassdensity
real(preal) function, public thermal_conduction_getmassdensity(ip, el)
returns homogenized mass density
Definition: thermal_conduction.f90:194
material
Parses material config file, either solverJobName.materialConfig or material.config.
Definition: material.f90:11
thermal_conduction::thermal_conduction_results
subroutine, public thermal_conduction_results(homog, group)
writes results to HDF5 output file
Definition: thermal_conduction.f90:244
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_conduction_label
character(len= *), parameter, public thermal_conduction_label
Definition: material.f90:25
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_conduction::thermal_conduction_getsourceanditstangent
subroutine, public thermal_conduction_getsourceanditstangent(Tdot, dTdot_dT, T, ip, el)
returns heat generation rate
Definition: thermal_conduction.f90:83
thermal_conduction::thermal_conduction_puttemperatureanditsrate
subroutine, public thermal_conduction_puttemperatureanditsrate(T, Tdot, ip, el)
updates thermal state with solution from heat conduction PDE
Definition: thermal_conduction.f90:221
material::thermalstate
type(tstate), dimension(:), allocatable, public thermalstate
Definition: material.f90:141
thermal_conduction::thermal_conduction_getspecificheat
real(preal) function, public thermal_conduction_getspecificheat(ip, el)
returns homogenized specific heat capacity
Definition: thermal_conduction.f90:168
thermal_conduction
material subroutine for temperature evolution from heat conduction
Definition: thermal_conduction.f90:9
thermal_conduction::param
type(tparameters), dimension(:), allocatable param
Definition: thermal_conduction.f90:27
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
thermal_conduction::thermal_conduction_getconductivity
real(preal) function, dimension(3, 3), public thermal_conduction_getconductivity(ip, el)
returns homogenized thermal conductivity in reference configuration
Definition: thermal_conduction.f90:142
material::material_homogenizationmemberat
integer, dimension(:,:), allocatable, target, public material_homogenizationmemberat
position of the element within its homogenization instance
Definition: material.f90:130
thermal_conduction::thermal_conduction_init
subroutine, public thermal_conduction_init
module initialization
Definition: thermal_conduction.f90:47
material::temperaturerate
type(group_float), dimension(:), allocatable, public temperaturerate
temperature change rate field
Definition: material.f90:174
thermal_conduction::tparameters
Definition: thermal_conduction.f90:22
config::config_homogenization
type(tpartitionedstringlist), dimension(:), allocatable, public, protected config_homogenization
Definition: config.f90:23
material::thermal_conduction_id
@, public thermal_conduction_id
Definition: material.f90:87