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