SPEED
MAKE_SDOF_SYSTEM.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine make_sdof_system (file_bldinfo, dtsite, id)
 Read input parameters for oscillators.
 

Function/Subroutine Documentation

◆ make_sdof_system()

subroutine make_sdof_system ( character*70  file_bldinfo,
real*8  dtsite,
integer*4  id 
)

Read input parameters for oscillators.

Author
Aline Herlin
Date
November, 2020
Version
1.0
Parameters
[in]file_sdoffile name (SDOF00000X.info)
[in]idmpi process id
[in]dtsitesoil time step

Definition at line 27 of file MAKE_SDOF_SYSTEM.f90.

28
29 use speed_sci
30
31 implicit none
32
33 integer*4 :: id, SDOFunit, i, j, dummy
34 character*70 :: file_bldinfo
35 real*8 :: dtsite, tempint
36 real*8, dimension(4,4) :: mat_temp
37
38 sdofunit=1+(id+1)*10
39
40 open(sdofunit,file=file_bldinfo) !!! open SDOFINFO.txt
41 read(sdofunit,*) n_bld !!! here is the number of SDOF oscillators
42
43 do i=1, n_bld
44 sys(i)%SFS = 0
45 sys(i)%ForceApplicationtype = 1
46
47 read (sdofunit,*) sys(i)%ID, sys(i)%StructType, sys(i)%const_law
48
49 if (sys(i)%StructType.eq.1) then
50
51 sys(i)%NDOF = 1
52 sys(i)%flag_Minv = 0;
53
54 read (sdofunit,*) sys(i)%SFS, sys(i)%ForceApplicationtype
55 sys(i)%SFS = 0
56
57 allocate(sys(i)%Ms(1,1), sys(i)%Ms_inv(1,1), sys(i)%Ks(1,1))
58 allocate(sys(i)%IDR(1,3), sys(i)%variIDR(1,3))
59 allocate(sys(i)%IntForce(1,3))
60 allocate(sys(i)%tempU0(1,3), sys(i)%tempU1(1,3), sys(i)%tempU(1,3))
61 allocate(sys(i)%tempA1(sys(i)%NDOF,3), sys(i)%tempRA1(sys(i)%NDOF,3))
62
63 if (sys(i)%const_law.eq.1) then
64 !read (SDOFunit,"(1I10)") sys(i)%SFS
65 if (sys(i)%SFS.eq.0) then
66 read(sdofunit,*) sys(i)%Ms, sys(i)%Ks, sys(i)%CSI, sys(i)%TN
67 write(*,*) 'SDOF No: ', i, sys(i)%Ms, sys(i)%Ks, sys(i)%CSI, sys(i)%TN
68 elseif (sys(i)%SFS.eq.1) then
69 read(sdofunit,*) sys(i)%Ms, sys(i)%Mf, sys(i)%J, sys(i)%height, &
70 sys(i)%Ks, sys(i)%K0, sys(i)%Kr, sys(i)%Kv, sys(i)%CSI, &
71 sys(i)%C0, sys(i)%Cr, sys(i)%Cv, sys(i)%TN
72
73 read(sdofunit,*) sys(i)%beta_newmark, sys(i)%gamma_newmark
74
75 sys(i)%MAT_M = 0; sys(i)%MAT_KS = 0; sys(i)%MAT_KI = 0; sys(i)%MAT_C = 0
76
77 sys(i)%MAT_M(1,1) = sys(i)%Ms(1,1); sys(i)%MAT_M(2,2) = sys(i)%Mf;
78 sys(i)%MAT_M(3,3) = sys(i)%J; sys(i)%MAT_M(4,4) = sys(i)%Ms(1,1) + sys(i)%Mf
79
80 sys(i)%MAT_KS(1,1) = sys(i)%Ks(1,1);
81 sys(i)%MAT_KS(1,2) = -sys(i)%Ks(1,1);
82 sys(i)%MAT_KS(1,3) = -sys(i)%Ks(1,1)*sys(i)%height;
83 sys(i)%MAT_KS(2,1) = -sys(i)%Ks(1,1);
84 sys(i)%MAT_KS(2,2) = sys(i)%Ks(1,1);
85 sys(i)%MAT_KS(2,3) = sys(i)%Ks(1,1)*sys(i)%height;
86 sys(i)%MAT_KS(3,1) = -sys(i)%Ks(1,1)*sys(i)%height;
87 sys(i)%MAT_KS(3,2) = sys(i)%Ks(1,1)*sys(i)%height;
88 sys(i)%MAT_KS(3,3) = sys(i)%Ks(1,1)*(sys(i)%height**2)
89
90 sys(i)%MAT_KI(2,2) = sys(i)%K0; sys(i)%MAT_KI(3,3) = sys(i)%Kr; sys(i)%MAT_KI(4,4) = sys(i)%Kv
91
92 sys(i)%Cs = datan(1.d0)*16.d0*sys(i)%Ms(1,1)*sys(i)%CSI/sys(i)%TN
93
94 sys(i)%MAT_C(1,1) = sys(i)%Cs; sys(i)%MAT_C(1,2) = -sys(i)%Cs; sys(i)%MAT_C(1,3) = -sys(i)%Cs*sys(i)%height
95 sys(i)%MAT_C(2,1) = -sys(i)%Cs; sys(i)%MAT_C(2,2) = sys(i)%Cs + sys(i)%C0; sys(i)%MAT_C(2,3) = sys(i)%Cs*sys(i)%height
96 sys(i)%MAT_C(3,1) = -sys(i)%Cs*sys(i)%height; sys(i)%MAT_C(3,2) = sys(i)%Cs*sys(i)%height;
97 sys(i)%MAT_C(3,3) = sys(i)%Cs*(sys(i)%height**2) + sys(i)%Cr
98 sys(i)%MAT_C(4,4) = sys(i)%Cv
99
100 sys(i)%u = 0; sys(i)%v = 0; sys(i)%a = 0; sys(i)%f = 0
101 endif
102
103 elseif (sys(i)%const_law.eq.2) then
104 read(sdofunit,"(5E15.7)") sys(i)%Ms, sys(i)%Ks, sys(i)%CSI, sys(i)%FY, sys(i)%TN
105
106 elseif (sys(i)%const_law.eq.3) then
107 ! Mass, Initial Stiffness, post-yield stiffness, Softening stiffness, dampingfactor, yield strength, displacement at peak strength, ultimate displacement, natural time period
108 read(sdofunit,"(9E15.7)") sys(i)%Ms, sys(i)%Ks, sys(i)%Hs, sys(i)%Ss, sys(i)%CSI, sys(i)%FY, sys(i)%EH, sys(i)%EU, sys(i)%TN
109 sys(i)%EY = sys(i)%FY/sys(i)%Ks(1,1) !Yield Displacement
110 sys(i)%FH = sys(i)%FY + sys(i)%Hs*(sys(i)%EH - sys(i)%EY) !Max. Peak strength
111 sys(i)%FU = sys(i)%FH + sys(i)%Ss*(sys(i)%EU - sys(i)%EH) !Ultimate failure strength after softening
112 sys(i)%branch = 1
113 sys(i)%damage = 0
114
115 write(*,*) 'SDOF No: ', i, sys(i)%Ms, sys(i)%Ks, sys(i)%CSI, sys(i)%TN
116
117 endif
118
119 sys(i)%Cs = datan(1.d0)*16.d0*sys(i)%Ms(1,1)*sys(i)%CSI/sys(i)%TN ! Damping Coefficient
120
121 elseif (sys(i)%StructType.eq.2) then
122
123
124 read(sdofunit,*) sys(i)%NDOF, sys(i)%Floor_h, sys(i)%Area, sys(i)%MasspArea, sys(i)%ForceApplicationtype
125 sys(i)%ForceApplicationtype = 1
126
127 sys(i)%Height=sys(i)%NDOF*sys(i)%Floor_h
128
129 allocate(sys(i)%props(sys(i)%NDOF,10,1))
130 allocate(sys(i)%dval(sys(i)%NDOF,4))
131 allocate(sys(i)%Ms(sys(i)%NDOF,sys(i)%NDOF))
132 allocate(sys(i)%Ks(sys(i)%NDOF,sys(i)%NDOF))
133 allocate(sys(i)%SysC(sys(i)%NDOF,sys(i)%NDOF))
134 allocate(sys(i)%Ms_Inv(sys(i)%NDOF,sys(i)%NDOF))
135 allocate(sys(i)%MDOFEt(sys(i)%NDOF,2))
136 allocate(sys(i)%MDOFstatev(sys(i)%NDOF,11,2))
137 allocate(sys(i)%MDOFspd(sys(i)%NDOF,2))
138 allocate(sys(i)%MDOFyield(sys(i)%NDOF,2))
139 allocate(sys(i)%IDR(sys(i)%NDOF,3))
140 allocate(sys(i)%variIDR(sys(i)%NDOF,3))
141 allocate(sys(i)%IntForce(sys(i)%NDOF,3))
142 allocate(sys(i)%tempU0(sys(i)%NDOF,3))
143 allocate(sys(i)%tempU1(sys(i)%NDOF,3))
144 allocate(sys(i)%tempU(sys(i)%NDOF,3))
145 allocate(sys(i)%MaxIDR(sys(i)%NDOF,2))
146 allocate(sys(i)%MDOFIDeath(sys(i)%NDOF,2))
147 allocate(sys(i)%MDOFstate(sys(i)%NDOF,2))
148 allocate(sys(i)%tempA1(sys(i)%NDOF,3))
149 allocate(sys(i)%tempRA1(sys(i)%NDOF,3))
150 sys(i)%Ms=0; sys(i)%Ks=0; sys(i)%SysC=0; sys(i)%Ms_Inv=0;
151 sys(i)%MDOFEt=0; sys(i)%MDOFstatev=0; sys(i)%MDOFspd=0; sys(i)%MDOFyield=0;
152 sys(i)%MaxIDR=0; sys(i)%MDOFIDeath=0; sys(i)%MDOFstate=0
153 sys(i)%flag_Minv = 0;
154
155 do j=1, sys(i)%NDOF
156 read(sdofunit,*)sys(i)%props(j,1:10,1)
157 enddo
158
159 endif
160
161 sys(i)%IDR=0
162 sys(i)%variIDR=0
163 sys(i)%IntForce=0
164 sys(i)%tempU0=0
165 sys(i)%tempU1=0
166 sys(i)%tempU=0
167
168 if (sys(i)%SFS.eq.1) then
169 sys(i)%MAT_F = 1/(sys(i)%beta_newmark*sys(i)%dt2) * sys(i)%MAT_M + &
170 sys(i)%gamma_newmark/(sys(i)%beta_newmark*sys(i)%dt)*sys(i)%MAT_C + sys(i)%MAT_KS
171 mat_temp = sys(i)%MAT_F + sys(i)%MAT_KI
172 call matinv(mat_temp, sys(i)%MAT_MCinv, 4)
173 endif
174 enddo
175
176 maxdof_loc = 0
177 do i=1, n_bld
178 if (maxdof_loc.lt.sys(i)%NDOF) then
179 maxdof_loc = sys(i)%NDOF
180 endif
181
182 if (sys(i)%StructType.eq.2) then
183 read(sdofunit,*)
184 do j=1, sys(i)%NDOF
185 read(sdofunit,*)sys(i)%dval(j,1:4)
186 enddo
187 endif
188 enddo
189
190 if (sys(1)%StructType.eq.2) read(sdofunit,*)
191 do i=1, n_bld
192 if (sys(i)%StructType.eq.2) then
193 read(sdofunit,*)dummy, sys(i)%TN, sys(i)%T2, sys(i)%Calpha, sys(i)%Cbeta
194 endif
195 enddo
196
197 close(sdofunit)
198
199
200 ! Making Mass, Damping Matrices
201 do i=1, n_bld
202 sys(i)%dt=sys(i)%TN/(datan(1.d0)*4.d0) !- Stable time step for explicit central difference method for SDOF
203 sys(i)%ndt=ceiling(dtsite/sys(i)%dt)
204 sys(i)%dt=dtsite/sys(i)%ndt
205 sys(i)%dt2=sys(i)%dt*sys(i)%dt
206
207 !!!!!!!!!!!MassMatrix!!!!!!!!!!!!
208 if (sys(i)%StructType.eq.2) then
209 do j=1,sys(i)%NDOF
210 sys(i)%Ms(j,j)=sys(i)%MasspArea*sys(i)%Area
211 enddo
212 endif
213
214 !!!!!!!!!!!StiffMatrix!!!!!!!!!!!!
215 if (sys(i)%StructType.eq.2) then
216 do j=1,sys(i)%NDOF
217 sys(i)%Ks(j,j)= sys(i)%Ks(j,j) + sys(i)%props(j,1,1)
218 if(j>1) then
219 sys(i)%Ks(j-1,j-1)= sys(i)%Ks(j-1,j-1) + sys(i)%props(j,1,1)
220 sys(i)%Ks(j,j-1) = sys(i)%Ks(j,j-1) - sys(i)%props(j,1,1)
221 sys(i)%Ks(j-1,j) = sys(i)%Ks(j-1,j) - sys(i)%props(j,1,1)
222 endif
223 enddo
224 endif
225
226 !!!!!!!!!!! DampingMatrix - Rayleigh !!!!!!!!!!!!
227 if (sys(i)%StructType.eq.2) then
228 sys(i)%SysC(1:sys(i)%NDOF,1:sys(i)%NDOF)= sys(i)%Calpha * sys(i)%Ms(1:sys(i)%NDOF,1:sys(i)%NDOF) &
229 + sys(i)%Cbeta * sys(i)%Ks(1:sys(i)%NDOF,1:sys(i)%NDOF)
230 deallocate(sys(i)%Ks)
231 endif
232 enddo
233
234 return
subroutine matinv(a, b, n)
Contains parameters for MDOF.
Definition MODULES.f90:725
integer *4 maxdof_loc
Definition MODULES.f90:775
integer *4 n_bld
Definition MODULES.f90:773
type(system), dimension(:), allocatable sys
SDOF system.
Definition MODULES.f90:771

References matinv(), speed_sci::maxdof_loc, speed_sci::n_bld, and speed_sci::sys.

Referenced by read_sdof_input_files().

Here is the call graph for this function:
Here is the caller graph for this function: