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

Go to the source code of this file.

Functions/Subroutines

subroutine s8zderict3 (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, nel)

Function/Subroutine Documentation

◆ s8zderict3()

subroutine s8zderict3 ( 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,
integer, intent(in) nel )

Definition at line 36 of file s8zderict3.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) :: NEL
75 double precision
76 . x1(*), x2(*), x3(*), x4(*), x5(*), x6(*), x7(*), x8(*),
77 . y1(*), y2(*), y3(*), y4(*), y5(*), y6(*), y7(*), y8(*),
78 . z1(*), z2(*), z3(*), z4(*), z5(*), z6(*), z7(*), z8(*)
79C REAL
81 . off(*),det(*),
82 . px1(*), px2(*), px3(*), px4(*),
83 . py1(*), py2(*), py3(*), py4(*),
84 . pz1(*), pz2(*), pz3(*), pz4(*),
85 . px1h1(*), px1h2(*), px1h3(*),px1h4(*),
86 . px2h1(*), px2h2(*), px2h3(*),px2h4(*),
87 . px3h1(*), px3h2(*), px3h3(*),px3h4(*),
88 . px4h1(*), px4h2(*), px4h3(*),px4h4(*),
89 . hx(mvsiz,4), hy(mvsiz,4), hz(mvsiz,4),
90 . jac1(*),jac2(*),jac3(*),
91 . jac4(*),jac5(*),jac6(*),
92 . jac7(*),jac8(*),jac9(*),smax(*)
93C-----------------------------------------------
94C L o c a l V a r i a b l e s
95C-----------------------------------------------
96 INTEGER NGL(*), I, J ,ICOR
97C REAL
98C 12
100 . dett ,
101 . jaci1, jaci2, jaci3,
102 . jaci4, jaci5, jaci6,
103 . jaci7, jaci8, jaci9,
104 . x17 , x28 , x35 , x46,
105 . y17 , y28 , y35 , y46,
106 . z17 , z28 , z35 , z46,
107 . jac_59_68(mvsiz), jac_67_49(mvsiz), jac_48_57(mvsiz),
108 . jac_38_29(mvsiz), jac_19_37(mvsiz), jac_27_18(mvsiz),
109 . jac_26_35(mvsiz), jac_34_16(mvsiz), jac_15_24(mvsiz),
110 . jaci12, jaci45, jaci78,
111 . x_17_46 , x_28_35 ,
112 . y_17_46 , y_28_35 ,
113 . z_17_46 , z_28_35
114C-----------------------------------------------
115C
116 DO i=1,nel
117 x17=x7(i)-x1(i)
118 x28=x8(i)-x2(i)
119 x35=x5(i)-x3(i)
120 x46=x6(i)-x4(i)
121 y17=y7(i)-y1(i)
122 y28=y8(i)-y2(i)
123 y35=y5(i)-y3(i)
124 y46=y6(i)-y4(i)
125 z17=z7(i)-z1(i)
126 z28=z8(i)-z2(i)
127 z35=z5(i)-z3(i)
128 z46=z6(i)-z4(i)
129C
130C Jacobian matrix
131 jac4(i)=x17+x28-x35-x46
132 jac5(i)=y17+y28-y35-y46
133 jac6(i)=z17+z28-z35-z46
134 x_17_46=x17+x46
135 x_28_35=x28+x35
136 y_17_46=y17+y46
137 y_28_35=y28+y35
138 z_17_46=z17+z46
139 z_28_35=z28+z35
140
141 jac7(i)=x_17_46+x_28_35
142 jac8(i)=y_17_46+y_28_35
143 jac9(i)=z_17_46+z_28_35
144 jac1(i)=x_17_46-x_28_35
145 jac2(i)=y_17_46-y_28_35
146 jac3(i)=z_17_46-z_28_35
147 ENDDO
148C
149 DO i=1,nel
150 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
151 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
152 jac_38_29(i)=(-jac2(i)*jac9(i)+jac3(i)*jac8(i))
153 jac_19_37(i)=( jac1(i)*jac9(i)-jac3(i)*jac7(i))
154 jac_27_18(i)=(-jac1(i)*jac8(i)+jac2(i)*jac7(i))
155 jac_26_35(i)=( jac2(i)*jac6(i)-jac3(i)*jac5(i))
156 jac_34_16(i)=(-jac1(i)*jac6(i)+jac3(i)*jac4(i))
157 jac_15_24(i)=( jac1(i)*jac5(i)-jac2(i)*jac4(i))
158 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
159 END DO
160C
161 DO i=1,nel
162 det(i)=one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
163 ENDDO
164C
165 CALL schkjab3(
166 1 off, det, ngl, nel)
167C
168C Jacobian matrix inverse
169C
170 DO i=1,nel
171 dett=one_over_64/det(i)
172 jaci1=dett*jac_59_68(i)
173 jaci4=dett*jac_67_49(i)
174 jaci7=dett*jac_48_57(i)
175 jaci2=dett*jac_38_29(i)
176 jaci5=dett*jac_19_37(i)
177 jaci8=dett*jac_27_18(i)
178 jaci3=dett*jac_26_35(i)
179 jaci6=dett*jac_34_16(i)
180 jaci9=dett*jac_15_24(i)
181!
182 jaci12=jaci1-jaci2
183 jaci45=jaci4-jaci5
184 jaci78=jaci7-jaci8
185!
186 px2(i)= jaci12-jaci3
187 py2(i)= jaci45-jaci6
188 pz2(i)= jaci78-jaci9
189 px4(i)=-jaci12-jaci3
190 py4(i)=-jaci45-jaci6
191 pz4(i)=-jaci78-jaci9
192!
193 jaci12=jaci1+jaci2
194 jaci45=jaci4+jaci5
195 jaci78=jaci7+jaci8
196!
197 px1(i)=-jaci12-jaci3
198 py1(i)=-jaci45-jaci6
199 pz1(i)=-jaci78-jaci9
200 px3(i)=jaci12-jaci3
201 py3(i)=jaci45-jaci6
202 pz3(i)=jaci78-jaci9
203 ENDDO
204C
205C mode 1
206C 1 1 -1 -1 -1 -1 1 1
207 DO i=1,nel
208 hx(i,1)=(x1(i)+x2(i)-x3(i)-x4(i)-x5(i)-x6(i)+x7(i)+x8(i))
209 hy(i,1)=(y1(i)+y2(i)-y3(i)-y4(i)-y5(i)-y6(i)+y7(i)+y8(i))
210 hz(i,1)=(z1(i)+z2(i)-z3(i)-z4(i)-z5(i)-z6(i)+z7(i)+z8(i))
211 px1h1(i)=px1(i)*hx(i,1)+ py1(i)*hy(i,1)+pz1(i)*hz(i,1)
212 px2h1(i)=px2(i)*hx(i,1)+ py2(i)*hy(i,1)+pz2(i)*hz(i,1)
213 px3h1(i)=px3(i)*hx(i,1)+ py3(i)*hy(i,1)+pz3(i)*hz(i,1)
214 px4h1(i)=px4(i)*hx(i,1)+ py4(i)*hy(i,1)+pz4(i)*hz(i,1)
215 ENDDO
216C mode 2
217C 1 -1 -1 1 -1 1 1 -1
218 DO i=1,nel
219 hx(i,2)=(x1(i)-x2(i)-x3(i)+x4(i)-x5(i)+x6(i)+x7(i)-x8(i))
220 hy(i,2)=(y1(i)-y2(i)-y3(i)+y4(i)-y5(i)+y6(i)+y7(i)-y8(i))
221 hz(i,2)=(z1(i)-z2(i)-z3(i)+z4(i)-z5(i)+z6(i)+z7(i)-z8(i))
222 px1h2(i)=px1(i)*hx(i,2)+ py1(i)*hy(i,2)+pz1(i)*hz(i,2)
223 px2h2(i)=px2(i)*hx(i,2)+ py2(i)*hy(i,2)+pz2(i)*hz(i,2)
224 px3h2(i)=px3(i)*hx(i,2)+ py3(i)*hy(i,2)+pz3(i)*hz(i,2)
225 px4h2(i)=px4(i)*hx(i,2)+ py4(i)*hy(i,2)+pz4(i)*hz(i,2)
226 ENDDO
227C mode 3
228C 1 -1 1 -1 1 -1 1 -1
229 DO i=1,nel
230 hx(i,3)=(x1(i)-x2(i)+x3(i)-x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
231 hy(i,3)=(y1(i)-y2(i)+y3(i)-y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
232 hz(i,3)=(z1(i)-z2(i)+z3(i)-z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
233 px1h3(i)=px1(i)*hx(i,3)+ py1(i)*hy(i,3)+pz1(i)*hz(i,3)
234 px2h3(i)=px2(i)*hx(i,3)+ py2(i)*hy(i,3)+pz2(i)*hz(i,3)
235 px3h3(i)=px3(i)*hx(i,3)+ py3(i)*hy(i,3)+pz3(i)*hz(i,3)
236 px4h3(i)=px4(i)*hx(i,3)+ py4(i)*hy(i,3)+pz4(i)*hz(i,3)
237 ENDDO
238C mode 4
239C -1 1 -1 1 1 -1 1 -1
240 DO i=1,nel
241 hx(i,4)=(-x1(i)+x2(i)-x3(i)+x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
242 hy(i,4)=(-y1(i)+y2(i)-y3(i)+y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
243 hz(i,4)=(-z1(i)+z2(i)-z3(i)+z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
244 px1h4(i)=px1(i)*hx(i,4)+ py1(i)*hy(i,4)+pz1(i)*hz(i,4)
245 px2h4(i)=px2(i)*hx(i,4)+ py2(i)*hy(i,4)+pz2(i)*hz(i,4)
246 px3h4(i)=px3(i)*hx(i,4)+ py3(i)*hy(i,4)+pz3(i)*hz(i,4)
247 px4h4(i)=px4(i)*hx(i,4)+ py4(i)*hy(i,4)+pz4(i)*hz(i,4)
248 ENDDO
249C----surface max mediane-- *16
250 DO i=1,nel
251 smax(i)= jac_59_68(i)*jac_59_68(i)+jac_67_49(i)*jac_67_49(i)
252 . +jac_48_57(i)*jac_48_57(i)
253 smax(i)= max(smax(i),jac_38_29(i)*jac_38_29(i)+jac_19_37(i)*jac_19_37(i)
254 . +jac_27_18(i)*jac_27_18(i))
255 smax(i)= max(smax(i),jac_26_35(i)*jac_26_35(i)+jac_34_16(i)*jac_34_16(i)
256 . +jac_15_24(i)*jac_15_24(i))
257 ENDDO
258 DO i=1,nel
259 IF(smax(i)<=zero)THEN
260 CALL ancmsg(msgid=173,anmode=aninfo,
261 . i1=ngl(i))
262 CALL arret(2)
263 ENDIF
264 smax(i)= one/sqrt(smax(i))
265 ENDDO
266 RETURN
267C
268 1000 FORMAT(/' ZERO OR NEGATIVE VOLUME : 3D-ELEMENT NB',i10/)
269 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 schkjab3(off, det, ngl, nel)
Definition schkjab3.F:39
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