OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
s10deri3.F File Reference
#include "implicit_f.inc"
#include "comlock.inc"
#include "mvsiz_p.inc"
#include "impl1_c.inc"
#include "scr17_c.inc"
#include "scr07_c.inc"
#include "scr18_c.inc"
#include "lockon.inc"
#include "lockoff.inc"
#include "vectorize.inc"

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ s10deri3()

subroutine s10deri3 ( off,
vol,
integer, dimension(*) ngl,
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, dimension(mvsiz,10) nc,
double precision, dimension(nel,30) sav,
offg,
integer nel,
integer, intent(in) npt,
integer, intent(in) ismstr,
integer, intent(in) jlag )

Definition at line 35 of file s10deri3.F.

45C-----------------------------------------------
46C M o d u l e s
47C-----------------------------------------------
48 USE message_mod
49C-----------------------------------------------
50C I m p l i c i t T y p e s
51C-----------------------------------------------
52#include "implicit_f.inc"
53#include "comlock.inc"
54C-----------------------------------------------
55C G l o b a l P a r a m e t e r s
56C-----------------------------------------------
57#include "mvsiz_p.inc"
58C-----------------------------------------------
59C C o m m o n B l o c k s
60C-----------------------------------------------
61#include "impl1_c.inc"
62#include "scr17_c.inc"
63#include "scr07_c.inc"
64#include "scr18_c.inc"
65C-----------------------------------------------
66C D u m m y A r g u m e n t s
67C-----------------------------------------------
68 INTEGER, INTENT(IN) :: NPT
69 INTEGER, INTENT(IN) :: ISMSTR
70 INTEGER, INTENT(IN) :: JLAG
71 INTEGER NGL(*), NC(MVSIZ,10), NEL
72C REAL
73 double precision
74 . xx(mvsiz,10), yy(mvsiz,10), zz(mvsiz,10),sav(nel,30),voldp(mvsiz,5)
76 . off(nel),vol(mvsiz,5),deltax(*),deltax2(*),
77 . rx(*),ry(*),rz(*), sx(*),sy(*),sz(*), tx(*),ty(*),tz(*),
78 . nx(mvsiz,10,5),voln(*),volg(mvsiz),
79 . px(mvsiz,10,5),py(mvsiz,10,5),pz(mvsiz,10,5),
80 . wip(5),alph(5),beta(5),offg(nel)
81C-----------------------------------------------
82C L o c a l V a r i a b l e s
83C-----------------------------------------------
84 INTEGER I,IP,N,K1,K2,K3,K4,K5,K6,K7,K8,K9,K10,
85 . IPERM(10,4),ICOR,NNEGA,INDEX(MVSIZ),J,NN
86 double precision
87 . xa(mvsiz,10),ya(mvsiz,10),za(mvsiz,10),
88 . xb(mvsiz,10),yb(mvsiz,10),zb(mvsiz,10),
89 . a4, b4, a4m1, b4m1
90 DATA iperm/
91 . 2, 4, 3, 1, 9,10, 6, 5, 8, 7,
92 . 4, 1, 3, 2, 8, 7,10, 9, 5, 6,
93 . 1, 4, 2, 3, 8, 9, 5, 7,10, 6,
94 . 1, 2, 3, 4, 5, 6, 7, 8, 9,10/
95C-----------------------------------------------
96 a4 = four * alph(1)
97 b4 = four * beta(1)
98 a4m1 = a4 - one
99 b4m1 = b4 - one
100C
101 DO i=1,nel
102 rx(i) = xx(i,1) - xx(i,4)
103 ry(i) = yy(i,1) - yy(i,4)
104 rz(i) = zz(i,1) - zz(i,4)
105 sx(i) = xx(i,2) - xx(i,4)
106 sy(i) = yy(i,2) - yy(i,4)
107 sz(i) = zz(i,2) - zz(i,4)
108 tx(i) = xx(i,3) - xx(i,4)
109 ty(i) = yy(i,3) - yy(i,4)
110 tz(i) = zz(i,3) - zz(i,4)
111 volg(i) =zero
112 ENDDO
113C
114 DO n=1,4
115 DO i=1,nel
116 xa(i,n) = a4m1*xx(i,n)
117 ya(i,n) = a4m1*yy(i,n)
118 za(i,n) = a4m1*zz(i,n)
119C
120 xb(i,n) = b4m1*xx(i,n)
121 yb(i,n) = b4m1*yy(i,n)
122 zb(i,n) = b4m1*zz(i,n)
123 ENDDO
124 ENDDO
125C
126 DO n=5,10
127 DO i=1,nel
128 xa(i,n) = a4*xx(i,n)
129 ya(i,n) = a4*yy(i,n)
130 za(i,n) = a4*zz(i,n)
131C
132 xb(i,n) = b4*xx(i,n)
133 yb(i,n) = b4*yy(i,n)
134 zb(i,n) = b4*zz(i,n)
135 ENDDO
136 ENDDO
137C
138 DO ip=1,4
139 k1 = iperm(1,ip)
140 k2 = iperm(2,ip)
141 k3 = iperm(3,ip)
142 k4 = iperm(4,ip)
143 k5 = iperm(5,ip)
144 k6 = iperm(6,ip)
145 k7 = iperm(7,ip)
146 k8 = iperm(8,ip)
147 k9 = iperm(9,ip)
148 k10= iperm(10,ip)
149 CALL s10jacob(
150 1 alph(ip), beta(ip), wip(ip), xb(1,k1),
151 2 xb(1,k2), xb(1,k3), xa(1,k4), xb(1,k5),
152 3 xb(1,k6), xb(1,k7), xb(1,k8), xb(1,k9),
153 4 xb(1,k10), xa(1,k8), xa(1,k9), xa(1,k10),
154 5 yb(1,k1), yb(1,k2), yb(1,k3), ya(1,k4),
155 6 yb(1,k5), yb(1,k6), yb(1,k7), yb(1,k8),
156 7 yb(1,k9), yb(1,k10), ya(1,k8), ya(1,k9),
157 8 ya(1,k10), zb(1,k1), zb(1,k2), zb(1,k3),
158 9 za(1,k4), zb(1,k5), zb(1,k6), zb(1,k7),
159 a zb(1,k8), zb(1,k9), zb(1,k10), za(1,k8),
160 b za(1,k9), za(1,k10), px(1,k1,ip), px(1,k2,ip),
161 c px(1,k3,ip), px(1,k4,ip), px(1,k5,ip), px(1,k6,ip),
162 d px(1,k7,ip), px(1,k8,ip), px(1,k9,ip), px(1,k10,ip),
163 e py(1,k1,ip), py(1,k2,ip), py(1,k3,ip), py(1,k4,ip),
164 f py(1,k5,ip), py(1,k6,ip), py(1,k7,ip), py(1,k8,ip),
165 g py(1,k9,ip), py(1,k10,ip),pz(1,k1,ip), pz(1,k2,ip),
166 h pz(1,k3,ip), pz(1,k4,ip), pz(1,k5,ip), pz(1,k6,ip),
167 i pz(1,k7,ip), pz(1,k8,ip), pz(1,k9,ip), pz(1,k10,ip),
168 j nx(1,k1,ip), nx(1,k2,ip), nx(1,k3,ip), nx(1,k4,ip),
169 k nx(1,k5,ip), nx(1,k6,ip), nx(1,k7,ip), nx(1,k8,ip),
170 l nx(1,k9,ip), nx(1,k10,ip),vol(1,ip), voldp(1,ip),
171 m nel, offg)
172c
173 ENDDO
174C
175C
176 IF(npt==5)THEN
177 ip = 5
178C
179 DO i=1,nel
180 xa(i,1) = zero
181 ENDDO
182C
183 CALL s10jacob(
184 1 alph(ip), beta(ip), wip(ip), xa(1,1),
185 2 xa(1,1), xa(1,1), xa(1,1), xx(1,k5),
186 3 xx(1,k6), xx(1,k7), xx(1,k8), xx(1,k9),
187 4 xx(1,k10), xx(1,k8), xx(1,k9), xx(1,k10),
188 5 xa(1,1), xa(1,1), xa(1,1), xa(1,1),
189 6 yy(1,k5), yy(1,k6), yy(1,k7), yy(1,k8),
190 7 yy(1,k9), yy(1,k10), yy(1,k8), yy(1,k9),
191 8 yy(1,k10), xa(1,1), xa(1,1), xa(1,1),
192 9 xa(1,1), zz(1,k5), zz(1,k6), zz(1,k7),
193 a zz(1,k8), zz(1,k9), zz(1,k10), zz(1,k8),
194 b zz(1,k9), zz(1,k10), px(1,k1,ip), px(1,k2,ip),
195 c px(1,k3,ip), px(1,k4,ip), px(1,k5,ip), px(1,k6,ip),
196 d px(1,k7,ip), px(1,k8,ip), px(1,k9,ip), px(1,k10,ip),
197 e py(1,k1,ip), py(1,k2,ip), py(1,k3,ip), py(1,k4,ip),
198 f py(1,k5,ip), py(1,k6,ip), py(1,k7,ip), py(1,k8,ip),
199 g py(1,k9,ip), py(1,k10,ip),pz(1,k1,ip), pz(1,k2,ip),
200 h pz(1,k3,ip), pz(1,k4,ip), pz(1,k5,ip), pz(1,k6,ip),
201 i pz(1,k7,ip), pz(1,k8,ip), pz(1,k9,ip), pz(1,k10,ip),
202 j nx(1,k1,ip), nx(1,k2,ip), nx(1,k3,ip), nx(1,k4,ip),
203 k nx(1,k5,ip), nx(1,k6,ip), nx(1,k7,ip), nx(1,k8,ip),
204 l nx(1,k9,ip), nx(1,k10,ip),vol(1,ip), voldp(1,ip),
205 m nel, offg)
206 ENDIF
207C
208C
209 nnega = 0
210 IF(jlag/=0.AND.(ismstr==10.OR.(ismstr==12.AND.idtmin(1)/=3))) THEN
211 DO i=1,nel
212 IF(offg(i) > one)THEN
213 nnega=nnega+1
214 index(nnega)=i
215 END IF
216 ENDDO
217 END IF
218 icor=0
219 DO ip=1,npt
220 DO i=1,nel
221 IF(off(i)==0.)THEN
222 vol(i,ip)=one
223 ELSEIF(off(i)> one)THEN
224 ELSEIF(vol(i,ip)<=zero)THEN
225 icor=1
226 ENDIF
227 ENDDO
228 ENDDO
229C-------------note: /STOP /DEL are not implemented
230 IF(jlag/=0)THEN
231 IF(icor/=0.AND.inconv==1)THEN
232 DO ip=1,npt
233 DO i=1,nel
234 IF(off(i) == zero.OR.offg(i) > one)THEN
235 ELSEIF(vol(i,ip)<=zero)THEN
236 nnega=nnega+1
237 index(nnega)=i
238#include "lockon.inc"
239 IF(ismstr<10) THEN
240 CALL ancmsg(msgid=260,anmode=aninfo,
241 . i1=ngl(i))
242 ELSE
243 CALL ancmsg(msgid=262,anmode=aninfo,
244 . i1=ngl(i))
245 END IF
246#include "lockoff.inc"
247 offg(i) = two
248 ENDIF
249 ENDDO
250 ENDDO
251 END if!(ICOR/=0.AND.INCONV==1)THEN
252 ENDIF !(JLAG/=0)
253C
254 IF(nnega >0 )THEN
255 DO n=1,10
256#include "vectorize.inc"
257 DO j=1,nnega
258 i = index(j)
259 nn = nc(i,n)
260 xx(i,n)=sav(i,n)
261 yy(i,n)=sav(i,n+10)
262 zz(i,n)=sav(i,n+20)
263 ENDDO
264 END DO
265#include "vectorize.inc"
266 DO j=1,nnega
267 i = index(j)
268 rx(i) = xx(i,1) - xx(i,4)
269 ry(i) = yy(i,1) - yy(i,4)
270 rz(i) = zz(i,1) - zz(i,4)
271 sx(i) = xx(i,2) - xx(i,4)
272 sy(i) = yy(i,2) - yy(i,4)
273 sz(i) = zz(i,2) - zz(i,4)
274 tx(i) = xx(i,3) - xx(i,4)
275 ty(i) = yy(i,3) - yy(i,4)
276 tz(i) = zz(i,3) - zz(i,4)
277 ENDDO
278
279 DO n=1,4
280#include "vectorize.inc"
281 DO j=1,nnega
282 i = index(j)
283 xa(i,n) = a4m1*xx(i,n)
284 ya(i,n) = a4m1*yy(i,n)
285 za(i,n) = a4m1*zz(i,n)
286
287 xb(i,n) = b4m1*xx(i,n)
288 yb(i,n) = b4m1*yy(i,n)
289 zb(i,n) = b4m1*zz(i,n)
290 ENDDO
291 ENDDO
292
293 DO n=5,10
294#include "vectorize.inc"
295 DO j=1,nnega
296 i = index(j)
297 xa(i,n) = a4*xx(i,n)
298 ya(i,n) = a4*yy(i,n)
299 za(i,n) = a4*zz(i,n)
300
301 xb(i,n) = b4*xx(i,n)
302 yb(i,n) = b4*yy(i,n)
303 zb(i,n) = b4*zz(i,n)
304 ENDDO
305 ENDDO
306 DO ip=1,4
307 k1 = iperm(1,ip)
308 k2 = iperm(2,ip)
309 k3 = iperm(3,ip)
310 k4 = iperm(4,ip)
311 k5 = iperm(5,ip)
312 k6 = iperm(6,ip)
313 k7 = iperm(7,ip)
314 k8 = iperm(8,ip)
315 k9 = iperm(9,ip)
316 k10= iperm(10,ip)
317 CALL s10jacob1(alph(ip),beta(ip),wip(ip),
318 . xb(1,k1),xb(1,k2),xb(1,k3),xa(1,k4),xb(1,k5),xb(1,k6),xb(1,k7),
319 . xb(1,k8),xb(1,k9),xb(1,k10),xa(1,k8),xa(1,k9),xa(1,k10),
320 . yb(1,k1),yb(1,k2),yb(1,k3),ya(1,k4),yb(1,k5),yb(1,k6),yb(1,k7),
321 . yb(1,k8),yb(1,k9),yb(1,k10),ya(1,k8),ya(1,k9),ya(1,k10),
322 . zb(1,k1),zb(1,k2),zb(1,k3),za(1,k4),zb(1,k5),zb(1,k6),zb(1,k7),
323 . zb(1,k8),zb(1,k9),zb(1,k10),za(1,k8),za(1,k9),za(1,k10),
324 . px(1,k1,ip) ,px(1,k2,ip),px(1,k3,ip),px(1,k4,ip),px(1,k5,ip),
325 . px(1,k6,ip) ,px(1,k7,ip),px(1,k8,ip),px(1,k9,ip),px(1,k10,ip),
326 . py(1,k1,ip) ,py(1,k2,ip),py(1,k3,ip),py(1,k4,ip),py(1,k5,ip),
327 . py(1,k6,ip) ,py(1,k7,ip),py(1,k8,ip),py(1,k9,ip),py(1,k10,ip),
328 . pz(1,k1,ip) ,pz(1,k2,ip),pz(1,k3,ip),pz(1,k4,ip),pz(1,k5,ip),
329 . pz(1,k6,ip) ,pz(1,k7,ip),pz(1,k8,ip),pz(1,k9,ip),pz(1,k10,ip),
330 . nx(1,k1,ip) ,nx(1,k2,ip),nx(1,k3,ip),nx(1,k4,ip),nx(1,k5,ip),
331 . nx(1,k6,ip) ,nx(1,k7,ip),nx(1,k8,ip),nx(1,k9,ip),nx(1,k10,ip),
332 . vol(1,ip) ,nnega, index ,voldp(1,ip))
333 ENDDO
334
335
336 IF(npt==5)THEN
337 ip = 5
338
339#include "vectorize.inc"
340 DO j=1,nnega
341 i = index(j)
342 xa(i,1) = zero
343 ENDDO
344
345 CALL s10jacob1(alph(ip),beta(ip),wip(ip),
346 . xa(1,1) ,xa(1,1) ,xa(1,1) ,xa(1,1) ,xx(1,k5),
347 . xx(1,k6),xx(1,k7),xx(1,k8),xx(1,k9),xx(1,k10),
348 . xx(1,k8),xx(1,k9),xx(1,k10),
349 . xa(1,1) ,xa(1,1) ,xa(1,1) ,xa(1,1) ,yy(1,k5),
350 . yy(1,k6),yy(1,k7),yy(1,k8),yy(1,k9),yy(1,k10),
351 . yy(1,k8),yy(1,k9),yy(1,k10),
352 . xa(1,1) ,xa(1,1) ,xa(1,1) ,xa(1,1) ,zz(1,k5),
353 . zz(1,k6),zz(1,k7),zz(1,k8),zz(1,k9),zz(1,k10),
354 . zz(1,k8),zz(1,k9),zz(1,k10),
355 . px(1,k1,ip) ,px(1,k2,ip),px(1,k3,ip),px(1,k4,ip),px(1,k5,ip),
356 . px(1,k6,ip) ,px(1,k7,ip),px(1,k8,ip),px(1,k9,ip),px(1,k10,ip),
357 . py(1,k1,ip) ,py(1,k2,ip),py(1,k3,ip),py(1,k4,ip),py(1,k5,ip),
358 . py(1,k6,ip) ,py(1,k7,ip),py(1,k8,ip),py(1,k9,ip),py(1,k10,ip),
359 . pz(1,k1,ip) ,pz(1,k2,ip),pz(1,k3,ip),pz(1,k4,ip),pz(1,k5,ip),
360 . pz(1,k6,ip) ,pz(1,k7,ip),pz(1,k8,ip),pz(1,k9,ip),pz(1,k10,ip),
361 . nx(1,k1,ip) ,nx(1,k2,ip),nx(1,k3,ip),nx(1,k4,ip),nx(1,k5,ip),
362 . nx(1,k6,ip) ,nx(1,k7,ip),nx(1,k8,ip),nx(1,k9,ip),nx(1,k10,ip),
363 . vol(1,ip) ,nnega, index ,voldp(1,ip))
364 ENDIF
365 IF(ineg_v==0)THEN
366 CALL ancmsg(msgid=280,anmode=aninfo)
367 mstop = 1
368 END IF !(INEG_V==0)THEN
369 ENDIF
370C
371 DO ip=1,npt
372 DO i=1,nel
373 volg(i) =volg(i) + vol(i,ip)
374 ENDDO
375 ENDDO
376C
377 DO i=1,nel
378 voln(i) =volg(i)/npt
379 ENDDO
380C-----------
381 RETURN
#define my_real
Definition cppsort.cpp:32
if(complex_arithmetic) id
subroutine s10jacob1(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, nnega, index, voldp)
Definition s10jacob1.F:37
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
subroutine ancmsg(msgid, msgtype, anmode, i1, i2, i3, i4, i5, i6, i7, i8, i9, i10, i11, i12, i13, i14, i15, i16, i17, i18, i19, i20, r1, r2, r3, r4, r5, r6, r7, r8, r9, c1, c2, c3, c4, c5, c6, c7, c8, c9, prmode)
Definition message.F:895