OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i20for3.F
Go to the documentation of this file.
1Copyright> OpenRadioss
2Copyright> Copyright (C) 1986-2025 Altair Engineering Inc.
3Copyright>
4Copyright> This program is free software: you can redistribute it and/or modify
5Copyright> it under the terms of the GNU Affero General Public License as published by
6Copyright> the Free Software Foundation, either version 3 of the License, or
7Copyright> (at your option) any later version.
8Copyright>
9Copyright> This program is distributed in the hope that it will be useful,
10Copyright> but WITHOUT ANY WARRANTY; without even the implied warranty of
11Copyright> MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12Copyright> GNU Affero General Public License for more details.
13Copyright>
14Copyright> You should have received a copy of the GNU Affero General Public License
15Copyright> along with this program. If not, see <https://www.gnu.org/licenses/>.
16Copyright>
17Copyright>
18Copyright> Commercial Alternative: Altair Radioss Software
19Copyright>
20Copyright> As an alternative to this open-source version, Altair also offers Altair Radioss
21Copyright> software under a commercial license. Contact Altair to discuss further if the
22Copyright> commercial version may interest you: https://www.altair.com/radioss/.
23!||====================================================================
24!|| i20for3 ../engine/source/interfaces/int20/i20for3.F
25!||--- called by ------------------------------------------------------
26!|| i20mainf ../engine/source/interfaces/int20/i20mainf.F
27!||--- calls -----------------------------------------------------
28!|| foat_to_6_float ../engine/source/system/parit.F
29!|| i7ass0 ../engine/source/interfaces/int07/i7ass3.F
30!|| i7ass05 ../engine/source/interfaces/int07/i7ass3.F
31!|| i7ass2 ../engine/source/interfaces/int07/i7ass3.F
32!|| i7ass25 ../engine/source/interfaces/int07/i7ass3.F
33!|| i7ass3 ../engine/source/interfaces/int07/i7ass3.F
34!|| i7ass35 ../engine/source/interfaces/int07/i7ass3.F
35!|| i7curv ../engine/source/interfaces/int07/i7curv.F
36!|| i7sms2 ../engine/source/interfaces/int07/i7sms2.F
37!|| ibcoff ../engine/source/interfaces/interf/ibcoff.F
38!||--- uses -----------------------------------------------------
39!|| anim_mod ../common_source/modules/output/anim_mod.F
40!|| h3d_mod ../engine/share/modules/h3d_mod.F
41!|| tri7box ../engine/share/modules/tri7box.F
42!||====================================================================
43 SUBROUTINE i20for3(JLT ,A ,VA ,IBCC ,ICODT ,
44 2 FSAV ,GAP ,FRIC ,MS ,VISC ,
45 3 VISCF ,NOINT ,STFA ,ITAB ,CN_LOC ,
46 4 STIGLO ,STIFN ,STIF ,FSKYI ,ISKY ,
47 6 NX1 ,NX2 ,NX3 ,NX4 ,NY1 ,
48 7 NY2 ,NY3 ,NY4 ,NZ1 ,NZ2 ,
49 8 NZ3 ,NZ4 ,LB1 ,LB2 ,LB3 ,
50 9 LB4 ,LC1 ,LC2 ,LC3 ,LC4 ,
51 A P1 ,P2 ,P3 ,P4 ,FCONT ,
52 B IX1L ,IX2L ,IX3L ,IX4L ,NSVG ,
53 C IVIS2 ,NELTST ,ITYPTST,DT2T ,
54 D GAPV ,INACTI ,INDEX ,NISKYFI,
55 E KINET ,NEWFRONT,ISECIN,NSTRF ,SECFCUM,
56 F X ,XA ,CE_LOC ,MFROT ,IFQ ,
57 G FROT_P ,CAND_FX,CAND_FY,CAND_FZ,ALPHA0,
58 H IFPEN ,GAPR ,DXANC ,NLN ,NLG ,
59 I IBAG ,ICONTACT,NSV ,PENIS ,PENIM ,
60 J VISCN ,VXI ,VYI ,VZI ,MSI ,
61 K KINI ,NIN ,NISUB ,LISUB ,ADDSUBS,
62 L ADDSUBM,LISUBS ,LISUBM ,FSAVSUB,CAND_N ,
63 M ILAGM ,ICURV ,NOD_NORMAL ,FNCONT ,FTCONT ,
64 N X1 ,X2 ,X3 ,X4 ,Y1 ,
65 O Y2 ,Y3 ,Y4 ,Z1 ,Z2 ,
66 P Z3 ,Z4 ,XI ,YI ,ZI ,
67 Q IADM ,RCURVI ,RCONTACT,ACONTACT,PCONTACT,
68 R ANGLMI ,PADM ,INTTH ,PHI , FTHE ,
69 S FTHESKYI,DAANC6,TEMP ,TEMPI ,RSTIF ,
70 T IFORM ,GAP_S ,IGAP ,ALPHAK ,MSKYI_SMS,
71 U ISKYI_SMS,NSMS ,CMAJ ,JTASK,ISENSINT,
72 V FSAVPARIT ,NFT ,H3D_DATA)
73C-----------------------------------------------
74C M o d u l e s
75C-----------------------------------------------
76 USE tri7box
77 USE h3d_mod
78 USE anim_mod
79C-----------------------------------------------
80C I m p l i c i t T y p e s
81C-----------------------------------------------
82#include "implicit_f.inc"
83#include "comlock.inc"
84C-----------------------------------------------
85C G l o b a l P a r a m e t e r s
86C-----------------------------------------------
87#include "mvsiz_p.inc"
88C-----------------------------------------------
89C C o m m o n B l o c k s
90C-----------------------------------------------
91#include "com01_c.inc"
92#include "com04_c.inc"
93#include "com06_c.inc"
94#include "com08_c.inc"
95#include "scr05_c.inc"
96#include "scr07_c.inc"
97#include "scr11_c.inc"
98#include "scr14_c.inc"
99#include "scr16_c.inc"
100#include "scr18_c.inc"
101#include "units_c.inc"
102#include "parit_c.inc"
103#include "param_c.inc"
104#include "impl1_c.inc"
105#include "sms_c.inc"
106#include "kincod_c.inc"
107C-----------------------------------------------
108C D u m m y A r g u m e n t s
109C-----------------------------------------------
110 INTEGER NELTST,ITYPTST,JLT,IBCC,IBCM,IBCS,IVIS2,INACTI,IBAG,NIN,
111 . NTY ,NLN,NLG(NLN),NSV(*),
112 . ICODT(*), ITAB(*), ISKY(*), KINET(*),
113 . MFROT, IFQ, NOINT,NEWFRONT,ISECIN, NSTRF(*),
114 . IFPEN(*) ,ICONTACT(*), CAND_N(*),
115 . KINI(*),
116 . ISET, NISKYFI,IADM,INTTH,IFORM, IGAP,JTASK
117 INTEGER IX1L(MVSIZ), IX2L(MVSIZ), IX3L(MVSIZ), IX4L(MVSIZ),
118 . CN_LOC(MVSIZ),CE_LOC(MVSIZ),INDEX(MVSIZ),NSVG(MVSIZ),
119 . NISUB, LISUB(*), ADDSUBS(*), ADDSUBM(*), LISUBS(*),
120 . LISUBM(*),ILAGM,ICURV(3),
121 . ISKYI_SMS(*), NSMS(*), ISENSINT(*),NFT
122 my_real
123 . STIGLO,FROT_P(*), X(3,*), XA(3,*),DXANC(3,*),
124 . A(3,*), MS(*), VA(3,*), FSAV(*),FCONT(3,*),
125 . CAND_FX(*),CAND_FY(*),CAND_FZ(*),ALPHA0,
126 . GAP, FRIC,VISC,VISCF,VIS,DT2T,STFA(*),STIFN(*),
127 . FSKYI(LSKYI,NFSKYI),FSAVSUB(NTHVKI,*), FNCONT(3,*),FTCONT(3,*),
128 . MSKYI_SMS(*)
129 my_real
130 . NX1(MVSIZ), NX2(MVSIZ), NX3(MVSIZ), NX4(MVSIZ),
131 . NY1(MVSIZ), NY2(MVSIZ), NY3(MVSIZ), NY4(MVSIZ),
132 . NZ1(MVSIZ), NZ2(MVSIZ), NZ3(MVSIZ), NZ4(MVSIZ),
133 . LB1(MVSIZ), LB2(MVSIZ), LB3(MVSIZ), LB4(MVSIZ),
134 . LC1(MVSIZ), LC2(MVSIZ), LC3(MVSIZ), LC4(MVSIZ),
135 . P1(MVSIZ), P2(MVSIZ), P3(MVSIZ), P4(MVSIZ), STIF(MVSIZ),
136 . GAPV(MVSIZ),GAPR(MVSIZ),SECFCUM(7,NUMNOD,NSECT),
137 . TMP(MVSIZ),STIFSAV(MVSIZ), VISCN(*),
138 . VXI(MVSIZ),VYI(MVSIZ),VZI(MVSIZ),MSI(MVSIZ),
139 . x1(mvsiz),y1(mvsiz),z1(mvsiz),
140 . x2(mvsiz),y2(mvsiz),z2(mvsiz),
141 . x3(mvsiz),y3(mvsiz),z3(mvsiz),
142 . x4(mvsiz),y4(mvsiz),z4(mvsiz),
143 . xi(mvsiz),yi(mvsiz),zi(mvsiz),penis(2,*),penim(2,*),
144 . phi(mvsiz), fthe(*),ftheskyi(lskyi),temp(*), tempi(mvsiz),
145 . rstif,fsavparit(nisub+1,11,*)
146 my_real
147 . nod_normal(3,*), rcurvi(*), rcontact(*), acontact(*),
148 . pcontact(*),padm, anglmi(*),gap_s(*),alphak(3,*),cmaj(mvsiz)
149 double precision
150 . daanc6(3,6,*)
151 TYPE(h3d_database) :: H3D_DATA
152C-----------------------------------------------
153C L o c a l V a r i a b l e s
154C-----------------------------------------------
155 INTEGER I,J1,IG,J,JG,IM,IS,K0,NBINTER,K1S,K,IL,IE,NN,NI,NA1,NA2,
156 . JSUB,KSUB,JJ,KK,IN,NSUB,ISIGN,IPROJ,IBID
157 INTEGER IX1G(MVSIZ), IX2G(MVSIZ), IX3G(MVSIZ), IX4G(MVSIZ)
158 my_real
159 . fxr(mvsiz), fyr(mvsiz), fzr(mvsiz),
160 . fxi(mvsiz), fyi(mvsiz), fzi(mvsiz), fni(mvsiz),
161 . fxt(mvsiz),fyt(mvsiz),fzt(mvsiz),
162 . fx1(mvsiz), fx2(mvsiz), fx3(mvsiz), fx4(mvsiz),
163 . fy1(mvsiz), fy2(mvsiz), fy3(mvsiz), fy4(mvsiz),
164 . fz1(mvsiz), fz2(mvsiz), fz3(mvsiz), fz4(mvsiz),
165 . n1(mvsiz), n2(mvsiz), n3(mvsiz), pene(mvsiz),
166 . vis2(mvsiz), dtmi(mvsiz), xmu(mvsiz),stif0(mvsiz),
167 . h1(mvsiz), h2(mvsiz), h3(mvsiz), h4(mvsiz),
168 . vx(mvsiz), vy(mvsiz), vz(mvsiz), vn(mvsiz),
169 . st1(mvsiz),st2(mvsiz),st3(mvsiz),st4(mvsiz),stv(mvsiz),
170 . kt(mvsiz),c(mvsiz),cf(mvsiz),
171 . ks(mvsiz),k1(mvsiz),k2(mvsiz),k3(mvsiz),k4(mvsiz),
172 . cs(mvsiz),c1(mvsiz),c2(mvsiz),c3(mvsiz),c4(mvsiz),
173 . p1s(mvsiz),p2s(mvsiz),p3s(mvsiz),p4s(mvsiz),
174 . phi1(mvsiz),phi2(mvsiz),phi3(mvsiz),phi4(mvsiz),
175 . fsavsub1(15,nisub),masm(mvsiz)
176 my_real
177 . vnx, vny, vnz, aa, crit,s2,dist,rdist,
178 . v2, fm2, dt1inv, visca, fac,ff,alphi,alpha,beta,
179 . fx, fy, fz, f2, mas2, m2sk, dtmi0,dti,ft,fn,fmax,ftn,
180 . facm1, econtt, econvt, h0, la1, la2, la3, la4,
181 . d1,d2,d3,d4,a1,a2,a3,a4,e10, h0d, s2d, sum,
182 . fsav1, fsav2, fsav3, fsav4, fsav5, fsav6, fsav7, fsav8,
183 . fsav9, fsav10, fsav11, fsav12, fsav13, fsav14, fsav15, ffo,
184 . la1d,la2d,la3d,la4d,t1,t1d,t2,t2d,ffd,visd,facd,d1d,
185 . d2d,d3d,d4d,vnxd,vnyd,vnzd,v2d,fm2d,f2d,aad,fxd,fyd,fzd,
186 . a1d,a2d,a3d,a4d,vv,ax1,ax2,ay1,ay2,az1,az2,ax,ay,az,
187 . area,p,vv1,vv2,v21,dmu, dti2,h00 ,a0x,a0y,a0z,rx,ry,rz,
188 . anx,any,anz,aan,aax,aay,aaz ,rr,rs,aaa,stfr,visr,
189 . prec,ps,xsa,pis,pplus,cx,cy,cfi,aux,tm,ts,impx,impy,impz,bb,
190 . nn1,nn2,nn3,nn4,xn1,yn1,zn1,xn2,yn2,zn2,xn3,yn3,zn3,xn4,yn4,
191 . zn4,dtmini,bid
192C
193 DOUBLE PRECISION FX6(6,MVSIZ), FY6(6,MVSIZ), FZ6(6,MVSIZ)
194C
195C-----------------------------------------------
196 IF (IRESP==1) then
197 prec = fiveem4
198 ELSE
199 prec = em10
200 ENDIF
201 IF(dt1>zero)THEN
202 dt1inv = one/dt1
203 ELSE
204 dt1inv =zero
205 ENDIF
206 econtt = zero
207 econvt = zero
208 DO i=1,jlt
209 stif0(i) = stif(i)
210 ix1g(i) = nlg(ix1l(i))
211 ix2g(i) = nlg(ix2l(i))
212 ix3g(i) = nlg(ix3l(i))
213 ix4g(i) = nlg(ix4l(i))
214 ENDDO
215C--------------------------------------------------------
216C UNIQUEMENT POUR PAQUET DE QUADRANGLE
217C--------------------------------------------------------
218C--------------------------------------------------------
219C CAS DES PAQUETS MIXTES
220C--------------------------------------------------------
221 IF(icurv(1) == 3) THEN
222 DO i=1,jlt
223C
224 bb = gapv(i)+cmaj(i)
225C
226 d1 = sqrt(p1(i))
227 p1(i) = max(zero, bb - d1)
228C
229 d2 = sqrt(p2(i))
230 p2(i) = max(zero, bb - d2)
231C
232 d3 = sqrt(p3(i))
233 p3(i) = max(zero, bb - d3)
234C
235 d4 = sqrt(p4(i))
236 p4(i) = max(zero, bb - d4)
237C
238 a1 = p1(i)/max(em20,d1)
239 a2 = p2(i)/max(em20,d2)
240 a3 = p3(i)/max(em20,d3)
241 a4 = p4(i)/max(em20,d4)
242 n1(i) = a1*nx1(i) + a2*nx2(i) + a3*nx3(i) + a4*nx4(i)
243 n2(i) = a1*ny1(i) + a2*ny2(i) + a3*ny3(i) + a4*ny4(i)
244 n3(i) = a1*nz1(i) + a2*nz2(i) + a3*nz3(i) + a4*nz4(i)
245 ENDDO
246 ELSE
247 DO i=1,jlt
248C
249 d1 = sqrt(p1(i))
250 p1(i) = max(zero, gapv(i) - d1)
251C
252 d2 = sqrt(p2(i))
253 p2(i) = max(zero, gapv(i) - d2)
254C
255 d3 = sqrt(p3(i))
256 p3(i) = max(zero, gapv(i) - d3)
257C
258 d4 = sqrt(p4(i))
259 p4(i) = max(zero, gapv(i) - d4)
260C
261 a1 = p1(i)/max(em20,d1)
262 a2 = p2(i)/max(em20,d2)
263 a3 = p3(i)/max(em20,d3)
264 a4 = p4(i)/max(em20,d4)
265 n1(i) = a1*nx1(i) + a2*nx2(i) + a3*nx3(i) + a4*nx4(i)
266 n2(i) = a1*ny1(i) + a2*ny2(i) + a3*ny3(i) + a4*ny4(i)
267 n3(i) = a1*nz1(i) + a2*nz2(i) + a3*nz3(i) + a4*nz4(i)
268 ENDDO
269 ENDIF
270C
271 DO i=1,jlt
272 IF(ix3g(i)/=ix4g(i))THEN
273 pene(i) = max(p1(i),p2(i),p3(i),p4(i))
274
275 la1 = one - lb1(i) - lc1(i)
276 la2 = one - lb2(i) - lc2(i)
277 la3 = one - lb3(i) - lc3(i)
278 la4 = one - lb4(i) - lc4(i)
279
280 h0 = fourth *
281 . (p1(i)*la1 + p2(i)*la2 + p3(i)*la3 + p4(i)*la4)
282 h1(i) = h0 + p1(i) * lb1(i) + p4(i) * lc4(i)
283 h2(i) = h0 + p2(i) * lb2(i) + p1(i) * lc1(i)
284 h3(i) = h0 + p3(i) * lb3(i) + p2(i) * lc2(i)
285 h4(i) = h0 + p4(i) * lb4(i) + p3(i) * lc3(i)
286
287 h00 = one/max(em20,h1(i) + h2(i) + h3(i) + h4(i))
288 h1(i) = h1(i) * h00
289 h2(i) = h2(i) * h00
290 h3(i) = h3(i) * h00
291 h4(i) = h4(i) * h00
292C
293 ELSE
294 pene(i) = p1(i)
295 n1(i) = nx1(i)
296 n2(i) = ny1(i)
297 n3(i) = nz1(i)
298 h1(i) = lb1(i)
299 h2(i) = lc1(i)
300 h3(i) = one - lb1(i) - lc1(i)
301 h4(i) = zero
302 ENDIF
303 ENDDO
304C---------------------
305C COURBURE FIXE
306C---------------------
307 IF(icurv(1)==1)THEN
308C spherique (que concave pour le moment)
309 na1 = icurv(2)
310 DO i=1,jlt
311 rr = 1.e30
312 a0x = xa(1,na1)
313 a0y = xa(2,na1)
314 a0z = xa(3,na1)
315C
316 rx = x1(i)-a0x
317 ry = y1(i)-a0y
318 rz = z1(i)-a0z
319 rr = min(rr , rx*rx + ry*ry + rz*rz)
320 rx = x2(i)-a0x
321 ry = y2(i)-a0y
322 rz = z2(i)-a0z
323 rr = min(rr , rx*rx + ry*ry + rz*rz)
324 rx = x3(i)-a0x
325 ry = y3(i)-a0y
326 rz = z3(i)-a0z
327 rr = min(rr , rx*rx + ry*ry + rz*rz)
328 IF(ix3g(i)/=ix4g(i))THEN
329 rx = x4(i)-a0x
330 ry = y4(i)-a0y
331 rz = z4(i)-a0z
332 rr = min(rr , rx*rx + ry*ry + rz*rz)
333 ENDIF
334 rx = xi(i)-a0x
335 ry = yi(i)-a0y
336 rz = zi(i)-a0z
337 rs = sqrt(rx*rx + ry*ry + rz*rz)
338 rr = sqrt(rr)
339 IF(rs-rr+gapv(i)<0.)THEN
340 stif(i) = 0.
341 pene(i) = 0.
342 ELSEIF(rs-rr+gapv(i)<pene(i))THEN
343 pene(i) = rs-rr+gapv(i)
344 ENDIF
345 n1(i) = -rx
346 n2(i) = -ry
347 n3(i) = -rz
348 ENDDO
349 ELSEIF(icurv(1)==2)THEN
350C cylindrique (que concave pour le moment)
351 na1 = icurv(2)
352 na2 = icurv(3)
353 DO i=1,jlt
354 rr = 1.e30
355 a0x = xa(1,na1)
356 a0y = xa(2,na1)
357 a0z = xa(3,na1)
358 anx = xa(1,na2)-a0x
359 any = xa(2,na2)-a0y
360 anz = xa(3,na2)-a0z
361 aan = 1. / (anx*anx + any*any + anz*anz)
362
363 aax = x1(i)-a0x
364 aay = y1(i)-a0y
365 aaz = z1(i)-a0z
366 aaa = (aax*anx + aay*any + aaz*anz) * aan
367 rx = aax - aaa * anx
368 ry = aay - aaa * any
369 rz = aaz - aaa * anz
370 rr = min(rr , rx*rx + ry*ry + rz*rz)
371
372 aax = x2(i)-a0x
373 aay = y2(i)-a0y
374 aaz = z2(i)-a0z
375 aaa = (aax*anx + aay*any + aaz*anz) * aan
376 rx = aax - aaa * anx
377 ry = aay - aaa * any
378 rz = aaz - aaa * anz
379 rr = min(rr , rx*rx + ry*ry + rz*rz)
380
381 aax = x3(i)-a0x
382 aay = y3(i)-a0y
383 aaz = z3(i)-a0z
384 aaa = (aax*anx + aay*any + aaz*anz) * aan
385 rx = aax - aaa * anx
386 ry = aay - aaa * any
387 rz = aaz - aaa * anz
388 rr = min(rr , rx*rx + ry*ry + rz*rz)
389 IF(ix3g(i)/=ix4g(i))THEN
390
391 aax = x4(i)-a0x
392 aay = y4(i)-a0y
393 aaz = z4(i)-a0z
394 aaa = (aax*anx + aay*any + aaz*anz) * aan
395 rx = aax - aaa * anx
396 ry = aay - aaa * any
397 rz = aaz - aaa * anz
398 rr = min(rr , rx*rx + ry*ry + rz*rz)
399 ENDIF
400 aax = xi(i)-a0x
401 aay = yi(i)-a0y
402 aaz = zi(i)-a0z
403
404 aaa = (aax*anx + aay*any + aaz*anz) * aan
405 rx = aax - aaa * anx
406 ry = aay - aaa * any
407 rz = aaz - aaa * anz
408 rs = sqrt(rx*rx + ry*ry + rz*rz)
409 rr = sqrt(rr)
410 IF(rs-rr+gapv(i)<0.)THEN
411 stif(i) = 0.
412 pene(i) = 0.
413 ELSEIF(rs-rr+gapv(i)<pene(i))THEN
414 pene(i) = rs-rr+gapv(i)
415 n1(i) = -rx
416 n2(i) = -ry
417 n3(i) = -rz
418 ELSEIF(rs-rr-gapv(i)>0.)THEN
419 stif(i) = 0.
420 pene(i) = 0.
421 ELSEIF(rs-rr-gapv(i) < pene(i))THEN
422 xn1 = x1(i) - xi(i)
423 yn1 = y1(i) - yi(i)
424 zn1 = z1(i) - zi(i)
425 xn2 = x2(i) - xi(i)
426 yn2 = y2(i) - yi(i)
427 zn2 = z2(i) - zi(i)
428 xn3 = x3(i) - xi(i)
429 yn3 = y3(i) - yi(i)
430 zn3 = z3(i) - zi(i)
431C --
432 nn1 = (yn1*zn2-yn2*zn1) * rx +
433 . (zn1*xn2-zn2*xn1) * ry +
434 . (xn1*yn2-xn2*yn1) * rz
435 nn2 = (yn2*zn3-yn3*zn2) * rx +
436 . (zn2*xn3-zn3*xn2) * ry +
437 . (xn2*yn3-xn3*yn2) * rz
438 nn3 = (yn3*zn4-yn4*zn3) * rx +
439 . (zn3*xn4-zn4*xn3) * ry +
440 . (xn3*yn4-xn4*yn3) * rz
441 IF(ix3l(i)/=ix4l(i))THEN
442 xn4 = x4(i) - xi(i)
443 yn4 = y4(i) - yi(i)
444 zn4 = z4(i) - zi(i)
445 nn4 = (yn4*zn1-yn1*zn4) * rx +
446 . (zn4*xn1-zn1*xn4) * ry +
447 . (xn4*yn1-xn1*yn4) * rz
448 ELSE
449 xn4 = zero
450 yn4 = zero
451 zn4 = zero
452 nn4=zero
453 ENDIF
454 IF( nn1>=zero .AND. nn2>=zero
455 . .AND. nn3>=zero .AND. nn4>=zero) THEN
456 iproj = 1
457 ELSEIF( nn1<=zero .AND. nn2<=zero
458 . .AND. nn3<=zero .AND. nn4<=zero) THEN
459 iproj = 1
460 ELSE
461 iproj = 0
462 ENDIF
463C --
464 IF(iproj == 1)THEN
465 pene(i) = -rs+rr+gapv(i)
466 n1(i) = rx
467 n2(i) = ry
468 n3(i) = rz
469 ENDIF
470 ENDIF
471 ENDDO
472
473 ELSEIF(icurv(1) == 3)THEN
474 CALL i7curv(jlt ,pene ,n1 ,n2 ,
475 1 n3 ,gapv ,xa ,nod_normal,
476 2 ix1l ,ix2l ,ix3l ,ix4l ,
477 3 h1 ,h2 ,h3 ,h4 ,
478 4 x1 ,x2 ,x3 ,x4 ,y1 ,
479 5 y2 ,y3 ,y4 ,z1 ,z2 ,
480 6 z3 ,z4 ,xi ,yi ,zi )
481
482 DO i=1,jlt
483 IF(pene(i)<zero)THEN
484 stif(i) =zero
485 pene(i) =zero
486 END IF
487 END DO
488 ENDIF
489
490 DO i=1,jlt
491 s2 = one/max(em30,sqrt(n1(i)**2 + n2(i)**2 + n3(i)**2))
492 n1(i) = n1(i)*s2
493 n2(i) = n2(i)*s2
494 n3(i) = n3(i)*s2
495 ENDDO
496C
497 DO i=1,jlt
498 vx(i) = vxi(i) - h1(i)*va(1,ix1l(i)) - h2(i)*va(1,ix2l(i))
499 . - h3(i)*va(1,ix3l(i)) - h4(i)*va(1,ix4l(i))
500 vy(i) = vyi(i) - h1(i)*va(2,ix1l(i)) - h2(i)*va(2,ix2l(i))
501 . - h3(i)*va(2,ix3l(i)) - h4(i)*va(2,ix4l(i))
502 vz(i) = vzi(i) - h1(i)*va(3,ix1l(i)) - h2(i)*va(3,ix2l(i))
503 . - h3(i)*va(3,ix3l(i)) - h4(i)*va(3,ix4l(i))
504 vn(i) = n1(i)*vx(i) + n2(i)*vy(i) + n3(i)*vz(i)
505 ENDDO
506
507 DO i=1,jlt
508C correction hourglass
509 h0 = -.25*(h1(i) - h2(i) + h3(i) - h4(i))
510 h0 = min(h0,h2(i),h4(i))
511 h0 = max(h0,-h1(i),-h3(i))
512 IF(ix3g(i)==ix4g(i))h0 = zero
513 h1(i) = h1(i) + h0
514 h2(i) = h2(i) - h0
515 h3(i) = h3(i) + h0
516 h4(i) = h4(i) - h0
517 ENDDO
518C---------------------
519C PENE INITIALE
520C---------------------
521 IF(inacti==5.or.inacti==6)THEN
522c DO I=1,JLT
523cC REDUCTION DE LA PENE INITIALE
524cC CAND_P(INDEX(I))=MIN(CAND_P(INDEX(I)),PENE(I))
525c CAND_P(INDEX(I))=MIN(CAND_P(INDEX(I)),
526c . ( (ONE-FIVEEM2)*CAND_P(INDEX(I))
527c . +FIVEEM2*(PENE(I)+FIVEEM2*(GAPV(I)-PENE(I)))) )
528cC SOUSTRACTION DE LA PENE INITIALE A LA PENE ET AU GAP
529c PENE(I)=MAX(ZERO,PENE(I)-CAND_P(INDEX(I)))
530c IF( PENE(I)==ZERO ) STIF(I) = ZERO
531c GAPV(I)=GAPV(I)-CAND_P(INDEX(I))
532c ENDDO
533#include "lockon.inc"
534C---
535 IF(igap > 0)THEN
536 DO i=1,jlt
537 is = cn_loc(i)
538 im = ce_loc(i)
539 nn = nsvg(i)
540 pplus = pene(i) + zep05*(gapv(i)-pene(i))
541 IF(nn > 0) THEN
542 IF (pplus < gap_s(is)) THEN
543 penis(2,is) = max(penis(2,is),pplus)
544 ELSE
545 penis(2,is) = max(penis(2,is),gap_s(is))
546 penim(2,im) = max(penim(2,im),pplus-gap_s(is))
547 END IF
548 ELSE
549 IF (pplus < gapfi(nin)%P(-nn)) THEN
550 penfi(nin)%P(2,-nn) = max(penfi(nin)%P(2,-nn),pplus)
551 ELSE
552 penfi(nin)%P(2,-nn) = max(penfi(nin)%P(2,-nn),
553 + gapfi(nin)%P(-nn))
554 penim(2,im) = max(penim(2,im),pplus-gapfi(nin)%P(-nn))
555 END IF
556 ENDIF
557 ENDDO
558 ELSE
559 DO i=1,jlt
560 im = ce_loc(i)
561 pplus = pene(i) + zep05*(gapv(i)-pene(i))
562 penim(2,im) = max(penim(2,im),pplus)
563 ENDDO
564 END IF
565C---
566c DO I=1,JLT
567c AAA = GAP_S(IS)/GAPV(I)
568c PPLUS=(PENE(I)+ZEP05*(GAPV(I)-PENE(I)))
569c NN = NSVG(I)
570c IF(NN > 0) THEN
571c PENIS(2,CN_LOC(I)) = MAX(PENIS(2,CN_LOC(I)),AAA*PPLUS)
572c ELSE
573c PENFI(NIN)%P(2,-NN) = MAX(PENFI(NIN)%P(2,-NN),AAA*PPLUS)
574c END IF
575c PENIM(2,CE_LOC(I)) = MAX(PENIM(2,CE_LOC(I)),(ONE-AAA)*PPLUS)
576c ENDDO
577C---
578#include "lockoff.inc"
579 DO i=1,jlt
580 is = cn_loc(i)
581 im = ce_loc(i)
582 nn = nsvg(i)
583 IF(nn > 0) THEN
584 pis = penis(1,is)
585 ELSE
586 pis = penfi(nin)%P(1,-nn)
587 END IF
588 pene(i) = pene(i) - pis - penim(1,im)
589 pene(i) = max(pene(i),zero)
590 IF (pene(i) == zero )stif(i)=zero
591 gapv(i) = gapv(i) - pis - penim(1,im)
592 END DO
593 ENDIF
594C---------------------
595C
596 dti = 1.e20
597C
598 DO 600 i=1,jlt
599 dist=gapv(i)-pene(i)
600 rdist = half*dist / max(em30,-vn(i))
601 dti = min(rdist,dti)
602 600 CONTINUE
603C intermediate variable coming from starter input deck interface cards
604C not read for this type of interface
605 dtmini=ep20
606C
607 IF(dti<=dtmin1(10))THEN
608 dti = 1.e20
609 DO i=1,jlt
610 dist=gapv(i)-pene(i)
611 dti2 = half*dist / max(em30,-vn(i))
612 IF(dti2<=dtmin1(10))THEN
613#include "lockon.inc"
614 WRITE(iout,'(A,E12.4,A,I10)')
615 . ' **WARNING MINIMUM TIME STEP ',dti2,
616 . ' IN INTERFACE ',noint
617 nn = nsvg(i)
618 IF(nn>0)THEN
619 ni = itab(nn)
620 ELSE
621 ni = itafi(nin)%P(-nn)
622 ENDIF
623#include "lockoff.inc"
624 IF(idtmin(10)==1)THEN
625#include "lockon.inc"
626 WRITE(iout,'(A,I10)') ' SECONDARY NODE : ',ni
627 WRITE(iout,'(A,4I10)')' MAIN NODES : ',
628 . itab(ix1g(i)),itab(ix2g(i)),itab(ix3g(i)),itab(ix4g(i))
629#include "lockoff.inc"
630 tstop = tt
631 ELSEIF(idtmin(10)==2)THEN
632#include "lockon.inc"
633 WRITE(iout,'(A,I10,A,I10)')' REMOVE SECONDARY NODE ',
634 . ni,' FROM INTERFACE ',noint
635 IF(nn>0) THEN
636 stfa(nsv(cn_loc(i))) = -abs(stfa(nsv(cn_loc(i))))
637 ELSE
638 stifi(nin)%P(-nn) = -abs(stifi(nin)%P(-nn))
639 ENDIF
640#include "lockoff.inc"
641 stif(i) = zero
642 newfront = -1
643 dti = dtmin1(10)
644 ELSEIF(idtmin(10)==5)THEN
645#include "lockon.inc"
646 WRITE(iout,'(A,I10)') ' SECONDARY NODE : ',ni
647 WRITE(iout,'(A,4I10)')' MAIN NODES : ',
648 . itab(ix1g(i)),itab(ix2g(i)),itab(ix3g(i)),itab(ix4g(i))
649#include "lockoff.inc"
650 mstop = 2
651 ELSEIF(idtmin(10)==6.AND.ilagm==2)THEN
652 ig=nsvg(i)
653 IF(kinet(ig)+kinet(ix1g(i))+kinet(ix2g(i))
654 . +kinet(ix3g(i))+kinet(ix4g(i))==0)THEN
655 cand_n(index(i)) = -iabs(cand_n(index(i)))
656 stif(i) = zero
657 dti2 = 1.e20
658#include "lockon.inc"
659 WRITE(iout,'(A,I10)') ' SECONDARY NODE : ',itab(nsvg(i))
660 WRITE(iout,'(A,4I10)')' MAIN NODES : ',
661 . itab(ix1g(i)),itab(ix2g(i)),itab(ix3g(i)),itab(ix4g(i))
662#include "lockoff.inc"
663 ENDIF
664 dti = min(dti2,dti)
665 ENDIF
666 ENDIF
667 ENDDO
668 ENDIF
669 IF(dti<dt2t)THEN
670 dt2t = dti
671 neltst = noint
672 ityptst = 10
673 ENDIF
674C-------------------------------------------
675 IF(impl_s>0)THEN
676 IF(imp_int7==2)THEN
677 DO i=1,jlt
678 IF(stiglo<=zero)THEN
679 stif(i) = half*stif(i)
680 ELSEIF(stif(i)/=zero)THEN
681 stif(i) = stiglo
682 ENDIF
683 fni(i)= -stif(i) * pene(i)
684 ENDDO
685 ELSE
686 DO i=1,jlt
687 fac = gapv(i)/max( em10,( gapv(i)-pene(i) ) )
688 facm1 = 1./fac
689 IF( (gapv(i)-pene(i))/gapv(i) <prec .AND.
690 . stif(i)>0. ) THEN
691 stif(i) = 0.
692 newfront = -1
693#include "lockon.inc"
694 nn = nsvg(i)
695 IF(nn>0)THEN
696 ni = itab(nn)
697 stfa(nsv(cn_loc(i))) = -abs(stfa(nsv(cn_loc(i))))
698 ELSE
699 ni = itafi(nin)%P(-nn)
700 stifi(nin)%P(-nn) = -abs(stifi(nin)%P(-nn))
701 ENDIF
702 WRITE(istdo,'(A,I10)')' WARNING INTERFACE ',noint
703 WRITE(istdo,'(A,I10,A)')' NODE ',ni,
704 . ' DE-ACTIVATED FROM INTERFACE'
705 WRITE(iout ,'(A,I10)')' WARNING INTERFACE ',noint
706 WRITE(iout ,'(A,I10,A)')' NODE ',ni,
707 . ' DE-ACTIVATED FROM INTERFACE'
708#include "lockoff.inc"
709 ENDIF
710 IF(stiglo<=zero)THEN
711 econtt = econtt + half*stif(i)*gapv(i)**2 *( facm1 -
712 . one -log(facm1) )
713 stif(i) = half*stif(i) * fac
714 ELSEIF(stif(i)/=zero)THEN
715 econtt = econtt + stiglo*gapv(i)**2 *( facm1 - one -
716 . log(facm1) )
717 stif(i) = stiglo * fac
718 ENDIF
719 fni(i)= -stif(i) * pene(i)
720 ENDDO
721 ENDIF
722 ELSE ! fin impl_s>0
723 DO 100 i=1,jlt
724 fac = gapv(i)/max( em10,( gapv(i)-pene(i) ) )
725 facm1 = 1./fac
726 IF( (gapv(i)-pene(i))/gapv(i) <prec .AND.
727 . stif(i)>0. ) THEN
728 stif(i) = 0.
729 newfront = -1
730#include "lockon.inc"
731 nn = nsvg(i)
732 IF(nn>0)THEN
733 ni = itab(nn)
734 stfa(nsv(cn_loc(i))) = -abs(stfa(nsv(cn_loc(i))))
735 ELSE
736 ni = itafi(nin)%P(-nn)
737 stifi(nin)%P(-nn) = -abs(stifi(nin)%P(-nn))
738 ENDIF
739 WRITE(istdo,'(A,I10)')' WARNING INTERFACE ',noint
740 WRITE(istdo,'(A,I10,A)')' NODE ',ni,
741 . ' DE-ACTIVATED FROM INTERFACE'
742 WRITE(iout ,'(A,I10)')' WARNING INTERFACE ',noint
743 WRITE(iout ,'(A,I10,A)')' NODE ',ni,
744 . ' DE-ACTIVATED FROM INTERFACE'
745#include "lockoff.inc"
746 ENDIF
747 IF(stiglo<=zero)THEN
748 econtt = econtt + half*stif(i)*gapv(i)**2 *( facm1 - one -
749 . log(facm1) )
750 stif(i) = half*stif(i) * fac
751 ELSEIF(stif(i)/=zero)THEN
752 econtt = econtt + stiglo*gapv(i)**2 *(facm1 - one - log(facm1))
753 stif(i) = stiglo * fac
754 ENDIF
755 fni(i)= -stif(i) * pene(i)
756 100 CONTINUE
757 ENDIF
758C---------------------------------
759C DAMPING + FRIC
760C---------------------------------
761 IF(visc/=zero.OR.viscf/=zero)THEN
762 IF(ivis2==0)THEN
763C---------------------------------
764C VISC QUAD TYPE V227
765C---------------------------------
766 DO i=1,jlt
767 vis2(i) = two * stif(i) * msi(i)
768 IF(vn(i)<zero) vis2(i) = vis2(i) /
769 . ( max(em10,(gapv(i)-pene(i))/gapv(i)) )
770 ENDDO
771C---------------------------------
772 visca = zep4
773 IF(kdtint==0.AND.idtmins/=2)THEN
774 DO i=1,jlt
775 fac = stif(i) / max(em30,stif(i))
776 vis = sqrt(vis2(i))
777 ff = fac * (
778 . visc * vis +
779 . visca**2 * two* msi(i) * max(zero,-vn(i)) /
780 . max((gapv(i) - pene(i)),em10) )
781 stif(i) = stif(i) * gapv(i) / max((gapv(i) - pene(i)),em10)
782 stif(i) = stif(i) + ff * dt1inv
783 stif(i) = max(stif(i) ,fac*sqrt(viscf)*vis*dt1inv)
784 ffo = ff
785 ff = ff * vn(i)
786 fni(i) = fni(i) + ff
787 ENDDO
788 ELSE
789 DO i=1,jlt
790 fac = stif(i) / max(em30,stif(i))
791 vis = sqrt(vis2(i))
792 c(i)= fac * (
793 . visc * vis +
794 . visca**2 * two * msi(i) * max(zero,-vn(i)) /
795 . max((gapv(i) - pene(i)),em10) )
796 stif(i) = stif(i) * gapv(i) / max((gapv(i) - pene(i)),em10)
797 kt(i)= stif(i)
798 stif(i) = stif(i) + c(i) * dt1inv
799 ff = c(i) * vn(i)
800 fni(i) = fni(i) + ff
801 cf(i) = fac*sqrt(viscf)*vis
802 stif(i) = max(stif(i) ,cf(i)*dt1inv)
803 ENDDO
804 ENDIF
805 ELSEIF(ivis2==1)THEN
806C---------------------------------
807C TEST
808C---------------------------------
809 DO i=1,jlt
810 masm(i) = ms(ix1g(i))*h1(i)
811 . + ms(ix2g(i))*h2(i)
812 . + ms(ix3g(i))*h3(i)
813 . + ms(ix4g(i))*h4(i)
814 masm(i) = msi(i) * masm(i) / max(em30,msi(i)+masm(i))
815 vis2(i) = two * stif(i) * masm(i)
816 IF(vn(i)<zero) vis2(i) = vis2(i) /
817 . ( max(em10,(gapv(i)-pene(i))/gapv(i)) )
818 ENDDO
819C---------------------------------
820 visca = zep4
821 IF(kdtint==0.AND.idtmins/=2)THEN
822 DO i=1,jlt
823 fac = stif(i) / max(em30,stif(i))
824 vis = sqrt(vis2(i))
825 ff = fac * (
826 . visc * vis +
827 . visca**2 * two* masm(i) * max(zero,-vn(i)) /
828 . max((gapv(i) - pene(i)),em10) )
829 stif(i) = stif(i) * gapv(i) / max((gapv(i) - pene(i)),em10)
830 stif(i) = stif(i) + ff * dt1inv
831 stif(i) = max(stif(i) ,fac*sqrt(viscf)*vis*dt1inv)
832 ffo = ff
833 ff = ff * vn(i)
834 fni(i) = fni(i) + ff
835 ENDDO
836 ELSE
837 DO i=1,jlt
838 fac = stif(i) / max(em30,stif(i))
839 vis = sqrt(vis2(i))
840 c(i)= fac * (
841 . visc * vis +
842 . visca**2 * two * masm(i) * max(zero,-vn(i)) /
843 . max((gapv(i) - pene(i)),em10) )
844 stif(i) = stif(i) * gapv(i) / max((gapv(i) - pene(i)),em10)
845 kt(i)= stif(i)
846 stif(i) = stif(i) + c(i) * dt1inv
847 ff = c(i) * vn(i)
848 fni(i) = fni(i) + ff
849 cf(i) = fac*sqrt(viscf)*vis
850 stif(i) = max(stif(i) ,cf(i)*dt1inv)
851 ENDDO
852 ENDIF
853 ELSEIF(ivis2==2)THEN
854C---------------------------------
855C VISC QUAD TYPE
856C---------------------------------
857 DO i=1,jlt
858 vis2(i) = two* stif(i) * msi(i)
859 vis2(i) = vis2(i) /
860 . ( max(em10,(gapv(i)-pene(i))/gapv(i)) )
861 ENDDO
862 visca = half
863 DO i=1,jlt
864 fac = stif(i) / max(em30,stif(i))
865 vis = sqrt(vis2(i))
866 ff = fac * (
867 . visc * vis +
868 . visca**2 * two * msi(i) * abs(vn(i)) /
869 . max((gapv(i) - pene(i)),em10) )
870 stif(i) = stif(i) * gapv(i) / max((gapv(i) - pene(i)),em10)
871 stif(i) = stif(i) + two * ff * dt1inv
872 stif(i) = max(stif(i) ,fac*sqrt(viscf)*vis*dt1inv)
873 ff = ff * vn(i)
874 fni(i) = fni(i) + ff
875 ENDDO
876 ELSEIF(ivis2==3)THEN
877C---------------------------------
878C VISC QUAD = 0
879C---------------------------------
880 DO i=1,jlt
881 vis2(i) = two * stif(i) * msi(i)
882 ENDDO
883C---------------------------------
884 DO i=1,jlt
885 fac = stif(i) / max(em30,stif(i))
886 vis = sqrt(vis2(i))
887 ff = fac * ( visc * vis ) /
888 . max((gapv(i) - pene(i)),em10)
889 stif(i) = stif(i) * gapv(i) /
890 . max((gapv(i) - pene(i)),em10)
891 stif(i) = stif(i) + two* ff * dt1inv
892 stif(i) = max(stif(i) ,fac*sqrt(viscf)*vis*dt1inv)
893 ff = ff * vn(i)
894 fni(i) = fni(i) + ff
895 ENDDO
896 ELSEIF(ivis2==4)THEN
897C---------------------------------
898C VISC = 0
899C---------------------------------
900 DO i=1,jlt
901 vis2(i) = two* stif(i) * msi(i)
902 vis = sqrt(vis2(i))
903 stif(i) = stif(i) * gapv(i) /
904 . max((gapv(i) - pene(i)),em10)
905 stif(i) = max(stif(i) ,fac*sqrt(viscf)*vis*dt1inv)
906 ENDDO
907 ELSEIF(ivis2==5)THEN
908C---------------------------------
909C VISC = 2M/dt => pour visc < 1 , stable : dt < 2M/visc = dt
910C M=M1*M2/M1+M2 pour visc = 1 , choc elastique
911C pour visc = 0.5 , choc mou
912C---------------------------------
913 DO i=1,jlt
914 mas2 = ms(ix1g(i))*h1(i)
915 . + ms(ix2g(i))*h2(i)
916 . + ms(ix3g(i))*h3(i)
917 . + ms(ix4g(i))*h4(i)
918 vis2(i) = two* stif(i) * msi(i)
919 vis = 2. * visc * dt1inv * msi(i) * mas2 /
920 . max(em30,msi(i)+mas2)
921 stif(i) = stif(i) * gapv(i) /
922 . max((gapv(i) - pene(i)),em10)
923 stif(i) = max(stif(i) ,fac*sqrt(viscf*vis2(i))*dt1inv)
924 ff = vis * vn(i)
925 econvt = econvt + min(zero,ff-fni(i)) * vn(i) * dt1
926 fni(i) = min(fni(i),ff)
927 ENDDO
928 ELSE
929 ENDIF
930 ELSE
931 DO i=1,jlt
932 vis2(i) = zero
933 stif(i) = stif(i) * gapv(i) /
934 . max((gapv(i) - pene(i)),em10)
935 ENDDO
936 ENDIF
937C---------------------------------
938C REDUCTION RIGIDITE ANCRAGE
939C---------------------------------
940#include "lockon.inc"
941 DO i=1,jlt
942 isign=1
943 aaa = one-pene(i)/gapv(i)
944 il = ix1l(i)
945 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
946 alphak(2,il)=isign*min(abs(alphak(2,il)),aaa)
947 il = ix2l(i)
948 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
949 alphak(2,il)=isign*min(abs(alphak(2,il)),aaa)
950 il = ix3l(i)
951 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
952 alphak(2,il)=isign*min(abs(alphak(2,il)),aaa)
953 il = ix4l(i)
954 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
955 alphak(2,il)=isign*min(abs(alphak(2,il)),aaa)
956 IF(nsvg(i)>0) THEN
957 il = nsv(cn_loc(i))
958 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
959 alphak(2,il)=isign*min(abs(alphak(2,il)),aaa)
960 ELSE
961C SPMD remote SECONDARYs
962 il = - nsvg(i)
963 IF(pene(i)>zero.OR.alphakfi(nin)%P(il)<zero)isign=-1
964 alphakfi(nin)%P(il)=isign*min(abs(alphakfi(nin)%P(il)),aaa)
965 ENDIF
966 ENDDO
967#include "lockoff.inc"
968C---------------------------------
969C SAUVEGARDE DE L'IMPULSION NORMALE
970C---------------------------------
971 fsav1 = zero
972 fsav2 = zero
973 fsav3 = zero
974
975 fsav8 = zero
976 fsav9 = zero
977 fsav10= zero
978 fsav11= zero
979 DO i=1,jlt
980 fxi(i)=n1(i)*fni(i)
981 fyi(i)=n2(i)*fni(i)
982 fzi(i)=n3(i)*fni(i)
983 impx=fxi(i)*dt12
984 impy=fyi(i)*dt12
985 impz=fzi(i)*dt12
986 fsav1 =fsav1 +impx
987 fsav2 =fsav2 +impy
988 fsav3 =fsav3 +impz
989 fsav8 =fsav8 +abs(impx)
990 fsav9 =fsav9 +abs(impy)
991 fsav10=fsav10+abs(impz)
992 fsav11=fsav11+fni(i)*dt12
993 ENDDO
994#include "lockon.inc"
995 fsav(1)=fsav(1)+fsav1
996 fsav(2)=fsav(2)+fsav2
997 fsav(3)=fsav(3)+fsav3
998
999 fsav(8)=fsav(8)+fsav8
1000 fsav(9)=fsav(9)+fsav9
1001 fsav(10)=fsav(10)+fsav10
1002 fsav(11)=fsav(11)+fsav11
1003#include "lockoff.inc"
1004C
1005 IF(isensint(1)/=0) THEN
1006 DO i=1,jlt
1007 fsavparit(1,1,i+nft) = fxi(i)
1008 fsavparit(1,2,i+nft) = fyi(i)
1009 fsavparit(1,3,i+nft) = fzi(i)
1010 ENDDO
1011 ENDIF
1012c
1013 IF((anim_v(12)+outp_v(12)+h3d_data%N_VECT_PCONT >0.AND.
1014 . ((tt>=tanim .AND. tt<=tanim_stop).OR.tt>=toutp.OR.(tt>=h3d_data%TH3D.AND.tt<=h3d_data%TH3D_STOP) .OR.
1015 . (manim>=4.AND.manim<=15).OR.h3d_data%MH3D/=0))
1016 . .OR.h3d_data%N_VECT_PCONT_MAX>0)THEN
1017#include "lockon.inc"
1018 DO i=1,jlt
1019 fncont(1,ix1g(i)) =fncont(1,ix1g(i)) + fxi(i)*h1(i)
1020 fncont(2,ix1g(i)) =fncont(2,ix1g(i)) + fyi(i)*h1(i)
1021 fncont(3,ix1g(i)) =fncont(3,ix1g(i)) + fzi(i)*h1(i)
1022 fncont(1,ix2g(i)) =fncont(1,ix2g(i)) + fxi(i)*h2(i)
1023 fncont(2,ix2g(i)) =fncont(2,ix2g(i)) + fyi(i)*h2(i)
1024 fncont(3,ix2g(i)) =fncont(3,ix2g(i)) + fzi(i)*h2(i)
1025 fncont(1,ix3g(i)) =fncont(1,ix3g(i)) + fxi(i)*h3(i)
1026 fncont(2,ix3g(i)) =fncont(2,ix3g(i)) + fyi(i)*h3(i)
1027 fncont(3,ix3g(i)) =fncont(3,ix3g(i)) + fzi(i)*h3(i)
1028 fncont(1,ix4g(i)) =fncont(1,ix4g(i)) + fxi(i)*h4(i)
1029 fncont(2,ix4g(i)) =fncont(2,ix4g(i)) + fyi(i)*h4(i)
1030 fncont(3,ix4g(i)) =fncont(3,ix4g(i)) + fzi(i)*h4(i)
1031 jg = nsvg(i)
1032 IF(jg>0) THEN
1033C en SPMD : traitement a refaire apres reception noeud remote si JG < 0
1034 fncont(1,jg)=fncont(1,jg)- fxi(i)
1035 fncont(2,jg)=fncont(2,jg)- fyi(i)
1036 fncont(3,jg)=fncont(3,jg)- fzi(i)
1037 ELSE ! cas noeud remote en SPMD
1038 jg = -jg
1039 fnconti(nin)%P(1,jg)=fnconti(nin)%P(1,jg)-fxi(i)
1040 fnconti(nin)%P(2,jg)=fnconti(nin)%P(2,jg)-fyi(i)
1041 fnconti(nin)%P(3,jg)=fnconti(nin)%P(3,jg)-fzi(i)
1042 ENDIF
1043 ENDDO
1044#include "lockoff.inc"
1045 ENDIF
1046C---------------------------------
1047C SORTIES TH PAR SOUS INTERFACE
1048C---------------------------------
1049 IF(nisub/=0)THEN
1050 DO jsub=1,nisub
1051 DO j=1,15
1052 fsavsub1(j,jsub)=zero
1053 END DO
1054 ENDDO
1055 DO i=1,jlt
1056 nn = nsvg(i)
1057 IF(nn>0)THEN
1058 in=cn_loc(i)
1059 ie=ce_loc(i)
1060 jj =addsubs(in)
1061 kk =addsubm(ie)
1062 DO WHILE(jj<addsubs(in+1))
1063 jsub=lisubs(jj)
1064 DO WHILE(kk<addsubm(ie+1))
1065 ksub=lisubm(kk)
1066 IF(ksub==jsub)THEN
1067 impx=fxi(i)*dt12
1068 impy=fyi(i)*dt12
1069 impz=fzi(i)*dt12
1070C MAIN side :
1071 fsavsub1(1,jsub)=fsavsub1(1,jsub)+impx
1072 fsavsub1(2,jsub)=fsavsub1(2,jsub)+impy
1073 fsavsub1(3,jsub)=fsavsub1(3,jsub)+impz
1074C
1075 fsavsub1(8,jsub) =fsavsub1(8,jsub) +abs(impx)
1076 fsavsub1(9,jsub) =fsavsub1(9,jsub) +abs(impy)
1077 fsavsub1(10,jsub)=fsavsub1(10,jsub)+abs(impz)
1078C
1079 fsavsub1(11,jsub)=fsavsub1(11,jsub)+fni(i)*dt12
1080 kk=kk+1
1081 GO TO 250
1082 ELSE IF(ksub<jsub)THEN
1083 kk=kk+1
1084 ELSE
1085 GO TO 250
1086 END IF
1087 END DO
1088 250 CONTINUE
1089 jj=jj+1
1090 END DO
1091 END IF
1092 END DO
1093
1094 IF(nspmd>1) THEN
1095C loop split because of a PGI bug
1096 DO i=1,jlt
1097 nn = nsvg(i)
1098 IF(nn<0)THEN
1099 nn = -nn
1100 ie=ce_loc(i)
1101 jj =addsubsfi(nin)%P(nn)
1102 kk =addsubm(ie)
1103 DO WHILE(jj<addsubsfi(nin)%P(nn+1))
1104 jsub=lisubsfi(nin)%P(jj)
1105 DO WHILE(kk<addsubm(ie+1))
1106 ksub=lisubm(kk)
1107 IF(ksub==jsub)THEN
1108 impx=fxi(i)*dt12
1109 impy=fyi(i)*dt12
1110 impz=fzi(i)*dt12
1111C MAIN side :
1112 fsavsub1(1,jsub)=fsavsub1(1,jsub)+impx
1113 fsavsub1(2,jsub)=fsavsub1(2,jsub)+impy
1114 fsavsub1(3,jsub)=fsavsub1(3,jsub)+impz
1115C
1116 fsavsub1(8,jsub) =fsavsub1(8,jsub) +abs(impx)
1117 fsavsub1(9,jsub) =fsavsub1(9,jsub) +abs(impy)
1118 fsavsub1(10,jsub)=fsavsub1(10,jsub)+abs(impz)
1119C
1120 fsavsub1(11,jsub)=fsavsub1(11,jsub)+fni(i)*dt12
1121 kk=kk+1
1122 GO TO 150
1123 ELSE IF(ksub<jsub)THEN
1124 kk=kk+1
1125 ELSE
1126 GO TO 150
1127 END IF
1128 END DO
1129 150 CONTINUE
1130 jj=jj+1
1131 END DO
1132 END IF
1133
1134 END DO
1135
1136 END IF
1137 END IF
1138
1139C---------------------------------
1140C NEW FRICTION MODELS
1141C---------------------------------
1142 IF (mfrot==0) THEN
1143C--- Coulomb friction
1144 DO i=1,jlt
1145 xmu(i) = fric
1146 ENDDO
1147 ELSEIF (mfrot==1) THEN
1148C--- Viscous friction
1149 DO i=1,jlt
1150 aa = n1(i)*vx(i) + n2(i)*vy(i) + n3(i)*vz(i)
1151 v2 = (vx(i) - n1(i)*aa)**2
1152 . + (vy(i) - n2(i)*aa)**2
1153 . + (vz(i) - n3(i)*aa)**2
1154 vv = sqrt(max(em30,v2))
1155 ax1 = x3(i) - x1(i)
1156 ay1 = y3(i) - y1(i)
1157 az1 = z3(i) - z1(i)
1158 ax2 = x4(i) - x2(i)
1159 ay2 = y4(i) - y2(i)
1160 az2 = z4(i) - z2(i)
1161 ax = ay1*az2 - az1*ay2
1162 ay = az1*ax2 - ax1*az2
1163 az = ax1*ay2 - ay1*ax2
1164 area = half*sqrt(ax*ax+ay*ay+az*az)
1165 p = -fni(i)/area
1166 xmu(i) = fric + (frot_p(1) + frot_p(4)*p ) * p
1167 . +(frot_p(2) + frot_p(3)*p) * vv + frot_p(5)*v2
1168 ENDDO
1169 ELSEIF(mfrot==2)THEN
1170C--- Loi Darmstad
1171 DO i=1,jlt
1172 aa = n1(i)*vx(i) + n2(i)*vy(i) + n3(i)*vz(i)
1173 v2 = (vx(i) - n1(i)*aa)**2
1174 . + (vy(i) - n2(i)*aa)**2
1175 . + (vz(i) - n3(i)*aa)**2
1176 vv = sqrt(max(em30,v2))
1177 ax1 = x3(i) - x1(i)
1178 ay1 = y3(i) - y1(i)
1179 az1 = z3(i) - z1(i)
1180 ax2 = x4(i) - x2(i)
1181 ay2 = y4(i) - y2(i)
1182 az2 = z4(i) - z2(i)
1183 ax = ay1*az2 - az1*ay2
1184 ay = az1*ax2 - ax1*az2
1185 az = ax1*ay2 - ay1*ax2
1186 area = half*sqrt(ax*ax+ay*ay+az*az)
1187 p = -fni(i)/area
1188 xmu(i) = fric
1189 . + frot_p(1)*exp(frot_p(2)*vv)*p*p
1190 . + frot_p(3)*exp(frot_p(4)*vv)*p
1191 . + frot_p(5)*exp(frot_p(6)*vv)
1192 ENDDO
1193 ELSEIF (mfrot==3) THEN
1194C--- Renard
1195 DO i=1,jlt
1196 aa = n1(i)*vx(i) + n2(i)*vy(i) + n3(i)*vz(i)
1197 v2 = (vx(i) - n1(i)*aa)**2
1198 . + (vy(i) - n2(i)*aa)**2
1199 . + (vz(i) - n3(i)*aa)**2
1200 vv = sqrt(max(em30,v2))
1201 IF(vv>=0.AND.vv<=frot_p(5)) THEN
1202 dmu = frot_p(3)-frot_p(1)
1203 vv1 = vv / frot_p(5)
1204 xmu(i) = frot_p(1)+ dmu*vv1*(two-vv1)
1205 ELSEIF(vv>frot_p(5).AND.vv<frot_p(6)) THEN
1206 dmu = frot_p(4)-frot_p(3)
1207 vv1 = (vv - frot_p(5))/(frot_p(6)-frot_p(5))
1208 xmu(i) = frot_p(3)+ dmu * (three-two*vv1)*vv1**2
1209 ELSE
1210 dmu = frot_p(2)-frot_p(4)
1211 vv2 = (vv - frot_p(6))**2
1212 xmu(i) = frot_p(2) - dmu / (one + dmu*vv2)
1213 ENDIF
1214 ENDDO
1215 ELSEIF(mfrot==4)THEN
1216C--- Exponential decay model
1217 DO i=1,jlt
1218 aa = n1(i)*vx(i) + n2(i)*vy(i) + n3(i)*vz(i)
1219 v2 = (vx(i) - n1(i)*aa)**2
1220 . + (vy(i) - n2(i)*aa)**2
1221 . + (vz(i) - n3(i)*aa)**2
1222 vv = sqrt(max(em30,v2))
1223 xmu(i) = frot_p(1)
1224 . + (fric-frot_p(1))*exp(-frot_p(2)*vv)
1225 xmu(i) = max(xmu(i),em30)
1226 ENDDO
1227 ENDIF
1228C------------------
1229C TANGENT FORCE CALCULATION
1230C------------------
1231 fsav4 = zero
1232 fsav5 = zero
1233 fsav6 = zero
1234
1235 fsav12= zero
1236 fsav13= zero
1237 fsav14= zero
1238 fsav15= zero
1239
1240 IF (ifq>=10) THEN
1241C---------------------------------
1242C INCREMENTAL (STIFFNESS) FORMULATION
1243C---------------------------------
1244 IF (ifq==13) THEN
1245 alpha = max(one,alpha0*dt12)
1246 ELSE
1247 alpha = alpha0
1248 ENDIF
1249 DO i=1,jlt
1250 fx = stif0(i)*vx(i)*dt12
1251 fy = stif0(i)*vy(i)*dt12
1252 fz = stif0(i)*vz(i)*dt12
1253
1254 fx = cand_fx(index(i)) + alpha*fx
1255 fy = cand_fy(index(i)) + alpha*fy
1256 fz = cand_fz(index(i)) + alpha*fz
1257
1258 ftn = fx*n1(i) + fy*n2(i) + fz*n3(i)
1259 fx = fx - ftn*n1(i)
1260 fy = fy - ftn*n2(i)
1261 fz = fz - ftn*n3(i)
1262 ft = fx*fx + fy*fy + fz*fz
1263 ft = max(ft,em30)
1264
1265 fn = fxi(i)**2+fyi(i)**2+fzi(i)**2
1266
1267 beta = min(one,xmu(i)*sqrt(fn/ft))
1268
1269 fxt(i) = fx * beta
1270 fyt(i) = fy * beta
1271 fzt(i) = fz * beta
1272
1273 cand_fx(index(i)) = fxt(i)
1274 cand_fy(index(i)) = fyt(i)
1275 cand_fz(index(i)) = fzt(i)
1276 ifpen(index(i)) = 1
1277
1278C------- total force
1279 fxi(i)=fxi(i) + fxt(i)
1280 fyi(i)=fyi(i) + fyt(i)
1281 fzi(i)=fzi(i) + fzt(i)
1282 econvt = econvt
1283 . + dt1*(vx(i)*fxt(i)+vy(i)*fyt(i)+vz(i)*fzt(i))
1284 ENDDO
1285C---------------------------------
1286C TOTAL (VISCOUS) FORMULATION + FRICTION FILTERING
1287C---------------------------------
1288 ELSEIF (ifq>0) THEN
1289
1290 IF (ifq==3) THEN
1291 alpha = max(one,alpha0*dt12)
1292 ELSE
1293 alpha = alpha0
1294 ENDIF
1295 alphi = one - alpha
1296 DO i=1,jlt
1297 vnx = n1(i)*vn(i)
1298 vny = n2(i)*vn(i)
1299 vnz = n3(i)*vn(i)
1300 vx(i) = vx(i) - vnx
1301 vy(i) = vy(i) - vny
1302 vz(i) = vz(i) - vnz
1303 v2 = vx(i)**2 + vy(i)**2 + vz(i)**2
1304 vis2(i) = viscf * vis2(i)
1305 fm2 = (xmu(i)*fni(i))**2
1306 f2 = vis2(i) * v2
1307 a2 = min(f2,fm2) / max(em30,f2)
1308 aa = sqrt(a2 * vis2(i))
1309 fx = aa * vx(i)
1310 fy = aa * vy(i)
1311 fz = aa * vz(i)
1312
1313 fxt(i) = alpha*fx + alphi*cand_fx(index(i))
1314 fyt(i) = alpha*fy + alphi*cand_fy(index(i))
1315 fzt(i) = alpha*fz + alphi*cand_fz(index(i))
1316 cand_fx(index(i)) = fxt(i)
1317 cand_fy(index(i)) = fyt(i)
1318 cand_fz(index(i)) = fzt(i)
1319 ifpen(index(i)) = 1
1320C------- total force
1321 fxi(i) = fxi(i) + fxt(i)
1322 fyi(i) = fyi(i) + fyt(i)
1323 fzi(i) = fzi(i) + fzt(i)
1324 econvt = econvt
1325 . + dt1*(vx(i)*fxt(i)+vy(i)*fyt(i)+vz(i)*fzt(i))
1326 ENDDO
1327 ELSE
1328C---------------------------------
1329C TOTAL (VISCOUS) FORMULATION / NO FRICTION FILTERING
1330C---------------------------------
1331 DO i=1,jlt
1332 vnx = n1(i)*vn(i)
1333 vny = n2(i)*vn(i)
1334 vnz = n3(i)*vn(i)
1335 vx(i) = vx(i) - vnx
1336 vy(i) = vy(i) - vny
1337 vz(i) = vz(i) - vnz
1338 v2 = vx(i)**2 + vy(i)**2 + vz(i)**2
1339 vis2(i) = viscf * vis2(i)
1340 fm2 = (xmu(i)*fni(i))**2
1341 f2 = vis2(i) * v2
1342 a2 = min(f2,fm2) / max(em30,f2)
1343 aa = sqrt(a2 * vis2(i))
1344 fxt(i) = aa * vx(i)
1345 fyt(i) = aa * vy(i)
1346 fzt(i) = aa * vz(i)
1347C------- total force
1348 fxi(i)=fxi(i) + fxt(i)
1349 fyi(i)=fyi(i) + fyt(i)
1350 fzi(i)=fzi(i) + fzt(i)
1351 econvt = econvt + aa * v2 * dt1
1352 ENDDO
1353 ENDIF
1354C---------------------------------
1355 IF((anim_v(12)+outp_v(12)+h3d_data%N_VECT_PCONT>0.AND.
1356 . ((tt>=tanim .AND. tt<=tanim_stop).OR.tt>=toutp.OR.(tt>=h3d_data%TH3D.AND.tt<=h3d_data%TH3D_STOP).OR.
1357 . (manim>=4.AND.manim<=15).OR.h3d_data%MH3D/=0))
1358 . .OR.h3d_data%N_VECT_PCONT_MAX>0)THEN
1359#include "lockon.inc"
1360 DO i=1,jlt
1361 ftcont(1,ix1g(i)) =ftcont(1,ix1g(i)) + fxt(i)*h1(i)
1362 ftcont(2,ix1g(i)) =ftcont(2,ix1g(i)) + fyt(i)*h1(i)
1363 ftcont(3,ix1g(i)) =ftcont(3,ix1g(i)) + fzt(i)*h1(i)
1364 ftcont(1,ix2g(i)) =ftcont(1,ix2g(i)) + fxt(i)*h2(i)
1365 ftcont(2,ix2g(i)) =ftcont(2,ix2g(i)) + fyt(i)*h2(i)
1366 ftcont(3,ix2g(i)) =ftcont(3,ix2g(i)) + fzt(i)*h2(i)
1367 ftcont(1,ix3g(i)) =ftcont(1,ix3g(i)) + fxt(i)*h3(i)
1368 ftcont(2,ix3g(i)) =ftcont(2,ix3g(i)) + fyt(i)*h3(i)
1369 ftcont(3,ix3g(i)) =ftcont(3,ix3g(i)) + fzt(i)*h3(i)
1370 ftcont(1,ix4g(i)) =ftcont(1,ix4g(i)) + fxt(i)*h4(i)
1371 ftcont(2,ix4g(i)) =ftcont(2,ix4g(i)) + fyt(i)*h4(i)
1372 ftcont(3,ix4g(i)) =ftcont(3,ix4g(i)) + fzt(i)*h4(i)
1373 jg = nsvg(i)
1374 IF(jg>0) THEN
1375C en SPMD : traitement a refaire apres reception noeud remote si JG < 0
1376 ftcont(1,jg)=ftcont(1,jg)- fxt(i)
1377 ftcont(2,jg)=ftcont(2,jg)- fyt(i)
1378 ftcont(3,jg)=ftcont(3,jg)- fzt(i)
1379 ELSE ! cas noeud remote en SPMD
1380 jg = -jg
1381 ftconti(nin)%P(1,jg)=ftconti(nin)%P(1,jg)-fxt(i)
1382 ftconti(nin)%P(2,jg)=ftconti(nin)%P(2,jg)-fyt(i)
1383 ftconti(nin)%P(3,jg)=ftconti(nin)%P(3,jg)-fzt(i)
1384 ENDIF
1385 ENDDO
1386#include "lockoff.inc"
1387 ENDIF
1388
1389C---------------------------------
1390 DO i=1,jlt
1391 impx=fxt(i)*dt12
1392 impy=fyt(i)*dt12
1393 impz=fzt(i)*dt12
1394 fsav4 =fsav4 +impx
1395 fsav5 =fsav5 +impy
1396 fsav6 =fsav6 +impz
1397 impx=fxi(i)*dt12
1398 impy=fyi(i)*dt12
1399 impz=fzi(i)*dt12
1400 fsav12=fsav12+abs(impx)
1401 fsav13=fsav13+abs(impy)
1402 fsav14=fsav14+abs(impz)
1403 fsav15=fsav15+sqrt(impx*impx+impy*impy+impz*impz)
1404 ENDDO
1405#include "lockon.inc"
1406 fsav(4) = fsav(4) + fsav4
1407 fsav(5) = fsav(5) + fsav5
1408 fsav(6) = fsav(6) + fsav6
1409
1410 fsav(12) = fsav(12) + fsav12
1411 fsav(13) = fsav(13) + fsav13
1412 fsav(14) = fsav(14) + fsav14
1413 fsav(15) = fsav(15) + fsav15
1414#include "lockoff.inc"
1415C
1416 IF(isensint(1)/=0) THEN
1417 DO i=1,jlt
1418 fsavparit(1,4,i+nft) = fxt(i)
1419 fsavparit(1,5,i+nft) = fyt(i)
1420 fsavparit(1,6,i+nft) = fzt(i)
1421 ENDDO
1422 ENDIF
1423C
1424C---------------------------------
1425C SORTIES TH PAR SOUS INTERFACE
1426C---------------------------------
1427 IF(nisub/=0)THEN
1428 DO i=1,jlt
1429 nn = nsvg(i)
1430 IF(nn>0)THEN
1431 in=cn_loc(i)
1432 ie=ce_loc(i)
1433 jj =addsubs(in)
1434 kk =addsubm(ie)
1435 DO WHILE(jj<addsubs(in+1))
1436 jsub=lisubs(jj)
1437 DO WHILE(kk<addsubm(ie+1))
1438 ksub=lisubm(kk)
1439 IF(ksub==jsub)THEN
1440 impx=fxt(i)*dt12
1441 impy=fyt(i)*dt12
1442 impz=fzt(i)*dt12
1443C MAIN side :
1444 fsavsub1(4,jsub)=fsavsub1(4,jsub)+impx
1445 fsavsub1(5,jsub)=fsavsub1(5,jsub)+impy
1446 fsavsub1(6,jsub)=fsavsub1(6,jsub)+impz
1447C
1448 impx=fxi(i)*dt12
1449 impy=fyi(i)*dt12
1450 impz=fzi(i)*dt12
1451 fsavsub1(12,jsub)=fsavsub1(12,jsub)+abs(impx)
1452 fsavsub1(13,jsub)=fsavsub1(13,jsub)+abs(impy)
1453 fsavsub1(14,jsub)=fsavsub1(14,jsub)+abs(impz)
1454C
1455 fsavsub1(15,jsub)= fsavsub1(15,jsub)
1456 . +sqrt(impx*impx+impy*impy+impz*impz)
1457 kk=kk+1
1458 GO TO 200
1459 ELSE IF(ksub<jsub)THEN
1460 kk=kk+1
1461 ELSE
1462 GO TO 200
1463 END IF
1464 END DO
1465 200 CONTINUE
1466 jj=jj+1
1467 END DO
1468 END IF
1469 END DO
1470
1471 IF(nspmd>1) THEN
1472
1473 DO i=1,jlt
1474 nn = nsvg(i)
1475 IF(nn<0)THEN
1476
1477 nn = -nn
1478 ie=ce_loc(i)
1479 jj =addsubsfi(nin)%P(nn)
1480 kk =addsubm(ie)
1481 DO WHILE(jj<addsubsfi(nin)%P(nn+1))
1482 jsub=lisubsfi(nin)%P(jj)
1483 DO WHILE(kk<addsubm(ie+1))
1484 ksub=lisubm(kk)
1485 IF(ksub==jsub)THEN
1486 impx=fxt(i)*dt12
1487 impy=fyt(i)*dt12
1488 impz=fzt(i)*dt12
1489C MAIN side :
1490 fsavsub1(4,jsub)=fsavsub1(4,jsub)+impx
1491 fsavsub1(5,jsub)=fsavsub1(5,jsub)+impy
1492 fsavsub1(6,jsub)=fsavsub1(6,jsub)+impz
1493C
1494 impx=fxi(i)*dt12
1495 impy=fyi(i)*dt12
1496 impz=fzi(i)*dt12
1497 fsavsub1(12,jsub)=fsavsub1(12,jsub)+abs(impx)
1498 fsavsub1(13,jsub)=fsavsub1(13,jsub)+abs(impy)
1499 fsavsub1(14,jsub)=fsavsub1(14,jsub)+abs(impz)
1500C
1501 fsavsub1(15,jsub)= fsavsub1(15,jsub)
1502 . +sqrt(impx*impx+impy*impy+impz*impz)
1503 kk=kk+1
1504 GO TO 300
1505 ELSE IF(ksub<jsub)THEN
1506 kk=kk+1
1507 ELSE
1508 GO TO 300
1509 END IF
1510 END DO
1511 300 CONTINUE
1512 jj=jj+1
1513 END DO
1514 END IF
1515
1516 END DO
1517
1518 END IF
1519#include "lockon.inc"
1520 DO jsub=1,nisub
1521 nsub=lisub(jsub)
1522 DO j=1,15
1523 fsavsub(j,nsub)=fsavsub(j,nsub)+fsavsub1(j,jsub)
1524 END DO
1525 END DO
1526#include "lockoff.inc"
1527 END IF
1528C---------------------------------
1529#include "lockon.inc"
1530 econtv = econtv + econvt
1531 econt = econt + econtt
1532#include "lockoff.inc"
1533C---------------------------------
1534 IF(kdtint==1)THEN
1535 IF( (visc/=zero.OR.viscf/=zero)
1536 . .AND.(ivis2==0.OR.ivis2==1))THEN
1537 DO i=1,jlt
1538C C(I)=2.*C(I)
1539 IF(msi(i)==zero)THEN
1540 ks(i) =zero
1541 cs(i) =zero
1542 stv(i)=zero
1543 ELSE
1544 cx = four*c(i)*c(i)
1545 cy = eight*msi(i)*kt(i)
1546 aux = sqrt(cx+cy)+two*c(i)
1547 stv(i)= kt(i)*aux*aux/max(cy,em30)
1548 aux = two*cf(i)*cf(i)/max(msi(i),em20)
1549 IF(aux>stv(i))THEN
1550 ks(i) =zero
1551 cs(i) =cf(i)
1552 stv(i)=aux
1553 ELSE
1554 ks(i)= kt(i)
1555 cs(i) =c(i)
1556 ENDIF
1557 ENDIF
1558C
1559 j1=ix1g(i)
1560 IF(ms(j1)==zero)THEN
1561 k1(i) =zero
1562 c1(i) =zero
1563 st1(i)=zero
1564 ELSE
1565 k1(i)=kt(i)*abs(h1(i))
1566 c1(i)=c(i)*abs(h1(i))
1567 cx =four*c1(i)*c1(i)
1568 cy =eight*ms(j1)*k1(i)
1569 aux = sqrt(cx+cy)+two*c1(i)
1570 st1(i)= k1(i)*aux*aux/max(cy,em30)
1571 cfi = cf(i)*abs(h1(i))
1572 aux = two*cfi*cfi/max(ms(j1),em20)
1573 IF(aux>st1(i))THEN
1574 k1(i) =zero
1575 c1(i) =cfi
1576 st1(i)=aux
1577 ENDIF
1578 ENDIF
1579C
1580 j1=ix2g(i)
1581 IF(ms(j1)==zero)THEN
1582 k2(i) =zero
1583 c2(i) =zero
1584 st2(i)=zero
1585 ELSE
1586 k2(i)=kt(i)*abs(h2(i))
1587 c2(i)=c(i)*abs(h2(i))
1588 cx =four*c2(i)*c2(i)
1589 cy =eight*ms(j1)*k2(i)
1590 aux = sqrt(cx+cy)+two*c2(i)
1591 st2(i)= k2(i)*aux*aux/max(cy,em30)
1592 cfi = cf(i)*abs(h2(i))
1593 aux = two*cfi*cfi/max(ms(j1),em20)
1594 IF(aux>st2(i))THEN
1595 k2(i) =zero
1596 c2(i) =cfi
1597 st2(i)=aux
1598 ENDIF
1599 ENDIF
1600C
1601 j1=ix3g(i)
1602 IF(ms(j1)==zero)THEN
1603 k3(i) =zero
1604 c3(i) =zero
1605 st3(i)=zero
1606 ELSE
1607 k3(i)=kt(i)*abs(h3(i))
1608 c3(i)=c(i)*abs(h3(i))
1609 cx =four*c3(i)*c3(i)
1610 cy =eight*ms(j1)*k3(i)
1611 aux = sqrt(cx+cy)+two*c3(i)
1612 st3(i)= k3(i)*aux*aux/max(cy,em30)
1613 cfi = cf(i)*abs(h3(i))
1614 aux = two*cfi*cfi/max(ms(j1),em20)
1615 IF(aux>st3(i))THEN
1616 k3(i) =zero
1617 c3(i) =cfi
1618 st3(i)=aux
1619 ENDIF
1620 ENDIF
1621C
1622 j1=ix4g(i)
1623 IF(ms(j1)==zero)THEN
1624 k4(i) =zero
1625 c4(i) =zero
1626 st4(i)=zero
1627 ELSE
1628 k4(i)=kt(i)*abs(h4(i))
1629 c4(i)=c(i)*abs(h4(i))
1630 cx =four*c4(i)*c4(i)
1631 cy =eight*ms(j1)*k4(i)
1632 aux = sqrt(cx+cy)+two*c4(i)
1633 st4(i)= k4(i)*aux*aux/max(cy,em30)
1634 cfi = cf(i)*abs(h4(i))
1635 aux = two*cfi*cfi/max(ms(j1),em20)
1636 IF(aux>st4(i))THEN
1637 k4(i) =zero
1638 c4(i) =cfi
1639 st4(i)=aux
1640 ENDIF
1641 ENDIF
1642 ENDDO
1643C
1644 ELSE
1645 DO i=1,jlt
1646 ks(i) =stif(i)
1647 cs(i) =zero
1648 stv(i)=ks(i)
1649 k1(i) =stif(i)*abs(h1(i))
1650 c1(i) =zero
1651 st1(i)=k1(i)
1652 k2(i) =stif(i)*abs(h2(i))
1653 c2(i) =zero
1654 st2(i)=k2(i)
1655 k3(i) =stif(i)*abs(h3(i))
1656 c3(i) =zero
1657 st3(i)=k3(i)
1658 k4(i) =stif(i)*abs(h4(i))
1659 c4(i) =zero
1660 st4(i)=k4(i)
1661 ENDDO
1662 ENDIF
1663 ENDIF
1664
1665 IF(idtmin(10)==1.OR.idtmin(10)==2.OR.
1666 . idtmin(10)==5.OR.idtmin(10)==6)THEN
1667
1668 dtmi0 = ep20
1669 IF(kdtint==0)THEN
1670 DO i=1,jlt
1671 dtmi(i) = ep20
1672 mas2 = two * msi(i)
1673 IF(mas2>zero.AND.stif(i)>zero.AND.
1674 . irb(kini(i))==0.AND.irb2(kini(i))==0)THEN
1675 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/stif(i)))
1676 ENDIF
1677 mas2 = two* ms(ix1g(i))
1678 IF(mas2>zero.AND.h1(i)*stif(i)>zero.AND.
1679 . irb(kinet(ix1g(i)))==0.AND.irb2(kinet(ix1g(i)))==0)THEN
1680 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(h1(i)*stif(i))))
1681 ENDIF
1682 mas2 = two * ms(ix2g(i))
1683 IF(mas2>zero.AND.h2(i)*stif(i)>zero.AND.
1684 . irb(kinet(ix2g(i)))==0.AND.irb2(kinet(ix2g(i)))==0)THEN
1685 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(h2(i)*stif(i))))
1686 ENDIF
1687 mas2 = two* ms(ix3g(i))
1688 IF(mas2>zero.AND.h3(i)*stif(i)>zero.AND.
1689 . irb(kinet(ix3g(i)))==0.AND.irb2(kinet(ix3g(i)))==0)THEN
1690 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(h3(i)*stif(i))))
1691 ENDIF
1692 mas2 = two * ms(ix4g(i))
1693 IF(mas2>zero.AND.h4(i)*stif(i)>zero.AND.
1694 . irb(kinet(ix4g(i)))==0.AND.irb2(kinet(ix4g(i)))==0)THEN
1695 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(h4(i)*stif(i))))
1696 ENDIF
1697 dtmi0 = min(dtmi0,dtmi(i))
1698 ENDDO
1699
1700 ELSE
1701 DO i=1,jlt
1702 dtmi(i) = ep20
1703 mas2 = two * msi(i)
1704 mas2 = two * msi(i)
1705 IF(mas2>zero.AND.stv(i)>zero.AND.
1706 . irb(kini(i))==0.AND.irb2(kini(i))==0)THEN
1707 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/stv(i)))
1708 ENDIF
1709 mas2 = two * ms(ix1g(i))
1710 IF(mas2>zero.AND.st1(i)>zero.AND.
1711 . irb(kinet(ix1g(i)))==0.AND.irb2(kinet(ix1g(i)))==0)THEN
1712 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(st1(i))))
1713 ENDIF
1714 mas2 = two * ms(ix2g(i))
1715 IF(mas2>zero.AND.st2(i)>zero.AND.
1716 . irb(kinet(ix2g(i)))==0.AND.irb2(kinet(ix2g(i)))==0)THEN
1717 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(st2(i))))
1718 ENDIF
1719 mas2 = two * ms(ix3g(i))
1720 IF(mas2>zero.AND.st3(i)>zero.AND.
1721 . irb(kinet(ix3g(i)))==0.AND.irb2(kinet(ix3g(i)))==0)THEN
1722 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(st3(i))))
1723 ENDIF
1724 mas2 = two * ms(ix4g(i))
1725 IF(mas2>zero.AND.st4(i)>zero.AND.
1726 . irb(kinet(ix4g(i)))==0.AND.irb2(kinet(ix4g(i)))==0)THEN
1727 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(st4(i))))
1728 ENDIF
1729 dtmi0 = min(dtmi0,dtmi(i))
1730 ENDDO
1731 ENDIF
1732 IF(dtmi0<=dtmin1(10))THEN
1733 DO i=1,jlt
1734 IF(dtmi(i)<=dtmin1(10))THEN
1735 jg = nsvg(i)
1736 IF(jg>0)THEN
1737 ni = itab(jg)
1738 ELSE
1739 ni = itafi(nin)%P(-jg)
1740 ENDIF
1741 IF(idtmin(10)==1)THEN
1742#include "lockon.inc"
1743 WRITE(iout,'(A,E12.4,A,I10)')
1744 . ' **WARNING MINIMUM TIME STEP ',dtmi(i),
1745 . ' IN INTERFACE ',noint
1746 WRITE(iout,'(A,I10)') ' SECONDARY NODE : ',ni
1747 WRITE(iout,'(A,4I10)')' MAIN NODES : ',
1748 . itab(ix1g(i)),itab(ix2g(i)),itab(ix3g(i)),itab(ix4g(i))
1749#include "lockoff.inc"
1750 tstop = tt
1751 ELSEIF(idtmin(10)==2)THEN
1752#include "lockon.inc"
1753 WRITE(iout,'(A,E12.4,A,I10)')
1754 . ' **WARNING MINIMUM TIME STEP ',dtmi(i),
1755 . ' IN INTERFACE ',noint
1756 WRITE(iout,'(A,I10,A,I10)')' DELETE SECONDARY NODE ',
1757 . ni,' FROM INTERFACE ',noint
1758 WRITE(iout,'(A,4I10)')' MAIN NODES : ',
1759 . itab(ix1g(i)),itab(ix2g(i)),itab(ix3g(i)),itab(ix4g(i))
1760 IF(jg>0) THEN
1761 stfa(nsv(cn_loc(i))) = -abs(stfa(nsv(cn_loc(i))))
1762 ELSE
1763 stifi(nin)%P(-jg) = -abs(stifi(nin)%P(-jg))
1764 ENDIF
1765#include "lockoff.inc"
1766 newfront = -1
1767 ELSEIF(idtmin(10)==5)THEN
1768#include "lockon.inc"
1769 WRITE(iout,'(A,E12.4,A,I10)')
1770 . ' **WARNING MINIMUM TIME STEP ',dtmi(i),
1771 . ' IN INTERFACE ',noint
1772 WRITE(iout,'(A,I10)') ' SECONDARY NODE : ',ni
1773 WRITE(iout,'(A,4I10)')' MAIN NODES : ',
1774 . itab(ix1g(i)),itab(ix2g(i)),itab(ix3g(i)),itab(ix4g(i))
1775#include "lockoff.inc"
1776 mstop = 2
1777 ELSEIF(idtmin(10)==6.AND.ilagm==2)THEN
1778 IF(kinet(jg)+kinet(ix1g(i))+kinet(ix2g(i))
1779 . +kinet(ix3g(i))+kinet(ix4g(i))==0 )THEN
1780 cand_n(index(i)) = -iabs(cand_n(index(i)))
1781 stif(i) = 0.
1782 fxi(i) = 0.
1783 fyi(i) = 0.
1784 fzi(i) = 0.
1785 ENDIF
1786 ENDIF
1787 ENDIF
1788 ENDDO
1789 ENDIF
1790 ENDIF
1791C=======================================================================
1792C FORCES sur noeuds mains
1793C=======================================================================
1794 DO i=1,jlt
1795 fx1(i)=fxi(i)*h1(i)
1796 fy1(i)=fyi(i)*h1(i)
1797 fz1(i)=fzi(i)*h1(i)
1798C
1799 fx2(i)=fxi(i)*h2(i)
1800 fy2(i)=fyi(i)*h2(i)
1801 fz2(i)=fzi(i)*h2(i)
1802C
1803 fx3(i)=fxi(i)*h3(i)
1804 fy3(i)=fyi(i)*h3(i)
1805 fz3(i)=fzi(i)*h3(i)
1806C
1807 fx4(i)=fxi(i)*h4(i)
1808 fy4(i)=fyi(i)*h4(i)
1809 fz4(i)=fzi(i)*h4(i)
1810 ENDDO
1811
1812C=======================================================================
1813C FORCES PARITH ON sur noeud d'ancrage SECONDARY
1814C=======================================================================
1815 CALL foat_to_6_float(1 ,jlt ,fxi, fx6)
1816 CALL foat_to_6_float(1 ,jlt ,fyi, fy6)
1817 CALL foat_to_6_float(1 ,jlt ,fzi, fz6)
1818#include "lockon.inc"
1819c noeuds second
1820 DO i = 1,jlt
1821 ig = nsvg(i)
1822 IF(ig > 0)THEN
1823 il = nsv(cn_loc(i))
1824 DO k = 1,6
1825 daanc6(1,k,il) = daanc6(1,k,il) - fx6(k,i)
1826 daanc6(2,k,il) = daanc6(2,k,il) - fy6(k,i)
1827 daanc6(3,k,il) = daanc6(3,k,il) - fz6(k,i)
1828 ENDDO
1829 ELSE
1830C
1831C SPMD remote SECONDARYs
1832C
1833 il = - ig
1834 DO k = 1,6
1835 daanc6fi(nin)%P(1,k,il) = daanc6fi(nin)%P(1,k,il) - fx6(k,i)
1836 daanc6fi(nin)%P(2,k,il) = daanc6fi(nin)%P(2,k,il) - fy6(k,i)
1837 daanc6fi(nin)%P(3,k,il) = daanc6fi(nin)%P(3,k,il) - fz6(k,i)
1838 ENDDO
1839 ENDIF
1840 ENDDO
1841
1842c noeuds matre
1843
1844 CALL foat_to_6_float(1 ,jlt ,fx1, fx6)
1845 CALL foat_to_6_float(1 ,jlt ,fy1, fy6)
1846 CALL foat_to_6_float(1 ,jlt ,fz1, fz6)
1847 DO i = 1,jlt
1848 il = ix1l(i)
1849 DO k = 1,6
1850 daanc6(1,k,il) = daanc6(1,k,il) + fx6(k,i)
1851 daanc6(2,k,il) = daanc6(2,k,il) + fy6(k,i)
1852 daanc6(3,k,il) = daanc6(3,k,il) + fz6(k,i)
1853 ENDDO
1854 ENDDO
1855
1856 CALL foat_to_6_float(1 ,jlt ,fx2, fx6)
1857 CALL foat_to_6_float(1 ,jlt ,fy2, fy6)
1858 CALL foat_to_6_float(1 ,jlt ,fz2, fz6)
1859 DO i = 1,jlt
1860 il = ix2l(i)
1861 DO k = 1,6
1862 daanc6(1,k,il) = daanc6(1,k,il) + fx6(k,i)
1863 daanc6(2,k,il) = daanc6(2,k,il) + fy6(k,i)
1864 daanc6(3,k,il) = daanc6(3,k,il) + fz6(k,i)
1865 ENDDO
1866 ENDDO
1867
1868 CALL foat_to_6_float(1 ,jlt ,fx3, fx6)
1869 CALL foat_to_6_float(1 ,jlt ,fy3, fy6)
1870 CALL foat_to_6_float(1 ,jlt ,fz3, fz6)
1871 DO i = 1,jlt
1872 il = ix3l(i)
1873 DO k = 1,6
1874 daanc6(1,k,il) = daanc6(1,k,il) + fx6(k,i)
1875 daanc6(2,k,il) = daanc6(2,k,il) + fy6(k,i)
1876 daanc6(3,k,il) = daanc6(3,k,il) + fz6(k,i)
1877 ENDDO
1878 ENDDO
1879
1880 CALL foat_to_6_float(1 ,jlt ,fx4, fx6)
1881 CALL foat_to_6_float(1 ,jlt ,fy4, fy6)
1882 CALL foat_to_6_float(1 ,jlt ,fz4, fz6)
1883 DO i = 1,jlt
1884 il = ix4l(i)
1885 DO k = 1,6
1886 daanc6(1,k,il) = daanc6(1,k,il) + fx6(k,i)
1887 daanc6(2,k,il) = daanc6(2,k,il) + fy6(k,i)
1888 daanc6(3,k,il) = daanc6(3,k,il) + fz6(k,i)
1889 ENDDO
1890 ENDDO
1891#include "lockoff.inc"
1892C=======================================================================
1893C=======================================================================
1894C mise a ZERO des FORCES sur noeuds mains et second
1895C si PENE (sur noeud second) < GAPR (gap reel)
1896C=======================================================================
1897 DO i = 1,jlt
1898 IF(gapv(i) > gapr(i))THEN
1899 ig = nsvg(i)
1900 IF(ig > 0)THEN
1901 il = nsv(cn_loc(i))
1902 xsa = n1(i)*(dxanc(1,il)-h1(i)*dxanc(1,ix1l(i))
1903 . -h2(i)*dxanc(1,ix2l(i))
1904 . -h3(i)*dxanc(1,ix3l(i))
1905 . -h4(i)*dxanc(1,ix4l(i)))
1906 . + n2(i)*(dxanc(2,il)-h1(i)*dxanc(2,ix1l(i))
1907 . -h2(i)*dxanc(2,ix2l(i))
1908 . -h3(i)*dxanc(2,ix3l(i))
1909 . -h4(i)*dxanc(2,ix4l(i)))
1910 . + n3(i)*(dxanc(3,il)-h1(i)*dxanc(3,ix1l(i))
1911 . -h2(i)*dxanc(3,ix2l(i))
1912 . -h3(i)*dxanc(3,ix3l(i))
1913 . -h4(i)*dxanc(3,ix4l(i)))
1914 ELSE
1915C
1916C SPMD remote SECONDARYs
1917C
1918C ******** ATTENTION DXANCFI a communiquer dans TRI20BOX ************
1919C
1920 il = - ig
1921 xsa = n1(i)*(dxancfi(nin)%P(1,il)-h1(i)*dxanc(1,ix1l(i))
1922 . -h2(i)*dxanc(1,ix2l(i))
1923 . -h3(i)*dxanc(1,ix3l(i))
1924 . -h4(i)*dxanc(1,ix4l(i)))
1925 . + n2(i)*(dxancfi(nin)%P(2,il)-h1(i)*dxanc(2,ix1l(i))
1926 . -h2(i)*dxanc(2,ix2l(i))
1927 . -h3(i)*dxanc(2,ix3l(i))
1928 . -h4(i)*dxanc(2,ix4l(i)))
1929 . + n3(i)*(dxancfi(nin)%P(3,il)-h1(i)*dxanc(3,ix1l(i))
1930 . -h2(i)*dxanc(3,ix2l(i))
1931 . -h3(i)*dxanc(3,ix3l(i))
1932 . -h4(i)*dxanc(3,ix4l(i)))
1933 END IF
1934 ps = pene(i) - xsa - gapv(i) + gapr(i)
1935 IF(ps <= zero)THEN
1936 stif(i) = zero
1937 fxi(i) = zero
1938 fyi(i) = zero
1939 fzi(i) = zero
1940 fx1(i) = zero
1941 fy1(i) = zero
1942 fz1(i) = zero
1943 fx2(i) = zero
1944 fy2(i) = zero
1945 fz2(i) = zero
1946 fx3(i) = zero
1947 fy3(i) = zero
1948 fz3(i) = zero
1949 fx4(i) = zero
1950 fy4(i) = zero
1951 fz4(i) = zero
1952 IF (ifq>0) THEN
1953 cand_fx(index(i)) = zero
1954 cand_fy(index(i)) = zero
1955 cand_fz(index(i)) = zero
1956C IFPEN(INDEX(I)) = 0
1957 ENDIF
1958 ENDIF
1959 ENDIF
1960 ENDDO
1961C=======================================================================
1962C FORCES sur noeuds maites et second
1963C=======================================================================
1964C---------------------------------
1965 IF(intth == 0 .OR. iform == 0) THEN
1966 DO i=1,jlt
1967 phi1(i) = zero
1968 phi2(i) = zero
1969 phi3(i) = zero
1970 phi4(i) = zero
1971C
1972 ENDDO
1973 ELSEIF(iform > 0) THEN
1974 DO i=1,jlt
1975 tm = h1(i)*temp(ix1g(i)) + h2(i)*temp(ix2g(i))
1976 . + h3(i)*temp(ix3g(i)) + h4(i)*temp(ix4g(i))
1977
1978 ts = tempi(i)
1979C
1980 ax1 = xa(1,ix3l(i)) - xa(1,ix1l(i))
1981 ay1 = xa(2,ix3l(i)) - xa(2,ix1l(i))
1982 az1 = xa(3,ix3l(i)) - xa(3,ix1l(i))
1983 ax2 = xa(1,ix4l(i)) - xa(1,ix2l(i))
1984 ay2 = xa(2,ix4l(i)) - xa(2,ix2l(i))
1985 az2 = xa(3,ix4l(i)) - xa(3,ix2l(i))
1986C
1987 ax = ay1*az2 - az1*ay2
1988 ay = az1*ax2 - ax1*az2
1989 az = ax1*ay2 - ay1*ax2
1990C
1991 area = one_over_8*sqrt(ax*ax+ay*ay+az*az)
1992 phi(i) = area* (tm - ts)*dt1 / rstif
1993 phi1(i) = -phi(i) *h1(i)
1994 phi2(i) = -phi(i) *h2(i)
1995 phi3(i) = -phi(i) *h3(i)
1996 phi4(i) = -phi(i) *h4(i)
1997 ENDDO
1998 ENDIF
1999C spmd : identification des noeuds interf. utiles a envoyer
2000 IF (nspmd>1) THEN
2001Ctmp+1 mic only
2002#include "mic_lockon.inc"
2003 DO i = 1,jlt
2004 nn = nsvg(i)
2005 IF(nn<0)THEN
2006C tag temporaire de NSVFI a -
2007 nsvfi(nin)%P(-nn) = -abs(nsvfi(nin)%P(-nn))
2008 ENDIF
2009 ENDDO
2010ctmp+1 mic only
2011#include "mic_lockoff.inc"
2012 ENDIF
2013C
2014 IF(idtmins==2.OR.idtmins_int/=0)THEN
2015 dti=dt2t
2016 CALL i7sms2(jlt ,ix1g ,ix2g ,ix3g ,ix4g ,
2017 2 nsvg ,h1 ,h2 ,h3 ,h4 ,stif ,
2018 3 nin ,noint ,mskyi_sms, iskyi_sms,nsms ,
2019 4 kt ,c ,cf ,dtmini,dti )
2020 IF(dti<dt2t)THEN
2021 dt2t = dti
2022 neltst = noint
2023 ityptst = 10
2024 ENDIF
2025 ENDIF
2026C
2027 IF(idtmins_int/=0)THEN
2028 stif(1:jlt)=zero
2029 END IF
2030C
2031 bid = zero
2032 ibid =0
2033 IF(iparit==3)THEN
2034 IF(kdtint==0)THEN
2035 CALL i7ass3(jlt ,ix1g ,ix2g ,ix3g ,ix4g ,
2036 2 nsvg ,h1 ,h2 ,h3 ,h4 ,stif ,
2037 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
2038 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
2039 5 fxi ,fyi ,fzi ,a ,stifn)
2040 ELSE
2041 CALL i7ass35(jlt ,ix1g ,ix2g ,ix3g ,ix4g ,
2042 2 nsvg ,h1 ,h2 ,h3 ,h4 ,stif ,
2043 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
2044 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
2045 5 fxi ,fyi ,fzi ,a ,stifn,viscn,
2046 6 ks ,k1 ,k2 ,k3 ,k4 ,cs ,
2047 7 c1 ,c2 ,c3 ,c4 )
2048 ENDIF
2049 ELSEIF(iparit==0)THEN
2050 IF(kdtint==0)THEN
2051 CALL i7ass0(jlt ,ix1g ,ix2g ,ix3g ,ix4g ,
2052 2 nsvg ,h1 ,h2 ,h3 ,h4 ,stif ,
2053 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
2054 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
2055 5 fxi ,fyi ,fzi ,a ,stifn ,nin ,
2056 6 intth ,phi ,fthe ,phi1 , phi2 ,phi3 ,
2057 7 phi4 ,bid ,bid ,jtask,ibid ,ibid )
2058
2059 ELSE
2060C
2061 CALL i7ass05(jlt ,ix1g ,ix2g ,ix3g ,ix4g ,
2062 2 nsvg ,h1 ,h2 ,h3 ,h4 ,
2063 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
2064 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
2065 5 fxi ,fyi ,fzi ,a ,stifn ,viscn ,
2066 6 ks ,k1 ,k2 ,k3 ,k4 ,cs ,
2067 7 c1 ,c2 ,c3 ,c4 ,nin ,intth ,
2068 8 phi ,fthe ,phi1 , phi2 ,phi3 , phi4,jtask,
2069 9 bid ,bid ,ibid ,ibid )
2070 ENDIF
2071C
2072 ELSE
2073 IF(kdtint==0)THEN
2074 CALL i7ass2(jlt ,ix1g ,ix2g ,ix3g ,ix4g ,
2075 2 nsvg ,h1 ,h2 ,h3 ,h4 ,stif ,
2076 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
2077 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
2078 5 fxi ,fyi ,fzi ,fskyi,isky ,niskyfi,
2079 6 nin ,noint ,intth,phi ,ftheskyi,phi1,
2080 7 phi2 ,phi3 , phi4,bid ,bid ,
2081 a ibid ,ibid )
2082 ELSE
2083 CALL i7ass25(jlt ,ix1g ,ix2g ,ix3g ,ix4g ,
2084 2 nsvg ,h1 ,h2 ,h3 ,h4 ,
2085 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
2086 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
2087 5 fxi ,fyi ,fzi ,fskyi,niskyfi,nin ,
2088 6 ks ,k1 ,k2 ,k3 ,k4 ,cs ,
2089 7 c1 ,c2 ,c3 ,c4 ,isky ,noint ,
2090 8 intth ,phi ,ftheskyi,phi1,phi2 ,phi3,
2091 9 phi4 ,bid ,bid ,ibid ,ibid )
2092 ENDIF
2093 ENDIF
2094C
2095 IF(anim_v(4)+outp_v(4)+h3d_data%N_VECT_CONT>0)THEN
2096#include "lockon.inc"
2097c goto 1234
2098 DO i=1,jlt
2099 fcont(1,ix1g(i)) =fcont(1,ix1g(i)) + fx1(i)
2100 fcont(2,ix1g(i)) =fcont(2,ix1g(i)) + fy1(i)
2101 fcont(3,ix1g(i)) =fcont(3,ix1g(i)) + fz1(i)
2102 fcont(1,ix2g(i)) =fcont(1,ix2g(i)) + fx2(i)
2103 fcont(2,ix2g(i)) =fcont(2,ix2g(i)) + fy2(i)
2104 fcont(3,ix2g(i)) =fcont(3,ix2g(i)) + fz2(i)
2105 fcont(1,ix3g(i)) =fcont(1,ix3g(i)) + fx3(i)
2106 fcont(2,ix3g(i)) =fcont(2,ix3g(i)) + fy3(i)
2107 fcont(3,ix3g(i)) =fcont(3,ix3g(i)) + fz3(i)
2108 fcont(1,ix4g(i)) =fcont(1,ix4g(i)) + fx4(i)
2109 fcont(2,ix4g(i)) =fcont(2,ix4g(i)) + fy4(i)
2110 fcont(3,ix4g(i)) =fcont(3,ix4g(i)) + fz4(i)
2111 jg = nsvg(i)
2112 IF(jg>0) THEN
2113C en SPMD : traitement a refaire apres reception noeud remote si JG < 0
2114 fcont(1,jg)=fcont(1,jg)- fxi(i)
2115 fcont(2,jg)=fcont(2,jg)- fyi(i)
2116 fcont(3,jg)=fcont(3,jg)- fzi(i)
2117 ENDIF
2118 ENDDO
2119c 1234 continue
2120#include "lockoff.inc"
2121 ENDIF
2122C-----------------------------------------------------
2123 IF(isecin>0)THEN
2124 k0=nstrf(25)
2125 IF(nstrf(1)+nstrf(2)/=0)THEN
2126 DO i=1,nsect
2127 nbinter=nstrf(k0+14)
2128 k1s=k0+30
2129 DO j=1,nbinter
2130 IF(nstrf(k1s)==noint)THEN
2131 IF(isecut/=0)THEN
2132#include "lockon.inc"
2133 DO k=1,jlt
2134C attention aux signes pour le cumul des efforts
2135C a rendre conforme avec CFORC3
2136 IF(secfcum(4,ix1g(k),i)==1.)THEN
2137 secfcum(1,ix1g(k),i)=secfcum(1,ix1g(k),i)-fx1(k)
2138 secfcum(2,ix1g(k),i)=secfcum(2,ix1g(k),i)-fy1(k)
2139 secfcum(3,ix1g(k),i)=secfcum(3,ix1g(k),i)-fz1(k)
2140 ENDIF
2141 IF(secfcum(4,ix2g(k),i)==1.)THEN
2142 secfcum(1,ix2g(k),i)=secfcum(1,ix2g(k),i)-fx2(k)
2143 secfcum(2,ix2g(k),i)=secfcum(2,ix2g(k),i)-fy2(k)
2144 secfcum(3,ix2g(k),i)=secfcum(3,ix2g(k),i)-fz2(k)
2145 ENDIF
2146 IF(secfcum(4,ix3g(k),i)==1.)THEN
2147 secfcum(1,ix3g(k),i)=secfcum(1,ix3g(k),i)-fx3(k)
2148 secfcum(2,ix3g(k),i)=secfcum(2,ix3g(k),i)-fy3(k)
2149 secfcum(3,ix3g(k),i)=secfcum(3,ix3g(k),i)-fz3(k)
2150 ENDIF
2151 IF(secfcum(4,ix4g(k),i)==1.)THEN
2152 secfcum(1,ix4g(k),i)=secfcum(1,ix4g(k),i)-fx4(k)
2153 secfcum(2,ix4g(k),i)=secfcum(2,ix4g(k),i)-fy4(k)
2154 secfcum(3,ix4g(k),i)=secfcum(3,ix4g(k),i)-fz4(k)
2155 ENDIF
2156
2157 jg = nsvg(k)
2158 IF(jg>0) THEN
2159C en SPMD : traitement a refaire apres reception noeud remote si JG < 0
2160 IF(secfcum(4,jg,i)==1.)THEN
2161 secfcum(1,jg,i)=secfcum(1,jg,i)+fxi(k)
2162 secfcum(2,jg,i)=secfcum(2,jg,i)+fyi(k)
2163 secfcum(3,jg,i)=secfcum(3,jg,i)+fzi(k)
2164 ENDIF
2165 ENDIF
2166
2167 ENDDO
2168#include "lockoff.inc"
2169 ENDIF
2170C +fsav(section)
2171 ENDIF
2172 k1s=k1s+1
2173 ENDDO
2174 k0=nstrf(k0+24)
2175 ENDDO
2176 ENDIF
2177 ENDIF
2178C-----------------------------------------------------
2179
2180 IF(ibag/=0.OR.iadm/=0)THEN
2181 DO i=1,jlt
2182
2183C IF(PENE(I)/=ZERO)THEN
2184C test modifie pour coherence avec communication SPMD (spmd_i7tools)
2185 IF(fxi(i)/=zero.OR.fyi(i)/=zero.OR.fzi(i)/=zero)THEN
2186
2187 jg = nsvg(i)
2188 IF(jg>0) THEN
2189C en SPMD : traitement a refaire apres reception noeud remote si JG < 0
2190 icontact(jg)=1
2191 ENDIF
2192
2193 icontact(ix1g(i))=1
2194 icontact(ix2g(i))=1
2195 icontact(ix3g(i))=1
2196 icontact(ix4g(i))=1
2197 ENDIF
2198 ENDDO
2199 ENDIF
2200
2201 IF(iadm/=0)THEN
2202 DO i=1,jlt
2203 jg = nsvg(i)
2204#include "lockon.inc"
2205 IF(jg>0) THEN
2206C en SPMD : traitement a refaire apres reception noeud remote si JG < 0
2207 rcontact(jg)=min(rcontact(jg),rcurvi(i))
2208 END IF
2209 rcontact(ix1g(i))=min(rcontact(ix1g(i)),rcurvi(i))
2210 rcontact(ix2g(i))=min(rcontact(ix2g(i)),rcurvi(i))
2211 rcontact(ix3g(i))=min(rcontact(ix3g(i)),rcurvi(i))
2212 rcontact(ix4g(i))=min(rcontact(ix4g(i)),rcurvi(i))
2213#include "lockoff.inc"
2214 END DO
2215 END IF
2216 IF(iadm>=2)THEN
2217 DO i=1,jlt
2218 jg = nsvg(i)
2219#include "lockon.inc"
2220 IF(jg>0) THEN
2221C en SPMD : traitement a refaire apres reception noeud remote si JG < 0
2222 pcontact(jg)=max(pcontact(jg),pene(i)/(padm*gapv(i)))
2223 acontact(jg)=min(acontact(jg),anglmi(i))
2224 END IF
2225#include "lockoff.inc"
2226 END DO
2227 END IF
2228
2229 IF(ibcc==0) RETURN
2230C
2231 DO 400 i=1,jlt
2232
2233 IF(pene(i)==zero)GOTO 400
2234 ibcm = ibcc / 8
2235 ibcs = ibcc - 8 * ibcm
2236 IF(ibcs>0) THEN
2237 ig=nsvg(i)
2238 IF(ig>0) THEN
2239C en SPMD : traitement a refaire apres reception noeud remote si JG < 0
2240 CALL ibcoff(ibcs,icodt(ig))
2241 ENDIF
2242 ENDIF
2243 IF(ibcm>0) THEN
2244 ig=ix1g(i)
2245 CALL ibcoff(ibcm,icodt(ig))
2246 ig=ix2g(i)
2247 CALL ibcoff(ibcm,icodt(ig))
2248 ig=ix3g(i)
2249 CALL ibcoff(ibcm,icodt(ig))
2250 ig=ix4g(i)
2251 CALL ibcoff(ibcm,icodt(ig))
2252 ENDIF
2253 400 CONTINUE
2254C
2255 RETURN
2256 END
2257!||====================================================================
2258!|| i20for3c ../engine/source/interfaces/int20/i20for3.F
2259!||--- called by ------------------------------------------------------
2260!|| i20mainf ../engine/source/interfaces/int20/i20mainf.F
2261!||--- uses -----------------------------------------------------
2262!|| icontact_mod ../engine/share/modules/icontact_mod.F
2263!|| tri7box ../engine/share/modules/tri7box.F
2264!||====================================================================
2265 SUBROUTINE i20for3c(NLN ,NLG ,MS ,DXANC ,
2266 2 DVANC ,STFA ,WEIGHT ,INACTI,
2267 3 DAANC6,STFAC,PENIA ,ALPHAK,
2268 4 DAANC ,KMIN )
2269C-----------------------------------------------
2270C M o d u l e s
2271C-----------------------------------------------
2272 USE tri7box
2273 USE icontact_mod
2274C-----------------------------------------------
2275C I m p l i c i t T y p e s
2276C-----------------------------------------------
2277#include "implicit_f.inc"
2278#include "comlock.inc"
2279C-----------------------------------------------
2280C C o m m o n B l o c k s
2281C-----------------------------------------------
2282#include "com06_c.inc"
2283#include "com08_c.inc"
2284#include "scr11_c.inc"
2285C-----------------------------------------------
2286C D u m m y A r g u m e n t s
2287C-----------------------------------------------
2288 INTEGER WEIGHT(*),NLN,INACTI,NLG(*)
2289 my_real
2290 . DVANC(3,*),DXANC(3,*),DAANC(3,*),STFA(*),PENIA(5,*),
2291 . STFAC,MS(*),ALPHAK(3,*),KMIN
2292 DOUBLE PRECISION
2293 . DAANC6(3,6,*)
2294C-----------------------------------------------
2295C L o c a l V a r i a b l e s
2296C-----------------------------------------------
2297 INTEGER I, J
2298
2299 my_real
2300 . DDX,STFR,VISR,UNSDT2,FX,FY,FZ,ECONTT, ECONVT,
2301 . FXR(NLN), FYR(NLN), FZR(NLN),DX(NLN),DY(NLN),DZ(NLN)
2302 DOUBLE PRECISION
2303 . FX6(6,NLN), FY6(6,NLN), FZ6(6,NLN)
2304
2305C-----------------------------------------------
2306
2307 unsdt2 = dt2/max(dt2*dt2,em30)
2308
2309C----------------------------------------------------------------------
2310C penetration initiale
2311C
2312C PENIA(1:3,I) vecteur directeur normee
2313C PENIA(4,I) module de la penetration initiale
2314C PENIA(4,I) module corrigee pour le cycle suivant
2315C----------------------------------------------------------------------
2316 IF(inacti >= 5)THEN
2317 DO i = 1,nln
2318 dx(i) = dxanc(1,i) - penia(1,i)*penia(4,i)
2319 dy(i) = dxanc(2,i) - penia(2,i)*penia(4,i)
2320 dz(i) = dxanc(3,i) - penia(3,i)*penia(4,i)
2321 ddx = dx(i)*penia(1,i) + dy(i)*penia(2,i) + dz(i)*penia(3,i)
2322 ddx = half*min(ddx,zero)
2323 penia(5,i) = max(zero,penia(5,i),penia(4,i)+ddx)
2324c a faire ici ou dans I20BUCE_CRIT IF(PENIA(5,I) /= ZERO)NACTI=NACTI+1
2325 ENDDO
2326 ELSE
2327 DO i = 1,nln
2328 dx(i) = dxanc(1,i)
2329 dy(i) = dxanc(2,i)
2330 dz(i) = dxanc(3,i)
2331 ENDDO
2332 ENDIF
2333C----------------------------------------------------------------------
2334C NOEUDS MAIN SECONDARY edge
2335C----------------------------------------------------------------------
2336
2337 econtt = zero
2338 econvt = zero
2339 econtv = zero
2340 DO i = 1,nln
2341 j = nlg(i)
2342 IF(stfac > zero)THEN
2343c STFR = HALF * STFAC * ABS(STFA(I))
2344 stfr = half * max(kmin,stfac*abs(stfa(i))) * alphak(1,i)
2345 ELSE
2346c STFR = HALF * ABS(STFAC)
2347 stfr = half * abs(stfac) * alphak(1,i)
2348 ENDIF
2349c essai viscosite critique
2350c VISR = 0.1 * TWO * SQRT(STFR * MS(J))
2351 visr = two * sqrt(stfr * ms(j))
2352
2353 fx = stfr * dx(i)
2354 fy = stfr * dy(i)
2355 fz = stfr * dz(i)
2356 fxr(i) = fx + visr * dvanc(1,i)
2357 fyr(i) = fy + visr * dvanc(2,i)
2358 fzr(i) = fz + visr * dvanc(3,i)
2359 IF(fx /= zero)alphak(3,i)=min(alphak(3,i),fxr(i)/fx)
2360 IF(fy /= zero)alphak(3,i)=min(alphak(3,i),fyr(i)/fy)
2361 IF(fz /= zero)alphak(3,i)=min(alphak(3,i),fzr(i)/fz)
2362 daanc(1,i) = - fxr(i)
2363 daanc(2,i) = - fyr(i)
2364 daanc(3,i) = - fzr(i)
2365
2366 IF(weight(j) == 1)THEN
2367 econtt = econtt + half*stfr*(dx(i)**2+dy(i)**2+dz(i)**2)
2368 econvt = econvt
2369 . + visr*(dvanc(1,i)**2+dvanc(2,i)**2+dvanc(3,i)**2)
2370 ENDIF
2371 ENDDO
2372
2373#include "lockon.inc"
2374 econtv = econtv + econvt*dt1
2375 econt = econt + econtt
2376#include "lockoff.inc"
2377
2378 RETURN
2379 END
2380!||====================================================================
2381!|| i20for3e ../engine/source/interfaces/int20/i20for3.F
2382!||--- called by ------------------------------------------------------
2383!|| i20mainf ../engine/source/interfaces/int20/i20mainf.F
2384!||--- calls -----------------------------------------------------
2385!|| foat_to_6_float ../engine/source/system/parit.F
2386!|| i20ass0 ../engine/source/interfaces/int20/i20for3.F
2387!|| i20ass05 ../engine/source/interfaces/int20/i20for3.F
2388!|| i20ass2 ../engine/source/interfaces/int20/i20for3.F
2389!|| i20ass25 ../engine/source/interfaces/int20/i20for3.F
2390!|| i20sms2e ../engine/source/interfaces/int20/i20sms2.F
2391!||--- uses -----------------------------------------------------
2392!|| h3d_mod ../engine/share/modules/h3d_mod.F
2393!|| tri7box ../engine/share/modules/tri7box.F
2394!||====================================================================
2395 SUBROUTINE i20for3e(
2396 1 JLT ,A ,V ,IBC ,ICODT ,
2397 2 FSAV ,GAP ,FRIC ,MS ,VISC ,
2398 3 VISCF ,NOINT ,ITAB ,CS_LOC ,CM_LOC ,
2399 4 STIGLO ,STIFN ,STIF ,FSKYI ,ISKY ,
2400 5 FCONT ,STFS ,STFM ,DT2T ,HS1 ,
2401 6 HS2 ,HM1 ,HM2 ,N1 ,N2 ,
2402 7 M1 ,M2 ,IVIS2 ,NELTST ,ITYPTST,
2403 8 NX ,NY ,NZ ,GAPV ,PENISE ,
2404 9 PENIME ,INACTI,NISKYFIE,NEWFRONT,ISECIN ,
2405 A NSTRF ,SECFCUM,VISCN ,NLINSA ,MS1 ,
2406 B MS2 ,MM1 ,MM2 ,VXS1 ,VYS1 ,
2407 C VZS1 ,VXS2 ,VYS2 ,VZS2 ,VXM1 ,
2408 D VYM1 ,VZM1 ,VXM2 ,VYM2 ,VZM2 ,
2409 E NIN ,N1L ,N2L ,M1L ,M2L ,
2410 F DAANC6 ,ALPHAK ,MSKYI_SMS,ISKYI_SMS,NSMS,
2411 G JTASK ,ISENSINT, FSAVPARIT ,NISUB ,NFT,
2412 H H3D_DATA)
2413C-----------------------------------------------
2414C M o d u l e s
2415C-----------------------------------------------
2416 USE tri7box
2417 USE h3d_mod
2418C-----------------------------------------------
2419C I m p l i c i t T y p e s
2420C-----------------------------------------------
2421#include "implicit_f.inc"
2422#include "comlock.inc"
2423C-----------------------------------------------
2424C G l o b a l P a r a m e t e r s
2425C-----------------------------------------------
2426#include "mvsiz_p.inc"
2427C-----------------------------------------------
2428C C o m m o n B l o c k s
2429C-----------------------------------------------
2430#include "com01_c.inc"
2431#include "com04_c.inc"
2432#include "com06_c.inc"
2433#include "com08_c.inc"
2434#include "scr05_c.inc"
2435#include "scr07_c.inc"
2436#include "scr11_c.inc"
2437#include "scr14_c.inc"
2438#include "scr16_c.inc"
2439#include "scr18_c.inc"
2440#include "units_c.inc"
2441#include "parit_c.inc"
2442#include "impl1_c.inc"
2443#include "sms_c.inc"
2444C-----------------------------------------------
2445C D u m m y A r g u m e n t s
2446C-----------------------------------------------
2447 INTEGER NELTST,ITYPTST,JLT,IBC,IVIS2,INACTI,NLINSA,NISKYFIE,NIN
2448 INTEGER ICODT(*), ITAB(*), ISKY(*),
2449 . NOINT,NEWFRONT,ISECIN, NSTRF(*), ISKYI_SMS(*)
2450 INTEGER N1(MVSIZ), N2(MVSIZ), M1(MVSIZ), M2(MVSIZ),
2451 . N1L(MVSIZ),N2L(MVSIZ),M1L(MVSIZ),M2L(MVSIZ),
2452 . CS_LOC(MVSIZ), CM_LOC(MVSIZ), NSMS(MVSIZ),JTASK,
2453 . ISENSINT(*),NISUB,NFT
2454 my_real
2455 . STIGLO,
2456 . A(3,*), MS(*), V(3,*), FSAV(*),FCONT(3,*),
2457 . STFS(*),STFM(*),STIFN(*),FSKYI(LSKYI,NFSKYI),GAPV(*),
2458 . PENISE(2,*), PENIME(2,*),ALPHAK(3,*), MSKYI_SMS(*),
2459 . GAP, FRIC,VISC,VISCF,VIS,DT2T
2460 my_real
2461 . HS1(MVSIZ), HS2(MVSIZ), HM1(MVSIZ), HM2(MVSIZ),
2462 . NX(MVSIZ), NY(MVSIZ), NZ(MVSIZ), STIF(MVSIZ),
2463 . SECFCUM(7,NUMNOD,NSECT), VISCN(*),
2464 . MS1(MVSIZ),MS2(MVSIZ),MM1(MVSIZ),MM2(MVSIZ),
2465 . vxs1(mvsiz),vys1(mvsiz),vzs1(mvsiz),vxs2(mvsiz),vys2(mvsiz),
2466 . vzs2(mvsiz),vxm1(mvsiz),vym1(mvsiz),vzm1(mvsiz),vxm2(mvsiz),
2467 . vym2(mvsiz),vzm2(mvsiz),fsavparit(nisub+1,11,*)
2468 DOUBLE PRECISION DAANC6(3,6,*)
2469 TYPE(H3D_DATABASE) :: H3D_DATA
2470C-----------------------------------------------
2471C L o c a l V a r i a b l e s
2472C-----------------------------------------------
2473 INTEGER I, J1, J , K0,NBINTER,K1S,K, NI, IL, IG
2474 INTEGER NISKYL,NISKYL1,ISIGN
2475 my_real
2476 . VX(MVSIZ), VY(MVSIZ), VZ(MVSIZ), VN(MVSIZ),
2477 . FXI(MVSIZ), FYI(MVSIZ), FZI(MVSIZ), FNI(MVSIZ),
2478 . FX1(MVSIZ), FX2(MVSIZ), FX3(MVSIZ), FX4(MVSIZ),
2479 . FY1(MVSIZ), FY2(MVSIZ), FY3(MVSIZ), FY4(MVSIZ),
2480 . FZ1(MVSIZ), FZ2(MVSIZ), FZ3(MVSIZ), FZ4(MVSIZ),
2481 . PENE(MVSIZ),MASMIN(MVSIZ),
2482 . VIS2(MVSIZ), DTMI(MVSIZ),
2483 . VNX, VNY, VNZ, AA, VMAX,S2,DIST,RDIST,
2484 . V2, FM2, DT1INV, VISCA, FAC, FF,
2485 . FX, FY, FZ, F2, MAS2, DTMI0,DTI,
2486 . FACM1, ECONTT, ECONVT, A2,MASM,
2487 . FSAV1, FSAV2, FSAV3, FSAV4, FSAV5, FSAV6,
2488 . DTI2, PPLUS
2489 my_real PREC
2490 my_real
2491 . ST1(MVSIZ),ST2(MVSIZ),ST3(MVSIZ),ST4(MVSIZ),
2492 . KT(MVSIZ),C(MVSIZ),CF(MVSIZ),
2493 . K1(MVSIZ),K2(MVSIZ),K3(MVSIZ),K4(MVSIZ),
2494 . C1(MVSIZ),C2(MVSIZ),C3(MVSIZ),C4(MVSIZ),
2495 . CX,CY,CFI,AUX,AAA
2496 double precision
2497 . fx6(6,mvsiz), fy6(6,mvsiz), fz6(6,mvsiz)
2498C-----------------------------------------------
2499 IF (iresp == 1) THEN
2500 prec = fiveem4
2501 ELSE
2502 prec = em10
2503 ENDIF
2504 IF(dt1>zero)THEN
2505 dt1inv = one/dt1
2506 ELSE
2507 dt1inv =zero
2508 ENDIF
2509 econtt = zero
2510 econvt = zero
2511C
2512 DO i=1,jlt
2513 s2 = sqrt(nx(i)**2 + ny(i)**2 + nz(i)**2)
2514 pene(i) = gapv(i) - s2
2515 s2 = one/max(em30,s2)
2516 nx(i) = nx(i)*s2
2517 ny(i) = ny(i)*s2
2518 nz(i) = nz(i)*s2
2519 ENDDO
2520C
2521 IF(inacti==5.or.inacti==6)THEN
2522#include "lockon.inc"
2523 DO i=1,jlt
2524 pplus=half*(pene(i)+fiveem2*(gapv(i)-pene(i)))
2525 IF(cs_loc(i)<=nlinsa) THEN
2526 penise(2,cs_loc(i)) = max(penise(2,cs_loc(i)),pplus)
2527 ELSE
2528 ni = cs_loc(i)-nlinsa
2529 penfie(nin)%P(2,ni) = max(penfie(nin)%P(2,ni),pplus)
2530 END IF
2531 penime(2,cm_loc(i)) = max(penime(2,cm_loc(i)),pplus)
2532 ENDDO
2533#include "lockoff.inc"
2534 DO i=1,jlt
2535 IF(cs_loc(i)<=nlinsa) THEN
2536 pene(i) = pene(i) - penise(1,cs_loc(i)) - penime(1,cm_loc(i))
2537 pene(i) = max(pene(i),zero)
2538 IF(pene(i)==zero)stif(i)=zero
2539 gapv(i) = gapv(i) - penise(1,cs_loc(i)) - penime(1,cm_loc(i))
2540 ELSE
2541 ni = cs_loc(i)-nlinsa
2542 pene(i) = pene(i) - penfie(nin)%P(1,ni) - penime(1,cm_loc(i))
2543 pene(i) = max(pene(i),zero)
2544 IF(pene(i)==zero)stif(i)=zero
2545 gapv(i) = gapv(i) - penfie(nin)%P(1,ni) - penime(1,cm_loc(i))
2546 END IF
2547 END DO
2548 ENDIF
2549
2550 vmax = zero
2551 DO i=1,jlt
2552 gapv(i) = zep9*gapv(i)
2553 vx(i) = hs1(i)*vxs1(i) + hs2(i)*vxs2(i)
2554 . - hm1(i)*vxm1(i) - hm2(i)*vxm2(i)
2555 vy(i) = hs1(i)*vys1(i) + hs2(i)*vys2(i)
2556 . - hm1(i)*vym1(i) - hm2(i)*vym2(i)
2557 vz(i) = hs1(i)*vzs1(i) + hs2(i)*vzs2(i)
2558 . - hm1(i)*vzm1(i) - hm2(i)*vzm2(i)
2559 vn(i) = nx(i)*vx(i) + ny(i)*vy(i) + nz(i)*vz(i)
2560 ENDDO
2561C-------------------------------------------
2562 DO i=1,jlt
2563 fac = gapv(i)/max( em10,( gapv(i)-pene(i) ) )
2564 facm1 = one/fac
2565 IF(( (gapv(i)-pene(i))/gapv(i) )<prec .AND.
2566 . stif(i)>zero ) THEN
2567 stif(i) = zero
2568 IF (impl_s==0) THEN
2569 newfront = -1
2570#include "lockon.inc"
2571 IF(cs_loc(i)<=nlinsa)THEN
2572 stfs(cs_loc(i)) = -abs(stfs(cs_loc(i)))
2573 WRITE(istdo,*)'WARNING INTERFACE NB',noint
2574 WRITE(istdo,*)'LINE CONNECTING NODES ',itab(n1(i)),
2575 . itab(n2(i)),'DE-ACTIVATED FROM INTERFACE'
2576 WRITE(istdo,*)'IMPACTED ON ',itab(m1(i)),itab(m2(i))
2577 WRITE(iout,*)'WARNING INTERFACE NB',noint
2578 WRITE(iout,*)'GAP=',gapv(i),'PENE=',pene(i)
2579 WRITE(iout,*)'LINE CONNECTING NODES ',itab(n1(i)),
2580 . itab(n2(i)),'DE-ACTIVATED FROM INTERFACE'
2581 WRITE(iout,*)'IMPACTED ON ',itab(m1(i)),itab(m2(i))
2582 ELSE
2583 ni = cs_loc(i)-nlinsa
2584 stifie(nin)%P(ni) = -abs(stifie(nin)%P(ni))
2585 WRITE(istdo,*)'WARNING INTERFACE NB',noint
2586 WRITE(istdo,*)'LINE CONNECTING NODES ',itafie(nin)%P(n1(i)),
2587 . itafie(nin)%P(n2(i)),'DE-ACTIVATED FROM INTERFACE'
2588 WRITE(iout,*)'WARNING INTERFACE NB',noint
2589 WRITE(iout,*)'GAP=',gapv(i),'PENE=',pene(i)
2590 WRITE(iout,*)'LINE CONNECTING NODES ',itafie(nin)%P(n1(i)),
2591 . itafie(nin)%P(n2(i)),'DE-ACTIVATED FROM INTERFACE'
2592 END IF
2593#include "lockoff.inc"
2594 ENDIF
2595 pene(i)= zero
2596 ENDIF
2597 econtt = econtt + half*stif(i)*gapv(i)**2 *( facm1 - one -
2598 . log(facm1) )
2599 stif(i) = half*stif(i) * fac
2600 fni(i)= -stif(i) * pene(i)
2601 ENDDO
2602
2603 dti = ep20
2604C
2605 DO i=1,jlt
2606 dist=gapv(i)-pene(i)
2607 rdist = half*dist / max(em30,-vn(i))
2608 dti = min(rdist,dti)
2609 ENDDO
2610C
2611 IF(dti<=dtmin1(10))THEN
2612 DO i=1,jlt
2613 dist=gapv(i)-pene(i)
2614 dti2 = half*dist / max(em30,-vn(i))
2615 IF(dti2<=dtmin1(10))THEN
2616#include "lockon.inc"
2617 IF(cs_loc(i)<=nlinsa)THEN
2618 WRITE(iout,*)
2619 . ' **WARNING MINIMUM TIME STEP ',dti2,
2620 . 'IN INTERFACE NB',noint
2621 WRITE(iout,*)'SECONDARY NODES NB',itab(n1(i)),
2622 . itab(n2(i))
2623 WRITE(iout,*)'MAIN NODES NB',itab(m1(i)),
2624 . itab(m2(i))
2625 ELSE
2626 WRITE(iout,*)
2627 . ' **WARNING MINIMUM TIME STEP ',dti2,
2628 . 'IN INTERFACE NB',noint
2629 WRITE(iout,*)'SECONDARY NODES NB',itafie(nin)%P(n1(i)),
2630 . itafie(nin)%P(n2(i))
2631 WRITE(iout,*)'MAIN NODES NB',itab(m1(i)),
2632 . itab(m2(i))
2633 END IF
2634#include "lockoff.inc"
2635 IF(idtmin(10)==1)THEN
2636 tstop = tt
2637 ELSEIF(idtmin(10)==2)THEN
2638#include "lockon.inc"
2639 WRITE(iout,*)'REMOVE SECONDARY LINE FROM INTERFACE'
2640 IF(cs_loc(i)<=nlinsa)THEN
2641 stfs(cs_loc(i)) = -abs(stfs(cs_loc(i)))
2642 ELSE
2643 ni = cs_loc(i)-nlinsa
2644 stifie(nin)%P(ni) = -abs(stifie(nin)%P(ni))
2645 END IF
2646#include "lockoff.inc"
2647 newfront = -1
2648 stif(i) = zero
2649 dti = dtmin1(10)
2650 ELSEIF(idtmin(10)==5)THEN
2651 mstop = 2
2652 ENDIF
2653 ENDIF
2654 ENDDO
2655 ENDIF
2656C
2657 IF(dti<dt2t)THEN
2658 dt2t = dti
2659 neltst = noint
2660 ityptst = 10
2661 ENDIF
2662C---------------------------------
2663C DAMPING + FRIC
2664C---------------------------------
2665 IF(visc/=zero.OR.viscf/=zero)THEN
2666 DO i=1,jlt
2667 mas2 = ms1(i)*hs1(i)
2668 . + ms2(i)*hs2(i)
2669 masm = mm1(i)*hm1(i)
2670 . + mm2(i)*hm2(i)
2671 masmin(i) = min(mas2,masm)
2672 vis2(i) = two * stif(i) * min(mas2,masm)
2673 ENDDO
2674 ENDIF
2675C---------------------------------
2676 IF(visc/=zero)THEN
2677 IF(ivis2==0.OR.ivis2==1)THEN
2678C---------------------------------
2679C VISC QUAD TYPE V227
2680C---------------------------------
2681 DO i=1,jlt
2682 IF(vn(i)<zero)
2683 . vis2(i) = vis2(i)/(max(em10,(gapv(i)-pene(i))/gapv(i)))
2684 ENDDO
2685C---------------------------------
2686 visca = zep4
2687 IF(kdtint==0.AND.idtmins/=2)THEN
2688 DO i=1,jlt
2689 fac = stif(i) / max(em30,stif(i))
2690 vis = sqrt(vis2(i))
2691 ff = fac * (
2692 . visc * vis +
2693 . visca**2 * two * masmin(i) * max(zero,-vn(i)) /
2694 . max((gapv(i) - pene(i)),em10) )
2695 stif(i) = stif(i) * gapv(i)/max((gapv(i)-pene(i)),em10)
2696 stif(i) = stif(i) + ff * dt1inv
2697 stif(i) = max(stif(i) ,fac*sqrt(viscf)*vis*dt1inv)
2698 ff = min(ff * vn(i),-fni(i))
2699c FF = MIN(FF * VN(I),ZERO)
2700 fni(i) = fni(i) + ff
2701cc ECONVT = ECONVT + FF * VN(I) * DT1
2702 ENDDO
2703
2704 ELSE
2705 DO i=1,jlt
2706 fac = stif(i) / max(em30,stif(i))
2707 vis = sqrt(vis2(i))
2708 c(i)= fac * (
2709 . visc * vis +
2710 . visca**2 * two * masmin(i) * max(zero,-vn(i)) /
2711 . max((gapv(i) - pene(i)),em10) )
2712 stif(i) = stif(i) * gapv(i) / max((gapv(i) - pene(i)),em10)
2713 kt(i) = stif(i)
2714 stif(i) = stif(i) + c(i) * dt1inv
2715 ff = min(c(i) * vn(i),-fni(i))
2716c FF = MIN(FF * VN(I),ZERO)
2717 fni(i) = fni(i) + ff
2718 cf(i) = fac*sqrt(viscf)*vis
2719 stif(i) = max(stif(i) ,cf(i)*dt1inv)
2720cc ECONVT = ECONVT + C(I) * VN(I) * DT1
2721 ENDDO
2722 ENDIF
2723
2724 ELSEIF(ivis2==2)THEN
2725C---------------------------------
2726C VISC QUAD TYPE
2727C---------------------------------
2728 DO i=1,jlt
2729 vis2(i) = vis2(i)/( max(em10,(gapv(i)-pene(i))/gapv(i)))
2730 ENDDO
2731C---------------------------------
2732 visca = half
2733 DO i=1,jlt
2734 fac = stif(i) / max(em30,stif(i))
2735 vis = sqrt(vis2(i))
2736 ff = fac * (
2737 . visc * vis +
2738 . visca**2 * two * masmin(i) * abs(vn(i)) /
2739 . max((gapv(i) - pene(i)),em10) )
2740 stif(i) = stif(i) * gapv(i) / max((gapv(i)-pene(i)),em10)
2741 stif(i) = stif(i) + two * ff * dt1inv
2742 stif(i) = max(stif(i) ,fac*sqrt(viscf)*vis*dt1inv)
2743 ff = min(ff * vn(i),-fni(i))
2744 fni(i) = fni(i) + ff
2745 ENDDO
2746 ELSEIF(ivis2==3)THEN
2747C---------------------------------
2748C VISC QUAD = 0
2749C---------------------------------
2750 DO i=1,jlt
2751 fac = stif(i) / max(em30,stif(i))
2752 vis = sqrt(vis2(i))
2753 ff = fac * ( visc * vis ) /
2754 . max((gapv(i) - pene(i)),em10)
2755 stif(i) = stif(i) * gapv(i) / max((gapv(i)-pene(i)),em10)
2756 stif(i) = stif(i) + two * ff * dt1inv
2757 stif(i) = max(stif(i) ,fac*sqrt(viscf)*vis*dt1inv)
2758 ff = min(ff * vn(i),-fni(i))
2759 fni(i) = fni(i) + ff
2760 ENDDO
2761 ELSEIF(ivis2==4)THEN
2762C---------------------------------
2763C VISC = 0
2764C---------------------------------
2765 DO i=1,jlt
2766 vis = sqrt(vis2(i))
2767 stif(i) = stif(i) * gapv(i) / max((gapv(i)-pene(i)),em10)
2768 stif(i) = max(stif(i) ,fac*sqrt(viscf)*vis*dt1inv)
2769 ENDDO
2770 ELSEIF(ivis2==5)THEN
2771C---------------------------------
2772C VISC = 2M/dt => pour visc < 1 , stable : dt < 2M/visc = dt
2773C M=M1*M2/M1+M2 pour visc = 1 , choc elastique
2774C pour visc = 0.5 , choc mou
2775C---------------------------------
2776 DO i=1,jlt
2777 mas2 = ms1(i)*hs1(i)
2778 . + ms2(i)*hs2(i)
2779 masm = mm1(i)*hm1(i)
2780 . + mm2(i)*hm2(i)
2781 vis = 2. * visc * dt1inv * masm * mas2 /
2782 . max(em30,masm+mas2)
2783 stif(i) = stif(i) * gapv(i) / max((gapv(i) -pene(i)),em10)
2784 stif(i) = max(stif(i) ,fac*sqrt(viscf*vis2(i))*dt1inv)
2785 ff = vis * vn(i)
2786 econvt = econvt + min(zero,ff-fni(i)) * vn(i) * dt1
2787 fni(i) = min(fni(i),ff)
2788 ENDDO
2789 ELSE
2790 ENDIF
2791 ELSE
2792 DO i=1,jlt
2793 stif(i) = stif(i) * gapv(i) / max((gapv(i) - pene(i)),em10)
2794 ENDDO
2795 ENDIF
2796C---------------------------------
2797C REDUCTION RIGIDITE ANCRAGE
2798C---------------------------------
2799#include "lockon.inc"
2800 DO i=1,jlt
2801 isign=1
2802 IF(pene(i)>zero)isign=-1
2803 aaa = one-pene(i)/gapv(i)
2804 il = m1l(i)
2805 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
2806 alphak(2,il)=isign*min(abs(alphak(2,il)),aaa)
2807 il = m2l(i)
2808 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
2809 alphak(2,il)=isign*min(abs(alphak(2,il)),aaa)
2810 IF(cs_loc(i) <= nlinsa)THEN
2811 il = n1l(i)
2812 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
2813 alphak(2,il)=isign*min(abs(alphak(2,il)),aaa)
2814 il = n2l(i)
2815 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
2816 alphak(2,il)=isign*min(abs(alphak(2,il)),aaa)
2817 ELSE
2818C SPMD remote SECONDARYs
2819 il = n1(i)
2820 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
2821 alphakfie(nin)%P(il)=isign*min(abs(alphakfie(nin)%P(il)),aaa)
2822 il = n2(i)
2823 IF(pene(i)>zero.OR.alphak(2,il)<zero)isign=-1
2824 alphakfie(nin)%P(il)=isign*min(abs(alphakfie(nin)%P(il)),aaa)
2825 ENDIF
2826 ENDDO
2827#include "lockoff.inc"
2828C---------------------------------
2829C SAUVEGARDE DE L'IMPULSION NORMALE
2830C---------------------------------
2831 fsav1 = zero
2832 fsav2 = zero
2833 fsav3 = zero
2834 DO i=1,jlt
2835 fxi(i)=nx(i)*fni(i)
2836 fyi(i)=ny(i)*fni(i)
2837 fzi(i)=nz(i)*fni(i)
2838 fsav1=fsav1+fxi(i)*dt12
2839 fsav2=fsav2+fyi(i)*dt12
2840 fsav3=fsav3+fzi(i)*dt12
2841 ENDDO
2842 IF (imconv==1) THEN
2843#include "lockon.inc"
2844 fsav(1)=fsav(1)+fsav1
2845 fsav(2)=fsav(2)+fsav2
2846 fsav(3)=fsav(3)+fsav3
2847#include "lockoff.inc"
2848 ENDIF
2849 IF(isensint(1)/=0) THEN
2850 DO i=1,jlt
2851 fsavparit(1,1,i+nft) = fxi(i)
2852 fsavparit(1,2,i+nft) = fyi(i)
2853 fsavparit(1,3,i+nft) = fzi(i)
2854 ENDDO
2855 ENDIF
2856C---------------------------------
2857C FRICTION
2858C---------------------------------
2859 IF(fric*viscf/=0.)THEN
2860 fsav4 = zero
2861 fsav5 = zero
2862 fsav6 = zero
2863 DO i=1,jlt
2864 vnx = nx(i)*vn(i)
2865 vny = ny(i)*vn(i)
2866 vnz = nz(i)*vn(i)
2867 vx(i) = vx(i) - vnx
2868 vy(i) = vy(i) - vny
2869 vz(i) = vz(i) - vnz
2870 v2 = vx(i)**2 + vy(i)**2 + vz(i)**2
2871 vis2(i) = viscf * vis2(i)
2872 fm2 = (fric*fni(i))**2
2873 f2 = vis2(i) * v2
2874 a2 = min(f2,fm2) / max(em30,f2)
2875 aa = sqrt(a2 * vis2(i))
2876 fx = aa * vx(i)
2877 fy = aa * vy(i)
2878 fz = aa * vz(i)
2879 fsav4 = fsav4 + fx*dt12
2880 fsav5 = fsav5 + fy*dt12
2881 fsav6 = fsav6 + fz*dt12
2882 fxi(i)=fxi(i) + fx
2883 fyi(i)=fyi(i) + fy
2884 fzi(i)=fzi(i) + fz
2885 econvt = econvt + aa * v2 * dt1
2886 ENDDO
2887 IF (imconv==1) THEN
2888#include "lockon.inc"
2889 fsav(4) = fsav(4) + fsav4
2890 fsav(5) = fsav(5) + fsav5
2891 fsav(6) = fsav(6) + fsav6
2892#include "lockoff.inc"
2893 ENDIF
2894 IF(isensint(1)/=0) THEN
2895 DO i=1,jlt
2896 fm2 = (fric*fni(i))**2
2897 f2 = vis2(i) * v2
2898 a2 = min(f2,fm2) / max(em30,f2)
2899 aa = sqrt(a2 * vis2(i))
2900 fsavparit(1,4,i+nft) = aa * vx(i)
2901 fsavparit(1,5,i+nft) = aa * vy(i)
2902 fsavparit(1,6,i+nft) = aa * vz(i)
2903 ENDDO
2904 ENDIF
2905 ENDIF
2906C
2907 IF (imconv==1) THEN
2908#include "lockon.inc"
2909 econtv = econtv + econvt
2910 econt = econt + econtt
2911#include "lockoff.inc"
2912 ENDIF
2913C---------------------------------
2914 DO i=1,jlt
2915 fx1(i)=-fxi(i)*hs1(i)
2916 fy1(i)=-fyi(i)*hs1(i)
2917 fz1(i)=-fzi(i)*hs1(i)
2918C
2919 fx2(i)=-fxi(i)*hs2(i)
2920 fy2(i)=-fyi(i)*hs2(i)
2921 fz2(i)=-fzi(i)*hs2(i)
2922C
2923 fx3(i)=fxi(i)*hm1(i)
2924 fy3(i)=fyi(i)*hm1(i)
2925 fz3(i)=fzi(i)*hm1(i)
2926C
2927 fx4(i)=fxi(i)*hm2(i)
2928 fy4(i)=fyi(i)*hm2(i)
2929 fz4(i)=fzi(i)*hm2(i)
2930C
2931 ENDDO
2932C
2933 IF (nspmd>1) THEN
2934Ctmp+1 mic only
2935#include "mic_lockon.inc"
2936 DO i = 1,jlt
2937 IF(cs_loc(i)>nlinsa)THEN
2938 ni = cs_loc(i)-nlinsa
2939C tag temporaire de NSVFI a -
2940 nsvfie(nin)%P(ni) = -abs(nsvfie(nin)%P(ni))
2941 ENDIF
2942 ENDDO
2943ctmp+1 mic only
2944#include "mic_lockoff.inc"
2945 ENDIF
2946C
2947 DO i=1,jlt
2948 stif(i) = two*stif(i)
2949 ENDDO
2950C
2951C---------------------------------
2952 IF(kdtint==1.OR.idtmins==2)THEN
2953 IF( (visc/=zero)
2954 . .AND.(ivis2==0.OR.ivis2==1))THEN
2955 DO i=1,jlt
2956 cx= c(i)*c(i)
2957C
2958 IF(ms1(i)==zero)THEN
2959 k1(i) =zero
2960 c1(i) =zero
2961 ELSE
2962 k1(i)=kt(i)*abs(hs1(i))
2963 c1(i)=c(i)*abs(hs1(i))
2964 cx =four*c1(i)*c1(i)
2965 cy =eight*ms1(i)*k1(i)
2966 aux = sqrt(cx+cy)+two*c1(i)
2967 st1(i)= k1(i)*aux*aux/max(cy,em30)
2968 cfi = cf(i)*abs(hs1(i))
2969 aux = two*cfi*cfi/max(ms1(i),em20)
2970 IF(aux>st1(i))THEN
2971 k1(i) =zero
2972 c1(i) =cfi
2973 ENDIF
2974 ENDIF
2975C
2976 IF(ms2(i)==zero)THEN
2977 k2(i) =zero
2978 c2(i) =zero
2979 ELSE
2980 k2(i)=kt(i)*abs(hs2(i))
2981 c2(i)=c(i)*abs(hs2(i))
2982 cx =four*c2(i)*c2(i)
2983 cy =eight*ms2(i)*k2(i)
2984 aux = sqrt(cx+cy)+two*c2(i)
2985 st2(i)= k2(i)*aux*aux/max(cy,em30)
2986 cfi = cf(i)*abs(hs2(i))
2987 aux = two*cfi*cfi/max(ms2(i),em20)
2988 IF(aux>st2(i))THEN
2989 k2(i) =zero
2990 c2(i) =cfi
2991 ENDIF
2992 ENDIF
2993C
2994 IF(mm1(i)==zero)THEN
2995 k3(i) =zero
2996 c3(i) =zero
2997 ELSE
2998 k3(i)=kt(i)*abs(hm1(i))
2999 c3(i)=c(i)*abs(hm1(i))
3000 cx =four*c3(i)*c3(i)
3001 cy =eight*mm1(i)*k3(i)
3002 aux = sqrt(cx+cy)+two*c3(i)
3003 st3(i)= k3(i)*aux*aux/max(cy,em30)
3004 cfi = cf(i)*abs(hm1(i))
3005 aux = two*cfi*cfi/max(mm1(i),em20)
3006 IF(aux>st3(i))THEN
3007 k3(i) =zero
3008 c3(i) =cfi
3009 ENDIF
3010 ENDIF
3011C
3012 IF(mm2(i)==zero)THEN
3013 k4(i) =zero
3014 c4(i) =zero
3015 ELSE
3016 k4(i)=kt(i)*abs(hm2(i))
3017 c4(i)=c(i)*abs(hm2(i))
3018 cx =four*c4(i)*c4(i)
3019 cy =eight*mm2(i)*k4(i)
3020 aux = sqrt(cx+cy)+two*c4(i)
3021 st4(i)= k4(i)*aux*aux/max(cy,em30)
3022 cfi = cf(i)*abs(hm2(i))
3023 aux = two*cfi*cfi/max(mm2(i),em20)
3024 IF(aux>st4(i))THEN
3025 k4(i) =zero
3026 c4(i) =cfi
3027 ENDIF
3028 ENDIF
3029 ENDDO
3030 ELSE
3031 DO i=1,jlt
3032 k1(i) =stif(i)*abs(hs1(i))
3033 c1(i) =zero
3034 k2(i) =stif(i)*abs(hs2(i))
3035 c2(i) =zero
3036 k3(i) =stif(i)*abs(hm1(i))
3037 c3(i) =zero
3038 k4(i) =stif(i)*abs(hm2(i))
3039 c4(i) =zero
3040 ENDDO
3041 ENDIF
3042 ENDIF
3043C=======================================================================
3044C FORCES PARITH ON sur noeud d'ancrage SECOND
3045C=======================================================================
3046 CALL foat_to_6_float(1 ,jlt ,fx1, fx6)
3047 CALL foat_to_6_float(1 ,jlt ,fy1, fy6)
3048 CALL foat_to_6_float(1 ,jlt ,fz1, fz6)
3049#include "lockon.inc"
3050 DO i = 1,jlt
3051 IF(cs_loc(i)<=nlinsa)THEN
3052 il = n1l(i)
3053C IG = N1(I)
3054C IF(IG > 0)THEN
3055 DO k = 1,6
3056 daanc6(1,k,il) = daanc6(1,k,il) + fx6(k,i)
3057 daanc6(2,k,il) = daanc6(2,k,il) + fy6(k,i)
3058 daanc6(3,k,il) = daanc6(3,k,il) + fz6(k,i)
3059 ENDDO
3060 ELSE
3061C SPMD remote SECONDARYs
3062C IL = - IG
3063 il = n1(i)
3064 DO k = 1,6
3065 daanc6fie(nin)%P(1,k,il) = daanc6fie(nin)%P(1,k,il)
3066 . + fx6(k,i)
3067 daanc6fie(nin)%P(2,k,il) = daanc6fie(nin)%P(2,k,il)
3068 . + fy6(k,i)
3069 daanc6fie(nin)%P(3,k,il) = daanc6fie(nin)%P(3,k,il)
3070 . + fz6(k,i)
3071 ENDDO
3072 ENDIF
3073 ENDDO
3074#include "lockoff.inc"
3075 CALL foat_to_6_float(1 ,jlt ,fx2, fx6)
3076 CALL foat_to_6_float(1 ,jlt ,fy2, fy6)
3077 CALL foat_to_6_float(1 ,jlt ,fz2, fz6)
3078#include "lockon.inc"
3079 DO i = 1,jlt
3080 IF(cs_loc(i)<=nlinsa)THEN
3081 il = n2l(i)
3082C IG = N2(I)
3083C IF(IG > 0)THEN
3084 DO k = 1,6
3085 daanc6(1,k,il) = daanc6(1,k,il) + fx6(k,i)
3086 daanc6(2,k,il) = daanc6(2,k,il) + fy6(k,i)
3087 daanc6(3,k,il) = daanc6(3,k,il) + fz6(k,i)
3088 ENDDO
3089 ELSE
3090C SPMD remote SECONDARYs
3091 il = n2(i)
3092C IL = - IG
3093 DO k = 1,6
3094 daanc6fie(nin)%P(1,k,il) = daanc6fie(nin)%P(1,k,il)
3095 . + fx6(k,i)
3096 daanc6fie(nin)%P(2,k,il) = daanc6fie(nin)%P(2,k,il)
3097 . + fy6(k,i)
3098 daanc6fie(nin)%P(3,k,il) = daanc6fie(nin)%P(3,k,il)
3099 . + fz6(k,i)
3100 ENDDO
3101 ENDIF
3102 ENDDO
3103#include "lockoff.inc"
3104C=======================================================================
3105C FORCES PARITH ON sur noeud d'ancrage MAIN
3106C=======================================================================
3107 CALL foat_to_6_float(1 ,jlt ,fx3, fx6)
3108 CALL foat_to_6_float(1 ,jlt ,fy3, fy6)
3109 CALL foat_to_6_float(1 ,jlt ,fz3, fz6)
3110#include "lockon.inc"
3111 DO i = 1,jlt
3112 il = m1l(i)
3113 DO k = 1,6
3114 daanc6(1,k,il) = daanc6(1,k,il) + fx6(k,i)
3115 daanc6(2,k,il) = daanc6(2,k,il) + fy6(k,i)
3116 daanc6(3,k,il) = daanc6(3,k,il) + fz6(k,i)
3117 ENDDO
3118 ENDDO
3119#include "lockoff.inc"
3120 CALL foat_to_6_float(1 ,jlt ,fx4, fx6)
3121 CALL foat_to_6_float(1 ,jlt ,fy4, fy6)
3122 CALL foat_to_6_float(1 ,jlt ,fz4, fz6)
3123#include "lockon.inc"
3124 DO i = 1,jlt
3125 il = m2l(i)
3126 DO k = 1,6
3127 daanc6(1,k,il) = daanc6(1,k,il) + fx6(k,i)
3128 daanc6(2,k,il) = daanc6(2,k,il) + fy6(k,i)
3129 daanc6(3,k,il) = daanc6(3,k,il) + fz6(k,i)
3130 ENDDO
3131 ENDDO
3132#include "lockoff.inc"
3133C=======================================================================
3134C mise a ZERO des FORCES sur noeuds maites et second
3135C si PENE (su noeud second) < GAPR (gap reel)
3136C=======================================================================
3137C=======================================================================
3138C FORCES sur noeuds maites et second
3139C=======================================================================
3140C---------------------------------
3141 IF(iparit==0)THEN
3142 IF(kdtint==0)THEN
3143 CALL i20ass0(jlt ,cs_loc,n1 ,n2 ,m1 ,
3144 2 m2 ,hs1 ,hs2 ,hm1 ,hm2 ,
3145 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,
3146 4 fz2 ,fx3 ,fy3 ,fz3 ,fx4 ,
3147 5 fy4 ,fz4 ,a ,stifn,stif ,
3148 6 nlinsa,nin ,jtask)
3149 ELSE
3150 CALL i20ass05(jlt ,cs_loc,n1 ,n2 ,m1 ,
3151 2 m2 ,hs1 ,hs2 ,hm1 ,hm2 ,
3152 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,
3153 4 fz2 ,fx3 ,fy3 ,fz3 ,fx4 ,
3154 5 fy4 ,fz4 ,a ,stifn,nlinsa,
3155 6 k1 ,k2 ,k3 ,k4 ,c1 ,
3156 7 c2 ,c3 ,c4 ,viscn,nin ,jtask )
3157 END IF
3158 ELSE
3159 IF(kdtint==0)THEN
3160 CALL i20ass2(jlt ,cs_loc ,n1 ,n2 ,m1 ,
3161 2 m2 ,hs1 ,hs2 ,hm1 ,hm2 ,
3162 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,
3163 4 fz2 ,fx3 ,fy3 ,fz3 ,fx4 ,
3164 5 fy4 ,fz4 ,fskyi ,isky ,niskyfie,
3165 6 stif ,nlinsa ,nin ,noint )
3166 ELSE
3167 CALL i20ass25(jlt ,cs_loc ,n1 ,n2 ,m1 ,
3168 2 m2 ,hs1 ,hs2 ,hm1 ,hm2 ,
3169 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,
3170 4 fz2 ,fx3 ,fy3 ,fz3 ,fx4 ,
3171 5 fy4 ,fz4 ,isky ,niskyfie,nlinsa ,
3172 6 k1 ,k2 ,k3 ,k4 ,c1 ,
3173 7 c2 ,c3 ,c4 ,nin , noint)
3174 END IF
3175 END IF
3176C
3177 IF(idtmins==2)
3178 . CALL i20sms2e(jlt ,cs_loc ,n1 ,n2 ,m1 ,
3179 2 m2 ,hs1 ,hs2 ,hm1 ,hm2 ,
3180 3 stif ,nin ,noint ,mskyi_sms ,iskyi_sms,
3181 4 nsms ,k1 ,k2 ,k3 ,k4 ,
3182 5 c1 ,c2 ,c3 ,c4 ,nlinsa )
3183C
3184 IF(idtmin(10)==1.OR.idtmin(10)==2)THEN
3185 dtmi0 = ep20
3186 DO i=1,jlt
3187 dtmi(i) = ep20
3188 mas2 = two * masmin(i)
3189 IF(mas2>zero.AND.stif(i)>zero)THEN
3190 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/stif(i)))
3191 ENDIF
3192 dtmi0 = min(dtmi0,dtmi(i))
3193 ENDDO
3194 IF(dtmi0<=dtmin1(10))THEN
3195 DO i=1,jlt
3196 IF(dtmi(i)<=dtmin1(10))THEN
3197 IF(idtmin(10)==1)THEN
3198#include "lockon.inc"
3199 IF(cs_loc(i)<=nlinsa) THEN
3200 WRITE(iout,*)
3201 . ' **WARNING MINIMUM TIME STEP ',dtmi(i),
3202 . ' IN INTERFACE NB',noint
3203 WRITE(iout,*)'SECONDARY NODES NB',itab(n1(i)),
3204 . itab(n2(i))
3205 WRITE(iout,*)'MAIN NODES NB',itab(m1(i)),
3206 . itab(m2(i))
3207 ELSE
3208 WRITE(iout,*)
3209 . ' **WARNING MINIMUM TIME STEP ',dtmi(i),
3210 . ' IN INTERFACE NB',noint
3211 WRITE(iout,*)'SECONDARY NODES NB',itafie(nin)%P(n1(i)),
3212 . itafie(nin)%P(n2(i))
3213 WRITE(iout,*)'MAIN NODES NB',itab(m1(i)),
3214 . itab(m2(i))
3215 END IF
3216#include "lockoff.inc"
3217 tstop = tt
3218 ELSEIF(idtmin(10)==2)THEN
3219#include "lockon.inc"
3220 IF(cs_loc(i)<=nlinsa) THEN
3221 WRITE(iout,*)
3222 . ' **WARNING MINIMUM TIME STEP ',dtmi(i),
3223 . ' IN INTERFACE NB',noint
3224 WRITE(iout,*)'SECONDARY NODES NB',itab(n1(i)),
3225 . itab(n2(i))
3226 WRITE(iout,*)'MAIN NODES NB',itab(m1(i)),
3227 . itab(m2(i))
3228 WRITE(iout,*)'DELETE SECONDARY LINE FROM INTERFACE'
3229 stfs(cs_loc(i)) = -abs(stfs(cs_loc(i)))
3230 ELSE
3231 ni = cs_loc(i)-nlinsa
3232 WRITE(iout,*)
3233 . ' **WARNING MINIMUM TIME STEP ',dtmi(i),
3234 . ' IN INTERFACE NB',noint
3235 WRITE(iout,*)'SECONDARY NODES NB',itafie(nin)%P(n1(i)),
3236 . itafie(nin)%P(n2(i))
3237 WRITE(iout,*)'MAIN NODES NB',itab(m1(i)),
3238 . itab(m2(i))
3239 WRITE(iout,*)'DELETE SECONDARY LINE FROM INTERFACE'
3240 stifie(nin)%P(ni) = -abs(stifie(nin)%P(ni))
3241 END IF
3242#include "lockoff.inc"
3243 newfront = -1
3244 ELSEIF(idtmin(10)==5)THEN
3245#include "lockon.inc"
3246 IF(cs_loc(i)<=nlinsa) THEN
3247 WRITE(iout,*)
3248 . ' **WARNING MINIMUM TIME STEP ',dtmi(i),
3249 . ' IN INTERFACE NB',noint
3250 WRITE(iout,*)'SECONDARY NODES NB',itab(n1(i)),
3251 . itab(n2(i))
3252 WRITE(iout,*)'MAIN NODES NB',itab(m1(i)),
3253 . itab(m2(i))
3254 ELSE
3255 WRITE(iout,*)
3256 . ' **WARNING MINIMUM TIME STEP ',dtmi(i),
3257 . ' IN INTERFACE NB',noint
3258 WRITE(iout,*)'SECONDARY NODES NB',itafie(nin)%P(n1(i)),
3259 . itafie(nin)%P(n2(i))
3260 WRITE(iout,*)'MAIN NODES NB',itab(m1(i)),
3261 . itab(m2(i))
3262 END IF
3263#include "lockoff.inc"
3264 mstop = 2
3265 ENDIF
3266 ENDIF
3267 ENDDO
3268 ENDIF
3269 ENDIF
3270C
3271 IF(anim_v(4)+outp_v(4)+h3d_data%N_VECT_CONT>0)THEN
3272#include "lockon.inc"
3273c goto 1234
3274 DO i=1,jlt
3275 IF(cs_loc(i)<=nlinsa) THEN
3276 fcont(1,n1(i)) =fcont(1,n1(i)) + fx1(i)
3277 fcont(2,n1(i)) =fcont(2,n1(i)) + fy1(i)
3278 fcont(3,n1(i)) =fcont(3,n1(i)) + fz1(i)
3279 fcont(1,n2(i)) =fcont(1,n2(i)) + fx2(i)
3280 fcont(2,n2(i)) =fcont(2,n2(i)) + fy2(i)
3281 fcont(3,n2(i)) =fcont(3,n2(i)) + fz2(i)
3282 END IF
3283 fcont(1,m1(i)) =fcont(1,m1(i)) + fx3(i)
3284 fcont(2,m1(i)) =fcont(2,m1(i)) + fy3(i)
3285 fcont(3,m1(i)) =fcont(3,m1(i)) + fz3(i)
3286 fcont(1,m2(i)) =fcont(1,m2(i)) + fx4(i)
3287 fcont(2,m2(i)) =fcont(2,m2(i)) + fy4(i)
3288 fcont(3,m2(i)) =fcont(3,m2(i)) + fz4(i)
3289 ENDDO
3290c 1234 continue
3291#include "lockoff.inc"
3292 ENDIF
3293C
3294 IF(isecin>0)THEN
3295 k0=nstrf(25)
3296 IF(nstrf(1)+nstrf(2)/=0)THEN
3297 DO i=1,nsect
3298 nbinter=nstrf(k0+14)
3299 k1s=k0+30
3300 DO j=1,nbinter
3301 IF(nstrf(k1s)==noint)THEN
3302 IF(isecut/=0)THEN
3303#include "lockon.inc"
3304 DO k=1,jlt
3305 IF(cs_loc(i)<=nlinsa) THEN
3306 IF(secfcum(4,n1(k),i)==1.)THEN
3307 secfcum(1,n1(k),i)=secfcum(1,n1(k),i)-fx1(k)
3308 secfcum(2,n1(k),i)=secfcum(2,n1(k),i)-fy1(k)
3309 secfcum(3,n1(k),i)=secfcum(3,n1(k),i)-fz1(k)
3310 ENDIF
3311 IF(secfcum(4,n2(k),i)==1.)THEN
3312 secfcum(1,n2(k),i)=secfcum(1,n2(k),i)-fx2(k)
3313 secfcum(2,n2(k),i)=secfcum(2,n2(k),i)-fy2(k)
3314 secfcum(3,n2(k),i)=secfcum(3,n2(k),i)-fz2(k)
3315 ENDIF
3316 END IF
3317 IF(secfcum(4,m1(k),i)==1.)THEN
3318 secfcum(1,m1(k),i)=secfcum(1,m1(k),i)-fx3(k)
3319 secfcum(2,m1(k),i)=secfcum(2,m1(k),i)-fy3(k)
3320 secfcum(3,m1(k),i)=secfcum(3,m1(k),i)-fz3(k)
3321 ENDIF
3322 IF(secfcum(4,m2(k),i)==1.)THEN
3323 secfcum(1,m2(k),i)=secfcum(1,m2(k),i)-fx4(k)
3324 secfcum(2,m2(k),i)=secfcum(2,m2(k),i)-fy4(k)
3325 secfcum(3,m2(k),i)=secfcum(3,m2(k),i)-fz4(k)
3326 ENDIF
3327 ENDDO
3328#include "lockoff.inc"
3329 ENDIF
3330C +fsav(section)
3331 ENDIF
3332 k1s=k1s+1
3333 ENDDO
3334 k0=nstrf(k0+24)
3335 ENDDO
3336 ENDIF
3337 ENDIF
3338C
3339 RETURN
3340 END
3341C
3342!||====================================================================
3343!|| i20ass0 ../engine/source/interfaces/int20/i20for3.F
3344!||--- called by ------------------------------------------------------
3345!|| i20for3e ../engine/source/interfaces/int20/i20for3.F
3346!||--- uses -----------------------------------------------------
3347!|| tri7box ../engine/share/modules/tri7box.F
3348!||====================================================================
3349 SUBROUTINE i20ass0(JLT ,CS_LOC,N1 ,N2 ,M1 ,
3350 2 M2 ,HS1 ,HS2 ,HM1 ,HM2 ,
3351 3 FX1 ,FY1 ,FZ1 ,FX2 ,FY2 ,
3352 4 FZ2 ,FX3 ,FY3 ,FZ3 ,FX4 ,
3353 5 FY4 ,FZ4 ,A ,STIFN,STIF ,
3354 6 NRTS ,NIN,JTASK)
3355C-----------------------------------------------
3356C M o d u l e s
3357C-----------------------------------------------
3358 USE tri7box
3359C-----------------------------------------------
3360C I m p l i c i t T y p e s
3361C-----------------------------------------------
3362#include "implicit_f.inc"
3363C-----------------------------------------------
3364C G l o b a l P a r a m e t e r s
3365C-----------------------------------------------
3366#include "mvsiz_p.inc"
3367C-----------------------------------------------
3368C D u m m y A r g u m e n t s
3369C-----------------------------------------------
3370 INTEGER JLT, NRTS, NIN,
3371 + cs_loc(*),
3372 + n1(mvsiz),n2(mvsiz),m1(mvsiz),m2(mvsiz),jtask
3373 my_real
3374 . hs1(mvsiz),hs2(mvsiz),hm1(mvsiz),hm2(mvsiz),
3375 . fx1(mvsiz),fy1(mvsiz),fz1(mvsiz),
3376 . fx2(mvsiz),fy2(mvsiz),fz2(mvsiz),
3377 . fx3(mvsiz),fy3(mvsiz),fz3(mvsiz),
3378 . fx4(mvsiz),fy4(mvsiz),fz4(mvsiz),
3379 . a(3,*), stifn(*), stif(*)
3380C-----------------------------------------------
3381C L o c a l V a r i a b l e s
3382C-----------------------------------------------
3383 INTEGER I, J1,NODFI,ISHIFT
3384C-----------------------------------------------
3385C
3386 nodfi = nlskyfi(nin)
3387 ishift = nodfi*(jtask-1)
3388C
3389 DO i=1,jlt
3390 IF(cs_loc(i)<=nrts) THEN
3391 j1=n1(i)
3392 a(1,j1)=a(1,j1)+fx1(i)
3393 a(2,j1)=a(2,j1)+fy1(i)
3394 a(3,j1)=a(3,j1)+fz1(i)
3395 stifn(j1) = stifn(j1) + stif(i)*abs(hs1(i))
3396C
3397 j1=n2(i)
3398 a(1,j1)=a(1,j1)+fx2(i)
3399 a(2,j1)=a(2,j1)+fy2(i)
3400 a(3,j1)=a(3,j1)+fz2(i)
3401 stifn(j1) = stifn(j1) + stif(i)*abs(hs2(i))
3402 ELSE
3403 j1=n1(i)
3404 afie(nin)%P(1,j1+ishift)=afie(nin)%P(1,j1+ishift)+fx1(i)
3405 afie(nin)%P(2,j1+ishift)=afie(nin)%P(2,j1+ishift)+fy1(i)
3406 afie(nin)%P(3,j1+ishift)=afie(nin)%P(3,j1+ishift)+fz1(i)
3407 stnfie(nin)%P(j1+ishift) = stnfie(nin)%P(j1+ishift) + stif(i)*abs(hs1(i))
3408C
3409 j1=n2(i)
3410 afie(nin)%P(1,j1+ishift)=afie(nin)%P(1,j1+ishift)+fx2(i)
3411 afie(nin)%P(2,j1+ishift)=afie(nin)%P(2,j1+ishift)+fy2(i)
3412 afie(nin)%P(3,j1+ishift)=afie(nin)%P(3,j1+ishift)+fz2(i)
3413 stnfie(nin)%P(j1+ishift) = stnfie(nin)%P(j1+ishift) + stif(i)*abs(hs2(i))
3414 END IF
3415 END DO
3416C
3417 DO i=1,jlt
3418 j1=m1(i)
3419 a(1,j1)=a(1,j1)+fx3(i)
3420 a(2,j1)=a(2,j1)+fy3(i)
3421 a(3,j1)=a(3,j1)+fz3(i)
3422 stifn(j1) = stifn(j1) + stif(i)*abs(hm1(i))
3423C
3424 j1=m2(i)
3425 a(1,j1)=a(1,j1)+fx4(i)
3426 a(2,j1)=a(2,j1)+fy4(i)
3427 a(3,j1)=a(3,j1)+fz4(i)
3428 stifn(j1) = stifn(j1) + stif(i)*abs(hm2(i))
3429 ENDDO
3430C
3431 RETURN
3432 END
3433C
3434!||====================================================================
3435!|| i20ass05 ../engine/source/interfaces/int20/i20for3.F
3436!||--- called by ------------------------------------------------------
3437!|| i20for3e ../engine/source/interfaces/int20/i20for3.F
3438!||--- uses -----------------------------------------------------
3439!|| tri7box ../engine/share/modules/tri7box.F
3440!||====================================================================
3441 SUBROUTINE i20ass05(JLT ,CS_LOC,N1 ,N2 ,M1 ,
3442 2 M2 ,HS1 ,HS2 ,HM1 ,HM2 ,
3443 3 FX1 ,FY1 ,FZ1 ,FX2 ,FY2 ,
3444 4 FZ2 ,FX3 ,FY3 ,FZ3 ,FX4 ,
3445 5 FY4 ,FZ4 ,A ,STIFN,NRTS ,
3446 6 K1 ,K2 ,K3 ,K4 ,C1 ,
3447 7 C2 ,C3 ,C4 ,VISCN,NIN ,JTASK )
3448C-----------------------------------------------
3449C M o d u l e s
3450C-----------------------------------------------
3451 USE tri7box
3452C-----------------------------------------------
3453C I m p l i c i t T y p e s
3454C-----------------------------------------------
3455#include "implicit_f.inc"
3456C-----------------------------------------------
3457C G l o b a l P a r a m e t e r s
3458C-----------------------------------------------
3459#include "mvsiz_p.inc"
3460C-----------------------------------------------
3461C D u m m y A r g u m e n t s
3462C-----------------------------------------------
3463 INTEGER JLT, NRTS, NIN,
3464 + cs_loc(*),
3465 + n1(mvsiz),n2(mvsiz),m1(mvsiz),m2(mvsiz),jtask
3466 my_real
3467 . hs1(mvsiz),hs2(mvsiz),hm1(mvsiz),hm2(mvsiz),
3468 . fx1(mvsiz),fy1(mvsiz),fz1(mvsiz),
3469 . fx2(mvsiz),fy2(mvsiz),fz2(mvsiz),
3470 . fx3(mvsiz),fy3(mvsiz),fz3(mvsiz),
3471 . fx4(mvsiz),fy4(mvsiz),fz4(mvsiz),
3472 . k1(mvsiz),k2(mvsiz),k3(mvsiz),k4(mvsiz),
3473 . c1(mvsiz),c2(mvsiz),c3(mvsiz),c4(mvsiz),
3474 . a(3,*), stifn(*), viscn(*)
3475C-----------------------------------------------
3476C L o c a l V a r i a b l e s
3477C-----------------------------------------------
3478 INTEGER I, J1,NODFI,ISHIFT
3479C-----------------------------------------------
3480C
3481 nodfi = nlskyfi(nin)
3482 ishift = nodfi*(jtask-1)
3483C
3484 DO i=1,jlt
3485 IF(cs_loc(i)<=nrts) THEN
3486 j1=n1(i)
3487 a(1,j1)=a(1,j1)+fx1(i)
3488 a(2,j1)=a(2,j1)+fy1(i)
3489 a(3,j1)=a(3,j1)+fz1(i)
3490 stifn(j1)=stifn(j1)+k1(i)
3491 viscn(j1)=viscn(j1)+c1(i)
3492C
3493 j1=n2(i)
3494 a(1,j1)=a(1,j1)+fx2(i)
3495 a(2,j1)=a(2,j1)+fy2(i)
3496 a(3,j1)=a(3,j1)+fz2(i)
3497 stifn(j1)=stifn(j1)+k2(i)
3498 viscn(j1)=viscn(j1)+c2(i)
3499 ELSE
3500 j1=n1(i)
3501 afie(nin)%P(1,j1+ishift)=afie(nin)%P(1,j1+ishift)+fx1(i)
3502 afie(nin)%P(2,j1+ishift)=afie(nin)%P(2,j1+ishift)+fy1(i)
3503 afie(nin)%P(3,j1+ishift)=afie(nin)%P(3,j1+ishift)+fz1(i)
3504 stnfie(nin)%P(j1+ishift)=stnfie(nin)%P(j1+ishift)+k1(i)
3505 vscfie(nin)%P(j1+ishift)=vscfie(nin)%P(j1+ishift)+c1(i)
3506C
3507 j1=n2(i)
3508 afie(nin)%P(1,j1+ishift)=afie(nin)%P(1,j1+ishift)+fx2(i)
3509 afie(nin)%P(2,j1+ishift)=afie(nin)%P(2,j1+ishift)+fy2(i)
3510 afie(nin)%P(3,j1+ishift)=afie(nin)%P(3,j1+ishift)+fz2(i)
3511 stnfie(nin)%P(j1+ishift)=stnfie(nin)%P(j1+ishift)+k2(i)
3512 vscfie(nin)%P(j1+ishift)=vscfie(nin)%P(j1+ishift)+c2(i)
3513 END IF
3514 END DO
3515C
3516 DO i=1,jlt
3517 j1=m1(i)
3518 a(1,j1)=a(1,j1)+fx3(i)
3519 a(2,j1)=a(2,j1)+fy3(i)
3520 a(3,j1)=a(3,j1)+fz3(i)
3521 stifn(j1)=stifn(j1)+k3(i)
3522 viscn(j1)=viscn(j1)+c3(i)
3523C
3524 j1=m2(i)
3525 a(1,j1)=a(1,j1)+fx4(i)
3526 a(2,j1)=a(2,j1)+fy4(i)
3527 a(3,j1)=a(3,j1)+fz4(i)
3528 stifn(j1)=stifn(j1)+k4(i)
3529 viscn(j1)=viscn(j1)+c4(i)
3530 ENDDO
3531C
3532 RETURN
3533 END
3534C
3535!||====================================================================
3536!|| i20ass2 ../engine/source/interfaces/int20/i20for3.F
3537!||--- called by ------------------------------------------------------
3538!|| i20for3e ../engine/source/interfaces/int20/i20for3.F
3539!||--- calls -----------------------------------------------------
3540!|| ancmsg ../engine/source/output/message/message.F
3541!|| arret ../engine/source/system/arret.F
3542!||--- uses -----------------------------------------------------
3543!|| message_mod ../engine/share/message_module/message_mod.F
3544!|| tri7box ../engine/share/modules/tri7box.F
3545!||====================================================================
3546 SUBROUTINE i20ass2(JLT ,CS_LOC ,N1 ,N2 ,M1 ,
3547 2 M2 ,HS1 ,HS2 ,HM1 ,HM2 ,
3548 3 FX1 ,FY1 ,FZ1 ,FX2 ,FY2 ,
3549 4 FZ2 ,FX3 ,FY3 ,FZ3 ,FX4 ,
3550 5 FY4 ,FZ4 ,FSKYI ,ISKY ,NISKYFIE,
3551 6 STIF ,NRTS ,NIN ,NOINT )
3552C-----------------------------------------------
3553C M o d u l e s
3554C-----------------------------------------------
3555 USE tri7box
3556 USE message_mod
3557C-----------------------------------------------
3558C I m p l i c i t T y p e s
3559C-----------------------------------------------
3560#include "implicit_f.inc"
3561#include "comlock.inc"
3562C-----------------------------------------------
3563C G l o b a l P a r a m e t e r s
3564C-----------------------------------------------
3565#include "mvsiz_p.inc"
3566C-----------------------------------------------
3567C C o m m o n B l o c k s
3568C-----------------------------------------------
3569#include "parit_c.inc"
3570C-----------------------------------------------
3571C D u m m y A r g u m e n t s
3572C-----------------------------------------------
3573 INTEGER JLT, NRTS,NISKYFIE,NIN,NOINT,
3574 + cs_loc(*),isky(*),
3575 + n1(mvsiz),n2(mvsiz),m1(mvsiz),m2(mvsiz)
3576 my_real
3577 . hs1(mvsiz),hs2(mvsiz),hm1(mvsiz),hm2(mvsiz),
3578 . fx1(mvsiz),fy1(mvsiz),fz1(mvsiz),
3579 . fx2(mvsiz),fy2(mvsiz),fz2(mvsiz),
3580 . fx3(mvsiz),fy3(mvsiz),fz3(mvsiz),
3581 . fx4(mvsiz),fy4(mvsiz),fz4(mvsiz),
3582 . fskyi(lskyi,nfskyi), stif(*)
3583C-----------------------------------------------
3584C L o c a l V a r i a b l e s
3585C-----------------------------------------------
3586 INTEGER I, J1, NISKYL1, NISKYL,IGP,IGM, NISKYFIEL
3587C
3588 niskyl1 = 0
3589 DO i = 1, jlt
3590 IF (hm1(i)/=zero) niskyl1 = niskyl1 + 1
3591 ENDDO
3592 DO i = 1, jlt
3593 IF (hm2(i)/=zero) niskyl1 = niskyl1 + 1
3594 ENDDO
3595
3596 igp = 0
3597 igm = 0
3598 DO i=1,jlt
3599 IF(cs_loc(i)<=nrts) THEN
3600 igp = igp+2
3601 ELSE
3602 igm = igm+1
3603 ENDIF
3604 ENDDO
3605
3606#include "lockon.inc"
3607 niskyl = nisky
3608 nisky = nisky + niskyl1 + igp
3609 niskyfiel = niskyfie
3610 niskyfie = niskyfie + igm
3611#include "lockoff.inc"
3612
3613 IF (niskyl+niskyl1+igp > lskyi) THEN
3614 CALL ancmsg(msgid=26,anmode=aninfo)
3615 CALL arret(2)
3616 ENDIF
3617 IF (niskyfiel+igm > nlskyfie(nin)) THEN
3618 CALL ancmsg(msgid=26,anmode=aninfo)
3619 CALL arret(2)
3620 ENDIF
3621C
3622 DO i=1,jlt
3623 IF(cs_loc(i)<=nrts) THEN
3624 niskyl = niskyl + 1
3625 fskyi(niskyl,1)=fx1(i)
3626 fskyi(niskyl,2)=fy1(i)
3627 fskyi(niskyl,3)=fz1(i)
3628 fskyi(niskyl,4)=stif(i)*abs(hs1(i))
3629 isky(niskyl) = n1(i)
3630C
3631 niskyl = niskyl + 1
3632 fskyi(niskyl,1)=fx2(i)
3633 fskyi(niskyl,2)=fy2(i)
3634 fskyi(niskyl,3)=fz2(i)
3635 fskyi(niskyl,4)=stif(i)*abs(hs2(i))
3636 isky(niskyl) = n2(i)
3637 ELSE
3638 niskyfiel = niskyfiel + 1
3639 fskyfie(nin)%P(1,niskyfiel)=fx1(i)
3640 fskyfie(nin)%P(2,niskyfiel)=fy1(i)
3641 fskyfie(nin)%P(3,niskyfiel)=fz1(i)
3642 fskyfie(nin)%P(4,niskyfiel)=stif(i)*abs(hs1(i))
3643 fskyfie(nin)%P(5,niskyfiel)=fx2(i)
3644 fskyfie(nin)%P(6,niskyfiel)=fy2(i)
3645 fskyfie(nin)%P(7,niskyfiel)=fz2(i)
3646 fskyfie(nin)%P(8,niskyfiel)=stif(i)*abs(hs2(i))
3647 iskyfie(nin)%P(niskyfiel) = cs_loc(i)-nrts
3648 END IF
3649 END DO
3650C
3651 DO i=1,jlt
3652 IF (hm1(i)/=zero) THEN
3653 niskyl = niskyl + 1
3654 fskyi(niskyl,1)=fx3(i)
3655 fskyi(niskyl,2)=fy3(i)
3656 fskyi(niskyl,3)=fz3(i)
3657 fskyi(niskyl,4)=stif(i)*abs(hm1(i))
3658 isky(niskyl) = m1(i)
3659 ENDIF
3660 ENDDO
3661 DO i=1,jlt
3662 IF (hm2(i)/=zero) THEN
3663 niskyl = niskyl + 1
3664 fskyi(niskyl,1)=fx4(i)
3665 fskyi(niskyl,2)=fy4(i)
3666 fskyi(niskyl,3)=fz4(i)
3667 fskyi(niskyl,4)=stif(i)*abs(hm2(i))
3668 isky(niskyl) = m2(i)
3669 ENDIF
3670 ENDDO
3671C
3672 RETURN
3673 END
3674C
3675!||====================================================================
3676!|| i20ass25 ../engine/source/interfaces/int20/i20for3.F
3677!||--- called by ------------------------------------------------------
3678!|| i20for3e ../engine/source/interfaces/int20/i20for3.F
3679!||--- calls -----------------------------------------------------
3680!|| ancmsg ../engine/source/output/message/message.F
3681!|| arret ../engine/source/system/arret.F
3682!||--- uses -----------------------------------------------------
3683!|| message_mod ../engine/share/message_module/message_mod.F
3684!|| tri7box ../engine/share/modules/tri7box.F
3685!||====================================================================
3686 SUBROUTINE i20ass25(JLT ,CS_LOC ,N1 ,N2 ,M1 ,
3687 2 M2 ,HS1 ,HS2 ,HM1 ,HM2 ,
3688 3 FX1 ,FY1 ,FZ1 ,FX2 ,FY2 ,
3689 4 FZ2 ,FX3 ,FY3 ,FZ3 ,FX4 ,
3690 5 FY4 ,FZ4 ,ISKY ,NISKYFIE,NRTS ,
3691 6 K1 ,K2 ,K3 ,K4 ,C1 ,
3692 7 C2 ,C3 ,C4 ,NIN ,NOINT )
3693C-----------------------------------------------
3694C M o d u l e s
3695C-----------------------------------------------
3696 USE tri7box
3697 USE message_mod
3698C-----------------------------------------------
3699C I m p l i c i t T y p e s
3700C-----------------------------------------------
3701#include "implicit_f.inc"
3702#include "comlock.inc"
3703C-----------------------------------------------
3704C G l o b a l P a r a m e t e r s
3705C-----------------------------------------------
3706#include "mvsiz_p.inc"
3707C-----------------------------------------------
3708C C o m m o n B l o c k s
3709C-----------------------------------------------
3710#include "parit_c.inc"
3711C-----------------------------------------------
3712C D u m m y A r g u m e n t s
3713C-----------------------------------------------
3714 INTEGER JLT, NRTS,NISKYFIE,NIN,NOINT,
3715 + cs_loc(*),isky(*),
3716 + n1(mvsiz),n2(mvsiz),m1(mvsiz),m2(mvsiz)
3717 my_real
3718 . hs1(mvsiz),hs2(mvsiz),hm1(mvsiz),hm2(mvsiz),
3719 . fx1(mvsiz),fy1(mvsiz),fz1(mvsiz),
3720 . fx2(mvsiz),fy2(mvsiz),fz2(mvsiz),
3721 . fx3(mvsiz),fy3(mvsiz),fz3(mvsiz),
3722 . fx4(mvsiz),fy4(mvsiz),fz4(mvsiz),
3723 . k1(mvsiz),k2(mvsiz),k3(mvsiz),k4(mvsiz),
3724 . c1(mvsiz),c2(mvsiz),c3(mvsiz),c4(mvsiz),
3725 . fskyi(lskyi,nfskyi)
3726C-----------------------------------------------
3727C L o c a l V a r i a b l e s
3728C-----------------------------------------------
3729 INTEGER I, J1, NISKYL1, NISKYL,IGP,IGM, NISKYFIEL
3730C
3731 niskyl1 = 0
3732 DO i = 1, jlt
3733 IF (hm1(i)/=zero) niskyl1 = niskyl1 + 1
3734 ENDDO
3735 DO i = 1, jlt
3736 IF (hm2(i)/=zero) niskyl1 = niskyl1 + 1
3737 ENDDO
3738
3739 igp = 0
3740 igm = 0
3741 DO i=1,jlt
3742 IF(cs_loc(i)<=nrts) THEN
3743 igp = igp+2
3744 ELSE
3745 igm = igm+1
3746 ENDIF
3747 ENDDO
3748
3749#include "lockon.inc"
3750 niskyl = nisky
3751 nisky = nisky + niskyl1 + igp
3752 niskyfiel = niskyfie
3753 niskyfie = niskyfie + igm
3754#include "lockoff.inc"
3755C
3756 IF (niskyl+niskyl1+igp > lskyi) THEN
3757 CALL ancmsg(msgid=26,anmode=aninfo)
3758 CALL arret(2)
3759 ENDIF
3760 IF (niskyfiel+igm > nlskyfie(nin)) THEN
3761 CALL ancmsg(msgid=26,anmode=aninfo)
3762 CALL arret(2)
3763 ENDIF
3764C
3765 DO i=1,jlt
3766 IF(cs_loc(i)<=nrts) THEN
3767 niskyl = niskyl + 1
3768 fskyi(niskyl,1)=fx1(i)
3769 fskyi(niskyl,2)=fy1(i)
3770 fskyi(niskyl,3)=fz1(i)
3771 fskyi(niskyl,4)=k1(i)
3772 fskyi(niskyl,5)=c1(i)
3773 isky(niskyl) = n1(i)
3774C
3775 niskyl = niskyl + 1
3776 fskyi(niskyl,1)=fx2(i)
3777 fskyi(niskyl,2)=fy2(i)
3778 fskyi(niskyl,3)=fz2(i)
3779 fskyi(niskyl,4)=k2(i)
3780 fskyi(niskyl,5)=c2(i)
3781 isky(niskyl) = n2(i)
3782 ELSE
3783 niskyfiel = niskyfiel + 1
3784 fskyfie(nin)%P(1,niskyfiel)=fx1(i)
3785 fskyfie(nin)%P(2,niskyfiel)=fy1(i)
3786 fskyfie(nin)%P(3,niskyfiel)=fz1(i)
3787 fskyfie(nin)%P(4,niskyfiel)=k1(i)
3788 fskyfie(nin)%P(5,niskyfiel)=c1(i)
3789 fskyfie(nin)%P(6,niskyfiel)=fx2(i)
3790 fskyfie(nin)%P(7,niskyfiel)=fy2(i)
3791 fskyfie(nin)%P(8,niskyfiel)=fz2(i)
3792 fskyfie(nin)%P(9,niskyfiel)=k2(i)
3793 fskyfie(nin)%P(10,niskyfiel)=c2(i)
3794 iskyfie(nin)%P(niskyfiel) = cs_loc(i)-nrts
3795 END IF
3796 END DO
3797C
3798 DO i=1,jlt
3799 IF (hm1(i)/=zero) THEN
3800 niskyl = niskyl + 1
3801 fskyi(niskyl,1)=fx3(i)
3802 fskyi(niskyl,2)=fy3(i)
3803 fskyi(niskyl,3)=fz3(i)
3804 fskyi(niskyl,4)=k3(i)
3805 fskyi(niskyl,5)=c3(i)
3806 isky(niskyl) = m1(i)
3807 ENDIF
3808 ENDDO
3809 DO i=1,jlt
3810 IF (hm2(i)/=zero) THEN
3811 niskyl = niskyl + 1
3812 fskyi(niskyl,1)=fx4(i)
3813 fskyi(niskyl,2)=fy4(i)
3814 fskyi(niskyl,3)=fz4(i)
3815 fskyi(niskyl,4)=k4(i)
3816 fskyi(niskyl,5)=c4(i)
3817 isky(niskyl) = m2(i)
3818 ENDIF
3819 ENDDO
3820C
3821 RETURN
3822 END
#define my_real
Definition cppsort.cpp:32
#define alpha
Definition eval.h:35
subroutine area(d1, x, x2, y, y2, eint, stif0)
subroutine i20for3e(jlt, a, v, ibc, icodt, fsav, gap, fric, ms, visc, viscf, noint, itab, cs_loc, cm_loc, stiglo, stifn, stif, fskyi, isky, fcont, stfs, stfm, dt2t, hs1, hs2, hm1, hm2, n1, n2, m1, m2, ivis2, neltst, ityptst, nx, ny, nz, gapv, penise, penime, inacti, niskyfie, newfront, isecin, nstrf, secfcum, viscn, nlinsa, ms1, ms2, mm1, mm2, vxs1, vys1, vzs1, vxs2, vys2, vzs2, vxm1, vym1, vzm1, vxm2, vym2, vzm2, nin, n1l, n2l, m1l, m2l, daanc6, alphak, mskyi_sms, iskyi_sms, nsms, jtask, isensint, fsavparit, nisub, nft, h3d_data)
Definition i20for3.F:2413
subroutine i20for3(jlt, a, va, ibcc, icodt, fsav, gap, fric, ms, visc, viscf, noint, stfa, itab, cn_loc, stiglo, stifn, stif, fskyi, isky, nx1, nx2, nx3, nx4, ny1, ny2, ny3, ny4, nz1, nz2, nz3, nz4, lb1, lb2, lb3, lb4, lc1, lc2, lc3, lc4, p1, p2, p3, p4, fcont, ix1l, ix2l, ix3l, ix4l, nsvg, ivis2, neltst, ityptst, dt2t, gapv, inacti, index, niskyfi, kinet, newfront, isecin, nstrf, secfcum, x, xa, ce_loc, mfrot, ifq, frot_p, cand_fx, cand_fy, cand_fz, alpha0, ifpen, gapr, dxanc, nln, nlg, ibag, icontact, nsv, penis, penim, viscn, vxi, vyi, vzi, msi, kini, nin, nisub, lisub, addsubs, addsubm, lisubs, lisubm, fsavsub, cand_n, ilagm, icurv, nod_normal, fncont, ftcont, x1, x2, x3, x4, y1, y2, y3, y4, z1, z2, z3, z4, xi, yi, zi, iadm, rcurvi, rcontact, acontact, pcontact, anglmi, padm, intth, phi, fthe, ftheskyi, daanc6, temp, tempi, rstif, iform, gap_s, igap, alphak, mskyi_sms, iskyi_sms, nsms, cmaj, jtask, isensint, fsavparit, nft, h3d_data)
Definition i20for3.F:73
subroutine i20for3c(nln, nlg, ms, dxanc, dvanc, stfa, weight, inacti, daanc6, stfac, penia, alphak, daanc, kmin)
Definition i20for3.F:2269
subroutine i20ass0(jlt, cs_loc, n1, n2, m1, m2, hs1, hs2, hm1, hm2, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, a, stifn, stif, nrts, nin, jtask)
Definition i20for3.F:3355
subroutine i20ass25(jlt, cs_loc, n1, n2, m1, m2, hs1, hs2, hm1, hm2, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, isky, niskyfie, nrts, k1, k2, k3, k4, c1, c2, c3, c4, nin, noint)
Definition i20for3.F:3693
subroutine i20ass05(jlt, cs_loc, n1, n2, m1, m2, hs1, hs2, hm1, hm2, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, a, stifn, nrts, k1, k2, k3, k4, c1, c2, c3, c4, viscn, nin, jtask)
Definition i20for3.F:3448
subroutine i20ass2(jlt, cs_loc, n1, n2, m1, m2, hs1, hs2, hm1, hm2, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, fskyi, isky, niskyfie, stif, nrts, nin, noint)
Definition i20for3.F:3552
subroutine i20sms2e(jlt, cs_loc, n1, n2, m1, m2, hs1, hs2, hm1, hm2, stif, nin, noint, mskyi_sms, iskyi_sms, nsms, k1, k2, k3, k4, c1, c2, c3, c4, nrts)
Definition i20sms2.F:39
subroutine i7ass0(jlt, ix1, ix2, ix3, ix4, nsvg, h1, h2, h3, h4, stif, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, fxi, fyi, fzi, a, stifn, nin, intth, phi, fthe, phi1, phi2, phi3, phi4, condn, condint, jtask, iform, nodadt_therm)
Definition i7ass3.F:718
subroutine i7ass35(jlt, ix1, ix2, ix3, ix4, nsvg, h1, h2, h3, h4, stif, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, fxi, fyi, fzi, i8a, i8stifn, i8viscn, ks, k1, k2, k3, k4, cs, c1, c2, c3, c4)
Definition i7ass3.F:507
subroutine i7ass3(jlt, ix1, ix2, ix3, ix4, nsvg, h1, h2, h3, h4, stif, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, fxi, fyi, fzi, i8a, i8stifn)
Definition i7ass3.F:332
subroutine i7ass2(jlt, ix1, ix2, ix3, ix4, nsvg, h1, h2, h3, h4, stif, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, fxi, fyi, fzi, fskyi, isky, niskyfi, nin, noint, intth, phi, ftheskyi, phi1, phi2, phi3, phi4, condnskyi, condint, iform, nodadt_therm)
Definition i7ass3.F:1150
subroutine i7ass05(jlt, ix1, ix2, ix3, ix4, nsvg, h1, h2, h3, h4, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, fxi, fyi, fzi, a, stifn, viscn, ks, k1, k2, k3, k4, cs, c1, c2, c3, c4, nin, intth, phi, fthe, phi1, phi2, phi3, phi4, jtask, condn, condint, iform, nodadt_therm)
Definition i7ass3.F:936
subroutine i7ass25(jlt, ix1, ix2, ix3, ix4, nsvg, h1, h2, h3, h4, fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, fxi, fyi, fzi, fskyi, niskyfi, nin, ks, k1, k2, k3, k4, cs, c1, c2, c3, c4, isky, noint, intth, phi, ftheskyi, phi1, phi2, phi3, phi4, condnskyi, condint, iform, nodadt_therm)
Definition i7ass3.F:1478
subroutine i7curv(jlt, pene, n1, n2, n3, gapv, x, nod_normal, ix1, ix2, ix3, ix4, h1, h2, h3, h4, x1, x2, x3, x4, y1, y2, y3, y4, z1, z2, z3, z4, xi, yi, zi)
Definition i7curv.F:443
subroutine i7sms2(jlt, ix1, ix2, ix3, ix4, nsvg, h1, h2, h3, h4, stif, nin, noint, mskyi_sms, iskyi_sms, nsms, kt, c, cf, dtmini, dti)
Definition i7sms2.F:40
subroutine ibcoff(ibc, icodt)
Definition ibcoff.F:44
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
type(real_pointer2), dimension(:), allocatable dxancfi
Definition tri7box.F:473
type(real_pointer2), dimension(:), allocatable fskyfie
Definition tri7box.F:459
type(real_pointer2), dimension(:), allocatable fnconti
Definition tri7box.F:510
type(real_pointer), dimension(:), allocatable stifi
Definition tri7box.F:449
type(int_pointer), dimension(:), allocatable iskyfie
Definition tri7box.F:480
type(r8_pointer3), dimension(:), allocatable daanc6fi
Definition tri7box.F:476
type(real_pointer2), dimension(:), allocatable penfi
Definition tri7box.F:459
type(real_pointer), dimension(:), allocatable alphakfi
Definition tri7box.F:449
type(int_pointer), dimension(:), allocatable lisubsfi
Definition tri7box.F:501
type(real_pointer), dimension(:), allocatable stifie
Definition tri7box.F:449
type(real_pointer), dimension(:), allocatable stnfie
Definition tri7box.F:449
type(real_pointer), dimension(:), allocatable gapfi
Definition tri7box.F:449
type(int_pointer), dimension(:), allocatable nsvfi
Definition tri7box.F:431
type(real_pointer2), dimension(:), allocatable afie
Definition tri7box.F:459
type(int_pointer), dimension(:), allocatable nsvfie
Definition tri7box.F:440
type(int_pointer), dimension(:), allocatable addsubsfi
Definition tri7box.F:509
type(r8_pointer3), dimension(:), allocatable daanc6fie
Definition tri7box.F:476
type(real_pointer), dimension(:), allocatable vscfie
Definition tri7box.F:449
type(real_pointer), dimension(:), allocatable alphakfie
Definition tri7box.F:449
type(int_pointer), dimension(:), allocatable itafie
Definition tri7box.F:440
integer, dimension(:), allocatable nlskyfie
Definition tri7box.F:512
type(real_pointer2), dimension(:), allocatable penfie
Definition tri7box.F:459
integer, dimension(:), allocatable nlskyfi
Definition tri7box.F:512
type(real_pointer2), dimension(:), allocatable ftconti
Definition tri7box.F:510
type(int_pointer), dimension(:), allocatable itafi
Definition tri7box.F:440
subroutine foat_to_6_float(jft, jlt, f, f6)
Definition parit.F:225
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
subroutine arret(nn)
Definition arret.F:87