OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
scderi3.F File Reference
#include "implicit_f.inc"
#include "comlock.inc"
#include "mvsiz_p.inc"
#include "com06_c.inc"
#include "units_c.inc"
#include "scr07_c.inc"
#include "scr17_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 scderi3 (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, jac1, jac2, jac3, jac4, jac5, jac6, rx0, ry0, sx0, sy0, vzl, volg, sav, offg, nel, ismstr)

Function/Subroutine Documentation

◆ scderi3()

subroutine scderi3 ( 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,
jac1,
jac2,
jac3,
jac4,
jac5,
jac6,
rx0,
ry0,
sx0,
sy0,
vzl,
volg,
double precision, dimension(nel,21) sav,
offg,
integer nel,
integer, intent(in) ismstr )

Definition at line 32 of file scderi3.F.

51C-----------------------------------------------
52C M o d u l e s
53C-----------------------------------------------
54 USE message_mod
55C-----------------------------------------------
56C I m p l i c i t T y p e s
57C-----------------------------------------------
58#include "implicit_f.inc"
59#include "comlock.inc"
60C-----------------------------------------------
61C G l o b a l P a r a m e t e r s
62C-----------------------------------------------
63#include "mvsiz_p.inc"
64C-----------------------------------------------
65C C o m m o n B l o c k s
66C-----------------------------------------------
67#include "com06_c.inc"
68#include "units_c.inc"
69#include "scr07_c.inc"
70#include "scr17_c.inc"
71#include "scr18_c.inc"
72C-----------------------------------------------
73C D u m m y A r g u m e n t s
74C-----------------------------------------------
75 INTEGER, INTENT(IN) :: ISMSTR
76 INTEGER :: NEL,NGL(*)
77C
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 . jac1(*),jac2(*),jac3(*),
94 . jac4(*),jac5(*),jac6(*),
95 . rx0(*),ry0(*),sx0(*),sy0(*),vzl(*),volg(*),offg(*)
96C-----------------------------------------------
97C L o c a l V a r i a b l e s
98C-----------------------------------------------
99 INTEGER :: I,J,NNEGA,INDEX(MVSIZ)
100C 12
101 my_real
102 . dett , jac7(mvsiz), jac8(mvsiz) ,jac9(mvsiz),
103 . jaci1, jaci2, jaci3,
104 . jaci4, jaci5, jaci6,
105 . jaci7, jaci8, jaci9,
106 . x17(mvsiz) , x28(mvsiz) , x35(mvsiz) , x46(mvsiz),
107 . y17(mvsiz) , y28(mvsiz) , y35(mvsiz) , y46(mvsiz),
108 . z17(mvsiz) , z28(mvsiz) , z35(mvsiz) , z46(mvsiz),
109 . jac_59_68(mvsiz), jac_67_49(mvsiz), jac_48_57(mvsiz),jac_19_37(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 ,area(mvsiz) ,
114 . hx,hy,hz,a_i,
115 . h1x(mvsiz),h1y(mvsiz),h2x(mvsiz),h2y(mvsiz),
116 . h1z(mvsiz),h2z(mvsiz),icor
117C=======================================================================
118C
119 nnega=0
120 DO i=1,nel
121 x17(i)=x7(i)-x1(i)
122 x28(i)=x8(i)-x2(i)
123 x35(i)=x5(i)-x3(i)
124 x46(i)=x6(i)-x4(i)
125 y17(i)=y7(i)-y1(i)
126 y28(i)=y8(i)-y2(i)
127 y35(i)=y5(i)-y3(i)
128 y46(i)=y6(i)-y4(i)
129 z17(i)=z7(i)-z1(i)
130 z28(i)=z8(i)-z2(i)
131 z35(i)=z5(i)-z3(i)
132 z46(i)=z6(i)-z4(i)
133 END DO
134C
135C Jacobian matrix
136 DO i=1,nel
137C -------ri.xi-----------
138 jac3(i)=x17(i)+x28(i)-x35(i)-x46(i)
139 jac1(i)=y17(i)+y28(i)-y35(i)-y46(i)
140 jac2(i)=z17(i)+z28(i)-z35(i)-z46(i)
141C-------
142 x_17_46=x17(i)+x46(i)
143 x_28_35=x28(i)+x35(i)
144 y_17_46=y17(i)+y46(i)
145 y_28_35=y28(i)+y35(i)
146 z_17_46=z17(i)+z46(i)
147 z_28_35=z28(i)+z35(i)
148C -------si.xi-----------
149 jac6(i)=x_17_46+x_28_35
150 jac4(i)=y_17_46+y_28_35
151 jac5(i)=z_17_46+z_28_35
152C -------ti.xi-----------
153 jac9(i)=x_17_46-x_28_35
154 jac7(i)=y_17_46-y_28_35
155 jac8(i)=z_17_46-z_28_35
156 ENDDO
157C
158 DO i=1,nel
159 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
160 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
161 jac_19_37(i)=jac1(i)*jac9(i)-jac3(i)*jac7(i)
162 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
163 ENDDO
164C
165 DO i=1,nel
166 area(i) = one_over_64*jac_19_37(i)
167 det(i)=one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
168 ENDDO
169C
170 IF(idtmin(1)==1)THEN
171 icor = 0
172 DO i=1,nel
173 IF(off(i) ==zero)THEN
174 area(i) = one
175 det(i)=one
176 ELSEIF((det(i)<=volmin).OR.(det(i)<=zero)
177 . .OR.(area(i)<=zero))THEN
178 icor = 1
179 ENDIF
180 ENDDO
181 IF (icor>0) THEN
182 DO i=1,nel
183 IF(off(i)/=zero)THEN
184 IF(det(i)<=volmin)THEN
185 area(i) = one
186 det(i)=one
187 off(i)=zero
188#include "lockon.inc"
189 WRITE(istdo,2000) ngl(i)
190 WRITE(iout ,2000) ngl(i)
191#include "lockoff.inc"
192 ELSEIF(det(i)<=zero.OR.(area(i)<=zero))THEN
193 CALL ancmsg(msgid=166,anmode=aninfo,
194 . i1=ngl(i))
195 mstop = 1
196 ENDIF
197 ENDIF
198 ENDDO
199 ENDIF
200 ELSEIF(idtmin(1)==2)THEN
201 icor = 0
202 DO i=1,nel
203 IF(off(i) ==zero)THEN
204 area(i) = one
205 det(i)=one
206 ELSEIF((det(i)<=volmin).OR.(det(i)<=zero)
207 . .OR.(area(i)<=zero))THEN
208 icor=1
209 ENDIF
210 ENDDO
211 IF (icor>0) THEN
212 DO i=1,nel
213 IF((off(i)/=0.).AND.
214 . (det(i)<=volmin.OR.det(i)<=zero)
215 . .OR.(area(i)<=zero))THEN
216 area(i) = one
217 det(i)=one
218 off(i)=zero
219#include "lockon.inc"
220 WRITE(istdo,2000) ngl(i)
221 WRITE(iout ,2000) ngl(i)
222#include "lockoff.inc"
223 ENDIF
224 ENDDO
225 ENDIF
226 ELSEIF (ismstr /=4 ) THEN
227 icor = 0
228 DO i=1,nel
229 IF(off(i) ==zero)THEN
230 det(i)=one
231 ELSEIF((det(i)<=volmin).OR.(area(i)<=zero))THEN
232 icor = 1
233 ENDIF
234 ENDDO
235 IF (icor>0) THEN
236 DO i=1,nel
237 IF(off(i) == zero)THEN
238 det(i)=one
239 area(i) = one
240 ELSEIF(offg(i) > one)THEN
241
242 ELSEIF((det(i)<=volmin).OR.(area(i)<=zero))THEN
243 nnega=nnega+1
244 index(nnega)=i
245#include "lockon.inc"
246 WRITE(istdo,3000) ngl(i)
247 WRITE(iout ,3000) ngl(i)
248#include "lockoff.inc"
249 ENDIF
250 ENDDO
251 IF (ineg_v==0) THEN
252 CALL ancmsg(msgid=280,anmode=aninfo)
253 mstop = 1
254 ENDIF
255 END IF !(ICOR>0) THEN
256 ELSE
257 icor = 0
258 DO i=1,nel
259 IF(off(i) ==zero)THEN
260 det(i)=one
261 area(i) = one
262 ELSEIF(det(i)<=zero)THEN
263 icor=1
264 ENDIF
265 ENDDO
266 IF (icor>0) THEN
267 DO i=1,nel
268 IF(off(i)/=zero)THEN
269 IF(det(i)<=zero)THEN
270 CALL ancmsg(msgid=166,anmode=aninfo,
271 . i1=ngl(i))
272 mstop = 1
273 ENDIF
274 ENDIF
275 ENDDO
276 ENDIF
277 ENDIF
278C
279 IF (nnega>0) THEN
280#include "vectorize.inc"
281 DO j=1,nnega
282 i = index(j)
283 x1(i)=sav(i,1)
284 y1(i)=sav(i,2)
285 z1(i)=sav(i,3)
286 x2(i)=sav(i,4)
287 y2(i)=sav(i,5)
288 z2(i)=sav(i,6)
289 x3(i)=sav(i,7)
290 y3(i)=sav(i,8)
291 z3(i)=sav(i,9)
292 x4(i)=sav(i,10)
293 y4(i)=sav(i,11)
294 z4(i)=sav(i,12)
295 x5(i)=sav(i,13)
296 y5(i)=sav(i,14)
297 z5(i)=sav(i,15)
298 x6(i)=sav(i,16)
299 y6(i)=sav(i,17)
300 z6(i)=sav(i,18)
301 x7(i)=sav(i,19)
302 y7(i)=sav(i,20)
303 z7(i)=sav(i,21)
304 x8(i)=zero
305 y8(i)=zero
306 z8(i)=zero
307C
308C---------
309 jac3(i)=x17(i)+x28(i)-x35(i)-x46(i)
310 jac1(i)=y17(i)+y28(i)-y35(i)-y46(i)
311 jac2(i)=z17(i)+z28(i)-z35(i)-z46(i)
312C-------
313 x_17_46=x17(i)+x46(i)
314 x_28_35=x28(i)+x35(i)
315 y_17_46=y17(i)+y46(i)
316 y_28_35=y28(i)+y35(i)
317 z_17_46=z17(i)+z46(i)
318 z_28_35=z28(i)+z35(i)
319C----------
320 jac6(i)=x_17_46+x_28_35
321 jac4(i)=y_17_46+y_28_35
322 jac5(i)=z_17_46+z_28_35
323 jac9(i)=x_17_46-x_28_35
324 jac7(i)=y_17_46-y_28_35
325 jac8(i)=z_17_46-z_28_35
326C
327 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
328 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
329 jac_19_37(i)=jac1(i)*jac9(i)-jac3(i)*jac7(i)
330 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
331C
332 area(i) = one_over_64*jac_19_37(i)
333 det(i)=one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
334 offg(i) = two
335 ENDDO
336 END IF
337C
338C Jacobian matrix inverse
339 DO i=1,nel
340 dett=one_over_64/det(i)
341 jaci1=dett*jac_59_68(i)
342 jaci4=dett*jac_67_49(i)
343 jaci7=dett*jac_48_57(i)
344 jaci2=dett*(-jac2(i)*jac9(i)+jac3(i)*jac8(i))
345 jaci5=dett*jac_19_37(i)
346 jaci8=dett*(-jac1(i)*jac8(i)+jac2(i)*jac7(i))
347 jaci3=dett*( jac2(i)*jac6(i)-jac3(i)*jac5(i))
348 jaci6=dett*(-jac1(i)*jac6(i)+jac3(i)*jac4(i))
349 jaci9=dett*( jac1(i)*jac5(i)-jac2(i)*jac4(i))
350C
351 jaci12=jaci1-jaci2
352 jaci45=jaci4-jaci5
353 jaci78=jaci7-jaci8
354C
355 py3(i)= jaci12+jaci3
356 pz3(i)= jaci45+jaci6
357 px3(i)= jaci78+jaci9
358 py4(i)= jaci12-jaci3
359 pz4(i)= jaci45-jaci6
360 px4(i)= jaci78-jaci9
361
362 jaci12=jaci1+jaci2
363 jaci45=jaci4+jaci5
364 jaci78=jaci7+jaci8
365
366 py1(i)=-jaci12-jaci3
367 pz1(i)=-jaci45-jaci6
368 px1(i)=-jaci78-jaci9
369 py2(i)=-jaci12+jaci3
370 pz2(i)=-jaci45+jaci6
371 px2(i)=-jaci78+jaci9
372 ENDDO
373C
374 DO i=1,nel
375C h3
376C 1 -1 1 -1 1 -1 1 -1
377 hx=one_over_8*(x1(i)-x2(i)+x3(i)-x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
378 hy=one_over_8*(y1(i)-y2(i)+y3(i)-y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
379 hz=one_over_8*(z1(i)-z2(i)+z3(i)-z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
380 px1h3(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
381 px2h3(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
382 px3h3(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
383 px4h3(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
384 END DO
385C h1
386C 1 1 -1 -1 -1 -1 1 1
387 DO i=1,nel
388 h1x(i)=x1(i)+x2(i)-x3(i)-x4(i)-x5(i)-x6(i)+x7(i)+x8(i)
389 h1y(i)=y1(i)+y2(i)-y3(i)-y4(i)-y5(i)-y6(i)+y7(i)+y8(i)
390 h1z(i)=z1(i)+z2(i)-z3(i)-z4(i)-z5(i)-z6(i)+z7(i)+z8(i)
391 hx=one_over_8*h1x(i)
392 hy=one_over_8*h1y(i)
393 hz=one_over_8*h1z(i)
394 px1h1(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
395 px2h1(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
396 px3h1(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
397 px4h1(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
398 END DO
399C h2
400C 1 -1 -1 1 -1 1 1 -1
401 DO i=1,nel
402 h2x(i)=x1(i)-x2(i)-x3(i)+x4(i)-x5(i)+x6(i)+x7(i)-x8(i)
403 h2y(i)=y1(i)-y2(i)-y3(i)+y4(i)-y5(i)+y6(i)+y7(i)-y8(i)
404 h2z(i)=z1(i)-z2(i)-z3(i)+z4(i)-z5(i)+z6(i)+z7(i)-z8(i)
405 hx=one_over_8*h2x(i)
406 hy=one_over_8*h2y(i)
407 hz=one_over_8*h2z(i)
408 px1h2(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
409 px2h2(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
410 px3h2(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
411 px4h2(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
412 END DO
413C h4
414C -1 1 -1 1 1 -1 1 -1
415 DO i=1,nel
416 hx=one_over_8*(-x1(i)+x2(i)-x3(i)+x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
417 hy=one_over_8*(-y1(i)+y2(i)-y3(i)+y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
418 hz=one_over_8*(-z1(i)+z2(i)-z3(i)+z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
419 px1h4(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
420 px2h4(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
421 px3h4(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
422 px4h4(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
423 END DO
424C
425 DO i=1,nel
426 a_i = one_over_8/area(i)
427 rx0(i) = jac9(i)*a_i
428 ry0(i) = jac7(i)*a_i
429 sx0(i) = jac3(i)*a_i
430 sy0(i) = jac1(i)*a_i
431 vzl(i) = one_over_64*jac5(i)*(
432 . jac9(i)*h1y(i)+jac1(i)*h2x(i)-jac3(i)*h2y(i)-jac7(i)*h1x(i))
433 . + one_over_64*jac4(i)*(
434 . jac3(i)*h2z(i)+jac8(i)*h1x(i)-jac9(i)*h1z(i)-jac2(i)*h2x(i))
435 . + one_over_64*jac6(i)*(
436 . jac7(i)*h1z(i)+jac2(i)*h2y(i)-jac1(i)*h2z(i)-jac8(i)*h1y(i))
437C
438 volg(i)= det(i)
439 ENDDO
440C-----------
441 RETURN
442 2000 FORMAT(/' ZERO OR NEGATIVE VOLUME : DELETE 3D-ELEMENT NB',i10/)
443 3000 FORMAT(/' ZERO OR NEGATIVE VOLUME : 3D-ELEMENT NB:',i10/,
444 + ' SOLID-SHELL ELEMENT IS SWITCHED TO SMALL STRAIN OPTION'/)
445C-----------
#define my_real
Definition cppsort.cpp:32
subroutine area(d1, x, x2, y, y2, eint, stif0)
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