FrontISTR 5.2.0
Large-scale structural analysis program with finit element method
Loading...
Searching...
No Matches
dynamic_output.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!-------------------------------------------------------------------------------
7 implicit none
8contains
9
11 !----------------------------------------------------------------------*
12 subroutine fstr_dynamic_output( hecMESH, fstrSOLID, fstrDYNAMIC, fstrPARAM )
13 !----------------------------------------------------------------------*
14 use m_fstr
18 type(hecmwst_local_mesh), intent(in) :: hecMESH
19 type(fstr_solid), intent(inout) :: fstrSOLID
20 type(fstr_dynamic), intent(in) :: fstrDYNAMIC
21 type(fstr_param), intent(in) :: fstrPARAM
22
23 type(hecmwst_result_data) :: fstrRESULT
24 integer(kind=kint) :: i, j, ndof, maxstep, interval, fnum, is, iE, gid, istep
25
26 ndof = hecmesh%n_dof
27
28 !C-- SET DISPLACEMENT etc.
29 istep = fstrdynamic%i_step
30
31 if( fstrsolid%TEMP_ngrp_tot>0 .or. fstrsolid%TEMP_irres>0 ) then
32 if( ndof==3 ) then
33 allocate( fstrsolid%tnstrain(hecmesh%n_node*6) )
34 allocate( fstrsolid%testrain(hecmesh%n_elem*6) )
35 else if( ndof==2 ) then
36 allocate( fstrsolid%tnstrain(hecmesh%n_node*3) )
37 allocate( fstrsolid%testrain(hecmesh%n_elem*3) )
38 else if( ndof== 4) then
39 write(*,*)'Error: This routine is not implemented'
40 stop
41 endif
42 endif
43
44 if( ndof==3 .or. ndof == 4 ) then
45 call fstr_nodalstress3d( hecmesh, fstrsolid )
46 else if( ndof==2 ) then
47 call fstr_nodalstress2d( hecmesh, fstrsolid )
48 else if( ndof==6 ) then
49 call fstr_nodalstress6d( hecmesh, fstrsolid )
50 endif
51
52 maxstep = fstrdynamic%n_step
53
54 if( (mod(istep,fstrsolid%output_ctrl(1)%freqency)==0 .or. istep==maxstep) ) then
55 fnum = fstrsolid%output_ctrl(1)%filenum
56 call fstr_dynamic_post( fnum, hecmesh, fstrsolid, fstrdynamic )
57 endif
58
59 if( fstrsolid%output_ctrl(2)%outinfo%grp_id>0 .and. &
60 (mod(istep,fstrsolid%output_ctrl(2)%freqency)==0 .or. istep==maxstep) ) then
61 is = fstrsolid%output_ctrl(2)%outinfo%grp_id
62 fnum = fstrsolid%output_ctrl(2)%filenum
63 do i = hecmesh%node_group%grp_index(is-1)+1, hecmesh%node_group%grp_index(is)
64 ie = hecmesh%node_group%grp_item(i)
65 gid = hecmesh%global_node_ID(ie)
66 write(fnum,'(2i10,1p6e15.7)') istep,gid,(fstrsolid%unode(ndof*(ie-1)+j),j=1,ndof)
67 enddo
68 endif
69
70 if( iresult==1 .and. &
71 (mod(istep,fstrsolid%output_ctrl(3)%freqency)==0 .or. istep==maxstep) ) then
72 call fstr_write_result( hecmesh, fstrsolid, fstrparam, istep, fstrdynamic%t_curr, 0, fstrdynamic )
73 endif
74
75 if( ivisual==1 .and. &
76 (mod(istep,fstrsolid%output_ctrl(4)%freqency)==0 .or. istep==maxstep) ) then
77 call fstr_make_result( hecmesh, fstrsolid, fstrresult, istep, fstrdynamic%t_curr, fstrdynamic )
78 call fstr2hecmw_mesh_conv( hecmesh )
79 call hecmw_visualize_init
80 call hecmw_visualize( hecmesh, fstrresult, istep )
81 call hecmw_visualize_finalize
82 call hecmw2fstr_mesh_conv( hecmesh )
83 call hecmw_result_free( fstrresult )
84 endif
85
86 end subroutine fstr_dynamic_output
87
89 !----------------------------------------------------------------------*
90 subroutine fstr_dynamic_post( fnum, hecMESH, fstrSOLID, fstrDYNAMIC )
91 !----------------------------------------------------------------------*
92 use m_fstr
93 integer, intent(in) :: fnum
94 type(hecmwst_local_mesh), intent(in) :: hecMESH
95 type(fstr_solid), intent(in) :: fstrSOLID
96 type(fstr_dynamic), intent(in) :: fstrDYNAMIC
97
98 real(kind=kreal) :: umax(6), umin(6), vmax(6), vmin(6), amax(6), amin(6)
99 real(kind=kreal) :: emax(14), emin(14), smax(14), smin(14)
100 real(kind=kreal) :: mmax(1), mmin(1), emmax(1), emmin(1)
101 real(kind=kreal) :: eemax(14), eemin(14), esmax(14), esmin(14)
102
103 real(kind=kreal) :: gumax(6), gumin(6), gvmax(6), gvmin(6), gamax(6), gamin(6)
104 real(kind=kreal) :: gemax(14), gemin(14), gsmax(14), gsmin(14)
105 real(kind=kreal) :: gmmax(1), gmmin(1), gemmax(1), gemmin(1)
106 real(kind=kreal) :: geemax(14), geemin(14), gesmax(14), gesmin(14)
107
108 integer(kind=kint) :: IUmax(6), IUmin(6), IVmax(6), IVmin(6), IAmax(6), IAmin(6)
109 integer(kind=kint) :: IEmax(14), IEmin(14), ISmax(14), ISmin(14)
110 integer(kind=kint) :: IMmax(1), IMmin(1), IEMmax(1), IEMmin(1)
111 integer(kind=kint) :: IEEmax(14), IEEmin(14), IESmax(14), IESmin(14)
112 integer(kind=kint) :: i, j, k, ndof, mdof, ID_area, idx
113 integer(kind=kint) :: label(6)
114
115 if( fstrdynamic%i_step==0 ) return
116 if( fstrdynamic%idx_eqa==1 .and. fstrdynamic%i_step>0 ) then
117 idx = 2
118 else
119 idx = 1
120 endif
121 ndof = hecmesh%n_dof
122
123 write( fnum, '(''#### Result step='',I6)') fstrdynamic%i_step
124 select case (ndof)
125 case (2)
126 mdof = 3
127 label(1)=11;label(2)=22
128 label(3)=12
129 case (3,4,6)
130 mdof = 6
131 label(1)=11;label(2)=22;label(3)=33;
132 label(4)=12;label(5)=23;label(6)=31;
133 end select
134
135 j = hecmesh%global_node_ID(1)
136 do k = 1, ndof
137 umax(k) = fstrdynamic%DISP(k,idx); umin(k) = fstrdynamic%DISP(k,idx)
138 vmax(k) = fstrdynamic%VEL(k,idx); vmin(k) = fstrdynamic%VEL(k,idx)
139 amax(k) = fstrdynamic%ACC(k,idx); amin(k) = fstrdynamic%ACC(k,idx)
140 iumax(k)= j; iumin(k)= j; ivmax(k)= j; ivmin(k)= j; iamax(k)= j; iamin(k)= j
141 enddo
142 do k = 1, mdof
143 emax(k) = fstrsolid%STRAIN(k); emin(k) = fstrsolid%STRAIN(k)
144 smax(k) = fstrsolid%STRESS(k); smin(k) = fstrsolid%STRESS(k)
145 eemax(k) = fstrsolid%ESTRAIN(k); eemin(k) = fstrsolid%ESTRAIN(k)
146 esmax(k) = fstrsolid%ESTRESS(k); esmin(k) = fstrsolid%ESTRESS(k)
147 iemax(k)= j; iemin(k)= j; ismax(k)= j; ismin(k)= j
148 ieemax(k)= j; ieemin(k)= j; iesmax(k)= j; iesmin(k)= j
149 enddo
150 mmax(1) = fstrsolid%MISES(1); mmin(1) = fstrsolid%MISES(1)
151 emmax(1) = fstrsolid%EMISES(1); emmin(1) = fstrsolid%EMISES(1)
152 immax(1)= j; immin(1)= j; iemmax(1)= j; iemmin(1)= j
153
154 !C*** Show Displacement / Velosity / Acc
155 do i = 1, hecmesh%nn_internal
156 if(fstrsolid%is_rot(i)==1)cycle
157 j = hecmesh%global_node_ID(i)
158 do k = 1, ndof
159 if ( fstrdynamic%DISP(ndof*(i-1)+k,idx) > umax(k) ) then
160 umax(k) = fstrdynamic%DISP(ndof*(i-1)+k,idx)
161 iumax(k)= j
162 else if( fstrdynamic%DISP(ndof*(i-1)+k,idx) < umin(k) ) then
163 umin(k) = fstrdynamic%DISP(ndof*(i-1)+k,idx)
164 iumin(k)= j
165 endif
166 if ( fstrdynamic%VEL(ndof*(i-1)+k,idx) > vmax(k) ) then
167 vmax(k) = fstrdynamic%VEL(ndof*(i-1)+k,idx)
168 ivmax(k)= j
169 else if( fstrdynamic%VEL(ndof*(i-1)+k,idx) < vmin(k) ) then
170 vmin(k) = fstrdynamic%VEL(ndof*(i-1)+k,idx)
171 ivmin(k)= j
172 endif
173 if ( fstrdynamic%ACC(ndof*(i-1)+k,idx) > amax(k) ) then
174 amax(k) = fstrdynamic%ACC(ndof*(i-1)+k,idx)
175 iamax(k)= j
176 else if( fstrdynamic%ACC(ndof*(i-1)+k,idx) < amin(k) ) then
177 amin(k) = fstrdynamic%ACC(ndof*(i-1)+k,idx)
178 iamin(k)= j
179 endif
180 enddo
181 enddo
182 !C*** Nodal Strain / Stress / MISES
183 !C @node
184 do i = 1, hecmesh%nn_internal
185 if(fstrsolid%is_rot(i)==1)cycle
186 j = hecmesh%global_node_ID(i)
187 do k = 1, mdof
188 if ( fstrsolid%STRAIN(mdof*(i-1)+k) > emax(k) ) then
189 emax(k) = fstrsolid%STRAIN(mdof*(i-1)+k)
190 iemax(k)= j
191 else if( fstrsolid%STRAIN(mdof*(i-1)+k) < emin(k) ) then
192 emin(k) = fstrsolid%STRAIN(mdof*(i-1)+k)
193 iemin(k)= j
194 endif
195 if ( fstrsolid%STRESS(mdof*(i-1)+k) > smax(k) ) then
196 smax(k) = fstrsolid%STRESS(mdof*(i-1)+k)
197 ismax(k)= j
198 else if( fstrsolid%STRESS(mdof*(i-1)+k) < smin(k) ) then
199 smin(k) = fstrsolid%STRESS(mdof*(i-1)+k)
200 ismin(k)= j
201 endif
202 enddo
203 if ( fstrsolid%MISES(i) > mmax(1) ) then
204 mmax(1) = fstrsolid%MISES(i)
205 immax(1)= j
206 else if( fstrsolid%MISES(i) < mmin(1) ) then
207 mmin(1) = fstrsolid%MISES(i)
208 immin(1)= j
209 endif
210 enddo
211 !C*** Elemental Strain / STRESS
212 do i = 1, hecmesh%n_elem
213 id_area = hecmesh%elem_ID(i*2)
214 if( id_area==hecmesh%my_rank ) then
215 j = hecmesh%global_elem_ID(i)
216 do k = 1, mdof
217 if( fstrsolid%ESTRAIN(mdof*(i-1)+k) > eemax(k) ) then
218 eemax(k) = fstrsolid%ESTRAIN(mdof*(i-1)+k)
219 ieemax(k)= j
220 else if( fstrsolid%ESTRAIN(mdof*(i-1)+k) < eemin(k) ) then
221 eemin(k) = fstrsolid%ESTRAIN(mdof*(i-1)+k)
222 ieemin(k)= j
223 endif
224 if( fstrsolid%ESTRESS(mdof*(i-1)+k) > esmax(k) ) then
225 esmax(k) = fstrsolid%ESTRESS(mdof*(i-1)+k)
226 iesmax(k)= j
227 else if( fstrsolid%ESTRESS(mdof*(i-1)+k) < esmin(k) ) then
228 esmin(k) = fstrsolid%ESTRESS(mdof*(i-1)+k)
229 iesmin(k)= j
230 endif
231 enddo
232 if( fstrsolid%EMISES(i) > emmax(1) ) then
233 emmax(1) = fstrsolid%EMISES(i)
234 iemmax(1)= j
235 else if( fstrsolid%EMISES(i) < emmin(1) ) then
236 emmin(1) = fstrsolid%EMISES(i)
237 iemmin(1)= j
238 endif
239 endif
240 enddo
241
242
243 write(ilog,*) '##### Local Summary @Node :Max/IdMax/Min/IdMin####'
244 do i = 1, ndof; write(ilog,1029) ' //U',i, ' ',umax(i),iumax(i),umin(i),iumin(i); end do
245 if (ndof == 4) write(ilog,1009) ' //P ' , umax(4),iumax(4),umin(4),iumin(4)
246 do i = 1, ndof; write(ilog,1029) ' //V',i, ' ',vmax(i),ivmax(i),vmin(i),ivmin(i); end do
247 do i = 1, ndof; write(ilog,1029) ' //A',i, ' ',amax(i),iamax(i),amin(i),iamin(i); end do
248 do i = 1, mdof; write(ilog,1029) ' //E',label(i),' ',emax(i),iemax(i),emin(i),iemin(i); end do
249 do i = 1, mdof; write(ilog,1029) ' //S',label(i),' ',smax(i),ismax(i),smin(i),ismin(i); end do
250 write(ilog,1009) '//SMS ' ,mmax(1),immax(1),mmin(1),immin(1)
251 write(ilog,*) '##### Local Summary @Element :Max/IdMax/Min/IdMin####'
252 do i = 1, mdof; write(ilog,1029) ' //E',label(i),' ',eemax(i),ieemax(i),eemin(i),ieemin(i); end do
253 do i = 1, mdof; write(ilog,1029) ' //S',label(i),' ',esmax(i),iesmax(i),esmin(i),iesmin(i); end do
254 write(ilog,1009) '//SMS ' ,emmax(1),iemmax(1),emmin(1),iemmin(1)
255
256 !C*** Show Summary
257 gumax = umax; gumin = umin; gvmax = vmax; gvmin = vmin; gamax = amax; gamin = amin;
258 gemax = emax; gemin = emin; geemax = eemax; geemin = eemin;
259 gsmax = smax; gsmin = smin; gesmax = esmax; gesmin = esmin;
260 gmmax = mmax; gmmin = mmin; gemmax = emmax; gemmin = emmin;
261
262 call hecmw_allreduce_r(hecmesh,gumax,ndof,hecmw_max)
263 call hecmw_allreduce_r(hecmesh,gumin,ndof,hecmw_min)
264 call hecmw_allreduce_r(hecmesh,gvmax,ndof,hecmw_max)
265 call hecmw_allreduce_r(hecmesh,gvmin,ndof,hecmw_min)
266 call hecmw_allreduce_r(hecmesh,gamax,ndof,hecmw_max)
267 call hecmw_allreduce_r(hecmesh,gamin,ndof,hecmw_min)
268 call hecmw_allreduce_r(hecmesh,gemax,mdof,hecmw_max)
269 call hecmw_allreduce_r(hecmesh,gemin,mdof,hecmw_min)
270 call hecmw_allreduce_r(hecmesh,gsmax,mdof,hecmw_max)
271 call hecmw_allreduce_r(hecmesh,gsmin,mdof,hecmw_min)
272 call hecmw_allreduce_r(hecmesh,gmmax,1,hecmw_max)
273 call hecmw_allreduce_r(hecmesh,gmmin,1,hecmw_min)
274 call hecmw_allreduce_r(hecmesh,geemax,mdof,hecmw_max)
275 call hecmw_allreduce_r(hecmesh,geemin,mdof,hecmw_min)
276 call hecmw_allreduce_r(hecmesh,gesmax,mdof,hecmw_max)
277 call hecmw_allreduce_r(hecmesh,gesmin,mdof,hecmw_min)
278 call hecmw_allreduce_r(hecmesh,gemmax,1,hecmw_max)
279 call hecmw_allreduce_r(hecmesh,gemmin,1,hecmw_min)
280
281 do i=1,ndof
282 if(gumax(i) > umax(i)) iumax(i) = 0
283 if(gvmax(i) > vmax(i)) ivmax(i) = 0
284 if(gamax(i) > amax(i)) iamax(i) = 0
285 if(gumin(i) < umin(i)) iumin(i) = 0
286 if(gvmin(i) < vmin(i)) ivmin(i) = 0
287 if(gamin(i) < amin(i)) iamin(i) = 0
288 enddo
289 do i=1,mdof
290 if(gemax(i) > emax(i)) iemax(i) = 0
291 if(gsmax(i) > smax(i)) ismax(i) = 0
292 if(geemax(i) > eemax(i)) ieemax(i) = 0
293 if(gesmax(i) > esmax(i)) iesmax(i) = 0
294 if(gemin(i) < emin(i)) iemin(i) = 0
295 if(gsmin(i) < smin(i)) ismin(i) = 0
296 if(geemin(i) < eemin(i)) ieemin(i) = 0
297 if(gesmin(i) < esmin(i)) iesmin(i) = 0
298 enddo
299 do i=1,1
300 if(gmmax(i) > mmax(i)) immax(i) = 0
301 if(gemmax(i) > emmax(i)) iemmax(i) = 0
302 if(gmmin(i) < mmin(i)) immin(i) = 0
303 if(gemmin(i) < emmin(i)) iemmin(i) = 0
304 enddo
305
306 call hecmw_allreduce_i(hecmesh,iumax,ndof,hecmw_max)
307 call hecmw_allreduce_i(hecmesh,iumin,ndof,hecmw_max)
308 call hecmw_allreduce_i(hecmesh,ivmax,ndof,hecmw_max)
309 call hecmw_allreduce_i(hecmesh,ivmin,ndof,hecmw_max)
310 call hecmw_allreduce_i(hecmesh,iamax,ndof,hecmw_max)
311 call hecmw_allreduce_i(hecmesh,iamin,ndof,hecmw_max)
312 call hecmw_allreduce_i(hecmesh,iemax,mdof,hecmw_max)
313 call hecmw_allreduce_i(hecmesh,iemin,mdof,hecmw_max)
314 call hecmw_allreduce_i(hecmesh,ismax,mdof,hecmw_max)
315 call hecmw_allreduce_i(hecmesh,ismin,mdof,hecmw_max)
316 call hecmw_allreduce_i(hecmesh,immax,1,hecmw_max)
317 call hecmw_allreduce_i(hecmesh,immin,1,hecmw_max)
318 call hecmw_allreduce_i(hecmesh,ieemax,mdof,hecmw_max)
319 call hecmw_allreduce_i(hecmesh,ieemin,mdof,hecmw_max)
320 call hecmw_allreduce_i(hecmesh,iesmax,mdof,hecmw_max)
321 call hecmw_allreduce_i(hecmesh,iesmin,mdof,hecmw_max)
322 call hecmw_allreduce_i(hecmesh,iemmax,1,hecmw_max)
323 call hecmw_allreduce_i(hecmesh,iemmin,1,hecmw_max)
324
325 if( hecmesh%my_rank==0 ) then
326 write(ilog,*) '##### Global Summary @Node :Max/IdMax/Min/IdMin####'
327 do i = 1, ndof; write(ilog,1029) ' //U',i, ' ',gumax(i),iumax(i),gumin(i),iumin(i); end do
328 if (ndof == 4) write(ilog,1009) ' //P ' ,gumax(4),iumax(4),gumin(4),iumin(4)
329 do i = 1, ndof; write(ilog,1029) ' //V',i, ' ',gvmax(i),ivmax(i),gvmin(i),ivmin(i); end do
330 do i = 1, ndof; write(ilog,1029) ' //A',i, ' ',gamax(i),iamax(i),gamin(i),iamin(i); end do
331 do i = 1, mdof; write(ilog,1029) ' //E',label(i),' ',gemax(i),iemax(i),gemin(i),iemin(i); end do
332 do i = 1, mdof; write(ilog,1029) ' //S',label(i),' ',gsmax(i),ismax(i),gsmin(i),ismin(i); end do
333 write(ilog,1009) '//SMS ' ,gmmax(1),immax(1),gmmin(1),immin(1)
334 write(ilog,*) '##### Global Summary @Element :Max/IdMax/Min/IdMin####'
335 do i = 1, mdof; write(ilog,1029) ' //E',label(i),' ',geemax(i),ieemax(i),geemin(i),ieemin(i); end do
336 do i = 1, mdof; write(ilog,1029) ' //S',label(i),' ',gesmax(i),iesmax(i),gesmin(i),iesmin(i); end do
337 write(ilog,1009) '//SMS ' ,gemmax(1),iemmax(1),gemmin(1),iemmin(1)
338 endif
339
340 1009 format(a7,1pe12.4,i10,1pe12.4,i10)
341 1029 format(a,i0,a,1pe12.4,i10,1pe12.4,i10)
342 end subroutine fstr_dynamic_post
343
344 !C================================================================C
345 !C-- subroutine dynamic_output_monit
346 !C================================================================C
347 subroutine dynamic_output_monit(hecMESH, fstrPARAM, fstrDYNAMIC, fstrEIG, fstrSOLID)
348 use m_fstr
349 type(hecmwst_local_mesh) :: hecMESH
350 type(fstr_param) :: fstrPARAM
351 type(fstr_dynamic) :: fstrDYNAMIC
352 type(fstr_eigen) :: fstrEIG
353 type(fstr_solid) :: fstrSOLID
354
355 integer(kind=kint) :: idx, ii, jj, ierr, ncmp
356 integer(kind=kint) :: num_monit, ig, is, iE, ik, iunitS, iunit
357 logical :: yes
358
359 if( mod(fstrdynamic%i_step,fstrdynamic%nout_monit)/=0 ) return
360
361 if( fstrdynamic%idx_eqa==1 .and. fstrdynamic%i_step>0 ) then
362 idx = 2
363 else
364 idx = 1
365 endif
366
367 num_monit = 0
368 ig = fstrdynamic%ngrp_monit
369 is = hecmesh%node_group%grp_index(ig-1)+1
370 ie = hecmesh%node_group%grp_index(ig)
371 do ik=is,ie
372 ii = hecmesh%node_group%grp_item(ik)
373 if (ii > hecmesh%nn_internal) cycle
374 num_monit = num_monit+1
375 jj = hecmesh%global_node_id(ii)
376 iunits = 10*(num_monit-1)
377
378 !C-- displacement
379 if( fstrdynamic%iout_list(1)==1 ) then
380 iunit = iunits + fstrdynamic%dynamic_IW4
381 write( iunit, '(i10,1pe13.4e3,i10,1p6e13.4e3)' ) &
382 fstrdynamic%i_step, fstrdynamic%t_curr, jj, &
383 fstrdynamic%DISP( hecmesh%n_dof*(ii-1)+1 : hecmesh%n_dof*ii , idx )
384 end if
385 !C-- velocity
386 if( fstrdynamic%iout_list(2)==1 ) then
387 iunit = iunits + fstrdynamic%dynamic_IW5
388 write( iunit, '(i10,1pe13.4e3,i10,1p6e13.4e3)' ) &
389 fstrdynamic%i_step, fstrdynamic%t_curr, jj, &
390 fstrdynamic%VEL( hecmesh%n_dof*(ii-1)+1 : hecmesh%n_dof*ii , idx )
391 end if
392 !C-- acceleration
393 if( fstrdynamic%iout_list(3)==1 ) then
394 iunit = iunits + fstrdynamic%dynamic_IW6
395 write( iunit, '(i10,1pe13.4e3,i10,1p6e13.4e3)' ) &
396 fstrdynamic%i_step, fstrdynamic%t_curr, jj, &
397 fstrdynamic%ACC( hecmesh%n_dof*(ii-1)+1 : hecmesh%n_dof*ii , idx )
398 end if
399 !C-- nodal force
400 if( fstrdynamic%iout_list(4)==1 ) then
401 iunit = iunits + fstrdynamic%dynamic_IW7
402 write( iunit, '(i10,1pe13.4e3,i10,1p6e13.4e3)' ) &
403 fstrdynamic%i_step, fstrdynamic%t_curr, jj, &
404 fstrsolid%QFORCE( hecmesh%n_dof*(ii-1)+1 : hecmesh%n_dof*ii )
405 end if
406 !C-- strain
407 if( fstrdynamic%iout_list(5) > 0 ) then
408 if (hecmesh%n_dof == 2 .or. hecmesh%n_dof == 3 .or. hecmesh%n_dof == 4) then
409 ncmp = 6
410 else
411 ncmp = 12
412 endif
413 iunit = iunits + fstrdynamic%dynamic_IW8
414 write( iunit, '(i10,1pe13.4e3,i10,1p6e13.4e3)') &
415 fstrdynamic%i_step, fstrdynamic%t_curr, jj, &
416 fstrsolid%STRAIN( ncmp*(ii-1)+1 : ncmp*ii )
417 end if
418 !C-- stress
419 if( fstrdynamic%iout_list(6) > 0 ) then
420 if (hecmesh%n_dof == 2 .or. hecmesh%n_dof == 3 .or. hecmesh%n_dof == 4) then
421 ncmp = 6
422 else
423 ncmp = 12
424 endif
425 iunit = iunits + fstrdynamic%dynamic_IW9
426 write( iunit, '(i10,1pe13.4e3,i10,1p7e13.4e3)') &
427 fstrdynamic%i_step, fstrdynamic%t_curr, jj, &
428 fstrsolid%STRESS( ncmp*(ii-1)+1 : ncmp*ii )
429 end if
430 enddo
431
432 if( hecmesh%my_rank==0 ) then
433 if( any(fstrdynamic%iout_list(1:3)==1) ) then
434 inquire( file='dyna_energy.txt', opened=yes )
435 if( .not. yes ) then
436 open( fstrdynamic%dynamic_IW10, file='dyna_energy.txt', status='replace', iostat=ierr )
437 if( ierr/=0 ) then
438 write(*,*) 'stop due to file opening error <dyna_enrgy.txt>'
439 call hecmw_abort( hecmw_comm_get_comm() )
440 endif
441 write( fstrdynamic%dynamic_IW10, * ) &
442 ' time step', ' time ', ' kinetic energy', ' strain energy', ' total energy'
443 endif
444 if(fstrdynamic%i_step==0) then
445 fstrdynamic%kineticEnergy = 0.0d0
446 do ii = 1, hecmesh%n_node*hecmesh%n_dof
447 fstrdynamic%kineticEnergy = fstrdynamic%kineticEnergy &
448 + 0.5d0 * fstreig%mass(ii) * fstrdynamic%VEL(ii,idx) * fstrdynamic%VEL(ii,idx)
449 enddo
450 endif
451 fstrdynamic%totalEnergy = fstrdynamic%kineticEnergy + fstrdynamic%strainEnergy
452 write( fstrdynamic%dynamic_IW10, '(i10,1pe13.4e3,1p3e16.4e3)' ) &
453 fstrdynamic%i_step, fstrdynamic%t_curr, &
454 fstrdynamic%kineticEnergy, fstrdynamic%strainEnergy, fstrdynamic%totalEnergy
455 endif
456 if( fstrdynamic%i_step==fstrdynamic%n_step ) close(fstrdynamic%dynamic_IW10)
457 endif
458 end subroutine dynamic_output_monit
459
460 !C================================================================C
461 !C-- subroutine matvec
462 !C================================================================C
463 subroutine matvec(y,x,hecMAT,ndof,D,AU,AL)
464 use m_fstr
465 type(hecmwst_matrix) :: hecMAT
466
467 integer(kind=kint) :: ndof, i, is, ie, j, icol, ki, kj, ix, iy, ip, nn
468 real(kind=kreal) :: d(ndof*ndof*hecmat%NP)
469 real(kind=kreal) :: au(ndof*ndof*hecmat%NPU)
470 real(kind=kreal) :: al(ndof*ndof*hecmat%NPL)
471 real(kind=kreal) :: x(ndof*hecmat%NP)
472 real(kind=kreal) :: y(ndof*hecmat%NP)
473
474 nn=ndof*ndof
475
476 y=0.0d0
477
478 do i=1,hecmat%NP
479 is=hecmat%indexU(i-1)+1
480 ie=hecmat%indexU(i)
481 do j=is,ie
482 icol=hecmat%itemU(j)
483 do ki=1,ndof
484 iy=ndof*(i-1)+ki
485 do kj=1,ndof
486 ix=ndof*(icol-1)+kj
487 ip=nn*(j-1)+ndof*(ki-1)+kj
488 y(iy)=y(iy)+au(ip)*x(ix)
489 enddo
490 enddo
491 enddo
492 enddo
493
494 do i=1,hecmat%NP
495 is=hecmat%indexL(i-1)+1
496 ie=hecmat%indexL(i)
497 do j=is,ie
498 icol=hecmat%itemL(j)
499 do ki=1,ndof
500 iy=ndof*(i-1)+ki
501 do kj=1,ndof
502 ix=ndof*(icol-1)+kj
503 ip=nn*(j-1)+ndof*(ki-1)+kj
504 y(iy)=y(iy)+al(ip)*x(ix)
505 enddo
506 enddo
507 enddo
508 enddo
509
510 do i=1,hecmat%NP
511 do ki=1,ndof
512 iy=ndof*(i-1)+ki
513 do kj=1,ndof
514 ix=ndof*(i-1)+kj
515 ip=nn*(i-1)+ndof*(ki-1)+kj
516 y(iy)=y(iy)+d(ip)*x(ix)
517 enddo
518 enddo
519 enddo
520 end subroutine matvec
521
522end module m_dynamic_output
This module provides functions to output result.
subroutine matvec(y, x, hecmat, ndof, d, au, al)
subroutine fstr_dynamic_post(fnum, hecmesh, fstrsolid, fstrdynamic)
Summarizer of output data which prints out max and min output values.
subroutine fstr_dynamic_output(hecmesh, fstrsolid, fstrdynamic, fstrparam)
Output result.
subroutine dynamic_output_monit(hecmesh, fstrparam, fstrdynamic, fstreig, fstrsolid)
This module provides functions to caluclation nodal stress.
subroutine fstr_nodalstress2d(hecmesh, fstrsolid)
Calculate NODAL STRESS of plane elements.
subroutine fstr_nodalstress6d(hecmesh, fstrsolid)
Calculate NODAL STRESS of shell elements.
subroutine fstr_nodalstress3d(hecmesh, fstrsolid)
Calculate NODAL STRESS of solid elements.
This module defined coomon data and basic structures for analysis.
Definition: m_fstr.f90:15
integer(kind=kint), pointer iresult
Definition: m_fstr.f90:106
integer(kind=kint), parameter ilog
FILE HANDLER.
Definition: m_fstr.f90:91
integer(kind=kint), pointer ivisual
Definition: m_fstr.f90:107
HECMW to FSTR Mesh Data Converter. Convering Conectivity of Element Type 232, 342 and 352.
subroutine hecmw2fstr_mesh_conv(hecmesh)
subroutine fstr2hecmw_mesh_conv(hecmesh)
This module provide a function to prepare output of static analysis.
Definition: make_result.f90:6
subroutine, public fstr_make_result(hecmesh, fstrsolid, fstrresult, istep, time, fstrdynamic)
MAKE RESULT for static and dynamic analysis (WITHOUT ELEMENTAL RESULTS) -----------------------------...
subroutine, public fstr_write_result(hecmesh, fstrsolid, fstrparam, istep, time, flag, fstrdynamic)
OUTPUT result file for static and dynamic analysis.
Definition: make_result.f90:23
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
Definition: m_fstr.f90:473
Package of data used by Lanczos eigenvalue solver.
Definition: m_fstr.f90:562
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.f90:138