OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
sderi3.F File Reference
#include "implicit_f.inc"
#include "mvsiz_p.inc"
#include "param_c.inc"
#include "com04_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine sderi3 (vol, veul, geo, igeo, xd1, xd2, xd3, xd4, xd5, xd6, xd7, xd8, yd1, yd2, yd3, yd4, yd5, yd6, yd7, yd8, zd1, zd2, zd3, zd4, zd5, zd6, zd7, zd8, jac1, jac2, jac3, jac4, jac5, jac6, ngl, ngeo, px1, px2, px3, px4, py1, py2, py3, py4, pz1, pz2, pz3, pz4, det, voldp, nel, jeul, nxref, imulti_fvm)
subroutine sjac_i (x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, jac_i, nel)

Function/Subroutine Documentation

◆ sderi3()

subroutine sderi3 ( vol,
veul,
geo,
integer, dimension(npropgi,numgeo) igeo,
double precision, dimension(mvsiz) xd1,
double precision, dimension(mvsiz) xd2,
double precision, dimension(mvsiz) xd3,
double precision, dimension(mvsiz) xd4,
double precision, dimension(mvsiz) xd5,
double precision, dimension(mvsiz) xd6,
double precision, dimension(mvsiz) xd7,
double precision, dimension(mvsiz) xd8,
double precision, dimension(mvsiz) yd1,
double precision, dimension(mvsiz) yd2,
double precision, dimension(mvsiz) yd3,
double precision, dimension(mvsiz) yd4,
double precision, dimension(mvsiz) yd5,
double precision, dimension(mvsiz) yd6,
double precision, dimension(mvsiz) yd7,
double precision, dimension(mvsiz) yd8,
double precision, dimension(mvsiz) zd1,
double precision, dimension(mvsiz) zd2,
double precision, dimension(mvsiz) zd3,
double precision, dimension(mvsiz) zd4,
double precision, dimension(mvsiz) zd5,
double precision, dimension(mvsiz) zd6,
double precision, dimension(mvsiz) zd7,
double precision, dimension(mvsiz) zd8,
jac1,
jac2,
jac3,
jac4,
jac5,
jac6,
integer, dimension(nel) ngl,
integer, dimension(nel) ngeo,
px1,
px2,
px3,
px4,
py1,
py2,
py3,
py4,
pz1,
pz2,
pz3,
pz4,
det,
double precision, dimension(nel) voldp,
integer nel,
integer jeul,
integer nxref,
integer imulti_fvm )

Definition at line 35 of file sderi3.F.

44C-----------------------------------------------
45C M o d u l e s
46C-----------------------------------------------
47 USE message_mod
48C----------------------------------------------------------
49C I m p l i c i t T y p e s
50C-----------------------------------------------
51#include "implicit_f.inc"
52C-----------------------------------------------
53C G l o b a l P a r a m e t e r s
54C-----------------------------------------------
55#include "mvsiz_p.inc"
56C-----------------------------------------------
57C C o m m o n B l o c k s
58C-----------------------------------------------
59#include "param_c.inc"
60#include "com04_c.inc"
61C-----------------------------------------------
62C D u m m y A r g u m e n t s
63C-----------------------------------------------
64 INTEGER :: IGEO(NPROPGI,NUMGEO), NGL(NEL), NGEO(NEL)
65 INTEGER :: NEL, JEUL, NXREF, IMULTI_FVM
66C
68 . vol(nel), veul(lveul,*) , geo(npropg,numgeo),
69 . jac1(nel), jac2(nel), jac3(nel), jac4(nel), jac5(nel), jac6(nel),
70 . jac12(nel), jac45(nel), jac78(nel),
71 . px1(nel), px2(nel), px3(nel), px4(nel),
72 . py1(nel), py2(nel), py3(nel), py4(nel),
73 . pz1(nel), pz2(nel), pz3(nel), pz4(nel), det(nel)
74 double precision
75 . xd1(mvsiz), xd2(mvsiz), xd3(mvsiz), xd4(mvsiz),
76 . xd5(mvsiz), xd6(mvsiz), xd7(mvsiz), xd8(mvsiz),
77 . yd1(mvsiz), yd2(mvsiz), yd3(mvsiz), yd4(mvsiz),
78 . yd5(mvsiz), yd6(mvsiz), yd7(mvsiz), yd8(mvsiz),
79 . zd1(mvsiz), zd2(mvsiz), zd3(mvsiz), zd4(mvsiz),
80 . zd5(mvsiz), zd6(mvsiz), zd7(mvsiz), zd8(mvsiz),
81 . voldp(nel)
82C-----------------------------------------------
83C L o c a l V a r i a b l e s
84C-----------------------------------------------
85 INTEGER :: I
86C
87 my_real :: jac7(mvsiz), jac8(mvsiz) , jac9(mvsiz)
89 . x_17_46, x_28_35,
90 . y_17_46, y_28_35,
91 . z_17_46, z_28_35
93 . dett(mvsiz),
94 . jaci1(mvsiz), jaci2(mvsiz), jaci3(mvsiz),
95 . jaci4(mvsiz), jaci5(mvsiz), jaci6(mvsiz),
96 . jaci7(mvsiz), jaci8(mvsiz), jaci9(mvsiz),
97 . jac_59_68(mvsiz), jac_67_49(mvsiz), jac_48_57(mvsiz)
98 double precision
99 . x17, x28, x35, x46,
100 . y17, y28, y35, y46,
101 . z17, z28, z35, z46
102C-----------------------------------------------
103C S o u r c e L i n e s
104C-----------------------------------------------
105C
106C Jacobian matrix
107 DO i=1,nel
108 x17=xd7(i)-xd1(i)
109 x28=xd8(i)-xd2(i)
110 x35=xd5(i)-xd3(i)
111 x46=xd6(i)-xd4(i)
112C
113 y17=yd7(i)-yd1(i)
114 y28=yd8(i)-yd2(i)
115 y35=yd5(i)-yd3(i)
116 y46=yd6(i)-yd4(i)
117C
118 z17=zd7(i)-zd1(i)
119 z28=zd8(i)-zd2(i)
120 z35=zd5(i)-zd3(i)
121 z46=zd6(i)-zd4(i)
122C
123 jac1(i)=x17+x28-x35-x46
124 jac2(i)=y17+y28-y35-y46
125 jac3(i)=z17+z28-z35-z46
126C
127 x_17_46=x17+x46
128 x_28_35=x28+x35
129 y_17_46=y17+y46
130 y_28_35=y28+y35
131 z_17_46=z17+z46
132 z_28_35=z28+z35
133C
134 jac4(i)=x_17_46+x_28_35
135 jac5(i)=y_17_46+y_28_35
136 jac6(i)=z_17_46+z_28_35
137C
138 jac7(i)=x_17_46-x_28_35
139 jac8(i)=y_17_46-y_28_35
140 jac9(i)=z_17_46-z_28_35
141 ENDDO
142C
143 DO i=1,nel
144 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
145 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
146 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
147 ENDDO
148C
149 DO i=1,nel
150 voldp(i)=one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
151 det(i)=voldp(i)
152 vol(i)=det(i)
153 ENDDO
154C
155 IF(jeul * (1 - imulti_fvm) /= 0)THEN
156 DO i=1,nel
157 veul(32,i) = vol(i)
158 ENDDO
159 ENDIF
160C
161 DO i=1,nel
162 IF (det(i) > zero) cycle
163 IF (igeo(11,ngeo(i)) /= 0 .AND. igeo(11,ngeo(i)) /= 43) THEN
164 CALL ancmsg(msgid=245,
165 . msgtype=msgerror,
166 . anmode=aninfo,
167 . i1=ngl(i))
168 ELSE
169 CALL ancmsg(msgid=635,
170 . msgtype=msgwarning,
171 . anmode=aninfo,
172 . i1=ngl(i))
173 ENDIF
174 ENDDO
175C
176 IF( jeul==0 .AND. nxref==0) RETURN
177C
178 DO i=1,nel
179 dett(i)=one_over_64/max(det(i),em20)
180 ENDDO
181C
182C Jacobian matrix inverse
183 DO i=1,nel
184 jaci1(i)=dett(i)*jac_59_68(i)
185 jaci4(i)=dett(i)*jac_67_49(i)
186 jaci7(i)=dett(i)*jac_48_57(i)
187 jaci2(i)=dett(i)*(-jac2(i)*jac9(i)+jac3(i)*jac8(i))
188 jaci5(i)=dett(i)*( jac1(i)*jac9(i)-jac3(i)*jac7(i))
189 jaci8(i)=dett(i)*(-jac1(i)*jac8(i)+jac2(i)*jac7(i))
190 jaci3(i)=dett(i)*( jac2(i)*jac6(i)-jac3(i)*jac5(i))
191 jaci6(i)=dett(i)*(-jac1(i)*jac6(i)+jac3(i)*jac4(i))
192 jaci9(i)=dett(i)*( jac1(i)*jac5(i)-jac2(i)*jac4(i))
193 ENDDO
194C
195C Shape functions partial derivatives
196 DO i=1,nel
197 jac12(i)=jaci1(i)-jaci2(i)
198 jac45(i)=jaci4(i)-jaci5(i)
199 jac78(i)=jaci7(i)-jaci8(i)
200 ENDDO
201C
202 DO i=1,nel
203 px3(i)= jac12(i)+jaci3(i)
204 py3(i)= jac45(i)+jaci6(i)
205 pz3(i)= jac78(i)+jaci9(i)
206 px4(i)= jac12(i)-jaci3(i)
207 py4(i)= jac45(i)-jaci6(i)
208 pz4(i)= jac78(i)-jaci9(i)
209 ENDDO
210C
211 DO i=1,nel
212 jac12(i)=jaci1(i)+jaci2(i)
213 jac45(i)=jaci4(i)+jaci5(i)
214 jac78(i)=jaci7(i)+jaci8(i)
215 ENDDO
216C
217 DO i=1,nel
218 px1(i)=-jac12(i)-jaci3(i)
219 py1(i)=-jac45(i)-jaci6(i)
220 pz1(i)=-jac78(i)-jaci9(i)
221 px2(i)=-jac12(i)+jaci3(i)
222 py2(i)=-jac45(i)+jaci6(i)
223 pz2(i)=-jac78(i)+jaci9(i)
224 ENDDO
225
226 IF(jeul * (1 - imulti_fvm) /= 0)THEN
227 DO i=1,nel
228 veul(1,i) = px1(i)
229 veul(2,i) = px2(i)
230 veul(3,i) = px3(i)
231 veul(4,i) = px4(i)
232 veul(5,i) = py1(i)
233 veul(6,i) = py2(i)
234 veul(7,i) = py3(i)
235 veul(8,i) = py4(i)
236 veul(9,i) = pz1(i)
237 veul(10,i)= pz2(i)
238 veul(11,i)= pz3(i)
239 veul(12,i)= pz4(i)
240 END DO
241C
242 IF (igeo(11,ngeo(1)) == 15) THEN
243 DO i=1,nel
244 vol(i)=vol(i)*geo(1,ngeo(i))
245 END DO
246 END IF
247 END IF
248C
249 RETURN
#define my_real
Definition cppsort.cpp:32
#define max(a, b)
Definition macros.h:21
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

◆ sjac_i()

subroutine sjac_i ( x1,
x2,
x3,
x4,
x5,
x6,
x7,
x8,
y1,
y2,
y3,
y4,
y5,
y6,
y7,
y8,
z1,
z2,
z3,
z4,
z5,
z6,
z7,
z8,
jac_i,
integer nel )

Definition at line 259 of file sderi3.F.

264C-----------------------------------------------
265C M o d u l e s
266C-----------------------------------------------
267 USE message_mod
268C-----------------------------------------------
269C I m p l i c i t T y p e s
270C-----------------------------------------------
271#include "implicit_f.inc"
272C-----------------------------------------------
273C G l o b a l P a r a m e t e r s
274C-----------------------------------------------
275#include "mvsiz_p.inc"
276C-----------------------------------------------
277C D u m m y A r g u m e n t s
278C-----------------------------------------------
279 INTEGER :: NEL
280C
281 my_real
282 . x1(*), x2(*), x3(*), x4(*), x5(*), x6(*),
283 . x7(*), x8(*), y1(*), y2(*), y3(*), y4(*), y5(*), y6(*), y7(*),
284 . y8(*), z1(*), z2(*), z3(*), z4(*), z5(*), z6(*), z7(*), z8(*),
285 . jac_i(10,mvsiz)
286C-----------------------------------------------
287C L o c a l V a r i a b l e s
288C-----------------------------------------------
289 INTEGER :: I
290C
291 double precision
292 . xd1(mvsiz), xd2(mvsiz), xd3(mvsiz), xd4(mvsiz),
293 . xd5(mvsiz), xd6(mvsiz), xd7(mvsiz), xd8(mvsiz),
294 . yd1(mvsiz), yd2(mvsiz), yd3(mvsiz), yd4(mvsiz),
295 . yd5(mvsiz), yd6(mvsiz), yd7(mvsiz), yd8(mvsiz),
296 . zd1(mvsiz), zd2(mvsiz), zd3(mvsiz), zd4(mvsiz),
297 . zd5(mvsiz), zd6(mvsiz), zd7(mvsiz), zd8(mvsiz)
298C
299 my_real
300 . jac1(mvsiz), jac2(mvsiz), jac3(mvsiz),
301 . jac4(mvsiz), jac5(mvsiz), jac6(mvsiz),
302 . jac7(mvsiz), jac8(mvsiz), jac9(mvsiz),
303 . jac_59_68(mvsiz), jac_67_49(mvsiz), jac_48_57(mvsiz)
304C
305 my_real
306 . det(mvsiz), dett(mvsiz),
307 . x17(mvsiz), x28(mvsiz), x35(mvsiz), x46(mvsiz),
308 . y17(mvsiz), y28(mvsiz), y35(mvsiz), y46(mvsiz),
309 . z17(mvsiz), z28(mvsiz), z35(mvsiz), z46(mvsiz)
310C=======================================================================
311 DO i=1,nel
312 xd1(i)=x1(i)
313 xd2(i)=x2(i)
314 xd3(i)=x3(i)
315 xd4(i)=x4(i)
316 xd5(i)=x5(i)
317 xd6(i)=x6(i)
318 xd7(i)=x7(i)
319 xd8(i)=x8(i)
320 yd1(i)=y1(i)
321 yd2(i)=y2(i)
322 yd3(i)=y3(i)
323 yd4(i)=y4(i)
324 yd5(i)=y5(i)
325 yd6(i)=y6(i)
326 yd7(i)=y7(i)
327 yd8(i)=y8(i)
328 zd1(i)=z1(i)
329 zd2(i)=z2(i)
330 zd3(i)=z3(i)
331 zd4(i)=z4(i)
332 zd5(i)=z5(i)
333 zd6(i)=z6(i)
334 zd7(i)=z7(i)
335 zd8(i)=z8(i)
336 ENDDO
337
338 DO i=1,nel
339 x17(i)=xd7(i)-xd1(i)
340 x28(i)=xd8(i)-xd2(i)
341 x35(i)=xd5(i)-xd3(i)
342 x46(i)=xd6(i)-xd4(i)
343C
344 y17(i)=yd7(i)-yd1(i)
345 y28(i)=yd8(i)-yd2(i)
346 y35(i)=yd5(i)-yd3(i)
347 y46(i)=yd6(i)-yd4(i)
348C
349 z17(i)=zd7(i)-zd1(i)
350 z28(i)=zd8(i)-zd2(i)
351 z35(i)=zd5(i)-zd3(i)
352 z46(i)=zd6(i)-zd4(i)
353 ENDDO
354C
355C Jacobian Matrix
356 DO i=1,nel
357 jac1(i)=x17(i)+x28(i)-x35(i)-x46(i)
358 jac2(i)=y17(i)+y28(i)-y35(i)-y46(i)
359 jac3(i)=z17(i)+z28(i)-z35(i)-z46(i)
360C
361 jac4(i)=x17(i)+x46(i)+x28(i)+x35(i)
362 jac5(i)=y17(i)+y46(i)+y28(i)+y35(i)
363 jac6(i)=z17(i)+z46(i)+z28(i)+z35(i)
364C
365 jac7(i)=x17(i)+x46(i)-x28(i)-x35(i)
366 jac8(i)=y17(i)+y46(i)-y28(i)-y35(i)
367 jac9(i)=z17(i)+z46(i)-z28(i)-z35(i)
368 ENDDO
369C
370 DO i=1,nel
371 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
372 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
373 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
374 ENDDO
375C
376 DO i=1,nel
377 det(i) =one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
378 dett(i)=one_over_64/max(det(i),em20)
379 ENDDO
380C
381C Jacobian matrix inverse
382 DO i=1,nel
383 jac_i(1,i)=dett(i)*jac_59_68(i)
384 jac_i(4,i)=dett(i)*jac_67_49(i)
385 jac_i(7,i)=dett(i)*jac_48_57(i)
386 jac_i(2,i)=dett(i)*(-jac2(i)*jac9(i)+jac3(i)*jac8(i))
387 jac_i(5,i)=dett(i)*( jac1(i)*jac9(i)-jac3(i)*jac7(i))
388 jac_i(8,i)=dett(i)*(-jac1(i)*jac8(i)+jac2(i)*jac7(i))
389 jac_i(3,i)=dett(i)*( jac2(i)*jac6(i)-jac3(i)*jac5(i))
390 jac_i(6,i)=dett(i)*(-jac1(i)*jac6(i)+jac3(i)*jac4(i))
391 jac_i(9,i)=dett(i)*( jac1(i)*jac5(i)-jac2(i)*jac4(i))
392 ENDDO
393C
394 DO i=1,nel
395 jac_i(10,i) = det(i)
396 ENDDO
397C
398 RETURN