OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i8lagm.F File Reference
#include "implicit_f.inc"
#include "mvsiz_p.inc"
#include "task_c.inc"
#include "com04_c.inc"
#include "com08_c.inc"
#include "comlock.inc"
#include "lockon.inc"
#include "lockoff.inc"
#include "constant.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine i8lagm (x, v, lll, jll, sll, xll, candn, cande, i_stok, ixs, iadll, eminx, nsv, nelem, n_mul_mx, itask, a, itied, nint, nkmax, comntag)
subroutine i8lll (llt, lll, jll, sll, xll, v, xx, yy, zz, iii, iadll, n_mul_mx, a, x, itied, nint, nkmax, comntag)
subroutine i8rst (llt, r, s, t, ni, nx, ny, nz, xx, yy, zz)
subroutine i8ni (llt, rr, ss, tt, ni)
subroutine i8rstn (llt, rr, ss, tt, ni, conv, drdx, drdy, drdz, dsdx, dsdy, dsdz, dtdx, dtdy, dtdz, xx, yy, zz, r, s, t)
subroutine i8deri (llt, rr, ss, tt, ni, drdx, drdy, drdz, dsdx, dsdy, dsdz, dtdx, dtdy, dtdz, dxdr, dydr, dzdr, dxdt, dydt, dzdt, xx, yy, zz)

Function/Subroutine Documentation

◆ i8deri()

subroutine i8deri ( integer llt,
rr,
ss,
tt,
ni,
drdx,
drdy,
drdz,
dsdx,
dsdy,
dsdz,
dtdx,
dtdy,
dtdz,
dxdr,
dydr,
dzdr,
dxdt,
dydt,
dzdt,
xx,
yy,
zz )

Definition at line 777 of file i8lagm.F.

781C-----------------------------------------------
782C I m p l i c i t T y p e s
783C-----------------------------------------------
784#include "implicit_f.inc"
785C-----------------------------------------------
786C G l o b a l P a r a m e t e r s
787C-----------------------------------------------
788#include "mvsiz_p.inc"
789C-----------------------------------------------
790C D u m m y A r g u m e n t s
791C-----------------------------------------------
792 INTEGER LLT
793C REAL
794 my_real
795 . dxdr(mvsiz), dydr(mvsiz), dzdr(mvsiz),
796 . dxdt(mvsiz), dydt(mvsiz), dzdt(mvsiz),
797 . drdx(mvsiz), dsdx(mvsiz), dtdx(mvsiz),
798 . drdy(mvsiz), dsdy(mvsiz), dtdy(mvsiz),
799 . drdz(mvsiz), dsdz(mvsiz), dtdz(mvsiz),
800 . xx(mvsiz,9) ,yy(mvsiz,9),zz(mvsiz,9),
801 . ni(mvsiz,9) ,rr(mvsiz) ,ss(mvsiz) ,tt(mvsiz)
802C-----------------------------------------------
803C L o c a l V a r i a b l e s
804C-----------------------------------------------
805 INTEGER I
806 my_real
807 . dxds(mvsiz), dyds(mvsiz), dzds(mvsiz),
808 . dnidr(8),dnids(8),dnidt(8),
809 . d, det(mvsiz)
810 my_real
811 . u_m_r,u_p_r,u_m_s,u_p_s,u_m_t,u_p_t,
812 . ums_umt,ums_upt,ups_umt,ups_upt,
813 .
814 .
815 . a,r05,s05,t05
816C-----------------------------------------------
817C/*
818C
819C
820C
821C ^ S _ T
822C | /|
823C | /
824C | /
825C ( 7)===========|=====================( 6)
826C //| | / //|
827C // | | / //||
828C // | | / // ||
829C // | (Is) | + // ||
830C // | | / // ||
831C // | + / // ||
832C // | | / // ||
833C // ( 3)-----------|--------------//-----( 2)
834C // / |/ // //
835C R <-------//- -+ - - - - - - - - -# - - - - - -//- -+ //
836C // / /| // //
837C ( 8)=================================( 5) //
838C || / / | || //
839C || / / + || //
840C || / / || //
841C || / + || //
842C || / || //
843C || / ||//
844C ||/ ||/
845C ( 4)=================================( 1)
846C
847C
848C*/
849C-----------------------------------------------
850C
851C-----------------------------------------------
852C _
853C \
854C x(r,s,t) = /_ (xi * Ni(r,s,t))
855C _
856C \
857C y(r,s,t) = /_ (yi * Ni(r,s,t))
858C _
859C \
860C z(r,s,t) = /_ zi * Ni(r,s,t))
861C
862C _
863C \
864C dx/dr = /_ (xi * dNi/dr)
865C ...
866C
867C [dx/dr dy/dr dz/dr]
868C [J] = |dx/ds dy/ds dz/ds|
869C [dx/dt dy/dt dz/dt]
870C
871C |r| |r| -1 |xs-x|
872C {s} = {s} + [J] {ys-y}
873C |t| |t| |zs-z|
874C
875C-----------------------------------------------------------------------
876C Ni; dNi/dr; dNi/ds; dNi/dt
877C-----------------------------------------------------------------------
878 DO i=1,llt
879 r05 = half*rr(i)
880 s05 = half*ss(i)
881 t05 = half*tt(i)
882C
883 u_m_r = half - r05
884 u_p_r = half + r05
885C
886 u_m_s = half - s05
887 u_p_s = half + s05
888C
889 u_m_t = half - t05
890 u_p_t = half + t05
891C
892 ums_umt = u_m_s * u_m_t
893 ums_upt = u_m_s * u_p_t
894 ups_umt = u_p_s * u_m_t
895 ups_upt = u_p_s * u_p_t
896C
897 ni(i,1) = u_m_r * ums_umt
898 ni(i,5) = u_m_r * ups_umt
899 ni(i,2) = u_m_r * ums_upt
900 ni(i,6) = u_m_r * ups_upt
901 ni(i,3) = u_p_r * ums_upt
902 ni(i,7) = u_p_r * ups_upt
903 ni(i,4) = u_p_r * ums_umt
904 ni(i,8) = u_p_r * ups_umt
905C
906 a = half
907 dnidr(1) = -ums_umt * a
908 dnidr(2) = -ums_upt * a
909 dnidr(3) = -dnidr(2)
910 dnidr(4) = -dnidr(1)
911 dnidr(5) = -ups_umt * a
912 dnidr(6) = -ups_upt * a
913 dnidr(7) = -dnidr(6)
914 dnidr(8) = -dnidr(5)
915C
916 dnids(1) = -u_m_t * u_m_r * a
917 dnids(2) = -u_p_t * u_m_r * a
918 dnids(3) = -u_p_t * u_p_r * a
919 dnids(4) = -u_m_t * u_p_r * a
920 dnids(5) = -dnids(1)
921 dnids(6) = -dnids(2)
922 dnids(7) = -dnids(3)
923 dnids(8) = -dnids(4)
924C
925 dnidt(1) = -u_m_r * u_m_s * a
926 dnidt(2) = -dnidt(1)
927 dnidt(3) = u_p_r * u_m_s * a
928 dnidt(4) = -dnidt(3)
929 dnidt(5) = -u_m_r * u_p_s * a
930 dnidt(6) = -dnidt(5)
931 dnidt(7) = u_p_r * u_p_s * a
932 dnidt(8) = -dnidt(7)
933C------------------------------------
934 ni(i,9) = -one
935C
936C-----------------------------------------------------------------------
937C dx/dr; dx/ds; dx/dt
938C-----------------------------------------------------------------------
939 dxdr(i) = dnidr(1)*xx(i,1) + dnidr(2)*xx(i,2) + dnidr(3)*xx(i,3)
940 + + dnidr(4)*xx(i,4) + dnidr(5)*xx(i,5) + dnidr(6)*xx(i,6)
941 + + dnidr(7)*xx(i,7) + dnidr(8)*xx(i,8)
942C
943 dxds(i) = dnids(1)*xx(i,1) + dnids(2)*xx(i,2) + dnids(3)*xx(i,3)
944 + + dnids(4)*xx(i,4) + dnids(5)*xx(i,5) + dnids(6)*xx(i,6)
945 + + dnids(7)*xx(i,7) + dnids(8)*xx(i,8)
946C
947 dxdt(i) = dnidt(1)*xx(i,1) + dnidt(2)*xx(i,2) + dnidt(3)*xx(i,3)
948 + + dnidt(4)*xx(i,4) + dnidt(5)*xx(i,5) + dnidt(6)*xx(i,6)
949 + + dnidt(7)*xx(i,7) + dnidt(8)*xx(i,8)
950C-----------------------------------------------------------------------
951C dy/dr; dy/ds; dy/dt
952C-----------------------------------------------------------------------
953 dydr(i) = dnidr(1)*yy(i,1) + dnidr(2)*yy(i,2) + dnidr(3)*yy(i,3)
954 + + dnidr(4)*yy(i,4) + dnidr(5)*yy(i,5) + dnidr(6)*yy(i,6)
955 + + dnidr(7)*yy(i,7) + dnidr(8)*yy(i,8)
956C
957 dyds(i) = dnids(1)*yy(i,1) + dnids(2)*yy(i,2) + dnids(3)*yy(i,3)
958 + + dnids(4)*yy(i,4) + dnids(5)*yy(i,5) + dnids(6)*yy(i,6)
959 + + dnids(7)*yy(i,7) + dnids(8)*yy(i,8)
960C
961 dydt(i) = dnidt(1)*yy(i,1) + dnidt(2)*yy(i,2) + dnidt(3)*yy(i,3)
962 + + dnidt(4)*yy(i,4) + dnidt(5)*yy(i,5) + dnidt(6)*yy(i,6)
963 + + dnidt(7)*yy(i,7) + dnidt(8)*yy(i,8)
964C-----------------------------------------------------------------------
965C dz/dr; dz/ds; dz/dt
966C-----------------------------------------------------------------------
967 dzdr(i) = dnidr(1)*zz(i,1) + dnidr(2)*zz(i,2) + dnidr(3)*zz(i,3)
968 + + dnidr(4)*zz(i,4) + dnidr(5)*zz(i,5) + dnidr(6)*zz(i,6)
969 + + dnidr(7)*zz(i,7) + dnidr(8)*zz(i,8)
970C
971 dzds(i) = dnids(1)*zz(i,1) + dnids(2)*zz(i,2) + dnids(3)*zz(i,3)
972 + + dnids(4)*zz(i,4) + dnids(5)*zz(i,5) + dnids(6)*zz(i,6)
973 + + dnids(7)*zz(i,7) + dnids(8)*zz(i,8)
974C
975 dzdt(i) = dnidt(1)*zz(i,1) + dnidt(2)*zz(i,2) + dnidt(3)*zz(i,3)
976 + + dnidt(4)*zz(i,4) + dnidt(5)*zz(i,5) + dnidt(6)*zz(i,6)
977 + + dnidt(7)*zz(i,7) + dnidt(8)*zz(i,8)
978C-----------------------------------------------------------------------
979C -1
980C [j] jacobian inversion
981C-----------------------------------------------------------------------
982 drdx(i)=dyds(i)*dzdt(i)-dzds(i)*dydt(i)
983 drdy(i)=dzds(i)*dxdt(i)-dxds(i)*dzdt(i)
984 drdz(i)=dxds(i)*dydt(i)-dyds(i)*dxdt(i)
985C
986 dsdz(i)=dxdt(i)*dydr(i)-dydt(i)*dxdr(i)
987 dsdy(i)=dzdt(i)*dxdr(i)-dxdt(i)*dzdr(i)
988 dsdx(i)=dydt(i)*dzdr(i)-dzdt(i)*dydr(i)
989C
990 dtdx(i)=dydr(i)*dzds(i)-dzdr(i)*dyds(i)
991 dtdy(i)=dzdr(i)*dxds(i)-dxdr(i)*dzds(i)
992 dtdz(i)=dxdr(i)*dyds(i)-dydr(i)*dxds(i)
993C
994 det(i) = dxdr(i) * drdx(i)
995 . + dydr(i) * drdy(i)
996 . + dzdr(i) * drdz(i)
997C
998c
999c print *, "Det",DET(I)
1000c print *, "DXDR(I),DYDR(I),DZDR(I)",DXDR(I),DYDR(I),DZDR(I)
1001c print *, "DXDs(I),DYDs(I),DZDs(I)",DXDs(I),DYDs(I),DZDs(I)
1002c print *, "DXDt(I),DYDt(I),DZDt(I)",DXDt(I),DYDt(I),DZDt(I)
1003c
1004 ENDDO
1005C
1006 DO i=1,llt
1007C-----------------------------------------------------------------------
1008C -1
1009C [J] Inversion of the Jacobian Suite
1010C-----------------------------------------------------------------------
1011 d = one/det(i)
1012 drdx(i)=d*drdx(i)
1013 dsdx(i)=d*dsdx(i)
1014 dtdx(i)=d*dtdx(i)
1015C
1016 drdy(i)=d*drdy(i)
1017 dsdy(i)=d*dsdy(i)
1018 dtdy(i)=d*dtdy(i)
1019C
1020 drdz(i)=d*drdz(i)
1021 dsdz(i)=d*dsdz(i)
1022 dtdz(i)=d*dtdz(i)
1023C
1024c
1025c print *, "DRDX(I),DRDY(I),DRDZ(I)",DRDX(I),DRDY(I),DRDZ(I)
1026c print *, "DSDX(I),DSDY(I),DSDZ(I)",DSDX(I),DSDY(I),DSDZ(I)
1027c print *, "DTDX(I),DTDY(I),DTDZ(I)",DTDX(I),DTDY(I),DTDZ(I)
1028c
1029 ENDDO
1030C-----------------------------------------------------------------------
1031 RETURN
#define my_real
Definition cppsort.cpp:32

◆ i8lagm()

subroutine i8lagm ( x,
v,
integer, dimension(*) lll,
integer, dimension(*) jll,
integer, dimension(*) sll,
xll,
integer, dimension(*) candn,
integer, dimension(*) cande,
integer i_stok,
integer, dimension(nixs,*) ixs,
integer, dimension(*) iadll,
eminx,
integer, dimension(*) nsv,
integer, dimension(*) nelem,
integer n_mul_mx,
integer itask,
a,
integer itied,
integer nint,
integer nkmax,
integer, dimension(*) comntag )

Definition at line 32 of file i8lagm.F.

37 use element_mod , only : nixs
38C-----------------------------------------------
39C I m p l i c i t T y p e s
40C-----------------------------------------------
41#include "implicit_f.inc"
42C-----------------------------------------------
43C G l o b a l P a r a m e t e r s
44C-----------------------------------------------
45#include "mvsiz_p.inc"
46C-----------------------------------------------
47C C o m m o n B l o c k s
48C-----------------------------------------------
49#include "task_c.inc"
50#include "com04_c.inc"
51#include "com08_c.inc"
52C-----------------------------------------------
53C D u m m y A r g u m e n t s
54C-----------------------------------------------
55 INTEGER I_STOK,N_MUL_MX,ITASK,ITIED,NINT,NKMAX ,
56 . LLL(*),JLL(*),SLL(*),CANDN(*),CANDE(*),COMNTAG(*),
57 . IXS(NIXS,*),IADLL(*),NSV(*) ,NELEM(*)
58C REAL
60 . x(3,*),v(3,*),xll(*),
61 . eminx(6,*),a(3,*)
62C-----------------------------------------------
63C L o c a l V a r i a b l e s
64C-----------------------------------------------
65 INTEGER I, K, IE, IS, IC, III(MVSIZ, 9), LLT, NFT, LE, FIRST, LAST
67 . xx(mvsiz,9),yy(mvsiz,9),zz(mvsiz,9),
68 . dist
69C-----------------------------------------------
70C
71C
72C | M | Lt| | a | M ao
73C |---+---| | = |
74C | l | 0 | | la | bo
75C
76C [M] a + [L]t la = [M] ao
77C [L] a = bo
78C
79C a = -[M]-1[L]t la + ao
80C [L][M]-1[L]t la = [L] ao - bo
81C
82C on pose:
83C [H] = [L][M]-1[L]t
84C b = [L] ao - bo
85C
86C [h] la = b
87C
88C a = ao - [M]-1[L]t la
89C-----------------------------------------------
90C
91C la : lambda(nc)
92C ao : A(NUMNOD)
93C L : XLL(NK,NC)
94C M : MAS(NUMNOD)
95C [L][M]-1[L]t la : HLA(NC)
96C [L] ao - b : B(NC)
97C [M]-1[L]t la : LTLA(NUMNOD)
98C
99C nc: number of contacts
100C NK: Number of node for contact (8+1.16+1.8+8.16+16)
101C
102C IC : contact number (1,NC)
103C IK : local node number for a contact (1,NK)
104C I : global node number (1,NUMNOD)
105C
106C IADLL(NC) : IAD = IADLL(IC)
107C LLL(NC*(9,27)) : I = LLL(IAD+1,2...IADNEXT-1)
108C-----------------------------------------------
109C evaluation of b:
110C
111C Vs = Somme(Ni Vi)
112C Vs_ + dt As = Somme(Ni Vi_) + Somme(dt Ni Ai)
113C Somme(dt Ni Ai) - dt As = Vs_ -Somme(Ni Vi_)
114C [L] = dt {N1,N2,..,N15,-1}
115C bo = [L] a = -[L]/dt v_
116C b = [L] ao - bo
117C b = [L] ao + [L]/dt v_ = [L] (v_ + ao dt)/dt
118C-----------------------------------------------
119C b = [L] vo+/dt + vout
120C-----------------------------------------------
121C-----------------------------------------------------------------------
122C loop over contact candidates
123C-----------------------------------------------------------------------
124 first = 1 + i_stok * itask / nthread
125 last = i_stok*(itask+1) / nthread
126 llt = 0
127 nft=llt+1
128 DO ic=first,last
129 le = cande(ic)
130 ie = nelem(le)
131C-----------------------------------------------------------------------
132C test si shell 16
133C-----------------------------------------------------------------------
134 IF(ie.ge .1.AND.ie.le .numels8)THEN
135 is = nsv(candn(ic))
136 dist = -ep30
137 dist = max(eminx(1,le)-x(1,is)-dt2*(v(1,is)+dt12*a(1,is)),
138 . x(1,is)+dt2*(v(1,is)+dt12*a(1,is))-eminx(4,le),dist)
139 dist = max(eminx(2,le)-x(2,is)-dt2*(v(2,is)+dt12*a(2,is)),
140 . x(2,is)+dt2*(v(2,is)+dt12*a(2,is))-eminx(5,le),dist)
141 dist = max(eminx(3,le)-x(3,is)-dt2*(v(3,is)+dt12*a(3,is)),
142 . x(3,is)+dt2*(v(3,is)+dt12*a(3,is))-eminx(6,le),dist)
143C-----------------------------------------------------------------------
144C test if inside the box
145C-----------------------------------------------------------------------
146 IF(dist.le .zero)THEN
147c
148c print *, "in la boite",XMIN,XMAX,YMIN,YMAX,ZMIN,ZMAX
149c
150 llt = llt+1
151 iii(llt,9)=is
152 xx(llt,9)=x(1,is)
153 yy(llt,9)=x(2,is)
154 zz(llt,9)=x(3,is)
155 DO k=1,8
156 iii(llt,k)=ixs(k+1,ie)
157 i = iii(llt,k)
158 xx(llt,k)=x(1,i)
159 yy(llt,k)=x(2,i)
160 zz(llt,k)=x(3,i)
161 ENDDO
162c
163c print *, "XX(1,1),XX(1,9)",XX(1,1),XX(1,9)
164c
165C-----------------------------------------------------------------------
166C calculation of [l] by mvsiz packet
167C-----------------------------------------------------------------------
168 IF(llt==mvsiz-1)THEN
169 CALL i8lll(
170 1 llt ,lll ,jll ,sll ,xll ,v ,
171 2 xx ,yy ,zz ,iii ,iadll ,
172 3 n_mul_mx ,a ,x ,itied ,nint ,nkmax ,
173 4 comntag )
174 nft=llt+1
175 llt = 0
176 ENDIF
177 ELSE
178c debug
179 k=0
180 ENDIF
181 ENDIF
182 ENDDO
183C-----------------------------------------------------------------------
184C calculation of [l] for the last packet
185C-----------------------------------------------------------------------
186 IF(llt/=0) CALL i8lll(
187 1 llt ,lll ,jll ,sll ,xll ,v ,
188 2 xx ,yy ,zz ,iii ,iadll ,
189 3 n_mul_mx ,a ,x ,itied ,nint ,nkmax ,
190 4 comntag )
191C
192C-----------------------------------------------
193 RETURN
subroutine i8lll(llt, lll, jll, sll, xll, v, xx, yy, zz, iii, iadll, n_mul_mx, a, x, itied, nint, nkmax, comntag)
Definition i8lagm.F:210
#define max(a, b)
Definition macros.h:21

◆ i8lll()

subroutine i8lll ( integer llt,
integer, dimension(*) lll,
integer, dimension(*) jll,
integer, dimension(*) sll,
xll,
v,
xx,
yy,
zz,
integer, dimension(mvsiz,9) iii,
integer, dimension(*) iadll,
integer n_mul_mx,
a,
x,
integer itied,
integer nint,
integer nkmax,
integer, dimension(*) comntag )

Definition at line 206 of file i8lagm.F.

210C-----------------------------------------------
211C M o d u l e s
212C-----------------------------------------------
213 USE message_mod
214C-----------------------------------------------
215C I m p l i c i t T y p e s
216C-----------------------------------------------
217#include "implicit_f.inc"
218#include "comlock.inc"
219C-----------------------------------------------
220C G l o b a l P a r a m e t e r s
221C-----------------------------------------------
222#include "mvsiz_p.inc"
223C-----------------------------------------------
224C C o m m o n B l o c k s
225C-----------------------------------------------
226#include "com08_c.inc"
227 COMMON /lagglob/n_mult
228 INTEGER N_MULT
229C-----------------------------------------------
230C D u m m y A r g u m e n t s
231C-----------------------------------------------
232 INTEGER LLT,N_MUL_MX,ITIED,NINT ,NKMAX
233 INTEGER LLL(*),JLL(*),SLL(*),COMNTAG(*),
234 . III(MVSIZ,9),IADLL(*)
235C REAL
236 my_real
237 . xll(*),v(3,*),a(3,*)
238 my_real
239 . xx(mvsiz,9),yy(mvsiz,9),zz(mvsiz,9),x(3,*)
240C-----------------------------------------------
241C L o c a l V a r i a b l e s
242C-----------------------------------------------
243 INTEGER I, IK, NK, IAD, NN
244 my_real
245 . vx,vy,vz,vn,aa
246 my_real
247 . r(mvsiz),s(mvsiz),t(mvsiz),
248 . nx(mvsiz), ny(mvsiz), nz(mvsiz),
249 . ni(mvsiz,9)
250C-----------------------------------------------
251C calculation de r,s,t
252C-----------------------------------------------
253c
254c print *, "XX(1,1),XX(1,9)",XX(1,1),XX(1,9)
255c
256 CALL i8rst(llt ,r ,s ,t ,ni ,
257 2 nx ,ny ,nz ,xx ,yy ,zz )
258C-----------------------------------------------
259C calculation of [l]
260C-----------------------------------------------
261 IF(itied==0)THEN
262 DO i=1,llt
263C-----------------------------------------------
264C Test if contact
265C-----------------------------------------------
266 IF(r(i)>=-one.AND.s(i)>=-one.AND.t(i)>=-one.AND.
267 . r(i)<= one.AND.s(i)<= one.AND.t(i)<= one)THEN
268C
269 nk = 9
270 vx = zero
271 vy = zero
272 vz = zero
273 DO ik=1,nk
274 vx = vx - (v(1,iii(i,ik))+dt12*a(1,iii(i,ik)))*ni(i,ik)
275 vy = vy - (v(2,iii(i,ik))+dt12*a(2,iii(i,ik)))*ni(i,ik)
276 vz = vz - (v(3,iii(i,ik))+dt12*a(3,iii(i,ik)))*ni(i,ik)
277 ENDDO
278c
279c print *, "vx,vy,vz s-m",vx,vy,vz
280c print *, "nx,ny,nz ", NX(I),NY(I),NZ(I)
281c
282 vn = nx(i)*vx + ny(i)*vy + nz(i)*vz
283C-----------------------------------------------
284C Test if incoming velocity in S
285C-----------------------------------------------
286 IF(s(i)*vn<=zero)THEN
287c
288c print *, "velocity entrante",vn
289c print *, "s = ",S(I)
290c
291c AA = DT12/SQRT(NX(I)*NX(I)+NY(I)*NY(I)+NZ(I)*NZ(I))
292 aa = one/sqrt(nx(i)*nx(i)+ny(i)*ny(i)+nz(i)*nz(i))
293 nx(i) = nx(i)*aa
294 ny(i) = ny(i)*aa
295 nz(i) = nz(i)*aa
296#include "lockon.inc"
297 n_mult=n_mult+1
298 IF(n_mult>n_mul_mx)THEN
299#include "lockoff.inc"
300 CALL ancmsg(msgid=84,anmode=aninfo)
301 CALL arret(2)
302 ENDIF
303 iadll(n_mult+1)=iadll(n_mult) + 27
304 IF(iadll(n_mult+1)-1>nkmax)THEN
305#include "lockoff.inc"
306 CALL ancmsg(msgid=84,anmode=aninfo)
307 CALL arret(2)
308 ENDIF
309 iad = iadll(n_mult) - 1
310 DO ik=1,9
311 lll(iad+ik) = iii(i,ik)
312 jll(iad+ik) = 1
313 sll(iad+ik) = 0
314 xll(iad+ik) = nx(i)*ni(i,ik)
315 lll(iad+ik+9) = iii(i,ik)
316 jll(iad+ik+9) = 2
317 sll(iad+ik+9) = 0
318 xll(iad+ik+9) = ny(i)*ni(i,ik)
319 lll(iad+ik+18) = iii(i,ik)
320 jll(iad+ik+18) = 3
321 sll(iad+ik+18) = 0
322 xll(iad+ik+18) = nz(i)*ni(i,ik)
323 nn = lll(iad+ik)
324 comntag(nn) = comntag(nn) + 1
325 ENDDO
326 sll(iad+9) = nint
327 sll(iad+18) = nint
328 sll(iad+27) = nint
329#include "lockoff.inc"
330 ENDIF
331 ENDIF
332 ENDDO
333 ELSEIF(itied==1)THEN
334C-----------------------------------------------
335C ITIED = 1
336C-----------------------------------------------
337 DO i=1,llt
338C-----------------------------------------------
339C Test if contact
340C-----------------------------------------------
341 IF(r(i)>=-one.AND.s(i)>=-one.AND.t(i)>=-one.AND.
342 . r(i)<= one.AND.s(i)<= one.AND.t(i)<= one)THEN
343C
344 nk = 9
345 vx = zero
346 vy = zero
347 vz = zero
348 DO ik=1,nk
349 vx = vx - (v(1,iii(i,ik))+dt12*a(1,iii(i,ik)))*ni(i,ik)
350 vy = vy - (v(2,iii(i,ik))+dt12*a(2,iii(i,ik)))*ni(i,ik)
351 vz = vz - (v(3,iii(i,ik))+dt12*a(3,iii(i,ik)))*ni(i,ik)
352 ENDDO
353c
354c print *, "vx,vy,vz s-m",vx,vy,vz
355c print *, "nx,ny,nz ", NX(I),NY(I),NZ(I)
356c
357 vn = nx(i)*vx + ny(i)*vy + nz(i)*vz
358C-----------------------------------------------
359C Test if incoming velocity in S
360C-----------------------------------------------
361 IF(s(i)*vn<=0.0)THEN
362c
363c print *, "velocity entrante",vn
364c print *, "s = ",S(I)
365c
366#include "lockon.inc"
367 IF(n_mult+3>n_mul_mx)THEN
368#include "lockoff.inc"
369 CALL ancmsg(msgid=84,anmode=aninfo)
370 CALL arret(2)
371 ENDIF
372 IF(iadll(n_mult+1)-1+9*3>nkmax)THEN
373#include "lockoff.inc"
374 CALL ancmsg(msgid=84,anmode=aninfo)
375 CALL arret(2)
376 ENDIF
377C
378 n_mult=n_mult+1
379 iadll(n_mult+1)=iadll(n_mult) + 9
380 iad = iadll(n_mult) - 1
381 DO ik=1,9
382 lll(iad+ik) = iii(i,ik)
383 jll(iad+ik) = 1
384 sll(iad+ik) = 0
385 xll(iad+ik) = ni(i,ik)
386 nn = lll(iad+ik)
387 comntag(nn) = comntag(nn) + 1
388 ENDDO
389 sll(iad+9) = nint
390C
391 n_mult=n_mult+1
392 iadll(n_mult+1)=iadll(n_mult) + 9
393 iad = iadll(n_mult) - 1
394 DO ik=1,9
395 lll(iad+ik) = iii(i,ik)
396 jll(iad+ik) = 2
397 sll(iad+ik) = 0
398 xll(iad+ik) = ni(i,ik)
399 nn = lll(iad+ik)
400 comntag(nn) = comntag(nn) + 1
401 ENDDO
402 sll(iad+9) = nint
403C
404 n_mult=n_mult+1
405 iadll(n_mult+1)=iadll(n_mult) + 9
406 iad = iadll(n_mult) - 1
407 DO ik=1,9
408 lll(iad+ik) = iii(i,ik)
409 jll(iad+ik) = 3
410 sll(iad+ik) = 0
411 xll(iad+ik) = ni(i,ik)
412 nn = lll(iad+ik)
413 comntag(nn) = comntag(nn) + 1
414 ENDDO
415 sll(iad+9) = nint
416#include "lockoff.inc"
417 ENDIF
418 ENDIF
419 ENDDO
420 ELSE
421C-----------------------------------------------
422C ITIED = 2
423C-----------------------------------------------
424 DO i=1,llt
425C-----------------------------------------------
426C Test if contact
427C-----------------------------------------------
428 IF(r(i)>=-one.AND.s(i)>=-one.AND.t(i)>=-one.AND.
429 . r(i)<=one.AND.s(i)<=one.AND.t(i)<=one)THEN
430C
431 nk = 9
432C-----------------------------------------------
433c print *, "s = ",S(I)
434c
435#include "lockon.inc"
436 IF(n_mult+3>n_mul_mx)THEN
437#include "lockoff.inc"
438 CALL ancmsg(msgid=84,anmode=aninfo)
439 CALL arret(2)
440 ENDIF
441 IF(iadll(n_mult+1)-1+9*3>nkmax)THEN
442#include "lockoff.inc"
443 CALL ancmsg(msgid=84,anmode=aninfo)
444 CALL arret(2)
445 ENDIF
446 n_mult=n_mult+1
447 iadll(n_mult+1)=iadll(n_mult) + 9
448 iad = iadll(n_mult) - 1
449 DO ik=1,9
450 lll(iad+ik) = iii(i,ik)
451 jll(iad+ik) = 1
452 sll(iad+ik) = 0
453 xll(iad+ik) = ni(i,ik)
454 nn = lll(iad+ik)
455 comntag(nn) = comntag(nn) + 1
456 ENDDO
457 sll(iad+9) = nint
458C
459 n_mult=n_mult+1
460 iadll(n_mult+1)=iadll(n_mult) + 9
461 iad = iadll(n_mult) - 1
462 DO ik=1,9
463 lll(iad+ik) = iii(i,ik)
464 jll(iad+ik) = 2
465 sll(iad+ik) = 0
466 xll(iad+ik) = ni(i,ik)
467 nn = lll(iad+ik)
468 comntag(nn) = comntag(nn) + 1
469 ENDDO
470 sll(iad+9) = nint
471C
472 n_mult=n_mult+1
473 iadll(n_mult+1)=iadll(n_mult) + 9
474 iad = iadll(n_mult) - 1
475 DO ik=1,9
476 lll(iad+ik) = iii(i,ik)
477 jll(iad+ik) = 3
478 sll(iad+ik) = 0
479 xll(iad+ik) = ni(i,ik)
480 nn = lll(iad+ik)
481 comntag(nn) = comntag(nn) + 1
482 ENDDO
483 sll(iad+9) = nint
484C
485#include "lockoff.inc"
486 ENDIF
487 ENDDO
488 ENDIF
489c
490c print *, "r,s,t",r(1),s(1),t(1)
491C
492 RETURN
subroutine i8rst(llt, r, s, t, ni, nx, ny, nz, xx, yy, zz)
Definition i8lagm.F:506
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:895
subroutine arret(nn)
Definition arret.F:86

◆ i8ni()

subroutine i8ni ( integer llt,
rr,
ss,
tt,
ni )

Definition at line 628 of file i8lagm.F.

629C-----------------------------------------------
630C I m p l i c i t T y p e s
631C-----------------------------------------------
632#include "implicit_f.inc"
633C-----------------------------------------------
634C G l o b a l P a r a m e t e r s
635C-----------------------------------------------
636#include "mvsiz_p.inc"
637C-----------------------------------------------
638C D u m m y A r g u m e n t s
639C-----------------------------------------------
640 INTEGER LLT
641 my_real
642 . rr(mvsiz),ss(mvsiz),tt(mvsiz),ni(mvsiz,9)
643C-----------------------------------------------
644C L o c a l V a r i a b l e s
645C-----------------------------------------------
646 INTEGER I
647 my_real
648 . u_m_r,u_p_r,u_m_s,u_p_s,u_m_t,u_p_t,
649 . ums_umt,ums_upt,ups_umt,ups_upt,
650 .
651 .
652 . r05,s05,t05
653C-----------------------------------------------------------------------
654C calculation of ni
655C-----------------------------------------------------------------------
656 DO i=1,llt
657C
658 r05 = half*rr(i)
659 s05 = half*ss(i)
660 t05 = half*tt(i)
661C
662 u_m_r = half - r05
663 u_p_r = half + r05
664C
665 u_m_s = half - s05
666 u_p_s = half + s05
667C
668 u_m_t = half - t05
669 u_p_t = half + t05
670C
671 ums_umt = u_m_s * u_m_t
672 ums_upt = u_m_s * u_p_t
673 ups_umt = u_p_s * u_m_t
674 ups_upt = u_p_s * u_p_t
675C
676 ni(i,1) = u_m_r * ums_umt
677 ni(i,5) = u_m_r * ups_umt
678 ni(i,2) = u_m_r * ums_upt
679 ni(i,6) = u_m_r * ups_upt
680 ni(i,3) = u_p_r * ums_upt
681 ni(i,7) = u_p_r * ups_upt
682 ni(i,4) = u_p_r * ums_umt
683 ni(i,8) = u_p_r * ups_umt
684C------------------------------------
685 ni(i,9) = -one
686C------------------------------------
687 ENDDO
688C-----------------------------------------------
689 RETURN

◆ i8rst()

subroutine i8rst ( integer llt,
r,
s,
t,
ni,
nx,
ny,
nz,
xx,
yy,
zz )

Definition at line 504 of file i8lagm.F.

506C-----------------------------------------------
507C I m p l i c i t T y p e s
508C-----------------------------------------------
509#include "implicit_f.inc"
510C-----------------------------------------------
511C G l o b a l P a r a m e t e r s
512C-----------------------------------------------
513#include "mvsiz_p.inc"
514C-----------------------------------------------
515C D u m m y A r g u m e n t s
516C-----------------------------------------------
517 INTEGER LLT
518C REAL
519 my_real
520 . xx(mvsiz,9),yy(mvsiz,9),zz(mvsiz,9)
521 my_real
522 . r(mvsiz),s(mvsiz),t(mvsiz),ni(mvsiz,9) ,
523 . nx(mvsiz),ny(mvsiz),nz(mvsiz)
524C-----------------------------------------------
525C L o c a l V a r i a b l e s
526C-----------------------------------------------
527 INTEGER I, ITER, NITERMAX, JTER, NJTERMAX, CONV
528 my_real
529 . drdx(mvsiz),drdy(mvsiz),drdz(mvsiz),
530 . dsdx(mvsiz),dsdy(mvsiz),dsdz(mvsiz),
531 . dtdx(mvsiz),dtdy(mvsiz),dtdz(mvsiz),
532 . dxdr(mvsiz),dydr(mvsiz),dzdr(mvsiz),
533 . dxdt(mvsiz),dydt(mvsiz),dzdt(mvsiz),
534 . rr(mvsiz),ss(mvsiz),tt(mvsiz)
535C-----------------------------------------------
536C
537C r=s=t=0
538C
539C +---> iter
540C |
541C | Ni(r,s,t) =
542C | dNi/dr =
543C | ... _
544C | \
545C | dx/dr = /_ (xi * dNi/dr)
546C | ...
547C |
548C | [dx/dr dy/dr dz/dr]
549C | [J] = |dx/ds dy/ds dz/ds|
550C | [dx/dt dy/dt dz/dt]
551C |
552C | +--> jter
553C | | _
554C | | \
555C | | x(r,s,t) = /_ (xi * Ni(r,s,t))
556C | | ...
557C | |
558C | | |r| |r| -1 |xs-x(r,s,t)|
559C | | {s} = {s} + [J] {ys-y(r,s,t)}
560C | | |t| |t| |zs-z(r,s,t)|
561C | |
562C | | Ni(r,s,t) =
563C +-+---
564C-----------------------------------------------
565 nitermax = 3
566 njtermax = 3
567 conv = 0
568C
569 DO i=1,llt
570 rr(i) = zero
571 ss(i) = zero
572 tt(i) = zero
573 ENDDO
574C-----------------------------------------------
575C calculation de r,s,t et Ni(r,s,t)
576C-----------------------------------------------
577 DO iter=1,nitermax
578c
579c print *, "iter",iter
580c
581C-----------------------------------------------
582C calculation de Ni(r,s,t); [J]; [J]-1
583C-----------------------------------------------
584 CALL i8deri(llt,rr ,ss ,tt ,ni ,
585 2 drdx ,drdy ,drdz ,dsdx ,dsdy ,dsdz ,
586 3 dtdx ,dtdy ,dtdz ,dxdr ,dydr ,dzdr ,
587 4 dxdt ,dydt ,dzdt ,xx ,yy ,zz )
588C
589 DO jter=1,njtermax
590c
591c print *, "jter",jter
592c
593C-----------------------------------------------
594C calculation de r,s,t new
595C-----------------------------------------------
596 CALL i8rstn(llt,rr,ss ,tt ,ni ,conv ,
597 2 drdx ,drdy ,drdz ,dsdx ,dsdy ,dsdz ,
598 3 dtdx ,dtdy ,dtdz ,xx ,yy ,zz ,
599 4 r ,s ,t )
600c
601c print *, "r,s,t",r(1),s(1),t(1)
602c print *, "rr,ss,tt",rr(1),ss(1),tt(1)
603c
604C-----------------------------------------------
605C calculation de Ni(-1<r<1 , -1<s<1 , -1<t<1)
606C-----------------------------------------------
607 CALL i8ni(llt,rr ,ss ,tt ,ni )
608C parith on problem if conv function of mvsiz !!!!!!
609C IF(CONV/=0)RETURN
610C
611 ENDDO
612 ENDDO
613C
614 DO i=1,llt
615 nx(i) = dydt(i)*dzdr(i) - dzdt(i)*dydr(i)
616 ny(i) = dzdt(i)*dxdr(i) - dxdt(i)*dzdr(i)
617 nz(i) = dxdt(i)*dydr(i) - dydt(i)*dxdr(i)
618 ENDDO
619C
620 RETURN
subroutine i8rstn(llt, rr, ss, tt, ni, conv, drdx, drdy, drdz, dsdx, dsdy, dsdz, dtdx, dtdy, dtdz, xx, yy, zz, r, s, t)
Definition i8lagm.F:700
subroutine i8deri(llt, rr, ss, tt, ni, drdx, drdy, drdz, dsdx, dsdy, dsdz, dtdx, dtdy, dtdz, dxdr, dydr, dzdr, dxdt, dydt, dzdt, xx, yy, zz)
Definition i8lagm.F:781
subroutine i8ni(llt, rr, ss, tt, ni)
Definition i8lagm.F:629

◆ i8rstn()

subroutine i8rstn ( integer llt,
rr,
ss,
tt,
ni,
integer conv,
drdx,
drdy,
drdz,
dsdx,
dsdy,
dsdz,
dtdx,
dtdy,
dtdz,
xx,
yy,
zz,
r,
s,
t )

Definition at line 696 of file i8lagm.F.

700C-----------------------------------------------
701C I m p l i c i t T y p e s
702C-----------------------------------------------
703c#include "implicit_f.inc"
704 implicit none
705C-----------------------------------------------
706C G l o b a l P a r a m e t e r s
707C-----------------------------------------------
708#include "mvsiz_p.inc"
709#include "constant.inc"
710C-----------------------------------------------
711C D u m m y A r g u m e n t s
712C-----------------------------------------------
713 INTEGER LLT,CONV
714 my_real
715 . r(mvsiz),s(mvsiz),t(mvsiz),ni(mvsiz,9) ,
716 . rr(mvsiz),ss(mvsiz),tt(mvsiz),
717 . xx(mvsiz,9) ,yy(mvsiz,9) ,zz(mvsiz,9) ,
718 . drdx(mvsiz),drdy(mvsiz),drdz(mvsiz),
719 . dsdx(mvsiz),dsdy(mvsiz),dsdz(mvsiz),
720 . dtdx(mvsiz),dtdy(mvsiz),dtdz(mvsiz)
721C-----------------------------------------------
722C L o c a l V a r i a b l e s
723C-----------------------------------------------
724 INTEGER I
725 my_real
726 . dx ,dy,dz,dr ,ds,dt,err
727C
728 err=zero
729C-----------------------------------------------
730 DO i=1,llt
731C
732 dx = xx(i,9)
733 + - ni(i, 1)*xx(i, 1) - ni(i, 2)*xx(i, 2) - ni(i, 3)*xx(i, 3)
734 + - ni(i, 4)*xx(i, 4) - ni(i, 5)*xx(i, 5) - ni(i, 6)*xx(i, 6)
735 + - ni(i, 7)*xx(i, 7) - ni(i, 8)*xx(i, 8)
736 dy = yy(i,9)
737 + - ni(i, 1)*yy(i, 1) - ni(i, 2)*yy(i, 2) - ni(i, 3)*yy(i, 3)
738 + - ni(i, 4)*yy(i, 4) - ni(i, 5)*yy(i, 5) - ni(i, 6)*yy(i, 6)
739 + - ni(i, 7)*yy(i, 7) - ni(i, 8)*yy(i, 8)
740 dz = zz(i,9)
741 + - ni(i, 1)*zz(i, 1) - ni(i, 2)*zz(i, 2) - ni(i, 3)*zz(i, 3)
742 + - ni(i, 4)*zz(i, 4) - ni(i, 5)*zz(i, 5) - ni(i, 6)*zz(i, 6)
743 + - ni(i, 7)*zz(i, 7) - ni(i, 8)*zz(i, 8)
744C
745 dr = drdx(i)*dx + drdy(i)*dy + drdz(i)*dz
746 ds = dsdx(i)*dx + dsdy(i)*dy + dsdz(i)*dz
747 dt = dtdx(i)*dx + dtdy(i)*dy + dtdz(i)*dz
748C
749 rr(i) = rr(i) + dr
750 ss(i) = ss(i) + ds
751 tt(i) = tt(i) + dt
752C
753 r(i) = rr(i)
754 s(i) = ss(i)
755 t(i) = tt(i)
756C
757 IF(r(i)>=-one.AND.s(i)>=-one.AND.t(i)>=-one.AND.
758 . r(i)<= one.AND.s(i)<= one.AND.t(i)<= one)THEN
759 err = max(err,abs(dr),abs(ds),abs(dt))
760 ELSE
761 rr(i) = max(min(rr(i),one),-one)
762 ss(i) = max(min(ss(i),one),-one)
763 tt(i) = max(min(tt(i),one),-one)
764 ENDIF
765c
766 ENDDO
767C
768 IF(err<em4) conv = 1
769C-----------------------------------------------
770 RETURN
#define min(a, b)
Definition macros.h:20