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