Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
v3d4.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(coor,matcstet,lmcstet,R_in,d,ci, &
54  s11,s22,s33,s12,s23,s13,numnp,nstart,nend,numcstet,numat_vol)
55 !________________________________________________________________________
56 !
57 ! V3D4 - Performs displacement based computations for Volumetric 3D
58 ! 4-node tetrahedra linear elastic elements with linear
59 ! interpolation functions. (constant strain tetrahedra).
60 ! Returns the internal force vector R_in.
61 !
62 ! DATE: 04.2000 AUTHOR: SCOT BREITENFELD
63 !________________________________________________________________________
64 
65  IMPLICIT NONE
66 !---- Global variables
67  INTEGER :: numnp ! number of nodal points
68  INTEGER :: numcstet ! number of CSTet elements
69  INTEGER :: numat_vol ! number of volumetric materials
70 !-- coordinate array
71  REAL*8, DIMENSION(1:3,1:numnp) :: coor
72 !-- elastic stiffness consts
73  REAL*8, DIMENSION(1:9,1:numat_vol) :: ci
74 !-- internal force
75  REAL*8, DIMENSION(1:3*numnp) :: r_in
76 !-- displacement vector
77  REAL*8, DIMENSION(1:3*numnp) :: d
78 !-- CSTet stress
79  REAL*8, DIMENSION(1:numcstet) :: s11, s22, s33, s12, s23, s13
80 !-- connectivity table for CSTet
81  INTEGER, DIMENSION(1:4,1:numcstet) :: lmcstet
82 !-- mat number for CSTet element
83  INTEGER, DIMENSION(1:numcstet) :: matcstet
84 !---- Local variables
85 !-- node numbers
86  INTEGER :: n1,n2,n3,n4
87 !-- x, y, and z displacements of nodes
88  REAL*8 :: u1,u2,u3,u4,v1,v2,v3,v4,w1,w2,w3,w4
89 !-- 6*volume, inverse(6*volume), and the volume
90  REAL*8 :: vx6, vx6inv, vol
91 !-- spacial derivatives
92  REAL*8 :: b1,b2,b3,b4,b5,b6,b7,b8,b9,b10,b11,b12
93 !-- strains
94  REAL*8 :: e11,e22,e33,e12,e23,e13
95 !-- coordinate holding variable
96  REAL*8 :: x1,x2,x3,x4,y1,y2,y3,y4,z1,z2,z3,z4
97 !-- Coordinate subtractions
98  REAL*8 :: x14, x24, x34, y14, y24, y34, z14, z24, z34
99 !-- Coordinate subtractions: These are to speed up B calculation
100  REAL*8 :: x12, x13, y12, y13, z12, z13
101 !-- Dummy
102  REAL*8 :: c11, c21, c31
103 !-- dummy and counters
104  INTEGER :: i,j,nstart,nend
105  INTEGER :: k1n1,k1n2,k1n3,k1n4,k2n1,k2n2,k2n3,k2n4
106  INTEGER :: k3n1,k3n2,k3n3,k3n4
107 
108  DO i = nstart, nend
109  j = matcstet(i)
110 
111  n1 = lmcstet(1,i)
112  n2 = lmcstet(2,i)
113  n3 = lmcstet(3,i)
114  n4 = lmcstet(4,i)
115 
116  k3n1 = 3*n1
117  k3n2 = 3*n2
118  k3n3 = 3*n3
119  k3n4 = 3*n4
120 
121  k2n1 = k3n1 - 1
122  k2n2 = k3n2 - 1
123  k2n3 = k3n3 - 1
124  k2n4 = k3n4 - 1
125 
126  k1n1 = k3n1 - 2
127  k1n2 = k3n2 - 2
128  k1n3 = k3n3 - 2
129  k1n4 = k3n4 - 2
130  ! k#n# dummy variables replaces:
131  u1 = d(k1n1) ! (3*n1-2)
132  u2 = d(k1n2) ! (3*n2-2)
133  u3 = d(k1n3) ! (3*n3-2)
134  u4 = d(k1n4) ! (3*n4-2)
135  v1 = d(k2n1) ! (3*n1-1)
136  v2 = d(k2n2) ! (3*n2-1)
137  v3 = d(k2n3) ! (3*n3-1)
138  v4 = d(k2n4) ! (3*n4-1)
139  w1 = d(k3n1) ! (3*n1)
140  w2 = d(k3n2) ! (3*n2)
141  w3 = d(k3n3) ! (3*n3)
142  w4 = d(k3n4) ! (3*n4)
143 
144  x1 = coor(1,n1) ! Node 1, x-coor
145  x2 = coor(1,n2) ! Node 2, x-coor
146  x3 = coor(1,n3) ! Node 3, x-coor
147  x4 = coor(1,n4) ! Node 4, x-coor
148  y1 = coor(2,n1) ! Node 1, y-coor
149  y2 = coor(2,n2) ! Node 2, y-coor
150  y3 = coor(2,n3) ! Node 3, y-coor
151  y4 = coor(2,n4) ! Node 4, y-coor
152  z1 = coor(3,n1) ! Node 1, z-coor
153  z2 = coor(3,n2) ! Node 2, z-coor
154  z3 = coor(3,n3) ! Node 3, z-coor
155  z4 = coor(3,n4) ! Node 4, z-coor
156 
157  x12 = x1 - x2 ! not used in vol. calc
158  x13 = x1 - x3 ! not used in vol. calc
159  x14 = x1 - x4
160  x24 = x2 - x4
161  x34 = x3 - x4
162  y12 = y1 - y2 ! not used in vol. calc
163  y13 = y1 - y3 ! not used in vol. calc
164  y14 = y1 - y4
165  y24 = y2 - y4
166  y34 = y3 - y4
167  z12 = z1 - z2 ! not used in vol. calc
168  z13 = z1 - z3 ! not used in vol. calc
169  z14 = z1 - z4
170  z24 = z2 - z4
171  z34 = z3 - z4
172 
173  c11 = y24*z34 - z24*y34
174  c21 = -( x24*z34 - z24*x34 )
175  c31 = x24*y34 - y24*x34
176 
177  vx6 = -( x14*c11 + y14*c21 + z14*c31 )
178 
179  vx6inv = 1.d0 / vx6
180 
181 ! See the maple worksheet 'V3D4.mws' for the derivation of [B]
182 ! NOTE: Factored for a more equivalent/compact form then maple's
183 
184  b1 = (y34*z24 - y24*z34) * vx6inv
185  b2 = (z34*x24 - z24*x34) * vx6inv
186  b3 = (x34*y24 - x24*y34) * vx6inv
187  b4 = (y13*z14 - y14*z13) * vx6inv
188  b5 = (z13*x14 - z14*x13) * vx6inv
189  b6 = (x13*y14 - x14*y13) * vx6inv
190  b7 = (y14*z12 - y12*z14) * vx6inv
191  b8 = (z14*x12 - z12*x14) * vx6inv
192  b9 = (x14*y12 - x12*y14) * vx6inv
193  b10 = (y12*z13 - y13*z12) * vx6inv
194  b11 = (z12*x13 - z13*x12) * vx6inv
195  b12 = (x12*y13 - x13*y12) * vx6inv
196 
197 ! calculate the strain
198 ! [E] = [B]{d}
199 !
200  e11 = b1*u1 + b4*u2 + b7*u3 + b10*u4
201  e22 = b2*v1 + b5*v2 + b8*v3 + b11*v4
202  e33 = b3*w1 + b6*w2 + b9*w3 + b12*w4
203  e12 = b2*u1 + b1*v1 + b5*u2 + b4*v2 + b8*u3 + b7*v3 + b11*u4 + b10*v4
204  e23 = b3*v1 + b2*w1 + b6*v2 + b5*w2 + b9*v3 + b8*w3 + b12*v4 + b11*w4
205  e13 = b3*u1 + b1*w1 + b6*u2 + b4*w2 + b9*u3 + b7*w3 + b12*u4 + b10*w4
206 
207 ! calculate the stress -1
208 ! [S] = [C] {E}
209 !
210  s11(i) = e11*ci(1,j) + e22*ci(2,j) + e33*ci(4,j)
211  s22(i) = e11*ci(2,j) + e22*ci(3,j) + e33*ci(5,j)
212  s33(i) = e11*ci(4,j) + e22*ci(5,j) + e33*ci(6,j)
213  s12(i) = e12*ci(7,j)
214  s23(i) = e23*ci(8,j)
215  s13(i) = e13*ci(9,j)
216 
217 ! calculate the volume
218 
219  vol = vx6 / 6.d0
220 
221 ! ASSEMBLE THE INTERNAL FORCE VECTOR
222 
223 ! local node 1
224  r_in(k1n1) = r_in(k1n1) - (s11(i)*b1 + s12(i)*b2 + s13(i)*b3)*vol
225  r_in(k2n1) = r_in(k2n1) - (s22(i)*b2 + s12(i)*b1 + s23(i)*b3)*vol
226  r_in(k3n1) = r_in(k3n1) - (s33(i)*b3 + s23(i)*b2 + s13(i)*b1)*vol
227 ! local node 2
228  r_in(k1n2) = r_in(k1n2) - (s11(i)*b4 + s12(i)*b5 + s13(i)*b6)*vol
229  r_in(k2n2) = r_in(k2n2) - (s22(i)*b5 + s12(i)*b4 + s23(i)*b6)*vol
230  r_in(k3n2) = r_in(k3n2) - (s33(i)*b6 + s23(i)*b5 + s13(i)*b4)*vol
231 ! local node 3
232  r_in(k1n3) = r_in(k1n3) - (s11(i)*b7 + s12(i)*b8 + s13(i)*b9 )*vol
233  r_in(k2n3) = r_in(k2n3) - (s22(i)*b8 + s12(i)*b7 + s23(i)*b9 )*vol
234  r_in(k3n3) = r_in(k3n3) - (s33(i)*b9 + s23(i)*b8 + s13(i)*b7 )*vol
235 ! local node 4
236  r_in(k1n4) = r_in(k1n4) - (s11(i)*b10 + s12(i)*b11 + s13(i)*b12 )*vol
237  r_in(k2n4) = r_in(k2n4) - (s22(i)*b11 + s12(i)*b10 + s23(i)*b12 )*vol
238  r_in(k3n4) = r_in(k3n4) - (s33(i)*b12 + s23(i)*b11 + s13(i)*b10 )*vol
239 
240  x1 = coor(1,n1) + u1
241  x2 = coor(1,n2) + u2
242  x3 = coor(1,n3) + u3
243  x4 = coor(1,n4) + u4
244  y1 = coor(2,n1) + v1
245  y2 = coor(2,n2) + v2
246  y3 = coor(2,n3) + v3
247  y4 = coor(2,n4) + v4
248  z1 = coor(3,n1) + w1
249  z2 = coor(3,n2) + w2
250  z3 = coor(3,n3) + w3
251  z4 = coor(3,n4) + w4
252 
253  x14 = x1 - x4
254  x24 = x2 - x4
255  x34 = x3 - x4
256  y14 = y1 - y4
257  y24 = y2 - y4
258  y34 = y3 - y4
259  z14 = z1 - z4
260  z24 = z2 - z4
261  z34 = z3 - z4
262 
263  c11 = y24*z34 - z24*y34
264  c21 = -( x24*z34 - z24*x34 )
265  c31 = x24*y34 - y24*x34
266 
267  vx6 = -( x14*c11 + y14*c21 + z14*c31 )
268 
269  IF(vx6.LT.0.d0) THEN
270  print*,'DEFORMED VOLUME TURNED NEGATIVE'
271  stop
272  ENDIF
273 
274  ENDDO
275  RETURN
276 END SUBROUTINE v3d4
277 
const NT & d
blockLoc i
Definition: read.cpp:79
subroutine v3d4(coor, matcstet, lmcstet, R_in, d, ci, S11, S22, S33, S12, S23, S13, numnp, nstart, nend, numcstet, numat_vol)
Definition: v3d4.f90:53
virtual std::ostream & print(std::ostream &os) const
j indices j
Definition: Indexing.h:6