FrontISTR  5.7.1
Large-scale structural analysis program with finit element method
precheck_LIB_shell.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  ! SHELL Element:
10  ! PRE_731( XX,YY,ZZ,thick,vol,almax,almin )
11  ! PRE_741( XX,YY,ZZ,thick,vol,almax,almin )
12  !----------------------------------------------------------------------*
13  subroutine pre_731( XX,YY,ZZ,thick,vol,almax,almin )
14  !----------------------------------------------------------------------*
15  !**
16  !** Precheck for 3nodes SHELL
17  !**
18  use hecmw
19  implicit none
20  ! I/F VARIABLES
21  real(kind=kreal) xx(*),yy(*),zz(*),thick,vol,almax,almin
22  ! LOCAL VARIABLES
23  real(kind=kreal) v1x,v1y,v1z
24  real(kind=kreal) v2x,v2y,v2z
25  real(kind=kreal) v3x,v3y,v3z
26  real(kind=kreal) area,a1,a2,a3
27  area = 0.0
28  vol = 0.0
29  !** FACE 1-2-3
30  v1x=xx(2)-xx(1)
31  v1y=yy(2)-yy(1)
32  v1z=zz(2)-zz(1)
33  v2x=xx(3)-xx(1)
34  v2y=yy(3)-yy(1)
35  v2z=zz(3)-zz(1)
36  v3x= v1y*v2z-v1z*v2y
37  v3y=-v1x*v2z+v1z*v2x
38  v3z= v1x*v2y-v1y*v2x
39  area=sqrt( v3x*v3x + v3y*v3y + v3z*v3z )*0.5
40  vol = area * thick
41  a1 = sqrt( (xx(2)-xx(1))**2+(yy(2)-yy(1))**2+(zz(2)-zz(1))**2 )
42  a2 = sqrt( (xx(3)-xx(2))**2+(yy(3)-yy(2))**2+(zz(3)-zz(2))**2 )
43  a3 = sqrt( (xx(1)-xx(3))**2+(yy(1)-yy(3))**2+(zz(1)-zz(3))**2 )
44  almax = dmax1( a1,a2,a3 )
45  almin = dmin1( a1,a2,a3 )
46 
47  end subroutine pre_731
48  !----------------------------------------------------------------------*
49  subroutine pre_741( XX,YY,ZZ,thick,vol,almax,almin )
50  !----------------------------------------------------------------------*
51  !**
52  !** Precheck for 3nodes SHELL
53  !**
54  use hecmw
55  use quadrature
56  implicit none
57  ! I/F VARIABLES
58  real(kind=kreal) xx(*),yy(*),zz(*),thick,vol,almax,almin
59  ! LOCAL VARIABLES
60  integer(kind=kint) NN
61  integer(kind=kint) NG
62  parameter(nn=8,ng=2)
63  real(kind=kreal) h(nn),hr(nn),hs(nn),ht(nn)
64  real(kind=kreal) ri,si,ti,rp,sp,tp,rm,sm,tm
65  real(kind=kreal) xj11,xj21,xj31,xj12,xj22,xj32,xj13,xj23,xj33,det,wg
66  integer(kind=kint) IG1,IG2,LX,LY,LZ,I
67  real(kind=kreal) vx,vy,vz,xcod,ycod,zcod
68  real(kind=kreal) ax,ay,az,rx,ry,rz,hx,hy,hz,val
69  real(kind=kreal) phx,phy,phz
70  real(kind=kreal) g1x,g1y,g1z
71  real(kind=kreal) g2x,g2y,g2z
72  real(kind=kreal) g3x,g3y,g3z
73  real(kind=kreal) xsum,coefx,coefy,coefz
74  real(kind=kreal) area,a1,a2,a3,a4
75  !
76  area = 0.0
77  vol = 0.0
78  ! INTEGRATION OVER SURFACE
79  do ig2=1,ng
80  si=gauss1d2(1,ig2)
81  do ig1=1,ng
82  ri=gauss1d2(1,ig1)
83  h(1)=0.25*(1.0-ri)*(1.0-si)
84  h(2)=0.25*(1.0+ri)*(1.0-si)
85  h(3)=0.25*(1.0+ri)*(1.0+si)
86  h(4)=0.25*(1.0-ri)*(1.0+si)
87  hr(1)=-.25*(1.0-si)
88  hr(2)= .25*(1.0-si)
89  hr(3)= .25*(1.0+si)
90  hr(4)=-.25*(1.0+si)
91  hs(1)=-.25*(1.0-ri)
92  hs(2)=-.25*(1.0+ri)
93  hs(3)= .25*(1.0+ri)
94  hs(4)= .25*(1.0-ri)
95  g1x=0.0
96  g1y=0.0
97  g1z=0.0
98  g2x=0.0
99  g2y=0.0
100  g2z=0.0
101  do i=1,nn
102  g1x=g1x+hr(i)*xx(i)
103  g1y=g1y+hr(i)*yy(i)
104  g1z=g1z+hr(i)*zz(i)
105  g2x=g2x+hs(i)*xx(i)
106  g2y=g2y+hs(i)*yy(i)
107  g2z=g2z+hs(i)*zz(i)
108  enddo
109  g3x=g1y*g2z-g1z*g2y
110  g3y=g1z*g2x-g1x*g2z
111  g3z=g1x*g2y-g1y*g2x
112  xsum=dsqrt(g3x**2+g3y**2+g3z**2)
113  g3x=g3x/xsum
114  g3y=g3y/xsum
115  g3z=g3z/xsum
116  !JACOBI MATRIX
117  xj11=g1x
118  xj12=g1y
119  xj13=g1z
120  xj21=g2x
121  xj22=g2y
122  xj23=g2z
123  xj31=g3x
124  xj32=g3y
125  xj33=g3z
126  !DETERMINANT OF JACOBIAN
127  det=xj11*xj22*xj33 &
128  +xj12*xj23*xj31 &
129  +xj13*xj21*xj32 &
130  -xj13*xj22*xj31 &
131  -xj12*xj21*xj33 &
132  -xj11*xj23*xj32
133  wg=weight1d2(ig1)*weight1d2(ig2)*det
134  do i = 1, nn
135  area = area + h(i)*wg
136  enddo
137  enddo
138  enddo
139 
140  vol = area*thick
141  a1 = sqrt( (xx(2)-xx(1))**2+(yy(2)-yy(1))**2+(zz(2)-zz(1))**2 )
142  a2 = sqrt( (xx(3)-xx(2))**2+(yy(3)-yy(2))**2+(zz(3)-zz(2))**2 )
143  a3 = sqrt( (xx(4)-xx(3))**2+(yy(4)-yy(3))**2+(zz(4)-zz(3))**2 )
144  a4 = sqrt( (xx(1)-xx(4))**2+(yy(1)-yy(4))**2+(zz(1)-zz(4))**2 )
145  almax = dmax1( a1,a2,a3,a4 )
146  almin = dmin1( a1,a2,a3,a4 )
147 
148  end subroutine pre_741
149 end module m_precheck_lib_shell
Definition: hecmw.f90:6
This module provides function to check input data of shell elements.
subroutine pre_741(XX, YY, ZZ, thick, vol, almax, almin)
subroutine pre_731(XX, YY, ZZ, thick, vol, almax, almin)
This module contains Gauss point information.
Definition: quadrature.f90:28
real(kind=kreal), dimension(2) weight1d2
Definition: quadrature.f90:32
real(kind=kreal), dimension(1, 2) gauss1d2
Definition: quadrature.f90:32