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