14 type(hecmwst_local_mesh) :: hecMESH
15 type(hecmwst_matrix) :: hecMAT
18 integer(kind=kint) :: N, NP, NDOF, NNDOF, NPNDOF
19 integer(kind=kint) :: i, j, k, in, jn, nget
20 real(kind=kreal) :: chk, gm
21 real(kind=kreal),
allocatable :: s(:), t(:), u(:), r(:)
22 real(kind=kreal),
pointer :: mass(:), eigval(:), eigvec(:,:)
32 eigval => fstreig%eigval
33 eigvec => fstreig%eigvec
35 allocate(fstreig%effmass(ndof*nget))
36 allocate(fstreig%partfactor(ndof*nget))
41 fstreig%effmass = 0.0d0
42 fstreig%partfactor = 0.0d0
50 r(k) = r(k) + mass(in)*eigvec(in,i)
51 gm = gm + mass(in)*eigvec(in,i)*eigvec(in,i)
55 call hecmw_allreduce_r(hecmesh, r, ndof, hecmw_sum)
56 call hecmw_allreduce_r1(hecmesh, gm, hecmw_sum)
61 fstreig%partfactor(in) = gm*r(j)
62 fstreig%effmass(in) = gm*r(j)*r(j)
66 call eglist(hecmesh, hecmat, fstreig)
70 write(
imsg,*)
'*----------------------------------------------*'
71 write(
imsg,*)
'*--E I G E N V A L U E C O N V E R G E N C E--*'
72 write(
imsg,*)
'Absolute residual = |(||Kx - lambda*Mx||)|'
73 write(
imsg,*)
' Iter.# Eigenvalue Abs. Residual '
74 write(
ilog,*)
' Iter.# Eigenvalue Abs. Residual '
75 write(
imsg,*)
' *-----* *---------* *--------------*'
78 do j = 1, fstreig%nget
86 s(i) = mass(i)*eigvec(i,j)
91 chk = chk + (t(i) - (eigval(j)-fstreig%sigma)*s(i))**2
93 call hecmw_allreduce_r1(hecmesh, chk, hecmw_sum)
97 write(
imsg,
'(2x,i5,2x,1p5e15.6)') j, eigval(j), chk
98 write(
ilog,
'(i5,1p5e12.4)') j, eigval(j), chk
103 write(
imsg,*)
'* ---END Eigenvalue listing--- *'
120 type(hecmwst_result_data) :: fstrRESULT
122 integer(kind=kint) :: i, istep, nget, NP, NDOF, NPNDOF, totalmpc, MPC_METHOD
123 real(kind=
kreal) :: t1
124 real(kind=
kreal),
pointer :: eigvec(:,:)
125 real(kind=
kreal),
allocatable :: x(:), egval(:)
126 character(len=HECMW_HEADER_LEN) :: header
127 character(len=HECMW_MSG_LEN) :: comment
128 character(len=HECMW_NAME_LEN) :: label
129 character(len=HECMW_NAME_LEN) :: nameID
133 npndof = hecmat%NP*hecmat%NDOF
137 eigvec => fstreig%eigvec
144 egval(1) = fstreig%eigval(istep)
146 x(i) = eigvec(i,istep)
160 call hecmw_update_r(hecmesh, x, hecmat%NP, ndof)
163 header =
"*fstrresult"
164 comment =
"eigen_result"
165 call hecmw_result_init(hecmesh,istep,header,comment)
167 call hecmw_result_add(hecmw_result_dtype_global,1,label,egval)
168 label =
"DISPLACEMENT"
169 call hecmw_result_add(hecmw_result_dtype_node,ndof,label,x)
171 call hecmw_result_write_by_name(nameid)
172 call hecmw_result_finalize
176 call hecmw_nullify_result_data(fstrresult)
177 fstrresult%ng_component = 1
178 fstrresult%nn_component = 1
179 fstrresult%ne_component = 0
180 allocate(fstrresult%ng_dof(1))
181 allocate(fstrresult%global_label(1))
182 allocate(fstrresult%global_val_item(1))
183 fstrresult%ng_dof(1) = 1
184 fstrresult%global_label(1) =
'EIGENVALUE'
185 fstrresult%global_val_item(1) = egval(1)
186 allocate(fstrresult%nn_dof(1))
187 allocate(fstrresult%node_label(1))
188 allocate(fstrresult%node_val_item(ndof*hecmat%NP))
189 fstrresult%nn_dof(1) = ndof
190 fstrresult%node_label(1) =
'DISPLACEMENT'
191 fstrresult%node_val_item = x
193 call hecmw_visualize_init
194 call hecmw_visualize( hecmesh, fstrresult, istep )
195 call hecmw_visualize_finalize
197 call hecmw_result_free(fstrresult)
206 subroutine eglist(hecMESH, hecMAT, fstrEIG)
215 integer(kind=kint) :: NDOF
216 integer(kind=kint) :: i, j, in, iter ,nget
217 real(kind=
kreal) :: pi, angle, freq, pf(3), em(3)
218 real(kind=
kreal),
pointer :: eigval(:)
223 pi = 4.0d0 * datan(1.0d0)
224 eigval => fstreig%eigval
228 write(
ilog,
"(a)")
"********************************"
229 write(
ilog,
"(a)")
"*RESULT OF EIGEN VALUE ANALYSIS*"
230 write(
ilog,
"(a)")
"********************************"
232 write(
ilog,
"(a,i8)")
"NUMBER OF ITERATIONS = ",iter
233 write(
ilog,
"(a,1pe12.4)")
"TOTAL MASS = ",fstreig%totalmass
235 write(
ilog,
"(3a)")
" ANGLE FREQUENCY ",&
236 "PARTICIPATION FACTOR EFFECTIVE MASS"
237 write(
ilog,
"(3a)")
" NO. EIGENVALUE FREQUENCY (HZ) ",&
239 write(
ilog,
"(3a)")
" --- ---------- ---------- ---------- ",&
240 "---------- ---------- ---------- ---------- ---------- ----------"
242 write(*,
"(a)")
"#----------------------------------#"
243 write(*,
"(a)")
"# RESULT OF EIGEN VALUE ANALYSIS #"
244 write(*,
"(a)")
"#----------------------------------#"
246 write(*,
"(a,i8)")
"### NUMBER OF ITERATIONS = ",iter
247 write(*,
"(a,1pe12.4)")
"### TOTAL MASS = ",fstreig%totalmass
249 write(*,
"(3a)")
" PERIOD FREQUENCY ",&
250 "PARTICIPATION FACTOR EFFECTIVE MASS"
251 write(*,
"(3a)")
" NO. [Sec] [HZ] ",&
253 write(*,
"(3a)")
" --- --------- --------- ",&
254 "--------- --------- --------- --------- --------- ---------"
259 if(eigval(i) < 0.0d0) eigval(i) = 0.0d0
260 angle = dsqrt(eigval(i))
261 freq = angle*0.5d0/pi
265 do j = 1, min(ndof, 3)
266 pf(j) = fstreig%partfactor(ndof*(i-1) + j)
267 em(j) = fstreig%effmass(ndof*(i-1) + j)
270 write(
ilog,
'(I5,1P9E12.4)') in, eigval(i), angle, freq, pf(1), pf(2), pf(3), em(1), em(2), em(3)
271 write(* ,
'(I5,1P8E11.3)') in, 1.0d0/freq, freq, pf(1), pf(2), pf(3), em(1), em(2), em(3)