FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
fstr_ctrl_common.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 !-------------------------------------------------------------------------------
6 
8  use m_fstr
9  use hecmw
10  use mcontact
11  use m_timepoint
12 
13  implicit none
14 
15  include 'fstr_ctrl_util_f.inc'
16 
17  private :: pc_strupr
18 
19 contains
20 
21  subroutine pc_strupr( s )
22  implicit none
23  character(*) :: s
24  integer :: i, n, a, da
25 
26  n = len_trim(s)
27  da = iachar('a') - iachar('A')
28  do i = 1, n
29  a = iachar(s(i:i))
30  if( a > iachar('Z')) then
31  a = a - da
32  s(i:i) = achar(a)
33  end if
34  end do
35  end subroutine pc_strupr
36 
37 
39  function fstr_ctrl_get_solution( ctrl, type, nlgeom )
40  integer(kind=kint) :: ctrl
41  integer(kind=kint) :: type
42  logical :: nlgeom
43  integer(kind=kint) :: fstr_ctrl_get_solution
44 
45  integer(kind=kint) :: ipt
46  character(len=80) :: s
47 
49 
50  s = 'ELEMCHECK,STATIC,EIGEN,HEAT,DYNAMIC,NLSTATIC,STATICEIGEN,NZPROF '
51  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', s, 1, 'P', type )/= 0) return
52  type = type -1
53 
54  ipt=0
55  if( fstr_ctrl_get_param_ex( ctrl, 'NONLINEAR ', '# ', 0, 'E', ipt )/= 0) return
56  if( ipt/=0 .and. ( type == kststatic .or. type == kstdynamic )) nlgeom = .true.
57 
58  if( type == 5 ) then !if type == NLSTATIC
59  type = kststatic
60  nlgeom = .true.
61  end if
62  if( type == kststaticeigen ) nlgeom = .true.
63 
65  end function fstr_ctrl_get_solution
66 
68  function fstr_ctrl_get_nonlinear_solver( ctrl, method )
69  integer(kind=kint) :: ctrl
70  integer(kind=kint) :: method
71  integer(kind=kint) :: fstr_ctrl_get_nonlinear_solver
72 
73  integer(kind=kint) :: ipt
74  character(len=80) :: s
75 
77 
78  s = 'NEWTON,QUASINEWTON '
79  if( fstr_ctrl_get_param_ex( ctrl, 'METHOD ', s, 1, 'P', method )/= 0) return
82 
84  function fstr_ctrl_get_solver( ctrl, method, precond, nset, iterlog, timelog, steplog, nier, &
85  iterpremax, nrest, nBFGS, scaling, &
86  dumptype, dumpexit, usejad, ncolor_in, mpc_method, estcond, method2, recyclepre, &
87  solver_opt, contact_elim, &
88  resid, singma_diag, sigma, thresh, filter )
89  integer(kind=kint) :: ctrl
90  integer(kind=kint) :: method
91  integer(kind=kint) :: precond
92  integer(kind=kint) :: nset
93  integer(kind=kint) :: iterlog
94  integer(kind=kint) :: timelog
95  integer(kind=kint) :: steplog
96  integer(kind=kint) :: nier
97  integer(kind=kint) :: iterpremax
98  integer(kind=kint) :: nrest
99  integer(kind=kint) :: nbfgs
100  integer(kind=kint) :: scaling
101  integer(kind=kint) :: dumptype
102  integer(kind=kint) :: dumpexit
103  integer(kind=kint) :: usejad
104  integer(kind=kint) :: ncolor_in
105  integer(kind=kint) :: mpc_method
106  integer(kind=kint) :: estcond
107  integer(kind=kint) :: method2
108  integer(kind=kint) :: recyclepre
109  integer(kind=kint) :: solver_opt(10)
110  integer(kind=kint) :: contact_elim
111  real(kind=kreal) :: resid
112  real(kind=kreal) :: singma_diag
113  real(kind=kreal) :: sigma
114  real(kind=kreal) :: thresh
115  real(kind=kreal) :: filter
116  integer(kind=kint) :: fstr_ctrl_get_solver
117 
118  character(92) :: mlist = '1,2,3,4,101,CG,BiCGSTAB,GMRES,GPBiCG,GMRESR,GMRESREN,DIRECT,DIRECTmkl,DIRECTlag,MUMPS,MKL '
119  !character(92) :: mlist = '1,2,3,4,5,101,CG,BiCGSTAB,GMRES,GPBiCG,DIRECT,DIRECTmkl,DIRECTlag,MUMPS,MKL '
120  character(24) :: dlist = '0,1,2,3,NONE,MM,CSR,BSR '
121 
122  integer(kind=kint) :: number_number = 5
123  integer(kind=kint) :: indirect_number = 6 ! GMRESR and GMRESREN need to be added
124  integer(kind=kint) :: iter, time, sclg, dmpt, dmpx, usjd, step
125 
127 
128  iter = iterlog+1
129  time = timelog+1
130  step = steplog+1
131  sclg = scaling+1
132  dmpt = dumptype+1
133  dmpx = dumpexit+1
134  usjd = usejad+1
135  !* parameter in header line -----------------------------------------------------------------*!
136 
137  ! JP-0
138  if( fstr_ctrl_get_param_ex( ctrl, 'METHOD ', mlist, 1, 'P', method ) /= 0) return
139  if( fstr_ctrl_get_param_ex( ctrl, 'PRECOND ', '1,2,3,4,5,6,7,8,9,10,11,12,20,21,30,31,32 ' ,0, 'I', precond ) /= 0) return
140  if( fstr_ctrl_get_param_ex( ctrl, 'NSET ', '0,-1,+1 ', 0, 'I', nset ) /= 0) return
141  if( fstr_ctrl_get_param_ex( ctrl, 'ITERLOG ', 'NO,YES ', 0, 'P', iter ) /= 0) return
142  if( fstr_ctrl_get_param_ex( ctrl, 'TIMELOG ', 'NO,YES,VERBOSE ', 0, 'P', time ) /= 0) return
143  if( fstr_ctrl_get_param_ex( ctrl, 'STEPLOG ', 'NO,YES ', 0, 'P', step ) /= 0) return
144  if( fstr_ctrl_get_param_ex( ctrl, 'SCALING ', 'NO,YES ', 0, 'P', sclg ) /= 0) return
145  if( fstr_ctrl_get_param_ex( ctrl, 'DUMPTYPE ', dlist, 0, 'P', dmpt ) /= 0) return
146  if( fstr_ctrl_get_param_ex( ctrl, 'DUMPEXIT ','NO,YES ', 0, 'P', dmpx ) /= 0) return
147  if( fstr_ctrl_get_param_ex( ctrl, 'USEJAD ' ,'NO,YES ', 0, 'P', usjd ) /= 0) return
148  if( fstr_ctrl_get_param_ex( ctrl, 'MPCMETHOD ','# ', 0, 'I',mpc_method) /= 0) return
149  if( fstr_ctrl_get_param_ex( ctrl, 'ESTCOND ' ,'# ', 0, 'I',estcond ) /= 0) return
150  if( fstr_ctrl_get_param_ex( ctrl, 'METHOD2 ', mlist, 0, 'P', method2 ) /= 0) return
151  if( fstr_ctrl_get_param_ex( ctrl, 'CONTACT_ELIM ','# ', 0, 'I',contact_elim ) /= 0) return
152  ! JP-1
153  if( method > number_number ) then ! JP-2
154  method = method - number_number
155  if( method > indirect_number ) then
156  ! JP-3
157  method = method - indirect_number + 100
158  if( method == 103 ) method = 101 ! DIRECTlag => DIRECT
159  if( method == 105 ) method = 102 ! MKL => DIRECTmkl
160  end if
161  end if
162  if( method2 > number_number ) then ! JP-2
163  method2 = method2 - number_number
164  if( method2 > indirect_number ) then
165  ! JP-3
166  method2 = method2 - indirect_number + 100
167  end if
168  end if
169 
170  dumptype = dmpt - 1
171  if( dumptype >= 4 ) then
172  dumptype = dumptype - 4
173  end if
174 
175  !* data --------------------------------------------------------------------------------------- *!
176  ! JP-4
177  if( fstr_ctrl_get_data_ex( ctrl, 1, 'iiiiii ', nier, iterpremax, nrest, ncolor_in, recyclepre, nbfgs )/= 0) return
178  if( fstr_ctrl_get_data_ex( ctrl, 2, 'rrr ', resid, singma_diag, sigma )/= 0) return
179 
180  if( precond == 20 .or. precond == 21) then
181  if( fstr_ctrl_get_data_ex( ctrl, 3, 'rr ', thresh, filter)/= 0) return
182  else if( precond == 5 ) then
183  if( fstr_ctrl_get_data_ex( ctrl, 3, 'iiiiiiiiii ', &
184  solver_opt(1), solver_opt(2), solver_opt(3), solver_opt(4), solver_opt(5), &
185  solver_opt(6), solver_opt(7), solver_opt(8), solver_opt(9), solver_opt(10) )/= 0) return
186  else if( method == 101 ) then
187  if( fstr_ctrl_get_data_ex( ctrl, 3, 'i ', solver_opt(1) )/= 0) return
188  end if
189 
190  iterlog = iter -1
191  timelog = time -1
192  steplog = step -1
193  scaling = sclg -1
194  dumpexit = dmpx -1
195  usejad = usjd -1
196 
198 
199  end function fstr_ctrl_get_solver
200 
201 
203  function fstr_ctrl_get_step( ctrl, amp, iproc )
204  integer(kind=kint) :: ctrl
205  character(len=HECMW_NAME_LEN) :: amp
206  integer(kind=kint) :: iproc
207  integer(kind=kint) :: fstr_ctrl_get_step
208 
209  integer(kind=kint) :: ipt = 0
210  integer(kind=kint) :: ip = 0
211 
212  fstr_ctrl_get_step = -1
213 
214  if( fstr_ctrl_get_param_ex( ctrl, 'AMP ', '# ', 0, 'S', amp )/= 0) return
215  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', 'STANDARD,NLGEOM ', 0, 'P', ipt )/= 0) return
216  if( fstr_ctrl_get_param_ex( ctrl, 'NLGEOM ', '# ', 0, 'E', ip )/= 0) return
217 
218  if( ipt == 2 .or. ip == 1 ) iproc = 1
219 
221 
222  end function fstr_ctrl_get_step
223 
225  logical function fstr_ctrl_get_istep( ctrl, hecMESH, steps, tpname, apname )
227  use m_step
228  integer(kind=kint), intent(in) :: ctrl
229  type (hecmwst_local_mesh), intent(in) :: hecmesh
230  type(step_info), intent(out) :: steps
231  character(len=*), intent(out) :: tpname
232  character(len=*), intent(out) :: apname
233 
234  character(len=HECMW_NAME_LEN) :: data_fmt,ss, data_fmt1
235  character(len=HECMW_NAME_LEN) :: amp
236  character(len=HECMW_NAME_LEN) :: header_name
237  integer(kind=kint) :: bcid
238  integer(kind=kint) :: i, n, sn, ierr
239  integer(kind=kint) :: bc_n, load_n, contact_n
240  real(kind=kreal) :: fn, f1, f2, f3
241 
242  fstr_ctrl_get_istep = .false.
243 
244  write(ss,*) hecmw_name_len
245  write( data_fmt, '(a,a,a)') 'S', trim(adjustl(ss)), 'I '
246  write( data_fmt1, '(a,a,a)') 'S', trim(adjustl(ss)),'rrr '
247 
248  call init_stepinfo(steps)
249  steps%solution = stepstatic
250  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', 'STATIC,VISCO ', 0, 'P', steps%solution )/= 0) return
251  steps%inc_type = stepfixedinc
252  if( fstr_ctrl_get_param_ex( ctrl, 'INC_TYPE ', 'FIXED,AUTO ', 0, 'P', steps%inc_type )/= 0) return
253  if( fstr_ctrl_get_param_ex( ctrl, 'SUBSTEPS ', '# ', 0, 'I', steps%num_substep )/= 0) return
254  steps%initdt = 1.d0/steps%num_substep
255  if( fstr_ctrl_get_param_ex( ctrl, 'ITMAX ', '# ', 0, 'I', steps%max_iter )/= 0) return
256  if( fstr_ctrl_get_param_ex( ctrl, 'MAXITER ', '# ', 0, 'I', steps%max_iter )/= 0) return
257  if( fstr_ctrl_get_param_ex( ctrl, 'MAXCONTITER ', '# ', 0, 'I', steps%max_contiter )/= 0) return
258  if( fstr_ctrl_get_param_ex( ctrl, 'CONVERG ', '# ', 0, 'R', steps%converg )/= 0) return
259  if( fstr_ctrl_get_param_ex( ctrl, 'CONVERG_LAG ', '# ', 0, 'R', steps%converg_lag )/= 0) return
260  if( fstr_ctrl_get_param_ex( ctrl, 'CONVERG_DDISP ', '# ', 0, 'R', steps%converg_ddisp )/= 0) return
261  if( fstr_ctrl_get_param_ex( ctrl, 'MAXRES ', '# ', 0, 'R', steps%maxres )/= 0) return
262  amp = ""
263  if( fstr_ctrl_get_param_ex( ctrl, 'AMP ', '# ', 0, 'S', amp )/= 0) return
264  if( len( trim(amp) )>0 ) then
265  call amp_name_to_id( hecmesh, '!STEP', amp, steps%amp_id )
266  endif
267  tpname=""
268  if( fstr_ctrl_get_param_ex( ctrl, 'TIMEPOINTS ', '# ', 0, 'S', tpname )/= 0) return
269  apname=""
270  if( fstr_ctrl_get_param_ex( ctrl, 'AUTOINCPARAM ', '# ', 0, 'S', apname )/= 0) return
271 
272  n = fstr_ctrl_get_data_line_n( ctrl )
273  if( n == 0 ) then
274  fstr_ctrl_get_istep = .true.; return
275  endif
276 
277  f2 = steps%mindt
278  f3 = steps%maxdt
279  if( fstr_ctrl_get_data_ex( ctrl, 1, data_fmt1, ss, f1, f2, f3 )/= 0) return
280  read( ss, * , iostat=ierr ) fn
281  sn=1
282  if( ierr==0 ) then
283  steps%initdt = fn
284  steps%elapsetime = f1
285  if( steps%inc_type == stepautoinc ) then
286  steps%mindt = min(f2,steps%initdt)
287  steps%maxdt = f3
288  endif
289  steps%num_substep = max(int((f1+0.999999999d0*fn)/fn),steps%num_substep)
290  !if( mod(f1,fn)/=0 ) steps%num_substep =steps%num_substep+1
291  sn = 2
292  endif
293 
294  bc_n = 0
295  load_n = 0
296  contact_n = 0
297  do i=sn,n
298  if( fstr_ctrl_get_data_ex( ctrl, i, data_fmt, header_name, bcid )/= 0) return
299  if( trim(header_name) == 'BOUNDARY' ) then
300  bc_n = bc_n + 1
301  else if( trim(header_name) == 'LOAD' ) then
302  load_n = load_n +1
303  else if( trim(header_name) == 'CONTACT' ) then
304  contact_n = contact_n+1
305  else if( trim(header_name) == 'TEMPERATURE' ) then
306  ! steps%Temperature = .true.
307  endif
308  end do
309 
310  if( bc_n>0 ) allocate( steps%Boundary(bc_n) )
311  if( load_n>0 ) allocate( steps%Load(load_n) )
312  if( contact_n>0 ) allocate( steps%Contact(contact_n) )
313 
314  bc_n = 0
315  load_n = 0
316  contact_n = 0
317  do i=sn,n
318  if( fstr_ctrl_get_data_ex( ctrl, i, data_fmt, header_name, bcid )/= 0) return
319  if( trim(header_name) == 'BOUNDARY' ) then
320  bc_n = bc_n + 1
321  steps%Boundary(bc_n) = bcid
322  else if( trim(header_name) == 'LOAD' ) then
323  load_n = load_n +1
324  steps%Load(load_n) = bcid
325  else if( trim(header_name) == 'CONTACT' ) then
326  contact_n = contact_n+1
327  steps%Contact(contact_n) = bcid
328  endif
329  end do
330 
331  fstr_ctrl_get_istep = .true.
332  end function fstr_ctrl_get_istep
333 
335  integer function fstr_ctrl_get_section( ctrl, hecMESH, sections )
337  integer(kind=kint), intent(in) :: ctrl
338  type (hecmwst_local_mesh), intent(inout) :: hecmesh
339  type (tsection), pointer, intent(inout) :: sections(:)
340 
341  integer(kind=kint) :: j, k, sect_id, ori_id, elemopt
342  integer(kind=kint),save :: cache = 1
343  character(len=HECMW_NAME_LEN) :: sect_orien
344  character(19) :: form341list = 'FI,SELECTIVE_ESNS '
345  character(16) :: form361list = 'FI,BBAR,IC,FBAR '
346 
348 
349  if( fstr_ctrl_get_param_ex( ctrl, 'SECNUM ', '# ', 1, 'I', sect_id )/= 0) return
350  if( sect_id > hecmesh%section%n_sect ) return
351 
352  elemopt = 0
353  if( fstr_ctrl_get_param_ex( ctrl, 'FORM341 ', form341list, 0, 'P', elemopt )/= 0) return
354  if( elemopt > 0 ) sections(sect_id)%elemopt341 = elemopt
355 
356  elemopt = 0
357  if( fstr_ctrl_get_param_ex( ctrl, 'FORM361 ', form361list, 0, 'P', elemopt )/= 0) return
358  if( elemopt > 0 ) sections(sect_id)%elemopt361 = elemopt
359 
360  ! sectional orientation ID
361  hecmesh%section%sect_orien_ID(sect_id) = -1
362  if( fstr_ctrl_get_param_ex( ctrl, 'ORIENTATION ', '# ', 0, 'S', sect_orien )/= 0) return
363 
364  if( associated(g_localcoordsys) ) then
365  call fstr_strupr(sect_orien)
366  k = size(g_localcoordsys)
367 
368  if(cache < k)then
369  if( sect_orien == g_localcoordsys(cache)%sys_name ) then
370  hecmesh%section%sect_orien_ID(sect_id) = cache
371  cache = cache + 1
373  return
374  endif
375  endif
376 
377  do j=1, k
378  if( sect_orien == g_localcoordsys(j)%sys_name ) then
379  hecmesh%section%sect_orien_ID(sect_id) = j
380  cache = j + 1
381  exit
382  endif
383  enddo
384  endif
385 
387 
388  end function fstr_ctrl_get_section
389 
390 
392  function fstr_ctrl_get_write( ctrl, res, visual, femap )
393  integer(kind=kint) :: ctrl
394  integer(kind=kint) :: res
395  integer(kind=kint) :: visual
396  integer(kind=kint) :: femap
397  integer(kind=kint) :: fstr_ctrl_get_write
398 
400 
401  ! JP-6
402  if( fstr_ctrl_get_param_ex( ctrl, 'RESULT ', '# ', 0, 'E', res )/= 0) return
403  if( fstr_ctrl_get_param_ex( ctrl, 'VISUAL ', '# ', 0, 'E', visual )/= 0) return
404  if( fstr_ctrl_get_param_ex( ctrl, 'FEMAP ', '# ', 0, 'E', femap )/= 0) return
405 
407 
408  end function fstr_ctrl_get_write
409 
411  function fstr_ctrl_get_echo( ctrl, echo )
412  integer(kind=kint) :: ctrl
413  integer(kind=kint) :: echo
414  integer(kind=kint) :: fstr_ctrl_get_echo
415 
416  echo = kon;
417 
419 
420  end function fstr_ctrl_get_echo
421 
423  function fstr_ctrl_get_couple( ctrl, fg_type, fg_first, fg_window, surf_id, surf_id_len )
424  integer(kind=kint) :: ctrl
425  integer(kind=kint) :: fg_type
426  integer(kind=kint) :: fg_first
427  integer(kind=kint) :: fg_window
428  character(len=HECMW_NAME_LEN),target :: surf_id(:)
429  character(len=HECMW_NAME_LEN),pointer :: surf_id_p
430  integer(kind=kint) :: surf_id_len
431  integer(kind=kint) :: fstr_ctrl_get_couple
432 
433  character(len=HECMW_NAME_LEN) :: data_fmt,ss
434  write(ss,*) surf_id_len
435  write(data_fmt,'(a,a,a)') 'S',trim(adjustl(ss)),' '
436 
438  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', '1,2,3,4,5,6 ', 0, 'I', fg_type )/= 0) return
439  if( fstr_ctrl_get_param_ex( ctrl, 'ISTEP ', '# ', 0, 'I', fg_first )/= 0) return
440  if( fstr_ctrl_get_param_ex( ctrl, 'WINDOW ', '# ', 0, 'I', fg_window )/= 0) return
441 
442  surf_id_p => surf_id(1)
444  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, surf_id_p )
445 
446  end function fstr_ctrl_get_couple
447 
449  function fstr_ctrl_get_mpc( ctrl, penalty )
450  integer(kind=kint), intent(in) :: ctrl
451  real(kind=kreal), intent(out) :: penalty
452  integer(kind=kint) :: fstr_ctrl_get_mpc
453 
454  fstr_ctrl_get_mpc = fstr_ctrl_get_data_ex( ctrl, 1, 'r ', penalty )
455  if( penalty <= 1.0 ) then
456  if (myrank == 0) then
457  write(imsg,*) "Warging : !MPC : too small penalty: ", penalty
458  write(*,*) "Warging : !MPC : too small penalty: ", penalty
459  endif
460  endif
461 
462  end function fstr_ctrl_get_mpc
463 
465  logical function fstr_ctrl_get_outitem( ctrl, hecMESH, outinfo )
467  use m_out
468  integer(kind=kint), intent(in) :: ctrl
469  type (hecmwst_local_mesh), intent(in) :: hecmesh
470  type( output_info ), intent(out) :: outinfo
471 
472  integer(kind=kint) :: rcode, ipos
473  integer(kind=kint) :: n, i, j
474  character(len=HECMW_NAME_LEN) :: data_fmt, ss
475  character(len=HECMW_NAME_LEN), allocatable :: header_name(:), onoff(:), vtype(:)
476 
477  write( ss, * ) hecmw_name_len
478  write( data_fmt, '(a,a,a,a,a)') 'S', trim(adjustl(ss)), 'S', trim(adjustl(ss)), ' '
479  ! write( data_fmt, '(a,a,a,a,a,a,a)') 'S', trim(adjustl(ss)), 'S', trim(adjustl(ss)), 'S', trim(adjustl(ss)), ' '
480 
481  fstr_ctrl_get_outitem = .false.
482 
483  outinfo%grp_id_name = "ALL"
484  rcode = fstr_ctrl_get_param_ex( ctrl, 'GROUP ', '# ', 0, 'S', outinfo%grp_id_name )
485  ipos = 0
486  rcode = fstr_ctrl_get_param_ex( ctrl, 'ACTION ', 'SUM ', 0, 'P', ipos )
487  outinfo%actn = ipos
488 
489  n = fstr_ctrl_get_data_line_n( ctrl )
490  if( n == 0 ) return
491  allocate( header_name(n), onoff(n), vtype(n) )
492  header_name(:) = ""; vtype(:) = ""; onoff(:) = ""
493  rcode = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, header_name, onoff )
494  ! rcode = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, header_name, onoff, vtype )
495 
496  do i = 1, n
497  do j = 1, outinfo%num_items
498  if( trim(header_name(i)) == outinfo%keyWord(j) ) then
499  outinfo%on(j) = .true.
500  if( trim(onoff(i)) == 'OFF' ) outinfo%on(j) = .false.
501  if( len( trim(vtype(i)) )>0 ) then
502  if( fstr_str2index( vtype(i), ipos ) ) then
503  outinfo%vtype(j) = ipos
504  else if( trim(vtype(i)) == "SCALER" ) then
505  outinfo%vtype(j) = -1
506  else if( trim(vtype(i)) == "VECTOR" ) then
507  outinfo%vtype(j) = -2
508  else if( trim(vtype(i)) == "SYMTENSOR" ) then
509  outinfo%vtype(j) = -3
510  else if( trim(vtype(i)) == "TENSOR" ) then
511  outinfo%vtype(j) = -4
512  endif
513  endif
514  endif
515  enddo
516  enddo
517 
518  deallocate( header_name, onoff, vtype )
519  fstr_ctrl_get_outitem = .true.
520 
521  end function fstr_ctrl_get_outitem
522 
524  function fstr_ctrl_get_contactalgo( ctrl, algo )
525  integer(kind=kint) :: ctrl
526  integer(kind=kint) :: algo
527  integer(kind=kint) :: fstr_ctrl_get_contactalgo
528 
529  integer(kind=kint) :: rcode
530  character(len=80) :: s
531  algo = kcaslagrange
532  s = 'SLAGRANGE,ALAGRANGE '
533  rcode = fstr_ctrl_get_param_ex( ctrl, 'TYPE ', s, 0, 'P', algo )
535  end function fstr_ctrl_get_contactalgo
536 
538  logical function fstr_ctrl_get_contact( ctrl, n, contact, np, tp, ntol, ttol, ctAlgo, cpname )
540  integer(kind=kint), intent(in) :: ctrl
541  integer(kind=kint), intent(in) :: n
542  integer(kind=kint), intent(in) :: ctalgo
543  type(tcontact), intent(out) :: contact(n)
544  real(kind=kreal), intent(out) :: np
545  real(kind=kreal), intent(out) :: tp
546  real(kind=kreal), intent(out) :: ntol
547  real(kind=kreal), intent(out) :: ttol
548  character(len=*), intent(out) :: cpname
549 
550  integer :: rcode, ipt
551  character(len=30) :: s1 = 'TIED,GLUED,SSLID,FSLID '
552  character(len=HECMW_NAME_LEN) :: data_fmt,ss
553  character(len=HECMW_NAME_LEN) :: cp_name(n)
554  real(kind=kreal) :: fcoeff(n),tpenalty(n)
555 
556  tpenalty = 1.0d6
557 
558  write(ss,*) hecmw_name_len
559 
560  fstr_ctrl_get_contact = .false.
561  contact(1)%ctype = 1 ! pure slave-master contact; default value
562  contact(1)%algtype = contactsslid ! small sliding contact; default value
563  rcode = fstr_ctrl_get_param_ex( ctrl, 'INTERACTION ', s1, 0, 'P', contact(1)%algtype )
564  if( contact(1)%algtype==contactglued ) contact(1)%algtype=contactfslid ! not complemented yet
565  if( fstr_ctrl_get_param_ex( ctrl, 'GRPID ', '# ', 1, 'I', contact(1)%group )/=0) return
566  do rcode=2,n
567  contact(rcode)%ctype = contact(1)%ctype
568  contact(rcode)%group = contact(1)%group
569  contact(rcode)%algtype = contact(1)%algtype
570  end do
571 
572  if( contact(1)%algtype==contactsslid .or. contact(1)%algtype==contactfslid ) then
573  write( data_fmt, '(a,a,a)') 'S', trim(adjustl(ss)),'Rr '
574  if( fstr_ctrl_get_data_array_ex( ctrl, data_fmt, cp_name, fcoeff, tpenalty ) /= 0 ) return
575  do rcode=1,n
576  call fstr_strupr(cp_name(rcode))
577  contact(rcode)%pair_name = cp_name(rcode)
578  contact(rcode)%fcoeff = fcoeff(rcode)
579  contact(rcode)%tPenalty = tpenalty(rcode)
580  enddo
581  else if( contact(1)%algtype==contacttied ) then
582  write( data_fmt, '(a,a)') 'S', trim(adjustl(ss))
583  if( fstr_ctrl_get_data_array_ex( ctrl, data_fmt, cp_name ) /= 0 ) return
584  do rcode=1,n
585  call fstr_strupr(cp_name(rcode))
586  contact(rcode)%pair_name = cp_name(rcode)
587  contact(rcode)%fcoeff = 0.d0
588  contact(rcode)%tPenalty = 1.d+4
589  enddo
590  endif
591 
592  np = 0.d0; tp=0.d0
593  ntol = 0.d0; ttol=0.d0
594  if( fstr_ctrl_get_param_ex( ctrl, 'NPENALTY ', '# ', 0, 'R', np ) /= 0 ) return
595  if( fstr_ctrl_get_param_ex( ctrl, 'TPENALTY ', '# ', 0, 'R', tp ) /= 0 ) return
596  if( fstr_ctrl_get_param_ex( ctrl, 'NTOL ', '# ', 0, 'R', ntol ) /= 0 ) return
597  if( fstr_ctrl_get_param_ex( ctrl, 'TTOL ', '# ', 0, 'R', ttol ) /= 0 ) return
598  cpname=""
599  if( fstr_ctrl_get_param_ex( ctrl, 'CONTACTPARAM ', '# ', 0, 'S', cpname )/= 0) return
600  fstr_ctrl_get_contact = .true.
601  end function fstr_ctrl_get_contact
602 
604  logical function fstr_ctrl_get_embed( ctrl, n, embed, cpname )
606  integer(kind=kint), intent(in) :: ctrl
607  integer(kind=kint), intent(in) :: n
608  type(tcontact), intent(out) :: embed(n)
609  character(len=*), intent(out) :: cpname
610 
611  integer :: rcode, ipt
612  character(len=30) :: s1 = 'TIED,GLUED,SSLID,FSLID '
613  character(len=HECMW_NAME_LEN) :: data_fmt,ss
614  character(len=HECMW_NAME_LEN) :: cp_name(n)
615  real(kind=kreal) :: fcoeff(n),tpenalty(n)
616 
617  tpenalty = 1.0d6
618 
619  write(ss,*) hecmw_name_len
620 
621  fstr_ctrl_get_embed = .false.
622  embed(1)%ctype = 1 ! pure slave-master contact; default value
623  embed(1)%algtype = contacttied ! small sliding contact; default value
624  if( fstr_ctrl_get_param_ex( ctrl, 'GRPID ', '# ', 1, 'I', embed(1)%group )/=0) return
625  do rcode=2,n
626  embed(rcode)%ctype = embed(1)%ctype
627  embed(rcode)%group = embed(1)%group
628  embed(rcode)%algtype = embed(1)%algtype
629  end do
630 
631  write( data_fmt, '(a,a)') 'S', trim(adjustl(ss))
632  if( fstr_ctrl_get_data_array_ex( ctrl, data_fmt, cp_name ) /= 0 ) return
633  do rcode=1,n
634  call fstr_strupr(cp_name(rcode))
635  embed(rcode)%pair_name = cp_name(rcode)
636  enddo
637 
638  cpname=""
639  if( fstr_ctrl_get_param_ex( ctrl, 'CONTACTPARAM ', '# ', 0, 'S', cpname )/= 0) return
640  fstr_ctrl_get_embed = .true.
641  end function
642 
644  function fstr_ctrl_get_contactparam( ctrl, contactparam )
645  implicit none
646  integer(kind=kint) :: ctrl
647  type( tcontactparam ) :: contactparam
648  integer(kind=kint) :: fstr_ctrl_get_contactparam
649 
650  integer(kind=kint) :: rcode
651  character(len=HECMW_NAME_LEN) :: data_fmt
652  character(len=128) :: msg
653  real(kind=kreal) :: clearance, clr_same_elem, clr_difflpos, clr_cal_norm
654  real(kind=kreal) :: distclr_init, distclr_free, distclr_nocheck, tensile_force
655  real(kind=kreal) :: box_exp_rate
656 
658 
659  !parameters
660  contactparam%name = ''
661  if( fstr_ctrl_get_param_ex( ctrl, 'NAME ', '# ', 1, 'S', contactparam%name ) /=0 ) return
662 
663  !read first line
664  data_fmt = 'rrrr '
665  rcode = fstr_ctrl_get_data_ex( ctrl, 1, data_fmt, &
666  & clearance, clr_same_elem, clr_difflpos, clr_cal_norm )
667  if( rcode /= 0 ) return
668  contactparam%CLEARANCE = clearance
669  contactparam%CLR_SAME_ELEM = clr_same_elem
670  contactparam%CLR_DIFFLPOS = clr_difflpos
671  contactparam%CLR_CAL_NORM = clr_cal_norm
672 
673  !read second line
674  data_fmt = 'rrrrr '
675  rcode = fstr_ctrl_get_data_ex( ctrl, 2, data_fmt, &
676  & distclr_init, distclr_free, distclr_nocheck, tensile_force, box_exp_rate )
677  if( rcode /= 0 ) return
678  contactparam%DISTCLR_INIT = distclr_init
679  contactparam%DISTCLR_FREE = distclr_free
680  contactparam%DISTCLR_NOCHECK = distclr_nocheck
681  contactparam%TENSILE_FORCE = tensile_force
682  contactparam%BOX_EXP_RATE = box_exp_rate
683 
684  !input check
685  rcode = 1
686  if( clearance<0.d0 .OR. 1.d0<clearance ) THEN
687  write(msg,*) 'fstr control file error : !CONTACT_PARAM : CLEARANCE must be 0 < CLEARANCE < 1.'
688  else if( clr_same_elem<0.d0 .or. 1.d0<clr_same_elem ) then
689  write(msg,*) 'fstr control file error : !CONTACT_PARAM : CLR_SAME_ELEM must be 0 < CLR_SAME_ELEM < 1.'
690  else if( clr_difflpos<0.d0 .or. 1.d0<clr_difflpos ) then
691  write(msg,*) 'fstr control file error : !CONTACT_PARAM : CLR_DIFFLPOS must be 0 < CLR_DIFFLPOS < 1.'
692  else if( clr_cal_norm<0.d0 .or. 1.d0<clr_cal_norm ) then
693  write(msg,*) 'fstr control file error : !CONTACT_PARAM : CLR_CAL_NORM must be 0 < CLR_CAL_NORM < 1.'
694  else if( distclr_init<0.d0 .or. 1.d0<distclr_init ) then
695  write(msg,*) 'fstr control file error : !CONTACT_PARAM : DISTCLR_INIT must be 0 < DISTCLR_INIT < 1.'
696  else if( distclr_free<-1.d0 .or. 1.d0<distclr_free ) then
697  write(msg,*) 'fstr control file error : !CONTACT_PARAM : DISTCLR_FREE must be -1 < DISTCLR_FREE < 1.'
698  else if( distclr_nocheck<0.5d0 ) then
699  write(msg,*) 'fstr control file error : !CONTACT_PARAM : DISTCLR_NOCHECK must be >= 0.5.'
700  else if( tensile_force>=0.d0 ) then
701  write(msg,*) 'fstr control file error : !CONTACT_PARAM : TENSILE_FORCE must be < 0.'
702  else if( box_exp_rate<=1.d0 .or. 2.0<box_exp_rate ) then
703  write(msg,*) 'fstr control file error : !CONTACT_PARAM : BOX_EXP_RATE must be 1 < BOX_EXP_RATE <= 2.'
704  else
705  rcode =0
706  end if
707  if( rcode /= 0 ) then
708  write(*,*) trim(msg)
709  write(ilog,*) trim(msg)
710  return
711  endif
712 
714  end function fstr_ctrl_get_contactparam
715 
717  function fstr_ctrl_get_elemopt( ctrl, elemopt361 )
718  integer(kind=kint) :: ctrl
719  integer(kind=kint) :: elemopt361
720  integer(kind=kint) :: fstr_ctrl_get_elemopt
721 
722  character(72) :: o361list = 'IC,Bbar '
723 
724  integer(kind=kint) :: o361
725 
727 
728  o361 = elemopt361 + 1
729 
730  !* parameter in header line -----------------------------------------------------------------*!
731  if( fstr_ctrl_get_param_ex( ctrl, '361 ', o361list, 0, 'P', o361 ) /= 0) return
732 
733  elemopt361 = o361 - 1
734 
736 
737  end function fstr_ctrl_get_elemopt
738 
739 
741  function fstr_get_autoinc( ctrl, aincparam )
742  implicit none
743  integer(kind=kint) :: ctrl
744  type( tparamautoinc ) :: aincparam
745  integer(kind=kint) :: fstr_get_autoinc
746 
747  integer(kind=kint) :: rcode
748  character(len=HECMW_NAME_LEN) :: data_fmt
749  character(len=128) :: msg
750  integer(kind=kint) :: bound_s(10), bound_l(10)
751  real(kind=kreal) :: rs, rl
752 
753  fstr_get_autoinc = -1
754 
755  bound_s(:) = 0
756  bound_l(:) = 0
757 
758  !parameters
759  aincparam%name = ''
760  if( fstr_ctrl_get_param_ex( ctrl, 'NAME ', '# ', 1, 'S', aincparam%name ) /=0 ) return
761 
762  !read first line ( decrease criteria )
763  data_fmt = 'riiii '
764  rcode = fstr_ctrl_get_data_ex( ctrl, 1, data_fmt, rs, &
765  & bound_s(1), bound_s(2), bound_s(3), aincparam%NRtimes_s )
766  if( rcode /= 0 ) return
767  aincparam%ainc_Rs = rs
768  aincparam%NRbound_s(knstmaxit) = bound_s(1)
769  aincparam%NRbound_s(knstsumit) = bound_s(2)
770  aincparam%NRbound_s(knstciter) = bound_s(3)
771 
772  !read second line ( increase criteria )
773  data_fmt = 'riiii '
774  rcode = fstr_ctrl_get_data_ex( ctrl, 2, data_fmt, rl, &
775  & bound_l(1), bound_l(2), bound_l(3), aincparam%NRtimes_l )
776  if( rcode /= 0 ) return
777  aincparam%ainc_Rl = rl
778  aincparam%NRbound_l(knstmaxit) = bound_l(1)
779  aincparam%NRbound_l(knstsumit) = bound_l(2)
780  aincparam%NRbound_l(knstciter) = bound_l(3)
781 
782  !read third line ( cutback criteria )
783  data_fmt = 'ri '
784  rcode = fstr_ctrl_get_data_ex( ctrl, 3, data_fmt, &
785  & aincparam%ainc_Rc, aincparam%CBbound )
786  if( rcode /= 0 ) return
787 
788  !input check
789  rcode = 1
790  if( rs<0.d0 .or. rs>1.d0 ) then
791  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : decrease ratio Rs must 0 < Rs < 1.'
792  else if( any(bound_s<0) ) then
793  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : decrease NR bound must >= 0.'
794  else if( aincparam%NRtimes_s < 1 ) then
795  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : # of times to decrease must > 0.'
796  else if( rl<1.d0 ) then
797  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : increase ratio Rl must > 1.'
798  else if( any(bound_l<0) ) then
799  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : increase NR bound must >= 0.'
800  else if( aincparam%NRtimes_l < 1 ) then
801  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : # of times to increase must > 0.'
802  elseif( aincparam%ainc_Rc<0.d0 .or. aincparam%ainc_Rc>1.d0 ) then
803  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : cutback decrease ratio Rc must 0 < Rc < 1.'
804  else if( aincparam%CBbound < 1 ) then
805  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : maximum # of cutback times must > 0.'
806  else
807  rcode =0
808  end if
809  if( rcode /= 0 ) then
810  write(*,*) trim(msg)
811  write(ilog,*) trim(msg)
812  return
813  endif
814 
815  fstr_get_autoinc = 0
816  end function fstr_get_autoinc
817 
819  function fstr_ctrl_get_timepoints( ctrl, tp )
820  integer(kind=kint) :: ctrl
821  type(time_points) :: tp
822  integer(kind=kint) :: fstr_ctrl_get_timepoints
823 
824  integer(kind=kint) :: i, n, rcode
825  logical :: generate
826  real(kind=kreal) :: stime, etime, interval
827 
829 
830  tp%name = ''
831  if( fstr_ctrl_get_param_ex( ctrl, 'NAME ', '# ', 1, 'S', tp%name ) /=0 ) return
832  tp%range_type = 1
833  if( fstr_ctrl_get_param_ex( ctrl, 'TIME ', 'STEP,TOTAL ', 0, 'P', tp%range_type ) /= 0 ) return
834  generate = .false.
835  if( fstr_ctrl_get_param_ex( ctrl, 'GENERATE ', '# ', 0, 'E', generate ) /= 0) return
836 
837  if( generate ) then
838  stime = 0.d0; etime = 0.d0; interval = 1.d0
839  if( fstr_ctrl_get_data_ex( ctrl, 1, 'rrr ', stime, etime, interval ) /= 0) return
840  tp%n_points = int((etime-stime)/interval)+1
841  allocate(tp%points(tp%n_points))
842  do i=1,tp%n_points
843  tp%points(i) = stime + dble(i-1)*interval
844  end do
845  else
846  n = fstr_ctrl_get_data_line_n( ctrl )
847  if( n == 0 ) return
848  tp%n_points = n
849  allocate(tp%points(tp%n_points))
850  if( fstr_ctrl_get_data_array_ex( ctrl, 'r ', tp%points ) /= 0 ) return
851  do i=1,tp%n_points-1
852  if( tp%points(i) < tp%points(i+1) ) cycle
853  write(*,*) 'Error in reading !TIME_POINT: time points must be given in ascending order.'
854  return
855  end do
856  end if
857 
859  end function fstr_ctrl_get_timepoints
860 
862  function fstr_ctrl_get_amplitude( ctrl, nline, name, type_def, type_time, type_val, n, val, table )
863  implicit none
864  integer(kind=kint), intent(in) :: ctrl
865  integer(kind=kint), intent(in) :: nline
866  character(len=HECMW_NAME_LEN), intent(out) :: name
867  integer(kind=kint), intent(out) :: type_def
868  integer(kind=kint), intent(out) :: type_time
869  integer(kind=kint), intent(out) :: type_val
870  integer(kind=kint), intent(out) :: n
871  real(kind=kreal), pointer :: val(:)
872  real(kind=kreal), pointer :: table(:)
873  integer(kind=kint) :: fstr_ctrl_get_amplitude
874 
875  integer(kind=kint) :: t_def, t_time, t_val
876  integer(kind=kint) :: i, j
877  real(kind=kreal) :: r(4), t(4)
878 
880 
881  name = ''
882  t_def = 1
883  t_time = 1
884  t_val = 1
885 
886  if( fstr_ctrl_get_param_ex( ctrl, 'NAME ', '# ', 1, 'S', name )/=0 ) return
887  if( fstr_ctrl_get_param_ex( ctrl, 'DEFINITION ', 'TABULAR ', 0, 'P', t_def )/=0 ) return
888  if( fstr_ctrl_get_param_ex( ctrl, 'TIME ', 'STEP ', 0, 'P', t_time )/=0 ) return
889  if( fstr_ctrl_get_param_ex( ctrl, 'VALUE ', 'RELATIVE,ABSOLUTE ', 0, 'P', t_val )/=0 ) return
890 
891  if( t_def==1 ) then
892  type_def = hecmw_amp_typedef_tabular
893  else
894  write(*,*) 'Error in reading !AMPLITUDE: invalid value for parameter DEFINITION.'
895  endif
896  if( t_time==1 ) then
897  type_time = hecmw_amp_typetime_step
898  else
899  write(*,*) 'Error in reading !AMPLITUDE: invalid value for parameter TIME.'
900  endif
901  if( t_val==1 ) then
902  type_val = hecmw_amp_typeval_relative
903  elseif( t_val==2 ) then
904  type_val = hecmw_amp_typeval_absolute
905  else
906  write(*,*) 'Error in reading !AMPLITUDE: invalid value for parameter VALUE.'
907  endif
908 
909  n = 0
910  do i = 1, nline
911  r(:)=huge(0.0d0); t(:)=huge(0.0d0)
912  if( fstr_ctrl_get_data_ex( ctrl, 1, 'RRrrrrrr ', r(1), t(1), r(2), t(2), r(3), t(3), r(4), t(4) ) /= 0) return
913  n = n+1
914  val(n) = r(1)
915  table(n) = t(1)
916  do j = 2, 4
917  if (r(j) < huge(0.0d0) .and. t(j) < huge(0.0d0)) then
918  n = n+1
919  val(n) = r(j)
920  table(n) = t(j)
921  else
922  exit
923  endif
924  enddo
925  enddo
927 
928  end function fstr_ctrl_get_amplitude
929 
930 end module fstr_ctrl_common
fstr_ctrl_common::fstr_ctrl_get_contactparam
integer(kind=kint) function fstr_ctrl_get_contactparam(ctrl, contactparam)
Read in !CONTACT_PARAM !
Definition: fstr_ctrl_common.f90:645
m_out
This module manages step information.
Definition: m_out.f90:6
fstr_ctrl_common::fstr_ctrl_get_amplitude
integer(kind=kint) function fstr_ctrl_get_amplitude(ctrl, nline, name, type_def, type_time, type_val, n, val, table)
Read in !AMPLITUDE.
Definition: fstr_ctrl_common.f90:863
m_out::output_info
output information
Definition: m_out.f90:17
m_fstr::kon
integer(kind=kint), parameter kon
Definition: m_fstr.f90:32
fstr_ctrl_common::fstr_ctrl_get_solver
integer(kind=kint) function fstr_ctrl_get_solver(ctrl, method, precond, nset, iterlog, timelog, steplog, nier, iterpremax, nrest, nBFGS, scaling, dumptype, dumpexit, usejad, ncolor_in, mpc_method, estcond, method2, recyclepre, solver_opt, contact_elim, resid, singma_diag, sigma, thresh, filter)
Read in !SOLVER.
Definition: fstr_ctrl_common.f90:89
fstr_ctrl_common
This module contains fstr control file data obtaining functions.
Definition: fstr_ctrl_common.f90:7
m_step::stepautoinc
integer, parameter stepautoinc
Definition: m_step.f90:15
fstr_ctrl_common::fstr_ctrl_get_mpc
integer(kind=kint) function fstr_ctrl_get_mpc(ctrl, penalty)
Read in !MPC.
Definition: fstr_ctrl_common.f90:450
m_fstr::myrank
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.f90:96
m_fstr::tsection
Data for section control.
Definition: m_fstr.f90:641
m_timepoint::time_points
Time points storage for output etc.
Definition: m_timepoint.f90:14
m_fstr::kststaticeigen
integer(kind=kint), parameter kststaticeigen
Definition: m_fstr.f90:42
fstr_setup_util::fstr_strupr
subroutine fstr_strupr(s)
Definition: fstr_setup_util.f90:58
m_step::stepfixedinc
integer, parameter stepfixedinc
Definition: m_step.f90:14
fstr_ctrl_common::fstr_ctrl_get_embed
logical function fstr_ctrl_get_embed(ctrl, n, embed, cpname)
Read in contact definition.
Definition: fstr_ctrl_common.f90:605
fstr_ctrl_common::fstr_ctrl_get_solution
integer(kind=kint) function fstr_ctrl_get_solution(ctrl, type, nlgeom)
Read in !SOLUTION.
Definition: fstr_ctrl_common.f90:40
m_fstr::kcaslagrange
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
Definition: m_fstr.f90:59
fstr_ctrl_common::fstr_ctrl_get_section
integer function fstr_ctrl_get_section(ctrl, hecMESH, sections)
Read in !SECTION.
Definition: fstr_ctrl_common.f90:336
m_fstr::kststatic
integer(kind=kint), parameter kststatic
Definition: m_fstr.f90:37
fstr_ctrl_common::fstr_ctrl_get_contactalgo
integer(kind=kint) function fstr_ctrl_get_contactalgo(ctrl, algo)
Read in !CONTACT.
Definition: fstr_ctrl_common.f90:525
m_timepoint
This module manages timepoint information.
Definition: m_timepoint.f90:6
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
fstr_ctrl_get_data_line_n
int fstr_ctrl_get_data_line_n(int *ctrl)
Definition: fstr_ctrl_util.c:1462
fstr_setup_util::amp_name_to_id
subroutine amp_name_to_id(hecMESH, header_name, aname, id)
Definition: fstr_setup_util.f90:1078
fstr_ctrl_common::fstr_ctrl_get_contact
logical function fstr_ctrl_get_contact(ctrl, n, contact, np, tp, ntol, ttol, ctAlgo, cpname)
Read in contact definition.
Definition: fstr_ctrl_common.f90:539
fstr_ctrl_common::fstr_ctrl_get_outitem
logical function fstr_ctrl_get_outitem(ctrl, hecMESH, outinfo)
Read in !OUTPUT_RES & !OUTPUT_VIS.
Definition: fstr_ctrl_common.f90:466
m_step::stepstatic
integer, parameter stepstatic
Definition: m_step.f90:12
fstr_ctrl_common::fstr_ctrl_get_timepoints
integer(kind=kint) function fstr_ctrl_get_timepoints(ctrl, tp)
Read in !TIME_POINTS.
Definition: fstr_ctrl_common.f90:820
fstr_ctrl_common::fstr_ctrl_get_step
integer(kind=kint) function fstr_ctrl_get_step(ctrl, amp, iproc)
Read in !STEP.
Definition: fstr_ctrl_common.f90:204
fstr_ctrl_common::fstr_ctrl_get_couple
integer(kind=kint) function fstr_ctrl_get_couple(ctrl, fg_type, fg_first, fg_window, surf_id, surf_id_len)
Read in !COUPLE.
Definition: fstr_ctrl_common.f90:424
hecmw
Definition: hecmw.f90:6
m_step::step_info
Step control such as active boundary condition, convergent condition etc.
Definition: m_step.f90:24
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_common::fstr_ctrl_get_echo
integer(kind=kint) function fstr_ctrl_get_echo(ctrl, echo)
Read in !ECHO.
Definition: fstr_ctrl_common.f90:412
fstr_ctrl_common::fstr_ctrl_get_istep
logical function fstr_ctrl_get_istep(ctrl, hecMESH, steps, tpname, apname)
Read in !STEP and !ISTEP.
Definition: fstr_ctrl_common.f90:226
fstr_setup_util::fstr_str2index
logical function fstr_str2index(s, x)
Definition: fstr_setup_util.f90:33
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_common::fstr_ctrl_get_elemopt
integer(kind=kint) function fstr_ctrl_get_elemopt(ctrl, elemopt361)
Read in !ELEMOPT.
Definition: fstr_ctrl_common.f90:718
m_step::init_stepinfo
subroutine init_stepinfo(stepinfo)
Initializer.
Definition: m_step.f90:67
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
fstr_ctrl_common::fstr_get_autoinc
integer(kind=kint) function fstr_get_autoinc(ctrl, aincparam)
Read in !AUTOINC_PARAM !
Definition: fstr_ctrl_common.f90:742
m_fstr::etime
real(kind=kreal) etime
Definition: m_fstr.f90:140
fstr_ctrl_common::fstr_ctrl_get_nonlinear_solver
integer(kind=kint) function fstr_ctrl_get_nonlinear_solver(ctrl, method)
Read in !NONLINEAR_SOLVER.
Definition: fstr_ctrl_common.f90:69
m_fstr::ilog
integer(kind=kint), parameter ilog
FILE HANDLER.
Definition: m_fstr.f90:107
fstr_ctrl_common::fstr_ctrl_get_write
integer(kind=kint) function fstr_ctrl_get_write(ctrl, res, visual, femap)
Read in !WRITE.
Definition: fstr_ctrl_common.f90:393
m_step
This module manages step information.
Definition: m_step.f90:6
fstr_setup_util
This module contains auxiliary functions in calculation setup.
Definition: fstr_setup_util.f90:7
m_fstr::kstdynamic
integer(kind=kint), parameter kstdynamic
Definition: m_fstr.f90:40
mcontact
This module provides functions to calculate contact stiff matrix.
Definition: fstr_contact.f90:6
m_fstr::imsg
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:110