OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
s10derit3.F File Reference
#include "implicit_f.inc"
#include "comlock.inc"
#include "mvsiz_p.inc"
#include "com01_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine s10derit3 (vol, deltax, deltax2, xx, yy, zz, px, py, pz, nx, rx, ry, rz, sx, sy, sz, tx, ty, tz, wip, alph, beta, voln, volg, voldp, nel, offg, npt)

Function/Subroutine Documentation

◆ s10derit3()

subroutine s10derit3 ( vol,
deltax,
deltax2,
double precision, dimension(mvsiz,10) xx,
double precision, dimension(mvsiz,10) yy,
double precision, dimension(mvsiz,10) zz,
px,
py,
pz,
nx,
rx,
ry,
rz,
sx,
sy,
sz,
tx,
ty,
tz,
wip,
alph,
beta,
voln,
volg,
double precision, dimension(mvsiz,5) voldp,
integer, intent(in) nel,
intent(in) offg,
integer, intent(in) npt )

Definition at line 33 of file s10derit3.F.

41C-----------------------------------------------
42C M o d u l e s
43C-----------------------------------------------
44 USE message_mod
45 USE elbufdef_mod
46C-----------------------------------------------
47C I m p l i c i t T y p e s
48C-----------------------------------------------
49#include "implicit_f.inc"
50#include "comlock.inc"
51C-----------------------------------------------
52C G l o b a l P a r a m e t e r s
53C-----------------------------------------------
54#include "mvsiz_p.inc"
55C-----------------------------------------------
56C C o m m o n B l o c k s
57C-----------------------------------------------
58#include "com01_c.inc"
59C-----------------------------------------------
60C D u m m y A r g u m e n t s
61C-----------------------------------------------
62 INTEGER, INTENT(IN) :: NPT
63C REAL
64 INTEGER, INTENT(IN) :: NEL ! number of element in the current group
65 my_real, DIMENSION(NEL), INTENT(IN) :: offg ! off array : 0 if the element is deleted
66 double precision
67 . xx(mvsiz,10), yy(mvsiz,10), zz(mvsiz,10),voldp(mvsiz,5)
69 . vol(mvsiz,5),deltax(*),deltax2(*),
70 . nx(mvsiz,10,5),voln(*),volg(mvsiz),
71 . rx(*),ry(*),rz(*), sx(*),sy(*),sz(*), tx(*),ty(*),tz(*),
72 . px(mvsiz,10,5),py(mvsiz,10,5),pz(mvsiz,10,5),
73 . wip(5),alph(5),beta(5)
74C-----------------------------------------------
75C L o c a l V a r i a b l e s
76C-----------------------------------------------
77 INTEGER I,IP,N,K1,K2,K3,K4,K5,K6,K7,K8,K9,K10,
78 . M,IPERM(10,4),ICOR,NNEGA,INDEX(MVSIZ),J,NN
79 INTEGER ITER,NITER
80 double precision
81 . xa(mvsiz,10),ya(mvsiz,10),za(mvsiz,10),
82 . xb(mvsiz,10),yb(mvsiz,10),zb(mvsiz,10),
83 . a4, b4, a4m1, b4m1
84 DATA iperm/
85 . 2, 4, 3, 1, 9,10, 6, 5, 8, 7,
86 . 4, 1, 3, 2, 8, 7,10, 9, 5, 6,
87 . 1, 4, 2, 3, 8, 9, 5, 7,10, 6,
88 . 1, 2, 3, 4, 5, 6, 7, 8, 9,10/
89C-----------------------------------------------
90 a4 = four * alph(1)
91 b4 = four * beta(1)
92 a4m1 = a4 - one
93 b4m1 = b4 - one
94C
95 DO i=1,nel
96 rx(i) = xx(i,1) - xx(i,4)
97 ry(i) = yy(i,1) - yy(i,4)
98 rz(i) = zz(i,1) - zz(i,4)
99 sx(i) = xx(i,2) - xx(i,4)
100 sy(i) = yy(i,2) - yy(i,4)
101 sz(i) = zz(i,2) - zz(i,4)
102 tx(i) = xx(i,3) - xx(i,4)
103 ty(i) = yy(i,3) - yy(i,4)
104 tz(i) = zz(i,3) - zz(i,4)
105 volg(i) =zero
106 ENDDO
107C
108 DO n=1,4
109 DO i=1,nel
110 xa(i,n) = a4m1*xx(i,n)
111 ya(i,n) = a4m1*yy(i,n)
112 za(i,n) = a4m1*zz(i,n)
113C
114 xb(i,n) = b4m1*xx(i,n)
115 yb(i,n) = b4m1*yy(i,n)
116 zb(i,n) = b4m1*zz(i,n)
117 ENDDO
118 ENDDO
119C
120 DO n=5,10
121 DO i=1,nel
122 xa(i,n) = a4*xx(i,n)
123 ya(i,n) = a4*yy(i,n)
124 za(i,n) = a4*zz(i,n)
125C
126 xb(i,n) = b4*xx(i,n)
127 yb(i,n) = b4*yy(i,n)
128 zb(i,n) = b4*zz(i,n)
129 ENDDO
130 ENDDO
131C
132 DO ip=1,4
133 k1 = iperm(1,ip)
134 k2 = iperm(2,ip)
135 k3 = iperm(3,ip)
136 k4 = iperm(4,ip)
137 k5 = iperm(5,ip)
138 k6 = iperm(6,ip)
139 k7 = iperm(7,ip)
140 k8 = iperm(8,ip)
141 k9 = iperm(9,ip)
142 k10= iperm(10,ip)
143 CALL s10jacob(
144 1 alph(ip), beta(ip), wip(ip), xb(1,k1),
145 2 xb(1,k2), xb(1,k3), xa(1,k4), xb(1,k5),
146 3 xb(1,k6), xb(1,k7), xb(1,k8), xb(1,k9),
147 4 xb(1,k10), xa(1,k8), xa(1,k9), xa(1,k10),
148 5 yb(1,k1), yb(1,k2), yb(1,k3), ya(1,k4),
149 6 yb(1,k5), yb(1,k6), yb(1,k7), yb(1,k8),
150 7 yb(1,k9), yb(1,k10), ya(1,k8), ya(1,k9),
151 8 ya(1,k10), zb(1,k1), zb(1,k2), zb(1,k3),
152 9 za(1,k4), zb(1,k5), zb(1,k6), zb(1,k7),
153 a zb(1,k8), zb(1,k9), zb(1,k10), za(1,k8),
154 b za(1,k9), za(1,k10), px(1,k1,ip), px(1,k2,ip),
155 c px(1,k3,ip), px(1,k4,ip), px(1,k5,ip), px(1,k6,ip),
156 d px(1,k7,ip), px(1,k8,ip), px(1,k9,ip), px(1,k10,ip),
157 e py(1,k1,ip), py(1,k2,ip), py(1,k3,ip), py(1,k4,ip),
158 f py(1,k5,ip), py(1,k6,ip), py(1,k7,ip), py(1,k8,ip),
159 g py(1,k9,ip), py(1,k10,ip),pz(1,k1,ip), pz(1,k2,ip),
160 h pz(1,k3,ip), pz(1,k4,ip), pz(1,k5,ip), pz(1,k6,ip),
161 i pz(1,k7,ip), pz(1,k8,ip), pz(1,k9,ip), pz(1,k10,ip),
162 j nx(1,k1,ip), nx(1,k2,ip), nx(1,k3,ip), nx(1,k4,ip),
163 k nx(1,k5,ip), nx(1,k6,ip), nx(1,k7,ip), nx(1,k8,ip),
164 l nx(1,k9,ip), nx(1,k10,ip),vol(1,ip), voldp(1,ip),
165 m nel, offg)
166c
167 ENDDO
168C
169C
170 IF(npt==5)THEN
171 ip = 5
172C
173 DO i=1,nel
174 xa(i,1) = zero
175 ENDDO
176CC
177 CALL s10jacob(
178 1 alph(ip), beta(ip), wip(ip), xa(1,1),
179 2 xa(1,1), xa(1,1), xa(1,1), xx(1,k5),
180 3 xx(1,k6), xx(1,k7), xx(1,k8), xx(1,k9),
181 4 xx(1,k10), xx(1,k8), xx(1,k9), xx(1,k10),
182 5 xa(1,1), xa(1,1), xa(1,1), xa(1,1),
183 6 yy(1,k5), yy(1,k6), yy(1,k7), yy(1,k8),
184 7 yy(1,k9), yy(1,k10), yy(1,k8), yy(1,k9),
185 8 yy(1,k10), xa(1,1), xa(1,1), xa(1,1),
186 9 xa(1,1), zz(1,k5), zz(1,k6), zz(1,k7),
187 a zz(1,k8), zz(1,k9), zz(1,k10), zz(1,k8),
188 b zz(1,k9), zz(1,k10), px(1,k1,ip), px(1,k2,ip),
189 c px(1,k3,ip), px(1,k4,ip), px(1,k5,ip), px(1,k6,ip),
190 d px(1,k7,ip), px(1,k8,ip), px(1,k9,ip), px(1,k10,ip),
191 e py(1,k1,ip), py(1,k2,ip), py(1,k3,ip), py(1,k4,ip),
192 f py(1,k5,ip), py(1,k6,ip), py(1,k7,ip), py(1,k8,ip),
193 g py(1,k9,ip), py(1,k10,ip),pz(1,k1,ip), pz(1,k2,ip),
194 h pz(1,k3,ip), pz(1,k4,ip), pz(1,k5,ip), pz(1,k6,ip),
195 i pz(1,k7,ip), pz(1,k8,ip), pz(1,k9,ip), pz(1,k10,ip),
196 j nx(1,k1,ip), nx(1,k2,ip), nx(1,k3,ip), nx(1,k4,ip),
197 k nx(1,k5,ip), nx(1,k6,ip), nx(1,k7,ip), nx(1,k8,ip),
198 l nx(1,k9,ip), nx(1,k10,ip),vol(1,ip), voldp(1,ip),
199 m nel, offg)
200 ENDIF
201
202 DO ip=1,npt
203 DO i=1,nel
204 volg(i) =volg(i) + vol(i,ip)
205 ENDDO
206 ENDDO
207C
208 DO i=1,nel
209 voln(i) =volg(i)/npt
210 ENDDO
211C-----------
212 RETURN
#define my_real
Definition cppsort.cpp:32
subroutine s10jacob(alph, beta, w, x1b, x2b, x3b, x4a, x5b, x6b, x7b, x8b, x9b, x10b, x8a, x9a, x10a, y1b, y2b, y3b, y4a, y5b, y6b, y7b, y8b, y9b, y10b, y8a, y9a, y10a, z1b, z2b, z3b, z4a, z5b, z6b, z7b, z8b, z9b, z10b, z8a, z9a, z10a, px1, px2, px3, px4, px5, px6, px7, px8, px9, px10, py1, py2, py3, py4, py5, py6, py7, py8, py9, py10, pz1, pz2, pz3, pz4, pz5, pz6, pz7, pz8, pz9, pz10, nx1, nx2, nx3, nx4, nx5, nx6, nx7, nx8, nx9, nx10, vol, voldp)
Definition s10deri3.F:271