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

Go to the source code of this file.

Functions/Subroutines

subroutine s8zderic3 (off, det, ngl, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, px1, px2, px3, px4, py1, py2, py3, py4, pz1, pz2, pz3, pz4, px1h1, px1h2, px1h3, px1h4, px2h1, px2h2, px2h3, px2h4, px3h1, px3h2, px3h3, px3h4, px4h1, px4h2, px4h3, px4h4, hx, hy, hz, jac1, jac2, jac3, jac4, jac5, jac6, jac7, jac8, jac9, smax, sav, offg, nel, ismstr, jlag)

Function/Subroutine Documentation

◆ s8zderic3()

subroutine s8zderic3 ( off,
det,
integer, dimension(*) ngl,
double precision, dimension(*) x1,
double precision, dimension(*) x2,
double precision, dimension(*) x3,
double precision, dimension(*) x4,
double precision, dimension(*) x5,
double precision, dimension(*) x6,
double precision, dimension(*) x7,
double precision, dimension(*) x8,
double precision, dimension(*) y1,
double precision, dimension(*) y2,
double precision, dimension(*) y3,
double precision, dimension(*) y4,
double precision, dimension(*) y5,
double precision, dimension(*) y6,
double precision, dimension(*) y7,
double precision, dimension(*) y8,
double precision, dimension(*) z1,
double precision, dimension(*) z2,
double precision, dimension(*) z3,
double precision, dimension(*) z4,
double precision, dimension(*) z5,
double precision, dimension(*) z6,
double precision, dimension(*) z7,
double precision, dimension(*) z8,
px1,
px2,
px3,
px4,
py1,
py2,
py3,
py4,
pz1,
pz2,
pz3,
pz4,
px1h1,
px1h2,
px1h3,
px1h4,
px2h1,
px2h2,
px2h3,
px2h4,
px3h1,
px3h2,
px3h3,
px3h4,
px4h1,
px4h2,
px4h3,
px4h4,
hx,
hy,
hz,
jac1,
jac2,
jac3,
jac4,
jac5,
jac6,
jac7,
jac8,
jac9,
smax,
double precision, dimension(nel,21) sav,
offg,
integer nel,
integer, intent(in) ismstr,
integer, intent(in) jlag )

Definition at line 35 of file s8zderic3.F.

55C-----------------------------------------------
56C M o d u l e s
57C-----------------------------------------------
58 USE message_mod
59C-----------------------------------------------
60C I m p l i c i t T y p e s
61C-----------------------------------------------
62#include "implicit_f.inc"
63#include "comlock.inc"
64C-----------------------------------------------
65C G l o b a l P a r a m e t e r s
66C-----------------------------------------------
67#include "mvsiz_p.inc"
68C-----------------------------------------------
69C C o m m o n B l o c k s
70C-----------------------------------------------
71C-----------------------------------------------
72C D u m m y A r g u m e n t s
73C-----------------------------------------------
74 INTEGER, INTENT(IN) :: JLAG
75 INTEGER, INTENT(IN) :: ISMSTR
76 INTEGER NEL
77C REAL
78 double precision
79 . x1(*), x2(*), x3(*), x4(*), x5(*), x6(*), x7(*), x8(*),
80 . y1(*), y2(*), y3(*), y4(*), y5(*), y6(*), y7(*), y8(*),
81 . z1(*), z2(*), z3(*), z4(*), z5(*), z6(*), z7(*), z8(*),
82 . sav(nel,21)
83
85 . off(*),det(*),
86 . px1(*), px2(*), px3(*), px4(*),
87 . py1(*), py2(*), py3(*), py4(*),
88 . pz1(*), pz2(*), pz3(*), pz4(*),
89 . px1h1(*), px1h2(*), px1h3(*),px1h4(*),
90 . px2h1(*), px2h2(*), px2h3(*),px2h4(*),
91 . px3h1(*), px3h2(*), px3h3(*),px3h4(*),
92 . px4h1(*), px4h2(*), px4h3(*),px4h4(*),
93 . hx(mvsiz,4), hy(mvsiz,4), hz(mvsiz,4),
94 . jac1(*),jac2(*),jac3(*),
95 . jac4(*),jac5(*),jac6(*),
96 . jac7(*),jac8(*),jac9(*),smax(*),offg(*)
97C-----------------------------------------------
98C L o c a l V a r i a b l e s
99C-----------------------------------------------
100 INTEGER NGL(*), I, J ,ICOR,NNEGA,INDEX(MVSIZ)
101C REAL
102C 12
103 my_real
104 . dett ,
105 . jaci1, jaci2, jaci3,
106 . jaci4, jaci5, jaci6,
107 . jaci7, jaci8, jaci9,
108 . x17 , x28 , x35 , x46,
109 . y17 , y28 , y35 , y46,
110 . z17 , z28 , z35 , z46,
111 . jac_59_68(mvsiz), jac_67_49(mvsiz), jac_48_57(mvsiz),
112 . jac_38_29(mvsiz), jac_19_37(mvsiz), jac_27_18(mvsiz),
113 . jac_26_35(mvsiz), jac_34_16(mvsiz), jac_15_24(mvsiz),
114 . jaci12, jaci45, jaci78,
115 . x_17_46 , x_28_35 ,
116 . y_17_46 , y_28_35 ,
117 . z_17_46 , z_28_35
118C-----------------------------------------------
119C
120 DO i=1,nel
121 x17=x7(i)-x1(i)
122 x28=x8(i)-x2(i)
123 x35=x5(i)-x3(i)
124 x46=x6(i)-x4(i)
125 y17=y7(i)-y1(i)
126 y28=y8(i)-y2(i)
127 y35=y5(i)-y3(i)
128 y46=y6(i)-y4(i)
129 z17=z7(i)-z1(i)
130 z28=z8(i)-z2(i)
131 z35=z5(i)-z3(i)
132 z46=z6(i)-z4(i)
133C
134 jac4(i)=x17+x28-x35-x46
135 jac5(i)=y17+y28-y35-y46
136 jac6(i)=z17+z28-z35-z46
137 x_17_46=x17+x46
138 x_28_35=x28+x35
139 y_17_46=y17+y46
140 y_28_35=y28+y35
141 z_17_46=z17+z46
142 z_28_35=z28+z35
143 jac7(i)=x_17_46+x_28_35
144 jac8(i)=y_17_46+y_28_35
145 jac9(i)=z_17_46+z_28_35
146 jac1(i)=x_17_46-x_28_35
147 jac2(i)=y_17_46-y_28_35
148 jac3(i)=z_17_46-z_28_35
149 ENDDO
150C
151C JACOBIAN
152C
153 DO i=1,nel
154 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
155 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
156 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
157 END DO
158C
159 DO i=1,nel
160 det(i)=one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
161 END DO
162C
163 CALL schkjabt3(
164 1 off, det, ngl, offg,
165 2 nnega, index, nel, ismstr,
166 3 jlag)
167 IF (nnega>0) THEN
168 IF (ismstr==10.OR.ismstr==12) THEN
169#include "vectorize.inc"
170 DO j=1,nnega
171 i = index(j)
172 x1(i)=sav(i,1)
173 y1(i)=sav(i,8)
174 z1(i)=sav(i,15)
175 x2(i)=sav(i,2)
176 y2(i)=sav(i,9)
177 z2(i)=sav(i,16)
178 x3(i)=sav(i,3)
179 y3(i)=sav(i,10)
180 z3(i)=sav(i,17)
181 x4(i)=sav(i,4)
182 y4(i)=sav(i,11)
183 z4(i)=sav(i,18)
184 x5(i)=sav(i,5)
185 y5(i)=sav(i,12)
186 z5(i)=sav(i,19)
187 x6(i)=sav(i,6)
188 y6(i)=sav(i,13)
189 z6(i)=sav(i,20)
190 x7(i)=sav(i,7)
191 y7(i)=sav(i,14)
192 z7(i)=sav(i,21)
193 x8(i)=zero
194 y8(i)=zero
195 z8(i)=zero
196 ENDDO
197 ELSE
198 DO j=1,nnega
199 i = index(j)
200 x1(i)=sav(i,1)
201 y1(i)=sav(i,2)
202 z1(i)=sav(i,3)
203 x2(i)=sav(i,4)
204 y2(i)=sav(i,5)
205 z2(i)=sav(i,6)
206 x3(i)=sav(i,7)
207 y3(i)=sav(i,8)
208 z3(i)=sav(i,9)
209 x4(i)=sav(i,10)
210 y4(i)=sav(i,11)
211 z4(i)=sav(i,12)
212 x5(i)=sav(i,13)
213 y5(i)=sav(i,14)
214 z5(i)=sav(i,15)
215 x6(i)=sav(i,16)
216 y6(i)=sav(i,17)
217 z6(i)=sav(i,18)
218 x7(i)=sav(i,19)
219 y7(i)=sav(i,20)
220 z7(i)=sav(i,21)
221 x8(i)=zero
222 y8(i)=zero
223 z8(i)=zero
224 ENDDO
225 END IF !(ISMSTR==10)
226#include "vectorize.inc"
227 DO j=1,nnega
228 i = index(j)
229C
230 x17=x7(i)-x1(i)
231 x28=x8(i)-x2(i)
232 x35=x5(i)-x3(i)
233 x46=x6(i)-x4(i)
234 y17=y7(i)-y1(i)
235 y28=y8(i)-y2(i)
236 y35=y5(i)-y3(i)
237 y46=y6(i)-y4(i)
238 z17=z7(i)-z1(i)
239 z28=z8(i)-z2(i)
240 z35=z5(i)-z3(i)
241 z46=z6(i)-z4(i)
242C
243 jac4(i)=x17+x28-x35-x46
244 jac5(i)=y17+y28-y35-y46
245 jac6(i)=z17+z28-z35-z46
246 x_17_46=x17+x46
247 x_28_35=x28+x35
248 y_17_46=y17+y46
249 y_28_35=y28+y35
250 z_17_46=z17+z46
251 z_28_35=z28+z35
252 jac7(i)=x_17_46+x_28_35
253 jac8(i)=y_17_46+y_28_35
254 jac9(i)=z_17_46+z_28_35
255 jac1(i)=x_17_46-x_28_35
256 jac2(i)=y_17_46-y_28_35
257 jac3(i)=z_17_46-z_28_35
258C
259 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
260 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
261 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
262C
263 det(i)=one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
264 offg(i) = two
265 ENDDO
266 END IF
267C
268 DO i=1,nel
269 jac_38_29(i)=(-jac2(i)*jac9(i)+jac3(i)*jac8(i))
270 jac_19_37(i)=( jac1(i)*jac9(i)-jac3(i)*jac7(i))
271 jac_27_18(i)=(-jac1(i)*jac8(i)+jac2(i)*jac7(i))
272 jac_26_35(i)=( jac2(i)*jac6(i)-jac3(i)*jac5(i))
273 jac_34_16(i)=(-jac1(i)*jac6(i)+jac3(i)*jac4(i))
274 jac_15_24(i)=( jac1(i)*jac5(i)-jac2(i)*jac4(i))
275 END DO
276C
277C Jacobian matrix inverse
278C
279 DO i=1,nel
280 dett=one_over_64/det(i)
281 jaci1=dett*jac_59_68(i)
282 jaci4=dett*jac_67_49(i)
283 jaci7=dett*jac_48_57(i)
284 jaci2=dett*jac_38_29(i)
285 jaci5=dett*jac_19_37(i)
286 jaci8=dett*jac_27_18(i)
287 jaci3=dett*jac_26_35(i)
288 jaci6=dett*jac_34_16(i)
289 jaci9=dett*jac_15_24(i)
290C
291 jaci12=jaci1-jaci2
292 jaci45=jaci4-jaci5
293 jaci78=jaci7-jaci8
294
295 px2(i)= jaci12-jaci3
296 py2(i)= jaci45-jaci6
297 pz2(i)= jaci78-jaci9
298 px4(i)=-jaci12-jaci3
299 py4(i)=-jaci45-jaci6
300 pz4(i)=-jaci78-jaci9
301
302 jaci12=jaci1+jaci2
303 jaci45=jaci4+jaci5
304 jaci78=jaci7+jaci8
305
306 px1(i)=-jaci12-jaci3
307 py1(i)=-jaci45-jaci6
308 pz1(i)=-jaci78-jaci9
309 px3(i)= jaci12-jaci3
310 py3(i)= jaci45-jaci6
311 pz3(i)= jaci78-jaci9
312 ENDDO
313C
314C mode 1
315C 1 1 -1 -1 -1 -1 1 1
316 DO i=1,nel
317 hx(i,1)=(x1(i)+x2(i)-x3(i)-x4(i)-x5(i)-x6(i)+x7(i)+x8(i))
318 hy(i,1)=(y1(i)+y2(i)-y3(i)-y4(i)-y5(i)-y6(i)+y7(i)+y8(i))
319 hz(i,1)=(z1(i)+z2(i)-z3(i)-z4(i)-z5(i)-z6(i)+z7(i)+z8(i))
320 px1h1(i)=px1(i)*hx(i,1)+ py1(i)*hy(i,1)+pz1(i)*hz(i,1)
321 px2h1(i)=px2(i)*hx(i,1)+ py2(i)*hy(i,1)+pz2(i)*hz(i,1)
322 px3h1(i)=px3(i)*hx(i,1)+ py3(i)*hy(i,1)+pz3(i)*hz(i,1)
323 px4h1(i)=px4(i)*hx(i,1)+ py4(i)*hy(i,1)+pz4(i)*hz(i,1)
324 ENDDO
325C mode 2
326C 1 -1 -1 1 -1 1 1 -1
327 DO i=1,nel
328 hx(i,2)=(x1(i)-x2(i)-x3(i)+x4(i)-x5(i)+x6(i)+x7(i)-x8(i))
329 hy(i,2)=(y1(i)-y2(i)-y3(i)+y4(i)-y5(i)+y6(i)+y7(i)-y8(i))
330 hz(i,2)=(z1(i)-z2(i)-z3(i)+z4(i)-z5(i)+z6(i)+z7(i)-z8(i))
331 px1h2(i)=px1(i)*hx(i,2)+ py1(i)*hy(i,2)+pz1(i)*hz(i,2)
332 px2h2(i)=px2(i)*hx(i,2)+ py2(i)*hy(i,2)+pz2(i)*hz(i,2)
333 px3h2(i)=px3(i)*hx(i,2)+ py3(i)*hy(i,2)+pz3(i)*hz(i,2)
334 px4h2(i)=px4(i)*hx(i,2)+ py4(i)*hy(i,2)+pz4(i)*hz(i,2)
335 ENDDO
336C mode 3
337C 1 -1 1 -1 1 -1 1 -1
338 DO i=1,nel
339 hx(i,3)=(x1(i)-x2(i)+x3(i)-x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
340 hy(i,3)=(y1(i)-y2(i)+y3(i)-y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
341 hz(i,3)=(z1(i)-z2(i)+z3(i)-z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
342 px1h3(i)=px1(i)*hx(i,3)+ py1(i)*hy(i,3)+pz1(i)*hz(i,3)
343 px2h3(i)=px2(i)*hx(i,3)+ py2(i)*hy(i,3)+pz2(i)*hz(i,3)
344 px3h3(i)=px3(i)*hx(i,3)+ py3(i)*hy(i,3)+pz3(i)*hz(i,3)
345 px4h3(i)=px4(i)*hx(i,3)+ py4(i)*hy(i,3)+pz4(i)*hz(i,3)
346 ENDDO
347C mode 4
348C -1 1 -1 1 1 -1 1 -1
349 DO i=1,nel
350 hx(i,4)=(-x1(i)+x2(i)-x3(i)+x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
351 hy(i,4)=(-y1(i)+y2(i)-y3(i)+y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
352 hz(i,4)=(-z1(i)+z2(i)-z3(i)+z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
353 px1h4(i)=px1(i)*hx(i,4)+ py1(i)*hy(i,4)+pz1(i)*hz(i,4)
354 px2h4(i)=px2(i)*hx(i,4)+ py2(i)*hy(i,4)+pz2(i)*hz(i,4)
355 px3h4(i)=px3(i)*hx(i,4)+ py3(i)*hy(i,4)+pz3(i)*hz(i,4)
356 px4h4(i)=px4(i)*hx(i,4)+ py4(i)*hy(i,4)+pz4(i)*hz(i,4)
357 ENDDO
358C----surface max mediane-- *16
359 DO i=1,nel
360 smax(i)= jac_59_68(i)*jac_59_68(i)+jac_67_49(i)*jac_67_49(i)
361 . +jac_48_57(i)*jac_48_57(i)
362 smax(i)= max(smax(i),jac_38_29(i)*jac_38_29(i)+jac_19_37(i)*jac_19_37(i)
363 . +jac_27_18(i)*jac_27_18(i))
364 smax(i)= max(smax(i),jac_26_35(i)*jac_26_35(i)+jac_34_16(i)*jac_34_16(i)
365 . +jac_15_24(i)*jac_15_24(i))
366 ENDDO
367 DO i=1,nel
368 IF(smax(i)<=zero)THEN
369 CALL ancmsg(msgid=173,anmode=aninfo,
370 . i1=ngl(i))
371 CALL arret(2)
372 ENDIF
373 smax(i)= one/sqrt(smax(i))
374 ENDDO
375 RETURN
376C
377 1000 FORMAT(/' ZERO OR NEGATIVE VOLUME : 3D-ELEMENT NB',i10/)
378 2000 FORMAT(/' ZERO OR NEGATIVE VOLUME : DELETE 3D-ELEMENT NB',i10/)
#define my_real
Definition cppsort.cpp:32
#define max(a, b)
Definition macros.h:21
subroutine schkjabt3(off, det, ngl, offg, nnega, index, nel, ismstr, jlag)
Definition schkjabt3.F:40
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:889
subroutine arret(nn)
Definition arret.F:87