11 private :: fstr_ctrl_get_nval
18 function fstr_ctrl_get_nval( ctrl, iType, amp, node_id, node_id_len, dof_ids, dof_ide, value )
20 integer(kind=kint) :: ctrl
21 integer(kind=kint) :: iType
22 character(len=HECMW_NAME_LEN) :: amp
23 character(len=HECMW_NAME_LEN) :: node_id(:)
24 integer(kind=kint) :: node_id_len
25 integer(kind=kint),
pointer :: dof_ids (:)
26 integer(kind=kint),
pointer :: dof_ide (:)
27 real(kind=kreal),
pointer :: value(:)
28 integer(kind=kint) :: fstr_ctrl_get_nval
29 integer(kind=kint) :: rcode
31 character(len=HECMW_NAME_LEN) :: data_fmt,ss
32 write(ss,*) node_id_len
33 write(data_fmt,
'(a,a,a)')
'S',trim(adjustl(ss)),
'IIr '
37 fstr_ctrl_get_nval = -1
41 fstr_ctrl_get_nval = &
44 end function fstr_ctrl_get_nval
49 integer(kind=kint) :: ctrl
50 integer(kind=kint) vtype
51 character(len=HECMW_NAME_LEN) :: amp
52 character(len=HECMW_NAME_LEN) :: node_id(:)
53 integer(kind=kint) :: node_id_len
54 integer(kind=kint),
pointer :: dof_ids (:)
55 integer(kind=kint),
pointer :: dof_ide (:)
56 real(kind=kreal),
pointer :: value(:)
60 fstr_ctrl_get_nval( ctrl, vtype, amp, node_id, node_id_len, dof_ids, dof_ide,
value )
67 integer(kind=kint) :: ctrl
68 integer(kind=kint) :: atype
69 character(len=HECMW_NAME_LEN) :: amp
70 character(len=HECMW_NAME_LEN) :: node_id(:)
71 integer(kind=kint) :: node_id_len
72 integer(kind=kint),
pointer :: dof_ids (:)
73 integer(kind=kint),
pointer :: dof_ide (:)
74 real(kind=kreal),
pointer :: value(:)
78 fstr_ctrl_get_nval( ctrl, atype, amp, node_id, node_id_len, dof_ids, dof_ide,
value )
84 idx_eqa, idx_resp, n_step, t_start, t_end, t_delta, &
85 gamma, beta, idx_mas, idx_dmp, ray_m, ray_k, &
86 nout, node_id, node_id_len, nout_monit, iout_list )
88 integer(kind=kint) :: ctrl
92 integer(kind=kint) :: idx_eqa
93 integer(kind=kint) :: idx_resp
96 integer(kind=kint) :: n_step
97 real(kind=kreal) :: t_start
98 real(kind=kreal) :: t_end
99 real(kind=kreal) :: t_delta
102 real(kind=kreal) :: gamma
103 real(kind=kreal) :: beta
106 integer(kind=kint) :: idx_mas
109 integer(kind=kint) :: idx_dmp
110 real(kind=kreal) :: ray_m
111 real(kind=kreal) :: ray_k
114 integer(kind=kint) :: nout
115 character(len=HECMW_NAME_LEN) :: node_id
116 integer(kind=kint) :: node_id_len
117 integer(kind=kint) :: nout_monit
118 integer(kind=kint) :: iout_list(6)
120 integer(kind=kint) :: rcode, nlflag
121 character(len=80) :: s
122 character(len=HECMW_NAME_LEN) :: data_fmt,ss
129 s =
'LINEAR,NONLINEAR '
132 if( nlflag/=0 ) nlgeom = (nlflag==2)
138 write(ss,*) node_id_len
139 write(data_fmt,
'(a,a,a)')
'iS',trim(adjustl(ss)),
'i '
142 iout_list(1), iout_list(2), iout_list(3), iout_list(4), iout_list(5), iout_list(6) )/=0 )
return
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 control file data obtaining functions for dynamic analysis.
integer(kind=kint) function fstr_ctrl_get_dynamic(ctrl, nlgeom, idx_eqa, idx_resp, n_step, t_start, t_end, t_delta, gamma, beta, idx_mas, idx_dmp, ray_m, ray_k, nout, node_id, node_id_len, nout_monit, iout_list)
Read in !DYNAMIC.
integer(kind=kint) function fstr_ctrl_get_velocity(ctrl, vType, amp, node_id, node_id_len, dof_ids, dof_ide, value)
Read in !VELOCITY.
integer(kind=kint) function fstr_ctrl_get_acceleration(ctrl, aType, amp, node_id, node_id_len, dof_ids, dof_ide, value)
Read in !ACCELERATION.
This module defines common data and basic structures for analysis.