20 subroutine pc_strupr( s )
23 integer :: i, n, a, da
26 da = iachar(
'a') - iachar(
'A')
29 if( a > iachar(
'Z'))
then
34 end subroutine pc_strupr
39 integer(kind=kint) :: ctrl
40 integer(kind=kint) :: type
44 integer(kind=kint) :: ipt
45 character(len=80) :: s
49 s =
'ELEMCHECK,STATIC,EIGEN,HEAT,DYNAMIC,NLSTATIC,STATICEIGEN,NZPROF '
68 integer(kind=kint) :: ctrl
69 integer(kind=kint) :: method
72 integer(kind=kint) :: ipt
73 character(len=80) :: s
77 s =
'NEWTON,QUASINEWTON '
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
117 character(92) :: mlist =
'1,2,3,4,101,CG,BiCGSTAB,GMRES,GPBiCG,GMRESR,GMRESREN,DIRECT,DIRECTmkl,DIRECTlag,MUMPS,MKL '
119 character(24) :: dlist =
'0,1,2,3,NONE,MM,CSR,BSR '
121 integer(kind=kint) :: number_number = 5
122 integer(kind=kint) :: indirect_number = 6
123 integer(kind=kint) :: iter, time, sclg, dmpt, dmpx, usjd, step
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
152 if( method > number_number )
then
153 method = method - number_number
154 if( method > indirect_number )
then
156 method = method - indirect_number + 100
157 if( method == 103 ) method = 101
158 if( method == 105 ) method = 102
161 if( method2 > number_number )
then
162 method2 = method2 - number_number
163 if( method2 > indirect_number )
then
165 method2 = method2 - indirect_number + 100
170 if( dumptype >= 4 )
then
171 dumptype = dumptype - 4
176 if(
fstr_ctrl_get_data_ex( ctrl, 1,
'iiiiii ', nier, iterpremax, nrest, ncolor_in, recyclepre, nbfgs )/= 0)
return
179 if( precond == 20 .or. precond == 21)
then
181 else if( precond == 5 )
then
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
203 integer(kind=kint) :: ctrl
204 character(len=HECMW_NAME_LEN) :: amp
205 integer(kind=kint) :: iproc
208 integer(kind=kint) :: ipt = 0
209 integer(kind=kint) :: ip = 0
217 if( ipt == 2 .or. ip == 1 ) iproc = 1
227 integer(kind=kint),
intent(in) :: ctrl
228 type (hecmwst_local_mesh),
intent(in) :: hecmesh
230 character(len=*),
intent(out) :: tpname
231 character(len=*),
intent(out) :: apname
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
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 '
253 steps%initdt = 1.d0/steps%num_substep
263 if( len( trim(amp) )>0 )
then
279 read( ss, * , iostat=ierr ) fn
283 steps%elapsetime = f1
285 steps%mindt = min(f2,steps%initdt)
288 steps%num_substep = max(int((f1+0.999999999d0*fn)/fn),steps%num_substep)
299 if( trim(header_name) ==
'BOUNDARY' )
then
301 else if( trim(header_name) ==
'LOAD' )
then
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
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) )
323 if( trim(header_name) ==
'BOUNDARY' )
then
325 steps%Boundary(bc_n) = bcid
326 else if( trim(header_name) ==
'LOAD' )
then
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
344 integer(kind=kint),
intent(in) :: ctrl
345 type (hecmwst_local_mesh),
intent(inout) :: hecmesh
346 type (
tsection),
pointer,
intent(inout) :: sections(:)
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 '
357 if( sect_id > hecmesh%section%n_sect )
return
361 if( elemopt > 0 ) sections(sect_id)%elemopt341 = elemopt
365 if( elemopt > 0 ) sections(sect_id)%elemopt361 = elemopt
368 hecmesh%section%sect_orien_ID(sect_id) = -1
371 if(
associated(g_localcoordsys) )
then
373 k =
size(g_localcoordsys)
376 if( sect_orien == g_localcoordsys(cache)%sys_name )
then
377 hecmesh%section%sect_orien_ID(sect_id) = cache
385 if( sect_orien == g_localcoordsys(j)%sys_name )
then
386 hecmesh%section%sect_orien_ID(sect_id) = j
400 integer(kind=kint) :: ctrl
401 integer(kind=kint) :: res
402 integer(kind=kint) :: visual
403 integer(kind=kint) :: femap
419 integer(kind=kint) :: ctrl
420 integer(kind=kint) :: echo
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
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)),
' '
444 if(
fstr_ctrl_get_param_ex( ctrl,
'TYPE ',
'1,2,3,4,5,6 ', 0,
'I', fg_type )/= 0)
return
455 integer(kind=kint),
intent(in) :: ctrl
456 real(kind=kreal),
intent(out) :: penalty
460 if( penalty <= 1.0 )
then
462 write(
imsg,*)
"Warging : !MPC : too small penalty: ", penalty
463 write(*,*)
"Warging : !MPC : too small penalty: ", penalty
473 integer(kind=kint),
intent(in) :: ctrl
474 type (hecmwst_local_mesh),
intent(in) :: hecmesh
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(:)
482 write( ss, * ) hecmw_name_len
483 write( data_fmt,
'(a,a,a,a,a)')
'S', trim(adjustl(ss)),
'S', trim(adjustl(ss)),
' '
488 outinfo%grp_id_name =
"ALL"
496 allocate( header_name(n), onoff(n), vtype(n) )
497 header_name(:) =
""; vtype(:) =
""; onoff(:) =
""
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
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
523 deallocate( header_name, onoff, vtype )
530 integer(kind=kint) :: ctrl
531 integer(kind=kint) :: algo
534 integer(kind=kint) :: rcode
535 character(len=80) :: s
537 s =
'SLAGRANGE,ALAGRANGE '
545 integer(kind=kint),
intent(in) :: ctrl
546 integer(kind=kint),
intent(in) :: n
547 integer(kind=kint),
intent(in) :: ctalgo
548 type(tcontact),
intent(out) :: contact(n)
549 real(kind=kreal),
intent(out) :: np
550 real(kind=kreal),
intent(out) :: tp
551 real(kind=kreal),
intent(out) :: ntol
552 real(kind=kreal),
intent(out) :: ttol
553 character(len=*),
intent(out) :: cpname
555 integer :: rcode, ipt
556 character(len=30) :: s1 =
'TIED,GLUED,SSLID,FSLID '
557 character(len=HECMW_NAME_LEN) :: data_fmt,ss
558 character(len=HECMW_NAME_LEN) :: cp_name(n)
559 real(kind=kreal) :: fcoeff(n),tpenalty(n)
563 write(ss,*) hecmw_name_len
567 contact(1)%algtype = contactsslid
569 if( contact(1)%algtype==contactglued ) contact(1)%algtype=contactfslid
572 contact(rcode)%ctype = contact(1)%ctype
573 contact(rcode)%group = contact(1)%group
574 contact(rcode)%algtype = contact(1)%algtype
577 if( contact(1)%algtype==contactsslid .or. contact(1)%algtype==contactfslid )
then
578 write( data_fmt,
'(a,a,a)')
'S', trim(adjustl(ss)),
'Rr '
582 contact(rcode)%pair_name = cp_name(rcode)
583 contact(rcode)%fcoeff = fcoeff(rcode)
584 contact(rcode)%tPenalty = tpenalty(rcode)
586 else if( contact(1)%algtype==contacttied )
then
587 write( data_fmt,
'(a,a)')
'S', trim(adjustl(ss))
591 contact(rcode)%pair_name = cp_name(rcode)
592 contact(rcode)%fcoeff = 0.d0
593 contact(rcode)%tPenalty = 1.d+4
598 ntol = 0.d0; ttol=0.d0
611 integer(kind=kint),
intent(in) :: ctrl
612 integer(kind=kint),
intent(in) :: n
613 type(tcontact),
intent(out) :: embed(n)
614 character(len=*),
intent(out) :: cpname
616 integer :: rcode, ipt
617 character(len=30) :: s1 =
'TIED,GLUED,SSLID,FSLID '
618 character(len=HECMW_NAME_LEN) :: data_fmt,ss
619 character(len=HECMW_NAME_LEN) :: cp_name(n)
620 real(kind=kreal) :: fcoeff(n),tpenalty(n)
624 write(ss,*) hecmw_name_len
628 embed(1)%algtype = contacttied
631 embed(rcode)%ctype = embed(1)%ctype
632 embed(rcode)%group = embed(1)%group
633 embed(rcode)%algtype = embed(1)%algtype
636 write( data_fmt,
'(a,a)')
'S', trim(adjustl(ss))
640 embed(rcode)%pair_name = cp_name(rcode)
651 integer(kind=kint) :: ctrl
652 type( tcontactparam ) :: contactparam
655 integer(kind=kint) :: rcode
656 character(len=HECMW_NAME_LEN) :: data_fmt
657 character(len=128) :: msg
658 real(kind=kreal) :: clearance, clr_same_elem, clr_difflpos, clr_cal_norm
659 real(kind=kreal) :: distclr_init, distclr_free, distclr_nocheck, tensile_force
660 real(kind=kreal) :: box_exp_rate
665 contactparam%name =
''
671 & clearance, clr_same_elem, clr_difflpos, clr_cal_norm )
672 if( rcode /= 0 )
return
673 contactparam%CLEARANCE = clearance
674 contactparam%CLR_SAME_ELEM = clr_same_elem
675 contactparam%CLR_DIFFLPOS = clr_difflpos
676 contactparam%CLR_CAL_NORM = clr_cal_norm
681 & distclr_init, distclr_free, distclr_nocheck, tensile_force, box_exp_rate )
682 if( rcode /= 0 )
return
683 contactparam%DISTCLR_INIT = distclr_init
684 contactparam%DISTCLR_FREE = distclr_free
685 contactparam%DISTCLR_NOCHECK = distclr_nocheck
686 contactparam%TENSILE_FORCE = tensile_force
687 contactparam%BOX_EXP_RATE = box_exp_rate
691 if( clearance<0.d0 .OR. 1.d0<clearance )
THEN
692 write(msg,*)
'fstr control file error : !CONTACT_PARAM : CLEARANCE must be 0 < CLEARANCE < 1.'
693 else if( clr_same_elem<0.d0 .or. 1.d0<clr_same_elem )
then
694 write(msg,*)
'fstr control file error : !CONTACT_PARAM : CLR_SAME_ELEM must be 0 < CLR_SAME_ELEM < 1.'
695 else if( clr_difflpos<0.d0 .or. 1.d0<clr_difflpos )
then
696 write(msg,*)
'fstr control file error : !CONTACT_PARAM : CLR_DIFFLPOS must be 0 < CLR_DIFFLPOS < 1.'
697 else if( clr_cal_norm<0.d0 .or. 1.d0<clr_cal_norm )
then
698 write(msg,*)
'fstr control file error : !CONTACT_PARAM : CLR_CAL_NORM must be 0 < CLR_CAL_NORM < 1.'
699 else if( distclr_init<0.d0 .or. 1.d0<distclr_init )
then
700 write(msg,*)
'fstr control file error : !CONTACT_PARAM : DISTCLR_INIT must be 0 < DISTCLR_INIT < 1.'
701 else if( distclr_free<-1.d0 .or. 1.d0<distclr_free )
then
702 write(msg,*)
'fstr control file error : !CONTACT_PARAM : DISTCLR_FREE must be -1 < DISTCLR_FREE < 1.'
703 else if( distclr_nocheck<0.5d0 )
then
704 write(msg,*)
'fstr control file error : !CONTACT_PARAM : DISTCLR_NOCHECK must be >= 0.5.'
705 else if( tensile_force>=0.d0 )
then
706 write(msg,*)
'fstr control file error : !CONTACT_PARAM : TENSILE_FORCE must be < 0.'
707 else if( box_exp_rate<=1.d0 .or. 2.0<box_exp_rate )
then
708 write(msg,*)
'fstr control file error : !CONTACT_PARAM : BOX_EXP_RATE must be 1 < BOX_EXP_RATE <= 2.'
712 if( rcode /= 0 )
then
714 write(
ilog,*) trim(msg)
724 integer(kind=kint),
intent(in) :: ctrl
725 integer(kind=kint),
intent(in) :: n
727 type(tcontactinterference),
intent(out) :: contact_if(n)
730 character(len=30) :: s1 =
'SLAVE,MASTER '
731 character(len=HECMW_NAME_LEN) :: data_fmt,ss
732 character(len=HECMW_NAME_LEN) :: cp_name(n)
733 real(kind=kreal) :: init_pos(n), end_pos(n)
737 write(ss,*) hecmw_name_len
740 write( data_fmt,
'(a,a,a)')
'S', trim(adjustl(ss)),
'rr '
741 init_pos = 0.d0; end_pos = 0.d0
744 contact_if(i)%if_type = contact_if(1)%if_type
745 contact_if(i)%etime = contact_if(1)%etime
747 contact_if(i)%cp_name = cp_name(i)
748 contact_if(i)%initial_pos = - init_pos(i)
749 contact_if(i)%end_pos = - end_pos(i)
750 if(contact_if(i)%if_type == c_if_slave .and. init_pos(i) /= 0.d0) contact_if(i)%initial_pos = 0.d0
758 integer(kind=kint) :: ctrl
759 integer(kind=kint) :: elemopt361
762 character(72) :: o361list =
'IC,Bbar '
764 integer(kind=kint) :: o361
768 o361 = elemopt361 + 1
773 elemopt361 = o361 - 1
783 integer(kind=kint) :: ctrl
784 type( tparamautoinc ) :: aincparam
787 integer(kind=kint) :: rcode
788 character(len=HECMW_NAME_LEN) :: data_fmt
789 character(len=128) :: msg
790 integer(kind=kint) :: bound_s(10), bound_l(10)
791 real(kind=kreal) :: rs, rl
805 & bound_s(1), bound_s(2), bound_s(3), aincparam%NRtimes_s )
806 if( rcode /= 0 )
return
807 aincparam%ainc_Rs = rs
808 aincparam%NRbound_s(knstmaxit) = bound_s(1)
809 aincparam%NRbound_s(knstsumit) = bound_s(2)
810 aincparam%NRbound_s(knstciter) = bound_s(3)
815 & bound_l(1), bound_l(2), bound_l(3), aincparam%NRtimes_l )
816 if( rcode /= 0 )
return
817 aincparam%ainc_Rl = rl
818 aincparam%NRbound_l(knstmaxit) = bound_l(1)
819 aincparam%NRbound_l(knstsumit) = bound_l(2)
820 aincparam%NRbound_l(knstciter) = bound_l(3)
825 & aincparam%ainc_Rc, aincparam%CBbound )
826 if( rcode /= 0 )
return
830 if( rs<0.d0 .or. rs>1.d0 )
then
831 write(msg,*)
'fstr control file error : !AUTOINC_PARAM : decrease ratio Rs must 0 < Rs < 1.'
832 else if( any(bound_s<0) )
then
833 write(msg,*)
'fstr control file error : !AUTOINC_PARAM : decrease NR bound must >= 0.'
834 else if( aincparam%NRtimes_s < 1 )
then
835 write(msg,*)
'fstr control file error : !AUTOINC_PARAM : # of times to decrease must > 0.'
836 else if( rl<1.d0 )
then
837 write(msg,*)
'fstr control file error : !AUTOINC_PARAM : increase ratio Rl must > 1.'
838 else if( any(bound_l<0) )
then
839 write(msg,*)
'fstr control file error : !AUTOINC_PARAM : increase NR bound must >= 0.'
840 else if( aincparam%NRtimes_l < 1 )
then
841 write(msg,*)
'fstr control file error : !AUTOINC_PARAM : # of times to increase must > 0.'
842 elseif( aincparam%ainc_Rc<0.d0 .or. aincparam%ainc_Rc>1.d0 )
then
843 write(msg,*)
'fstr control file error : !AUTOINC_PARAM : cutback decrease ratio Rc must 0 < Rc < 1.'
844 else if( aincparam%CBbound < 1 )
then
845 write(msg,*)
'fstr control file error : !AUTOINC_PARAM : maximum # of cutback times must > 0.'
849 if( rcode /= 0 )
then
851 write(
ilog,*) trim(msg)
860 integer(kind=kint) :: ctrl
864 integer(kind=kint) :: i, n, rcode
866 real(kind=kreal) :: stime,
etime, interval
878 stime = 0.d0;
etime = 0.d0; interval = 1.d0
880 tp%n_points = int((
etime-stime)/interval)+1
881 allocate(tp%points(tp%n_points))
883 tp%points(i) = stime + dble(i-1)*interval
889 allocate(tp%points(tp%n_points))
892 if( tp%points(i) < tp%points(i+1) ) cycle
893 write(*,*)
'Error in reading !TIME_POINT: time points must be given in ascending order.'
904 integer(kind=kint),
intent(in) :: ctrl
905 integer(kind=kint),
intent(in) :: nline
906 character(len=HECMW_NAME_LEN),
intent(out) :: name
907 integer(kind=kint),
intent(out) :: type_def
908 integer(kind=kint),
intent(out) :: type_time
909 integer(kind=kint),
intent(out) :: type_val
910 integer(kind=kint),
intent(out) :: n
911 real(kind=kreal),
pointer :: val(:)
912 real(kind=kreal),
pointer :: table(:)
915 integer(kind=kint) :: t_def, t_time, t_val
916 integer(kind=kint) :: i, j
917 real(kind=kreal) :: r(4), t(4)
932 type_def = hecmw_amp_typedef_tabular
934 write(*,*)
'Error in reading !AMPLITUDE: invalid value for parameter DEFINITION.'
937 type_time = hecmw_amp_typetime_step
939 write(*,*)
'Error in reading !AMPLITUDE: invalid value for parameter TIME.'
942 type_val = hecmw_amp_typeval_relative
943 elseif( t_val==2 )
then
944 type_val = hecmw_amp_typeval_absolute
946 write(*,*)
'Error in reading !AMPLITUDE: invalid value for parameter VALUE.'
951 r(:)=huge(0.0d0); t(:)=huge(0.0d0)
952 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
957 if (r(j) < huge(0.0d0) .and. t(j) < huge(0.0d0))
then
976 integer(kind=kint) :: ctrl
977 character(len=HECMW_NAME_LEN) :: amp
978 real(kind=kreal) ::
eps
979 character(len=HECMW_NAME_LEN),
target :: grp_id_name(:)
980 integer(kind=kint) :: dtype
981 integer(kind=kint) :: state
982 real(kind=kreal),
target :: thlow(:), thup(:)
985 character(len=HECMW_NAME_LEN),
pointer :: element_id_p
986 real(kind=kreal),
pointer :: thlow_p(:), thup_p(:)
988 integer(kind=kint) :: i, n
989 integer(kind=kint) :: rcode
990 character(len=HECMW_NAME_LEN) :: data_fmt,s1
991 integer(kind=kint) :: lid
992 character(len=HECMW_NAME_LEN) :: state_str
1000 state_str =
'ON,OFF'
1005 s1 =
'NONE,STRESS,STRAIN '
1008 write(s1,*) hecmw_name_len
1013 element_id_p => grp_id_name(1)
1017 if( dtype == 1 )
then
1018 write( data_fmt,
'(a,a)')
'S', trim(adjustl(s1))
1021 write( data_fmt,
'(a,a,a)')
'S', trim(adjustl(s1)),
'RR'
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_contactalgo(ctrl, algo)
Read in !CONTACT.
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_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.
integer(kind=kint) function fstr_ctrl_get_element_activation(ctrl, amp, eps, grp_id_name, dtype, state, thlow, thup)
Read in !ELEMENT_ACTIVATION.
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)
This module defines common data and basic structures for analysis.
integer(kind=kint) myrank
PARALLEL EXECUTION.
integer(kind=kint), parameter imsg
integer(kind=kint), parameter kstdynamic
integer(kind=kint), parameter kon
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
integer(kind=kint), parameter ilog
FILE HANDLER.
integer(kind=kint), parameter kststatic
integer(kind=kint), parameter kststaticeigen
This module manages step information.
This module manages step information.
subroutine init_stepinfo(stepinfo)
Initializer.
integer, parameter stepfixedinc
integer, parameter stepautoinc
integer, parameter stepstatic
This module manages timepoint information.
Data for section control.
Step control such as active boundary condition, convergent condition etc.
Time points storage for output etc.