Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
v3d4_NeoHookeanInCompressPrin.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_neohookeanincompressprin(coor,matcstet,lmcstet,R_in,d,ci, &
54  s11,s22,s33,s12,s23,s13,numnp,nstart,nend,numcstet,numat_vol, &
55  xmu,xkappa,xlambda)
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 holding variable
93  REAL*8 :: x1d,x2d,x3d,x4d,y1d,y2d,y3d,y4d,z1d,z2d,z3d,z4d
94 !-- Coordinate subtractions
95  REAL*8 :: x14, x24, x34, y14, y24, y34, z14, z24, z34
96 !-- Coordinate subtractions: These are to speed up B calculation
97  REAL*8 :: x12, x13, y12, y13, z12, z13
98 !-- 6*volume, inverse(6*volume), and the volume
99 ! REAL*8 :: Vx6, vol
100 !-- spacial derivatives
101  REAL*8 :: b1,b2,b3,b4,b5,b6,b7,b8,b9,b10,b11,b12
102 !-- partial derivatives of the displacement
103  REAL*8 :: dudx,dvdy,dwdz,dudy,dvdx,dvdz,dwdy,dudz,dwdx
104 !-- strains
105  REAL*8 :: e11,e22,e33,e12,e23,e13
106 ! -- connectivity table for CSTet
107  INTEGER, DIMENSION(1:4,1:numcstet) :: lmcstet
108 ! -- mat number for CST element
109  INTEGER, DIMENSION(1:numcstet) :: matcstet
110 ! -- node numbers
111  INTEGER :: n1,n2,n3,n4
112 ! -- dummy and counters
113  INTEGER :: i,j,nstart,nend
114  INTEGER :: k1n1,k1n2,k1n3,k1n4,k2n1,k2n2,k2n3,k2n4
115  INTEGER :: k3n1,k3n2,k3n3,k3n4
116 ! --
117  REAL*8 :: f11, f12, f13, f21, f22, f23, f31, f32, f33
118  REAL*8 :: j1,jsqrd,jsqrdinv
119  REAL*8 :: c11, c12, c13, c21, c22, c23, c31, c32, c33
120  REAL*8 :: xmu, xlambda
121  REAL*8, DIMENSION(1:4,1:4) :: eledb
122  REAL*8,DIMENSION(1:3,1:3) :: xj
123  REAL*8,DIMENSION(1:3,1:3) :: xjv0
124  REAL*8 :: weight = 1.d0/6.d0
125  INTEGER :: id, jd, in, ip
126  REAL*8 :: xcoor
127  REAL*8 :: detjb,detjbv0,evol,theta,xkappa,press
128  REAL*8 :: ic, iiic,vx6old0,vx6old
129  REAL*8 :: onethird = 1.d0/3.d0
130  REAL*8 :: val11, val21, val31
131 
132  REAL*8 :: v0x6, v0x6inv
133  REAL*8 :: vx6, vx6inv
134  REAL*8 :: vol, vol0
135  REAL*8 :: term1, term2, cinv
136  REAL*8, DIMENSION(1:3,1:3) :: btens,princ,sigma
137  REAL*8, DIMENSION(1:3) :: stret,sprin
138  REAL*8 :: detf
139  REAL*8 :: sh1, sh2, sh3, sh4, sh5, sh6, sh7, sh8, sh9, sh10, sh11, sh12
140  REAL*8 :: voldef,vx6invdef,vx6def
141 
142  DO i = nstart, nend
143  j = matcstet(i)
144 
145  n1 = lmcstet(1,i)
146  n2 = lmcstet(2,i)
147  n3 = lmcstet(3,i)
148  n4 = lmcstet(4,i)
149 
150  k3n1 = 3*n1
151  k3n2 = 3*n2
152  k3n3 = 3*n3
153  k3n4 = 3*n4
154 
155  k2n1 = k3n1 - 1
156  k2n2 = k3n2 - 1
157  k2n3 = k3n3 - 1
158  k2n4 = k3n4 - 1
159 
160  k1n1 = k3n1 - 2
161  k1n2 = k3n2 - 2
162  k1n3 = k3n3 - 2
163  k1n4 = k3n4 - 2
164  ! k#n# dummy variables replaces:
165  u1 = d(k1n1) ! (3*n1-2)
166  u2 = d(k1n2) ! (3*n2-2)
167  u3 = d(k1n3) ! (3*n3-2)
168  u4 = d(k1n4) ! (3*n4-2)
169  v1 = d(k2n1) ! (3*n1-1)
170  v2 = d(k2n2) ! (3*n2-1)
171  v3 = d(k2n3) ! (3*n3-1)
172  v4 = d(k2n4) ! (3*n4-1)
173  w1 = d(k3n1) ! (3*n1)
174  w2 = d(k3n2) ! (3*n2)
175  w3 = d(k3n3) ! (3*n3)
176  w4 = d(k3n4) ! (3*n4)
177 
178  x1 = coor(1,n1)
179  x2 = coor(1,n2)
180  x3 = coor(1,n3)
181  x4 = coor(1,n4)
182  y1 = coor(2,n1)
183  y2 = coor(2,n2)
184  y3 = coor(2,n3)
185  y4 = coor(2,n4)
186  z1 = coor(3,n1)
187  z2 = coor(3,n2)
188  z3 = coor(3,n3)
189  z4 = coor(3,n4)
190 
191  x12 = x1 - x2 ! not used in vol. calc
192  x13 = x1 - x3 ! not used in vol. calc
193  x14 = x1 - x4
194  x24 = x2 - x4
195  x34 = x3 - x4
196  y12 = y1 - y2 ! not used in vol. calc
197  y13 = y1 - y3 ! not used in vol. calc
198  y14 = y1 - y4
199  y24 = y2 - y4
200  y34 = y3 - y4
201  z12 = z1 - z2 ! not used in vol. calc
202  z13 = z1 - z3 ! not used in vol. calc
203  z14 = z1 - z4
204  z24 = z2 - z4
205  z34 = z3 - z4
206 
207  val11 = y24*z34 - z24*y34
208  val21 = -( x24*z34 - z24*x34 )
209  val31 = x24*y34 - y24*x34
210 
211  v0x6 = -( x14*val11 + y14*val21 + z14*val31 )
212 
213  vx6inv = 1.d0/v0x6
214 
215 ! See the maple worksheet 'V3D4.mws' for the derivation of Vx6
216 
217 ! See the maple worksheet 'V3D4.mws' for the derivation of [B]
218 ! NOTE: Factored for a more equivalent/compact form then maple's
219 
220  b1 = (y34*z24 - y24*z34) * vx6inv
221  b2 = (z34*x24 - z24*x34) * vx6inv
222  b3 = (x34*y24 - x24*y34) * vx6inv
223  b4 = (y13*z14 - y14*z13) * vx6inv
224  b5 = (z13*x14 - z14*x13) * vx6inv
225  b6 = (x13*y14 - x14*y13) * vx6inv
226  b7 = (y14*z12 - y12*z14) * vx6inv
227  b8 = (z14*x12 - z12*x14) * vx6inv
228  b9 = (x14*y12 - x12*y14) * vx6inv
229  b10 = (y12*z13 - y13*z12) * vx6inv
230  b11 = (z12*x13 - z13*x12) * vx6inv
231  b12 = (x12*y13 - x13*y12) * vx6inv
232 !
233 ! deformation gradients F
234 !
235  f11 = 1.d0 + ( b1*u1 + b4*u2 + b7*u3 + b10*u4 ) ! 1 + ( dudx )
236  f22 = 1.d0 + ( b2*v1 + b5*v2 + b8*v3 + b11*v4 ) ! 1 + ( dvdy )
237  f33 = 1.d0 + ( b3*w1 + b6*w2 + b9*w3 + b12*w4 ) ! 1 + ( dwdz )
238  f12 = b2*u1 + b5*u2 + b8*u3 + b11*u4 ! dudy
239  f21 = b1*v1 + b4*v2 + b7*v3 + b10*v4 ! dvdx
240  f23 = b3*v1 + b6*v2 + b9*v3 + b12*v4 ! dvdz
241  f32 = b2*w1 + b5*w2 + b8*w3 + b11*w4 ! dwdy
242  f13 = b3*u1 + b6*u2 + b9*u3 + b12*u4 ! dudz
243  f31 = b1*w1 + b4*w2 + b7*w3 + b10*w4 ! dwdx
244 
245  detf = f11*(f22*f33-f23*f32)+f12*(f31*f23-f21*f33)+f13*(-f31*f22+f21*f32)
246  IF(detf.LE.0.d0) THEN
247  WRITE(*,*)'Jacobian has become zero for element',i
248  stop
249  ENDIF
250 
251 !
252 ! obtains the tensor b, left cauchy-green tensor using equation ???
253 !
254  btens(1,1) = f11**2+f12**2+f13**2
255  btens(1,2) = f21*f11+f12*f22+f13*f23
256  btens(1,3) = f31*f11+f12*f32+f13*f33
257  btens(2,1) = btens(1,2)
258  btens(2,2) = f21**2+f22**2+f23**2
259  btens(2,3) = f21*f31+f32*f22+f23*f33
260  btens(3,1) = btens(1,3)
261  btens(3,2) = btens(2,3)
262  btens(3,3) = f31**2+f32**2+f33**2
263 
264  CALL jacobi(btens,stret,princ)
265 
266  CALL cauchystressprinc(3,xmu,xlambda,detf,stret,princ,sigma,sprin)
267 
268  s11(i) = sigma(1,1)
269  s12(i) = sigma(1,2)
270  s13(i) = sigma(1,3)
271  s22(i) = sigma(2,2)
272  s23(i) = sigma(2,3)
273  s33(i) = sigma(3,3)
274 
275  x1 = x1 + u1
276  x2 = x2 + u2
277  x3 = x3 + u3
278  x4 = x4 + u4
279  y1 = y1 + v1
280  y2 = y2 + v2
281  y3 = y3 + v3
282  y4 = y4 + v4
283  z1 = z1 + w1
284  z2 = z2 + w2
285  z3 = z3 + w3
286  z4 = z4 + w4
287 
288  x12 = x1 - x2 ! not used in vol. calc
289  x13 = x1 - x3 ! not used in vol. calc
290  x14 = x1 - x4
291  x24 = x2 - x4
292  x34 = x3 - x4
293  y12 = y1 - y2 ! not used in vol. calc
294  y13 = y1 - y3 ! not used in vol. calc
295  y14 = y1 - y4
296  y24 = y2 - y4
297  y34 = y3 - y4
298  z12 = z1 - z2 ! not used in vol. calc
299  z13 = z1 - z3 ! not used in vol. calc
300  z14 = z1 - z4
301  z24 = z2 - z4
302  z34 = z3 - z4
303 
304  val11 = y24*z34 - z24*y34
305  val21 = -( x24*z34 - z24*x34 )
306  val31 = x24*y34 - y24*x34
307 
308  vx6def = -( x14*val11 + y14*val21 + z14*val31 )
309 
310  vx6invdef = 1.d0 / vx6def
311 
312 ! calculate the volume
313  voldef = vx6def/6.d0
314 
315  IF(vx6def.LE.0.d0) THEN
316  WRITE(*,100) i
317  stop
318  ENDIF
319 
320  sh1 = (y34*z24 - y24*z34) * vx6invdef
321  sh2 = (z34*x24 - z24*x34) * vx6invdef
322  sh3 = (x34*y24 - x24*y34) * vx6invdef
323  sh4 = (y13*z14 - y14*z13) * vx6invdef
324  sh5 = (z13*x14 - z14*x13) * vx6invdef
325  sh6 = (x13*y14 - x14*y13) * vx6invdef
326  sh7 = (y14*z12 - y12*z14) * vx6invdef
327  sh8 = (z14*x12 - z12*x14) * vx6invdef
328  sh9 = (x14*y12 - x12*y14) * vx6invdef
329  sh10 = (y12*z13 - y13*z12) * vx6invdef
330  sh11 = (z12*x13 - z13*x12) * vx6invdef
331  sh12 = (x12*y13 - x13*y12) * vx6invdef
332 
333 ! ASSEMBLE THE INTERNAL FORCE VECTOR (uses Cauchy Stress)
334 !
335 ! local node 1
336  r_in(k1n1) = r_in(k1n1) - voldef* &
337  ( s11(i)*sh1 + s12(i)*sh2 + s13(i)*sh3 )
338  r_in(k2n1) = r_in(k2n1) - voldef* &
339  ( s12(i)*sh1 + s22(i)*sh2 + s23(i)*sh3 )
340  r_in(k3n1) = r_in(k3n1) - voldef* &
341  ( s13(i)*sh1 + s23(i)*sh2 + s33(i)*sh3 )
342 ! local node 2
343  r_in(k1n2) = r_in(k1n2) - voldef* &
344  ( s11(i)*sh4 + s12(i)*sh5 + s13(i)*sh6 )
345  r_in(k2n2) = r_in(k2n2) - voldef* &
346  ( s12(i)*sh4 + s22(i)*sh5 + s23(i)*sh6 )
347  r_in(k3n2) = r_in(k3n2) - voldef* &
348  ( s13(i)*sh4 + s23(i)*sh5 + s33(i)*sh6 )
349 ! local node 3
350  r_in(k1n3) = r_in(k1n3) - voldef* &
351  ( s11(i)*sh7 + s12(i)*sh8 + s13(i)*sh9 )
352  r_in(k2n3) = r_in(k2n3) - voldef* &
353  ( s12(i)*sh7 + s22(i)*sh8 + s23(i)*sh9 )
354  r_in(k3n3) = r_in(k3n3) - voldef* &
355  ( s13(i)*sh7 + s23(i)*sh8 + s33(i)*sh9 )
356 ! local node 4
357  r_in(k1n4) = r_in(k1n4) - voldef* &
358  ( s11(i)*sh10 + s12(i)*sh11 + s13(i)*sh12 )
359  r_in(k2n4) = r_in(k2n4) - voldef* &
360  ( s12(i)*sh10 + s22(i)*sh11 + s23(i)*sh12 )
361  r_in(k3n4) = r_in(k3n4) - voldef* &
362  ( s13(i)*sh10 + s23(i)*sh11 + s33(i)*sh12 )
363 
364  ENDDO
365 
366  RETURN
367 
368 100 FORMAT(' Negative Jacobian for element: ',i10)
369 
370 END SUBROUTINE v3d4_neohookeanincompressprin
371 
subroutine cauchystressprinc(ndime, xmu, xlamb, detf, stret, princ, sigma, sprin)
const NT & d
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE ic
blockLoc i
Definition: read.cpp:79
subroutine jacobi(btens, stret, princ)
Definition: jacobi.f90:53
subroutine v3d4_neohookeanincompressprin(coor, matcstet, lmcstet, R_in, d, ci, S11, S22, S33, S12, S23, S13, numnp, nstart, nend, numcstet, numat_vol, xmu, xkappa, xlambda)
j indices j
Definition: Indexing.h:6
unsigned long id(const Leda_like_handle &x)
Definition: Handle.h:107