FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
precheck_LIB_2d.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  ! 2D Element: PreCheck
10  ! PRE_231( xx,yy,thick,vol,almax,almin )
11  ! PRE_241( xx,yy,thick,vol,almax,almin )
12  ! PRE_232( xx,yy,thick,vol,almax,almin )
13  ! PRE_242( xx,yy,thick,vol,almax,almin )
14  !***********************************************************************
15  !----------------------------------------------------------------------*
16  subroutine pre_231( XX,YY,thick,vol,almax,almin )
17  !----------------------------------------------------------------------*
18  !
19  ! CALCULATION 2D 3 NODE PLANE ELEMENT
20  !
21  use hecmw
23  implicit none
24  ! I/F VARIABLES
25  real(kind=kreal) xx(*),yy(*),thick,vol,almax,almin
26  ! LOCAL VARIABLES
27  integer(kind=kint) NN
28  integer(kind=kint) NG
29  parameter(nn=3,ng=2)
30  real(kind=kreal) h(nn),hl1(nn),hl2(nn),hl3(nn)
31  integer(kind=kint) L1,L2,I
32  real(kind=kreal) x1,x2,x3,xl1,xl2
33  real(kind=kreal) ri,si,rp,sp,rm,sm
34  real(kind=kreal) xj11,xj21,xj12,xj22,det,wg,area
35  real(kind=kreal) a1,a2,a3
36  !C
37  area = 0.0
38  vol = 0.0
39  !* LOOP OVER ALL INTEGRATION POINTS
40  do l2=1,ng
41  xl2=xg(ng,l2)
42  x2 =(xl2+1.0)*0.5
43  do l1=1,ng
44  xl1=xg(ng,l1)
45  x1=0.5*(1.0-x2)*(xl1+1.0)
46  ! INTERPOLATION FUNCTION
47  x3=1.0-x1-x2
48  h(1)=x1
49  h(2)=x2
50  h(3)=x3
51  ! DERIVATIVE OF INTERPOLATION FUNCTION
52  ! FOR L1-COORDINATE
53  hl1(1)=1.0
54  hl1(2)=0.0
55  hl1(3)=0.0
56  ! FOR L2-COORDINATE
57  hl2(1)=0.0
58  hl2(2)=1.0
59  hl2(3)=0.0
60  ! FOR L3-COORDINATE
61  hl3(1)=0.0
62  hl3(2)=0.0
63  hl3(3)=1.0
64  ! JACOBI MATRIX
65  xj11=0.0
66  xj21=0.0
67  xj12=0.0
68  xj22=0.0
69  do i=1,nn
70  xj11=xj11+(hl1(i)-hl3(i))*xx(i)
71  xj21=xj21+(hl2(i)-hl3(i))*xx(i)
72  xj12=xj12+(hl1(i)-hl3(i))*yy(i)
73  xj22=xj22+(hl2(i)-hl3(i))*yy(i)
74  enddo
75  det=xj11*xj22-xj21*xj12
76  wg=wgt(ng,l1)*wgt(ng,l2)*det*(1.0-x2)*0.25
77  do i = 1, nn
78  area = area + h(i)*wg
79  enddo
80  enddo
81  enddo
82 
83  vol = area * thick
84  a1 = sqrt( (xx(2)-xx(1))**2+(yy(2)-yy(1))**2 )
85  a2 = sqrt( (xx(3)-xx(2))**2+(yy(3)-yy(2))**2 )
86  a3 = sqrt( (xx(1)-xx(3))**2+(yy(1)-yy(3))**2 )
87  almax = dmax1( a1,a2,a3 )
88  almin = dmin1( a1,a2,a3 )
89 
90  end subroutine pre_231
91  !----------------------------------------------------------------------*
92  subroutine pre_241( XX,YY,thick,vol,almax,almin )
93  !----------------------------------------------------------------------*
94  !
95  ! CALCULATION 2D 4 NODE PLANE ELEMENT
96  !
97  use hecmw
99  implicit none
100  ! I/F VARIABLES
101  real(kind=kreal) xx(*),yy(*),thick,vol,almax,almin
102  ! LOCAL VARIABLES
103  integer(kind=kint) NN
104  integer(kind=kint) NG
105  parameter(nn=4,ng=2)
106  real(kind=kreal) h(nn),hr(nn),hs(nn)
107  integer(kind=kint) LX,LY,I
108  real(kind=kreal) ri,si,rp,sp,rm,sm
109  real(kind=kreal) xj11,xj21,xj12,xj22,det,wg,area
110  real(kind=kreal) a1,a2,a3,a4
111  !C
112  area = 0.0
113  vol = 0.0
114  !* LOOP OVER ALL INTEGRATION POINTS
115  do lx=1,ng
116  ri=xg(ng,lx)
117  do ly=1,ng
118  si=xg(ng,ly)
119  rp=1.0+ri
120  sp=1.0+si
121  rm=1.0-ri
122  sm=1.0-si
123  h(1)=0.25*rm*sm
124  h(2)=0.25*rp*sm
125  h(3)=0.25*rp*sp
126  h(4)=0.25*rm*sp
127  hr(1)=-.25*sm
128  hr(2)= .25*sm
129  hr(3)= .25*sp
130  hr(4)=-.25*sp
131  hs(1)=-.25*rm
132  hs(2)=-.25*rp
133  hs(3)= .25*rp
134  hs(4)= .25*rm
135  xj11=0.0
136  xj21=0.0
137  xj12=0.0
138  xj22=0.0
139  do i=1,nn
140  xj11=xj11+hr(i)*xx(i)
141  xj21=xj21+hs(i)*xx(i)
142  xj12=xj12+hr(i)*yy(i)
143  xj22=xj22+hs(i)*yy(i)
144  enddo
145  det=xj11*xj22-xj21*xj12
146  wg=wgt(ng,lx)*wgt(ng,ly)*det
147  do i = 1, nn
148  area = area + h(i)*wg
149  enddo
150  enddo
151  enddo
152 
153  vol = area * thick
154  a1 = sqrt( (xx(2)-xx(1))**2+(yy(2)-yy(1))**2 )
155  a2 = sqrt( (xx(3)-xx(2))**2+(yy(3)-yy(2))**2 )
156  a3 = sqrt( (xx(4)-xx(3))**2+(yy(4)-yy(3))**2 )
157  a4 = sqrt( (xx(1)-xx(4))**2+(yy(1)-yy(4))**2 )
158  almax = dmax1( a1,a2,a3,a4 )
159  almin = dmin1( a1,a2,a3,a4 )
160 
161  end subroutine pre_241
162  !----------------------------------------------------------------------*
163  subroutine pre_232( XX,YY,thick,vol,almax,almin )
164  !----------------------------------------------------------------------*
165  !
166  ! CALCULATION 2D 6 NODE PLANE ELEMENT
167  !
168  use hecmw
170  implicit none
171  ! I/F VARIABLES
172  real(kind=kreal) xx(*),yy(*),thick,vol,tline,almax,almin
173  ! LOCAL VARIABLES
174  integer(kind=kint) NN
175  integer(kind=kint) NG
176  parameter(nn=6,ng=3)
177  real(kind=kreal) h(nn),hl1(nn),hl2(nn),hl3(nn)
178  integer(kind=kint) L1,L2,I
179  real(kind=kreal) x1,x2,x3,xl1,xl2
180  real(kind=kreal) ri,si,rp,sp,rm,sm
181  real(kind=kreal) xj11,xj21,xj12,xj22,det,wg,area
182  real(kind=kreal) a1,a2,al1,al2,al3
183 
184  !*************************
185  area = 0.0
186  vol = 0.0
187  !* LOOP OVER ALL INTEGRATION POINTS
188  do l2=1,ng
189  xl2=xg(ng,l2)
190  x2 =(xl2+1.0)*0.5
191  do l1=1,ng
192  xl1=xg(ng,l1)
193  x1=0.5*(1.0-x2)*(xl1+1.0)
194  ! INTERPOLATION FUNCTION
195  x3=1.0-x1-x2
196  h(1)= x1*(2.0*x1-1.)
197  h(2)= x2*(2.0*x2-1.)
198  h(3)= x3*(2.0*x3-1.)
199  h(4)= 4.0*x1*x2
200  h(5)= 4.0*x2*x3
201  h(6)= 4.0*x1*x3
202  ! DERIVATIVE OF INTERPOLATION FUNCTION
203  ! FOR L1-COORDINATE
204  hl1(1)=4.0*x1-1.0
205  hl1(2)= 0.0
206  hl1(3)= 0.0
207  hl1(4)= 4.0*x2
208  hl1(5)= 0.0
209  hl1(6)= 4.0*x3
210  ! FOR L2-COORDINATE
211  hl2(1)= 0.0
212  hl2(2)= 4.0*x2-1.0
213  hl2(3)= 0.0
214  hl2(4)= 4.0*x1
215  hl2(5)= 4.0*x3
216  hl2(6)= 0.0
217  ! FOR L3-COORDINATE
218  hl3(1)= 0.0
219  hl3(2)= 0.0
220  hl3(3)= 4.0*x3-1.0
221  hl3(4)= 0.0
222  hl3(5)= 4.0*x2
223  hl3(6)= 4.0*x1
224  ! JACOBI MATRIX
225  xj11=0.0
226  xj21=0.0
227  xj12=0.0
228  xj22=0.0
229  do i=1,nn
230  xj11=xj11+(hl1(i)-hl3(i))*xx(i)
231  xj21=xj21+(hl2(i)-hl3(i))*xx(i)
232  xj12=xj12+(hl1(i)-hl3(i))*yy(i)
233  xj22=xj22+(hl2(i)-hl3(i))*yy(i)
234  enddo
235  det=xj11*xj22-xj21*xj12
236  wg=wgt(ng,l1)*wgt(ng,l2)*det*(1.0-x2)*0.25
237  do i = 1, nn
238  area = area + h(i)*wg
239  enddo
240  enddo
241  enddo
242 
243  vol = area * thick
244  a1 = sqrt( (xx(4)-xx(1))**2+(yy(4)-yy(1))**2 )
245  a2 = sqrt( (xx(2)-xx(4))**2+(yy(2)-yy(4))**2 )
246  al1 = a1 + a2
247  a1 = sqrt( (xx(5)-xx(2))**2+(yy(5)-yy(2))**2 )
248  a2 = sqrt( (xx(3)-xx(5))**2+(yy(3)-yy(5))**2 )
249  al2 = a1 + a2
250  a1 = sqrt( (xx(6)-xx(3))**2+(yy(6)-yy(3))**2 )
251  a2 = sqrt( (xx(1)-xx(6))**2+(yy(1)-yy(6))**2 )
252  al3 = a1 + a2
253  almax = dmax1( al1,al2,al3 )
254  almin = dmin1( al1,al2,al3 )
255 
256  end subroutine pre_232
257  !----------------------------------------------------------------------*
258  subroutine pre_242( XX,YY,thick,vol,almax,almin )
259  !----------------------------------------------------------------------*
260  !
261  ! CALCULATION 2D 8 NODE PLANE ELEMENT
262  !
263  use hecmw
265  implicit none
266  ! I/F VARIABLES
267  real(kind=kreal) xx(*),yy(*),thick,vol,almax,almin
268  ! LOCAL VARIABLES
269  integer(kind=kint) NN
270  integer(kind=kint) NG
271  parameter(nn=8,ng=3)
272  real(kind=kreal) h(nn),hr(nn),hs(nn)
273  integer(kind=kint) LX,LY,I
274  real(kind=kreal) ri,si,rp,sp,rm,sm
275  real(kind=kreal) xj11,xj21,xj12,xj22,det,wg,area
276  !REAL(kind=kreal) RX(8),SX(8)
277  real(kind=kreal) a1,a2,al1,al2,al3,al4
278  !DATA RX/-1., 1.,1.,-1., 0., 1., 0., -1./
279  !DATA SX/-1.,-1.,1., 1.,-1., 0., 1., 0./
280  !*************************
281  area = 0.0
282  vol = 0.0
283  !* LOOP OVER ALL INTEGRATION POINTS
284  do lx=1,ng
285  ri=xg(ng,lx)
286  do ly=1,ng
287  si=xg(ng,ly)
288  rp=1.0+ri
289  sp=1.0+si
290  rm=1.0-ri
291  sm=1.0-si
292  h(1)=0.25*rm*sm*(-1.0-ri-si)
293  h(2)=0.25*rp*sm*(-1.0+ri-si)
294  h(3)=0.25*rp*sp*(-1.0+ri+si)
295  h(4)=0.25*rm*sp*(-1.0-ri+si)
296  h(5)=0.5*(1.0-ri*ri)*(1.0-si)
297  h(6)=0.5*(1.0-si*si)*(1.0+ri)
298  h(7)=0.5*(1.0-ri*ri)*(1.0+si)
299  h(8)=0.5*(1.0-si*si)*(1.0-ri)
300  hr(1)=-.25*sm*(-1.0-ri-si)-0.25*rm*sm
301  hr(2)= .25*sm*(-1.0+ri-si)+0.25*rp*sm
302  hr(3)= .25*sp*(-1.0+ri+si)+0.25*rp*sp
303  hr(4)=-.25*sp*(-1.0-ri+si)-0.25*rm*sp
304  hr(5)=-ri*(1.0-si)
305  hr(6)= 0.5*(1.0-si*si)
306  hr(7)=-ri*(1.0+si)
307  hr(8)=-0.5*(1.0-si*si)
308  hs(1)=-.25*rm*(-1.0-ri-si)-0.25*rm*sm
309  hs(2)=-.25*rp*(-1.0+ri-si)-0.25*rp*sm
310  hs(3)= .25*rp*(-1.0+ri+si)+0.25*rp*sp
311  hs(4)= .25*rm*(-1.0-ri+si)+0.25*rm*sp
312  hs(5)=-0.5*(1.0-ri*ri)
313  hs(6)=-si *(1.0+ri)
314  hs(7)= 0.5*(1.0-ri*ri)
315  hs(8)= -si*(1.0-ri)
316  xj11=0.0
317  xj21=0.0
318  xj12=0.0
319  xj22=0.0
320  do i=1,nn
321  xj11=xj11+hr(i)*xx(i)
322  xj21=xj21+hs(i)*xx(i)
323  xj12=xj12+hr(i)*yy(i)
324  xj22=xj22+hs(i)*yy(i)
325  enddo
326  det=xj11*xj22-xj21*xj12
327  wg=wgt(ng,lx)*wgt(ng,ly)*det
328  do i = 1, nn
329  area = area + h(i)*wg
330  enddo
331  enddo
332  enddo
333 
334  vol = area * thick
335  a1 = sqrt( (xx(5)-xx(1))**2+(yy(5)-yy(1))**2 )
336  a2 = sqrt( (xx(2)-xx(5))**2+(yy(2)-yy(5))**2 )
337  al1 = a1 + a2
338  a1 = sqrt( (xx(6)-xx(2))**2+(yy(6)-yy(2))**2 )
339  a2 = sqrt( (xx(3)-xx(6))**2+(yy(3)-yy(6))**2 )
340  al2 = a1 + a2
341  a1 = sqrt( (xx(7)-xx(3))**2+(yy(7)-yy(3))**2 )
342  a2 = sqrt( (xx(4)-xx(7))**2+(yy(4)-yy(7))**2 )
343  al3 = a1 + a2
344  a1 = sqrt( (xx(8)-xx(4))**2+(yy(8)-yy(4))**2 )
345  a2 = sqrt( (xx(1)-xx(8))**2+(yy(1)-yy(8))**2 )
346  al4 = a1 + a2
347  almax = dmax1( al1,al2,al3,al4 )
348  almin = dmin1( al1,al2,al3,al4 )
349 
350  end subroutine pre_242
351 end module m_precheck_lib_2d
gauss_integration
This module provides data for gauss quadrature.
Definition: GaussM.f90:6
m_precheck_lib_2d::pre_231
subroutine pre_231(XX, YY, thick, vol, almax, almin)
Definition: precheck_LIB_2d.f90:17
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
hecmw
Definition: hecmw.f90:6
gauss_integration::xg
real(kind=kreal), dimension(3, 3) xg
abscissa of gauss points
Definition: GaussM.f90:9
m_precheck_lib_2d::pre_232
subroutine pre_232(XX, YY, thick, vol, almax, almin)
Definition: precheck_LIB_2d.f90:164
gauss_integration::wgt
real(kind=kreal), dimension(3, 3) wgt
weight of gauss points
Definition: GaussM.f90:10
m_precheck_lib_2d
This module provides function to check input data of 2d static analysis.
Definition: precheck_LIB_2d.f90:6