OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
szderi3.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 szderi3 (off, det, ngl, ismstr, 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, jac1, jac2, jac3, jac4, jac5, jac6, jac9, sav, offg, nel, voldp, jlag)
subroutine szderit3 (off, det, ngl, nel, 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, jac1, jac2, jac3, jac4, jac5, jac6, jac9, jlag)
subroutine szderito3 (off, det, px1, px2, px3, px4, py1, py2, py3, py4, pz1, pz2, pz3, pz4, jac_i, nel, jlag)

Function/Subroutine Documentation

◆ szderi3()

subroutine szderi3 ( off,
det,
integer, dimension(*) ngl,
integer ismstr,
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,
jac1,
jac2,
jac3,
jac4,
jac5,
jac6,
jac9,
double precision, dimension(nel,21) sav,
offg,
integer nel,
double precision, dimension(*) voldp,
integer, intent(in) jlag )

Definition at line 30 of file szderi3.F.

48C-----------------------------------------------
49C I m p l i c i t T y p e s
50C-----------------------------------------------
51#include "implicit_f.inc"
52#include "comlock.inc"
53C-----------------------------------------------
54C G l o b a l P a r a m e t e r s
55C-----------------------------------------------
56#include "mvsiz_p.inc"
57C-----------------------------------------------
58C D u m m y A r g u m e n t s
59C-----------------------------------------------
60 INTEGER, INTENT(IN) :: JLAG
61 INTEGER :: ISMSTR,NEL,NGL(*)
62 double precision
63 . x1(*), x2(*), x3(*), x4(*), x5(*), x6(*), x7(*), x8(*),
64 . y1(*), y2(*), y3(*), y4(*), y5(*), y6(*), y7(*), y8(*),
65 . z1(*), z2(*), z3(*), z4(*), z5(*), z6(*), z7(*), z8(*),
66 . sav(nel,21),voldp(*)
67C
69 . off(*),det(*),
70 . px1(*), px2(*), px3(*), px4(*),
71 . py1(*), py2(*), py3(*), py4(*),
72 . pz1(*), pz2(*), pz3(*), pz4(*),
73 . px1h1(*), px1h2(*), px1h3(*),px1h4(*),
74 . px2h1(*), px2h2(*), px2h3(*),px2h4(*),
75 . px3h1(*), px3h2(*), px3h3(*),px3h4(*),
76 . px4h1(*), px4h2(*), px4h3(*),px4h4(*),
77 . jac1(*),jac2(*),jac3(*),
78 . jac4(*),jac5(*),jac6(*),jac9(*),offg(*)
79C-----------------------------------------------
80C L o c a l V a r i a b l e s
81C-----------------------------------------------
82 INTEGER :: I,J,NNEGA,INDEX(MVSIZ)
83C
85 . dett , jac7(mvsiz), jac8(mvsiz) ,
86 . jaci1, jaci2, jaci3,
87 . jaci4, jaci5, jaci6,
88 . jaci7, jaci8, jaci9,
89 . jac_59_68(mvsiz), jac_67_49(mvsiz), jac_48_57(mvsiz),
90 . jaci12, jaci45, jaci78,
91 . x_17_46 , x_28_35 ,
92 . y_17_46 , y_28_35 ,
93 . z_17_46 , z_28_35 ,
94 . hx,hy,hz
95 double precision
96 . x17 , x28 , x35 , x46,
97 . y17 , y28 , y35 , y46,
98 . z17 , z28 , z35 , z46
99C-----------------------------------------------
100 DO i=1,nel
101 x17=x7(i)-x1(i)
102 x28=x8(i)-x2(i)
103 x35=x5(i)-x3(i)
104 x46=x6(i)-x4(i)
105 y17=y7(i)-y1(i)
106 y28=y8(i)-y2(i)
107 y35=y5(i)-y3(i)
108 y46=y6(i)-y4(i)
109 z17=z7(i)-z1(i)
110 z28=z8(i)-z2(i)
111 z35=z5(i)-z3(i)
112 z46=z6(i)-z4(i)
113C
114 jac4(i)=x17+x28-x35-x46
115 jac5(i)=y17+y28-y35-y46
116 jac6(i)=z17+z28-z35-z46
117 x_17_46=x17+x46
118 x_28_35=x28+x35
119 y_17_46=y17+y46
120 y_28_35=y28+y35
121 z_17_46=z17+z46
122 z_28_35=z28+z35
123 jac7(i)=x_17_46+x_28_35
124 jac8(i)=y_17_46+y_28_35
125 jac9(i)=z_17_46+z_28_35
126 jac1(i)=x_17_46-x_28_35
127 jac2(i)=y_17_46-y_28_35
128 jac3(i)=z_17_46-z_28_35
129 ENDDO
130C
131 DO i=1,nel
132 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
133 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
134 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
135 END DO
136C
137 DO i=1,nel
138 voldp(i)=one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
139 det(i)=voldp(i)
140 END DO
141C
142 CALL schkjabt3(
143 1 off, det, ngl, offg,
144 2 nnega, index, nel, ismstr,
145 3 jlag)
146 IF (nnega>0) THEN
147 IF (ismstr==10.OR.ismstr==12) THEN
148#include "vectorize.inc"
149 DO j=1,nnega
150 i = index(j)
151 x1(i)=sav(i,1)
152 y1(i)=sav(i,8)
153 z1(i)=sav(i,15)
154 x2(i)=sav(i,2)
155 y2(i)=sav(i,9)
156 z2(i)=sav(i,16)
157 x3(i)=sav(i,3)
158 y3(i)=sav(i,10)
159 z3(i)=sav(i,17)
160 x4(i)=sav(i,4)
161 y4(i)=sav(i,11)
162 z4(i)=sav(i,18)
163 x5(i)=sav(i,5)
164 y5(i)=sav(i,12)
165 z5(i)=sav(i,19)
166 x6(i)=sav(i,6)
167 y6(i)=sav(i,13)
168 z6(i)=sav(i,20)
169 x7(i)=sav(i,7)
170 y7(i)=sav(i,14)
171 z7(i)=sav(i,21)
172 x8(i)=zero
173 y8(i)=zero
174 z8(i)=zero
175 ENDDO
176 ELSE
177 DO j=1,nnega
178 i = index(j)
179 x1(i)=sav(i,1)
180 y1(i)=sav(i,2)
181 z1(i)=sav(i,3)
182 x2(i)=sav(i,4)
183 y2(i)=sav(i,5)
184 z2(i)=sav(i,6)
185 x3(i)=sav(i,7)
186 y3(i)=sav(i,8)
187 z3(i)=sav(i,9)
188 x4(i)=sav(i,10)
189 y4(i)=sav(i,11)
190 z4(i)=sav(i,12)
191 x5(i)=sav(i,13)
192 y5(i)=sav(i,14)
193 z5(i)=sav(i,15)
194 x6(i)=sav(i,16)
195 y6(i)=sav(i,17)
196 z6(i)=sav(i,18)
197 x7(i)=sav(i,19)
198 y7(i)=sav(i,20)
199 z7(i)=sav(i,21)
200 x8(i)=zero
201 y8(i)=zero
202 z8(i)=zero
203 ENDDO
204 END IF
205#include "vectorize.inc"
206 DO j=1,nnega
207 i = index(j)
208C
209 x17=x7(i)-x1(i)
210 x28=x8(i)-x2(i)
211 x35=x5(i)-x3(i)
212 x46=x6(i)-x4(i)
213 y17=y7(i)-y1(i)
214 y28=y8(i)-y2(i)
215 y35=y5(i)-y3(i)
216 y46=y6(i)-y4(i)
217 z17=z7(i)-z1(i)
218 z28=z8(i)-z2(i)
219 z35=z5(i)-z3(i)
220 z46=z6(i)-z4(i)
221C
222 jac4(i)=x17+x28-x35-x46
223 jac5(i)=y17+y28-y35-y46
224 jac6(i)=z17+z28-z35-z46
225 jac7(i)=x17+x46+x28+x35
226 jac8(i)=y17+y46+y28+y35
227 jac9(i)=z17+z46+z28+z35
228 jac1(i)=x17+x46-x28-x35
229 jac2(i)=y17+y46-y28-y35
230 jac3(i)=z17+z46-z28-z35
231C
232 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
233 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
234 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
235C
236 voldp(i)=one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
237 det(i)=voldp(i)
238 offg(i) = two
239 ENDDO
240 END IF
241C
242C Jacobian matrix inverse and shape functions derivatives Pij
243 DO i=1,nel
244 dett=one_over_64/det(i)
245 jaci1=dett*jac_59_68(i)
246 jaci4=dett*jac_67_49(i)
247 jaci7=dett*jac_48_57(i)
248 jaci2=dett*(-jac2(i)*jac9(i)+jac3(i)*jac8(i))
249 jaci5=dett*( jac1(i)*jac9(i)-jac3(i)*jac7(i))
250 jaci8=dett*(-jac1(i)*jac8(i)+jac2(i)*jac7(i))
251 jaci3=dett*( jac2(i)*jac6(i)-jac3(i)*jac5(i))
252 jaci6=dett*(-jac1(i)*jac6(i)+jac3(i)*jac4(i))
253 jaci9=dett*( jac1(i)*jac5(i)-jac2(i)*jac4(i))
254
255 jaci12=jaci1-jaci2
256 jaci45=jaci4-jaci5
257 jaci78=jaci7-jaci8
258 px2(i)= jaci12-jaci3
259 py2(i)= jaci45-jaci6
260 pz2(i)= jaci78-jaci9
261 px4(i)=-jaci12-jaci3
262 py4(i)=-jaci45-jaci6
263 pz4(i)=-jaci78-jaci9
264
265 jaci12=jaci1+jaci2
266 jaci45=jaci4+jaci5
267 jaci78=jaci7+jaci8
268 px1(i)=-jaci12-jaci3
269 py1(i)=-jaci45-jaci6
270 pz1(i)=-jaci78-jaci9
271 px3(i)=jaci12-jaci3
272 py3(i)=jaci45-jaci6
273 pz3(i)=jaci78-jaci9
274 ENDDO
275C
276 DO i=1,nel
277C 1 -1 1 -1 1 -1 1 -1
278 hx=one_over_8*(x1(i)-x2(i)+x3(i)-x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
279 hy=one_over_8*(y1(i)-y2(i)+y3(i)-y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
280 hz=one_over_8*(z1(i)-z2(i)+z3(i)-z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
281 px1h3(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
282 px2h3(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
283 px3h3(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
284 px4h3(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
285 END DO
286C 1 1 -1 -1 -1 -1 1 1
287 DO i=1,nel
288 hx=one_over_8*(x1(i)+x2(i)-x3(i)-x4(i)-x5(i)-x6(i)+x7(i)+x8(i))
289 hy=one_over_8*(y1(i)+y2(i)-y3(i)-y4(i)-y5(i)-y6(i)+y7(i)+y8(i))
290 hz=one_over_8*(z1(i)+z2(i)-z3(i)-z4(i)-z5(i)-z6(i)+z7(i)+z8(i))
291 px1h1(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
292 px2h1(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
293 px3h1(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
294 px4h1(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
295 END DO
296C 1 -1 -1 1 -1 1 1 -1
297 DO i=1,nel
298 hx=one_over_8*(x1(i)-x2(i)-x3(i)+x4(i)-x5(i)+x6(i)+x7(i)-x8(i))
299 hy=one_over_8*(y1(i)-y2(i)-y3(i)+y4(i)-y5(i)+y6(i)+y7(i)-y8(i))
300 hz=one_over_8*(z1(i)-z2(i)-z3(i)+z4(i)-z5(i)+z6(i)+z7(i)-z8(i))
301 px1h2(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
302 px2h2(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
303 px3h2(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
304 px4h2(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
305 END DO
306C -1 1 -1 1 1 -1 1 -1
307 DO i=1,nel
308 hx=one_over_8*(-x1(i)+x2(i)-x3(i)+x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
309 hy=one_over_8*(-y1(i)+y2(i)-y3(i)+y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
310 hz=one_over_8*(-z1(i)+z2(i)-z3(i)+z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
311 px1h4(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
312 px2h4(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
313 px3h4(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
314 px4h4(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
315 END DO
316 RETURN
317C
#define my_real
Definition cppsort.cpp:32
subroutine schkjabt3(off, det, ngl, offg, nnega, index, nel, ismstr, jlag)
Definition schkjabt3.F:40

◆ szderit3()

subroutine szderit3 ( off,
det,
integer, dimension(*) ngl,
integer nel,
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,
jac1,
jac2,
jac3,
jac4,
jac5,
jac6,
jac9,
integer, intent(in) jlag )

Definition at line 326 of file szderi3.F.

343C-----------------------------------------------
344C I m p l i c i t T y p e s
345C-----------------------------------------------
346#include "implicit_f.inc"
347#include "comlock.inc"
348C-----------------------------------------------
349C G l o b a l P a r a m e t e r s
350C-----------------------------------------------
351#include "mvsiz_p.inc"
352C-----------------------------------------------
353C D u m m y A r g u m e n t s
354C-----------------------------------------------
355 INTEGER, INTENT(IN) :: JLAG
356 INTEGER :: NEL,NGL(*)
357 double precision
358 . x1(*), x2(*), x3(*), x4(*), x5(*), x6(*), x7(*), x8(*),
359 . y1(*), y2(*), y3(*), y4(*), y5(*), y6(*), y7(*), y8(*),
360 . z1(*), z2(*), z3(*), z4(*), z5(*), z6(*), z7(*), z8(*)
361 my_real
362 . off(*),det(*),
363 . px1(*), px2(*), px3(*), px4(*),
364 . py1(*), py2(*), py3(*), py4(*),
365 . pz1(*), pz2(*), pz3(*), pz4(*),
366 . px1h1(*), px1h2(*), px1h3(*),px1h4(*),
367 . px2h1(*), px2h2(*), px2h3(*),px2h4(*),
368 . px3h1(*), px3h2(*), px3h3(*),px3h4(*),
369 . px4h1(*), px4h2(*), px4h3(*),px4h4(*),
370 . jac1(*),jac2(*),jac3(*),
371 . jac4(*),jac5(*),jac6(*),jac9(*)
372C-----------------------------------------------
373C L o c a l V a r i a b l e s
374C-----------------------------------------------
375 INTEGER :: I
376C 12
377 my_real
378 . dett , jac7(mvsiz), jac8(mvsiz) ,
379 . jaci1, jaci2, jaci3,
380 . jaci4, jaci5, jaci6,
381 . jaci7, jaci8, jaci9,
382 . x17 , x28 , x35 , x46,
383 . y17 , y28 , y35 , y46,
384 . z17 , z28 , z35 , z46,
385 . jac_59_68(mvsiz), jac_67_49(mvsiz), jac_48_57(mvsiz),
386 . jaci12, jaci45, jaci78,
387 . x_17_46 , x_28_35 ,
388 . y_17_46 , y_28_35 ,
389 . z_17_46 , z_28_35 ,
390 . hx,hy,hz
391C-----------------------------------------------
392C
393C Jacobian matrix
394 DO i=1,nel
395 x17=x7(i)-x1(i)
396 x28=x8(i)-x2(i)
397 x35=x5(i)-x3(i)
398 x46=x6(i)-x4(i)
399 y17=y7(i)-y1(i)
400 y28=y8(i)-y2(i)
401 y35=y5(i)-y3(i)
402 y46=y6(i)-y4(i)
403 z17=z7(i)-z1(i)
404 z28=z8(i)-z2(i)
405 z35=z5(i)-z3(i)
406 z46=z6(i)-z4(i)
407C
408 jac4(i)=x17+x28-x35-x46
409 jac5(i)=y17+y28-y35-y46
410 jac6(i)=z17+z28-z35-z46
411C
412 x_17_46=x17+x46
413 x_28_35=x28+x35
414 y_17_46=y17+y46
415 y_28_35=y28+y35
416 z_17_46=z17+z46
417 z_28_35=z28+z35
418C
419 jac7(i)=x_17_46+x_28_35
420 jac8(i)=y_17_46+y_28_35
421 jac9(i)=z_17_46+z_28_35
422 jac1(i)=x_17_46-x_28_35
423 jac2(i)=y_17_46-y_28_35
424 jac3(i)=z_17_46-z_28_35
425 ENDDO
426C
427 DO i=1,nel
428 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
429 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
430 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
431 END DO
432C
433 DO i=1,nel
434 det(i)=one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
435 END DO
436C
437 CALL schkjab3(
438 1 off, det, ngl, nel)
439C
440C Jacobian matrix inverse and shape functions derivatives Pij
441 DO i=1,nel
442 dett=one_over_64/det(i)
443 jaci1=dett*jac_59_68(i)
444 jaci4=dett*jac_67_49(i)
445 jaci7=dett*jac_48_57(i)
446 jaci2=dett*(-jac2(i)*jac9(i)+jac3(i)*jac8(i))
447 jaci5=dett*( jac1(i)*jac9(i)-jac3(i)*jac7(i))
448 jaci8=dett*(-jac1(i)*jac8(i)+jac2(i)*jac7(i))
449 jaci3=dett*( jac2(i)*jac6(i)-jac3(i)*jac5(i))
450 jaci6=dett*(-jac1(i)*jac6(i)+jac3(i)*jac4(i))
451 jaci9=dett*( jac1(i)*jac5(i)-jac2(i)*jac4(i))
452C
453 jaci12=jaci1-jaci2
454 jaci45=jaci4-jaci5
455 jaci78=jaci7-jaci8
456
457 px2(i)= jaci12-jaci3
458 py2(i)= jaci45-jaci6
459 pz2(i)= jaci78-jaci9
460 px4(i)=-jaci12-jaci3
461 py4(i)=-jaci45-jaci6
462 pz4(i)=-jaci78-jaci9
463
464 jaci12=jaci1+jaci2
465 jaci45=jaci4+jaci5
466 jaci78=jaci7+jaci8
467
468 px1(i)=-jaci12-jaci3
469 py1(i)=-jaci45-jaci6
470 pz1(i)=-jaci78-jaci9
471 px3(i)=jaci12-jaci3
472 py3(i)=jaci45-jaci6
473 pz3(i)=jaci78-jaci9
474 ENDDO
475C
476 DO i=1,nel
477C 1 -1 1 -1 1 -1 1 -1
478 hx=one_over_8*(x1(i)-x2(i)+x3(i)-x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
479 hy=one_over_8*(y1(i)-y2(i)+y3(i)-y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
480 hz=one_over_8*(z1(i)-z2(i)+z3(i)-z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
481 px1h3(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
482 px2h3(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
483 px3h3(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
484 px4h3(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
485 END DO
486C 1 1 -1 -1 -1 -1 1 1
487 DO i=1,nel
488 hx=one_over_8*(x1(i)+x2(i)-x3(i)-x4(i)-x5(i)-x6(i)+x7(i)+x8(i))
489 hy=one_over_8*(y1(i)+y2(i)-y3(i)-y4(i)-y5(i)-y6(i)+y7(i)+y8(i))
490 hz=one_over_8*(z1(i)+z2(i)-z3(i)-z4(i)-z5(i)-z6(i)+z7(i)+z8(i))
491 px1h1(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
492 px2h1(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
493 px3h1(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
494 px4h1(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
495 END DO
496C 1 -1 -1 1 -1 1 1 -1
497 DO i=1,nel
498 hx=one_over_8*(x1(i)-x2(i)-x3(i)+x4(i)-x5(i)+x6(i)+x7(i)-x8(i))
499 hy=one_over_8*(y1(i)-y2(i)-y3(i)+y4(i)-y5(i)+y6(i)+y7(i)-y8(i))
500 hz=one_over_8*(z1(i)-z2(i)-z3(i)+z4(i)-z5(i)+z6(i)+z7(i)-z8(i))
501 px1h2(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
502 px2h2(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
503 px3h2(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
504 px4h2(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
505 END DO
506C -1 1 -1 1 1 -1 1 -1
507 DO i=1,nel
508 hx=one_over_8*(-x1(i)+x2(i)-x3(i)+x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
509 hy=one_over_8*(-y1(i)+y2(i)-y3(i)+y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
510 hz=one_over_8*(-z1(i)+z2(i)-z3(i)+z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
511 px1h4(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
512 px2h4(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
513 px3h4(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
514 px4h4(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
515 END DO
516 RETURN
517C
subroutine schkjab3(off, det, ngl, nel)
Definition schkjab3.F:39

◆ szderito3()

subroutine szderito3 ( off,
det,
px1,
px2,
px3,
px4,
py1,
py2,
py3,
py4,
pz1,
pz2,
pz3,
pz4,
jac_i,
integer nel,
integer, intent(in) jlag )

Definition at line 524 of file szderi3.F.

530C-----------------------------------------------
531C I m p l i c i t T y p e s
532C-----------------------------------------------
533#include "implicit_f.inc"
534C-----------------------------------------------
535C G l o b a l P a r a m e t e r s
536C-----------------------------------------------
537#include "mvsiz_p.inc"
538C-----------------------------------------------
539C D u m m y A r g u m e n t s
540C-----------------------------------------------
541 INTEGER, INTENT(IN) :: JLAG
542 INTEGER :: NEL
543 my_real
544 . off(*),det(*),
545 . px1(*), px2(*), px3(*), px4(*),
546 . py1(*), py2(*), py3(*), py4(*),
547 . pz1(*), pz2(*), pz3(*), pz4(*),
548 . jac_i(10,mvsiz)
549C-----------------------------------------------
550C L o c a l V a r i a b l e s
551C-----------------------------------------------
552 INTEGER I
553C 12
554 my_real
555 . jaci1, jaci2, jaci3,
556 . jaci4, jaci5, jaci6,
557 . jaci7, jaci8, jaci9,
558 . jaci12, jaci45, jaci78
559C-----------------------------------------------
560C----JAC_I is calculated based on H8, permutation :[t,r,s]:[r,s,t]
561 DO i=1,nel
562 jaci1=jac_i(3,i)
563 jaci4=jac_i(6,i)
564 jaci7=jac_i(9,i)
565 jaci2=jac_i(1,i)
566 jaci5=jac_i(4,i)
567 jaci8=jac_i(7,i)
568 jaci3=jac_i(2,i)
569 jaci6=jac_i(5,i)
570 jaci9=jac_i(8,i)
571 det(i) =jac_i(10,i)
572C
573 jaci12=jaci1-jaci2
574 jaci45=jaci4-jaci5
575 jaci78=jaci7-jaci8
576 px2(i)= jaci12-jaci3
577 py2(i)= jaci45-jaci6
578 pz2(i)= jaci78-jaci9
579 px4(i)=-jaci12-jaci3
580 py4(i)=-jaci45-jaci6
581 pz4(i)=-jaci78-jaci9
582C
583 jaci12=jaci1+jaci2
584 jaci45=jaci4+jaci5
585 jaci78=jaci7+jaci8
586 px1(i)=-jaci12-jaci3
587 py1(i)=-jaci45-jaci6
588 pz1(i)=-jaci78-jaci9
589 px3(i)= jaci12-jaci3
590 py3(i)= jaci45-jaci6
591 pz3(i)= jaci78-jaci9
592 ENDDO
593C
594 RETURN
595C