Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MP/Source/v3d4_NeoHookeanCompress.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 v3d4_neohookeancompress(coor,matcstet,lmcstet,R_in,d,ci, &
54  s11,s22,s33,s12,s23,s13,numnp,nstart,nend,numcstet,numat_vol, &
55  xmu,xlamda)
56 
57 !________________________________________________________________________
58 !
59 ! V3D4 - Performs displacement based computations for Volumetric 3D
60 ! 4-node tetrahedra large deformation (i.e. nonlinear) elastic
61 ! elements with linear interpolation functions. (constant strain
62 ! tetrahedra). Returns the internal force vector R_in.
63 !
64 ! DATE: 04.2000 AUTHOR: SCOT BREITENFELD
65 !
66 ! Source:
67 ! "Nonlinear continuum mechanics for finite element analysis"
68 ! Javier Bonet and Richard D. Wood
69 ! ISBN 0-521-57272-X
70 
71 !________________________________________________________________________
72 
73  IMPLICIT NONE
74 !---- Global variables
75  INTEGER :: numnp ! number of nodal points
76  INTEGER :: numcstet ! number of CSTet elements
77  INTEGER :: numat_vol ! number of volumetric materials
78 !-- coordinate array
79  REAL*8, DIMENSION(1:3,1:numnp) :: coor
80 !-- elastic stiffness consts
81  REAL*8, DIMENSION(1:9,1:numat_vol) :: ci
82 !-- internal force
83  REAL*8, DIMENSION(1:3*numnp) :: r_in
84 !-- displacement vector
85  REAL*8, DIMENSION(1:3*numnp) :: d
86 !-- CSTet stress
87  REAL*8, DIMENSION(1:numcstet) :: s11, s22, s33, s12, s23, s13
88 !-- x, y and z displacements of nodes
89  REAL*8 :: u1,u2,u3,u4,v1,v2,v3,v4,w1,w2,w3,w4
90 !-- coordinate holding variable
91  REAL*8 :: x1,x2,x3,x4,y1,y2,y3,y4,z1,z2,z3,z4
92 !-- Coordinate subtractions
93  REAL*8 :: x14, x24, x34, y14, y24, y34, z14, z24, z34
94 !--
95  REAL*8 :: c11, c21, c31
96 !-- 6*volume, inverse(6*volume), and the volume
97  REAL*8 :: vx6, vx6inv, vol
98 !-- spacial derivatives
99  REAL*8 :: b1,b2,b3,b4,b5,b6,b7,b8,b9,b10,b11,b12
100 !-- partial derivatives of the displacement
101  REAL*8 :: dudx,dvdy,dwdz,dudy,dvdx,dvdz,dwdy,dudz,dwdx
102 !-- strains
103  REAL*8 :: e11,e22,e33,e12,e23,e13
104 ! -- connectivity table for CSTet
105  INTEGER, DIMENSION(1:4,1:numcstet) :: lmcstet
106 ! -- mat number for CST element
107  INTEGER, DIMENSION(1:numcstet) :: matcstet
108 ! -- node numbers
109  INTEGER :: n1,n2,n3,n4
110 ! -- dummy and counters
111  INTEGER :: i,j,nstart,nend
112  INTEGER :: k1n1,k1n2,k1n3,k1n4,k2n1,k2n2,k2n3,k2n4
113  INTEGER :: k3n1,k3n2,k3n3,k3n4
114 ! --
115  REAL*8 :: f11, f12, f13, f21, f22, f23, f31, f32, f33
116  REAL*8 :: j1,j2,j2inv
117  REAL*8 :: c11, c12, c13, c21, c22, c23, c31, c32, c33
118  REAL*8 :: xmu, xlambda
119 
120  DO i = nstart, nend
121  j = matcstet(i)
122 
123  n1 = lmcstet(1,i)
124  n2 = lmcstet(2,i)
125  n3 = lmcstet(3,i)
126  n4 = lmcstet(4,i)
127 
128  k3n1 = 3*n1
129  k3n2 = 3*n2
130  k3n3 = 3*n3
131  k3n4 = 3*n4
132 
133  k2n1 = k3n1 - 1
134  k2n2 = k3n2 - 1
135  k2n3 = k3n3 - 1
136  k2n4 = k3n4 - 1
137 
138  k1n1 = k3n1 - 2
139  k1n2 = k3n2 - 2
140  k1n3 = k3n3 - 2
141  k1n4 = k3n4 - 2
142  ! k#n# dummy variables replaces:
143  u1 = d(k1n1) ! (3*n1-2)
144  u2 = d(k1n2) ! (3*n2-2)
145  u3 = d(k1n3) ! (3*n3-2)
146  u4 = d(k1n4) ! (3*n4-2)
147  v1 = d(k2n1) ! (3*n1-1)
148  v2 = d(k2n2) ! (3*n2-1)
149  v3 = d(k2n3) ! (3*n3-1)
150  v4 = d(k2n4) ! (3*n4-1)
151  w1 = d(k3n1) ! (3*n1)
152  w2 = d(k3n2) ! (3*n2)
153  w3 = d(k3n3) ! (3*n3)
154  w4 = d(k3n4) ! (3*n4)
155 
156  x1 = coor(1,n1)
157  x2 = coor(1,n2)
158  x3 = coor(1,n3)
159  x4 = coor(1,n4)
160  y1 = coor(2,n1)
161  y2 = coor(2,n2)
162  y3 = coor(2,n3)
163  y4 = coor(2,n4)
164  z1 = coor(3,n1)
165  z2 = coor(3,n2)
166  z3 = coor(3,n3)
167  z4 = coor(3,n4)
168 
169 
170  x14 = x1 - x4
171  x24 = x2 - x4
172  x34 = x3 - x4
173  y14 = y1 - y4
174  y24 = y2 - y4
175  y34 = y3 - y4
176  z14 = z1 - z4
177  z24 = z2 - z4
178  z34 = z3 - z4
179 
180  c11 = y24*z34 - z24*y34
181  c21 = -( x24*z34 - z24*x34 )
182  c31 = x24*y34 - y24*x34
183 
184  vx6 = -( x14*c11 + y14*c21 + z14*c31 )
185 
186 ! See the maple worksheet 'V3D4.mws' for the derivation of Vx6
187 
188 ! Vx6 = x2*y3*z4 - x2*y4*z3 - y2*x3*z4 + y2*x4*z3 + z2*x3*y4
189 ! $ - z2*x4*y3 - x1*y3*z4 + x1*y4*z3 + x1*y2*z4 - x1*y2*z3
190 ! $ - x1*z2*y4 + x1*z2*y3 + y1*x3*z4 - y1*x4*z3 - y1*x2*z4
191 ! $ + y1*x2*z3 + y1*z2*x4 - y1*z2*x3 - z1*x3*y4 + z1*x4*y3
192 ! $ + z1*x2*y4 - z1*x2*y3 - z1*y2*x4 + z1*y2*x3
193 
194  vx6inv = 1.d0 / vx6
195 
196 ! See the maple worksheet 'V3D4.mws' for the derivation of [B]
197 
198  b1 = -(y3*z4 - y4*z3 - y2*z4 + y2*z3 + z2*y4 - z2*y3) * vx6inv
199  b2 = (x3*z4 - x4*z3 - x2*z4 + x2*z3 + z2*x4 - z2*x3) * vx6inv
200  b3 = -(x3*y4 - x4*y3 - x2*y4 + x2*y3 + y2*x4 - y2*x3) * vx6inv
201  b4 = (y3*z4 - y4*z3 - y1*z4 + y1*z3 + z1*y4 - z1*y3) * vx6inv
202  b5 = -(x3*z4 - x4*z3 - x1*z4 + x1*z3 + z1*x4 - z1*x3) * vx6inv
203  b6 = (x3*y4 - x4*y3 - x1*y4 + x1*y3 + y1*x4 - y1*x3) * vx6inv
204  b7 = -(y2*z4 - z2*y4 - y1*z4 + y1*z2 + z1*y4 - z1*y2) * vx6inv
205  b8 = (x2*z4 - z2*x4 - x1*z4 + x1*z2 + z1*x4 - z1*x2) * vx6inv
206  b9 = -(x2*y4 - y2*x4 - x1*y4 + x1*y2 + y1*x4 - y1*x2) * vx6inv
207  b10 = (y2*z3 - z2*y3 - y1*z3 + y1*z2 + z1*y3 - z1*y2) * vx6inv
208  b11 = -(x2*z3 - z2*x3 - x1*z3 + x1*z2 + z1*x3 - z1*x2) * vx6inv
209  b12 = (x2*y3 - y2*x3 - x1*y3 + x1*y2 + y1*x3 - y1*x2) * vx6inv
210 !
211 ! Calculate displacement gradient (H)
212 !
213  dudx = b1*u1 + b4*u2 + b7*u3 + b10*u4
214  dvdy = b2*v1 + b5*v2 + b8*v3 + b11*v4
215  dwdz = b3*w1 + b6*w2 + b9*w3 + b12*w4
216  dudy = b2*u1 + b5*u2 + b8*u3 + b11*u4
217  dvdx = b1*v1 + b4*v2 + b7*v3 + b10*v4
218  dvdz = b3*v1 + b6*v2 + b9*v3 + b12*v4
219  dwdy = b2*w1 + b5*w2 + b8*w3 + b11*w4
220  dudz = b3*u1 + b6*u2 + b9*u3 + b12*u4
221  dwdx = b1*w1 + b4*w2 + b7*w3 + b10*w4
222 !
223 ! deformation gradients F
224 !
225  f11 = 1.d0 + dudx
226  f12 = dudy
227  f13 = dudz
228  f21 = dvdx
229  f22 = 1.d0 + dvdy
230  f23 = dvdz
231  f31 = dwdx
232  f32 = dwdy
233  f33 = 1.d0 + dwdz
234 
235  j1 = f11*(f22*f33-f23*f32)+f12*(f31*f23-f21*f33)+f13*(-f31*f22+f21*f32)
236  IF(j1.LE.0.d0) THEN
237  WRITE(*,*)'area has become zero for element',i
238  stop
239  ENDIF
240  j2 = j1*j1
241 
242 ! T
243 ! C = F F = right Cauchy-Green deformation tensor
244 !
245 ! Divided by third invariant squared
246 
247  j2inv = 1.d0/j2
248 
249  c11 = (f11*f11+f21*f21+f31*f31)*j2inv
250  c12 = (f11*f12+f21*f22+f31*f32)*j2inv
251  c13 = (f11*f13+f21*f23+f31*f33)*j2inv
252  c21 = c12
253  c22 = (f12*f12+f22*f22+f32*f32)*j2inv
254  c23 = (f12*f13+f22*f23+f32*f33)*j2inv
255  c31 = c13
256  c32 = c23
257  c33 = (f13*f13+f23*f23+f33*f33)*j2inv
258 
259 !
260 ! Second Piola-Kirchoff tensor
261 ! Eq. (5.28), pg. 124
262 
263  s11(i) = xmu*(1.d0-(c22*c33-c23*c32)) + xlambda*log(j1)*(c22*c33-c23*c32)
264  s22(i) = xmu*(1.d0-(c11*c33-c31*c13)) + xlambda*log(j1)*(c11*c33-c31*c13)
265  s33(i) = xmu*(1.d0-(c11*c22-c12*c21)) + xlambda*log(j1)*(c11*c22-c12*c21)
266  s12(i) = (-c12*c33+c13*c32)*(xmu - xlambda*log(j1))
267  s23(i) = (-c11*c23+c13*c21)*(xmu - xlambda*log(j1))
268  s13(i) = (c12*c23-c13*c22)*(xmu - xlambda*log(j1))
269 
270 
271 ! calculate the volume
272 
273  vol = vx6 / 6.d0
274 
275 ! ASSEMBLE THE INTERNAL FORCE VECTOR
276 !
277 ! local node 1
278  r_in(k1n1) = r_in(k1n1) - vol* &
279  ( s11(i)*b1*(1.d0+dudx) + s22(i)*b2*dudy + s33(i)*b3*dudz &
280  + s12(i)*( b2*(1.d0+dudx) + b1*dudy ) &
281  + s23(i)*( b3*dudy + b2*dudz ) &
282  + s13(i)*( b3*(1.d0+dudx) + b1*dudz ) )
283  r_in(k2n1) = r_in(k2n1) - vol* &
284  ( s11(i)*b1*dvdx + s22(i)*b2*(1.d0+dvdy) + s33(i)*b3*dvdz &
285  + s12(i)*( b1*(1.d0+dvdy) + b2*dvdx ) &
286  + s23(i)*( b3*(1.d0+dvdy) + b2*dvdz ) &
287  + s13(i)*( b3*dvdx + b1*dvdz ) )
288  r_in(k3n1) = r_in(k3n1) - vol* &
289  ( s11(i)*b1*dwdx + s22(i)*b2*dwdy + s33(i)*b3*(1.d0+dwdz) &
290  + s12(i)*( b2*dwdx + b1*dwdy ) &
291  + s23(i)*( b3*dwdy + b2*(1.d0 + dwdz) ) &
292  + s13(i)*( b3*dwdx + b1*(1.d0 + dwdz) ) )
293 ! local node 2
294  r_in(k1n2) = r_in(k1n2) - vol* &
295  ( s11(i)*b4*(1.d0+dudx) + s22(i)*b5*dudy + s33(i)*b6*dudz &
296  + s12(i)*( b5*(1.d0+dudx) + b4*dudy ) &
297  + s23(i)*( b6*dudy + b5*dudz ) &
298  + s13(i)*( b6*(1.d0+dudx) + b4*dudz ) )
299  r_in(k2n2) = r_in(k2n2) - vol* &
300  ( s11(i)*b4*dvdx + s22(i)*b5*(1.d0+dvdy) + s33(i)*b6*dvdz &
301  + s12(i)*( b4*(1.d0+dvdy) + b5*dvdx ) &
302  + s23(i)*( b6*(1.d0+dvdy) + b5*dvdz ) &
303  + s13(i)*( b6*dvdx + b4*dvdz ) )
304  r_in(k3n2) = r_in(k3n2) - vol* &
305  ( s11(i)*b4*dwdx + s22(i)*b5*dwdy + s33(i)*b6*(1.d0+dwdz) &
306  + s12(i)*( b5*dwdx + b4*dwdy ) &
307  + s23(i)*( b6*dwdy + b5*(1.d0 + dwdz) ) &
308  + s13(i)*( b6*dwdx + b4*(1.d0 + dwdz) ) )
309 ! local node 3
310  r_in(k1n3) = r_in(k1n3) - vol* &
311  ( s11(i)*b7*(1.d0+dudx) + s22(i)*b8*dudy + s33(i)*b9*dudz &
312  + s12(i)*( b8*(1.d0+dudx) + b7*dudy ) &
313  + s23(i)*( b9*dudy + b8*dudz ) &
314  + s13(i)*( b9*(1.d0+dudx) + b7*dudz ) )
315  r_in(k2n3) = r_in(k2n3) - vol* &
316  ( s11(i)*b7*dvdx + s22(i)*b8*(1.d0+dvdy) + s33(i)*b9*dvdz &
317  + s12(i)*( b7*(1.d0+dvdy) + b8*dvdx ) &
318  + s23(i)*( b9*(1.d0+dvdy) + b8*dvdz ) &
319  + s13(i)*( b9*dvdx + b7*dvdz ) )
320  r_in(k3n3) = r_in(k3n3) - vol* &
321  ( s11(i)*b7*dwdx + s22(i)*b8*dwdy + s33(i)*b9*(1.d0+dwdz) &
322  + s12(i)*( b8*dwdx + b7*dwdy ) &
323  + s23(i)*( b9*dwdy + b8*(1.d0 + dwdz) ) &
324  + s13(i)*( b9*dwdx + b7*(1.d0 + dwdz) ) )
325 ! local node 4
326  r_in(k1n4) = r_in(k1n4) - vol* &
327  ( s11(i)*b10*(1.d0+dudx) + s22(i)*b11*dudy+s33(i)*b12*dudz &
328  + s12(i)*( b11*(1.d0+dudx) + b10*dudy ) &
329  + s23(i)*( b12*dudy + b11*dudz ) &
330  + s13(i)*( b12*(1.d0+dudx) + b10*dudz ) )
331  r_in(k2n4) = r_in(k2n4) - vol* &
332  ( s11(i)*b10*dvdx + s22(i)*b11*(1.d0+dvdy)+s33(i)*b12*dvdz &
333  + s12(i)*( b10*(1.d0+dvdy) + b11*dvdx ) &
334  + s23(i)*( b12*(1.d0+dvdy) + b11*dvdz ) &
335  + s13(i)*( b12*dvdx + b10*dvdz ) )
336  r_in(k3n4) = r_in(k3n4) - vol* &
337  ( s11(i)*b10*dwdx + s22(i)*b11*dwdy+s33(i)*b12*(1.d0+dwdz) &
338  + s12(i)*( b11*dwdx + b10*dwdy ) &
339  + s23(i)*( b12*dwdy + b11*(1.d0 + dwdz) ) &
340  + s13(i)*( b12*dwdx + b10*(1.d0 + dwdz) ) )
341 
342 
343  ENDDO
344 
345  RETURN
346 END SUBROUTINE v3d4_neohookeancompress
347 
const NT & d
blockLoc i
Definition: read.cpp:79
j indices j
Definition: Indexing.h:6
subroutine v3d4_neohookeancompress(coor, matcstet, lmcstet, R_in, d, ci, S11, S22, S33, S12, S23, S13, numnp, nstart, nend, numcstet, numat_vol, xmu, xlamda)