FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
static_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
8 contains
9 
11  !----------------------------------------------------------------------*
12  subroutine fstr_static_output( cstep, istep, time, hecMESH, fstrSOLID, fstrPARAM, flag, outflag )
13  !----------------------------------------------------------------------*
14  use m_fstr
16  use m_make_result
18  use mcontact
19  integer, intent(in) :: cstep
20  integer, intent(in) :: istep
21  real(kind=kreal), intent(in) :: time
22  type (hecmwST_local_mesh), intent(in) :: hecMESH
23  type (fstr_solid), intent(inout) :: fstrSOLID
24  type (fstr_param ) :: fstrPARAM
25  integer, intent(in) :: flag
26  logical, intent(in) :: outflag
27 
28  type ( hecmwST_result_data ) :: fstrRESULT
29  integer(kind=kint) :: i, j, ndof, fnum, is, iE, gid
30 
31  ndof = hecmesh%n_dof
32 
33  if( fstrsolid%TEMP_ngrp_tot>0 .or. fstrsolid%TEMP_irres>0 ) then
34  if( ndof==3 ) then
35  if( .not. associated(fstrsolid%tnstrain) ) allocate( fstrsolid%tnstrain(hecmesh%n_node*6) )
36  if( .not. associated(fstrsolid%testrain) ) allocate( fstrsolid%testrain(hecmesh%n_elem*6) )
37  else if( ndof==2 ) then
38  if( .not. associated(fstrsolid%tnstrain) ) allocate( fstrsolid%tnstrain(hecmesh%n_node*3) )
39  if( .not. associated(fstrsolid%testrain) ) allocate( fstrsolid%testrain(hecmesh%n_elem*3) )
40  endif
41  endif
42 
43  if( ndof==3 ) then
44  call fstr_nodalstress3d( hecmesh, fstrsolid )
45  else if( ndof==2 ) then
46  call fstr_nodalstress2d( hecmesh, fstrsolid )
47  else if( ndof==6 ) then
48  call fstr_nodalstress6d( hecmesh, fstrsolid )
49  endif
50 
51  if( associated( fstrsolid%contacts ) ) then
52  call setup_contact_output_variables( hecmesh, fstrsolid, -1 )
53  endif
54 
55  if( flag==kststaticeigen ) then
56  if( iresult==1 .and. &
57  (mod(istep,fstrsolid%output_ctrl(3)%frequency)==0 .or. outflag) ) then
58  call fstr_write_result( hecmesh, fstrsolid, fstrparam, istep, time, 1 )
59  endif
60  return
61  endif
62 
63  if( iresult==1 .and. &
64  (mod(istep,fstrsolid%output_ctrl(3)%frequency)==0 .or. outflag) ) then
65  if( associated( fstrsolid%contacts ) ) &
66  & call setup_contact_output_variables( hecmesh, fstrsolid, 3 )
67  call fstr_write_result( hecmesh, fstrsolid, fstrparam, istep, time, 0 )
68  endif
69 
70  if( ivisual==1 .and. &
71  (mod(istep,fstrsolid%output_ctrl(4)%frequency)==0 .or. outflag) ) then
72 
73  if( associated( fstrsolid%contacts ) ) &
74  & call setup_contact_output_variables( hecmesh, fstrsolid, 4 )
75  call fstr_make_result( hecmesh, fstrsolid, fstrresult, istep, time )
76  call fstr2hecmw_mesh_conv( hecmesh )
77  call hecmw_visualize_init
78  call hecmw_visualize( hecmesh, fstrresult, istep )
79  call hecmw_visualize_finalize
80  call hecmw2fstr_mesh_conv( hecmesh )
81  call hecmw_result_free( fstrresult )
82  endif
83 
84  if( (mod(istep,fstrsolid%output_ctrl(1)%frequency)==0 .or. outflag) ) then
85  fnum = fstrsolid%output_ctrl(1)%filenum
86  call fstr_static_post( fnum, hecmesh, fstrsolid, istep )
87  endif
88 
89  if( fstrsolid%output_ctrl(2)%outinfo%grp_id>0 .and. &
90  (mod(istep,fstrsolid%output_ctrl(2)%frequency)==0 .or. outflag) ) then
91  is = fstrsolid%output_ctrl(2)%outinfo%grp_id
92  fnum = fstrsolid%output_ctrl(2)%filenum
93  do i = hecmesh%node_group%grp_index(is-1)+1, hecmesh%node_group%grp_index(is)
94  ie = hecmesh%node_group%grp_item(i)
95  gid = hecmesh%global_node_ID(ie)
96  write(fnum,'(2i10,1p6e15.7)') istep,gid,(fstrsolid%unode(ndof*(ie-1)+j),j=1,ndof)
97  enddo
98  endif
99 
100  end subroutine fstr_static_output
101 
103  !----------------------------------------------------------------------*
104  subroutine fstr_static_post( fnum, hecMESH, fstrSOLID, istep )
105  !----------------------------------------------------------------------*
106  use m_fstr
107  integer, intent(in) :: fnum, istep
108  type (hecmwST_local_mesh), intent(in) :: hecMESH
109  type (fstr_solid), intent(in) :: fstrSOLID
110 
111  real(kind=kreal) :: umax(6), umin(6), emax(14), emin(14), smax(14), smin(14)
112  real(kind=kreal) :: mmax(1), mmin(1), emmax(1), emmin(1)
113  real(kind=kreal) :: eemax(14), eemin(14), esmax(14), esmin(14)
114 
115  real(kind=kreal) :: gumax(6), gumin(6), gemax(14), gemin(14), gsmax(14), gsmin(14)
116  real(kind=kreal) :: gmmax(1), gmmin(1), gemmax(1), gemmin(1)
117  real(kind=kreal) :: geemax(14), geemin(14), gesmax(14), gesmin(14)
118 
119  integer(kind=kint) :: IUmax(6), IUmin(6), IEmax(14), IEmin(14), ISmax(14), ISmin(14)
120  integer(kind=kint) :: IMmax(1), IMmin(1), IEMmax(1), IEMmin(1)
121  integer(kind=kint) :: IEEmax(14), IEEmin(14), IESmax(14), IESmin(14)
122  integer(kind=kint) :: i, j, k, ndof, mdof, ID_area
123  integer(kind=kint) :: label(6)
124 
125  ndof = hecmesh%n_dof
126 
127  write( fnum, '(''#### Result step='',I6)') istep
128  select case (ndof)
129  case (2)
130  mdof = 3
131  label(1)=11;label(2)=22
132  label(3)=12
133  case (3,6)
134  mdof = 6
135  label(1)=11;label(2)=22;label(3)=33;
136  label(4)=12;label(5)=23;label(6)=31;
137  end select
138 
139  j = hecmesh%global_node_ID(1)
140  do k = 1, ndof
141  umax(k) = fstrsolid%unode(k); umin(k) = fstrsolid%unode(k)
142  iumax(k)= j; iumin(k)= j
143  enddo
144  do k = 1, mdof
145  emax(k) = fstrsolid%STRAIN(k); emin(k) = fstrsolid%STRAIN(k)
146  smax(k) = fstrsolid%STRESS(k); smin(k) = fstrsolid%STRESS(k)
147  eemax(k) = fstrsolid%ESTRAIN(k); eemin(k) = fstrsolid%ESTRAIN(k)
148  esmax(k) = fstrsolid%ESTRESS(k); esmin(k) = fstrsolid%ESTRESS(k)
149  iemax(k)= j; iemin(k)= j; ismax(k)= j; ismin(k)= j
150  ieemax(k)= j; ieemin(k)= j; iesmax(k)= j; iesmin(k)= j
151  enddo
152  mmax(1) = fstrsolid%MISES(1); mmin(1) = fstrsolid%MISES(1)
153  emmax(1) = fstrsolid%EMISES(1); emmin(1) = fstrsolid%EMISES(1)
154  immax(1)= j; immin(1)= j; iemmax(1)= j; iemmin(1)= j
155 
156  !C*** Show Displacement
157  do i = 1, hecmesh%nn_internal
158  if(fstrsolid%is_rot(i)==1)cycle
159  j = hecmesh%global_node_ID(i)
160  do k = 1, ndof
161  if ( fstrsolid%unode(ndof*(i-1)+k) > umax(k) ) then
162  umax(k) = fstrsolid%unode(ndof*(i-1)+k)
163  iumax(k)= j
164  else if( fstrsolid%unode(ndof*(i-1)+k) < umin(k) ) then
165  umin(k) = fstrsolid%unode(ndof*(i-1)+k)
166  iumin(k)= j
167  endif
168  enddo
169  enddo
170  !C*** Nodal Strain / Stress / MISES
171  !C @node
172  do i = 1, hecmesh%nn_internal
173  if(fstrsolid%is_rot(i)==1)cycle
174  j = hecmesh%global_node_ID(i)
175  do k = 1, mdof
176  if ( fstrsolid%STRAIN(mdof*(i-1)+k) > emax(k) ) then
177  emax(k) = fstrsolid%STRAIN(mdof*(i-1)+k)
178  iemax(k)= j
179  else if( fstrsolid%STRAIN(mdof*(i-1)+k) < emin(k) ) then
180  emin(k) = fstrsolid%STRAIN(mdof*(i-1)+k)
181  iemin(k)= j
182  endif
183  if ( fstrsolid%STRESS(mdof*(i-1)+k) > smax(k) ) then
184  smax(k) = fstrsolid%STRESS(mdof*(i-1)+k)
185  ismax(k)= j
186  else if( fstrsolid%STRESS(mdof*(i-1)+k) < smin(k) ) then
187  smin(k) = fstrsolid%STRESS(mdof*(i-1)+k)
188  ismin(k)= j
189  endif
190  enddo
191  if ( fstrsolid%MISES(i) > mmax(1) ) then
192  mmax(1) = fstrsolid%MISES(i)
193  immax(1)= j
194  else if( fstrsolid%MISES(i) < mmin(1) ) then
195  mmin(1) = fstrsolid%MISES(i)
196  immin(1)= j
197  endif
198  enddo
199  !C*** Elemental Strain / STRESS
200  do i = 1, hecmesh%n_elem
201  id_area = hecmesh%elem_ID(i*2)
202  if( id_area==hecmesh%my_rank ) then
203  j = hecmesh%global_elem_ID(i)
204  do k = 1, mdof
205  if( fstrsolid%ESTRAIN(mdof*(i-1)+k) > eemax(k) ) then
206  eemax(k) = fstrsolid%ESTRAIN(mdof*(i-1)+k)
207  ieemax(k)= j
208  else if( fstrsolid%ESTRAIN(mdof*(i-1)+k) < eemin(k) ) then
209  eemin(k) = fstrsolid%ESTRAIN(mdof*(i-1)+k)
210  ieemin(k)= j
211  endif
212  if( fstrsolid%ESTRESS(mdof*(i-1)+k) > esmax(k) ) then
213  esmax(k) = fstrsolid%ESTRESS(mdof*(i-1)+k)
214  iesmax(k)= j
215  else if( fstrsolid%ESTRESS(mdof*(i-1)+k) < esmin(k) ) then
216  esmin(k) = fstrsolid%ESTRESS(mdof*(i-1)+k)
217  iesmin(k)= j
218  endif
219  enddo
220  if( fstrsolid%EMISES(i) > emmax(1) ) then
221  emmax(1) = fstrsolid%EMISES(i)
222  iemmax(1)= j
223  else if( fstrsolid%EMISES(i) < emmin(1) ) then
224  emmin(1) = fstrsolid%EMISES(i)
225  iemmin(1)= j
226  endif
227  endif
228  enddo
229 
230 
231  write(ilog,*) '##### Local Summary @Node :Max/IdMax/Min/IdMin####'
232  do i = 1, ndof; write(ilog,1029) ' //U',i, ' ',umax(i),iumax(i),umin(i),iumin(i); end do
233  do i = 1, mdof; write(ilog,1029) ' //E',label(i),' ',emax(i),iemax(i),emin(i),iemin(i); end do
234  do i = 1, mdof; write(ilog,1029) ' //S',label(i),' ',smax(i),ismax(i),smin(i),ismin(i); end do
235  write(ilog,1009) '//SMS ' ,mmax(1),immax(1),mmin(1),immin(1)
236  write(ilog,*) '##### Local Summary @Element :Max/IdMax/Min/IdMin####'
237  do i = 1, mdof; write(ilog,1029) ' //E',label(i),' ',eemax(i),ieemax(i),eemin(i),ieemin(i); end do
238  do i = 1, mdof; write(ilog,1029) ' //S',label(i),' ',esmax(i),iesmax(i),esmin(i),iesmin(i); end do
239  write(ilog,1009) '//SMS ' ,emmax(1),iemmax(1),emmin(1),iemmin(1)
240 
241  !C*** Show Summary
242  gumax = umax; gumin = umin;
243  gemax = emax; gemin = emin; geemax = eemax; geemin = eemin;
244  gsmax = smax; gsmin = smin; gesmax = esmax; gesmin = esmin;
245  gmmax = mmax; gmmin = mmin; gemmax = emmax; gemmin = emmin;
246 
247  call hecmw_allreduce_r(hecmesh,gumax,ndof,hecmw_max)
248  call hecmw_allreduce_r(hecmesh,gumin,ndof,hecmw_min)
249  call hecmw_allreduce_r(hecmesh,gemax,mdof,hecmw_max)
250  call hecmw_allreduce_r(hecmesh,gemin,mdof,hecmw_min)
251  call hecmw_allreduce_r(hecmesh,gsmax,mdof,hecmw_max)
252  call hecmw_allreduce_r(hecmesh,gsmin,mdof,hecmw_min)
253  call hecmw_allreduce_r(hecmesh,gmmax,1,hecmw_max)
254  call hecmw_allreduce_r(hecmesh,gmmin,1,hecmw_min)
255  call hecmw_allreduce_r(hecmesh,geemax,mdof,hecmw_max)
256  call hecmw_allreduce_r(hecmesh,geemin,mdof,hecmw_min)
257  call hecmw_allreduce_r(hecmesh,gesmax,mdof,hecmw_max)
258  call hecmw_allreduce_r(hecmesh,gesmin,mdof,hecmw_min)
259  call hecmw_allreduce_r(hecmesh,gemmax,1,hecmw_max)
260  call hecmw_allreduce_r(hecmesh,gemmin,1,hecmw_min)
261 
262  do i=1,ndof
263  if(gumax(i) > umax(i)) iumax(i) = 0
264  if(gumin(i) < umin(i)) iumin(i) = 0
265  enddo
266  do i=1,mdof
267  if(gemax(i) > emax(i)) iemax(i) = 0
268  if(gsmax(i) > smax(i)) ismax(i) = 0
269  if(geemax(i) > eemax(i)) ieemax(i) = 0
270  if(gesmax(i) > esmax(i)) iesmax(i) = 0
271  if(gemin(i) < emin(i)) iemin(i) = 0
272  if(gsmin(i) < smin(i)) ismin(i) = 0
273  if(geemin(i) < eemin(i)) ieemin(i) = 0
274  if(gesmin(i) < esmin(i)) iesmin(i) = 0
275  enddo
276  do i=1,1
277  if(gmmax(i) > mmax(i)) immax(i) = 0
278  if(gemmax(i) > emmax(i)) iemmax(i) = 0
279  if(gmmin(i) < mmin(i)) immin(i) = 0
280  if(gemmin(i) < emmin(i)) iemmin(i) = 0
281  enddo
282 
283  call hecmw_allreduce_i(hecmesh,iumax,ndof,hecmw_max)
284  call hecmw_allreduce_i(hecmesh,iumin,ndof,hecmw_max)
285  call hecmw_allreduce_i(hecmesh,iemax,mdof,hecmw_max)
286  call hecmw_allreduce_i(hecmesh,iemin,mdof,hecmw_max)
287  call hecmw_allreduce_i(hecmesh,ismax,mdof,hecmw_max)
288  call hecmw_allreduce_i(hecmesh,ismin,mdof,hecmw_max)
289  call hecmw_allreduce_i(hecmesh,immax,1,hecmw_max)
290  call hecmw_allreduce_i(hecmesh,immin,1,hecmw_max)
291  call hecmw_allreduce_i(hecmesh,ieemax,mdof,hecmw_max)
292  call hecmw_allreduce_i(hecmesh,ieemin,mdof,hecmw_max)
293  call hecmw_allreduce_i(hecmesh,iesmax,mdof,hecmw_max)
294  call hecmw_allreduce_i(hecmesh,iesmin,mdof,hecmw_max)
295  call hecmw_allreduce_i(hecmesh,iemmax,1,hecmw_max)
296  call hecmw_allreduce_i(hecmesh,iemmin,1,hecmw_max)
297 
298  if( hecmesh%my_rank==0 ) then
299  write(ilog,*) '##### Global Summary @Node :Max/IdMax/Min/IdMin####'
300  do i = 1, ndof; write(ilog,1029) ' //U',i, ' ',gumax(i),iumax(i),gumin(i),iumin(i); end do
301  do i = 1, mdof; write(ilog,1029) ' //E',label(i),' ',gemax(i),iemax(i),gemin(i),iemin(i); end do
302  do i = 1, mdof; write(ilog,1029) ' //S',label(i),' ',gsmax(i),ismax(i),gsmin(i),ismin(i); end do
303  write(ilog,1009) '//SMS ' ,gmmax(1),immax(1),gmmin(1),immin(1)
304  write(ilog,*) '##### Global Summary @Element :Max/IdMax/Min/IdMin####'
305  do i = 1, mdof; write(ilog,1029) ' //E',label(i),' ',geemax(i),ieemax(i),geemin(i),ieemin(i); end do
306  do i = 1, mdof; write(ilog,1029) ' //S',label(i),' ',gesmax(i),iesmax(i),gesmin(i),iesmin(i); end do
307  write(ilog,1009) '//SMS ' ,gemmax(1),iemmax(1),gemmin(1),iemmin(1)
308  endif
309 
310  1009 format(a7,1pe12.4,i10,1pe12.4,i10)
311  1029 format(a,i0,a,1pe12.4,i10,1pe12.4,i10)
312  end subroutine fstr_static_post
313 
314 end module m_static_output
m_make_result::setup_contact_output_variables
subroutine, public setup_contact_output_variables(hecMESH, fstrSOLID, phase)
Definition: make_result.f90:1612
m_make_result::fstr_write_result
subroutine, public fstr_write_result(hecMESH, fstrSOLID, fstrPARAM, istep, time, flag, fstrDYNAMIC)
OUTPUT result file for static and dynamic analysis.
Definition: make_result.f90:22
m_fstr_nodalstress::fstr_nodalstress6d
subroutine fstr_nodalstress6d(hecMESH, fstrSOLID)
Calculate NODAL STRESS of shell elements.
Definition: fstr_NodalStress.f90:952
m_fstr::kststaticeigen
integer(kind=kint), parameter kststaticeigen
Definition: m_fstr.f90:42
m_fstr_nodalstress::fstr_nodalstress2d
subroutine fstr_nodalstress2d(hecMESH, fstrSOLID)
Calculate NODAL STRESS of plane elements.
Definition: fstr_NodalStress.f90:706
m_make_result
This module provide a function to prepare output of static analysis.
Definition: make_result.f90:6
m_static_output
This module provides functions to output result.
Definition: static_output.f90:6
m_hecmw2fstr_mesh_conv::fstr2hecmw_mesh_conv
subroutine fstr2hecmw_mesh_conv(hecMESH)
Definition: hecmw2fstr_mesh_conv.f90:27
m_make_result::fstr_make_result
subroutine, public fstr_make_result(hecMESH, fstrSOLID, fstrRESULT, istep, time, fstrDYNAMIC)
MAKE RESULT for static and dynamic analysis (WITHOUT ELEMENTAL RESULTS) -----------------------------...
Definition: make_result.f90:608
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_hecmw2fstr_mesh_conv::hecmw2fstr_mesh_conv
subroutine hecmw2fstr_mesh_conv(hecMESH)
Definition: hecmw2fstr_mesh_conv.f90:15
m_fstr::ivisual
integer(kind=kint), pointer ivisual
Definition: m_fstr.f90:123
m_fstr::iresult
integer(kind=kint), pointer iresult
Definition: m_fstr.f90:122
m_hecmw2fstr_mesh_conv
HECMW to FSTR Mesh Data Converter. Converting Connectivity of Element Type 232, 342 and 352.
Definition: hecmw2fstr_mesh_conv.f90:8
m_static_output::fstr_static_output
subroutine fstr_static_output(cstep, istep, time, hecMESH, fstrSOLID, fstrPARAM, flag, outflag)
Output result.
Definition: static_output.f90:13
m_fstr::ilog
integer(kind=kint), parameter ilog
FILE HANDLER.
Definition: m_fstr.f90:107
m_fstr_nodalstress
This module provides functions to calculation nodal stress.
Definition: fstr_NodalStress.f90:6
m_static_output::fstr_static_post
subroutine fstr_static_post(fnum, hecMESH, fstrSOLID, istep)
Summarizer of output data which prints out max and min output values.
Definition: static_output.f90:105
m_fstr_nodalstress::fstr_nodalstress3d
subroutine fstr_nodalstress3d(hecMESH, fstrSOLID)
Calculate NODAL STRESS of solid elements.
Definition: fstr_NodalStress.f90:16
mcontact
This module provides functions to calculate contact stiff matrix.
Definition: fstr_contact.f90:6