OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i18for3_mod Module Reference

Functions/Subroutines

subroutine i18for3 (jlt, a, v, ibcc, icodt, fsav, gap, fric, ms, visc, viscf, noint, stfn, itab, cn_loc, stfval, 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, ix1, ix2, ix3, ix4, nsvg, ivis2, neltst, ityptst, dt2t, ixs, gapv, cand_p, index, niskyfi, kinet, newfront, isecin, nstrf, secfcum, x, irect, ce_loc, mfrot, ifq, frot_p, cand_fx, cand_fy, cand_fz, ifpen, icontact, igroups, iparg, viscn, vxi, vyi, vzi, msi, kini, nin, nisub, lisub, addsubs, addsubm, lisubs, lisubm, fsavsub, cand_n, ilagm, icurv, fncont, ms0, jtask, isensint, fsavparit, nft, multi_fvm, h3d_data, elbuf_tab, idir)

Function/Subroutine Documentation

◆ i18for3()

subroutine i18for3_mod::i18for3 ( integer jlt,
a,
v,
integer ibcc,
integer, dimension(*) icodt,
fsav,
gap,
fric,
ms,
visc,
viscf,
integer noint,
stfn,
integer, dimension(*) itab,
integer, dimension(mvsiz) cn_loc,
stfval,
stifn,
stif,
fskyi,
integer, dimension(*) 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,
integer, dimension(mvsiz) ix1,
integer, dimension(mvsiz) ix2,
integer, dimension(mvsiz) ix3,
integer, dimension(mvsiz) ix4,
integer, dimension(mvsiz) nsvg,
integer ivis2,
integer neltst,
integer ityptst,
dt2t,
integer, dimension(nixs,numels) ixs,
gapv,
cand_p,
integer, dimension(mvsiz) index,
integer niskyfi,
integer, dimension(*) kinet,
integer newfront,
integer isecin,
integer, dimension(*) nstrf,
secfcum,
x,
integer, dimension(4,*) irect,
integer, dimension(mvsiz) ce_loc,
integer mfrot,
integer ifq,
frot_p,
cand_fx,
cand_fy,
cand_fz,
integer, dimension(*) ifpen,
integer, dimension(*) icontact,
integer, dimension(numels) igroups,
integer, dimension(nparg) iparg,
viscn,
vxi,
vyi,
vzi,
msi,
integer, dimension(*) kini,
integer nin,
integer nisub,
integer, dimension(*) lisub,
integer, dimension(*) addsubs,
integer, dimension(*) addsubm,
integer, dimension(*) lisubs,
integer, dimension(*) lisubm,
fsavsub,
integer, dimension(*) cand_n,
integer ilagm,
integer, dimension(3) icurv,
fncont,
ms0,
integer jtask,
integer, dimension(*) isensint,
fsavparit,
integer nft,
type(multi_fvm_struct), intent(inout) multi_fvm,
type(h3d_database) h3d_data,
type (elbuf_struct_), dimension(ngroup) elbuf_tab,
integer, intent(in) idir )

Definition at line 48 of file i18for3.F.

70C-----------------------------------------------
71C D e s c r i p t i o n
72C-----------------------------------------------
73C This subroutine is computing reaction forces
74C for fluid structure interaction /INTER/TYPE18
75C It also outputs contour and time histories.
76C
77C I=1:JLT : local loop (JLT<=128)
78C INDEX(I) is Node id (internal)
79C
80C-----------------------------------------------
81C M o d u l e s
82C-----------------------------------------------
83 USE tri7box
84 USE multi_fvm_mod
85 USE h3d_mod
86 USE aleanim_mod !FANI_CELL
87 USE elbufdef_mod
88 USE anim_mod
89C-----------------------------------------------
90C I m p l i c i t T y p e s
91C-----------------------------------------------
92#include "implicit_f.inc"
93#include "comlock.inc"
94C-----------------------------------------------
95C G l o b a l P a r a m e t e r s
96C-----------------------------------------------
97#include "mvsiz_p.inc"
98C-----------------------------------------------
99C C o m m o n B l o c k s
100C-----------------------------------------------
101#include "com01_c.inc"
102#include "com04_c.inc"
103#include "com06_c.inc"
104#include "com08_c.inc"
105#include "scr07_c.inc"
106#include "scr14_c.inc"
107#include "scr16_c.inc"
108#include "scr18_c.inc"
109#include "units_c.inc"
110#include "parit_c.inc"
111#include "param_c.inc"
112#include "kincod_c.inc"
113C-----------------------------------------------
114C D u m m y A r g u m e n t s
115C-----------------------------------------------
116 INTEGER, INTENT(IN) :: IDIR
117 TYPE (ELBUF_STRUCT_), DIMENSION(NGROUP) :: ELBUF_TAB
118 INTEGER NELTST,ITYPTST,JLT,IBCC,IVIS2,NIN,
119 . ICODT(*), ITAB(*), ISKY(*), KINET(*),
120 . MFROT, IFQ, NOINT,NEWFRONT,ISECIN, NSTRF(*),
121 . IRECT(4,*),IFPEN(*) ,ICONTACT(*), CAND_N(*),
122 . KINI(*),IGROUPS(NUMELS),
123 . ISET, NISKYFI,INTTH,IFORM,JTASK,NFT,IPARG(NPARG)
124 INTEGER IX1(MVSIZ), IX2(MVSIZ), IX3(MVSIZ), IX4(MVSIZ),
125 . CN_LOC(MVSIZ),CE_LOC(MVSIZ),INDEX(MVSIZ),NSVG(MVSIZ),
126 . NISUB, LISUB(*), ADDSUBS(*), ADDSUBM(*), LISUBS(*),
127 . LISUBM(*),ILAGM,ICURV(3),ISENSINT(*),IXS(NIXS,NUMELS)
128 my_real
129 . stfval,cand_p(*),frot_p(*), x(3,*),ms0(*),
130 . a(3,*), ms(*), v(3,*), fsav(*),fcont(3,*),
131 . cand_fx(*),cand_fy(*),cand_fz(*),
132 . gap, fric,visc,viscf,vis,dt2t,stfn(*),stifn(*),
133 . fskyi(lskyi,nfskyi),fsavsub(nthvki,*),fncont(3,*),
134 . fsavparit(nisub+1,11,*)
135 my_real
136 . nx1(mvsiz), nx2(mvsiz), nx3(mvsiz), nx4(mvsiz),
137 . ny1(mvsiz), ny2(mvsiz), ny3(mvsiz), ny4(mvsiz),
138 . nz1(mvsiz), nz2(mvsiz), nz3(mvsiz), nz4(mvsiz),
139 . lb1(mvsiz), lb2(mvsiz), lb3(mvsiz), lb4(mvsiz),
140 . lc1(mvsiz), lc2(mvsiz), lc3(mvsiz), lc4(mvsiz),
141 . p1(mvsiz), p2(mvsiz), p3(mvsiz), p4(mvsiz), stif(mvsiz),
142 . gapv(mvsiz),
143 . secfcum(7,numnod,nsect), tmp(mvsiz),
144 . stifsav(mvsiz), viscn(*),
145 . vxi(mvsiz),vyi(mvsiz),vzi(mvsiz),msi(2*mvsiz),
146 . rstif
147 TYPE(MULTI_FVM_STRUCT), INTENT(INOUT) :: MULTI_FVM
148 TYPE(H3D_DATABASE) :: H3D_DATA
149C-----------------------------------------------
150C L o c a l V a r i a b l e s
151C-----------------------------------------------
152 INTEGER I, J1, IG, J, JG , K0,NBINTER,K1S,K,IL,IE, NN, NI,
153 . NA1,NA2
154 my_real
155 . fxi(mvsiz), fyi(mvsiz), fzi(mvsiz), fni(mvsiz),
156 . fxt(mvsiz),fyt(mvsiz),fzt(mvsiz),
157 . fx1(mvsiz), fx2(mvsiz), fx3(mvsiz), fx4(mvsiz),
158 . fy1(mvsiz), fy2(mvsiz), fy3(mvsiz), fy4(mvsiz),
159 . fz1(mvsiz), fz2(mvsiz), fz3(mvsiz), fz4(mvsiz),
160 . n1(mvsiz), n2(mvsiz), n3(mvsiz), pene(mvsiz),
161 . vis2(mvsiz), dtmi(mvsiz), xmu(mvsiz),
162 . h1(mvsiz), h2(mvsiz), h3(mvsiz), h4(mvsiz),
163 . vx(mvsiz), vy(mvsiz), vz(mvsiz), vn(mvsiz),dist(mvsiz),
164 . vnx, vny, vnz, aa, crit,s2,rdist,
165 . v2, fm2, dt1inv, visca, fac,ff,alphi,alpha,beta,
166 . fx, fy, fz, f2, mas2, m2sk, dtmi0,ft,fn,fmax,ftn,
167 . facm1, econtt, econvt, h0, la1, la2, la3, la4,
168 . d1,d2,d3,d4,a1,a2,a3,a4,econtdt,
169 . fsav1, fsav2, fsav3, fsav4, fsav5, fsav6, fsav7, fsav8,
170 . fsav9, fsav10, fsav11, fsav12, fsav13, fsav14, fsav15, ffo,
171 . e10, h0d, s2d, sum,
172 . la1d,la2d,la3d,la4d,t1,t1d,t2,t2d,ffd,visd,facd,d1d,
173 . p1s(mvsiz),p2s(mvsiz),p3s(mvsiz),p4s(mvsiz),
174 . d2d,d3d,d4d,vnxd,vnyd,vnzd,v2d,fm2d,f2d,aad,fxd,fyd,fzd,
175 . a1d,a2d,a3d,a4d,vv,ax1,ax2,ay1,ay2,az1,az2,ax,ay,az,
176 . area,p,vv1,vv2,v21,dmu, h00 ,a0x,a0y,a0z,rx,ry,rz,
177 . anx,any,anz,aan,aax,aay,aaz ,rr,rs,aaa ,tm,ts
178 my_real
179 . surfx,surfy,surfz,surf
180 my_real
181 . st1(mvsiz),st2(mvsiz),st3(mvsiz),st4(mvsiz),stv(mvsiz),
182 . kt(mvsiz),c(mvsiz),cf(mvsiz),
183 . ks(mvsiz),k1(mvsiz),k2(mvsiz),k3(mvsiz),k4(mvsiz),
184 . cs(mvsiz),c1(mvsiz),c2(mvsiz),c3(mvsiz),c4(mvsiz),
185 . cx,cy,cfi,aux,phi1(mvsiz), phi2(mvsiz), phi3(mvsiz),
186 . phi4(mvsiz),dx, dti
187 INTEGER JSUB,KSUB,JJ,KK,IN,NSUB,IBID,ITASK,NELFT,NELLT
188 my_real fsavsub1(15,nisub),impx,impy,impz,pp1,pp2,pp3,pp4,bid
189C-----------------------------------------------
190C S o u r c e L i n e s
191C-----------------------------------------------
192
193 ibid = 0
194 bid = zero
195 econtt = zero
196 econtdt = zero
197
198 dt1inv =zero
199 IF(dt1 > zero)dt1inv = one/dt1
200
201 itask = jtask-1
202
203 !---------------------
204 ! PENETRATION CALCULATION (PENE)
205 ! & LOCAL COORDINATES (H1,H2,H3,H4) ON FACE FOR VELOCITY INTERPOLATION ON FACE
206 !---------------------
207 DO i=1,jlt
208 IF(ix3(i) /= ix4(i))THEN
209 !distance from triangles
210 d1 = sqrt(p1(i))
211 d2 = sqrt(p2(i))
212 d3 = sqrt(p3(i))
213 d4 = sqrt(p4(i))
214 !penetration into each triangle gap
215 pp1 = max(zero, gap - d1)
216 pp2 = max(zero, gap - d2)
217 pp3 = max(zero, gap - d3)
218 pp4 = max(zero, gap - d4)
219 !MAIN face penetration
220 pene(i) = max(pp1,pp2,pp3,pp4)
221 !ratios
222 a1 = pp1/max(em20,d1)
223 a2 = pp2/max(em20,d2)
224 a3 = pp3/max(em20,d3)
225 a4 = pp4/max(em20,d4)
226 !sum(ratio[k]*2Sn[k],k=1..4)
227 n1(i) = a1*nx1(i) + a2*nx2(i) + a3*nx3(i) + a4*nx4(i)
228 n2(i) = a1*ny1(i) + a2*ny2(i) + a3*ny3(i) + a4*ny4(i)
229 n3(i) = a1*nz1(i) + a2*nz2(i) + a3*nz3(i) + a4*nz4(i)
230 la1 = one - lb1(i) - lc1(i)
231 la2 = one - lb2(i) - lc2(i)
232 la3 = one - lb3(i) - lc3(i)
233 la4 = one - lb4(i) - lc4(i)
234 h0 = fourth * (pp1*la1 + pp2*la2 + pp3*la3 + pp4*la4)
235 h1(i) = h0 + pp1 * lb1(i) + pp4 * lc4(i)
236 h2(i) = h0 + pp2 * lb2(i) + pp1 * lc1(i)
237 h3(i) = h0 + pp3 * lb3(i) + pp2 * lc2(i)
238 h4(i) = h0 + pp4 * lb4(i) + pp3 * lc3(i)
239 h00 = one/max(em20,h1(i) + h2(i) + h3(i) + h4(i))
240 h1(i) = h1(i) * h00
241 h2(i) = h2(i) * h00
242 h3(i) = h3(i) * h00
243 h4(i) = h4(i) * h00
244 ELSE
245 d1 = sqrt(p1(i))
246 pp1 = max(zero, gap - d1)
247 pene(i) = pp1
248 n1(i) = nx1(i)
249 n2(i) = ny1(i)
250 n3(i) = nz1(i)
251 h1(i) = lb1(i)
252 h2(i) = lc1(i)
253 h3(i) = one - lb1(i) - lc1(i)
254 h4(i) = zero
255 ENDIF
256 ENDDO
257 !---------------------
258 ! UNITARY NORMAL VECTOR
259 !---------------------
260 DO i=1,jlt
261 s2 = one/max(em30,sqrt(n1(i)*n1(i) + n2(i)*n2(i) + n3(i)*n3(i)))
262 n1(i) = n1(i)*s2
263 n2(i) = n2(i)*s2
264 n3(i) = n3(i)*s2
265 ENDDO
266 !---------------------
267 ! RELATIVE VELOCITY (VX,VY,VZ)
268 ! & ITS NORMAL COMPONENT (VN)
269 !---------------------
270 DO i=1,jlt
271 vx(i) = vxi(i) - h1(i)*v(1,ix1(i)) - h2(i)*v(1,ix2(i)) - h3(i)*v(1,ix3(i)) - h4(i)*v(1,ix4(i))
272 vy(i) = vyi(i) - h1(i)*v(2,ix1(i)) - h2(i)*v(2,ix2(i)) - h3(i)*v(2,ix3(i)) - h4(i)*v(2,ix4(i))
273 vz(i) = vzi(i) - h1(i)*v(3,ix1(i)) - h2(i)*v(3,ix2(i)) - h3(i)*v(3,ix3(i)) - h4(i)*v(3,ix4(i))
274 vn(i) = n1(i)*vx(i) + n2(i)*vy(i) + n3(i)*vz(i)
275 !hourglass correction
276 h0 = -fourth*(h1(i) - h2(i) + h3(i) - h4(i))
277 h0 = min(h0,h2(i),h4(i))
278 h0 = max(h0,-h1(i),-h3(i))
279 IF(ix3(i) == ix4(i))h0 = zero
280 h1(i) = h1(i) + h0
281 h2(i) = h2(i) - h0
282 h3(i) = h3(i) + h0
283 h4(i) = h4(i) - h0
284 ENDDO
285
286 !----------------------------------------------!
287 ! CONSTANT STIFFNESS !
288 !----------------------------------------------!
289 DO i=1,jlt
290 IF(pene(i) > zero )THEN
291 dist(i)=gap-pene(i)
292 cand_p(index(i)) = cand_p(index(i)) + vn(i)*dt1
293 ! IF(CAND_P(INDEX(I)) > ZERO)CAND_P(INDEX(I))=ZERO !contact lost
294 stif(i) = stfval * pene(i)/gap !STVAL:constant user param
295 ELSE
296 !remove old candidates
297 cand_p(index(i)) = zero
298 stif(i) = zero
299 ENDIF
300 ENDDO
301
302 !----------------------------------------------!
303 ! REACTION FORCE !
304 !----------------------------------------------!
305 DO i=1,jlt
306 IF(pene(i) > zero )THEN
307 fni(i) = stif(i) * cand_p(index(i))
308 ELSE
309 fni(i)=zero
310 ENDIF
311 ENDDO
312
313 !!KINEMATIC TIME STEP
314 !DTI = EP20
315 !IF(INTER18_AUTOPARAM == 1)THEN
316 ! DO I=1,JLT
317 ! DX =GAP/TEN
318 ! RDIST = DX / MAX(EM20,ABS(VN(I)))
319 ! DTI = MIN(RDIST,DTI)
320 ! ENDDO
321 !ENDIF
322 !IF(DTI<DT2T)THEN
323 ! DT2T = DTI
324 ! NELTST = NOINT
325 ! ITYPTST = 10
326 !ENDIF
327
328 !---------------------------------
329 ! EXPERIMENTAL SURFACE ORIENTATION
330 !---------------------------------
331 IF(idir == -1)THEN
332 DO i=1,jlt
333 fni(i) = min(fni(i),zero)
334 ENDDO
335 ELSEIF(idir == 1)THEN
336 DO i=1,jlt
337 fni(i) = max(fni(i),zero)
338 ENDDO
339 ENDIF
340
341 !---------------------------------
342 ! DAMPING
343 !---------------------------------
344 IF(visc /= zero)THEN
345 DO i=1,jlt
346 IF(vn(i) > zero)THEN
347 fac = stif(i) / max(em30,stif(i))
348 ff = fac * visc * pene(i)/gap
349 stif(i) = stif(i) + two * ff * dt1inv
350 ff = ff * vn(i)
351 econtdt = econtdt + ff * vn(i) * dt1 ! Damping Energy
352 fni(i) = fni(i) + ff
353 ENDIF
354 ENDDO
355 ENDIF
356
357 !---------------------------------
358 ! NORMAL IMPULSE OUTPUT
359 !---------------------------------
360 fsav1 = zero
361 fsav2 = zero
362 fsav3 = zero
363 fsav8 = zero
364 fsav9 = zero
365 fsav10= zero
366 fsav11= zero
367 DO i=1,jlt
368 econtt = econtt + dt1*vn(i)*fni(i) ! Elastic energy : it is cumulated energy per cycle
369 fxi(i)=n1(i)*fni(i)
370 fyi(i)=n2(i)*fni(i)
371 fzi(i)=n3(i)*fni(i)
372 impx=fxi(i)*dt12
373 impy=fyi(i)*dt12
374 impz=fzi(i)*dt12
375 fsav1 =fsav1 +impx
376 fsav2 =fsav2 +impy
377 fsav3 =fsav3 +impz
378 fsav8 =fsav8 +abs(impx)
379 fsav9 =fsav9 +abs(impy)
380 fsav10=fsav10+abs(impz)
381 fsav11=fsav11+fni(i)*dt12
382 ENDDO
383#include "lockon.inc"
384 fsav(1)=fsav(1)+fsav1
385 fsav(2)=fsav(2)+fsav2
386 fsav(3)=fsav(3)+fsav3
387 fsav(8)=fsav(8)+fsav8
388 fsav(9)=fsav(9)+fsav9
389 fsav(10)=fsav(10)+fsav10
390 fsav(11)=fsav(11)+fsav11
391#include "lockoff.inc"
392 IF(isensint(1)/=0) THEN
393 DO i=1,jlt
394 fsavparit(1,1,i+nft) = fsavparit(1,1,i+nft) + sqrt((fxi(i)**2)+(fyi(i)**2)+(fzi(i)**2))
395 ENDDO
396 ENDIF
397 !---------------------------------
398 ! TH OUTPUT FOR SUB INTERFACES
399 !---------------------------------
400 IF(nisub /= 0)THEN
401 DO jsub=1,nisub
402 DO j=1,15
403 fsavsub1(j,jsub)=zero
404 END DO
405 END DO
406 DO i=1,jlt
407 nn = nsvg(i)
408 IF(nn > 0)THEN
409 in=cn_loc(i)
410 ie=ce_loc(i)
411 jj =addsubs(in)
412 kk =addsubm(ie)
413 DO WHILE(jj<addsubs(in+1))
414 jsub=lisubs(jj)
415 DO WHILE(kk<addsubm(ie+1))
416 ksub=lisubm(kk)
417 IF(ksub == jsub)THEN
418 impx=fxi(i)*dt12
419 impy=fyi(i)*dt12
420 impz=fzi(i)*dt12
421 !MAIN side
422 fsavsub1(1,jsub)=fsavsub1(1,jsub)+impx
423 fsavsub1(2,jsub)=fsavsub1(2,jsub)+impy
424 fsavsub1(3,jsub)=fsavsub1(3,jsub)+impz
425 IF(isensint(jsub)/=0) THEN
426 fsavparit(jsub+1,1,i+nft) = fsavparit(jsub+1,1,i+nft) + sqrt((fxi(i)**2)+(fyi(i)**2)+(fzi(i)**2))
427 ENDIF
428 fsavsub1(8,jsub) =fsavsub1(8,jsub) +abs(impx)
429 fsavsub1(9,jsub) =fsavsub1(9,jsub) +abs(impy)
430 fsavsub1(10,jsub)=fsavsub1(10,jsub)+abs(impz)
431 fsavsub1(11,jsub)=fsavsub1(11,jsub)+fni(i)*dt12
432 kk=kk+1
433 EXIT
434 ELSE IF(ksub<jsub)THEN
435 kk=kk+1
436 ELSE
437 EXIT
438 END IF
439 END DO
440 jj=jj+1
441 END DO !WHILE
442 ELSE
443 nn = -nn
444 ie=ce_loc(i)
445 jj =addsubsfi(nin)%P(nn)
446 kk =addsubm(ie)
447 DO WHILE(jj<addsubsfi(nin)%P(nn+1))
448 jsub=lisubsfi(nin)%P(jj)
449 DO WHILE(kk<addsubm(ie+1))
450 ksub=lisubm(kk)
451 IF(ksub == jsub)THEN
452 impx=fxi(i)*dt12
453 impy=fyi(i)*dt12
454 impz=fzi(i)*dt12
455 ! MAIN side
456 fsavsub1(1,jsub)=fsavsub1(1,jsub)+impx
457 fsavsub1(2,jsub)=fsavsub1(2,jsub)+impy
458 fsavsub1(3,jsub)=fsavsub1(3,jsub)+impz
459 fsavsub1(8,jsub) =fsavsub1(8,jsub) +abs(impx)
460 fsavsub1(9,jsub) =fsavsub1(9,jsub) +abs(impy)
461 fsavsub1(10,jsub)=fsavsub1(10,jsub)+abs(impz)
462 fsavsub1(11,jsub)=fsavsub1(11,jsub)+fni(i)*dt12
463 kk=kk+1
464 EXIT
465 ELSE IF(ksub<jsub)THEN
466 kk=kk+1
467 ELSE
468 EXIT
469 END IF
470 END do! WHILE (KK<ADDSUBM(IE+1))
471 jj=jj+1
472 END do! WHILE (JJ<ADDSUBSFI(NIN)%P(NN+1))
473 END IF
474 END DO !NEXT I=1,JLT
475 END IF
476 !---------------------------------
477 ! SORTIES TH PAR SOUS INTERFACE
478 !---------------------------------
479 IF(nisub /= 0)THEN
480 DO i=1,jlt
481 nn = nsvg(i)
482 IF(nn > 0)THEN
483 in=cn_loc(i)
484 ie=ce_loc(i)
485 jj =addsubs(in)
486 kk =addsubm(ie)
487 DO WHILE(jj<addsubs(in+1))
488 jsub=lisubs(jj)
489 DO WHILE(kk<addsubm(ie+1))
490 ksub=lisubm(kk)
491 IF(ksub == jsub)THEN
492 !CURRENTLY NO TANGENTIAL FORCE IN INTER18
493 !IMPX=FXT(I)*DT12
494 !IMPY=FYT(I)*DT12
495 !IMPZ=FZT(I)*DT12
496 !FSAVSUB1(4,JSUB)=FSAVSUB1(4,JSUB)+IMPX
497 !FSAVSUB1(5,JSUB)=FSAVSUB1(5,JSUB)+IMPY
498 !FSAVSUB1(6,JSUB)=FSAVSUB1(6,JSUB)+IMPZ
499 impx=fxi(i)*dt12
500 impy=fyi(i)*dt12
501 impz=fzi(i)*dt12
502 fsavsub1(12,jsub)=fsavsub1(12,jsub)+abs(impx)
503 fsavsub1(13,jsub)=fsavsub1(13,jsub)+abs(impy)
504 fsavsub1(14,jsub)=fsavsub1(14,jsub)+abs(impz)
505 fsavsub1(15,jsub)= fsavsub1(15,jsub) + sqrt(impx*impx+impy*impy+impz*impz)
506 kk=kk+1
507 EXIT
508 ELSE IF(ksub<jsub)THEN
509 kk=kk+1
510 ELSE
511 EXIT
512 END IF
513 END do!WHILE (KK<ADDSUBM(IE+1))
514 jj=jj+1
515 END do! WHILE (JJ<ADDSUBS(IN+1))
516 ELSE
517 nn = -nn
518 ie=ce_loc(i)
519 jj =addsubsfi(nin)%P(nn)
520 kk =addsubm(ie)
521 DO WHILE(jj<addsubsfi(nin)%P(nn+1))
522 jsub=lisubsfi(nin)%P(jj)
523 DO WHILE(kk<addsubm(ie+1))
524 ksub=lisubm(kk)
525 IF(ksub == jsub)THEN
526 !CURRENTLY NO TANGENTIAL FORCE IN INTER18
527 !IMPX=FXT(I)*DT12
528 !IMPY=FYT(I)*DT12
529 !IMPZ=FZT(I)*DT12
530 !!MAIN side :
531 !FSAVSUB1(4,JSUB)=FSAVSUB1(4,JSUB)+IMPX
532 !FSAVSUB1(5,JSUB)=FSAVSUB1(5,JSUB)+IMPY
533 !FSAVSUB1(6,JSUB)=FSAVSUB1(6,JSUB)+IMPZ
534 impx=fxi(i)*dt12
535 impy=fyi(i)*dt12
536 impz=fzi(i)*dt12
537 fsavsub1(12,jsub)=fsavsub1(12,jsub)+abs(impx)
538 fsavsub1(13,jsub)=fsavsub1(13,jsub)+abs(impy)
539 fsavsub1(14,jsub)=fsavsub1(14,jsub)+abs(impz)
540 fsavsub1(15,jsub)= fsavsub1(15,jsub)+sqrt(impx*impx+impy*impy+impz*impz)
541 kk=kk+1
542 EXIT
543 ELSE IF(ksub<jsub)THEN
544 kk=kk+1
545 ELSE
546 EXIT
547 END IF
548 END do!DO WHILE(KK<ADDSUBM(IE+1))
549 jj=jj+1
550 END do!WHILE(JJ<ADDSUBSFI(NIN)%P(NN+1))
551 END IF
552 END DO !NEXT I=1,JLT
553#include "lockon.inc"
554 DO jsub=1,nisub
555 nsub=lisub(jsub)
556 DO j=1,15
557 fsavsub(j,nsub)=fsavsub(j,nsub)+fsavsub1(j,jsub)
558 END DO
559 END DO
560#include "lockoff.inc"
561 END IF
562C---------------------------------
563#include "lockon.inc"
564 econtd = econtd + econtdt ! Damping Energy
565 econt_cumu = econt_cumu + econtt ! ELASTIC ENERGY
566 fsav(26) = fsav(26) + econtt
567 fsav(28) = fsav(28) + econtdt
568#include "lockoff.inc"
569C---------------------------------
570
571 IF(idtmin(10) == 1.OR.idtmin(10) == 2.OR.idtmin(10) == 5.OR.idtmin(10) == 6)THEN
572 dtmi0 = ep20
573 DO i=1,jlt
574 dtmi(i) = ep20
575 mas2 = two * msi(i)
576 IF(mas2>zero.AND.stif(i)>zero .AND. irb(kini(i))==0.AND.irb2(kini(i))==0)THEN
577 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/stif(i)))
578 ENDIF
579 mas2 = two* ms(ix1(i))
580 IF(mas2>zero.AND.h1(i)*stif(i)>zero .AND. irb(kinet(ix1(i)))==0.AND.irb2(kinet(ix1(i)))==0)THEN
581 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(h1(i)*stif(i))))
582 ENDIF
583 mas2 = two * ms(ix2(i))
584 IF(mas2>zero.AND.h2(i)*stif(i)>zero .AND. irb(kinet(ix2(i)))==0.AND.irb2(kinet(ix2(i)))==0)THEN
585 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(h2(i)*stif(i))))
586 ENDIF
587 mas2 = two* ms(ix3(i))
588 IF(mas2 > zero.AND.h3(i)*stif(i) > zero .AND. irb(kinet(ix3(i))) == 0.AND.irb2(kinet(ix3(i))) == 0)THEN
589 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(h3(i)*stif(i))))
590 ENDIF
591 mas2 = two * ms(ix4(i))
592 IF(mas2 > zero.AND.h4(i)*stif(i) > zero .AND. irb(kinet(ix4(i))) == 0.AND.irb2(kinet(ix4(i))) == 0)THEN
593 dtmi(i) = min(dtmi(i),dtfac1(10)*sqrt(mas2/(h4(i)*stif(i))))
594 ENDIF
595 dtmi0 = min(dtmi0,dtmi(i))
596 ENDDO
597 IF(dtmi0<=dtmin1(10))THEN
598 DO i=1,jlt
599 IF(dtmi(i)<=dtmin1(10))THEN
600 jg = nsvg(i)
601 IF(jg > 0)THEN
602 ni = itab(jg)
603 ELSE
604 ni = itafi(nin)%P(-jg)
605 ENDIF
606 IF(idtmin(10) == 1)THEN
607#include "lockon.inc"
608 WRITE(iout,'(A,E12.4,A,I10,A,E12.4,A)')
609 . ' **WARNING MINIMUM TIME STEP ',dtmi(i),
610 . ' IN INTERFACE ',noint
611 WRITE(iout,'(A,I10)') ' SECONDARY NODE : ',ni
612 WRITE(iout,'(A,4I10)')' MAIN NODES : ',
613 . itab(ix1(i)),itab(ix2(i)),itab(ix3(i)),itab(ix4(i))
614#include "lockoff.inc"
615 tstop = tt
616 IF ( istamping == 1) THEN
617 WRITE(istdo,'(A)')'The run encountered a problem in an interface Type 7.'
618 WRITE(istdo,'(A)')'You may need to check if there is enou gh clearance between the tools,'
619 WRITE(istdo,'(A)')'and that they do not penetrate each other during their travel'
620 WRITE(iout, '(A)')'The run encountered a problem in an interface Type 7.'
621 WRITE(iout, '(A)')'You may need to check if there is enough clearance between the tools,'
622 WRITE(iout, '(A)')'and that they do not penetrate each other during their travel'
623 ENDIF
624 ELSEIF(idtmin(10) == 2)THEN
625#include "lockon.inc"
626 WRITE(iout,'(A,E12.4,A,I10,A,E12.4,A)')' **WARNING MINIMUM TIME STEP ',dtmi(i),' IN INTERFACE ',noint
627 WRITE(iout,'(A,I10,A,I10)')' DELETE SECONDARY NODE ',ni,' FROM INTERFACE ',noint
628 WRITE(iout,'(A,4I10)')' MAIN NODES : ',itab(ix1(i)),itab(ix2(i)),itab(ix3(i)),itab(ix4(i))
629 IF(jg > 0) THEN
630 stfn(cn_loc(i)) = -abs(stfn(cn_loc(i)))
631 ELSE
632 stifi(nin)%P(-jg) = -abs(stifi(nin)%P(-jg))
633 ENDIF
634#include "lockoff.inc"
635 IF ( istamping == 1) THEN
636 WRITE(istdo,'(A)')'The run encountered a problem in an interface Type 7.'
637 WRITE(istdo,'(A)')'You may need to check if there is enou gh clearance between the tools,'
638 WRITE(istdo,'(A)')'and that they do not penetrate each other during their travel'
639 WRITE(iout, '(A)')'The run encountered a problem in an interface Type 7.'
640 WRITE(iout, '(A)')'You may need to check if there is enough clearance between the tools,'
641 WRITE(iout, '(A)')'and that they do not penetrate each other during their travel'
642 ENDIF
643 newfront = -1
644 ELSEIF(idtmin(10) == 5)THEN
645#include "lockon.inc"
646 WRITE(iout,'(A,E12.4,A,I10,A,E12.4,A)')' **WARNING MINIMUM TIME STEP ',dtmi(i),' IN INTERFACE ',noint
647 WRITE(iout,'(A,I10)') ' SECONDARY NODE : ',ni
648 WRITE(iout,'(A,4I10)')' MAIN NODES : ',
649 . itab(ix1(i)),itab(ix2(i)),itab(ix3(i)),itab(ix4(i))
650#include "lockoff.inc"
651 mstop = 2
652 IF ( istamping == 1) THEN
653 WRITE(istdo,'(A)')'The run encountered a problem in an interface Type 7.'
654 WRITE(istdo,'(A)')'You may need to check if there is enou gh clearance between the tools,'
655 WRITE(istdo,'(A)')'and that they do not penetrate each other during their travel'
656 WRITE(iout, '(A)')'The run encountered a problem in an interface Type 7.'
657 WRITE(iout, '(A)')'You may need to check if there is enough clearance between the tools,'
658 WRITE(iout, '(A)')'and that they do not penetrate each other during their travel'
659 ENDIF
660 ELSEIF(idtmin(10) == 6.AND.ilagm == 2)THEN
661 IF(kinet(jg)+kinet(ix1(i))+kinet(ix2(i))+kinet(ix3(i))+kinet(ix4(i)) == 0 )THEN
662 cand_n(index(i)) = -iabs(cand_n(index(i)))
663 stif(i) = zero
664 fxi(i) = zero
665 fyi(i) = zero
666 fzi(i) = zero
667 ENDIF
668 ENDIF
669 ENDIF
670 ENDDO
671 ENDIF
672 ENDIF
673C---------------------------------
674 DO i=1,jlt
675 fx1(i)=fxi(i)*h1(i)
676 fy1(i)=fyi(i)*h1(i)
677 fz1(i)=fzi(i)*h1(i)
678 !
679 fx2(i)=fxi(i)*h2(i)
680 fy2(i)=fyi(i)*h2(i)
681 fz2(i)=fzi(i)*h2(i)
682 !
683 fx3(i)=fxi(i)*h3(i)
684 fy3(i)=fyi(i)*h3(i)
685 fz3(i)=fzi(i)*h3(i)
686 !
687 fx4(i)=fxi(i)*h4(i)
688 fy4(i)=fyi(i)*h4(i)
689 fz4(i)=fzi(i)*h4(i)
690 !case law151 with single lagrangian fluid brick and constant velocity
691 !FXI(I)=ZERO
692 !FYI(I)=ZERO
693 !FZI(I)=ZERO
694 ENDDO
695
696 !SPMD : identification des noeuds interf. utiles a envoyer
697 IF (nspmd > 1) THEN
698 DO i = 1,jlt
699 nn = nsvg(i)
700 IF(nn<0)THEN
701 ! tag temporaire de NSVFI a -
702 nsvfi(nin)%P(-nn) = -abs(nsvfi(nin)%P(-nn))
703 ENDIF
704 ENDDO
705 ENDIF
706 IF (multi_fvm%IS_USED) THEN
707 IF (iparit == 0) THEN
709 1 dt1 ,jlt ,ix1 ,ix2 ,ix3 ,ix4 ,
710 2 nsvg ,h1 ,h2 ,h3 ,h4 ,stif ,
711 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
712 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
713 5 fxi ,fyi ,fzi ,a ,stifn ,nin ,
714 7 jtask ,multi_fvm ,x ,ixs ,v ,
715 8 elbuf_tab,igroups ,iparg,msi)
716 ELSE
718 1 jlt ,ix1 ,ix2 ,ix3 ,ix4 ,
719 2 nsvg ,h1 ,h2 ,h3 ,h4 ,stif ,
720 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
721 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
722 5 fxi ,fyi ,fzi ,fskyi,isky ,niskyfi,
723 6 nin ,noint ,multi_fvm,dt1,jtask)
724 ENDIF
725 ELSE
726 IF(iparit == 3)THEN
727 CALL i7ass3( jlt ,ix1 ,ix2 ,ix3 ,ix4 ,
728 2 nsvg ,h1 ,h2 ,h3 ,h4 ,stif ,
729 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
730 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
731 5 fxi ,fyi ,fzi ,a ,stifn)
732
733 ELSEIF(iparit == 0)THEN
734 CALL i7ass0(jlt ,ix1 ,ix2 ,ix3 ,ix4 ,
735 2 nsvg ,h1 ,h2 ,h3 ,h4 ,stif ,
736 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
737 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
738 5 fxi ,fyi ,fzi ,a ,stifn ,nin ,
739 6 ibid ,bid ,bid ,bid ,bid ,bid ,
740 7 bid ,bid ,bid ,jtask,ibid ,ibid )
741 ELSE
742 CALL i7ass2(jlt ,ix1 ,ix2 ,ix3 ,ix4 ,
743 2 nsvg ,h1 ,h2 ,h3 ,h4 ,stif ,
744 3 fx1 ,fy1 ,fz1 ,fx2 ,fy2 ,fz2 ,
745 4 fx3 ,fy3 ,fz3 ,fx4 ,fy4 ,fz4 ,
746 5 fxi ,fyi ,fzi ,fskyi,isky ,niskyfi,
747 6 nin ,noint ,ibid ,bid ,bid ,bid ,
748 7 bid ,bid ,bid ,bid ,bid ,
749 a ibid ,ibid )
750 ENDIF
751 ENDIF
752 IF(anim_v(4)+outp_v(4)+h3d_data%N_VECT_CONT > 0)THEN
753#include "lockon.inc"
754 !---REACTION FORCE ON THE SOLID SURFACE
755 ! /H3D/NODA/CONT (STAGGERED SCHEME)
756 DO i=1,jlt
757 fcont(1,ix1(i)) =fcont(1,ix1(i)) + fx1(i)
758 fcont(2,ix1(i)) =fcont(2,ix1(i)) + fy1(i)
759 fcont(3,ix1(i)) =fcont(3,ix1(i)) + fz1(i)
760 fcont(1,ix2(i)) =fcont(1,ix2(i)) + fx2(i)
761 fcont(2,ix2(i)) =fcont(2,ix2(i)) + fy2(i)
762 fcont(3,ix2(i)) =fcont(3,ix2(i)) + fz2(i)
763 fcont(1,ix3(i)) =fcont(1,ix3(i)) + fx3(i)
764 fcont(2,ix3(i)) =fcont(2,ix3(i)) + fy3(i)
765 fcont(3,ix3(i)) =fcont(3,ix3(i)) + fz3(i)
766 fcont(1,ix4(i)) =fcont(1,ix4(i)) + fx4(i)
767 fcont(2,ix4(i)) =fcont(2,ix4(i)) + fy4(i)
768 fcont(3,ix4(i)) =fcont(3,ix4(i)) + fz4(i)
769 ENDDO
770 !---REACTION FORCE ON THE FLUID SIDE
771 ! /H3D/ELEM/VECT/CONT WITH LAW151 (COLOCATED SCHEME)
772 IF (multi_fvm%IS_USED) THEN
773 IF(h3d_data%N_VECT_CONT > 0)THEN
774 DO i=1,jlt
775 jg=nsvg(i)-numnod
776 IF(jg>0) THEN ! local cell
777 fani_cell%F18(1,jg)=fani_cell%F18(1,jg)-fxi(i)
778 fani_cell%F18(2,jg)=fani_cell%F18(2,jg)-fyi(i)
779 fani_cell%F18(3,jg)=fani_cell%F18(3,jg)-fzi(i)
780 ENDIF
781 ENDDO
782 ENDIF
783 ELSE
784 !/H3D/NODA/CONT (STAGGERED SCHEME)
785 DO i=1,jlt
786 jg = nsvg(i)
787 IF(jg > 0) THEN
788 ! SPMD : reprocess required after receiving a remote node if JG < 0
789 fcont(1,jg)=fcont(1,jg)- fxi(i)
790 fcont(2,jg)=fcont(2,jg)- fyi(i)
791 fcont(3,jg)=fcont(3,jg)- fzi(i)
792 ENDIF
793 ENDDO
794 ENDIF
795#include "lockoff.inc"
796 ENDIF
797 IF((anim_v(12)+outp_v(12)+h3d_data%N_VECT_PCONT > 0 .AND.
798 . ((tt>=tanim .AND. tt<=tanim_stop).OR.tt>=toutp.OR.(tt>=h3d_data%TH3D.AND.tt<=h3d_data%TH3D_STOP) .OR.
799 . (manim>=4.AND.manim<=15).OR.h3d_data%MH3D /= 0))
800 . .OR.h3d_data%N_VECT_PCONT_MAX>0)THEN
801#include "lockon.inc"
802 DO i=1,jlt
803 fncont(1,ix1(i)) =fncont(1,ix1(i)) + fxi(i)*h1(i)
804 fncont(2,ix1(i)) =fncont(2,ix1(i)) + fyi(i)*h1(i)
805 fncont(3,ix1(i)) =fncont(3,ix1(i)) + fzi(i)*h1(i)
806 fncont(1,ix2(i)) =fncont(1,ix2(i)) + fxi(i)*h2(i)
807 fncont(2,ix2(i)) =fncont(2,ix2(i)) + fyi(i)*h2(i)
808 fncont(3,ix2(i)) =fncont(3,ix2(i)) + fzi(i)*h2(i)
809 fncont(1,ix3(i)) =fncont(1,ix3(i)) + fxi(i)*h3(i)
810 fncont(2,ix3(i)) =fncont(2,ix3(i)) + fyi(i)*h3(i)
811 fncont(3,ix3(i)) =fncont(3,ix3(i)) + fzi(i)*h3(i)
812 fncont(1,ix4(i)) =fncont(1,ix4(i)) + fxi(i)*h4(i)
813 fncont(2,ix4(i)) =fncont(2,ix4(i)) + fyi(i)*h4(i)
814 fncont(3,ix4(i)) =fncont(3,ix4(i)) + fzi(i)*h4(i)
815 ENDDO
816#include "lockoff.inc"
817 ENDIF
818C-----------------------------------------------------
819 IF(isecin > 0)THEN
820 k0=nstrf(25)
821 IF(nstrf(1)+nstrf(2) /= 0)THEN
822 DO i=1,nsect
823 nbinter=nstrf(k0+14)
824 k1s=k0+30
825 DO j=1,nbinter
826 IF(nstrf(k1s) == noint)THEN
827 IF(isecut /= 0)THEN
828#include "lockon.inc"
829 DO k=1,jlt
830 ! Be careful with sign conventions during force accumulation
831 ! Ensure consistency with CFORC3 implementation
832 IF(secfcum(4,ix1(k),i) == one)THEN
833 secfcum(1,ix1(k),i)=secfcum(1,ix1(k),i)-fx1(k)
834 secfcum(2,ix1(k),i)=secfcum(2,ix1(k),i)-fy1(k)
835 secfcum(3,ix1(k),i)=secfcum(3,ix1(k),i)-fz1(k)
836 ENDIF
837 IF(secfcum(4,ix2(k),i) == one)THEN
838 secfcum(1,ix2(k),i)=secfcum(1,ix2(k),i)-fx2(k)
839 secfcum(2,ix2(k),i)=secfcum(2,ix2(k),i)-fy2(k)
840 secfcum(3,ix2(k),i)=secfcum(3,ix2(k),i)-fz2(k)
841 ENDIF
842 IF(secfcum(4,ix3(k),i) == one)THEN
843 secfcum(1,ix3(k),i)=secfcum(1,ix3(k),i)-fx3(k)
844 secfcum(2,ix3(k),i)=secfcum(2,ix3(k),i)-fy3(k)
845 secfcum(3,ix3(k),i)=secfcum(3,ix3(k),i)-fz3(k)
846 ENDIF
847 IF(secfcum(4,ix4(k),i) == one)THEN
848 secfcum(1,ix4(k),i)=secfcum(1,ix4(k),i)-fx4(k)
849 secfcum(2,ix4(k),i)=secfcum(2,ix4(k),i)-fy4(k)
850 secfcum(3,ix4(k),i)=secfcum(3,ix4(k),i)-fz4(k)
851 ENDIF
852 jg = nsvg(k)
853 IF(jg > 0) THEN
854 ! SPMD : reprocess required after receiving a remote node if JG < 0
855 IF(secfcum(4,jg,i) == 1.)THEN
856 secfcum(1,jg,i)=secfcum(1,jg,i)+fxi(k)
857 secfcum(2,jg,i)=secfcum(2,jg,i)+fyi(k)
858 secfcum(3,jg,i)=secfcum(3,jg,i)+fzi(k)
859 ENDIF
860 ENDIF
861 ENDDO
862#include "lockoff.inc"
863 ENDIF
864 ENDIF
865 k1s=k1s+1
866 enddo!NEXT J
867 k0=nstrf(k0+24)
868 enddo!NEXT I
869 endif!(NSTRF(1)+NSTRF(2) /= 0)
870 endif!(ISECIN > 0)
871C-----------------------------------------------------
872 RETURN
#define my_real
Definition cppsort.cpp:32
#define alpha
Definition eval.h:35
subroutine area(d1, x, x2, y, y2, eint, stif0)
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 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
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine multi_i18_force_poff(dt, 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, jtask, multi_fvm, x, ixs, v, elbuf_tab, igroups, iparg, msi)
subroutine multi_i18_force_pon(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, multi_fvm, dt, jtask)
type(fani_cell_) fani_cell
Definition aleanim_mod.F:55
type(real_pointer), dimension(:), allocatable stifi
Definition tri7box.F:449
type(int_pointer), dimension(:), allocatable lisubsfi
Definition tri7box.F:501
type(int_pointer), dimension(:), allocatable nsvfi
Definition tri7box.F:431
type(int_pointer), dimension(:), allocatable addsubsfi
Definition tri7box.F:509
type(int_pointer), dimension(:), allocatable itafi
Definition tri7box.F:440