FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
fstr_precheck.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 contains
8 
9  !=============================================================================!
11  !=============================================================================!
12 
13  subroutine fstr_input_precheck( hecMESH, hecMAT, fstrSOLID )
14  use m_fstr
15  implicit none
16  type (hecmwST_local_mesh), intent(in) :: hecMESH
17  type (hecmwST_matrix), intent(in) :: hecMAT
18  type(fstr_solid), intent(in) :: fstrSOLID
19  integer(kint) :: i, cid
20 
21  if(fstrpr%solution_type == ksteigen .or. &
22  fstrpr%solution_type == kststaticeigen .or. &
23  fstrpr%solution_type == kstdynamic)then
24  do i = 1, hecmesh%section%n_sect
25  if(hecmesh%section%sect_type(i) == 4) cycle
26  cid = hecmesh%section%sect_mat_ID_item(i)
27 
28  if(fstrsolid%materials(cid)%variables(m_density) == 0.0d0)then
29  write(*,*) "*** error: density is not assigned or set to zero"
30  call hecmw_abort(hecmw_comm_get_comm())
31  endif
32  enddo
33  endif
34 
35  if(fstrpr%solution_type == kstprecheck)then
36  if(myrank == 0)then
37  write(imsg,*)
38  write(imsg,*) " **** STAGE Precheck **"
39  endif
40  call fstr_precheck_elem(hecmesh, hecmat)
41  write(idbg,*) "fstr_precheck_elem: OK"
42  endif
43 
44  if(fstrpr%solution_type == kstnzprof)then
45  call hecmw_nonzero_profile(hecmesh, hecmat)
46  endif
47  end subroutine fstr_input_precheck
48 
49  subroutine fstr_get_thickness(hecMESH,mid,thick)
50  use hecmw
51  use m_fstr
52  implicit none
53  type (hecmwST_local_mesh) :: hecMESH
54  integer(kind=kint) :: mid, ihead
55  real(kind=kreal) :: thick
56 
57  ihead = hecmesh%section%sect_R_index(mid-1)
58  thick = hecmesh%section%sect_R_item(ihead+1)
59  ! if(thick.LE.0.0) then
60  ! write(*,*) "Zero thickness <= 0 is illegal"
61  ! call hecmw_abort( hecmw_comm_get_comm())
62  ! endif
63  end subroutine fstr_get_thickness
64  !C
65  !C***
66  !C*** Pre Check for FSTR solver
67  !C***
68  !C
69  subroutine fstr_precheck ( hecMESH, hecMAT, soltype )
70  use m_fstr
71  implicit none
72  type (hecmwST_local_mesh) :: hecMESH
73  type (hecmwST_matrix ) :: hecMAT
74  integer(kind=kint) :: soltype
75 
76  if(myrank .EQ. 0) then
77  write(imsg,*)
78  write(imsg,*) ' **** STAGE PreCheck **'
79  endif
80 
81  call fstr_precheck_elem ( hecmesh, hecmat )
82  write(idbg,*) 'fstr_precheck_elem: OK'
83 
84  !C
85  !C Output sparsity pattern
86  !C
87  if( soltype == kstnzprof )then
88  call hecmw_nonzero_profile(hecmesh, hecmat)
89  endif
90 
91  end subroutine fstr_precheck
92  !C
93  !C
94  subroutine fstr_precheck_elem ( hecMESH, hecMAT )
95  use m_fstr
99 
100  implicit none
101 
102  type (hecmwST_matrix) :: hecMAT
103  type (hecmwST_local_mesh) :: hecMESH
104  type (fstr_solid) :: fstrSOLID
105 
106  !** Local variables
107  integer(kind=kint) :: nelem, mid, j, isect, nline, tline, icel, iiS
108  integer(kind=kint) :: ndof2, nelem_wo_mpc
109  integer(kind=kint) :: ie, ia, jelem, ic_type, nn, jS, jE, itype
110  integer(kind=kint) :: nodLOCAL(20),NTOTsum(1)
111  integer(kind=kint) :: nonzero
112  real(kind=kreal) :: ntdof2
113  real(kind=kreal) :: al, almin, almax, aa, thick, vol, avvol
114  real(kind=kreal) :: tvol, tvmax, tvmin, tlmax, tlmin, asp, aspmax
115  real(kind=kreal) :: xx(20),yy(20),zz(20)
116  real(kind=kreal) :: totsum(1),totmax(3),totmin(2)
117  !C
118  !C INIT
119  !C
120  nelem = 0
121  tvol = 0.0
122  tvmax = 0.0
123  tvmin = 1.0e+20
124  tlmax = 0.0
125  tlmin = 1.0e+20
126  aspmax = 0.0
127 
128  !C
129  !C Mesh Summary
130  !C
131  write(ilog,"(a)") '### Mesh Summary ###'
132  write(ilog,"(a,i12)") ' Num of node:',hecmesh%n_node
133  write(* ,"(a,i12)") ' Num of node:',hecmesh%n_node
134  write(ilog,"(a,i12)") ' Num of DOF :',hecmesh%n_node*hecmesh%n_dof
135  write(* ,"(a,i12)") ' Num of DOF :',hecmesh%n_node*hecmesh%n_dof
136  ndof2 = hecmesh%n_dof**2
137  ntdof2 = dble(hecmesh%n_node*hecmesh%n_dof)**2
138  write(ilog,"(a,i12)") ' Num of elem:',hecmesh%n_elem
139  write(* ,"(a,i12)") ' Num of elem:',hecmesh%n_elem
140  nelem_wo_mpc = 0
141  do itype = 1, hecmesh%n_elem_type
142  js = hecmesh%elem_type_index(itype-1)
143  je = hecmesh%elem_type_index(itype )
144  ic_type = hecmesh%elem_type_item(itype)
145  if(.not. (hecmw_is_etype_link(ic_type) .or. hecmw_is_etype_patch(ic_type))) &
146  nelem_wo_mpc = nelem_wo_mpc + je-js
147  enddo
148  write(ilog,"(a,i12)") ' ** w/o MPC/Patch:',nelem_wo_mpc
149  write(* ,"(a,i12)") ' ** w/o MPC/Patch:',nelem_wo_mpc
150  do itype = 1, hecmesh%n_elem_type
151  js = hecmesh%elem_type_index(itype-1)
152  je = hecmesh%elem_type_index(itype )
153  ic_type = hecmesh%elem_type_item(itype)
154  write(ilog,"(a,i4,a,i12)") ' Num of ',ic_type,':',je-js
155  write(* ,"(a,i4,a,i12)") ' Num of ',ic_type,':',je-js
156  enddo
157  nonzero = ndof2*(hecmat%NP + hecmat%NPU + hecmat%NPL)
158  write(ilog,"(a,i12)") ' Num of NZ :',nonzero
159  write(* ,"(a,i12)") ' Num of NZ :',nonzero
160  write(ilog,"(a,1pe12.5,a)") ' Sparsity :',100.0d0*dble(nonzero)/dble(ntdof2),"[%]"
161  write(* ,"(a,1pe12.5,a)") ' Sparsity :',100.0d0*dble(nonzero)/dble(ntdof2),"[%]"
162 
163  !C
164  !C 3D
165  !C
166  if( hecmesh%n_dof .eq. 3 ) then
167  do ie = 1, hecmesh%n_elem
168  ia = hecmesh%elem_ID(ie*2)
169  if( ia.ne.hecmesh%my_rank ) cycle
170  ! je = hecMESH%elem_ID(ie*2-1)
171  jelem = hecmesh%global_elem_ID(ie)
172  !C
173  ic_type = hecmesh%elem_type(ie)
174  !C
175  if (.not. (hecmw_is_etype_rod(ic_type) .or. hecmw_is_etype_solid(ic_type) &
176  & .or. hecmw_is_etype_beam(ic_type) .or. hecmw_is_etype_shell(ic_type))) then
177  write(ilog,*) jelem, ' This Element cannot be checked. Type=',ic_type
178  cycle
179  endif
180  nn = hecmw_get_max_node(ic_type)
181  !C
182  js = hecmesh%elem_node_index(ie-1)
183  je = hecmesh%elem_node_index(ie)
184 
185  do j = 1, nn
186  nodlocal(j) = hecmesh%elem_node_item (js+j)
187  xx(j) = hecmesh%node(3*nodlocal(j)-2)
188  yy(j) = hecmesh%node(3*nodlocal(j)-1)
189  zz(j) = hecmesh%node(3*nodlocal(j) )
190  enddo
191  !C
192  if ( ic_type.eq.111 ) then
193  isect = hecmesh%section_ID(ie)
194  mid = hecmesh%section%sect_mat_ID_item(isect)
195  call fstr_get_thickness( hecmesh,mid,aa )
196  al = sqrt( (xx(2)-xx(1))**2+(yy(2)-yy(1))**2+(zz(2)-zz(1))**2 )
197  nline = 1
198  tline = al
199  vol = aa*al
200  almax = al
201  almin = al
202  elseif( ic_type.eq.341 ) then
203  call pre_341 ( xx,yy,zz,vol,almax,almin )
204  elseif( ic_type.eq.351 ) then
205  call pre_351 ( xx,yy,zz,vol,almax,almin )
206  elseif( ic_type.eq.361 ) then
207  call pre_361 ( xx,yy,zz,vol,almax,almin )
208  elseif( ic_type.eq.342 ) then
209  call pre_342 ( xx,yy,zz,vol,almax,almin )
210  elseif( ic_type.eq.352 ) then
211  call pre_352 ( xx,yy,zz,vol,almax,almin )
212  elseif( ic_type.eq.362 ) then
213  call pre_362 ( xx,yy,zz,vol,almax,almin )
214  elseif( ic_type.eq.641 ) then
215  vol = 1.0d-12
216  elseif( ic_type.eq.761 ) then
217  vol = 1.0d-12
218  elseif( ic_type.eq.781 ) then
219  vol = 1.0d-12
220  endif
221  !C
222  if( vol.le.0.0 ) then
223  write(ilog,*) ' %%% ERROR %%% Volume of Element no.=',jelem,' is zero or negative.'
224  endif
225  nelem = nelem + 1
226  tvol = tvol + vol
227  if( vol.gt.tvmax ) tvmax = vol
228  if( vol.lt.tvmin ) tvmin = vol
229  if( almax.gt.tlmax ) tlmax = almax
230  if( almin.lt.tlmin ) tlmin = almin
231  asp = almax/almin
232  if( asp.gt.aspmax ) aspmax = asp
233  if( asp.gt.50 ) then
234  write(ilog,*) ' %%% WARNING %%% Aspect ratio of Element no.=',jelem,' exceeds 50.'
235  write(ilog,*) ' Maximum length =',almax
236  write(ilog,*) ' Minimum length =',almin
237  endif
238  enddo
239  !C
240  !C 2D
241  !C
242  elseif( hecmesh%n_dof .eq. 2 ) then
243  do ie = 1, hecmesh%n_elem
244  ia = hecmesh%elem_ID(ie*2)
245  if( ia.ne.hecmesh%my_rank ) cycle
246  ! je = hecMESH%elem_ID(ie*2-1)
247  jelem = hecmesh%global_elem_ID(ie)
248  !C
249  ic_type = hecmesh%elem_type(ie)
250  !C
251  if (.not. (hecmw_is_etype_rod(ic_type) .or. hecmw_is_etype_surface(ic_type))) then
252  write(ilog,*) jelem, ' This Element cannot be checked. Type=',ic_type
253  cycle
254  endif
255  nn = hecmw_get_max_node(ic_type)
256  !C
257  js = hecmesh%elem_node_index(ie-1)
258  do j = 1, nn
259  nodlocal(j) = hecmesh%elem_node_item (js+j)
260  xx(j) = hecmesh%node(3*nodlocal(j)-2)
261  yy(j) = hecmesh%node(3*nodlocal(j)-1)
262  enddo
263  !C
264  isect = hecmesh%section_ID(ie)
265  mid = hecmesh%section%sect_mat_ID_item(isect)
266  call fstr_get_thickness( hecmesh,mid,aa )
267  !C
268  if ( ic_type.eq.111 ) then
269  al = sqrt( (xx(2)-xx(1))**2+(yy(2)-yy(1))**2 )
270  vol = aa*al
271  if( al.gt.tlmax ) tlmax = al
272  if( al.lt.tlmin ) tlmin = al
273  aspmax = 1.0
274  elseif( ic_type.eq.231 ) then
275  call pre_231 ( xx,yy,aa,vol,almax,almin )
276  elseif( ic_type.eq.241 ) then
277  call pre_241 ( xx,yy,aa,vol,almax,almin )
278  elseif( ic_type.eq.232 ) then
279  call pre_232 ( xx,yy,aa,vol,almax,almin )
280  elseif( ic_type.eq.242 ) then
281  call pre_242 ( xx,yy,aa,vol,almax,almin )
282  else
283  vol = 0.0
284  endif
285  !C
286  if( vol.le.0.0 ) then
287  write(ilog,*) ' %%% ERROR %%% Volume of Element no.=',jelem,' is zero or negative.'
288  endif
289  nelem = nelem + 1
290  tvol = tvol + vol
291  if( vol.gt.tvmax ) tvmax = vol
292  if( vol.lt.tvmin ) tvmin = vol
293  if( almax.gt.tlmax ) tlmax = almax
294  if( almin.lt.tlmin ) tlmin = almin
295  asp = almax/almin
296  if( asp.gt.aspmax ) aspmax = asp
297  if( asp.gt.50 ) then
298  write(ilog,*) ' %%% WARNING %%% Aspect ratio of Element no.=',jelem,' exceeds 50.'
299  write(ilog,*) ' Maximum length =',almax
300  write(ilog,*) ' Minimum length =',almin
301  endif
302  enddo
303  !C
304  !C SHELL
305  !C
306  elseif( hecmesh%n_dof .eq. 6 ) then
307  do ie = 1, hecmesh%n_elem
308  ia = hecmesh%elem_ID(ie*2)
309  if( ia.ne.hecmesh%my_rank ) cycle
310  ! je = hecMESH%elem_ID(ie*2-1)
311  jelem = hecmesh%global_elem_ID(ie)
312  !C
313  ic_type = hecmesh%elem_type(ie)
314  !C
315  if (.not. (hecmw_is_etype_beam(ic_type) .or. hecmw_is_etype_shell(ic_type))) then
316  write(ilog,*) jelem, ' This Element cannot be checked. Type=',ic_type
317  cycle
318  endif
319  nn = hecmw_get_max_node(ic_type)
320  !C
321  js = hecmesh%elem_node_index(ie-1)
322  do j = 1, nn
323  nodlocal(j) = hecmesh%elem_node_item (js+j)
324  xx(j) = hecmesh%node(3*nodlocal(j)-2)
325  yy(j) = hecmesh%node(3*nodlocal(j)-1)
326  zz(j) = hecmesh%node(3*nodlocal(j) )
327  enddo
328  !C
329  isect = hecmesh%section_ID(ie)
330  mid = hecmesh%section%sect_mat_ID_item(isect)
331  call fstr_get_thickness( hecmesh,mid,aa )
332  !C
333  if ( ic_type.eq.111 ) then
334  al = sqrt( (xx(2)-xx(1))**2+(yy(2)-yy(1))**2+(zz(2)-zz(1))**2 )
335  nline = nline + 1
336  tline = tline + al
337  vol = aa*al
338  if( al.gt.tlmax ) tlmax = al
339  if( al.lt.tlmin ) tlmin = al
340  aspmax = 1.0
341  elseif( ic_type.eq.731 ) then
342  call pre_731 ( xx,yy,zz,aa,vol,almax,almin )
343  elseif( ic_type.eq.741 ) then
344  call pre_741 ( xx,yy,zz,aa,vol,almax,almin )
345  endif
346  !C
347  if( vol.le.0.0 ) then
348  write(ilog,*) ' %%% ERROR %%% Volume of Element no.=',jelem,' is zero or negative.'
349  endif
350  nelem = nelem + 1
351  tvol = tvol + vol
352  if( vol.gt.tvmax ) tvmax = vol
353  if( vol.lt.tvmin ) tvmin = vol
354  if( almax.gt.tlmax ) tlmax = almax
355  if( almin.lt.tlmin ) tlmin = almin
356  asp = almax/almin
357  if( asp.gt.aspmax ) aspmax = asp
358  if( asp.gt.50 ) then
359  write(ilog,*) ' %%% WARNING %%% Aspect ratio of Element no.=',jelem,' exceeds 50.'
360  write(ilog,*) ' Maximum length =',almax
361  write(ilog,*) ' Minimum length =',almin
362  endif
363  enddo
364  endif
365  !C
366  avvol = tvol / nelem
367  write(ilog,*) '### Sub Summary ###'
368  write(ilog,*) ' Total Volumes in this region = ',tvol
369  write(ilog,*) ' Average Volume of elements = ',avvol
370  write(ilog,*) ' Maximum Volume of elements = ',tvmax
371  write(ilog,*) ' Minimum Volume of elements = ',tvmin
372  write(ilog,*) ' Maximum length of element edges = ',tlmax
373  write(ilog,*) ' Minimum length of element edges = ',tlmin
374 
375  write(ilog,*) ' Maximum aspect ratio in this region = ',aspmax
376  totsum(1) = tvol
377  call hecmw_allreduce_r(hecmesh,totsum,1,hecmw_sum)
378  ntotsum(1) = nelem
379  call hecmw_allreduce_i(hecmesh,ntotsum,1,hecmw_sum)
380  totmax(1) = tvmax
381  totmax(2) = tlmax
382  totmax(3) = aspmax
383  call hecmw_allreduce_r(hecmesh,totmax,3,hecmw_max)
384  totmin(1) = tvmin
385  totmin(2) = tlmin
386  call hecmw_allreduce_r(hecmesh,totmin,2,hecmw_min)
387  if( hecmesh%my_rank .eq. 0 ) then
388  avvol = totsum(1) / ntotsum(1)
389  write(ilog,*) '### Global Summary ###'
390  write(ilog,*) ' TOTAL VOLUME = ',totsum(1)
391  write(*,*) ' TOTAL VOLUME = ',totsum(1)
392  write(ilog,*) ' AVERAGE VOLUME OF ELEMENTS = ',avvol
393  write(*,*) ' AVERAGE VOLUME OF ELEMENTS = ',avvol
394  write(ilog,*) ' MAXIMUM VOLUME OF ELEMENTS = ',totmax(1)
395  write(*,*) ' MAXIMUM VOLUME OF ELEMENTS = ',totmax(1)
396  write(ilog,*) ' MINIMUM VOLUME OF ELEMENTS = ',totmin(1)
397  write(*,*) ' MINIMUM VOLUME OF ELEMENTS = ',totmin(1)
398  write(ilog,*) ' MAXIMUM LENGTH OF ELEMENT EDGES = ',totmax(2)
399  write(*,*) ' MAXIMUM LENGTH OF ELEMENT EDGES = ',totmax(2)
400  write(ilog,*) ' MINIMUM LENGTH OF ELEMENT EDGES = ',totmin(2)
401  write(*,*) ' MINIMUM LENGTH OF ELEMENT EDGES = ',totmin(2)
402  write(ilog,*) ' MAXIMUM ASPECT RATIO = ',totmax(3)
403  write(*,*) ' MAXIMUM ASPECT RATIO = ',totmax(3)
404  endif
405  !C
406  end subroutine fstr_precheck_elem
407 
408  subroutine hecmw_nonzero_profile(hecMESH, hecMAT)
410  implicit none
411  type (hecmwST_local_mesh) :: hecMESH
412  type (hecmwST_matrix) :: hecMAT
413 
414  integer(kind=kint) :: i, j, in, jS, jE, ftype, n, ndof, nnz, fio
415  real(kind=kreal) :: rnum, dens, cond
416  character :: fileid*3
417 
418  fio = 70 + hecmesh%my_rank
419  write(fileid,"(i3.3)")hecmesh%my_rank
420 
421  !ftype = 2: eps
422  !ftype = 4: png
423  ftype = 4
424 
425  n = hecmat%N
426  ndof = 3*n
427  nnz = 9*n + 9*2*hecmat%indexL(n)
428  dens = 100*dble(nnz)/dble(9*n*n)
429  rnum = (7.21d0+0.01*dlog10(dble(hecmat%N)))*10.0d0/dble(hecmat%N)
430  cond = 1.0d0
431  !rnum = (7.25d0)*10.0d0/dble(hecMAT%N)
432 
433  open(fio,file='nonzero.dat.'//fileid, status='replace')
434  !write(fio,"(a,f12.5,i0)")"##magic number 10 : 7.2, ",rnum,hecMAT%N
435  do i= 1, n
436  js= hecmat%indexL(i-1) + 1
437  je= hecmat%indexL(i )
438  write(fio,"(i0,a,i0)")i," ",i
439  do j= js, je
440  in = hecmat%itemL(j)
441  write(fio,"(i0,a,i0)")i, " ",in
442  write(fio,"(i0,a,i0)")in," ",i
443  enddo
444  enddo
445  close(fio)
446 
447  open(fio,file='nonzero.plt.'//fileid, status='replace')
448  if(ftype == 4)then
449  write(fio,"(a)")'set terminal png size 1500,1500'
450  else
451  write(fio,"(a)")'set terminal postscript eps enhanced color solid "TimesNewRomanPSMT" 20'
452  endif
453 
454  write(fio,"(a)")'unset key'
455  write(fio,"(a)")'unset xtics'
456  write(fio,"(a)")'unset ytics'
457  write(fio,"(a)")'set size ratio 1.0'
458  write(fio,"(a)")'set border lw 1.0'
459  write(fio,"(a,i0,a)")'set xrange[0.5:',n,'.5]'
460  write(fio,"(a,i0,a)")'set yrange[0.5:',n,'.5] reverse '
461 
462  if(ftype == 4)then
463  write(fio,"(a)")'set out "image.'//fileid//'.png"'
464  else
465  write(fio,"(a)")'set out "image.'//fileid//'.eps"'
466  write(fio,"(a)" )'set label 1 "Name" at graph 1.1,0.9'
467  write(fio,"(a)")'set label 2 "N" at graph 1.1,0.85'
468  write(fio,"(a)")'set label 3 "Non-Zero Elem." at graph 1.1,0.8'
469  write(fio,"(a)")'set label 4 "Density [%]" at graph 1.1,0.75'
470  write(fio,"(a)")'set label 9 "Condition Num." at graph 1.1,0.7'
471  write(fio,"(a)" )'set label 5 ": matrix" at graph 1.4,0.9'
472  write(fio,"(a,i0,a)")'set label 6 ": ',ndof,'" at graph 1.4,0.85'
473  write(fio,"(a,i0,a)")'set label 7 ": ',nnz,'" at graph 1.4,0.8'
474  write(fio,"(a,1pe9.2,a)")'set label 8 ": ',dens,'" at graph 1.4,0.75'
475  write(fio,"(a,1pe9.2,a)")'set label 10 ": ',cond,'" at graph 1.4,0.7'
476  endif
477 
478  write(fio,"(a,f12.5,a)")'plot "nonzero.dat.'//fileid//'" pointtype 5 pointsize ',rnum,' linecolor rgb "#F96566"'
479  close(fio)
480 
481  write(*,*)''
482  write(*,*)' ### Command recommendation'
483  write(*,*)' gnuplot -persist "nonzero.plt"'
484 
485  !call system('gnuplot -persist "nonzero.plt"')
486  !open(fio,file='nonzero.dat',status='old')
487  !close(fio,status="delete")
488  end subroutine hecmw_nonzero_profile
489 end module m_fstr_precheck
m_fstr_precheck::hecmw_nonzero_profile
subroutine hecmw_nonzero_profile(hecMESH, hecMAT)
Definition: fstr_precheck.f90:409
m_fstr_precheck::fstr_precheck
subroutine fstr_precheck(hecMESH, hecMAT, soltype)
Definition: fstr_precheck.f90:70
m_fstr::myrank
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.f90:96
m_fstr::kststaticeigen
integer(kind=kint), parameter kststaticeigen
Definition: m_fstr.f90:42
m_fstr::fstr_solid
Definition: m_fstr.f90:238
m_precheck_lib_shell::pre_741
subroutine pre_741(XX, YY, ZZ, thick, vol, almax, almin)
Definition: precheck_LIB_shell.f90:51
m_fstr::idbg
integer(kind=kint), parameter idbg
Definition: m_fstr.f90:111
m_precheck_lib_2d::pre_231
subroutine pre_231(XX, YY, thick, vol, almax, almin)
Definition: precheck_LIB_2d.f90:17
m_fstr_precheck::fstr_get_thickness
subroutine fstr_get_thickness(hecMESH, mid, thick)
Definition: fstr_precheck.f90:50
m_fstr_precheck::fstr_precheck_elem
subroutine fstr_precheck_elem(hecMESH, hecMAT)
Definition: fstr_precheck.f90:95
m_fstr::fstrpr
type(fstr_param), target fstrpr
GLOBAL VARIABLE INITIALIZED IN FSTR_SETUP.
Definition: m_fstr.f90:208
m_precheck_lib_3d::pre_361
subroutine pre_361(XX, YY, ZZ, vol, almax, almin)
Definition: precheck_LIB_3d.f90:245
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
m_precheck_lib_3d
This module provides function to check input data of 3d static analysis.
Definition: precheck_LIB_3d.f90:6
m_precheck_lib_2d::pre_241
subroutine pre_241(XX, YY, thick, vol, almax, almin)
Definition: precheck_LIB_2d.f90:93
m_precheck_lib_2d::pre_242
subroutine pre_242(XX, YY, thick, vol, almax, almin)
Definition: precheck_LIB_2d.f90:259
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
m_fstr_precheck::fstr_input_precheck
subroutine fstr_input_precheck(hecMESH, hecMAT, fstrSOLID)
fstr_input_precheck !
Definition: fstr_precheck.f90:14
m_precheck_lib_3d::pre_362
subroutine pre_362(XX, YY, ZZ, vol, almax, almin)
Definition: precheck_LIB_3d.f90:705
m_fstr::ksteigen
integer(kind=kint), parameter ksteigen
Definition: m_fstr.f90:38
hecmw
Definition: hecmw.f90:6
m_precheck_lib_3d::pre_341
subroutine pre_341(XX, YY, ZZ, vol, almax, almin)
Definition: precheck_LIB_3d.f90:18
m_fstr::kstprecheck
integer(kind=kint), parameter kstprecheck
solution type (st)
Definition: m_fstr.f90:36
m_fstr_precheck
This module provides function to check input data of IFSTR solver.
Definition: fstr_precheck.f90:6
m_precheck_lib_2d::pre_232
subroutine pre_232(XX, YY, thick, vol, almax, almin)
Definition: precheck_LIB_2d.f90:164
m_precheck_lib_3d::pre_342
subroutine pre_342(XX, YY, ZZ, vol, almax, almin)
Definition: precheck_LIB_3d.f90:372
m_precheck_lib_shell
This module provides function to check input data of shell elements.
Definition: precheck_LIB_shell.f90:6
m_precheck_lib_3d::pre_351
subroutine pre_351(XX, YY, ZZ, vol, almax, almin)
Definition: precheck_LIB_3d.f90:125
m_fstr::ilog
integer(kind=kint), parameter ilog
FILE HANDLER.
Definition: m_fstr.f90:107
m_precheck_lib_3d::pre_352
subroutine pre_352(XX, YY, ZZ, vol, almax, almin)
Definition: precheck_LIB_3d.f90:521
m_precheck_lib_shell::pre_731
subroutine pre_731(XX, YY, ZZ, thick, vol, almax, almin)
Definition: precheck_LIB_shell.f90:14
m_precheck_lib_2d
This module provides function to check input data of 2d static analysis.
Definition: precheck_LIB_2d.f90:6
m_fstr::kstdynamic
integer(kind=kint), parameter kstdynamic
Definition: m_fstr.f90:40
m_fstr::kstnzprof
integer(kind=kint), parameter kstnzprof
Definition: m_fstr.f90:43
m_fstr::imsg
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:110