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(inout) :: 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  steps%solution = stepstatic
248  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', 'STATIC,VISCO ', 0, 'P', steps%solution )/= 0) return
249  steps%inc_type = stepfixedinc
250  if( fstr_ctrl_get_param_ex( ctrl, 'INC_TYPE ', 'FIXED,AUTO ', 0, 'P', steps%inc_type )/= 0) return
251  if( fstr_ctrl_get_param_ex( ctrl, 'SUBSTEPS ', '# ', 0, 'I', steps%num_substep )/= 0) return
252  steps%initdt = 1.d0/steps%num_substep
253  if( fstr_ctrl_get_param_ex( ctrl, 'ITMAX ', '# ', 0, 'I', steps%max_iter )/= 0) return
254  if( fstr_ctrl_get_param_ex( ctrl, 'MAXITER ', '# ', 0, 'I', steps%max_iter )/= 0) return
255  if( fstr_ctrl_get_param_ex( ctrl, 'MAXCONTITER ', '# ', 0, 'I', steps%max_contiter )/= 0) return
256  if( fstr_ctrl_get_param_ex( ctrl, 'CONVERG ', '# ', 0, 'R', steps%converg )/= 0) return
257  if( fstr_ctrl_get_param_ex( ctrl, 'CONVERG_LAG ', '# ', 0, 'R', steps%converg_lag )/= 0) return
258  if( fstr_ctrl_get_param_ex( ctrl, 'CONVERG_DDISP ', '# ', 0, 'R', steps%converg_ddisp )/= 0) return
259  if( fstr_ctrl_get_param_ex( ctrl, 'MAXRES ', '# ', 0, 'R', steps%maxres )/= 0) return
260  amp = ""
261  if( fstr_ctrl_get_param_ex( ctrl, 'AMP ', '# ', 0, 'S', amp )/= 0) return
262  if( len( trim(amp) )>0 ) then
263  call amp_name_to_id( hecmesh, '!STEP', amp, steps%amp_id )
264  endif
265  tpname=""
266  if( fstr_ctrl_get_param_ex( ctrl, 'TIMEPOINTS ', '# ', 0, 'S', tpname )/= 0) return
267  apname=""
268  if( fstr_ctrl_get_param_ex( ctrl, 'AUTOINCPARAM ', '# ', 0, 'S', apname )/= 0) return
269 
270  n = fstr_ctrl_get_data_line_n( ctrl )
271  if( n == 0 ) then
272  fstr_ctrl_get_istep = .true.; return
273  endif
274 
275  f2 = steps%mindt
276  f3 = steps%maxdt
277  if( fstr_ctrl_get_data_ex( ctrl, 1, data_fmt1, ss, f1, f2, f3 )/= 0) return
278  read( ss, * , iostat=ierr ) fn
279  sn=1
280  if( ierr==0 ) then
281  steps%initdt = fn
282  steps%elapsetime = f1
283  if( steps%inc_type == stepautoinc ) then
284  steps%mindt = min(f2,steps%initdt)
285  steps%maxdt = f3
286  endif
287  steps%num_substep = max(int((f1+0.999999999d0*fn)/fn),steps%num_substep)
288  !if( mod(f1,fn)/=0 ) steps%num_substep =steps%num_substep+1
289  sn = 2
290  endif
291 
292  bc_n = 0
293  load_n = 0
294  contact_n = 0
295  elemact_n = 0
296  do i=sn,n
297  if( fstr_ctrl_get_data_ex( ctrl, i, data_fmt, header_name, bcid )/= 0) return
298  if( trim(header_name) == 'BOUNDARY' ) then
299  bc_n = bc_n + 1
300  else if( trim(header_name) == 'LOAD' ) then
301  load_n = load_n +1
302  else if( trim(header_name) == 'CONTACT' ) then
303  contact_n = contact_n+1
304  else if( trim(header_name) == 'ELEMACT' ) then
305  elemact_n = elemact_n+1
306  else if( trim(header_name) == 'TEMPERATURE' ) then
307  ! steps%Temperature = .true.
308  endif
309  end do
310 
311  if( bc_n>0 ) allocate( steps%Boundary(bc_n) )
312  if( load_n>0 ) allocate( steps%Load(load_n) )
313  if( contact_n>0 ) allocate( steps%Contact(contact_n) )
314  if( elemact_n>0 ) allocate( steps%ElemActivation(elemact_n) )
315 
316  bc_n = 0
317  load_n = 0
318  contact_n = 0
319  elemact_n = 0
320  do i=sn,n
321  if( fstr_ctrl_get_data_ex( ctrl, i, data_fmt, header_name, bcid )/= 0) return
322  if( trim(header_name) == 'BOUNDARY' ) then
323  bc_n = bc_n + 1
324  steps%Boundary(bc_n) = bcid
325  else if( trim(header_name) == 'LOAD' ) then
326  load_n = load_n +1
327  steps%Load(load_n) = bcid
328  else if( trim(header_name) == 'CONTACT' ) then
329  contact_n = contact_n+1
330  steps%Contact(contact_n) = bcid
331  else if( trim(header_name) == 'ELEMACT' ) then
332  elemact_n = elemact_n+1
333  steps%ElemActivation(elemact_n) = bcid
334  endif
335  end do
336 
337  fstr_ctrl_get_istep = .true.
338  end function fstr_ctrl_get_istep
339 
341  integer function fstr_ctrl_get_section( ctrl, hecMESH, sections )
342  use fstr_setup_util
343  integer(kind=kint), intent(in) :: ctrl
344  type (hecmwst_local_mesh), intent(inout) :: hecmesh
345  type (tsection), pointer, intent(inout) :: sections(:)
346 
347  integer(kind=kint) :: j, k, sect_id, ori_id, elemopt
348  integer(kind=kint),save :: cache = 1
349  character(len=HECMW_NAME_LEN) :: sect_orien
350  character(19) :: form341list = 'FI,SELECTIVE_ESNS '
351  character(16) :: form361list = 'FI,BBAR,IC,FBAR '
352 
354 
355  if( fstr_ctrl_get_param_ex( ctrl, 'SECNUM ', '# ', 1, 'I', sect_id )/= 0) return
356  if( sect_id > hecmesh%section%n_sect ) return
357 
358  elemopt = 0
359  if( fstr_ctrl_get_param_ex( ctrl, 'FORM341 ', form341list, 0, 'P', elemopt )/= 0) return
360  if( elemopt > 0 ) sections(sect_id)%elemopt341 = elemopt
361 
362  elemopt = 0
363  if( fstr_ctrl_get_param_ex( ctrl, 'FORM361 ', form361list, 0, 'P', elemopt )/= 0) return
364  if( elemopt > 0 ) sections(sect_id)%elemopt361 = elemopt
365 
366  ! sectional orientation ID
367  hecmesh%section%sect_orien_ID(sect_id) = -1
368  if( fstr_ctrl_get_param_ex( ctrl, 'ORIENTATION ', '# ', 0, 'S', sect_orien )/= 0) return
369 
370  if( associated(g_localcoordsys) ) then
371  call fstr_strupr(sect_orien)
372  k = size(g_localcoordsys)
373 
374  if(cache < k)then
375  if( sect_orien == g_localcoordsys(cache)%sys_name ) then
376  hecmesh%section%sect_orien_ID(sect_id) = cache
377  cache = cache + 1
379  return
380  endif
381  endif
382 
383  do j=1, k
384  if( sect_orien == g_localcoordsys(j)%sys_name ) then
385  hecmesh%section%sect_orien_ID(sect_id) = j
386  cache = j + 1
387  exit
388  endif
389  enddo
390  endif
391 
393 
394  end function fstr_ctrl_get_section
395 
396 
398  function fstr_ctrl_get_write( ctrl, res, visual, femap )
399  integer(kind=kint) :: ctrl
400  integer(kind=kint) :: res
401  integer(kind=kint) :: visual
402  integer(kind=kint) :: femap
403  integer(kind=kint) :: fstr_ctrl_get_write
404 
406 
407  ! JP-6
408  if( fstr_ctrl_get_param_ex( ctrl, 'RESULT ', '# ', 0, 'E', res )/= 0) return
409  if( fstr_ctrl_get_param_ex( ctrl, 'VISUAL ', '# ', 0, 'E', visual )/= 0) return
410  if( fstr_ctrl_get_param_ex( ctrl, 'FEMAP ', '# ', 0, 'E', femap )/= 0) return
411 
413 
414  end function fstr_ctrl_get_write
415 
417  function fstr_ctrl_get_echo( ctrl, echo )
418  integer(kind=kint) :: ctrl
419  integer(kind=kint) :: echo
420  integer(kind=kint) :: fstr_ctrl_get_echo
421 
422  echo = kon;
423 
425 
426  end function fstr_ctrl_get_echo
427 
429  function fstr_ctrl_get_couple( ctrl, fg_type, fg_first, fg_window, surf_id, surf_id_len )
430  integer(kind=kint) :: ctrl
431  integer(kind=kint) :: fg_type
432  integer(kind=kint) :: fg_first
433  integer(kind=kint) :: fg_window
434  character(len=HECMW_NAME_LEN) :: surf_id(:)
435  integer(kind=kint) :: surf_id_len
436  integer(kind=kint) :: fstr_ctrl_get_couple
437 
438  character(len=HECMW_NAME_LEN) :: data_fmt,ss
439  write(ss,*) surf_id_len
440  write(data_fmt,'(a,a,a)') 'S',trim(adjustl(ss)),' '
441 
443  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', '1,2,3,4,5,6 ', 0, 'I', fg_type )/= 0) return
444  if( fstr_ctrl_get_param_ex( ctrl, 'ISTEP ', '# ', 0, 'I', fg_first )/= 0) return
445  if( fstr_ctrl_get_param_ex( ctrl, 'WINDOW ', '# ', 0, 'I', fg_window )/= 0) return
446 
448  fstr_ctrl_get_data_array_ex( ctrl, data_fmt, surf_id )
449 
450  end function fstr_ctrl_get_couple
451 
453  function fstr_ctrl_get_mpc( ctrl, penalty )
454  integer(kind=kint), intent(in) :: ctrl
455  real(kind=kreal), intent(out) :: penalty
456  integer(kind=kint) :: fstr_ctrl_get_mpc
457 
458  fstr_ctrl_get_mpc = fstr_ctrl_get_data_ex( ctrl, 1, 'r ', penalty )
459  if( penalty <= 1.0 ) then
460  if (myrank == 0) then
461  write(imsg,*) "Warging : !MPC : too small penalty: ", penalty
462  write(*,*) "Warging : !MPC : too small penalty: ", penalty
463  endif
464  endif
465 
466  end function fstr_ctrl_get_mpc
467 
469  logical function fstr_ctrl_get_outitem( ctrl, hecMESH, outinfo )
470  use fstr_setup_util
471  use m_out
472  integer(kind=kint), intent(in) :: ctrl
473  type (hecmwst_local_mesh), intent(in) :: hecmesh
474  type( output_info ), intent(out) :: outinfo
475 
476  integer(kind=kint) :: rcode, ipos
477  integer(kind=kint) :: n, i, j
478  character(len=HECMW_NAME_LEN) :: data_fmt, ss
479  character(len=HECMW_NAME_LEN), allocatable :: header_name(:), onoff(:), vtype(:)
480 
481  write( ss, * ) hecmw_name_len
482  write( data_fmt, '(a,a,a,a,a)') 'S', trim(adjustl(ss)), 'S', trim(adjustl(ss)), ' '
483  ! write( data_fmt, '(a,a,a,a,a,a,a)') 'S', trim(adjustl(ss)), 'S', trim(adjustl(ss)), 'S', trim(adjustl(ss)), ' '
484 
485  fstr_ctrl_get_outitem = .false.
486 
487  outinfo%grp_id_name = "ALL"
488  rcode = fstr_ctrl_get_param_ex( ctrl, 'GROUP ', '# ', 0, 'S', outinfo%grp_id_name )
489  ipos = 0
490  rcode = fstr_ctrl_get_param_ex( ctrl, 'ACTION ', 'SUM ', 0, 'P', ipos )
491  outinfo%actn = ipos
492 
493  n = fstr_ctrl_get_data_line_n( ctrl )
494  if( n == 0 ) return
495  allocate( header_name(n), onoff(n), vtype(n) )
496  header_name(:) = ""; vtype(:) = ""; onoff(:) = ""
497  rcode = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, header_name, onoff )
498  ! rcode = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, header_name, onoff, vtype )
499 
500  do i = 1, n
501  do j = 1, outinfo%num_items
502  if( trim(header_name(i)) == outinfo%keyWord(j) ) then
503  outinfo%on(j) = .true.
504  if( trim(onoff(i)) == 'OFF' ) outinfo%on(j) = .false.
505  if( len( trim(vtype(i)) )>0 ) then
506  if( fstr_str2index( vtype(i), ipos ) ) then
507  outinfo%vtype(j) = ipos
508  else if( trim(vtype(i)) == "SCALER" ) then
509  outinfo%vtype(j) = -1
510  else if( trim(vtype(i)) == "VECTOR" ) then
511  outinfo%vtype(j) = -2
512  else if( trim(vtype(i)) == "SYMTENSOR" ) then
513  outinfo%vtype(j) = -3
514  else if( trim(vtype(i)) == "TENSOR" ) then
515  outinfo%vtype(j) = -4
516  endif
517  endif
518  endif
519  enddo
520  enddo
521 
522  deallocate( header_name, onoff, vtype )
523  fstr_ctrl_get_outitem = .true.
524 
525  end function fstr_ctrl_get_outitem
526 
528  function fstr_ctrl_get_contactalgo( ctrl, algo, augiter )
529  integer(kind=kint) :: ctrl
530  integer(kind=kint) :: algo
531  integer(kind=kint) :: augiter
532  integer(kind=kint) :: fstr_ctrl_get_contactalgo
533 
534  integer(kind=kint) :: rcode
535  character(len=80) :: s
536  s = 'SLAGRANGE,ALAGRANGE '
537  rcode = fstr_ctrl_get_param_ex( ctrl, 'TYPE ', s, 0, 'P', algo )
538  if( rcode /= 0 ) then
540  return
541  endif
542  rcode = fstr_ctrl_get_param_ex( ctrl, 'AUGITER ', '# ', 0, 'I', augiter )
544  end function fstr_ctrl_get_contactalgo
545 
547  logical function fstr_ctrl_get_contact( ctrl, n, contact, np, tp, ntol, ttol, ctAlgo, cpname, smoothing )
548  use fstr_setup_util
549  integer(kind=kint), intent(in) :: ctrl
550  integer(kind=kint), intent(in) :: n
551  integer(kind=kint), intent(in) :: ctalgo
552  type(tcontact), intent(out) :: contact(n)
553  real(kind=kreal), intent(out) :: np
554  real(kind=kreal), intent(out) :: tp
555  real(kind=kreal), intent(out) :: ntol
556  real(kind=kreal), intent(out) :: ttol
557  character(len=*), intent(out) :: cpname
558  integer(kind=kint), intent(out) :: smoothing
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  real(kind=kreal) :: damp_alpha, damp_gact
566 
567  write(ss,*) hecmw_name_len
568 
569  fstr_ctrl_get_contact = .false.
570  contact(1)%ctype = 1 ! pure slave-master contact; default value
571  contact(1)%algtype = contactsslid ! small sliding contact; default value
572  rcode = fstr_ctrl_get_param_ex( ctrl, 'INTERACTION ', s1, 0, 'P', contact(1)%algtype )
573  if( contact(1)%algtype==contactglued ) contact(1)%algtype=contactfslid ! not complemented yet
574  if( fstr_ctrl_get_param_ex( ctrl, 'GRPID ', '# ', 1, 'I', contact(1)%group )/=0) return
575  smoothing = kcsnone + 1
576  if( fstr_ctrl_get_param_ex( ctrl, 'SMOOTHING ', 'NONE,NAGATA ', 0, 'P', smoothing ) /= 0 ) return
577  smoothing = smoothing - 1
578  do rcode=2,n
579  contact(rcode)%ctype = contact(1)%ctype
580  contact(rcode)%group = contact(1)%group
581  contact(rcode)%algtype = contact(1)%algtype
582  end do
583 
584  tpenalty = 0.5d0
585 
586  if( contact(1)%algtype==contactsslid .or. contact(1)%algtype==contactfslid ) then
587  write( data_fmt, '(a,a,a)') 'S', trim(adjustl(ss)),'Rr '
588  if( fstr_ctrl_get_data_array_ex( ctrl, data_fmt, cp_name, fcoeff, tpenalty ) /= 0 ) return
589  do rcode=1,n
590  call fstr_strupr(cp_name(rcode))
591  contact(rcode)%pair_name = cp_name(rcode)
592  contact(rcode)%fcoeff = fcoeff(rcode)
593  contact(rcode)%nPenalty = 5.0d0
594  contact(rcode)%tPenalty = tpenalty(rcode)
595  contact(rcode)%refStiff = 1.d0
596  contact(rcode)%damp_alpha = 0.0d0
597  contact(rcode)%damp_gact = 0.0d0
598  enddo
599  else if( contact(1)%algtype==contacttied ) then
600  write( data_fmt, '(a,a)') 'S', trim(adjustl(ss))
601  if( fstr_ctrl_get_data_array_ex( ctrl, data_fmt, cp_name ) /= 0 ) return
602  do rcode=1,n
603  call fstr_strupr(cp_name(rcode))
604  contact(rcode)%pair_name = cp_name(rcode)
605  contact(rcode)%nPenalty = 5.0d0
606  contact(rcode)%fcoeff = 0.d0
607  contact(rcode)%tPenalty = 1.d0
608  contact(rcode)%damp_alpha = 0.0d0
609  contact(rcode)%damp_gact = 0.0d0
610  enddo
611  endif
612 
613  np = 0.d0; tp=0.d0
614  ntol = 0.d0; ttol=0.d0
615  damp_alpha = 0.0d0; damp_gact = 0.0d0
616  if( fstr_ctrl_get_param_ex( ctrl, 'NPENALTY ', '# ', 0, 'R', np ) /= 0 ) return
617  if( fstr_ctrl_get_param_ex( ctrl, 'TPENALTY ', '# ', 0, 'R', tp ) /= 0 ) return
618  if( fstr_ctrl_get_param_ex( ctrl, 'NTOL ', '# ', 0, 'R', ntol ) /= 0 ) return
619  if( fstr_ctrl_get_param_ex( ctrl, 'TTOL ', '# ', 0, 'R', ttol ) /= 0 ) return
620  if( fstr_ctrl_get_param_ex( ctrl, 'DAMP_ALPHA ', '# ', 0, 'R', damp_alpha ) /= 0 ) return
621  if( fstr_ctrl_get_param_ex( ctrl, 'DAMP_GACT ', '# ', 0, 'R', damp_gact ) /= 0 ) return
622  cpname=""
623  if( fstr_ctrl_get_param_ex( ctrl, 'CONTACTPARAM ', '# ', 0, 'S', cpname )/= 0) return
624 
625  ! Set penalty coefficients to contact structure if specified (for ALagrange method)
626  if( np > 0.d0 ) then
627  do rcode=1,n
628  contact(rcode)%nPenalty = np
629  enddo
630  endif
631  if( tp > 0.d0 ) then
632  do rcode=1,n
633  contact(rcode)%tPenalty = tp
634  enddo
635  endif
636  if( damp_alpha > 0.0d0 ) then
637  do rcode=1,n
638  contact(rcode)%damp_alpha = damp_alpha
639  enddo
640  endif
641  if( damp_gact > 0.0d0 ) then
642  do rcode=1,n
643  contact(rcode)%damp_gact = damp_gact
644  enddo
645  endif
646 
647  fstr_ctrl_get_contact = .true.
648  end function fstr_ctrl_get_contact
649 
651  logical function fstr_ctrl_get_embed( ctrl, n, embed, cpname, smoothing )
652  use fstr_setup_util
653  integer(kind=kint), intent(in) :: ctrl
654  integer(kind=kint), intent(in) :: n
655  type(tcontact), intent(out) :: embed(n)
656  character(len=*), intent(out) :: cpname
657  integer(kind=kint), intent(out) :: smoothing
658 
659  integer :: rcode, ipt
660  character(len=30) :: s1 = 'TIED,GLUED,SSLID,FSLID '
661  character(len=HECMW_NAME_LEN) :: data_fmt,ss
662  character(len=HECMW_NAME_LEN) :: cp_name(n)
663  real(kind=kreal) :: fcoeff(n),tpenalty(n)
664 
665  tpenalty = 1.0d6
666 
667  write(ss,*) hecmw_name_len
668 
669  fstr_ctrl_get_embed = .false.
670  embed(1)%ctype = 1 ! pure slave-master contact; default value
671  embed(1)%algtype = contacttied ! small sliding contact; default value
672  if( fstr_ctrl_get_param_ex( ctrl, 'GRPID ', '# ', 1, 'I', embed(1)%group )/=0) return
673  smoothing = kcsnone + 1
674  if( fstr_ctrl_get_param_ex( ctrl, 'SMOOTHING ', 'NONE,NAGATA ', 0, 'P', smoothing ) /= 0 ) return
675  smoothing = smoothing - 1
676  do rcode=2,n
677  embed(rcode)%ctype = embed(1)%ctype
678  embed(rcode)%group = embed(1)%group
679  embed(rcode)%algtype = embed(1)%algtype
680  end do
681 
682  write( data_fmt, '(a,a)') 'S', trim(adjustl(ss))
683  if( fstr_ctrl_get_data_array_ex( ctrl, data_fmt, cp_name ) /= 0 ) return
684  do rcode=1,n
685  call fstr_strupr(cp_name(rcode))
686  embed(rcode)%pair_name = cp_name(rcode)
687  enddo
688 
689  cpname=""
690  if( fstr_ctrl_get_param_ex( ctrl, 'CONTACTPARAM ', '# ', 0, 'S', cpname )/= 0) return
691  fstr_ctrl_get_embed = .true.
692  end function
693 
695  function fstr_ctrl_get_contactparam( ctrl, contactparam )
696  implicit none
697  integer(kind=kint) :: ctrl
698  type( tcontactparam ) :: contactparam
699  integer(kind=kint) :: fstr_ctrl_get_contactparam
700 
701  integer(kind=kint) :: rcode
702  character(len=HECMW_NAME_LEN) :: data_fmt
703  character(len=128) :: msg
704  real(kind=kreal) :: clearance, clr_same_elem, clr_difflpos, clr_cal_norm
705  real(kind=kreal) :: distclr_init, distclr_free, distclr_nocheck, tensile_force
706  real(kind=kreal) :: box_exp_rate
707 
709 
710  !parameters
711  contactparam%name = ''
712  if( fstr_ctrl_get_param_ex( ctrl, 'NAME ', '# ', 1, 'S', contactparam%name ) /=0 ) return
713 
714  !read first line
715  data_fmt = 'rrrr '
716  rcode = fstr_ctrl_get_data_ex( ctrl, 1, data_fmt, &
717  & clearance, clr_same_elem, clr_difflpos, clr_cal_norm )
718  if( rcode /= 0 ) return
719  contactparam%CLEARANCE = clearance
720  contactparam%CLR_SAME_ELEM = clr_same_elem
721  contactparam%CLR_DIFFLPOS = clr_difflpos
722  contactparam%CLR_CAL_NORM = clr_cal_norm
723 
724  !read second line
725  data_fmt = 'rrrrr '
726  rcode = fstr_ctrl_get_data_ex( ctrl, 2, data_fmt, &
727  & distclr_init, distclr_free, distclr_nocheck, tensile_force, box_exp_rate )
728  if( rcode /= 0 ) return
729  contactparam%DISTCLR_INIT = distclr_init
730  contactparam%DISTCLR_FREE = distclr_free
731  contactparam%DISTCLR_NOCHECK = distclr_nocheck
732  contactparam%TENSILE_FORCE = tensile_force
733  contactparam%BOX_EXP_RATE = box_exp_rate
734 
735  !input check
736  rcode = 1
737  if( clearance<0.d0 .OR. 1.d0<clearance ) THEN
738  write(msg,*) 'fstr control file error : !CONTACT_PARAM : CLEARANCE must be 0 < CLEARANCE < 1.'
739  else if( clr_same_elem<0.d0 .or. 1.d0<clr_same_elem ) then
740  write(msg,*) 'fstr control file error : !CONTACT_PARAM : CLR_SAME_ELEM must be 0 < CLR_SAME_ELEM < 1.'
741  else if( clr_difflpos<0.d0 .or. 1.d0<clr_difflpos ) then
742  write(msg,*) 'fstr control file error : !CONTACT_PARAM : CLR_DIFFLPOS must be 0 < CLR_DIFFLPOS < 1.'
743  else if( clr_cal_norm<0.d0 .or. 1.d0<clr_cal_norm ) then
744  write(msg,*) 'fstr control file error : !CONTACT_PARAM : CLR_CAL_NORM must be 0 < CLR_CAL_NORM < 1.'
745  else if( distclr_init<0.d0 .or. 1.d0<distclr_init ) then
746  write(msg,*) 'fstr control file error : !CONTACT_PARAM : DISTCLR_INIT must be 0 < DISTCLR_INIT < 1.'
747  else if( distclr_free<-1.d0 .or. 1.d0<distclr_free ) then
748  write(msg,*) 'fstr control file error : !CONTACT_PARAM : DISTCLR_FREE must be -1 < DISTCLR_FREE < 1.'
749  else if( distclr_nocheck<0.5d0 ) then
750  write(msg,*) 'fstr control file error : !CONTACT_PARAM : DISTCLR_NOCHECK must be >= 0.5.'
751  else if( tensile_force>=0.d0 ) then
752  write(msg,*) 'fstr control file error : !CONTACT_PARAM : TENSILE_FORCE must be < 0.'
753  else if( box_exp_rate<=1.d0 .or. 2.0<box_exp_rate ) then
754  write(msg,*) 'fstr control file error : !CONTACT_PARAM : BOX_EXP_RATE must be 1 < BOX_EXP_RATE <= 2.'
755  else
756  rcode =0
757  end if
758  if( rcode /= 0 ) then
759  write(*,*) trim(msg)
760  write(ilog,*) trim(msg)
761  return
762  endif
763 
765  end function fstr_ctrl_get_contactparam
766 
768  function fstr_ctrl_get_contact_if( ctrl, n, contact_if )
769  use fstr_setup_util
770  integer(kind=kint), intent(in) :: ctrl
771  integer(kind=kint), intent(in) :: n
772  !
773  type(tcontactinterference), intent(out) :: contact_if(n)
774 
775  integer :: rcode, i
776  character(len=30) :: s1 = 'SLAVE,MASTER '
777  character(len=HECMW_NAME_LEN) :: data_fmt,ss
778  character(len=HECMW_NAME_LEN) :: cp_name(n)
779  real(kind=kreal) :: init_pos(n), end_pos(n)
780  integer(kind=kint) :: fstr_ctrl_get_contact_if
781 
783  write(ss,*) hecmw_name_len
784  if( fstr_ctrl_get_param_ex( ctrl, 'TYPE ', s1, 0, 'P', contact_if(1)%if_type ) /= 0 ) return
785  if( fstr_ctrl_get_param_ex( ctrl, 'END ', '# ', 0, 'R', contact_if(1)%etime ) /= 0 ) return
786  write( data_fmt, '(a,a,a)') 'S', trim(adjustl(ss)),'rr '
787  init_pos = 0.d0; end_pos = 0.d0
788  if( fstr_ctrl_get_data_array_ex( ctrl, data_fmt, cp_name, init_pos, end_pos ) /= 0 ) return
789  do i = 1, n
790  contact_if(i)%if_type = contact_if(1)%if_type
791  contact_if(i)%etime = contact_if(1)%etime
792 
793  contact_if(i)%cp_name = cp_name(i)
794  contact_if(i)%initial_pos = - init_pos(i)
795  contact_if(i)%end_pos = - end_pos(i)
796  if(contact_if(i)%if_type == c_if_slave .and. init_pos(i) /= 0.d0) contact_if(i)%initial_pos = 0.d0
797  end do
799 
800  end function fstr_ctrl_get_contact_if
801 
803  function fstr_ctrl_get_elemopt( ctrl, elemopt361 )
804  integer(kind=kint) :: ctrl
805  integer(kind=kint) :: elemopt361
806  integer(kind=kint) :: fstr_ctrl_get_elemopt
807 
808  character(72) :: o361list = 'IC,Bbar '
809 
810  integer(kind=kint) :: o361
811 
813 
814  o361 = elemopt361 + 1
815 
816  !* parameter in header line -----------------------------------------------------------------*!
817  if( fstr_ctrl_get_param_ex( ctrl, '361 ', o361list, 0, 'P', o361 ) /= 0) return
818 
819  elemopt361 = o361 - 1
820 
822 
823  end function fstr_ctrl_get_elemopt
824 
825 
827  function fstr_get_autoinc( ctrl, aincparam )
828  implicit none
829  integer(kind=kint) :: ctrl
830  type( tparamautoinc ) :: aincparam
831  integer(kind=kint) :: fstr_get_autoinc
832 
833  integer(kind=kint) :: rcode
834  character(len=HECMW_NAME_LEN) :: data_fmt
835  character(len=128) :: msg
836  integer(kind=kint) :: bound_s(10), bound_l(10)
837  real(kind=kreal) :: rs, rl
838 
839  fstr_get_autoinc = -1
840 
841  bound_s(:) = 0
842  bound_l(:) = 0
843 
844  !parameters
845  aincparam%name = ''
846  if( fstr_ctrl_get_param_ex( ctrl, 'NAME ', '# ', 1, 'S', aincparam%name ) /=0 ) return
847 
848  !read first line ( decrease criteria )
849  data_fmt = 'riiii '
850  rcode = fstr_ctrl_get_data_ex( ctrl, 1, data_fmt, rs, &
851  & bound_s(1), bound_s(2), bound_s(3), aincparam%NRtimes_s )
852  if( rcode /= 0 ) return
853  aincparam%ainc_Rs = rs
854  aincparam%NRbound_s(knstmaxit) = bound_s(1)
855  aincparam%NRbound_s(knstsumit) = bound_s(2)
856  aincparam%NRbound_s(knstciter) = bound_s(3)
857 
858  !read second line ( increase criteria )
859  data_fmt = 'riiii '
860  rcode = fstr_ctrl_get_data_ex( ctrl, 2, data_fmt, rl, &
861  & bound_l(1), bound_l(2), bound_l(3), aincparam%NRtimes_l )
862  if( rcode /= 0 ) return
863  aincparam%ainc_Rl = rl
864  aincparam%NRbound_l(knstmaxit) = bound_l(1)
865  aincparam%NRbound_l(knstsumit) = bound_l(2)
866  aincparam%NRbound_l(knstciter) = bound_l(3)
867 
868  !read third line ( cutback criteria )
869  data_fmt = 'ri '
870  rcode = fstr_ctrl_get_data_ex( ctrl, 3, data_fmt, &
871  & aincparam%ainc_Rc, aincparam%CBbound )
872  if( rcode /= 0 ) return
873 
874  !input check
875  rcode = 1
876  if( rs<0.d0 .or. rs>1.d0 ) then
877  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : decrease ratio Rs must 0 < Rs < 1.'
878  else if( any(bound_s<0) ) then
879  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : decrease NR bound must >= 0.'
880  else if( aincparam%NRtimes_s < 1 ) then
881  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : # of times to decrease must > 0.'
882  else if( rl<1.d0 ) then
883  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : increase ratio Rl must > 1.'
884  else if( any(bound_l<0) ) then
885  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : increase NR bound must >= 0.'
886  else if( aincparam%NRtimes_l < 1 ) then
887  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : # of times to increase must > 0.'
888  elseif( aincparam%ainc_Rc<0.d0 .or. aincparam%ainc_Rc>1.d0 ) then
889  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : cutback decrease ratio Rc must 0 < Rc < 1.'
890  else if( aincparam%CBbound < 1 ) then
891  write(msg,*) 'fstr control file error : !AUTOINC_PARAM : maximum # of cutback times must > 0.'
892  else
893  rcode =0
894  end if
895  if( rcode /= 0 ) then
896  write(*,*) trim(msg)
897  write(ilog,*) trim(msg)
898  return
899  endif
900 
901  fstr_get_autoinc = 0
902  end function fstr_get_autoinc
903 
905  function fstr_ctrl_get_timepoints( ctrl, tp )
906  integer(kind=kint) :: ctrl
907  type(time_points) :: tp
908  integer(kind=kint) :: fstr_ctrl_get_timepoints
909 
910  integer(kind=kint) :: i, n, rcode
911  logical :: generate
912  real(kind=kreal) :: stime, etime, interval
913 
915 
916  tp%name = ''
917  if( fstr_ctrl_get_param_ex( ctrl, 'NAME ', '# ', 1, 'S', tp%name ) /=0 ) return
918  tp%range_type = 1
919  if( fstr_ctrl_get_param_ex( ctrl, 'TIME ', 'STEP,TOTAL ', 0, 'P', tp%range_type ) /= 0 ) return
920  generate = .false.
921  if( fstr_ctrl_get_param_ex( ctrl, 'GENERATE ', '# ', 0, 'E', generate ) /= 0) return
922 
923  if( generate ) then
924  stime = 0.d0; etime = 0.d0; interval = 1.d0
925  if( fstr_ctrl_get_data_ex( ctrl, 1, 'rrr ', stime, etime, interval ) /= 0) return
926  tp%n_points = int((etime-stime)/interval)+1
927  allocate(tp%points(tp%n_points))
928  do i=1,tp%n_points
929  tp%points(i) = stime + dble(i-1)*interval
930  end do
931  else
932  n = fstr_ctrl_get_data_line_n( ctrl )
933  if( n == 0 ) return
934  tp%n_points = n
935  allocate(tp%points(tp%n_points))
936  if( fstr_ctrl_get_data_array_ex( ctrl, 'r ', tp%points ) /= 0 ) return
937  do i=1,tp%n_points-1
938  if( tp%points(i) < tp%points(i+1) ) cycle
939  write(*,*) 'Error in reading !TIME_POINT: time points must be given in ascending order.'
940  return
941  end do
942  end if
943 
945  end function fstr_ctrl_get_timepoints
946 
948  function fstr_ctrl_get_amplitude( ctrl, nline, name, type_def, type_time, type_val, n, val, table )
949  implicit none
950  integer(kind=kint), intent(in) :: ctrl
951  integer(kind=kint), intent(in) :: nline
952  character(len=HECMW_NAME_LEN), intent(out) :: name
953  integer(kind=kint), intent(out) :: type_def
954  integer(kind=kint), intent(out) :: type_time
955  integer(kind=kint), intent(out) :: type_val
956  integer(kind=kint), intent(out) :: n
957  real(kind=kreal), pointer :: val(:)
958  real(kind=kreal), pointer :: table(:)
959  integer(kind=kint) :: fstr_ctrl_get_amplitude
960 
961  integer(kind=kint) :: t_def, t_time, t_val
962  integer(kind=kint) :: i, j
963  real(kind=kreal) :: r(4), t(4)
964 
966 
967  name = ''
968  t_def = 1
969  t_time = 1
970  t_val = 1
971 
972  if( fstr_ctrl_get_param_ex( ctrl, 'NAME ', '# ', 1, 'S', name )/=0 ) return
973  if( fstr_ctrl_get_param_ex( ctrl, 'DEFINITION ', 'TABULAR ', 0, 'P', t_def )/=0 ) return
974  if( fstr_ctrl_get_param_ex( ctrl, 'TIME ', 'STEP ', 0, 'P', t_time )/=0 ) return
975  if( fstr_ctrl_get_param_ex( ctrl, 'VALUE ', 'RELATIVE,ABSOLUTE ', 0, 'P', t_val )/=0 ) return
976 
977  if( t_def==1 ) then
978  type_def = hecmw_amp_typedef_tabular
979  else
980  write(*,*) 'Error in reading !AMPLITUDE: invalid value for parameter DEFINITION.'
981  endif
982  if( t_time==1 ) then
983  type_time = hecmw_amp_typetime_step
984  else
985  write(*,*) 'Error in reading !AMPLITUDE: invalid value for parameter TIME.'
986  endif
987  if( t_val==1 ) then
988  type_val = hecmw_amp_typeval_relative
989  elseif( t_val==2 ) then
990  type_val = hecmw_amp_typeval_absolute
991  else
992  write(*,*) 'Error in reading !AMPLITUDE: invalid value for parameter VALUE.'
993  endif
994 
995  n = 0
996  do i = 1, nline
997  r(:)=huge(0.0d0); t(:)=huge(0.0d0)
998  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
999  n = n+1
1000  val(n) = r(1)
1001  table(n) = t(1)
1002  do j = 2, 4
1003  if (r(j) < huge(0.0d0) .and. t(j) < huge(0.0d0)) then
1004  n = n+1
1005  val(n) = r(j)
1006  table(n) = t(j)
1007  else
1008  exit
1009  endif
1010  enddo
1011  enddo
1013 
1014  end function fstr_ctrl_get_amplitude
1015 
1017  function fstr_ctrl_get_element_activation( ctrl, amp, eps, grp_id_name, mode, measure, state, thlow, thup )
1018  implicit none
1019  integer(kind=kint) :: ctrl
1020  character(len=HECMW_NAME_LEN) :: amp
1021  real(kind=kreal) :: eps
1022  character(len=HECMW_NAME_LEN),target :: grp_id_name(:)
1023  integer(kind=kint) :: mode ! 1=FIXED, 2=AMPLITUDE, 3=DAMAGE
1024  integer(kind=kint) :: measure ! 1=NONE, 2=STRESS, 3=STRAIN
1025  integer(kind=kint) :: state ! 0=ACTIVE, 1=INACTIVE
1026  real(kind=kreal), target :: thlow(:), thup(:)
1027  integer(kind=kint) :: fstr_ctrl_get_element_activation
1028 
1029  character(len=HECMW_NAME_LEN),pointer :: element_id_p
1030  real(kind=kreal),pointer :: thlow_p(:), thup_p(:)
1031  integer(kind=kint) :: rcode, n
1032  character(len=HECMW_NAME_LEN) :: data_fmt, s1
1033 
1035 
1036  ! MODE (required)
1037  if( fstr_ctrl_get_param_ex( ctrl, 'MODE ', 'FIXED,AMPLITUDE,DAMAGE ', 1, 'P', mode ) /= 0 ) return
1038 
1039  ! Mode-specific parameters
1040  if( mode == 1 ) then
1041  ! FIXED: STATE required
1042  if( fstr_ctrl_get_param_ex( ctrl, 'STATE ', 'ON,OFF ', 1, 'P', state ) /= 0 ) return
1043  state = state - 1
1044  measure = 1
1045  amp = ''
1046  elseif( mode == 2 ) then
1047  ! AMPLITUDE: AMP required
1048  if( fstr_ctrl_get_param_ex( ctrl, 'AMP ', '# ', 1, 'S', amp ) /= 0 ) return
1049  state = 0
1050  measure = 1
1051  elseif( mode == 3 ) then
1052  ! DAMAGE: MEASURE required
1053  if( fstr_ctrl_get_param_ex( ctrl, 'MEASURE ', 'STRESS,STRAIN ', 1, 'P', measure ) /= 0 ) return
1054  measure = measure + 1
1055  state = 0
1056  amp = ''
1057  endif
1058 
1059  ! EPSILON (optional)
1060  eps = 1.0d-6
1061  if( fstr_ctrl_get_param_ex( ctrl, 'EPSILON ', '# ', 0, 'R', eps ) /= 0 ) return
1062 
1063  write(s1,*) hecmw_name_len
1064  n = fstr_ctrl_get_data_line_n(ctrl)
1065  element_id_p => grp_id_name(1)
1066  thlow_p => thlow
1067  thup_p => thup
1068 
1069  if( mode == 3 ) then
1070  write( data_fmt, '(a,a,a)') 'S', trim(adjustl(s1)), 'RR'
1071  rcode = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, element_id_p, thlow_p, thup_p )
1072  else
1073  write( data_fmt, '(a,a)') 'S', trim(adjustl(s1))
1074  rcode = fstr_ctrl_get_data_array_ex( ctrl, data_fmt, element_id_p )
1075  endif
1076 
1078 
1080 
1081 
1082 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.
logical function fstr_ctrl_get_contact(ctrl, n, contact, np, tp, ntol, ttol, ctAlgo, cpname, smoothing)
Read in contact definition.
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_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, smoothing)
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:144
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.F90:98
integer(kind=kint), parameter imsg
Definition: m_fstr.F90:112
integer(kind=kint), parameter kstdynamic
Definition: m_fstr.F90:42
real(kind=kreal) etime
Definition: m_fstr.F90:142
integer(kind=kint), parameter kon
Definition: m_fstr.F90:34
integer(kind=kint), parameter ilog
FILE HANDLER.
Definition: m_fstr.F90:109
integer(kind=kint), parameter kststatic
Definition: m_fstr.F90:39
integer(kind=kint), parameter kststaticeigen
Definition: m_fstr.F90:44
This module manages step information.
Definition: m_out.f90:6
This module manages step information.
Definition: m_step.f90:6
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:653
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