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