Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MP/Source/v3d10.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(coor,matcstet,lmcstet,R_in,d,ci, &
54  s11l,s22l,s33l,s12l,s23l,s13l,&
55  numnp,nstart,nend,numcstet,numat_vol)
56 !-----Performs displacement based computations for 3D, 10-node
57 !-----Quadratic Linear Elastic Tetrahedron with quadratic interpolation
58 !-----functions. (linear strain tetrahedra)
59 
60 !----- 4 Guass points
61  IMPLICIT NONE
62 !-----Global variables
63  INTEGER :: numnp ! number of nodes
64  INTEGER :: numat_vol ! number of volumetric materials
65  INTEGER :: numcstet ! number of LSTets
66 !-- coordinate array
67  REAL*8, DIMENSION(1:3,1:numnp) :: coor
68 !-- elastic stiffness consts
69  REAL*8, DIMENSION(1:9,1:numat_vol) :: ci
70 !-- internal force
71  REAL*8, DIMENSION(1:3*numnp) :: r_in
72 !-- displacement vector
73  REAL*8, DIMENSION(1:3*numnp) :: d
74 !-- CSTet stress
75  REAL*8, DIMENSION(1:4,1:numcstet) :: s11l, s22l, s33l, s12l, s23l, s13l
76 !-- connectivity table for CSTet
77  INTEGER, DIMENSION(1:10,1:numcstet) :: lmcstet
78 !-- mat number for CSTet element
79  INTEGER, DIMENSION(1:numcstet) :: matcstet
80 !---- Local variables
81 !-- node numbers
82  INTEGER :: n1,n2,n3,n4,n5,n6,n7,n8,n9,n10
83 !-- x, y, and z displacements of nodes
84  REAL*8 :: u1,u2,u3,u4,u5,u6,u7,u8,u9,u10
85  REAL*8 :: v1,v2,v3,v4,v5,v6,v7,v8,v9,v10
86  REAL*8 :: w1,w2,w3,w4,w5,w6,w7,w8,w9,w10
87 !-- 6*volume and inverse of 6*volume
88  REAL*8 :: vx6, vx6inv
89 !-- spacial derivatives
90  REAL*8 :: b1,b2,b3,b4,b5,b6,b7,b8,b9,b10
91  REAL*8 :: b11,b12,b13,b14,b15,b16,b17,b18,b19,b20
92  REAL*8 :: b21,b22,b23,b24,b25,b26,b27,b28,b29,b30
93 !-- strains
94  REAL*8 :: e11,e22,e33,e12,e23,e13
95 !-- coordinate holding variable
96  REAL*8 :: x1,x2,x3,x4,x5,x6,x7,x8,x9,x10
97  REAL*8 :: y1,y2,y3,y4,y5,y6,y7,y8,y9,y10
98  REAL*8 :: z1,z2,z3,z4,z5,z6,z7,z8,z9,z10
99 !-- dummy and counters
100  INTEGER :: i,j,nstart,nend
101  REAL*8 :: aux1,aux2,aux3,aux4,aux5,aux6,aux7,aux8,aux9,aux10,aux11,aux12
102 !-- partial internal force
103  REAL*8 :: r1,r2,r3,r4,r5,r6,r7,r8,r9,r10,r11,r12,r13,r14,r15,r16,r17,r18
104  REAL*8 :: r19,r20,r21,r22,r23,r24,r25,r26,r27,r28,r29,r30
105  REAL*8 :: g1, g2, g3, g4
106  REAL*8 :: xn1, xn2, xn3, xn4
107 !-- Coordinate subtractions
108  REAL*8 :: x14, x24, x34, y14, y24, y34, z14, z24, z34
109 !--
110  REAL*8 :: c11, c21, c31
111  INTEGER :: k1n1,k1n2,k1n3,k1n4,k1n5,k1n6,k1n7,k1n8,k1n9,k1n10
112  INTEGER :: k2n1,k2n2,k2n3,k2n4,k2n5,k2n6,k2n7,k2n8,k2n9,k2n10
113  INTEGER :: k3n1,k3n2,k3n3,k3n4,k3n5,k3n6,k3n7,k3n8,k3n9,k3n10
114 
115 
116  DO i = nstart, nend
117  j = matcstet(i)
118 
119  n1 = lmcstet(1,i)
120  n2 = lmcstet(2,i)
121  n3 = lmcstet(3,i)
122  n4 = lmcstet(4,i)
123  n5 = lmcstet(5,i)
124  n6 = lmcstet(6,i)
125  n7 = lmcstet(7,i)
126  n8 = lmcstet(8,i)
127  n9 = lmcstet(9,i)
128  n10 = lmcstet(10,i)
129 
130  k3n1 = 3*n1
131  k3n2 = 3*n2
132  k3n3 = 3*n3
133  k3n4 = 3*n4
134  k3n5 = 3*n5
135  k3n6 = 3*n6
136  k3n7 = 3*n7
137  k3n8 = 3*n8
138  k3n9 = 3*n9
139  k3n10 = 3*n10
140 
141  k2n1 = k3n1 - 1
142  k2n2 = k3n2 - 1
143  k2n3 = k3n3 - 1
144  k2n4 = k3n4 - 1
145  k2n5 = k3n5 - 1
146  k2n6 = k3n6 - 1
147  k2n7 = k3n7 - 1
148  k2n8 = k3n8 - 1
149  k2n9 = k3n9 - 1
150  k2n10 = k3n10 - 1
151 
152  k1n1 = k3n1 - 2
153  k1n2 = k3n2 - 2
154  k1n3 = k3n3 - 2
155  k1n4 = k3n4 - 2
156  k1n5 = k3n5 - 2
157  k1n6 = k3n6 - 2
158  k1n7 = k3n7 - 2
159  k1n8 = k3n8 - 2
160  k1n9 = k3n9 - 2
161  k1n10 = k3n10 - 2
162  ! k#n# dummy variables replaces:
163  u1 = d(k1n1) ! (3*n1 -2)
164  u2 = d(k1n2) ! (3*n2 -2)
165  u3 = d(k1n3) ! (3*n3 -2)
166  u4 = d(k1n4) ! (3*n4 -2)
167  u5 = d(k1n5) ! (3*n5 -2)
168  u6 = d(k1n6) ! (3*n6 -2)
169  u7 = d(k1n7) ! (3*n7 -2)
170  u8 = d(k1n8) ! (3*n8 -2)
171  u9 = d(k1n9) ! (3*n9 -2)
172  u10 = d(k1n10) ! (3*n10-2)
173  v1 = d(k2n1) ! (3*n1 -1)
174  v2 = d(k2n2) ! (3*n2 -1)
175  v3 = d(k2n3) ! (3*n3 -1)
176  v4 = d(k2n4) ! (3*n4 -1)
177  v5 = d(k2n5) ! (3*n5 -1)
178  v6 = d(k2n6) ! (3*n6 -1)
179  v7 = d(k2n7) ! (3*n7 -1)
180  v8 = d(k2n8) ! (3*n8 -1)
181  v9 = d(k2n9) ! (3*n9 -1)
182  v10 = d(k2n10) ! (3*n10-1)
183  w1 = d(k3n1) ! (3*n1)
184  w2 = d(k3n2) ! (3*n2)
185  w3 = d(k3n3) ! (3*n3)
186  w4 = d(k3n4) ! (3*n4)
187  w5 = d(k3n5) ! (3*n5)
188  w6 = d(k3n6) ! (3*n6)
189  w7 = d(k3n7) ! (3*n7)
190  w8 = d(k3n8) ! (3*n8)
191  w9 = d(k3n9) ! (3*n9)
192  w10 = d(k3n10) ! (3*n10)
193 
194  x1 = coor(1,n1)
195  x2 = coor(1,n2)
196  x3 = coor(1,n3)
197  x4 = coor(1,n4)
198  y1 = coor(2,n1)
199  y2 = coor(2,n2)
200  y3 = coor(2,n3)
201  y4 = coor(2,n4)
202  z1 = coor(3,n1)
203  z2 = coor(3,n2)
204  z3 = coor(3,n3)
205  z4 = coor(3,n4)
206 
207  x14 = x1 - x4
208  x24 = x2 - x4
209  x34 = x3 - x4
210  y14 = y1 - y4
211  y24 = y2 - y4
212  y34 = y3 - y4
213  z14 = z1 - z4
214  z24 = z2 - z4
215  z34 = z3 - z4
216 
217  c11 = y24*z34 - z24*y34
218  c21 = -( x24*z34 - z24*x34 )
219  c31 = x24*y34 - y24*x34
220 
221  vx6 = -( x14*c11 + y14*c21 + z14*c31 )
222 
223 ! Vx6 = x2*y3*z4 - x2*y4*z3 - y2*x3*z4 + y2*x4*z3 + z2*x3*y4
224 ! $ - z2*x4*y3 - x1*y3*z4 + x1*y4*z3 + x1*y2*z4 - x1*y2*z3
225 ! $ - x1*z2*y4 + x1*z2*y3 + y1*x3*z4 - y1*x4*z3 - y1*x2*z4
226 ! $ + y1*x2*z3 + y1*z2*x4 - y1*z2*x3 - z1*x3*y4 + z1*x4*y3
227 ! $ + z1*x2*y4 - z1*x2*y3 - z1*y2*x4 + z1*y2*x3
228 
229  vx6inv = 1.d0 / vx6
230 
231  aux1 = -(y3*z4 - y4*z3 - y2*z4 + y2*z3 + z2*y4 - z2*y3)
232  aux2 = (x3*z4 - x4*z3 - x2*z4 + x2*z3 + z2*x4 - z2*x3)
233  aux3 = -(x3*y4 - x4*y3 - x2*y4 + x2*y3 + y2*x4 - y2*x3)
234  aux4 = (y3*z4 - y4*z3 - y1*z4 + y1*z3 + z1*y4 - z1*y3)
235  aux5 = -(x3*z4 - x4*z3 - x1*z4 + x1*z3 + z1*x4 - z1*x3)
236  aux6 = (x3*y4 - x4*y3 - x1*y4 + x1*y3 + y1*x4 - y1*x3)
237  aux7 = -(y2*z4 - z2*y4 - y1*z4 + y1*z2 + z1*y4 - z1*y2)
238  aux8 = (x2*z4 - z2*x4 - x1*z4 + x1*z2 + z1*x4 - z1*x2)
239  aux9 = -(x2*y4 - y2*x4 - x1*y4 + x1*y2 + y1*x4 - y1*x2)
240  aux10 = (y2*z3 - z2*y3 - y1*z3 + y1*z2 + z1*y3 - z1*y2)
241  aux11 =-(x2*z3 - z2*x3 - x1*z3 + x1*z2 + z1*x3 - z1*x2)
242  aux12 = (x2*y3 - y2*x3 - x1*y3 + x1*y2 + y1*x3 - y1*x2)
243 
244 !-- Gauss Point 1
245 
246  g1 = 0.58541020d0
247  g2 = 0.13819660d0
248  g3 = 0.13819660d0
249  g4 = 0.13819660d0
250 
251  xn1 = (4.d0*g1-1.d0) ! derivative of shape function
252  xn2 = (4.d0*g2-1.d0) ! dN_i/dzeta_i
253  xn3 = (4.d0*g3-1.d0)
254  xn4 = (4.d0*g4-1.d0)
255 ! xN5 = 4.d0*g1*g2
256 ! xN6 = 4.d0*g2*g3
257 ! xN7 = 4.d0*g3*g1
258 ! xN8 = 4.d0*g1*g4
259 ! xN9 = 4.d0*g2*g4
260 ! xN10= 4.d0*g3*g4
261 
262  b1 = aux1*xn1
263  b2 = aux2*xn1
264  b3 = aux3*xn1
265  b4 = aux4*xn2
266  b5 = aux5*xn2
267  b6 = aux6*xn2
268  b7 = aux7*xn3
269  b8 = aux8*xn3
270  b9 = aux9*xn3
271  b10 = aux10*xn4
272  b11 = aux11*xn4
273  b12 = aux12*xn4
274 
275  b13 = 4.d0*(g2*aux1 + g1*aux4)
276  b14 = 4.d0*(g2*aux2 + g1*aux5)
277  b15 = 4.d0*(g2*aux3 + g1*aux6)
278 
279  b16 = 4.d0*(g3*aux4 + g2*aux7)
280  b17 = 4.d0*(g3*aux5 + g2*aux8)
281  b18 = 4.d0*(g3*aux6 + g2*aux9)
282 
283  b19 = 4.d0*(g1*aux7 + g3*aux1)
284  b20 = 4.d0*(g1*aux8 + g3*aux2)
285  b21 = 4.d0*(g1*aux9 + g3*aux3)
286 
287  b22 = 4.d0*(g4*aux1 + g1*aux10)
288  b23 = 4.d0*(g4*aux2 + g1*aux11)
289  b24 = 4.d0*(g4*aux3 + g1*aux12)
290 
291  b25 = 4.d0*(g4*aux4 + g2*aux10)
292  b26 = 4.d0*(g4*aux5 + g2*aux11)
293  b27 = 4.d0*(g4*aux6 + g2*aux12)
294 
295  b28 = 4.d0*(g4*aux7 + g3*aux10)
296  b29 = 4.d0*(g4*aux8 + g3*aux11)
297  b30 = 4.d0*(g4*aux9 + g3*aux12)
298 
299  e11 = (b1*u1 + b4*u2 + b7*u3 + b10*u4 + b13*u5 + &
300  b16*u6 + b19*u7 + b22*u8 + b25*u9 + b28*u10) * vx6inv
301  e22 = (b2*v1 + b5*v2 + b8*v3 + b11*v4 + b14*v5 + &
302  b17*v6 + b20*v7 + b23*v8 + b26*v9 + b29*v10) * vx6inv
303  e33 = (b3*w1 + b6*w2 + b9*w3 + b12*w4 + b15*w5 + &
304  b18*w6 + b21*w7 + b24*w8 + b27*w9 + b30*w10) * vx6inv
305  e12 = (b2*u1 + b1*v1 + b5*u2 + b4*v2 &
306  + b8*u3 + b7*v3 + b11*u4 + b10*v4 &
307  + b14*u5 + b13*v5 + b17*u6 + b16*v6 &
308  + b20*u7 + b19*v7 + b23*u8 + b22*v8 &
309  + b26*u9 + b25*v9 + b29*u10 + b28*v10) * vx6inv
310  e23 = (b3*v1 + b2*w1 + b6*v2 + b5*w2 &
311  + b9*v3 + b8*w3 + b12*v4 + b11*w4 &
312  + b15*v5 + b14*w5 + b18*v6 + b17*w6 &
313  + b21*v7 + b20*w7 + b24*v8 + b23*w8 &
314  + b27*v9 + b26*w9 + b30*v10 + b29*w10) * vx6inv
315  e13 = (b3*u1 + b1*w1 + b6*u2 + b4*w2 &
316  + b9*u3 + b7*w3 + b12*u4 + b10*w4 &
317  + b15*u5 + b13*w5 + b18*u6 + b16*w6 &
318  + b21*u7 + b19*w7 + b24*u8 + b22*w8 &
319  + b27*u9 + b25*w9 + b30*u10 + b28*w10) * vx6inv
320 
321  s11l(1,i) = e11*ci(1,j) + e22*ci(2,j) + e33*ci(4,j)
322  s22l(1,i) = e11*ci(2,j) + e22*ci(3,j) + e33*ci(5,j)
323  s33l(1,i) = e11*ci(4,j) + e22*ci(5,j) + e33*ci(6,j)
324  s12l(1,i) = e12*ci(7,j)
325  s23l(1,i) = e23*ci(8,j)
326  s13l(1,i) = e13*ci(9,j)
327 
328  r1 = s11l(1,i)*b1 + s12l(1,i)*b2 + s13l(1,i)*b3
329  r2 = s22l(1,i)*b2 + s12l(1,i)*b1 + s23l(1,i)*b3
330  r3 = s33l(1,i)*b3 + s23l(1,i)*b2 + s13l(1,i)*b1
331 
332  r4 = s11l(1,i)*b4 + s12l(1,i)*b5 + s13l(1,i)*b6
333  r5 = s22l(1,i)*b5 + s12l(1,i)*b4 + s23l(1,i)*b6
334  r6 = s33l(1,i)*b6 + s23l(1,i)*b5 + s13l(1,i)*b4
335 
336  r7 = s11l(1,i)*b7 + s12l(1,i)*b8 + s13l(1,i)*b9
337  r8 = s22l(1,i)*b8 + s12l(1,i)*b7 + s23l(1,i)*b9
338  r9 = s33l(1,i)*b9 + s23l(1,i)*b8 + s13l(1,i)*b7
339 
340  r10 = s11l(1,i)*b10 + s12l(1,i)*b11 + s13l(1,i)*b12
341  r11 = s22l(1,i)*b11 + s12l(1,i)*b10 + s23l(1,i)*b12
342  r12 = s33l(1,i)*b12 + s23l(1,i)*b11 + s13l(1,i)*b10
343 
344  r13 = s11l(1,i)*b13 + s12l(1,i)*b14 + s13l(1,i)*b15
345  r14 = s22l(1,i)*b14 + s12l(1,i)*b13 + s23l(1,i)*b15
346  r15 = s33l(1,i)*b15 + s23l(1,i)*b14 + s13l(1,i)*b13
347 
348  r16 = s11l(1,i)*b16 + s12l(1,i)*b17 + s13l(1,i)*b18
349  r17 = s22l(1,i)*b17 + s12l(1,i)*b16 + s23l(1,i)*b18
350  r18 = s33l(1,i)*b18 + s23l(1,i)*b17 + s13l(1,i)*b16
351 
352  r19 = s11l(1,i)*b19 + s12l(1,i)*b20 + s13l(1,i)*b21
353  r20 = s22l(1,i)*b20 + s12l(1,i)*b19 + s23l(1,i)*b21
354  r21 = s33l(1,i)*b21 + s23l(1,i)*b20 + s13l(1,i)*b19
355 
356  r22 = s11l(1,i)*b22 + s12l(1,i)*b23 + s13l(1,i)*b24
357  r23 = s22l(1,i)*b23 + s12l(1,i)*b22 + s23l(1,i)*b24
358  r24 = s33l(1,i)*b24 + s23l(1,i)*b23 + s13l(1,i)*b22
359 
360  r25 = s11l(1,i)*b25 + s12l(1,i)*b26 + s13l(1,i)*b27
361  r26 = s22l(1,i)*b26 + s12l(1,i)*b25 + s23l(1,i)*b27
362  r27 = s33l(1,i)*b27 + s23l(1,i)*b26 + s13l(1,i)*b25
363 
364  r28 = s11l(1,i)*b28 + s12l(1,i)*b29 + s13l(1,i)*b30
365  r29 = s22l(1,i)*b29 + s12l(1,i)*b28 + s23l(1,i)*b30
366  r30 = s33l(1,i)*b30 + s23l(1,i)*b29 + s13l(1,i)*b28
367 !-- Gauss Point 2
368 
369  g1 = 0.13819660d0
370  g2 = 0.58541020d0
371  g3 = 0.13819660d0
372  g4 = 0.13819660d0
373 
374  xn1 = (4.d0*g1-1.d0) ! derivative of shape function
375  xn2 = (4.d0*g2-1.d0) ! dN_i/dzeta_i
376  xn3 = (4.d0*g3-1.d0)
377  xn4 = (4.d0*g4-1.d0)
378 ! xN5 = 4.d0*g1*g2
379 ! xN6 = 4.d0*g2*g3
380 ! xN7 = 4.d0*g3*g1
381 ! xN8 = 4.d0*g1*g4
382 ! xN9 = 4.d0*g2*g4
383 ! xN10= 4.d0*g3*g4
384 
385  b1 = aux1*xn1
386  b2 = aux2*xn1
387  b3 = aux3*xn1
388  b4 = aux4*xn2
389  b5 = aux5*xn2
390  b6 = aux6*xn2
391  b7 = aux7*xn3
392  b8 = aux8*xn3
393  b9 = aux9*xn3
394  b10 = aux10*xn4
395  b11 = aux11*xn4
396  b12 = aux12*xn4
397  b13 = 4.d0*(g2*aux1 + g1*aux4)
398  b14 = 4.d0*(g2*aux2 + g1*aux5)
399  b15 = 4.d0*(g2*aux3 + g1*aux6)
400  b16 = 4.d0*(g3*aux4 + g2*aux7)
401  b17 = 4.d0*(g3*aux5 + g2*aux8)
402  b18 = 4.d0*(g3*aux6 + g2*aux9)
403  b19 = 4.d0*(g1*aux7 + g3*aux1)
404  b20 = 4.d0*(g1*aux8 + g3*aux2)
405  b21 = 4.d0*(g1*aux9 + g3*aux3)
406  b22 = 4.d0*(g4*aux1 + g1*aux10)
407  b23 = 4.d0*(g4*aux2 + g1*aux11)
408  b24 = 4.d0*(g4*aux3 + g1*aux12)
409  b25 = 4.d0*(g4*aux4 + g2*aux10)
410  b26 = 4.d0*(g4*aux5 + g2*aux11)
411  b27 = 4.d0*(g4*aux6 + g2*aux12)
412  b28 = 4.d0*(g4*aux7 + g3*aux10)
413  b29 = 4.d0*(g4*aux8 + g3*aux11)
414  b30 = 4.d0*(g4*aux9 + g3*aux12)
415 
416  e11 = (b1*u1 + b4*u2 + b7*u3 + b10*u4 + b13*u5 + &
417  b16*u6 + b19*u7 + b22*u8 + b25*u9 + b28*u10) * vx6inv
418  e22 = (b2*v1 + b5*v2 + b8*v3 + b11*v4 + b14*v5 + &
419  b17*v6 + b20*v7 + b23*v8 + b26*v9 + b29*v10) * vx6inv
420  e33 = (b3*w1 + b6*w2 + b9*w3 + b12*w4 + b15*w5 + &
421  b18*w6 + b21*w7 + b24*w8 + b27*w9 + b30*w10) * vx6inv
422  e12 = (b2*u1 + b1*v1 + b5*u2 + b4*v2 &
423  + b8*u3 + b7*v3 + b11*u4 + b10*v4 &
424  + b14*u5 + b13*v5 + b17*u6 + b16*v6 &
425  + b20*u7 + b19*v7 + b23*u8 + b22*v8 &
426  + b26*u9 + b25*v9 + b29*u10 + b28*v10) * vx6inv
427  e23 = (b3*v1 + b2*w1 + b6*v2 + b5*w2 &
428  + b9*v3 + b8*w3 + b12*v4 + b11*w4 &
429  + b15*v5 + b14*w5 + b18*v6 + b17*w6 &
430  + b21*v7 + b20*w7 + b24*v8 + b23*w8 &
431  + b27*v9 + b26*w9 + b30*v10 + b29*w10) * vx6inv
432  e13 = (b3*u1 + b1*w1 + b6*u2 + b4*w2 &
433  + b9*u3 + b7*w3 + b12*u4 + b10*w4 &
434  + b15*u5 + b13*w5 + b18*u6 + b16*w6 &
435  + b21*u7 + b19*w7 + b24*u8 + b22*w8 &
436  + b27*u9 + b25*w9 + b30*u10 + b28*w10) * vx6inv
437 
438  s11l(2,i) = e11*ci(1,j) + e22*ci(2,j) + e33*ci(4,j)
439  s22l(2,i) = e11*ci(2,j) + e22*ci(3,j) + e33*ci(5,j)
440  s33l(2,i) = e11*ci(4,j) + e22*ci(5,j) + e33*ci(6,j)
441  s12l(2,i) = e12*ci(7,j)
442  s23l(2,i) = e23*ci(8,j)
443  s13l(2,i) = e13*ci(9,j)
444 
445  r1 = r1 +s11l(2,i)*b1 + s12l(2,i)*b2 + s13l(2,i)*b3
446  r2 = r2 +s22l(2,i)*b2 + s12l(2,i)*b1 + s23l(2,i)*b3
447  r3 = r3 +s33l(2,i)*b3 + s23l(2,i)*b2 + s13l(2,i)*b1
448  r4 = r4 +s11l(2,i)*b4 + s12l(2,i)*b5 + s13l(2,i)*b6
449  r5 = r5 +s22l(2,i)*b5 + s12l(2,i)*b4 + s23l(2,i)*b6
450  r6 = r6 +s33l(2,i)*b6 + s23l(2,i)*b5 + s13l(2,i)*b4
451  r7 = r7 +s11l(2,i)*b7 + s12l(2,i)*b8 + s13l(2,i)*b9
452  r8 = r8 +s22l(2,i)*b8 + s12l(2,i)*b7 + s23l(2,i)*b9
453  r9 = r9 +s33l(2,i)*b9 + s23l(2,i)*b8 + s13l(2,i)*b7
454  r10 = r10 +s11l(2,i)*b10 + s12l(2,i)*b11 + s13l(2,i)*b12
455  r11 = r11 +s22l(2,i)*b11 + s12l(2,i)*b10 + s23l(2,i)*b12
456  r12 = r12 +s33l(2,i)*b12 + s23l(2,i)*b11 + s13l(2,i)*b10
457  r13 = r13 +s11l(2,i)*b13 + s12l(2,i)*b14 + s13l(2,i)*b15
458  r14 = r14 +s22l(2,i)*b14 + s12l(2,i)*b13 + s23l(2,i)*b15
459  r15 = r15 +s33l(2,i)*b15 + s23l(2,i)*b14 + s13l(2,i)*b13
460  r16 = r16 +s11l(2,i)*b16 + s12l(2,i)*b17 + s13l(2,i)*b18
461  r17 = r17 +s22l(2,i)*b17 + s12l(2,i)*b16 + s23l(2,i)*b18
462  r18 = r18 +s33l(2,i)*b18 + s23l(2,i)*b17 + s13l(2,i)*b16
463  r19 = r19 +s11l(2,i)*b19 + s12l(2,i)*b20 + s13l(2,i)*b21
464  r20 = r20 +s22l(2,i)*b20 + s12l(2,i)*b19 + s23l(2,i)*b21
465  r21 = r21 +s33l(2,i)*b21 + s23l(2,i)*b20 + s13l(2,i)*b19
466  r22 = r22 +s11l(2,i)*b22 + s12l(2,i)*b23 + s13l(2,i)*b24
467  r23 = r23 +s22l(2,i)*b23 + s12l(2,i)*b22 + s23l(2,i)*b24
468  r24 = r24 +s33l(2,i)*b24 + s23l(2,i)*b23 + s13l(2,i)*b22
469  r25 = r25 +s11l(2,i)*b25 + s12l(2,i)*b26 + s13l(2,i)*b27
470  r26 = r26 +s22l(2,i)*b26 + s12l(2,i)*b25 + s23l(2,i)*b27
471  r27 = r27 +s33l(2,i)*b27 + s23l(2,i)*b26 + s13l(2,i)*b25
472  r28 = r28 +s11l(2,i)*b28 + s12l(2,i)*b29 + s13l(2,i)*b30
473  r29 = r29 +s22l(2,i)*b29 + s12l(2,i)*b28 + s23l(2,i)*b30
474  r30 = r30 +s33l(2,i)*b30 + s23l(2,i)*b29 + s13l(2,i)*b28
475 !-- Gauss Point 3
476 
477  g1 = 0.13819660d0
478  g2 = 0.13819660d0
479  g3 = 0.58541020d0
480  g4 = 0.13819660d0
481 
482  xn1 = (4.d0*g1-1.d0) ! derivative of shape function
483  xn2 = (4.d0*g2-1.d0) ! dN_i/dzeta_i
484  xn3 = (4.d0*g3-1.d0)
485  xn4 = (4.d0*g4-1.d0)
486 ! xN5 = 4.d0*g1*g2
487 ! xN6 = 4.d0*g2*g3
488 ! xN7 = 4.d0*g3*g1
489 ! xN8 = 4.d0*g1*g4
490 ! xN9 = 4.d0*g2*g4
491 ! xN10= 4.d0*g3*g4
492 
493  b1 = aux1*xn1
494  b2 = aux2*xn1
495  b3 = aux3*xn1
496  b4 = aux4*xn2
497  b5 = aux5*xn2
498  b6 = aux6*xn2
499  b7 = aux7*xn3
500  b8 = aux8*xn3
501  b9 = aux9*xn3
502  b10 = aux10*xn4
503  b11 = aux11*xn4
504  b12 = aux12*xn4
505  b13 = 4.d0*(g2*aux1 + g1*aux4)
506  b14 = 4.d0*(g2*aux2 + g1*aux5)
507  b15 = 4.d0*(g2*aux3 + g1*aux6)
508  b16 = 4.d0*(g3*aux4 + g2*aux7)
509  b17 = 4.d0*(g3*aux5 + g2*aux8)
510  b18 = 4.d0*(g3*aux6 + g2*aux9)
511  b19 = 4.d0*(g1*aux7 + g3*aux1)
512  b20 = 4.d0*(g1*aux8 + g3*aux2)
513  b21 = 4.d0*(g1*aux9 + g3*aux3)
514  b22 = 4.d0*(g4*aux1 + g1*aux10)
515  b23 = 4.d0*(g4*aux2 + g1*aux11)
516  b24 = 4.d0*(g4*aux3 + g1*aux12)
517  b25 = 4.d0*(g4*aux4 + g2*aux10)
518  b26 = 4.d0*(g4*aux5 + g2*aux11)
519  b27 = 4.d0*(g4*aux6 + g2*aux12)
520  b28 = 4.d0*(g4*aux7 + g3*aux10)
521  b29 = 4.d0*(g4*aux8 + g3*aux11)
522  b30 = 4.d0*(g4*aux9 + g3*aux12)
523 
524  e11 = (b1*u1 + b4*u2 + b7*u3 + b10*u4 + b13*u5 + &
525  b16*u6 + b19*u7 + b22*u8 + b25*u9 + b28*u10) * vx6inv
526  e22 = (b2*v1 + b5*v2 + b8*v3 + b11*v4 + b14*v5 + &
527  b17*v6 + b20*v7 + b23*v8 + b26*v9 + b29*v10) * vx6inv
528  e33 = (b3*w1 + b6*w2 + b9*w3 + b12*w4 + b15*w5 + &
529  b18*w6 + b21*w7 + b24*w8 + b27*w9 + b30*w10) * vx6inv
530  e12 = (b2*u1 + b1*v1 + b5*u2 + b4*v2 &
531  + b8*u3 + b7*v3 + b11*u4 + b10*v4 &
532  + b14*u5 + b13*v5 + b17*u6 + b16*v6 &
533  + b20*u7 + b19*v7 + b23*u8 + b22*v8 &
534  + b26*u9 + b25*v9 + b29*u10 + b28*v10) * vx6inv
535  e23 = (b3*v1 + b2*w1 + b6*v2 + b5*w2 &
536  + b9*v3 + b8*w3 + b12*v4 + b11*w4 &
537  + b15*v5 + b14*w5 + b18*v6 + b17*w6 &
538  + b21*v7 + b20*w7 + b24*v8 + b23*w8 &
539  + b27*v9 + b26*w9 + b30*v10 + b29*w10) * vx6inv
540  e13 = (b3*u1 + b1*w1 + b6*u2 + b4*w2 &
541  + b9*u3 + b7*w3 + b12*u4 + b10*w4 &
542  + b15*u5 + b13*w5 + b18*u6 + b16*w6 &
543  + b21*u7 + b19*w7 + b24*u8 + b22*w8 &
544  + b27*u9 + b25*w9 + b30*u10 + b28*w10) * vx6inv
545 
546  s11l(3,i) = e11*ci(1,j) + e22*ci(2,j) + e33*ci(4,j)
547  s22l(3,i) = e11*ci(2,j) + e22*ci(3,j) + e33*ci(5,j)
548  s33l(3,i) = e11*ci(4,j) + e22*ci(5,j) + e33*ci(6,j)
549  s12l(3,i) = e12*ci(7,j)
550  s23l(3,i) = e23*ci(8,j)
551  s13l(3,i) = e13*ci(9,j)
552 
553  r1 = r1 + s11l(3,i)*b1 + s12l(3,i)*b2 + s13l(3,i)*b3
554  r2 = r2 + s22l(3,i)*b2 + s12l(3,i)*b1 + s23l(3,i)*b3
555  r3 = r3 +s33l(3,i)*b3 + s23l(3,i)*b2 + s13l(3,i)*b1
556  r4 = r4 +s11l(3,i)*b4 + s12l(3,i)*b5 + s13l(3,i)*b6
557  r5 = r5 +s22l(3,i)*b5 + s12l(3,i)*b4 + s23l(3,i)*b6
558  r6 = r6 +s33l(3,i)*b6 + s23l(3,i)*b5 + s13l(3,i)*b4
559  r7 = r7 +s11l(3,i)*b7 + s12l(3,i)*b8 + s13l(3,i)*b9
560  r8 = r8 +s22l(3,i)*b8 + s12l(3,i)*b7 + s23l(3,i)*b9
561  r9 = r9 +s33l(3,i)*b9 + s23l(3,i)*b8 + s13l(3,i)*b7
562  r10 = r10 +s11l(3,i)*b10 + s12l(3,i)*b11 + s13l(3,i)*b12
563  r11 = r11 +s22l(3,i)*b11 + s12l(3,i)*b10 + s23l(3,i)*b12
564  r12 = r12 +s33l(3,i)*b12 + s23l(3,i)*b11 + s13l(3,i)*b10
565  r13 = r13 +s11l(3,i)*b13 + s12l(3,i)*b14 + s13l(3,i)*b15
566  r14 = r14 +s22l(3,i)*b14 + s12l(3,i)*b13 + s23l(3,i)*b15
567  r15 = r15 +s33l(3,i)*b15 + s23l(3,i)*b14 + s13l(3,i)*b13
568  r16 = r16 +s11l(3,i)*b16 + s12l(3,i)*b17 + s13l(3,i)*b18
569  r17 = r17 +s22l(3,i)*b17 + s12l(3,i)*b16 + s23l(3,i)*b18
570  r18 = r18 +s33l(3,i)*b18 + s23l(3,i)*b17 + s13l(3,i)*b16
571  r19 = r19 +s11l(3,i)*b19 + s12l(3,i)*b20 + s13l(3,i)*b21
572  r20 = r20 +s22l(3,i)*b20 + s12l(3,i)*b19 + s23l(3,i)*b21
573  r21 = r21 +s33l(3,i)*b21 + s23l(3,i)*b20 + s13l(3,i)*b19
574  r22 = r22 +s11l(3,i)*b22 + s12l(3,i)*b23 + s13l(3,i)*b24
575  r23 = r23 +s22l(3,i)*b23 + s12l(3,i)*b22 + s23l(3,i)*b24
576  r24 = r24 +s33l(3,i)*b24 + s23l(3,i)*b23 + s13l(3,i)*b22
577  r25 = r25 +s11l(3,i)*b25 + s12l(3,i)*b26 + s13l(3,i)*b27
578  r26 = r26 +s22l(3,i)*b26 + s12l(3,i)*b25 + s23l(3,i)*b27
579  r27 = r27 +s33l(3,i)*b27 + s23l(3,i)*b26 + s13l(3,i)*b25
580  r28 = r28 +s11l(3,i)*b28 + s12l(3,i)*b29 + s13l(3,i)*b30
581  r29 = r29 +s22l(3,i)*b29 + s12l(3,i)*b28 + s23l(3,i)*b30
582  r30 = r30 +s33l(3,i)*b30 + s23l(3,i)*b29 + s13l(3,i)*b28
583 !-- Gauss Point 4
584 
585  g1 = 0.13819660d0
586  g2 = 0.13819660d0
587  g3 = 0.13819660d0
588  g4 = 0.58541020d0
589 
590  xn1 = (4.d0*g1-1.d0) ! derivative of shape function
591  xn2 = (4.d0*g2-1.d0) ! dN_i/dzeta_i
592  xn3 = (4.d0*g3-1.d0)
593  xn4 = (4.d0*g4-1.d0)
594 ! xN5 = 4.d0*g1*g2
595 ! xN6 = 4.d0*g2*g3
596 ! xN7 = 4.d0*g3*g1
597 ! xN8 = 4.d0*g1*g4
598 ! xN9 = 4.d0*g2*g4
599 ! xN10= 4.d0*g3*g4
600 
601  b1 = aux1*xn1
602  b2 = aux2*xn1
603  b3 = aux3*xn1
604  b4 = aux4*xn2
605  b5 = aux5*xn2
606  b6 = aux6*xn2
607  b7 = aux7*xn3
608  b8 = aux8*xn3
609  b9 = aux9*xn3
610  b10 = aux10*xn4
611  b11 = aux11*xn4
612  b12 = aux12*xn4
613  b13 = 4.d0*(g2*aux1 + g1*aux4)
614  b14 = 4.d0*(g2*aux2 + g1*aux5)
615  b15 = 4.d0*(g2*aux3 + g1*aux6)
616  b16 = 4.d0*(g3*aux4 + g2*aux7)
617  b17 = 4.d0*(g3*aux5 + g2*aux8)
618  b18 = 4.d0*(g3*aux6 + g2*aux9)
619  b19 = 4.d0*(g1*aux7 + g3*aux1)
620  b20 = 4.d0*(g1*aux8 + g3*aux2)
621  b21 = 4.d0*(g1*aux9 + g3*aux3)
622  b22 = 4.d0*(g4*aux1 + g1*aux10)
623  b23 = 4.d0*(g4*aux2 + g1*aux11)
624  b24 = 4.d0*(g4*aux3 + g1*aux12)
625  b25 = 4.d0*(g4*aux4 + g2*aux10)
626  b26 = 4.d0*(g4*aux5 + g2*aux11)
627  b27 = 4.d0*(g4*aux6 + g2*aux12)
628  b28 = 4.d0*(g4*aux7 + g3*aux10)
629  b29 = 4.d0*(g4*aux8 + g3*aux11)
630  b30 = 4.d0*(g4*aux9 + g3*aux12)
631 
632  e11 = (b1*u1 + b4*u2 + b7*u3 + b10*u4 + b13*u5 + &
633  b16*u6 + b19*u7 + b22*u8 + b25*u9 + b28*u10) * vx6inv
634  e22 = (b2*v1 + b5*v2 + b8*v3 + b11*v4 + b14*v5 + &
635  b17*v6 + b20*v7 + b23*v8 + b26*v9 + b29*v10) * vx6inv
636  e33 = (b3*w1 + b6*w2 + b9*w3 + b12*w4 + b15*w5 + &
637  b18*w6 + b21*w7 + b24*w8 + b27*w9 + b30*w10) * vx6inv
638  e12 = (b2*u1 + b1*v1 + b5*u2 + b4*v2 &
639  + b8*u3 + b7*v3 + b11*u4 + b10*v4 &
640  + b14*u5 + b13*v5 + b17*u6 + b16*v6 &
641  + b20*u7 + b19*v7 + b23*u8 + b22*v8 &
642  + b26*u9 + b25*v9 + b29*u10 + b28*v10) * vx6inv
643  e23 = (b3*v1 + b2*w1 + b6*v2 + b5*w2 &
644  + b9*v3 + b8*w3 + b12*v4 + b11*w4 &
645  + b15*v5 + b14*w5 + b18*v6 + b17*w6 &
646  + b21*v7 + b20*w7 + b24*v8 + b23*w8 &
647  + b27*v9 + b26*w9 + b30*v10 + b29*w10) * vx6inv
648  e13 = (b3*u1 + b1*w1 + b6*u2 + b4*w2 &
649  + b9*u3 + b7*w3 + b12*u4 + b10*w4 &
650  + b15*u5 + b13*w5 + b18*u6 + b16*w6 &
651  + b21*u7 + b19*w7 + b24*u8 + b22*w8 &
652  + b27*u9 + b25*w9 + b30*u10 + b28*w10) * vx6inv
653 
654  s11l(4,i) = e11*ci(1,j) + e22*ci(2,j) + e33*ci(4,j)
655  s22l(4,i) = e11*ci(2,j) + e22*ci(3,j) + e33*ci(5,j)
656  s33l(4,i) = e11*ci(4,j) + e22*ci(5,j) + e33*ci(6,j)
657  s12l(4,i) = e12*ci(7,j)
658  s23l(4,i) = e23*ci(8,j)
659  s13l(4,i) = e13*ci(9,j)
660 
661  r1 = r1 + s11l(4,i)*b1 + s12l(4,i)*b2 + s13l(4,i)*b3
662  r2 = r2 + s22l(4,i)*b2 + s12l(4,i)*b1 + s23l(4,i)*b3
663  r3 = r3 + s33l(4,i)*b3 + s23l(4,i)*b2 + s13l(4,i)*b1
664  r4 = r4 + s11l(4,i)*b4 + s12l(4,i)*b5 + s13l(4,i)*b6
665  r5 = r5 + s22l(4,i)*b5 + s12l(4,i)*b4 + s23l(4,i)*b6
666  r6 = r6 + s33l(4,i)*b6 + s23l(4,i)*b5 + s13l(4,i)*b4
667  r7 = r7 + s11l(4,i)*b7 + s12l(4,i)*b8 + s13l(4,i)*b9
668  r8 = r8 + s22l(4,i)*b8 + s12l(4,i)*b7 + s23l(4,i)*b9
669  r9 = r9 + s33l(4,i)*b9 + s23l(4,i)*b8 + s13l(4,i)*b7
670  r10 = r10 +s11l(4,i)*b10 + s12l(4,i)*b11 + s13l(4,i)*b12
671  r11 = r11 +s22l(4,i)*b11 + s12l(4,i)*b10 + s23l(4,i)*b12
672  r12 = r12 +s33l(4,i)*b12 + s23l(4,i)*b11 + s13l(4,i)*b10
673  r13 = r13 +s11l(4,i)*b13 + s12l(4,i)*b14 + s13l(4,i)*b15
674  r14 = r14 +s22l(4,i)*b14 + s12l(4,i)*b13 + s23l(4,i)*b15
675  r15 = r15 +s33l(4,i)*b15 + s23l(4,i)*b14 + s13l(4,i)*b13
676  r16 = r16 +s11l(4,i)*b16 + s12l(4,i)*b17 + s13l(4,i)*b18
677  r17 = r17 +s22l(4,i)*b17 + s12l(4,i)*b16 + s23l(4,i)*b18
678  r18 = r18 +s33l(4,i)*b18 + s23l(4,i)*b17 + s13l(4,i)*b16
679  r19 = r19 +s11l(4,i)*b19 + s12l(4,i)*b20 + s13l(4,i)*b21
680  r20 = r20 +s22l(4,i)*b20 + s12l(4,i)*b19 + s23l(4,i)*b21
681  r21 = r21 +s33l(4,i)*b21 + s23l(4,i)*b20 + s13l(4,i)*b19
682  r22 = r22 +s11l(4,i)*b22 + s12l(4,i)*b23 + s13l(4,i)*b24
683  r23 = r23 +s22l(4,i)*b23 + s12l(4,i)*b22 + s23l(4,i)*b24
684  r24 = r24 +s33l(4,i)*b24 + s23l(4,i)*b23 + s13l(4,i)*b22
685  r25 = r25 +s11l(4,i)*b25 + s12l(4,i)*b26 + s13l(4,i)*b27
686  r26 = r26 +s22l(4,i)*b26 + s12l(4,i)*b25 + s23l(4,i)*b27
687  r27 = r27 +s33l(4,i)*b27 + s23l(4,i)*b26 + s13l(4,i)*b25
688  r28 = r28 +s11l(4,i)*b28 + s12l(4,i)*b29 + s13l(4,i)*b30
689  r29 = r29 +s22l(4,i)*b29 + s12l(4,i)*b28 + s23l(4,i)*b30
690  r30 = r30 +s33l(4,i)*b30 + s23l(4,i)*b29 + s13l(4,i)*b28
691 
692 ! Wi (i.e. weight) for 4 guass integration is 1/4
693 
694 ! Wi * 1/6 because the volume of a reference tetrahedra in
695 ! volume coordinates is 1/6
696 
697 ! ASSEMBLE THE INTERNAL FORCE VECTOR
698 !
699 ! local node 1
700  r_in(k1n1) = r_in(k1n1) - r1*0.04166666666666667d0
701  r_in(k2n1) = r_in(k2n1) - r2*0.04166666666666667d0
702  r_in(k3n1) = r_in(k3n1) - r3*0.04166666666666667d0
703 ! local node 2
704  r_in(k1n2) = r_in(k1n2) - r4*0.04166666666666667d0
705  r_in(k2n2) = r_in(k2n2) - r5*0.04166666666666667d0
706  r_in(k3n2) = r_in(k3n2) - r6*0.04166666666666667d0
707 ! local node 3
708  r_in(k1n3) = r_in(k1n3) - r7*0.04166666666666667d0
709  r_in(k2n3) = r_in(k2n3) - r8*0.04166666666666667d0
710  r_in(k3n3) = r_in(k3n3) - r9*0.04166666666666667d0
711 ! local node 4
712  r_in(k1n4) = r_in(k1n4) - r10*0.04166666666666667d0
713  r_in(k2n4) = r_in(k2n4) - r11*0.04166666666666667d0
714  r_in(k3n4) = r_in(k3n4) - r12*0.04166666666666667d0
715 ! local node 5
716  r_in(k1n5) = r_in(k1n5) - r13*0.04166666666666667d0
717  r_in(k2n5) = r_in(k2n5) - r14*0.04166666666666667d0
718  r_in(k3n5) = r_in(k3n5) - r15*0.04166666666666667d0
719 ! local node 6
720  r_in(k1n6) = r_in(k1n6) - r16*0.04166666666666667d0
721  r_in(k2n6) = r_in(k2n6) - r17*0.04166666666666667d0
722  r_in(k3n6) = r_in(k3n6) - r18*0.04166666666666667d0
723 ! local node 7
724  r_in(k1n7) = r_in(k1n7) - r19*0.04166666666666667d0
725  r_in(k2n7) = r_in(k2n7) - r20*0.04166666666666667d0
726  r_in(k3n7) = r_in(k3n7) - r21*0.04166666666666667d0
727 ! local node 8
728  r_in(k1n8) = r_in(k1n8) - r22*0.04166666666666667d0
729  r_in(k2n8) = r_in(k2n8) - r23*0.04166666666666667d0
730  r_in(k3n8) = r_in(k3n8) - r24*0.04166666666666667d0
731 ! local node 9
732  r_in(k1n9) = r_in(k1n9) - r25*0.04166666666666667d0
733  r_in(k2n9) = r_in(k2n9) - r26*0.04166666666666667d0
734  r_in(k3n9) = r_in(k3n9) - r27*0.04166666666666667d0
735 ! local node 10
736  r_in(k1n10) = r_in(k1n10) - r28*0.04166666666666667d0
737  r_in(k2n10) = r_in(k2n10) - r29*0.04166666666666667d0
738  r_in(k3n10) = r_in(k3n10) - r30*0.04166666666666667d0
739  ENDDO
740 
741  RETURN
742 END SUBROUTINE v3d10
743 
744 
745 
746 
subroutine v3d10(coor, matcstet, lmcstet, R_in, d, ci, S11l, S22l, S33l, S12l, S23l, S13l, numnp, nstart, nend, numcstet, numat_vol)
Definition: v3d10.f90:53
const NT & d
blockLoc i
Definition: read.cpp:79
j indices j
Definition: Indexing.h:6