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

Go to the source code of this file.

Functions/Subroutines

subroutine s8zjac_ic (xd1, xd2, xd3, xd4, xd5, xd6, xd7, xd8, yd1, yd2, yd3, yd4, yd5, yd6, yd7, yd8, zd1, zd2, zd3, zd4, zd5, zd6, zd7, zd8, aj1, aj2, aj3, aj4, aj5, aj6, aj7, aj8, aj9, hx, hy, hz, pxc1, pxc2, pxc3, pxc4, pyc1, pyc2, pyc3, pyc4, pzc1, pzc2, pzc3, pzc4, jac_i, index, nel)

Function/Subroutine Documentation

◆ s8zjac_ic()

subroutine s8zjac_ic ( 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,
aj1,
aj2,
aj3,
aj4,
aj5,
aj6,
aj7,
aj8,
aj9,
hx,
hy,
hz,
pxc1,
pxc2,
pxc3,
pxc4,
pyc1,
pyc2,
pyc3,
pyc4,
pzc1,
pzc2,
pzc3,
pzc4,
jac_i,
integer, dimension(*) index,
integer nel )

Definition at line 28 of file s8zjac_ic.F.

40C-----------------------------------------------
41C I m p l i c i t T y p e s
42C-----------------------------------------------
43#include "implicit_f.inc"
44C-----------------------------------------------
45C G l o b a l P a r a m e t e r s
46C-----------------------------------------------
47#include "mvsiz_p.inc"
48C-----------------------------------------------
49C D u m m y A r g u m e n t s
50C-----------------------------------------------
51 INTEGER NEL,INDEX(*)
53 . hx(4,*), hy(4,*), hz(4,*),
54 . aj1(*),aj2(*),aj3(*),
55 . aj4(*),aj5(*),aj6(*),
56 . aj7(*),aj8(*),aj9(*),jac_i(10,*),
57 . pxc1(*), pxc2(*), pxc3(*), pxc4(*),
58 . pyc1(*), pyc2(*), pyc3(*), pyc4(*),
59 . pzc1(*), pzc2(*), pzc3(*), pzc4(*)
60 double precision
61 . xd1(mvsiz), xd2(mvsiz), xd3(mvsiz), xd4(mvsiz),
62 . xd5(mvsiz), xd6(mvsiz), xd7(mvsiz), xd8(mvsiz),
63 . yd1(mvsiz), yd2(mvsiz), yd3(mvsiz), yd4(mvsiz),
64 . yd5(mvsiz), yd6(mvsiz), yd7(mvsiz), yd8(mvsiz),
65 . zd1(mvsiz), zd2(mvsiz), zd3(mvsiz), zd4(mvsiz),
66 . zd5(mvsiz), zd6(mvsiz), zd7(mvsiz), zd8(mvsiz)
67C-----------------------------------------------
68C L o c a l V a r i a b l e s
69C-----------------------------------------------
70 INTEGER I, J
71
73 . det(mvsiz) ,dett(mvsiz)
74C 12
76 . x17(mvsiz) , x28(mvsiz) , x35(mvsiz) , x46(mvsiz),
77 . y17(mvsiz) , y28(mvsiz) , y35(mvsiz) , y46(mvsiz),
78 . z17(mvsiz) , z28(mvsiz) , z35(mvsiz) , z46(mvsiz),
79 . jac_59_68(mvsiz), jac_67_49(mvsiz), jac_48_57(mvsiz),
80 . jac_38_29(mvsiz), jac_19_37(mvsiz), jac_27_18(mvsiz),
81 . jac_26_35(mvsiz), jac_34_16(mvsiz), jac_15_24(mvsiz),
82 . aji1, aji2, aji3,aji4, aji5, aji6,aji7, aji8, aji9,
83 . aj12, aj45, aj78,aj12p, aj45p, aj78p,
84 . a17(mvsiz) , a28(mvsiz) ,
85 . b17(mvsiz) , b28(mvsiz) ,
86 . c17(mvsiz) , c28(mvsiz) ,jac_1(10,mvsiz)
87C-----------------------------------------------
88 DO i=1,nel
89 x17(i)=xd7(i)-xd1(i)
90 x28(i)=xd8(i)-xd2(i)
91 x35(i)=xd5(i)-xd3(i)
92 x46(i)=xd6(i)-xd4(i)
93 y17(i)=yd7(i)-yd1(i)
94 y28(i)=yd8(i)-yd2(i)
95 y35(i)=yd5(i)-yd3(i)
96 y46(i)=yd6(i)-yd4(i)
97 z17(i)=zd7(i)-zd1(i)
98 z28(i)=zd8(i)-zd2(i)
99 z35(i)=zd5(i)-zd3(i)
100 z46(i)=zd6(i)-zd4(i)
101 END DO
102C
103 DO i=1,nel
104 aj4(i)=x17(i)+x28(i)-x35(i)-x46(i)
105 aj5(i)=y17(i)+y28(i)-y35(i)-y46(i)
106 aj6(i)=z17(i)+z28(i)-z35(i)-z46(i)
107 a17(i)=x17(i)+x46(i)
108 a28(i)=x28(i)+x35(i)
109 b17(i)=y17(i)+y46(i)
110 b28(i)=y28(i)+y35(i)
111 c17(i)=z17(i)+z46(i)
112 c28(i)=z28(i)+z35(i)
113C
114 aj7(i)=a17(i)+a28(i)
115 aj8(i)=b17(i)+b28(i)
116 aj9(i)=c17(i)+c28(i)
117 aj1(i)=a17(i)-a28(i)
118 aj2(i)=b17(i)-b28(i)
119 aj3(i)=c17(i)-c28(i)
120 END DO
121C
122C JACOBIAN
123C
124 DO i=1,nel
125 jac_59_68(i)=aj5(i)*aj9(i)-aj6(i)*aj8(i)
126 jac_67_49(i)=aj6(i)*aj7(i)-aj4(i)*aj9(i)
127 jac_38_29(i)=(-aj2(i)*aj9(i)+aj3(i)*aj8(i))
128 jac_19_37(i)=( aj1(i)*aj9(i)-aj3(i)*aj7(i))
129 jac_27_18(i)=(-aj1(i)*aj8(i)+aj2(i)*aj7(i))
130 jac_26_35(i)=( aj2(i)*aj6(i)-aj3(i)*aj5(i))
131 jac_34_16(i)=(-aj1(i)*aj6(i)+aj3(i)*aj4(i))
132 jac_15_24(i)=( aj1(i)*aj5(i)-aj2(i)*aj4(i))
133 jac_48_57(i)=aj4(i)*aj8(i)-aj5(i)*aj7(i)
134C
135 det(i)=one_over_64*(aj1(i)*jac_59_68(i)+aj2(i)*jac_67_49(i)+aj3(i)*jac_48_57(i))
136 dett(i)=one_over_64/det(i)
137 ENDDO
138C
139 DO i=1,nel
140 aji1=dett(i)*jac_59_68(i)
141 aji4=dett(i)*jac_67_49(i)
142 aji7=dett(i)*jac_48_57(i)
143 aji2=dett(i)*jac_38_29(i)
144 aji5=dett(i)*jac_19_37(i)
145 aji8=dett(i)*jac_27_18(i)
146 aji3=dett(i)*jac_26_35(i)
147 aji6=dett(i)*jac_34_16(i)
148 aji9=dett(i)*jac_15_24(i)
149 aj12=aji1-aji2
150 aj45=aji4-aji5
151 aj78=aji7-aji8
152 pxc2(i)= aj12-aji3
153 pyc2(i)= aj45-aji6
154 pzc2(i)= aj78-aji9
155 pxc4(i)=-aj12-aji3
156 pyc4(i)=-aj45-aji6
157 pzc4(i)=-aj78-aji9
158 aj12p=aji1+aji2
159 aj45p=aji4+aji5
160 aj78p=aji7+aji8
161 pxc1(i)=-aj12p-aji3
162 pyc1(i)=-aj45p-aji6
163 pzc1(i)=-aj78p-aji9
164 pxc3(i)=aj12p-aji3
165 pyc3(i)=aj45p-aji6
166 pzc3(i)=aj78p-aji9
167 jac_1(1,i)=aji1
168 jac_1(4,i)=aji4
169 jac_1(7,i)=aji7
170 jac_1(2,i)=aji2
171 jac_1(5,i)=aji5
172 jac_1(8,i)=aji8
173 jac_1(3,i)=aji3
174 jac_1(6,i)=aji6
175 jac_1(9,i)=aji9
176 jac_1(10,i)=det(i)
177 ENDDO
178C mode 1
179C 1 1 -1 -1 -1 -1 1 1
180 DO i=1,nel
181 hx(1,i)=(xd1(i)+xd2(i)-xd3(i)-xd4(i)-xd5(i)-xd6(i)+xd7(i)+xd8(i))
182 hy(1,i)=(yd1(i)+yd2(i)-yd3(i)-yd4(i)-yd5(i)-yd6(i)+yd7(i)+yd8(i))
183 hz(1,i)=(zd1(i)+zd2(i)-zd3(i)-zd4(i)-zd5(i)-zd6(i)+zd7(i)+zd8(i))
184 ENDDO
185C mode 2
186C 1 -1 -1 1 -1 1 1 -1
187 DO i=1,nel
188 hx(2,i)=(xd1(i)-xd2(i)-xd3(i)+xd4(i)-xd5(i)+xd6(i)+xd7(i)-xd8(i))
189 hy(2,i)=(yd1(i)-yd2(i)-yd3(i)+yd4(i)-yd5(i)+yd6(i)+yd7(i)-yd8(i))
190 hz(2,i)=(zd1(i)-zd2(i)-zd3(i)+zd4(i)-zd5(i)+zd6(i)+zd7(i)-zd8(i))
191 ENDDO
192C mode 3
193C 1 -1 1 -1 1 -1 1 -1
194 DO i=1,nel
195 hx(3,i)=(xd1(i)-xd2(i)+xd3(i)-xd4(i)+xd5(i)-xd6(i)+xd7(i)-xd8(i))
196 hy(3,i)=(yd1(i)-yd2(i)+yd3(i)-yd4(i)+yd5(i)-yd6(i)+yd7(i)-yd8(i))
197 hz(3,i)=(zd1(i)-zd2(i)+zd3(i)-zd4(i)+zd5(i)-zd6(i)+zd7(i)-zd8(i))
198 ENDDO
199C mode 4
200C -1 1 -1 1 1 -1 1 -1
201 DO i=1,nel
202 hx(4,i)=(-xd1(i)+xd2(i)-xd3(i)+xd4(i)+xd5(i)-xd6(i)+xd7(i)-xd8(i))
203 hy(4,i)=(-yd1(i)+yd2(i)-yd3(i)+yd4(i)+yd5(i)-yd6(i)+yd7(i)-yd8(i))
204 hz(4,i)=(-zd1(i)+zd2(i)-zd3(i)+zd4(i)+zd5(i)-zd6(i)+zd7(i)-zd8(i))
205 ENDDO
206#include "vectorize.inc"
207 DO j=1,nel
208 i = index(j)
209 jac_i(1:10,i) = jac_1(1:10,j)
210 ENDDO
211C
212 RETURN
213C
#define my_real
Definition cppsort.cpp:32