OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
sderit3.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 sderit3 (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, px2h1, px2h2, px2h3, px3h1, px3h2, px3h3, px4h1, px4h2, px4h3, jac1, jac2, jac3, jac4, jac5, jac6, nel, jhbe)

Function/Subroutine Documentation

◆ sderit3()

subroutine sderit3 ( off,
det,
integer, dimension(*), intent(in) 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,
px2h1,
px2h2,
px2h3,
px3h1,
px3h2,
px3h3,
px4h1,
px4h2,
px4h3,
jac1,
jac2,
jac3,
jac4,
jac5,
jac6,
integer, intent(in) nel,
integer, intent(in) jhbe )

Definition at line 30 of file sderit3.F.

44C-----------------------------------------------
45C I m p l i c i t T y p e s
46C-----------------------------------------------
47#include "implicit_f.inc"
48#include "comlock.inc"
49C-----------------------------------------------
50C G l o b a l P a r a m e t e r s
51C-----------------------------------------------
52#include "mvsiz_p.inc"
53C-----------------------------------------------
54C D u m m y A r g u m e n t s
55C-----------------------------------------------
56 double precision
57 . x1(*), x2(*), x3(*), x4(*), x5(*), x6(*), x7(*), x8(*),
58 . y1(*), y2(*), y3(*), y4(*), y5(*), y6(*), y7(*), y8(*),
59 . z1(*), z2(*), z3(*), z4(*), z5(*), z6(*), z7(*), z8(*)
60 INTEGER, INTENT(IN) :: NGL(*), NEL, JHBE
62 . off(*),det(*),
63 . px1(*), px2(*), px3(*), px4(*),
64 . py1(*), py2(*), py3(*), py4(*),
65 . pz1(*), pz2(*), pz3(*), pz4(*),
66 . px1h1(*), px1h2(*), px1h3(*),
67 . px2h1(*), px2h2(*), px2h3(*),
68 . px3h1(*), px3h2(*), px3h3(*),
69 . px4h1(*), px4h2(*), px4h3(*),
70 . jac1(*), jac2(*), jac3(*),
71 . jac4(*), jac5(*), jac6(*)
72C-----------------------------------------------
73C L o c a l V a r i a b l e s
74C-----------------------------------------------
75 INTEGER I, J
77 . dett , jac7(mvsiz), jac8(mvsiz) , jac9(mvsiz),
78 . jaci1, jaci2, jaci3,
79 . jaci4, jaci5, jaci6,
80 . jaci7, jaci8, jaci9,
81 . x17 , x28 , x35 , x46,
82 . y17 , y28 , y35 , y46,
83 . z17 , z28 , z35 , z46,
84 . jac_59_68(mvsiz), jac_67_49(mvsiz), jac_48_57(mvsiz),
85 . x_17_46, x_28_35,
86 . y_17_46, y_28_35,
87 . z_17_46, z_28_35,
88 . hx, hy, hz
89C
90C Jacobian matrix
91 DO i=1,nel
92 x17=x7(i)-x1(i)
93 x28=x8(i)-x2(i)
94 x35=x5(i)-x3(i)
95 x46=x6(i)-x4(i)
96C
97 y17=y7(i)-y1(i)
98 y28=y8(i)-y2(i)
99 y35=y5(i)-y3(i)
100 y46=y6(i)-y4(i)
101C
102 z17=z7(i)-z1(i)
103 z28=z8(i)-z2(i)
104 z35=z5(i)-z3(i)
105 z46=z6(i)-z4(i)
106C
107 jac1(i)=x17+x28-x35-x46
108 jac2(i)=y17+y28-y35-y46
109 jac3(i)=z17+z28-z35-z46
110C
111 x_17_46=x17+x46
112 x_28_35=x28+x35
113 y_17_46=y17+y46
114 y_28_35=y28+y35
115 z_17_46=z17+z46
116 z_28_35=z28+z35
117C
118 jac4(i)=x_17_46+x_28_35
119 jac5(i)=y_17_46+y_28_35
120 jac6(i)=z_17_46+z_28_35
121C
122 jac7(i)=x_17_46-x_28_35
123 jac8(i)=y_17_46-y_28_35
124 jac9(i)=z_17_46-z_28_35
125 ENDDO
126C
127 DO i=1,nel
128 jac_59_68(i)=jac5(i)*jac9(i)-jac6(i)*jac8(i)
129 jac_67_49(i)=jac6(i)*jac7(i)-jac4(i)*jac9(i)
130 jac_48_57(i)=jac4(i)*jac8(i)-jac5(i)*jac7(i)
131 ENDDO
132C
133 DO i=1,nel
134 det(i)=one_over_64*(jac1(i)*jac_59_68(i)+jac2(i)*jac_67_49(i)+jac3(i)*jac_48_57(i))
135 ENDDO
136C
137C Check volume
138 CALL schkjab3(
139 1 off, det, ngl, nel)
140C
141C Jacobian matrix inverse
142 DO i=1,nel
143 dett=one_over_64/det(i)
144 jaci1=dett*jac_59_68(i)
145 jaci4=dett*jac_67_49(i)
146 jaci7=dett*jac_48_57(i)
147 jaci2=dett*(jac3(i)*jac8(i)-jac2(i)*jac9(i))
148 jaci5=dett*(jac1(i)*jac9(i)-jac3(i)*jac7(i))
149 jaci8=dett*(jac2(i)*jac7(i)-jac1(i)*jac8(i))
150 jaci3=dett*(jac2(i)*jac6(i)-jac3(i)*jac5(i))
151 jaci6=dett*(jac3(i)*jac4(i)-jac1(i)*jac6(i))
152 jaci9=dett*(jac1(i)*jac5(i)-jac2(i)*jac4(i))
153C
154 px1(i)=-jaci1-jaci2-jaci3
155 py1(i)=-jaci4-jaci5-jaci6
156 pz1(i)=-jaci7-jaci8-jaci9
157C
158 px2(i)=-jaci1-jaci2+jaci3
159 py2(i)=-jaci4-jaci5+jaci6
160 pz2(i)=-jaci7-jaci8+jaci9
161C
162 px3(i)= jaci1-jaci2+jaci3
163 py3(i)= jaci4-jaci5+jaci6
164 pz3(i)= jaci7-jaci8+jaci9
165C
166 px4(i)= jaci1-jaci2-jaci3
167 py4(i)= jaci4-jaci5-jaci6
168 pz4(i)= jaci7-jaci8-jaci9
169 ENDDO
170C
171 IF(jhbe/=0)THEN
172 DO i=1,nel
173 hx=(x1(i)-x2(i)+x3(i)-x4(i)+x5(i)-x6(i)+x7(i)-x8(i))
174 hy=(y1(i)-y2(i)+y3(i)-y4(i)+y5(i)-y6(i)+y7(i)-y8(i))
175 hz=(z1(i)-z2(i)+z3(i)-z4(i)+z5(i)-z6(i)+z7(i)-z8(i))
176 px1h1(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
177 px2h1(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
178 px3h1(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
179 px4h1(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
180 ENDDO
181C
182 DO i=1,nel
183 hx=(x1(i)+x2(i)-x3(i)-x4(i)-x5(i)-x6(i)+x7(i)+x8(i))
184 hy=(y1(i)+y2(i)-y3(i)-y4(i)-y5(i)-y6(i)+y7(i)+y8(i))
185 hz=(z1(i)+z2(i)-z3(i)-z4(i)-z5(i)-z6(i)+z7(i)+z8(i))
186 px1h2(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
187 px2h2(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
188 px3h2(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
189 px4h2(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
190 ENDDO
191C
192 DO i=1,nel
193 hx=(x1(i)-x2(i)-x3(i)+x4(i)-x5(i)+x6(i)+x7(i)-x8(i))
194 hy=(y1(i)-y2(i)-y3(i)+y4(i)-y5(i)+y6(i)+y7(i)-y8(i))
195 hz=(z1(i)-z2(i)-z3(i)+z4(i)-z5(i)+z6(i)+z7(i)-z8(i))
196 px1h3(i)=px1(i)*hx+ py1(i)*hy+pz1(i)*hz
197 px2h3(i)=px2(i)*hx+ py2(i)*hy+pz2(i)*hz
198 px3h3(i)=px3(i)*hx+ py3(i)*hy+pz3(i)*hz
199 px4h3(i)=px4(i)*hx+ py4(i)*hy+pz4(i)*hz
200 ENDDO
201 ENDIF
202 RETURN
203C-----------
#define my_real
Definition cppsort.cpp:32
subroutine schkjab3(off, det, ngl, nel)
Definition schkjab3.F:39