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