FrontISTR  5.7.1
Large-scale structural analysis program with finit element method
fstr_ctrl_heat.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
7  use m_fstr
8  use hecmw
10 
11  private :: pc_strupr
12 
13 contains
14 
15  subroutine pc_strupr( s )
16  implicit none
17  character(*) :: s
18  integer :: i, n, a, da
19 
20  n = len_trim(s)
21  da = iachar('a') - iachar('A')
22  do i = 1, n
23  a = iachar(s(i:i))
24  if( a > iachar('Z')) then
25  a = a - da
26  s(i:i) = achar(a)
27  end if
28  end do
29  end subroutine pc_strupr
30 
31 
32 
34  function fstr_ctrl_get_heat( ctrl, dt, etime, dtmin, deltmx, itmax, eps, tpname, beta )
35  implicit none
36  integer(kind=kint) :: ctrl
37  real(kind=kreal),pointer :: dt(:)
38  real(kind=kreal),pointer :: etime(:)
39  real(kind=kreal),pointer :: dtmin(:)
40  real(kind=kreal),pointer :: deltmx(:)
41  integer(kind=kint),pointer :: itmax(:)
42  real(kind=kreal),pointer :: eps(:)
43  character(len=*), intent(out) :: tpname
44  integer(kind=kint) :: fstr_ctrl_get_heat
45  integer(kind=kint) :: result
46  real(kind=kreal) :: beta, t_beta
47 
49 
50  tpname=""
51  if( fstr_ctrl_get_param_ex( ctrl, 'TIMEPOINTS ', '# ', 0, 'S', tpname )/= 0) return
52 
53  ! JP-7
54  if( fstr_ctrl_get_data_array_ex( ctrl, 'rrrrir ', dt, etime, dtmin, deltmx, itmax, eps )/= 0) return
55 
56  beta = -1.0d0
57  t_beta = -1.0d0
58  if( fstr_ctrl_get_param_ex( ctrl, 'BETA ', '# ', 0, 'R', t_beta)/=0 ) return
59  if(0.0d0 <= t_beta .and. t_beta <= 1.0d0) beta = t_beta
60 
62  end function fstr_ctrl_get_heat
63 
65  function fstr_ctrl_get_fixtemp( ctrl, amp, node_grp_name, node_grp_name_len, value )
66  implicit none
67  integer(kind=kint) :: ctrl
68  character(len=HECMW_NAME_LEN) :: amp
69  character(len=HECMW_NAME_LEN) :: node_grp_name(:)
70  integer(kind=kint) :: node_grp_name_len
71  real(kind=kreal), pointer :: value(:)
72  integer(kind=kint) :: fstr_ctrl_get_fixtemp
73 
74  character(len=HECMW_NAME_LEN) :: data_fmt,ss
75 
77 
78  ! JP-8
79  if( fstr_ctrl_get_param_ex( ctrl, 'AMP ', '# ', 0, 'S', amp )/= 0) return
80 
81  write(ss,*) node_grp_name_len
82  write(data_fmt,'(a,a,a)') 'S',trim(adjustl(ss)),'r '
83 
84  fstr_ctrl_get_fixtemp = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, node_grp_name, value )
85 
86  end function fstr_ctrl_get_fixtemp
87 
88 
90  function fstr_ctrl_get_cflux( ctrl, amp, node_grp_name, node_grp_name_len, value )
91  implicit none
92  integer(kind=kint) :: ctrl
93  character(len=HECMW_NAME_LEN) :: amp
94  character(len=HECMW_NAME_LEN) :: node_grp_name(:)
95  integer(kind=kint) :: node_grp_name_len
96  real(kind=kreal),pointer :: value(:)
97  integer(kind=kint) :: fstr_ctrl_get_cflux
98 
99  character(len=HECMW_NAME_LEN) :: data_fmt,ss
100 
102 
103  ! JP-9
104  if( fstr_ctrl_get_param_ex( ctrl, 'AMP ', '# ', 0, 'S', amp )/= 0) return
105 
106  write(ss,*) node_grp_name_len
107  write(data_fmt,'(a,a,a)') 'S',trim(adjustl(ss)),'r '
108 
109  fstr_ctrl_get_cflux = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, node_grp_name, value )
110 
111  end function fstr_ctrl_get_cflux
112 
114  function fstr_ctrl_get_dflux( ctrl, amp, elem_grp_name, elem_grp_name_len, load_type, value )
115  implicit none
116  integer(kind=kint) :: ctrl
117  character(len=HECMW_NAME_LEN) :: amp
118  character(len=HECMW_NAME_LEN) :: elem_grp_name(:)
119  integer(kind=kint) :: elem_grp_name_len
120  integer(kind=kint),pointer :: load_type(:)
121  real(kind=kreal),pointer :: value(:)
122  integer(kind=kint) :: fstr_ctrl_get_dflux
123 
124  integer(kind=kint), parameter :: type_name_size = 5
125  integer(kind=kint) :: i, n
126  character(len=HECMW_NAME_LEN) :: data_fmt,s1,s2
127  character(len=type_name_size),pointer :: type_name_list(:)
128  integer(kind=kint) :: rcode
129  integer(kind=kint) :: lid = -1
130 
132 
133  ! JP-10
134  if( fstr_ctrl_get_param_ex( ctrl, 'AMP ', '# ', 0, 'S', amp )/= 0) return
135 
136  write(s1,*) elem_grp_name_len
137  write(s2,*) type_name_size
138  write(data_fmt,'(a,a,a,a,a)') 'S',trim(adjustl(s1)),'S',trim(adjustl(s2)),'r '
139 
140  n = fstr_ctrl_get_data_line_n(ctrl)
141  allocate( type_name_list(n) )
142 
143  rcode = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, elem_grp_name, type_name_list, value)
144 
145  if( rcode /= 0 ) then
146  deallocate( type_name_list )
147  return
148  end if
149 
150  do i = 1, n
151  lid = -1;
152  call pc_strupr( type_name_list(i) )
153  if( type_name_list(i)(1:2) == 'BF' ) then; lid = 0
154  else if( type_name_list(i)(1:2) == 'S0' ) then; lid = 1
155  else if( type_name_list(i)(1:2) == 'S1' ) then; lid = 1
156  else if( type_name_list(i)(1:2) == 'S2' ) then; lid = 2
157  else if( type_name_list(i)(1:2) == 'S3' ) then; lid = 3
158  else if( type_name_list(i)(1:2) == 'S4' ) then; lid = 4
159  else if( type_name_list(i)(1:2) == 'S5' ) then; lid = 5
160  else if( type_name_list(i)(1:2) == 'S6' ) then; lid = 6
161  end if
162  if( lid < 0 ) then
163  write(ilog,*) 'Error : !DFLUX : Load type ',type_name_list(i),' is unknown'
164  deallocate( type_name_list )
165  return
166  end if
167  load_type(i) = lid
168  end do
169 
170  deallocate( type_name_list )
172  end function fstr_ctrl_get_dflux
173 
174  !* ----------------------------------------------------------------------------------------------- *!
176  !* ----------------------------------------------------------------------------------------------- *!
177 
178  function fstr_ctrl_get_sflux( ctrl, amp, surface_grp_name, surface_grp_name_len, value )
179  implicit none
180  integer(kind=kint) :: ctrl
181  character(len=HECMW_NAME_LEN) :: amp
182  character(len=HECMW_NAME_LEN) :: surface_grp_name(:)
183  integer(kind=kint) :: surface_grp_name_len
184  real(kind=kreal),pointer :: value(:)
185  integer(kind=kint) :: fstr_ctrl_get_sflux
186 
187  character(len=HECMW_NAME_LEN) :: data_fmt,ss
188 
190 
191  ! JP-11
192  if( fstr_ctrl_get_param_ex( ctrl, 'AMP ', '# ', 0, 'S', amp )/= 0) return
193 
194  write(ss,*) surface_grp_name_len
195  write(data_fmt,'(a,a,a)') 'S',trim(adjustl(ss)),'r '
196 
197  fstr_ctrl_get_sflux = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, surface_grp_name, value )
198  end function fstr_ctrl_get_sflux
199 
200 
201  !* ----------------------------------------------------------------------------------------------- *!
203  !* ----------------------------------------------------------------------------------------------- *!
204 
205  function fstr_ctrl_get_film( ctrl, amp1, amp2, elem_grp_name, elem_grp_name_len, load_type, value, sink)
206  implicit none
207  integer(kind=kint) :: ctrl
208  character(len=HECMW_NAME_LEN) :: amp1
209  character(len=HECMW_NAME_LEN) :: amp2
210  character(len=HECMW_NAME_LEN) :: elem_grp_name(:)
211  integer(kind=kint) :: elem_grp_name_len
212  integer(kind=kint),pointer :: load_type(:)
213  real(kind=kreal),pointer :: value(:)
214  real(kind=kreal),pointer :: sink(:)
215  integer(kind=kint) :: fstr_ctrl_get_film
216 
217  integer(kind=kint),parameter :: type_name_size = 5
218  integer(kind=kint) :: i, n
219  character(len=HECMW_NAME_LEN) :: data_fmt,s1,s2
220  character(len=type_name_size),pointer :: type_name_list(:)
221  integer(kind=kint) :: lid
222  integer(kind=kint) :: rcode
223 
224  fstr_ctrl_get_film = -1
225 
226  ! JP-12
227  if( fstr_ctrl_get_param_ex( ctrl, 'AMP1 ', '# ', 0, 'S', amp1 )/= 0) return
228  if( fstr_ctrl_get_param_ex( ctrl, 'AMP2 ', '# ', 0, 'S', amp2 )/= 0) return
229 
230  write(s1,*) elem_grp_name_len
231  write(s2,*) type_name_size
232  write(data_fmt,'(a,a,a,a,a)') 'S',trim(adjustl(s1)),'S',trim(adjustl(s2)),'Rr '
233 
234  n = fstr_ctrl_get_data_line_n(ctrl)
235  allocate( type_name_list(n) )
236 
237  rcode = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, elem_grp_name, type_name_list, value, sink)
238 
239  if( rcode /= 0 ) then
240  deallocate( type_name_list )
241  return
242  end if
243 
244  do i = 1, n
245  lid = -1;
246  call pc_strupr( type_name_list(i) )
247  if( type_name_list(i)(1:2) == 'F0' ) then; lid = 1
248  else if( type_name_list(i)(1:2) == 'F1' ) then; lid = 1
249  else if( type_name_list(i)(1:2) == 'F2' ) then; lid = 2
250  else if( type_name_list(i)(1:2) == 'F3' ) then; lid = 3
251  else if( type_name_list(i)(1:2) == 'F4' ) then; lid = 4
252  else if( type_name_list(i)(1:2) == 'F5' ) then; lid = 5
253  else if( type_name_list(i)(1:2) == 'F6' ) then; lid = 6
254  end if
255  if( lid < 0 ) then
256  write(ilog,*) 'Error : !FILM : Load type ',type_name_list(i),' is unknown'
257  deallocate( type_name_list )
258  return
259  end if
260  load_type(i) = lid
261  end do
262 
263  deallocate( type_name_list )
264 
266  end function fstr_ctrl_get_film
267 
268  !* ----------------------------------------------------------------------------------------------- *!
270  !* ----------------------------------------------------------------------------------------------- *!
271 
272  function fstr_ctrl_get_sfilm( ctrl, amp1, amp2, surface_grp_name, surface_grp_name_len, value, sink)
273  implicit none
274  integer(kind=kint) :: ctrl
275  character(len=HECMW_NAME_LEN) :: amp1
276  character(len=HECMW_NAME_LEN) :: amp2
277  character(len=HECMW_NAME_LEN) :: surface_grp_name(:)
278  integer(kind=kint) :: surface_grp_name_len
279  real(kind=kreal),pointer :: value(:)
280  real(kind=kreal),pointer :: sink(:)
281  integer(kind=kint) :: fstr_ctrl_get_sfilm
282 
283  character(len=HECMW_NAME_LEN) :: data_fmt,ss
284 
286 
287  ! JP-13
288  if( fstr_ctrl_get_param_ex( ctrl, 'AMP1 ', '# ', 0, 'S', amp1 )/= 0) return
289  if( fstr_ctrl_get_param_ex( ctrl, 'AMP2 ', '# ', 0, 'S', amp2 )/= 0) return
290 
291  write(ss,*) surface_grp_name_len
292  write(data_fmt,'(a,a,a)') 'S',trim(adjustl(ss)),'Rr '
293 
294  fstr_ctrl_get_sfilm = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, surface_grp_name, value, sink )
295  end function fstr_ctrl_get_sfilm
296 
297 
298  !* ----------------------------------------------------------------------------------------------- *!
300  !* ----------------------------------------------------------------------------------------------- *!
301 
302  function fstr_ctrl_get_radiate( ctrl, amp1, amp2, elem_grp_name, elem_grp_name_len, load_type, value, sink)
303  implicit none
304  integer(kind=kint) :: ctrl
305  character(len=HECMW_NAME_LEN) :: amp1
306  character(len=HECMW_NAME_LEN) :: amp2
307  character(len=HECMW_NAME_LEN) :: elem_grp_name(:)
308  integer(kind=kint) :: elem_grp_name_len
309  integer(kind=kint),pointer :: load_type(:)
310  real(kind=kreal),pointer :: value(:)
311  real(kind=kreal),pointer :: sink(:)
312  integer(kind=kint) :: fstr_ctrl_get_radiate
313 
314  integer(kind=kint),parameter :: type_name_size = 5
315  integer(kind=kint) :: i, n
316  character(len=HECMW_NAME_LEN) :: data_fmt,s1,s2
317  character(len=type_name_size),pointer :: type_name_list(:)
318  integer(kind=kint) :: lid
319  integer(kind=kint) :: rcode
320 
322 
323  ! JP-14
324  if( fstr_ctrl_get_param_ex( ctrl, 'AMP1 ', '# ', 0, 'S', amp1 )/= 0) return
325  if( fstr_ctrl_get_param_ex( ctrl, 'AMP2 ', '# ', 0, 'S', amp2 )/= 0) return
326 
327  write(s1,*) elem_grp_name_len
328  write(s2,*) type_name_size
329  write(data_fmt,'(a,a,a,a,a)') 'S',trim(adjustl(s1)),'S',trim(adjustl(s2)),'Rr '
330 
331  n = fstr_ctrl_get_data_line_n(ctrl)
332  allocate( type_name_list(n) )
333 
334  rcode = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, elem_grp_name ,type_name_list, value, sink)
335 
336  if( rcode /= 0 ) then
337  deallocate( type_name_list )
338  return
339  end if
340 
341  do i = 1, n
342  lid = -1;
343  call pc_strupr( type_name_list(i) )
344  if( type_name_list(i)(1:2) == 'R0' ) then; lid = 1
345  else if( type_name_list(i)(1:2) == 'R1' ) then; lid = 1
346  else if( type_name_list(i)(1:2) == 'R2' ) then; lid = 2
347  else if( type_name_list(i)(1:2) == 'R3' ) then; lid = 3
348  else if( type_name_list(i)(1:2) == 'R4' ) then; lid = 4
349  else if( type_name_list(i)(1:2) == 'R5' ) then; lid = 5
350  else if( type_name_list(i)(1:2) == 'R6' ) then; lid = 6
351  end if
352  if( lid < 0 ) then
353  write(ilog,*) 'Error : !RADIATE : Load type ',type_name_list(i),' is unknown'
354  deallocate( type_name_list )
355  return
356  end if
357  load_type(i) = lid
358  end do
359 
360  deallocate( type_name_list )
362  end function fstr_ctrl_get_radiate
363 
364 
365  !* ----------------------------------------------------------------------------------------------- *!
367  !* ----------------------------------------------------------------------------------------------- *!
368 
369  function fstr_ctrl_get_sradiate( ctrl, amp1, amp2, surface_grp_name, surface_grp_name_len, value, sink)
370  implicit none
371  integer(kind=kint) :: ctrl
372  character(len=HECMW_NAME_LEN) :: amp1
373  character(len=HECMW_NAME_LEN) :: amp2
374  character(len=HECMW_NAME_LEN) :: surface_grp_name(:)
375  integer(kind=kint) :: surface_grp_name_len
376  real(kind=kreal),pointer :: value(:)
377  real(kind=kreal),pointer :: sink(:)
378  integer(kind=kint) :: fstr_ctrl_get_sradiate
379 
380  character(len=HECMW_NAME_LEN) :: data_fmt
381  character(len=HECMW_NAME_LEN) :: s1
382 
384 
385  ! JP-15
386  if( fstr_ctrl_get_param_ex( ctrl, 'AMP1 ', '# ', 0, 'S', amp1 )/= 0) return
387  if( fstr_ctrl_get_param_ex( ctrl, 'AMP2 ', '# ', 0, 'S', amp2 )/= 0) return
388 
389  write(s1,*) surface_grp_name_len;
390  write(data_fmt,'(a,a,a)') 'S', trim(adjustl(s1)), 'Rr '
391 
392  fstr_ctrl_get_sradiate = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, surface_grp_name, value, sink )
393 
394  end function fstr_ctrl_get_sradiate
395 
396  !* ----------------------------------------------------------------------------------------------- *!
398  !* ----------------------------------------------------------------------------------------------- *!
399 
400  function fstr_ctrl_get_weldline( ctrl, hecMESH, grp_name_len, weldline )
401  use fstr_setup_util
402  implicit none
403  integer(kind=kint), intent(in) :: ctrl
404  type(hecmwst_local_mesh), intent(in) :: hecmesh
405  integer(kind=kint), intent(in) :: grp_name_len
406  type(tweldline), intent(inout) :: weldline
407  integer(kind=kint) :: fstr_ctrl_get_weldline
408 
409  character(len=HECMW_NAME_LEN) :: data_fmt
410  character(len=HECMW_NAME_LEN) :: s1, grp_id_name(1)
411  integer :: grp_id(1)
412 
414  if( fstr_ctrl_get_data_ex( ctrl, 1, 'RRRR ', weldline%I, weldline%U, weldline%coe, weldline%v )/=0 ) return
415  write(s1,*) grp_name_len
416  write(data_fmt,'(a,a,a)') 'S', trim(adjustl(s1)), 'IRRRR '
417  if( fstr_ctrl_get_data_ex( ctrl, 2, data_fmt, grp_id_name(1), weldline%xyz, weldline%n1, &
418  weldline%n2, weldline%distol, weldline%tstart )/=0 ) return
419  call elem_grp_name_to_id( hecmesh, 'WELD_LINE ', 1, grp_id_name, grp_id )
420  weldline%egrpid = grp_id(1)
421 
423  end function fstr_ctrl_get_weldline
424 
425  !* ----------------------------------------------------------------------------------------------- *!
426 end module fstr_ctrl_heat
int fstr_ctrl_get_param_ex(int *ctrl, const char *param_name, const char *value_list, int *necessity, char *type, void *val)
int fstr_ctrl_get_data_array_ex(int *ctrl, const char *format,...)
int fstr_ctrl_get_data_ex(int *ctrl, int *line_no, const char *format,...)
This module contains control file data obtaining functions for heat conductive analysis.
integer(kind=kint) function fstr_ctrl_get_dflux(ctrl, amp, elem_grp_name, elem_grp_name_len, load_type, value)
Read in !DFLUX (heat)
integer(kind=kint) function fstr_ctrl_get_sflux(ctrl, amp, surface_grp_name, surface_grp_name_len, value)
Read in !SFLUX (heat)
integer(kind=kint) function fstr_ctrl_get_weldline(ctrl, hecMESH, grp_name_len, weldline)
Read in !WELD_LINE (heat)
integer(kind=kint) function fstr_ctrl_get_heat(ctrl, dt, etime, dtmin, deltmx, itmax, eps, tpname, beta)
Read in !HEAT.
integer(kind=kint) function fstr_ctrl_get_film(ctrl, amp1, amp2, elem_grp_name, elem_grp_name_len, load_type, value, sink)
Read in !FILM (heat)
integer(kind=kint) function fstr_ctrl_get_radiate(ctrl, amp1, amp2, elem_grp_name, elem_grp_name_len, load_type, value, sink)
Read in !RADIATE (heat)
integer(kind=kint) function fstr_ctrl_get_cflux(ctrl, amp, node_grp_name, node_grp_name_len, value)
Read in !CFLUX (heat)
integer(kind=kint) function fstr_ctrl_get_fixtemp(ctrl, amp, node_grp_name, node_grp_name_len, value)
Read in !FIXTEMP.
integer(kind=kint) function fstr_ctrl_get_sfilm(ctrl, amp1, amp2, surface_grp_name, surface_grp_name_len, value, sink)
Read in !SFILM (heat)
integer(kind=kint) function fstr_ctrl_get_sradiate(ctrl, amp1, amp2, surface_grp_name, surface_grp_name_len, value, sink)
Read in !SRADIATE (heat)
This module contains auxiliary functions in calculation setup.
subroutine elem_grp_name_to_id(hecMESH, header_name, n, grp_id_name, grp_ID)
Definition: hecmw.f90:6
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
real(kind=kreal) eps
Definition: m_fstr.f90:142
real(kind=kreal) etime
Definition: m_fstr.f90:140
integer(kind=kint) itmax
Definition: m_fstr.f90:141
integer(kind=kint), parameter ilog
FILE HANDLER.
Definition: m_fstr.f90:107
real(kind=kreal) dt
ANALYSIS CONTROL for NLGEOM and HEAT.
Definition: m_fstr.f90:139
-1:not relation, >1:index of coupled_node
Definition: m_fstr.f90:629