Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
v3d10_nl.f90
Go to the documentation of this file.
1 !*********************************************************************
2 !* Illinois Open Source License *
3 !* *
4 !* University of Illinois/NCSA *
5 !* Open Source License *
6 !* *
7 !* Copyright@2008, University of Illinois. All rights reserved. *
8 !* *
9 !* Developed by: *
10 !* *
11 !* Center for Simulation of Advanced Rockets *
12 !* *
13 !* University of Illinois *
14 !* *
15 !* www.csar.uiuc.edu *
16 !* *
17 !* Permission is hereby granted, free of charge, to any person *
18 !* obtaining a copy of this software and associated documentation *
19 !* files (the "Software"), to deal with the Software without *
20 !* restriction, including without limitation the rights to use, *
21 !* copy, modify, merge, publish, distribute, sublicense, and/or *
22 !* sell copies of the Software, and to permit persons to whom the *
23 !* Software is furnished to do so, subject to the following *
24 !* conditions: *
25 !* *
26 !* *
27 !* @ Redistributions of source code must retain the above copyright *
28 !* notice, this list of conditions and the following disclaimers. *
29 !* *
30 !* @ Redistributions in binary form must reproduce the above *
31 !* copyright notice, this list of conditions and the following *
32 !* disclaimers in the documentation and/or other materials *
33 !* provided with the distribution. *
34 !* *
35 !* @ Neither the names of the Center for Simulation of Advanced *
36 !* Rockets, the University of Illinois, nor the names of its *
37 !* contributors may be used to endorse or promote products derived *
38 !* from this Software without specific prior written permission. *
39 !* *
40 !* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
41 !* EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
42 !* OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
43 !* NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
44 !* COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
45 !* LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
46 !* ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
47 !* USE OR OTHER DEALINGS WITH THE SOFTWARE. *
48 !*********************************************************************
49 !* Please acknowledge The University of Illinois Center for *
50 !* Simulation of Advanced Rockets in works and publications *
51 !* resulting from this software or its derivatives. *
52 !*********************************************************************
53 SUBROUTINE v3d10_nl(coor,matcstet,lmcstet,R_in,d,ci,&
54  s11,s22,s33,s12,s23,s13, &
55  numnp,nstart,nend,numcstet,numat_vol)
56 
57 !!****f* Rocfrac/Rocfrac/Source/v3d10_nl.f90
58 !!
59 !! NAME
60 !! v3d10_nl
61 !!
62 !! FUNCTION
63 !!
64 !! Performs displacement based computations for Volumetric 3D,
65 !! 10-node Quadratic Linear Elastic Tetrahedron with quadratic
66 !! interpolation functions. (linear strain tetrahedra).
67 !! Large Deformation.
68 !!
69 !! INPUTS
70 !!
71 !! NumNP -- Number of nodes
72 !! numcstet -- Number of elements
73 !! Coor -- number of coordinates
74 !! Matcstet -- Material id
75 !! d -- nodal displacement
76 !! ci -- compliance matrix
77 !!
78 !! stored in the follow format:
79 !!
80 !! [ ci(1) ci(2) ci(4) ]
81 !! [ ci(3) ci(5) ]
82 !! [ ci(6) ]
83 !! [ ci(7) ]
84 !! [ ci(8) ]
85 !! [ ci(9) ]
86 !!
87 !! lmcstet -- Nodal connectivity
88 !! nstart, nend -- element beginning and end loop counter
89 !! numat_vol -- number of materials
90 !!
91 !! OUTPUT
92 !!
93 !! R_in -- internal force vector
94 !! S11,S22,S33,S12,S23,S13 -- Stresses at a guass point
95 !!
96 !!****
97 
98 
99  IMPLICIT NONE
100 !-----Global variables
101  INTEGER :: numnp ! number of nodes
102  INTEGER :: numat_vol ! number of volumetric materials
103  INTEGER :: numcstet ! number of LSTets
104 !-- coordinate array
105  REAL*8, DIMENSION(1:3,1:numnp) :: coor
106 !-- elastic stiffness consts
107  REAL*8, DIMENSION(1:9,1:numat_vol) :: ci
108 !-- internal force
109  REAL*8, DIMENSION(1:3*numnp) :: r_in
110 !-- displacement vector
111  REAL*8, DIMENSION(1:3*numnp) :: d
112 !-- CSTet stress
113  REAL*8, DIMENSION(1:4,1:numcstet) :: s11, s22, s33, s12, s23, s13
114 !-- connectivity table for CSTet
115  INTEGER, DIMENSION(1:10,1:numcstet) :: lmcstet
116 !-- mat number for CSTet element
117  INTEGER, DIMENSION(1:numcstet) :: matcstet
118 !---- Local variables
119 !-- node numbers
120  INTEGER :: n1,n2,n3,n4,n5,n6,n7,n8,n9,n10
121 !-- x, y, and z displacements of nodes
122  REAL*8 :: u1,u2,u3,u4,u5,u6,u7,u8,u9,u10
123  REAL*8 :: v1,v2,v3,v4,v5,v6,v7,v8,v9,v10
124  REAL*8 :: w1,w2,w3,w4,w5,w6,w7,w8,w9,w10
125 !-- 6*volume and the volume
126  REAL*8 :: vx6,vx6inv
127 !-- spacial derivatives
128  REAL*8 :: b1,b2,b3,b4,b5,b6,b7,b8,b9,b10
129  REAL*8 :: b11,b12,b13,b14,b15,b16,b17,b18,b19,b20
130  REAL*8 :: b21,b22,b23,b24,b25,b26,b27,b28,b29,b30
131 !-- partial derivatives of the displacement
132  REAL*8 :: dudx,dvdy,dwdz,dudy,dvdx,dvdz,dwdy,dudz,dwdx
133 !-- strains
134  REAL*8 :: e11,e22,e33,e12,e23,e13
135 !-- coordinate holding variable
136  REAL*8 :: x1,x2,x3,x4,x5,x6,x7,x8,x9,x10
137  REAL*8 :: y1,y2,y3,y4,y5,y6,y7,y8,y9,y10
138  REAL*8 :: z1,z2,z3,z4,z5,z6,z7,z8,z9,z10
139 !-- dummy and counters
140  INTEGER :: i,j,nstart,nend
141  REAL*8 :: aux1,aux2,aux3,aux4,aux5,aux6,aux7,aux8,aux9,aux10,aux11,aux12
142 !-- partial internal force
143  REAL*8 :: r1,r2,r3,r4,r5,r6,r7,r8,r9,r10,r11,r12,r13,r14,r15,r16,r17,r18
144  REAL*8 :: r19,r20,r21,r22,r23,r24,r25,r26,r27,r28,r29,r30
145  REAL*8 :: g1, g2, g3, g4
146  REAL*8 :: xn1, xn2, xn3, xn4
147 !-- Coordinate subtractions
148  REAL*8 :: x14, x24, x34, y14, y24, y34, z14, z24, z34
149 !--
150  REAL*8 :: c11, c21, c31
151  INTEGER :: k1n1,k1n2,k1n3,k1n4,k1n5,k1n6,k1n7,k1n8,k1n9,k1n10
152  INTEGER :: k2n1,k2n2,k2n3,k2n4,k2n5,k2n6,k2n7,k2n8,k2n9,k2n10
153  INTEGER :: k3n1,k3n2,k3n3,k3n4,k3n5,k3n6,k3n7,k3n8,k3n9,k3n10
154 
155  DO i = nstart, nend
156 
157  j = matcstet(i)
158 
159  n1 = lmcstet(1,i)
160  n2 = lmcstet(2,i)
161  n3 = lmcstet(3,i)
162  n4 = lmcstet(4,i)
163  n5 = lmcstet(5,i)
164  n6 = lmcstet(6,i)
165  n7 = lmcstet(7,i)
166  n8 = lmcstet(8,i)
167  n9 = lmcstet(9,i)
168  n10 = lmcstet(10,i)
169 
170  k3n1 = 3*n1
171  k3n2 = 3*n2
172  k3n3 = 3*n3
173  k3n4 = 3*n4
174  k3n5 = 3*n5
175  k3n6 = 3*n6
176  k3n7 = 3*n7
177  k3n8 = 3*n8
178  k3n9 = 3*n9
179  k3n10 = 3*n10
180 
181  k2n1 = k3n1 - 1
182  k2n2 = k3n2 - 1
183  k2n3 = k3n3 - 1
184  k2n4 = k3n4 - 1
185  k2n5 = k3n5 - 1
186  k2n6 = k3n6 - 1
187  k2n7 = k3n7 - 1
188  k2n8 = k3n8 - 1
189  k2n9 = k3n9 - 1
190  k2n10 = k3n10 - 1
191 
192  k1n1 = k3n1 - 2
193  k1n2 = k3n2 - 2
194  k1n3 = k3n3 - 2
195  k1n4 = k3n4 - 2
196  k1n5 = k3n5 - 2
197  k1n6 = k3n6 - 2
198  k1n7 = k3n7 - 2
199  k1n8 = k3n8 - 2
200  k1n9 = k3n9 - 2
201  k1n10 = k3n10 - 2
202  ! k#n# dummy variables replaces:
203  u1 = d(k1n1) ! (3*n1 -2)
204  u2 = d(k1n2) ! (3*n2 -2)
205  u3 = d(k1n3) ! (3*n3 -2)
206  u4 = d(k1n4) ! (3*n4 -2)
207  u5 = d(k1n5) ! (3*n5 -2)
208  u6 = d(k1n6) ! (3*n6 -2)
209  u7 = d(k1n7) ! (3*n7 -2)
210  u8 = d(k1n8) ! (3*n8 -2)
211  u9 = d(k1n9) ! (3*n9 -2)
212  u10 = d(k1n10) ! (3*n10-2)
213  v1 = d(k2n1) ! (3*n1 -1)
214  v2 = d(k2n2) ! (3*n2 -1)
215  v3 = d(k2n3) ! (3*n3 -1)
216  v4 = d(k2n4) ! (3*n4 -1)
217  v5 = d(k2n5) ! (3*n5 -1)
218  v6 = d(k2n6) ! (3*n6 -1)
219  v7 = d(k2n7) ! (3*n7 -1)
220  v8 = d(k2n8) ! (3*n8 -1)
221  v9 = d(k2n9) ! (3*n9 -1)
222  v10 = d(k2n10) ! (3*n10-1)
223  w1 = d(k3n1) ! (3*n1)
224  w2 = d(k3n2) ! (3*n2)
225  w3 = d(k3n3) ! (3*n3)
226  w4 = d(k3n4) ! (3*n4)
227  w5 = d(k3n5) ! (3*n5)
228  w6 = d(k3n6) ! (3*n6)
229  w7 = d(k3n7) ! (3*n7)
230  w8 = d(k3n8) ! (3*n8)
231  w9 = d(k3n9) ! (3*n9)
232  w10 = d(k3n10) ! (3*n10)
233 
234  x1 = coor(1,n1)
235  x2 = coor(1,n2)
236  x3 = coor(1,n3)
237  x4 = coor(1,n4)
238  y1 = coor(2,n1)
239  y2 = coor(2,n2)
240  y3 = coor(2,n3)
241  y4 = coor(2,n4)
242  z1 = coor(3,n1)
243  z2 = coor(3,n2)
244  z3 = coor(3,n3)
245  z4 = coor(3,n4)
246 
247  vx6 = x2*y3*z4 - x2*y4*z3 - y2*x3*z4 + y2*x4*z3 + z2*x3*y4 &
248  - z2*x4*y3 - x1*y3*z4 + x1*y4*z3 + x1*y2*z4 - x1*y2*z3 &
249  - x1*z2*y4 + x1*z2*y3 + y1*x3*z4 - y1*x4*z3 - y1*x2*z4 &
250  + y1*x2*z3 + y1*z2*x4 - y1*z2*x3 - z1*x3*y4 + z1*x4*y3 &
251  + z1*x2*y4 - z1*x2*y3 - z1*y2*x4 + z1*y2*x3
252 
253  vx6inv = 1.d0 / vx6
254 
255  aux1 = -(y3*z4 - y4*z3 - y2*z4 + y2*z3 + z2*y4 - z2*y3)
256  aux2 = (x3*z4 - x4*z3 - x2*z4 + x2*z3 + z2*x4 - z2*x3)
257  aux3 = -(x3*y4 - x4*y3 - x2*y4 + x2*y3 + y2*x4 - y2*x3)
258  aux4 = (y3*z4 - y4*z3 - y1*z4 + y1*z3 + z1*y4 - z1*y3)
259  aux5 = -(x3*z4 - x4*z3 - x1*z4 + x1*z3 + z1*x4 - z1*x3)
260  aux6 = (x3*y4 - x4*y3 - x1*y4 + x1*y3 + y1*x4 - y1*x3)
261  aux7 = -(y2*z4 - z2*y4 - y1*z4 + y1*z2 + z1*y4 - z1*y2)
262  aux8 = (x2*z4 - z2*x4 - x1*z4 + x1*z2 + z1*x4 - z1*x2)
263  aux9 = -(x2*y4 - y2*x4 - x1*y4 + x1*y2 + y1*x4 - y1*x2)
264  aux10 = (y2*z3 - z2*y3 - y1*z3 + y1*z2 + z1*y3 - z1*y2)
265  aux11 =-(x2*z3 - z2*x3 - x1*z3 + x1*z2 + z1*x3 - z1*x2)
266  aux12 = (x2*y3 - y2*x3 - x1*y3 + x1*y2 + y1*x3 - y1*x2)
267 
268 !-- Gauss Point 1
269 
270  g1 = 0.58541020d0
271  g2 = 0.13819660d0
272  g3 = 0.13819660d0
273  g4 = 0.13819660d0
274 
275  xn1 = (4.d0*g1-1.d0) ! derivative of shape function
276  xn2 = (4.d0*g2-1.d0) ! dN_i/dzeta_i
277  xn3 = (4.d0*g3-1.d0)
278  xn4 = (4.d0*g4-1.d0)
279 ! xN5 = 4.d0*g1*g2
280 ! xN6 = 4.d0*g2*g3
281 ! xN7 = 4.d0*g3*g1
282 ! xN8 = 4.d0*g1*g4
283 ! xN9 = 4.d0*g2*g4
284 ! xN10= 4.d0*g3*g4
285 
286  b1 = aux1*xn1
287  b2 = aux2*xn1
288  b3 = aux3*xn1
289  b4 = aux4*xn2
290  b5 = aux5*xn2
291  b6 = aux6*xn2
292  b7 = aux7*xn3
293  b8 = aux8*xn3
294  b9 = aux9*xn3
295  b10 = aux10*xn4
296  b11 = aux11*xn4
297  b12 = aux12*xn4
298 
299  b13 = 4.d0*(g2*aux1 + g1*aux4)
300  b14 = 4.d0*(g2*aux2 + g1*aux5)
301  b15 = 4.d0*(g2*aux3 + g1*aux6)
302 
303  b16 = 4.d0*(g3*aux4 + g2*aux7)
304  b17 = 4.d0*(g3*aux5 + g2*aux8)
305  b18 = 4.d0*(g3*aux6 + g2*aux9)
306 
307  b19 = 4.d0*(g1*aux7 + g3*aux1)
308  b20 = 4.d0*(g1*aux8 + g3*aux2)
309  b21 = 4.d0*(g1*aux9 + g3*aux3)
310 
311  b22 = 4.d0*(g4*aux1 + g1*aux10)
312  b23 = 4.d0*(g4*aux2 + g1*aux11)
313  b24 = 4.d0*(g4*aux3 + g1*aux12)
314 
315  b25 = 4.d0*(g4*aux4 + g2*aux10)
316  b26 = 4.d0*(g4*aux5 + g2*aux11)
317  b27 = 4.d0*(g4*aux6 + g2*aux12)
318 
319  b28 = 4.d0*(g4*aux7 + g3*aux10)
320  b29 = 4.d0*(g4*aux8 + g3*aux11)
321  b30 = 4.d0*(g4*aux9 + g3*aux12)
322 
323 
324 !-----Calculate displacement gradient (H)
325  dudx = (b1*u1 + b4*u2 + b7*u3 + b10*u4 + b13*u5 + b16*u6 + b19*u7 + b22*u8 + b25*u9 + b28*u10)*vx6inv
326  dvdy = (b2*v1 + b5*v2 + b8*v3 + b11*v4 + b14*v5 + b17*v6 + b20*v7 + b23*v8 + b26*v9 + b29*v10)*vx6inv
327  dwdz = (b3*w1 + b6*w2 + b9*w3 + b12*w4 + b15*w5 + b18*w6 + b21*w7 + b24*w8 + b27*w9 + b30*w10)*vx6inv
328  dudy = (b2*u1 + b5*u2 + b8*u3 + b11*u4 + b14*u5 + b17*u6 + b20*u7 + b23*u8 + b26*u9 + b29*u10)*vx6inv
329  dvdx = (b1*v1 + b4*v2 + b7*v3 + b10*v4 + b13*v5 + b16*v6 + b19*v7 + b22*v8 + b25*v9 + b28*v10)*vx6inv
330  dvdz = (b3*v1 + b6*v2 + b9*v3 + b12*v4 + b15*v5 + b18*v6 + b21*v7 + b24*v8 + b27*v9 + b30*v10)*vx6inv
331  dwdy = (b2*w1 + b5*w2 + b8*w3 + b11*w4 + b14*w5 + b17*w6 + b20*w7 + b23*w8 + b26*w9 + b29*w10)*vx6inv
332  dudz = (b3*u1 + b6*u2 + b9*u3 + b12*u4 + b15*u5 + b18*u6 + b21*u7 + b24*u8 + b27*u9 + b30*u10)*vx6inv
333  dwdx = (b1*w1 + b4*w2 + b7*w3 + b10*w4 + b13*w5 + b16*w6 + b19*w7 + b22*w8 + b25*w9 + b28*w10)*vx6inv
334 !
335 ! Green-Lagrange Strain (E)
336 !
337  e11 = dudx + 0.5d0*(dudx*dudx + dvdx*dvdx + dwdx*dwdx)
338  e22 = dvdy + 0.5d0*(dudy*dudy + dvdy*dvdy + dwdy*dwdy)
339  e33 = dwdz + 0.5d0*(dudz*dudz + dvdz*dvdz + dwdz*dwdz)
340  e12 = dudy + dvdx + ( dudx*dudy + dvdx*dvdy + dwdx*dwdy )
341  e23 = dvdz + dwdy + ( dudy*dudz + dvdy*dvdz + dwdy*dwdz )
342  e13 = dwdx + dudz + ( dudz*dudx + dvdz*dvdx + dwdz*dwdx )
343 !
344 ! Second Piola-Kirchhoff stress (S)
345 !
346  s11(1,i) = e11*ci(1,j) + e22*ci(2,j) + e33*ci(4,j)
347  s22(1,i) = e11*ci(2,j) + e22*ci(3,j) + e33*ci(5,j)
348  s33(1,i) = e11*ci(4,j) + e22*ci(5,j) + e33*ci(6,j)
349  s12(1,i) = e12*ci(7,j)
350  s23(1,i) = e23*ci(8,j)
351  s13(1,i) = e13*ci(9,j)
352 
353  r1 = &
354  ( s11(1,i)*b1*(1.d0+dudx) + s22(1,i)*b2*dudy + s33(1,i)*b3*dudz &
355  + s12(1,i)*( b2*(1.d0+dudx) + b1*dudy ) &
356  + s23(1,i)*( b3*dudy + b2*dudz ) &
357  + s13(1,i)*( b3*(1.d0+dudx) + b1*dudz ) )
358  r2 = &
359  ( s11(1,i)*b1*dvdx + s22(1,i)*b2*(1.d0+dvdy) + s33(1,i)*b3*dvdz &
360  + s12(1,i)*( b1*(1.d0+dvdy) + b2*dvdx ) &
361  + s23(1,i)*( b3*(1.d0+dvdy) + b2*dvdz ) &
362  + s13(1,i)*( b3*dvdx + b1*dvdz ) )
363  r3 = &
364  ( s11(1,i)*b1*dwdx + s22(1,i)*b2*dwdy + s33(1,i)*b3*(1.d0+dwdz) &
365  + s12(1,i)*( b2*dwdx + b1*dwdy ) &
366  + s23(1,i)*( b3*dwdy + b2*(1.d0 + dwdz) ) &
367  + s13(1,i)*( b3*dwdx + b1*(1.d0 + dwdz) ) )
368 
369  r4 = &
370  ( s11(1,i)*b4*(1.d0+dudx) + s22(1,i)*b5*dudy + s33(1,i)*b6*dudz &
371  + s12(1,i)*( b5*(1.d0+dudx) + b4*dudy ) &
372  + s23(1,i)*( b6*dudy + b5*dudz ) &
373  + s13(1,i)*( b6*(1.d0+dudx) + b4*dudz ) )
374  r5 = &
375  ( s11(1,i)*b4*dvdx + s22(1,i)*b5*(1.d0+dvdy) + s33(1,i)*b6*dvdz &
376  + s12(1,i)*( b4*(1.d0+dvdy) + b5*dvdx ) &
377  + s23(1,i)*( b6*(1.d0+dvdy) + b5*dvdz ) &
378  + s13(1,i)*( b6*dvdx + b4*dvdz ) )
379  r6 = &
380  ( s11(1,i)*b4*dwdx + s22(1,i)*b5*dwdy + s33(1,i)*b6*(1.d0+dwdz) &
381  + s12(1,i)*( b5*dwdx + b4*dwdy ) &
382  + s23(1,i)*( b6*dwdy + b5*(1.d0 + dwdz) ) &
383  + s13(1,i)*( b6*dwdx + b4*(1.d0 + dwdz) ) )
384 
385  r7 = &
386  ( s11(1,i)*b7*(1.d0+dudx) + s22(1,i)*b8*dudy + s33(1,i)*b9*dudz &
387  + s12(1,i)*( b8*(1.d0+dudx) + b7*dudy ) &
388  + s23(1,i)*( b9*dudy + b8*dudz ) &
389  + s13(1,i)*( b9*(1.d0+dudx) + b7*dudz ) )
390  r8 = &
391  ( s11(1,i)*b7*dvdx + s22(1,i)*b8*(1.d0+dvdy) + s33(1,i)*b9*dvdz &
392  + s12(1,i)*( b7*(1.d0+dvdy) + b8*dvdx ) &
393  + s23(1,i)*( b9*(1.d0+dvdy) + b8*dvdz ) &
394  + s13(1,i)*( b9*dvdx + b7*dvdz ) )
395  r9 = &
396  ( s11(1,i)*b7*dwdx + s22(1,i)*b8*dwdy + s33(1,i)*b9*(1.d0+dwdz) &
397  + s12(1,i)*( b8*dwdx + b7*dwdy ) &
398  + s23(1,i)*( b9*dwdy + b8*(1.d0 + dwdz) ) &
399  + s13(1,i)*( b9*dwdx + b7*(1.d0 + dwdz) ) )
400 
401  r10 = &
402  ( s11(1,i)*b10*(1.d0+dudx) + s22(1,i)*b11*dudy+s33(1,i)*b12*dudz &
403  + s12(1,i)*( b11*(1.d0+dudx) + b10*dudy ) &
404  + s23(1,i)*( b12*dudy + b11*dudz ) &
405  + s13(1,i)*( b12*(1.d0+dudx) + b10*dudz ) )
406  r11 = &
407  ( s11(1,i)*b10*dvdx + s22(1,i)*b11*(1.d0+dvdy)+s33(1,i)*b12*dvdz &
408  + s12(1,i)*( b10*(1.d0+dvdy) + b11*dvdx ) &
409  + s23(1,i)*( b12*(1.d0+dvdy) + b11*dvdz ) &
410  + s13(1,i)*( b12*dvdx + b10*dvdz ) )
411  r12 = &
412  ( s11(1,i)*b10*dwdx + s22(1,i)*b11*dwdy+s33(1,i)*b12*(1.d0+dwdz) &
413  + s12(1,i)*( b11*dwdx + b10*dwdy ) &
414  + s23(1,i)*( b12*dwdy + b11*(1.d0 + dwdz) ) &
415  + s13(1,i)*( b12*dwdx + b10*(1.d0 + dwdz) ) )
416 
417  r13 = &
418  ( s11(1,i)*b13*(1.d0+dudx) + s22(1,i)*b14*dudy + s33(1,i)*b15*dudz &
419  + s12(1,i)*( b14*(1.d0+dudx) + b13*dudy ) &
420  + s23(1,i)*( b15*dudy + b14*dudz ) &
421  + s13(1,i)*( b15*(1.d0+dudx) + b13*dudz ) )
422  r14 = &
423  ( s11(1,i)*b13*dvdx + s22(1,i)*b14*(1.d0+dvdy) + s33(1,i)*b15*dvdz &
424  + s12(1,i)*( b13*(1.d0+dvdy) + b14*dvdx ) &
425  + s23(1,i)*( b15*(1.d0+dvdy) + b14*dvdz ) &
426  + s13(1,i)*( b15*dvdx + b13*dvdz ) )
427  r15 = &
428  ( s11(1,i)*b13*dwdx + s22(1,i)*b14*dwdy + s33(1,i)*b15*(1.d0+dwdz) &
429  + s12(1,i)*( b14*dwdx + b13*dwdy ) &
430  + s23(1,i)*( b15*dwdy + b14*(1.d0 + dwdz) ) &
431  + s13(1,i)*( b15*dwdx + b13*(1.d0 + dwdz) ) )
432 
433  r16 = &
434  ( s11(1,i)*b16*(1.d0+dudx) + s22(1,i)*b17*dudy + s33(1,i)*b18*dudz &
435  + s12(1,i)*( b17*(1.d0+dudx) + b16*dudy ) &
436  + s23(1,i)*( b18*dudy + b17*dudz ) &
437  + s13(1,i)*( b18*(1.d0+dudx) + b16*dudz ) )
438  r17 = &
439  ( s11(1,i)*b16*dvdx + s22(1,i)*b17*(1.d0+dvdy) + s33(1,i)*b18*dvdz &
440  + s12(1,i)*( b16*(1.d0+dvdy) + b17*dvdx ) &
441  + s23(1,i)*( b18*(1.d0+dvdy) + b17*dvdz ) &
442  + s13(1,i)*( b18*dvdx + b16*dvdz ) )
443  r18 = &
444  ( s11(1,i)*b16*dwdx + s22(1,i)*b17*dwdy + s33(1,i)*b18*(1.d0+dwdz) &
445  + s12(1,i)*( b17*dwdx + b16*dwdy ) &
446  + s23(1,i)*( b18*dwdy + b17*(1.d0 + dwdz) ) &
447  + s13(1,i)*( b18*dwdx + b16*(1.d0 + dwdz) ) )
448 
449  r19 = &
450  ( s11(1,i)*b19*(1.d0+dudx) + s22(1,i)*b20*dudy + s33(1,i)*b21*dudz &
451  + s12(1,i)*( b20*(1.d0+dudx) + b19*dudy ) &
452  + s23(1,i)*( b21*dudy + b20*dudz ) &
453  + s13(1,i)*( b21*(1.d0+dudx) + b19*dudz ) )
454  r20 = &
455  ( s11(1,i)*b19*dvdx + s22(1,i)*b20*(1.d0+dvdy) + s33(1,i)*b21*dvdz &
456  + s12(1,i)*( b19*(1.d0+dvdy) + b20*dvdx ) &
457  + s23(1,i)*( b21*(1.d0+dvdy) + b20*dvdz ) &
458  + s13(1,i)*( b21*dvdx + b19*dvdz ) )
459  r21 = &
460  ( s11(1,i)*b19*dwdx + s22(1,i)*b20*dwdy + s33(1,i)*b21*(1.d0+dwdz) &
461  + s12(1,i)*( b20*dwdx + b19*dwdy ) &
462  + s23(1,i)*( b21*dwdy + b20*(1.d0 + dwdz) ) &
463  + s13(1,i)*( b21*dwdx + b19*(1.d0 + dwdz) ) )
464 
465  r22 = &
466  ( s11(1,i)*b22*(1.d0+dudx) + s22(1,i)*b23*dudy+s33(1,i)*b24*dudz &
467  + s12(1,i)*( b23*(1.d0+dudx) + b22*dudy ) &
468  + s23(1,i)*( b24*dudy + b23*dudz ) &
469  + s13(1,i)*( b24*(1.d0+dudx) + b22*dudz ) )
470  r23 = &
471  ( s11(1,i)*b22*dvdx + s22(1,i)*b23*(1.d0+dvdy)+s33(1,i)*b24*dvdz &
472  + s12(1,i)*( b22*(1.d0+dvdy) + b23*dvdx ) &
473  + s23(1,i)*( b24*(1.d0+dvdy) + b23*dvdz ) &
474  + s13(1,i)*( b24*dvdx + b22*dvdz ) )
475  r24 = &
476  ( s11(1,i)*b22*dwdx + s22(1,i)*b23*dwdy+s33(1,i)*b24*(1.d0+dwdz) &
477  + s12(1,i)*( b23*dwdx + b22*dwdy ) &
478  + s23(1,i)*( b24*dwdy + b23*(1.d0 + dwdz) ) &
479  + s13(1,i)*( b24*dwdx + b22*(1.d0 + dwdz) ) )
480 
481  r25 = &
482  ( s11(1,i)*b25*(1.d0+dudx) + s22(1,i)*b26*dudy+s33(1,i)*b27*dudz &
483  + s12(1,i)*( b26*(1.d0+dudx) + b25*dudy ) &
484  + s23(1,i)*( b27*dudy + b26*dudz ) &
485  + s13(1,i)*( b27*(1.d0+dudx) + b25*dudz ) )
486  r26 = &
487  ( s11(1,i)*b25*dvdx + s22(1,i)*b26*(1.d0+dvdy)+s33(1,i)*b27*dvdz &
488  + s12(1,i)*( b25*(1.d0+dvdy) + b26*dvdx ) &
489  + s23(1,i)*( b27*(1.d0+dvdy) + b26*dvdz ) &
490  + s13(1,i)*( b27*dvdx + b25*dvdz ) )
491  r27 = &
492  ( s11(1,i)*b25*dwdx + s22(1,i)*b26*dwdy+s33(1,i)*b27*(1.d0+dwdz) &
493  + s12(1,i)*( b26*dwdx + b25*dwdy ) &
494  + s23(1,i)*( b27*dwdy + b26*(1.d0 + dwdz) ) &
495  + s13(1,i)*( b27*dwdx + b25*(1.d0 + dwdz) ) )
496 
497  r28 = &
498  ( s11(1,i)*b28*(1.d0+dudx) + s22(1,i)*b29*dudy+s33(1,i)*b30*dudz &
499  + s12(1,i)*( b29*(1.d0+dudx) + b28*dudy ) &
500  + s23(1,i)*( b30*dudy + b29*dudz ) &
501  + s13(1,i)*( b30*(1.d0+dudx) + b28*dudz ) )
502  r29 = &
503  ( s11(1,i)*b28*dvdx + s22(1,i)*b29*(1.d0+dvdy)+s33(1,i)*b30*dvdz &
504  + s12(1,i)*( b28*(1.d0+dvdy) + b29*dvdx ) &
505  + s23(1,i)*( b30*(1.d0+dvdy) + b29*dvdz ) &
506  + s13(1,i)*( b30*dvdx + b28*dvdz ) )
507  r30 = &
508  ( s11(1,i)*b28*dwdx + s22(1,i)*b29*dwdy+s33(1,i)*b30*(1.d0+dwdz) &
509  + s12(1,i)*( b29*dwdx + b28*dwdy ) &
510  + s23(1,i)*( b30*dwdy + b29*(1.d0 + dwdz) ) &
511  + s13(1,i)*( b30*dwdx + b28*(1.d0 + dwdz) ) )
512 
513 !-- Gauss Point 2
514 
515  g1 = 0.13819660d0
516  g2 = 0.58541020d0
517  g3 = 0.13819660d0
518  g4 = 0.13819660d0
519 
520  xn1 = (4.d0*g1-1.d0) ! derivative of shape function
521  xn2 = (4.d0*g2-1.d0) ! dN_i/dzeta_i
522  xn3 = (4.d0*g3-1.d0)
523  xn4 = (4.d0*g4-1.d0)
524 ! xN5 = 4.d0*g1*g2
525 ! xN6 = 4.d0*g2*g3
526 ! xN7 = 4.d0*g3*g1
527 ! xN8 = 4.d0*g1*g4
528 ! xN9 = 4.d0*g2*g4
529 ! xN10= 4.d0*g3*g4
530 
531  b1 = aux1*xn1
532  b2 = aux2*xn1
533  b3 = aux3*xn1
534  b4 = aux4*xn2
535  b5 = aux5*xn2
536  b6 = aux6*xn2
537  b7 = aux7*xn3
538  b8 = aux8*xn3
539  b9 = aux9*xn3
540  b10 = aux10*xn4
541  b11 = aux11*xn4
542  b12 = aux12*xn4
543 
544  b13 = 4.d0*(g2*aux1 + g1*aux4)
545  b14 = 4.d0*(g2*aux2 + g1*aux5)
546  b15 = 4.d0*(g2*aux3 + g1*aux6)
547 
548  b16 = 4.d0*(g3*aux4 + g2*aux7)
549  b17 = 4.d0*(g3*aux5 + g2*aux8)
550  b18 = 4.d0*(g3*aux6 + g2*aux9)
551 
552  b19 = 4.d0*(g1*aux7 + g3*aux1)
553  b20 = 4.d0*(g1*aux8 + g3*aux2)
554  b21 = 4.d0*(g1*aux9 + g3*aux3)
555 
556  b22 = 4.d0*(g4*aux1 + g1*aux10)
557  b23 = 4.d0*(g4*aux2 + g1*aux11)
558  b24 = 4.d0*(g4*aux3 + g1*aux12)
559 
560  b25 = 4.d0*(g4*aux4 + g2*aux10)
561  b26 = 4.d0*(g4*aux5 + g2*aux11)
562  b27 = 4.d0*(g4*aux6 + g2*aux12)
563 
564  b28 = 4.d0*(g4*aux7 + g3*aux10)
565  b29 = 4.d0*(g4*aux8 + g3*aux11)
566  b30 = 4.d0*(g4*aux9 + g3*aux12)
567 !-----Calculate displacement gradient (H)
568  dudx = (b1*u1 + b4*u2 + b7*u3 + b10*u4 + b13*u5 + b16*u6 + b19*u7 + b22*u8 + b25*u9 + b28*u10)*vx6inv
569  dvdy = (b2*v1 + b5*v2 + b8*v3 + b11*v4 + b14*v5 + b17*v6 + b20*v7 + b23*v8 + b26*v9 + b29*v10)*vx6inv
570  dwdz = (b3*w1 + b6*w2 + b9*w3 + b12*w4 + b15*w5 + b18*w6 + b21*w7 + b24*w8 + b27*w9 + b30*w10)*vx6inv
571  dudy = (b2*u1 + b5*u2 + b8*u3 + b11*u4 + b14*u5 + b17*u6 + b20*u7 + b23*u8 + b26*u9 + b29*u10)*vx6inv
572  dvdx = (b1*v1 + b4*v2 + b7*v3 + b10*v4 + b13*v5 + b16*v6 + b19*v7 + b22*v8 + b25*v9 + b28*v10)*vx6inv
573  dvdz = (b3*v1 + b6*v2 + b9*v3 + b12*v4 + b15*v5 + b18*v6 + b21*v7 + b24*v8 + b27*v9 + b30*v10)*vx6inv
574  dwdy = (b2*w1 + b5*w2 + b8*w3 + b11*w4 + b14*w5 + b17*w6 + b20*w7 + b23*w8 + b26*w9 + b29*w10)*vx6inv
575  dudz = (b3*u1 + b6*u2 + b9*u3 + b12*u4 + b15*u5 + b18*u6 + b21*u7 + b24*u8 + b27*u9 + b30*u10)*vx6inv
576  dwdx = (b1*w1 + b4*w2 + b7*w3 + b10*w4 + b13*w5 + b16*w6 + b19*w7 + b22*w8 + b25*w9 + b28*w10)*vx6inv
577 !
578 ! Green-Lagrange Strain (E)
579 !
580  e11 = dudx + 0.5d0*(dudx*dudx + dvdx*dvdx + dwdx*dwdx)
581  e22 = dvdy + 0.5d0*(dudy*dudy + dvdy*dvdy + dwdy*dwdy)
582  e33 = dwdz + 0.5d0*(dudz*dudz + dvdz*dvdz + dwdz*dwdz)
583  e12 = dudy + dvdx + ( dudx*dudy + dvdx*dvdy + dwdx*dwdy )
584  e23 = dvdz + dwdy + ( dudy*dudz + dvdy*dvdz + dwdy*dwdz )
585  e13 = dwdx + dudz + ( dudz*dudx + dvdz*dvdx + dwdz*dwdx )
586 !
587 ! Second Piola-Kirchhoff stress (S)
588 !
589  s11(2,i) = e11*ci(1,j) + e22*ci(2,j) + e33*ci(4,j)
590  s22(2,i) = e11*ci(2,j) + e22*ci(3,j) + e33*ci(5,j)
591  s33(2,i) = e11*ci(4,j) + e22*ci(5,j) + e33*ci(6,j)
592  s12(2,i) = e12*ci(7,j)
593  s23(2,i) = e23*ci(8,j)
594  s13(2,i) = e13*ci(9,j)
595 
596  r1 = r1 + &
597  ( s11(2,i)*b1*(1.d0+dudx) + s22(2,i)*b2*dudy + s33(2,i)*b3*dudz &
598  + s12(2,i)*( b2*(1.d0+dudx) + b1*dudy ) &
599  + s23(2,i)*( b3*dudy + b2*dudz ) &
600  + s13(2,i)*( b3*(1.d0+dudx) + b1*dudz ) )
601  r2 = r2 + &
602  ( s11(2,i)*b1*dvdx + s22(2,i)*b2*(1.d0+dvdy) + s33(2,i)*b3*dvdz &
603  + s12(2,i)*( b1*(1.d0+dvdy) + b2*dvdx ) &
604  + s23(2,i)*( b3*(1.d0+dvdy) + b2*dvdz ) &
605  + s13(2,i)*( b3*dvdx + b1*dvdz ) )
606  r3 = r3 + &
607  ( s11(2,i)*b1*dwdx + s22(2,i)*b2*dwdy + s33(2,i)*b3*(1.d0+dwdz) &
608  + s12(2,i)*( b2*dwdx + b1*dwdy ) &
609  + s23(2,i)*( b3*dwdy + b2*(1.d0 + dwdz) ) &
610  + s13(2,i)*( b3*dwdx + b1*(1.d0 + dwdz) ) )
611 
612  r4 = r4 + &
613  ( s11(2,i)*b4*(1.d0+dudx) + s22(2,i)*b5*dudy + s33(2,i)*b6*dudz &
614  + s12(2,i)*( b5*(1.d0+dudx) + b4*dudy ) &
615  + s23(2,i)*( b6*dudy + b5*dudz ) &
616  + s13(2,i)*( b6*(1.d0+dudx) + b4*dudz ) )
617  r5 = r5 + &
618  ( s11(2,i)*b4*dvdx + s22(2,i)*b5*(1.d0+dvdy) + s33(2,i)*b6*dvdz &
619  + s12(2,i)*( b4*(1.d0+dvdy) + b5*dvdx ) &
620  + s23(2,i)*( b6*(1.d0+dvdy) + b5*dvdz ) &
621  + s13(2,i)*( b6*dvdx + b4*dvdz ) )
622  r6 = r6 + &
623  ( s11(2,i)*b4*dwdx + s22(2,i)*b5*dwdy + s33(2,i)*b6*(1.d0+dwdz) &
624  + s12(2,i)*( b5*dwdx + b4*dwdy ) &
625  + s23(2,i)*( b6*dwdy + b5*(1.d0 + dwdz) ) &
626  + s13(2,i)*( b6*dwdx + b4*(1.d0 + dwdz) ) )
627 
628  r7 = r7 + &
629  ( s11(2,i)*b7*(1.d0+dudx) + s22(2,i)*b8*dudy + s33(2,i)*b9*dudz &
630  + s12(2,i)*( b8*(1.d0+dudx) + b7*dudy ) &
631  + s23(2,i)*( b9*dudy + b8*dudz ) &
632  + s13(2,i)*( b9*(1.d0+dudx) + b7*dudz ) )
633  r8 = r8 + &
634  ( s11(2,i)*b7*dvdx + s22(2,i)*b8*(1.d0+dvdy) + s33(2,i)*b9*dvdz &
635  + s12(2,i)*( b7*(1.d0+dvdy) + b8*dvdx ) &
636  + s23(2,i)*( b9*(1.d0+dvdy) + b8*dvdz ) &
637  + s13(2,i)*( b9*dvdx + b7*dvdz ) )
638  r9 = r9 + &
639  ( s11(2,i)*b7*dwdx + s22(2,i)*b8*dwdy + s33(2,i)*b9*(1.d0+dwdz) &
640  + s12(2,i)*( b8*dwdx + b7*dwdy ) &
641  + s23(2,i)*( b9*dwdy + b8*(1.d0 + dwdz) ) &
642  + s13(2,i)*( b9*dwdx + b7*(1.d0 + dwdz) ) )
643 
644  r10 = r10 + &
645  ( s11(2,i)*b10*(1.d0+dudx) + s22(2,i)*b11*dudy+s33(2,i)*b12*dudz &
646  + s12(2,i)*( b11*(1.d0+dudx) + b10*dudy ) &
647  + s23(2,i)*( b12*dudy + b11*dudz ) &
648  + s13(2,i)*( b12*(1.d0+dudx) + b10*dudz ) )
649  r11 = r11 + &
650  ( s11(2,i)*b10*dvdx + s22(2,i)*b11*(1.d0+dvdy)+s33(2,i)*b12*dvdz &
651  + s12(2,i)*( b10*(1.d0+dvdy) + b11*dvdx ) &
652  + s23(2,i)*( b12*(1.d0+dvdy) + b11*dvdz ) &
653  + s13(2,i)*( b12*dvdx + b10*dvdz ) )
654  r12 = r12 + &
655  ( s11(2,i)*b10*dwdx + s22(2,i)*b11*dwdy+s33(2,i)*b12*(1.d0+dwdz) &
656  + s12(2,i)*( b11*dwdx + b10*dwdy ) &
657  + s23(2,i)*( b12*dwdy + b11*(1.d0 + dwdz) ) &
658  + s13(2,i)*( b12*dwdx + b10*(1.d0 + dwdz) ) )
659 
660  r13 = r13 + &
661  ( s11(2,i)*b13*(1.d0+dudx) + s22(2,i)*b14*dudy + s33(2,i)*b15*dudz &
662  + s12(2,i)*( b14*(1.d0+dudx) + b13*dudy ) &
663  + s23(2,i)*( b15*dudy + b14*dudz ) &
664  + s13(2,i)*( b15*(1.d0+dudx) + b13*dudz ) )
665  r14 = r14 + &
666  ( s11(2,i)*b13*dvdx + s22(2,i)*b14*(1.d0+dvdy) + s33(2,i)*b15*dvdz &
667  + s12(2,i)*( b13*(1.d0+dvdy) + b14*dvdx ) &
668  + s23(2,i)*( b15*(1.d0+dvdy) + b14*dvdz ) &
669  + s13(2,i)*( b15*dvdx + b13*dvdz ) )
670  r15 = r15 + &
671  ( s11(2,i)*b13*dwdx + s22(2,i)*b14*dwdy + s33(2,i)*b15*(1.d0+dwdz) &
672  + s12(2,i)*( b14*dwdx + b13*dwdy ) &
673  + s23(2,i)*( b15*dwdy + b14*(1.d0 + dwdz) ) &
674  + s13(2,i)*( b15*dwdx + b13*(1.d0 + dwdz) ) )
675 
676  r16 = r16 + &
677  ( s11(2,i)*b16*(1.d0+dudx) + s22(2,i)*b17*dudy + s33(2,i)*b18*dudz &
678  + s12(2,i)*( b17*(1.d0+dudx) + b16*dudy ) &
679  + s23(2,i)*( b18*dudy + b17*dudz ) &
680  + s13(2,i)*( b18*(1.d0+dudx) + b16*dudz ) )
681  r17 = r17 + &
682  ( s11(2,i)*b16*dvdx + s22(2,i)*b17*(1.d0+dvdy) + s33(2,i)*b18*dvdz &
683  + s12(2,i)*( b16*(1.d0+dvdy) + b17*dvdx ) &
684  + s23(2,i)*( b18*(1.d0+dvdy) + b17*dvdz ) &
685  + s13(2,i)*( b18*dvdx + b16*dvdz ) )
686  r18 = r18 + &
687  ( s11(2,i)*b16*dwdx + s22(2,i)*b17*dwdy + s33(2,i)*b18*(1.d0+dwdz) &
688  + s12(2,i)*( b17*dwdx + b16*dwdy ) &
689  + s23(2,i)*( b18*dwdy + b17*(1.d0 + dwdz) ) &
690  + s13(2,i)*( b18*dwdx + b16*(1.d0 + dwdz) ) )
691 
692  r19 = r19 + &
693  ( s11(2,i)*b19*(1.d0+dudx) + s22(2,i)*b20*dudy + s33(2,i)*b21*dudz &
694  + s12(2,i)*( b20*(1.d0+dudx) + b19*dudy ) &
695  + s23(2,i)*( b21*dudy + b20*dudz ) &
696  + s13(2,i)*( b21*(1.d0+dudx) + b19*dudz ) )
697  r20 = r20 + &
698  ( s11(2,i)*b19*dvdx + s22(2,i)*b20*(1.d0+dvdy) + s33(2,i)*b21*dvdz &
699  + s12(2,i)*( b19*(1.d0+dvdy) + b20*dvdx ) &
700  + s23(2,i)*( b21*(1.d0+dvdy) + b20*dvdz ) &
701  + s13(2,i)*( b21*dvdx + b19*dvdz ) )
702  r21 = r21 + &
703  ( s11(2,i)*b19*dwdx + s22(2,i)*b20*dwdy + s33(2,i)*b21*(1.d0+dwdz) &
704  + s12(2,i)*( b20*dwdx + b19*dwdy ) &
705  + s23(2,i)*( b21*dwdy + b20*(1.d0 + dwdz) ) &
706  + s13(2,i)*( b21*dwdx + b19*(1.d0 + dwdz) ) )
707 
708  r22 = r22 + &
709  ( s11(2,i)*b22*(1.d0+dudx) + s22(2,i)*b23*dudy+s33(2,i)*b24*dudz &
710  + s12(2,i)*( b23*(1.d0+dudx) + b22*dudy ) &
711  + s23(2,i)*( b24*dudy + b23*dudz ) &
712  + s13(2,i)*( b24*(1.d0+dudx) + b22*dudz ) )
713  r23 = r23 + &
714  ( s11(2,i)*b22*dvdx + s22(2,i)*b23*(1.d0+dvdy)+s33(2,i)*b24*dvdz &
715  + s12(2,i)*( b22*(1.d0+dvdy) + b23*dvdx ) &
716  + s23(2,i)*( b24*(1.d0+dvdy) + b23*dvdz ) &
717  + s13(2,i)*( b24*dvdx + b22*dvdz ) )
718  r24 = r24 + &
719  ( s11(2,i)*b22*dwdx + s22(2,i)*b23*dwdy+s33(2,i)*b24*(1.d0+dwdz) &
720  + s12(2,i)*( b23*dwdx + b22*dwdy ) &
721  + s23(2,i)*( b24*dwdy + b23*(1.d0 + dwdz) ) &
722  + s13(2,i)*( b24*dwdx + b22*(1.d0 + dwdz) ) )
723 
724  r25 = r25 + &
725  ( s11(2,i)*b25*(1.d0+dudx) + s22(2,i)*b26*dudy+s33(2,i)*b27*dudz &
726  + s12(2,i)*( b26*(1.d0+dudx) + b25*dudy ) &
727  + s23(2,i)*( b27*dudy + b26*dudz ) &
728  + s13(2,i)*( b27*(1.d0+dudx) + b25*dudz ) )
729  r26 = r26 + &
730  ( s11(2,i)*b25*dvdx + s22(2,i)*b26*(1.d0+dvdy)+s33(2,i)*b27*dvdz &
731  + s12(2,i)*( b25*(1.d0+dvdy) + b26*dvdx ) &
732  + s23(2,i)*( b27*(1.d0+dvdy) + b26*dvdz ) &
733  + s13(2,i)*( b27*dvdx + b25*dvdz ) )
734  r27 = r27 + &
735  ( s11(2,i)*b25*dwdx + s22(2,i)*b26*dwdy+s33(2,i)*b27*(1.d0+dwdz) &
736  + s12(2,i)*( b26*dwdx + b25*dwdy ) &
737  + s23(2,i)*( b27*dwdy + b26*(1.d0 + dwdz) ) &
738  + s13(2,i)*( b27*dwdx + b25*(1.d0 + dwdz) ) )
739 
740  r28 = r28 + &
741  ( s11(2,i)*b28*(1.d0+dudx) + s22(2,i)*b29*dudy+s33(2,i)*b30*dudz &
742  + s12(2,i)*( b29*(1.d0+dudx) + b28*dudy ) &
743  + s23(2,i)*( b30*dudy + b29*dudz ) &
744  + s13(2,i)*( b30*(1.d0+dudx) + b28*dudz ) )
745  r29 = r29 + &
746  ( s11(2,i)*b28*dvdx + s22(2,i)*b29*(1.d0+dvdy)+s33(2,i)*b30*dvdz &
747  + s12(2,i)*( b28*(1.d0+dvdy) + b29*dvdx ) &
748  + s23(2,i)*( b30*(1.d0+dvdy) + b29*dvdz ) &
749  + s13(2,i)*( b30*dvdx + b28*dvdz ) )
750  r30 = r30 + &
751  ( s11(2,i)*b28*dwdx + s22(2,i)*b29*dwdy+s33(2,i)*b30*(1.d0+dwdz) &
752  + s12(2,i)*( b29*dwdx + b28*dwdy ) &
753  + s23(2,i)*( b30*dwdy + b29*(1.d0 + dwdz) ) &
754  + s13(2,i)*( b30*dwdx + b28*(1.d0 + dwdz) ) )
755 !-- Gauss Point 3
756 
757  g1 = 0.13819660d0
758  g2 = 0.13819660d0
759  g3 = 0.58541020d0
760  g4 = 0.13819660d0
761 
762  xn1 = (4.d0*g1-1.d0) ! derivative of shape function
763  xn2 = (4.d0*g2-1.d0) ! dN_i/dzeta_i
764  xn3 = (4.d0*g3-1.d0)
765  xn4 = (4.d0*g4-1.d0)
766  ! xN5 = 4.d0*g1*g2
767  ! xN6 = 4.d0*g2*g3
768  ! xN7 = 4.d0*g3*g1
769 ! xN8 = 4.d0*g1*g4
770 ! xN9 = 4.d0*g2*g4
771 ! xN10= 4.d0*g3*g4
772 
773  b1 = aux1*xn1
774  b2 = aux2*xn1
775  b3 = aux3*xn1
776  b4 = aux4*xn2
777  b5 = aux5*xn2
778  b6 = aux6*xn2
779  b7 = aux7*xn3
780  b8 = aux8*xn3
781  b9 = aux9*xn3
782  b10 = aux10*xn4
783  b11 = aux11*xn4
784  b12 = aux12*xn4
785 
786  b13 = 4.d0*(g2*aux1 + g1*aux4)
787  b14 = 4.d0*(g2*aux2 + g1*aux5)
788  b15 = 4.d0*(g2*aux3 + g1*aux6)
789 
790  b16 = 4.d0*(g3*aux4 + g2*aux7)
791  b17 = 4.d0*(g3*aux5 + g2*aux8)
792  b18 = 4.d0*(g3*aux6 + g2*aux9)
793 
794  b19 = 4.d0*(g1*aux7 + g3*aux1)
795  b20 = 4.d0*(g1*aux8 + g3*aux2)
796  b21 = 4.d0*(g1*aux9 + g3*aux3)
797 
798  b22 = 4.d0*(g4*aux1 + g1*aux10)
799  b23 = 4.d0*(g4*aux2 + g1*aux11)
800  b24 = 4.d0*(g4*aux3 + g1*aux12)
801 
802  b25 = 4.d0*(g4*aux4 + g2*aux10)
803  b26 = 4.d0*(g4*aux5 + g2*aux11)
804  b27 = 4.d0*(g4*aux6 + g2*aux12)
805 
806  b28 = 4.d0*(g4*aux7 + g3*aux10)
807  b29 = 4.d0*(g4*aux8 + g3*aux11)
808  b30 = 4.d0*(g4*aux9 + g3*aux12)
809 !-----Calculate displacement gradient (H)
810  dudx = (b1*u1 + b4*u2 + b7*u3 + b10*u4 + b13*u5 + b16*u6 + b19*u7 + b22*u8 + b25*u9 + b28*u10)*vx6inv
811  dvdy = (b2*v1 + b5*v2 + b8*v3 + b11*v4 + b14*v5 + b17*v6 + b20*v7 + b23*v8 + b26*v9 + b29*v10)*vx6inv
812  dwdz = (b3*w1 + b6*w2 + b9*w3 + b12*w4 + b15*w5 + b18*w6 + b21*w7 + b24*w8 + b27*w9 + b30*w10)*vx6inv
813  dudy = (b2*u1 + b5*u2 + b8*u3 + b11*u4 + b14*u5 + b17*u6 + b20*u7 + b23*u8 + b26*u9 + b29*u10)*vx6inv
814  dvdx = (b1*v1 + b4*v2 + b7*v3 + b10*v4 + b13*v5 + b16*v6 + b19*v7 + b22*v8 + b25*v9 + b28*v10)*vx6inv
815  dvdz = (b3*v1 + b6*v2 + b9*v3 + b12*v4 + b15*v5 + b18*v6 + b21*v7 + b24*v8 + b27*v9 + b30*v10)*vx6inv
816  dwdy = (b2*w1 + b5*w2 + b8*w3 + b11*w4 + b14*w5 + b17*w6 + b20*w7 + b23*w8 + b26*w9 + b29*w10)*vx6inv
817  dudz = (b3*u1 + b6*u2 + b9*u3 + b12*u4 + b15*u5 + b18*u6 + b21*u7 + b24*u8 + b27*u9 + b30*u10)*vx6inv
818  dwdx = (b1*w1 + b4*w2 + b7*w3 + b10*w4 + b13*w5 + b16*w6 + b19*w7 + b22*w8 + b25*w9 + b28*w10)*vx6inv
819 ! Green-Lagrange Strain (E)
820 !
821  e11 = dudx + 0.5d0*(dudx*dudx + dvdx*dvdx + dwdx*dwdx)
822  e22 = dvdy + 0.5d0*(dudy*dudy + dvdy*dvdy + dwdy*dwdy)
823  e33 = dwdz + 0.5d0*(dudz*dudz + dvdz*dvdz + dwdz*dwdz)
824  e12 = dudy + dvdx + ( dudx*dudy + dvdx*dvdy + dwdx*dwdy )
825  e23 = dvdz + dwdy + ( dudy*dudz + dvdy*dvdz + dwdy*dwdz )
826  e13 = dwdx + dudz + ( dudz*dudx + dvdz*dvdx + dwdz*dwdx )
827 !
828 ! Second Piola-Kirchhoff stress (S)
829 !
830  s11(3,i) = e11*ci(1,j) + e22*ci(2,j) + e33*ci(4,j)
831  s22(3,i) = e11*ci(2,j) + e22*ci(3,j) + e33*ci(5,j)
832  s33(3,i) = e11*ci(4,j) + e22*ci(5,j) + e33*ci(6,j)
833  s12(3,i) = e12*ci(7,j)
834  s23(3,i) = e23*ci(8,j)
835  s13(3,i) = e13*ci(9,j)
836 
837  r1 = r1 + &
838  ( s11(3,i)*b1*(1.d0+dudx) + s22(3,i)*b2*dudy + s33(3,i)*b3*dudz &
839  + s12(3,i)*( b2*(1.d0+dudx) + b1*dudy ) &
840  + s23(3,i)*( b3*dudy + b2*dudz ) &
841  + s13(3,i)*( b3*(1.d0+dudx) + b1*dudz ) )
842  r2 = r2 + &
843  ( s11(3,i)*b1*dvdx + s22(3,i)*b2*(1.d0+dvdy) + s33(3,i)*b3*dvdz &
844  + s12(3,i)*( b1*(1.d0+dvdy) + b2*dvdx ) &
845  + s23(3,i)*( b3*(1.d0+dvdy) + b2*dvdz ) &
846  + s13(3,i)*( b3*dvdx + b1*dvdz ) )
847  r3 = r3 + &
848  ( s11(3,i)*b1*dwdx + s22(3,i)*b2*dwdy + s33(3,i)*b3*(1.d0+dwdz) &
849  + s12(3,i)*( b2*dwdx + b1*dwdy ) &
850  + s23(3,i)*( b3*dwdy + b2*(1.d0 + dwdz) ) &
851  + s13(3,i)*( b3*dwdx + b1*(1.d0 + dwdz) ) )
852 
853  r4 = r4 + &
854  ( s11(3,i)*b4*(1.d0+dudx) + s22(3,i)*b5*dudy + s33(3,i)*b6*dudz &
855  + s12(3,i)*( b5*(1.d0+dudx) + b4*dudy ) &
856  + s23(3,i)*( b6*dudy + b5*dudz ) &
857  + s13(3,i)*( b6*(1.d0+dudx) + b4*dudz ) )
858  r5 = r5 + &
859  ( s11(3,i)*b4*dvdx + s22(3,i)*b5*(1.d0+dvdy) + s33(3,i)*b6*dvdz &
860  + s12(3,i)*( b4*(1.d0+dvdy) + b5*dvdx ) &
861  + s23(3,i)*( b6*(1.d0+dvdy) + b5*dvdz ) &
862  + s13(3,i)*( b6*dvdx + b4*dvdz ) )
863  r6 = r6 + &
864  ( s11(3,i)*b4*dwdx + s22(3,i)*b5*dwdy + s33(3,i)*b6*(1.d0+dwdz) &
865  + s12(3,i)*( b5*dwdx + b4*dwdy ) &
866  + s23(3,i)*( b6*dwdy + b5*(1.d0 + dwdz) ) &
867  + s13(3,i)*( b6*dwdx + b4*(1.d0 + dwdz) ) )
868 
869  r7 = r7 + &
870  ( s11(3,i)*b7*(1.d0+dudx) + s22(3,i)*b8*dudy + s33(3,i)*b9*dudz &
871  + s12(3,i)*( b8*(1.d0+dudx) + b7*dudy ) &
872  + s23(3,i)*( b9*dudy + b8*dudz ) &
873  + s13(3,i)*( b9*(1.d0+dudx) + b7*dudz ) )
874  r8 = r8 + &
875  ( s11(3,i)*b7*dvdx + s22(3,i)*b8*(1.d0+dvdy) + s33(3,i)*b9*dvdz &
876  + s12(3,i)*( b7*(1.d0+dvdy) + b8*dvdx ) &
877  + s23(3,i)*( b9*(1.d0+dvdy) + b8*dvdz ) &
878  + s13(3,i)*( b9*dvdx + b7*dvdz ) )
879  r9 = r9 + &
880  ( s11(3,i)*b7*dwdx + s22(3,i)*b8*dwdy + s33(3,i)*b9*(1.d0+dwdz) &
881  + s12(3,i)*( b8*dwdx + b7*dwdy ) &
882  + s23(3,i)*( b9*dwdy + b8*(1.d0 + dwdz) ) &
883  + s13(3,i)*( b9*dwdx + b7*(1.d0 + dwdz) ) )
884 
885  r10 = r10 + &
886  ( s11(3,i)*b10*(1.d0+dudx) + s22(3,i)*b11*dudy+s33(3,i)*b12*dudz &
887  + s12(3,i)*( b11*(1.d0+dudx) + b10*dudy ) &
888  + s23(3,i)*( b12*dudy + b11*dudz ) &
889  + s13(3,i)*( b12*(1.d0+dudx) + b10*dudz ) )
890  r11 = r11 + &
891  ( s11(3,i)*b10*dvdx + s22(3,i)*b11*(1.d0+dvdy)+s33(3,i)*b12*dvdz &
892  + s12(3,i)*( b10*(1.d0+dvdy) + b11*dvdx ) &
893  + s23(3,i)*( b12*(1.d0+dvdy) + b11*dvdz ) &
894  + s13(3,i)*( b12*dvdx + b10*dvdz ) )
895  r12 = r12 + &
896  ( s11(3,i)*b10*dwdx + s22(3,i)*b11*dwdy+s33(3,i)*b12*(1.d0+dwdz) &
897  + s12(3,i)*( b11*dwdx + b10*dwdy ) &
898  + s23(3,i)*( b12*dwdy + b11*(1.d0 + dwdz) ) &
899  + s13(3,i)*( b12*dwdx + b10*(1.d0 + dwdz) ) )
900 
901  r13 = r13 + &
902  ( s11(3,i)*b13*(1.d0+dudx) + s22(3,i)*b14*dudy + s33(3,i)*b15*dudz &
903  + s12(3,i)*( b14*(1.d0+dudx) + b13*dudy ) &
904  + s23(3,i)*( b15*dudy + b14*dudz ) &
905  + s13(3,i)*( b15*(1.d0+dudx) + b13*dudz ) )
906  r14 = r14 + &
907  ( s11(3,i)*b13*dvdx + s22(3,i)*b14*(1.d0+dvdy) + s33(3,i)*b15*dvdz &
908  + s12(3,i)*( b13*(1.d0+dvdy) + b14*dvdx ) &
909  + s23(3,i)*( b15*(1.d0+dvdy) + b14*dvdz ) &
910  + s13(3,i)*( b15*dvdx + b13*dvdz ) )
911  r15 = r15 + &
912  ( s11(3,i)*b13*dwdx + s22(3,i)*b14*dwdy + s33(3,i)*b15*(1.d0+dwdz) &
913  + s12(3,i)*( b14*dwdx + b13*dwdy ) &
914  + s23(3,i)*( b15*dwdy + b14*(1.d0 + dwdz) ) &
915  + s13(3,i)*( b15*dwdx + b13*(1.d0 + dwdz) ) )
916 
917  r16 = r16 + &
918  ( s11(3,i)*b16*(1.d0+dudx) + s22(3,i)*b17*dudy + s33(3,i)*b18*dudz &
919  + s12(3,i)*( b17*(1.d0+dudx) + b16*dudy ) &
920  + s23(3,i)*( b18*dudy + b17*dudz ) &
921  + s13(3,i)*( b18*(1.d0+dudx) + b16*dudz ) )
922  r17 = r17 + &
923  ( s11(3,i)*b16*dvdx + s22(3,i)*b17*(1.d0+dvdy) + s33(3,i)*b18*dvdz &
924  + s12(3,i)*( b16*(1.d0+dvdy) + b17*dvdx ) &
925  + s23(3,i)*( b18*(1.d0+dvdy) + b17*dvdz ) &
926  + s13(3,i)*( b18*dvdx + b16*dvdz ) )
927  r18 = r18 + &
928  ( s11(3,i)*b16*dwdx + s22(3,i)*b17*dwdy + s33(3,i)*b18*(1.d0+dwdz) &
929  + s12(3,i)*( b17*dwdx + b16*dwdy ) &
930  + s23(3,i)*( b18*dwdy + b17*(1.d0 + dwdz) ) &
931  + s13(3,i)*( b18*dwdx + b16*(1.d0 + dwdz) ) )
932 
933  r19 = r19 + &
934  ( s11(3,i)*b19*(1.d0+dudx) + s22(3,i)*b20*dudy + s33(3,i)*b21*dudz &
935  + s12(3,i)*( b20*(1.d0+dudx) + b19*dudy ) &
936  + s23(3,i)*( b21*dudy + b20*dudz ) &
937  + s13(3,i)*( b21*(1.d0+dudx) + b19*dudz ) )
938  r20 = r20 + &
939  ( s11(3,i)*b19*dvdx + s22(3,i)*b20*(1.d0+dvdy) + s33(3,i)*b21*dvdz &
940  + s12(3,i)*( b19*(1.d0+dvdy) + b20*dvdx ) &
941  + s23(3,i)*( b21*(1.d0+dvdy) + b20*dvdz ) &
942  + s13(3,i)*( b21*dvdx + b19*dvdz ) )
943  r21 = r21 + &
944  ( s11(3,i)*b19*dwdx + s22(3,i)*b20*dwdy + s33(3,i)*b21*(1.d0+dwdz) &
945  + s12(3,i)*( b20*dwdx + b19*dwdy ) &
946  + s23(3,i)*( b21*dwdy + b20*(1.d0 + dwdz) ) &
947  + s13(3,i)*( b21*dwdx + b19*(1.d0 + dwdz) ) )
948 
949  r22 = r22 + &
950  ( s11(3,i)*b22*(1.d0+dudx) + s22(3,i)*b23*dudy+s33(3,i)*b24*dudz &
951  + s12(3,i)*( b23*(1.d0+dudx) + b22*dudy ) &
952  + s23(3,i)*( b24*dudy + b23*dudz ) &
953  + s13(3,i)*( b24*(1.d0+dudx) + b22*dudz ) )
954  r23 = r23 + &
955  ( s11(3,i)*b22*dvdx + s22(3,i)*b23*(1.d0+dvdy)+s33(3,i)*b24*dvdz &
956  + s12(3,i)*( b22*(1.d0+dvdy) + b23*dvdx ) &
957  + s23(3,i)*( b24*(1.d0+dvdy) + b23*dvdz ) &
958  + s13(3,i)*( b24*dvdx + b22*dvdz ) )
959  r24 = r24 + &
960  ( s11(3,i)*b22*dwdx + s22(3,i)*b23*dwdy+s33(3,i)*b24*(1.d0+dwdz) &
961  + s12(3,i)*( b23*dwdx + b22*dwdy ) &
962  + s23(3,i)*( b24*dwdy + b23*(1.d0 + dwdz) ) &
963  + s13(3,i)*( b24*dwdx + b22*(1.d0 + dwdz) ) )
964 
965  r25 = r25 + &
966  ( s11(3,i)*b25*(1.d0+dudx) + s22(3,i)*b26*dudy+s33(3,i)*b27*dudz &
967  + s12(3,i)*( b26*(1.d0+dudx) + b25*dudy ) &
968  + s23(3,i)*( b27*dudy + b26*dudz ) &
969  + s13(3,i)*( b27*(1.d0+dudx) + b25*dudz ) )
970  r26 = r26 + &
971  ( s11(3,i)*b25*dvdx + s22(3,i)*b26*(1.d0+dvdy)+s33(3,i)*b27*dvdz &
972  + s12(3,i)*( b25*(1.d0+dvdy) + b26*dvdx ) &
973  + s23(3,i)*( b27*(1.d0+dvdy) + b26*dvdz ) &
974  + s13(3,i)*( b27*dvdx + b25*dvdz ) )
975  r27 = r27 + &
976  ( s11(3,i)*b25*dwdx + s22(3,i)*b26*dwdy+s33(3,i)*b27*(1.d0+dwdz) &
977  + s12(3,i)*( b26*dwdx + b25*dwdy ) &
978  + s23(3,i)*( b27*dwdy + b26*(1.d0 + dwdz) ) &
979  + s13(3,i)*( b27*dwdx + b25*(1.d0 + dwdz) ) )
980 
981  r28 = r28 + &
982  ( s11(3,i)*b28*(1.d0+dudx) + s22(3,i)*b29*dudy+s33(3,i)*b30*dudz &
983  + s12(3,i)*( b29*(1.d0+dudx) + b28*dudy ) &
984  + s23(3,i)*( b30*dudy + b29*dudz ) &
985  + s13(3,i)*( b30*(1.d0+dudx) + b28*dudz ) )
986  r29 = r29 + &
987  ( s11(3,i)*b28*dvdx + s22(3,i)*b29*(1.d0+dvdy)+s33(3,i)*b30*dvdz &
988  + s12(3,i)*( b28*(1.d0+dvdy) + b29*dvdx ) &
989  + s23(3,i)*( b30*(1.d0+dvdy) + b29*dvdz ) &
990  + s13(3,i)*( b30*dvdx + b28*dvdz ) )
991  r30 = r30 + &
992  ( s11(3,i)*b28*dwdx + s22(3,i)*b29*dwdy+s33(3,i)*b30*(1.d0+dwdz) &
993  + s12(3,i)*( b29*dwdx + b28*dwdy ) &
994  + s23(3,i)*( b30*dwdy + b29*(1.d0 + dwdz) ) &
995  + s13(3,i)*( b30*dwdx + b28*(1.d0 + dwdz) ) )
996 
997 !-- Gauss Point 4
998 
999  g1 = 0.13819660d0
1000  g2 = 0.13819660d0
1001  g3 = 0.13819660d0
1002  g4 = 0.58541020d0
1003  xn1 = (4.d0*g1-1.d0) ! derivative of shape function
1004  xn2 = (4.d0*g2-1.d0) ! dN_i/dzeta_i
1005  xn3 = (4.d0*g3-1.d0)
1006  xn4 = (4.d0*g4-1.d0)
1007 ! xN5 = 4.d0*g1*g2
1008 ! xN6 = 4.d0*g2*g3
1009 ! xN7 = 4.d0*g3*g1
1010 ! xN8 = 4.d0*g1*g4
1011 ! xN9 = 4.d0*g2*g4
1012 ! xN10= 4.d0*g3*g4
1013 
1014  b1 = aux1*xn1
1015  b2 = aux2*xn1
1016  b3 = aux3*xn1
1017  b4 = aux4*xn2
1018  b5 = aux5*xn2
1019  b6 = aux6*xn2
1020  b7 = aux7*xn3
1021  b8 = aux8*xn3
1022  b9 = aux9*xn3
1023  b10 = aux10*xn4
1024  b11 = aux11*xn4
1025  b12 = aux12*xn4
1026 
1027  b13 = 4.d0*(g2*aux1 + g1*aux4)
1028  b14 = 4.d0*(g2*aux2 + g1*aux5)
1029  b15 = 4.d0*(g2*aux3 + g1*aux6)
1030 
1031  b16 = 4.d0*(g3*aux4 + g2*aux7)
1032  b17 = 4.d0*(g3*aux5 + g2*aux8)
1033  b18 = 4.d0*(g3*aux6 + g2*aux9)
1034 
1035  b19 = 4.d0*(g1*aux7 + g3*aux1)
1036  b20 = 4.d0*(g1*aux8 + g3*aux2)
1037  b21 = 4.d0*(g1*aux9 + g3*aux3)
1038 
1039  b22 = 4.d0*(g4*aux1 + g1*aux10)
1040  b23 = 4.d0*(g4*aux2 + g1*aux11)
1041  b24 = 4.d0*(g4*aux3 + g1*aux12)
1042 
1043  b25 = 4.d0*(g4*aux4 + g2*aux10)
1044  b26 = 4.d0*(g4*aux5 + g2*aux11)
1045  b27 = 4.d0*(g4*aux6 + g2*aux12)
1046 
1047  b28 = 4.d0*(g4*aux7 + g3*aux10)
1048  b29 = 4.d0*(g4*aux8 + g3*aux11)
1049  b30 = 4.d0*(g4*aux9 + g3*aux12)
1050 !-----Calculate displacement gradient (H)
1051  dudx = (b1*u1 + b4*u2 + b7*u3 + b10*u4 + b13*u5 + b16*u6 + b19*u7 + b22*u8 + b25*u9 + b28*u10)*vx6inv
1052  dvdy = (b2*v1 + b5*v2 + b8*v3 + b11*v4 + b14*v5 + b17*v6 + b20*v7 + b23*v8 + b26*v9 + b29*v10)*vx6inv
1053  dwdz = (b3*w1 + b6*w2 + b9*w3 + b12*w4 + b15*w5 + b18*w6 + b21*w7 + b24*w8 + b27*w9 + b30*w10)*vx6inv
1054  dudy = (b2*u1 + b5*u2 + b8*u3 + b11*u4 + b14*u5 + b17*u6 + b20*u7 + b23*u8 + b26*u9 + b29*u10)*vx6inv
1055  dvdx = (b1*v1 + b4*v2 + b7*v3 + b10*v4 + b13*v5 + b16*v6 + b19*v7 + b22*v8 + b25*v9 + b28*v10)*vx6inv
1056  dvdz = (b3*v1 + b6*v2 + b9*v3 + b12*v4 + b15*v5 + b18*v6 + b21*v7 + b24*v8 + b27*v9 + b30*v10)*vx6inv
1057  dwdy = (b2*w1 + b5*w2 + b8*w3 + b11*w4 + b14*w5 + b17*w6 + b20*w7 + b23*w8 + b26*w9 + b29*w10)*vx6inv
1058  dudz = (b3*u1 + b6*u2 + b9*u3 + b12*u4 + b15*u5 + b18*u6 + b21*u7 + b24*u8 + b27*u9 + b30*u10)*vx6inv
1059  dwdx = (b1*w1 + b4*w2 + b7*w3 + b10*w4 + b13*w5 + b16*w6 + b19*w7 + b22*w8 + b25*w9 + b28*w10)*vx6inv
1060 !
1061 ! Green-Lagrange Strain (E)
1062 !
1063  e11 = dudx + 0.5d0*(dudx*dudx + dvdx*dvdx + dwdx*dwdx)
1064  e22 = dvdy + 0.5d0*(dudy*dudy + dvdy*dvdy + dwdy*dwdy)
1065  e33 = dwdz + 0.5d0*(dudz*dudz + dvdz*dvdz + dwdz*dwdz)
1066  e12 = dudy + dvdx + ( dudx*dudy + dvdx*dvdy + dwdx*dwdy )
1067  e23 = dvdz + dwdy + ( dudy*dudz + dvdy*dvdz + dwdy*dwdz )
1068  e13 = dwdx + dudz + ( dudz*dudx + dvdz*dvdx + dwdz*dwdx )
1069 !
1070 ! Second Piola-Kirchhoff stress (S)
1071 !
1072  s11(4,i) = e11*ci(1,j) + e22*ci(2,j) + e33*ci(4,j)
1073  s22(4,i) = e11*ci(2,j) + e22*ci(3,j) + e33*ci(5,j)
1074  s33(4,i) = e11*ci(4,j) + e22*ci(5,j) + e33*ci(6,j)
1075  s12(4,i) = e12*ci(7,j)
1076  s23(4,i) = e23*ci(8,j)
1077  s13(4,i) = e13*ci(9,j)
1078 
1079  r1 = r1 + &
1080  ( s11(4,i)*b1*(1.d0+dudx) + s22(4,i)*b2*dudy + s33(4,i)*b3*dudz &
1081  + s12(4,i)*( b2*(1.d0+dudx) + b1*dudy ) &
1082  + s23(4,i)*( b3*dudy + b2*dudz ) &
1083  + s13(4,i)*( b3*(1.d0+dudx) + b1*dudz ) )
1084  r2 = r2 + &
1085  ( s11(4,i)*b1*dvdx + s22(4,i)*b2*(1.d0+dvdy) + s33(4,i)*b3*dvdz &
1086  + s12(4,i)*( b1*(1.d0+dvdy) + b2*dvdx ) &
1087  + s23(4,i)*( b3*(1.d0+dvdy) + b2*dvdz ) &
1088  + s13(4,i)*( b3*dvdx + b1*dvdz ) )
1089  r3 = r3 + &
1090  ( s11(4,i)*b1*dwdx + s22(4,i)*b2*dwdy + s33(4,i)*b3*(1.d0+dwdz) &
1091  + s12(4,i)*( b2*dwdx + b1*dwdy ) &
1092  + s23(4,i)*( b3*dwdy + b2*(1.d0 + dwdz) ) &
1093  + s13(4,i)*( b3*dwdx + b1*(1.d0 + dwdz) ) )
1094 
1095  r4 = r4 + &
1096  ( s11(4,i)*b4*(1.d0+dudx) + s22(4,i)*b5*dudy + s33(4,i)*b6*dudz &
1097  + s12(4,i)*( b5*(1.d0+dudx) + b4*dudy ) &
1098  + s23(4,i)*( b6*dudy + b5*dudz ) &
1099  + s13(4,i)*( b6*(1.d0+dudx) + b4*dudz ) )
1100  r5 = r5 + &
1101  ( s11(4,i)*b4*dvdx + s22(4,i)*b5*(1.d0+dvdy) + s33(4,i)*b6*dvdz &
1102  + s12(4,i)*( b4*(1.d0+dvdy) + b5*dvdx ) &
1103  + s23(4,i)*( b6*(1.d0+dvdy) + b5*dvdz ) &
1104  + s13(4,i)*( b6*dvdx + b4*dvdz ) )
1105  r6 = r6 + &
1106  ( s11(4,i)*b4*dwdx + s22(4,i)*b5*dwdy + s33(4,i)*b6*(1.d0+dwdz) &
1107  + s12(4,i)*( b5*dwdx + b4*dwdy ) &
1108  + s23(4,i)*( b6*dwdy + b5*(1.d0 + dwdz) ) &
1109  + s13(4,i)*( b6*dwdx + b4*(1.d0 + dwdz) ) )
1110 
1111  r7 = r7 + &
1112  ( s11(4,i)*b7*(1.d0+dudx) + s22(4,i)*b8*dudy + s33(4,i)*b9*dudz &
1113  + s12(4,i)*( b8*(1.d0+dudx) + b7*dudy ) &
1114  + s23(4,i)*( b9*dudy + b8*dudz ) &
1115  + s13(4,i)*( b9*(1.d0+dudx) + b7*dudz ) )
1116  r8 = r8 + &
1117  ( s11(4,i)*b7*dvdx + s22(4,i)*b8*(1.d0+dvdy) + s33(4,i)*b9*dvdz &
1118  + s12(4,i)*( b7*(1.d0+dvdy) + b8*dvdx ) &
1119  + s23(4,i)*( b9*(1.d0+dvdy) + b8*dvdz ) &
1120  + s13(4,i)*( b9*dvdx + b7*dvdz ) )
1121  r9 = r9 + &
1122  ( s11(4,i)*b7*dwdx + s22(4,i)*b8*dwdy + s33(4,i)*b9*(1.d0+dwdz) &
1123  + s12(4,i)*( b8*dwdx + b7*dwdy ) &
1124  + s23(4,i)*( b9*dwdy + b8*(1.d0 + dwdz) ) &
1125  + s13(4,i)*( b9*dwdx + b7*(1.d0 + dwdz) ) )
1126 
1127  r10 = r10 + &
1128  ( s11(4,i)*b10*(1.d0+dudx) + s22(4,i)*b11*dudy+s33(4,i)*b12*dudz &
1129  + s12(4,i)*( b11*(1.d0+dudx) + b10*dudy ) &
1130  + s23(4,i)*( b12*dudy + b11*dudz ) &
1131  + s13(4,i)*( b12*(1.d0+dudx) + b10*dudz ) )
1132  r11 = r11 + &
1133  ( s11(4,i)*b10*dvdx + s22(4,i)*b11*(1.d0+dvdy)+s33(4,i)*b12*dvdz &
1134  + s12(4,i)*( b10*(1.d0+dvdy) + b11*dvdx ) &
1135  + s23(4,i)*( b12*(1.d0+dvdy) + b11*dvdz ) &
1136  + s13(4,i)*( b12*dvdx + b10*dvdz ) )
1137  r12 = r12 + &
1138  ( s11(4,i)*b10*dwdx + s22(4,i)*b11*dwdy+s33(4,i)*b12*(1.d0+dwdz) &
1139  + s12(4,i)*( b11*dwdx + b10*dwdy ) &
1140  + s23(4,i)*( b12*dwdy + b11*(1.d0 + dwdz) ) &
1141  + s13(4,i)*( b12*dwdx + b10*(1.d0 + dwdz) ) )
1142 
1143  r13 = r13 + &
1144  ( s11(4,i)*b13*(1.d0+dudx) + s22(4,i)*b14*dudy + s33(4,i)*b15*dudz &
1145  + s12(4,i)*( b14*(1.d0+dudx) + b13*dudy ) &
1146  + s23(4,i)*( b15*dudy + b14*dudz ) &
1147  + s13(4,i)*( b15*(1.d0+dudx) + b13*dudz ) )
1148  r14 = r14 + &
1149  ( s11(4,i)*b13*dvdx + s22(4,i)*b14*(1.d0+dvdy) + s33(4,i)*b15*dvdz &
1150  + s12(4,i)*( b13*(1.d0+dvdy) + b14*dvdx ) &
1151  + s23(4,i)*( b15*(1.d0+dvdy) + b14*dvdz ) &
1152  + s13(4,i)*( b15*dvdx + b13*dvdz ) )
1153  r15 = r15 + &
1154  ( s11(4,i)*b13*dwdx + s22(4,i)*b14*dwdy + s33(4,i)*b15*(1.d0+dwdz) &
1155  + s12(4,i)*( b14*dwdx + b13*dwdy ) &
1156  + s23(4,i)*( b15*dwdy + b14*(1.d0 + dwdz) ) &
1157  + s13(4,i)*( b15*dwdx + b13*(1.d0 + dwdz) ) )
1158 
1159  r16 = r16 + &
1160  ( s11(4,i)*b16*(1.d0+dudx) + s22(4,i)*b17*dudy + s33(4,i)*b18*dudz &
1161  + s12(4,i)*( b17*(1.d0+dudx) + b16*dudy ) &
1162  + s23(4,i)*( b18*dudy + b17*dudz ) &
1163  + s13(4,i)*( b18*(1.d0+dudx) + b16*dudz ) )
1164  r17 = r17 + &
1165  ( s11(4,i)*b16*dvdx + s22(4,i)*b17*(1.d0+dvdy) + s33(4,i)*b18*dvdz &
1166  + s12(4,i)*( b16*(1.d0+dvdy) + b17*dvdx ) &
1167  + s23(4,i)*( b18*(1.d0+dvdy) + b17*dvdz ) &
1168  + s13(4,i)*( b18*dvdx + b16*dvdz ) )
1169  r18 = r18 + &
1170  ( s11(4,i)*b16*dwdx + s22(4,i)*b17*dwdy + s33(4,i)*b18*(1.d0+dwdz) &
1171  + s12(4,i)*( b17*dwdx + b16*dwdy ) &
1172  + s23(4,i)*( b18*dwdy + b17*(1.d0 + dwdz) ) &
1173  + s13(4,i)*( b18*dwdx + b16*(1.d0 + dwdz) ) )
1174 
1175  r19 = r19 + &
1176  ( s11(4,i)*b19*(1.d0+dudx) + s22(4,i)*b20*dudy + s33(4,i)*b21*dudz &
1177  + s12(4,i)*( b20*(1.d0+dudx) + b19*dudy ) &
1178  + s23(4,i)*( b21*dudy + b20*dudz ) &
1179  + s13(4,i)*( b21*(1.d0+dudx) + b19*dudz ) )
1180  r20 = r20 + &
1181  ( s11(4,i)*b19*dvdx + s22(4,i)*b20*(1.d0+dvdy) + s33(4,i)*b21*dvdz &
1182  + s12(4,i)*( b19*(1.d0+dvdy) + b20*dvdx ) &
1183  + s23(4,i)*( b21*(1.d0+dvdy) + b20*dvdz ) &
1184  + s13(4,i)*( b21*dvdx + b19*dvdz ) )
1185  r21 = r21 + &
1186  ( s11(4,i)*b19*dwdx + s22(4,i)*b20*dwdy + s33(4,i)*b21*(1.d0+dwdz) &
1187  + s12(4,i)*( b20*dwdx + b19*dwdy ) &
1188  + s23(4,i)*( b21*dwdy + b20*(1.d0 + dwdz) ) &
1189  + s13(4,i)*( b21*dwdx + b19*(1.d0 + dwdz) ) )
1190 
1191  r22 = r22 + &
1192  ( s11(4,i)*b22*(1.d0+dudx) + s22(4,i)*b23*dudy+s33(4,i)*b24*dudz &
1193  + s12(4,i)*( b23*(1.d0+dudx) + b22*dudy ) &
1194  + s23(4,i)*( b24*dudy + b23*dudz ) &
1195  + s13(4,i)*( b24*(1.d0+dudx) + b22*dudz ) )
1196  r23 = r23 + &
1197  ( s11(4,i)*b22*dvdx + s22(4,i)*b23*(1.d0+dvdy)+s33(4,i)*b24*dvdz &
1198  + s12(4,i)*( b22*(1.d0+dvdy) + b23*dvdx ) &
1199  + s23(4,i)*( b24*(1.d0+dvdy) + b23*dvdz ) &
1200  + s13(4,i)*( b24*dvdx + b22*dvdz ) )
1201  r24 = r24 + &
1202  ( s11(4,i)*b22*dwdx + s22(4,i)*b23*dwdy+s33(4,i)*b24*(1.d0+dwdz) &
1203  + s12(4,i)*( b23*dwdx + b22*dwdy ) &
1204  + s23(4,i)*( b24*dwdy + b23*(1.d0 + dwdz) ) &
1205  + s13(4,i)*( b24*dwdx + b22*(1.d0 + dwdz) ) )
1206 
1207  r25 = r25 + &
1208  ( s11(4,i)*b25*(1.d0+dudx) + s22(4,i)*b26*dudy+s33(4,i)*b27*dudz &
1209  + s12(4,i)*( b26*(1.d0+dudx) + b25*dudy ) &
1210  + s23(4,i)*( b27*dudy + b26*dudz ) &
1211  + s13(4,i)*( b27*(1.d0+dudx) + b25*dudz ) )
1212  r26 = r26 + &
1213  ( s11(4,i)*b25*dvdx + s22(4,i)*b26*(1.d0+dvdy)+s33(4,i)*b27*dvdz &
1214  + s12(4,i)*( b25*(1.d0+dvdy) + b26*dvdx ) &
1215  + s23(4,i)*( b27*(1.d0+dvdy) + b26*dvdz ) &
1216  + s13(4,i)*( b27*dvdx + b25*dvdz ) )
1217  r27 = r27 + &
1218  ( s11(4,i)*b25*dwdx + s22(4,i)*b26*dwdy+s33(4,i)*b27*(1.d0+dwdz) &
1219  + s12(4,i)*( b26*dwdx + b25*dwdy ) &
1220  + s23(4,i)*( b27*dwdy + b26*(1.d0 + dwdz) ) &
1221  + s13(4,i)*( b27*dwdx + b25*(1.d0 + dwdz) ) )
1222 
1223  r28 = r28 + &
1224  ( s11(4,i)*b28*(1.d0+dudx) + s22(4,i)*b29*dudy+s33(4,i)*b30*dudz &
1225  + s12(4,i)*( b29*(1.d0+dudx) + b28*dudy ) &
1226  + s23(4,i)*( b30*dudy + b29*dudz ) &
1227  + s13(4,i)*( b30*(1.d0+dudx) + b28*dudz ) )
1228  r29 = r29 + &
1229  ( s11(4,i)*b28*dvdx + s22(4,i)*b29*(1.d0+dvdy)+s33(4,i)*b30*dvdz &
1230  + s12(4,i)*( b28*(1.d0+dvdy) + b29*dvdx ) &
1231  + s23(4,i)*( b30*(1.d0+dvdy) + b29*dvdz ) &
1232  + s13(4,i)*( b30*dvdx + b28*dvdz ) )
1233  r30 = r30 + &
1234  ( s11(4,i)*b28*dwdx + s22(4,i)*b29*dwdy+s33(4,i)*b30*(1.d0+dwdz) &
1235  + s12(4,i)*( b29*dwdx + b28*dwdy ) &
1236  + s23(4,i)*( b30*dwdy + b29*(1.d0 + dwdz) ) &
1237  + s13(4,i)*( b30*dwdx + b28*(1.d0 + dwdz) ) )
1238 
1239 ! Wi (i.e. weight) for 4 guass point integration is 1/4
1240 
1241 ! Wi * 1/6 because the volume of a reference tetrahedra in
1242 ! volume coordinates is 1/6
1243 
1244 ! ASSEMBLE THE INTERNAL FORCE VECTOR
1245 !
1246 ! local node 1
1247  r_in(k1n1) = r_in(k1n1) - r1*0.04166666666666667d0
1248  r_in(k2n1) = r_in(k2n1) - r2*0.04166666666666667d0
1249  r_in(k3n1) = r_in(k3n1) - r3*0.04166666666666667d0
1250 ! local node 2
1251  r_in(k1n2) = r_in(k1n2) - r4*0.04166666666666667d0
1252  r_in(k2n2) = r_in(k2n2) - r5*0.04166666666666667d0
1253  r_in(k3n2) = r_in(k3n2) - r6*0.04166666666666667d0
1254 ! local node 3
1255  r_in(k1n3) = r_in(k1n3) - r7*0.04166666666666667d0
1256  r_in(k2n3) = r_in(k2n3) - r8*0.04166666666666667d0
1257  r_in(k3n3) = r_in(k3n3) - r9*0.04166666666666667d0
1258 ! local node 4
1259  r_in(k1n4) = r_in(k1n4) - r10*0.04166666666666667d0
1260  r_in(k2n4) = r_in(k2n4) - r11*0.04166666666666667d0
1261  r_in(k3n4) = r_in(k3n4) - r12*0.04166666666666667d0
1262 ! local node 5
1263  r_in(k1n5) = r_in(k1n5) - r13*0.04166666666666667d0
1264  r_in(k2n5) = r_in(k2n5) - r14*0.04166666666666667d0
1265  r_in(k3n5) = r_in(k3n5) - r15*0.04166666666666667d0
1266 ! local node 6
1267  r_in(k1n6) = r_in(k1n6) - r16*0.04166666666666667d0
1268  r_in(k2n6) = r_in(k2n6) - r17*0.04166666666666667d0
1269  r_in(k3n6) = r_in(k3n6) - r18*0.04166666666666667d0
1270 ! local node 7
1271  r_in(k1n7) = r_in(k1n7) - r19*0.04166666666666667d0
1272  r_in(k2n7) = r_in(k2n7) - r20*0.04166666666666667d0
1273  r_in(k3n7) = r_in(k3n7) - r21*0.04166666666666667d0
1274 ! local node 8
1275  r_in(k1n8) = r_in(k1n8) - r22*0.04166666666666667d0
1276  r_in(k2n8) = r_in(k2n8) - r23*0.04166666666666667d0
1277  r_in(k3n8) = r_in(k3n8) - r24*0.04166666666666667d0
1278 ! local node 9
1279  r_in(k1n9) = r_in(k1n9) - r25*0.04166666666666667d0
1280  r_in(k2n9) = r_in(k2n9) - r26*0.04166666666666667d0
1281  r_in(k3n9) = r_in(k3n9) - r27*0.04166666666666667d0
1282 ! local node 10
1283  r_in(k1n10) = r_in(k1n10) - r28*0.04166666666666667d0
1284  r_in(k2n10) = r_in(k2n10) - r29*0.04166666666666667d0
1285  r_in(k3n10) = r_in(k3n10) - r30*0.04166666666666667d0
1286 
1287  ENDDO
1288  RETURN
1289 END SUBROUTINE v3d10_nl
1290 
subroutine v3d10_nl(coor, matcstet, lmcstet, R_in, d, ci, S11, S22, S33, S12, S23, S13, numnp, nstart, nend, numcstet, numat_vol)
Definition: v3d10_nl.f90:53
const NT & d
blockLoc i
Definition: read.cpp:79
j indices j
Definition: Indexing.h:6