OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i7tri.F File Reference
#include "implicit_f.inc"
#include "mvsiz_p.inc"
#include "com01_c.inc"
#include "com04_c.inc"
#include "param_c.inc"
#include "ige3d_c.inc"
#include "vectorize.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine i7tri (add, nsn, renum, nsnr, isznsnr, irect, x, stf, stfn, xyzm, i_add, nsv, maxsiz, ii_stok, cand_n, cand_e, mulnsn, noint, tzinf, maxbox, minbox, i_mem, nb_n_b, i_add_max, eshift, inacti, ifq, cand_a, cand_p, ifpen, nrtm, nsnrold, igap, gap, gap_s, gap_m, gapmin, gapmax, marge, curv_max, nin, gap_s_l, gap_m_l, intth, drad, itied, cand_f, kremnod, remnod, flagremnode, dgapload, intheat, idt_therm, nodadt_therm)

Function/Subroutine Documentation

◆ i7tri()

subroutine i7tri ( integer, dimension(2,*) add,
integer nsn,
integer, dimension(*) renum,
integer nsnr,
integer isznsnr,
integer, dimension(4,*) irect,
x,
stf,
stfn,
xyzm,
integer i_add,
integer, dimension(*) nsv,
integer maxsiz,
integer ii_stok,
integer, dimension(*) cand_n,
integer, dimension(*) cand_e,
integer mulnsn,
integer noint,
tzinf,
maxbox,
minbox,
integer i_mem,
integer nb_n_b,
integer i_add_max,
integer eshift,
integer inacti,
integer ifq,
integer, dimension(*) cand_a,
cand_p,
integer, dimension(*) ifpen,
integer nrtm,
integer nsnrold,
integer igap,
gap,
gap_s,
gap_m,
gapmin,
gapmax,
marge,
curv_max,
integer nin,
gap_s_l,
gap_m_l,
integer intth,
intent(in) drad,
integer itied,
cand_f,
integer, dimension(*) kremnod,
integer, dimension(*) remnod,
integer flagremnode,
intent(in) dgapload,
integer, intent(in) intheat,
integer, intent(in) idt_therm,
integer, intent(in) nodadt_therm )

Definition at line 34 of file i7tri.F.

46C=======================================================================
47C M o d u l e s
48C-----------------------------------------------
49 USE tri7box
50C-----------------------------------------------
51C I m p l i c i t T y p e s
52C-----------------------------------------------
53#include "implicit_f.inc"
54C-----------------------------------------------
55C G l o b a l P a r a m e t e r s
56C-----------------------------------------------
57#include "mvsiz_p.inc"
58c parameter setting the size for the vector (orig version is 128)
59 INTEGER NVECSZ
60 parameter(nvecsz = mvsiz)
61C-----------------------------------------------
62C C o m m o n B l o c k s
63C-----------------------------------------------
64#include "com01_c.inc"
65#include "com04_c.inc"
66#include "param_c.inc"
67#include "ige3d_c.inc"
68C-----------------------------------------------
69C role of the routine:
70C ===================
71C classify BPE elements and BPN nodes in two zones
72C > or < to a boundary determined here and output all
73C in bpe,hpe, and bpn,hpn
74C-----------------------------------------------
75C D u m m y A r g u m e n t s
76C
77C NOM DESCRIPTION E/S
78C
79C BPE ARRAY OF FACETTES TO SORT => Local
80C and from the result max side
81C PE ARRAY OF FACETTES => Local
82C RESULTAT COTE MIN
83C BPN SORTED NODES ARRAY => Local
84C and from the result max side
85C PN NODES ARRAY => Local
86C RESULTAT COTE MIN
87C ADD(2,*) ARRAY OF ADRESSES E/S
88C 1.......ADRESSES NODES C 2.......ADRESSES ELEMENTS
89C ZYZM(6,*) ARRAY OF XYZMIN E/S
90C 1.......XMIN BOITE
91C 2.......YMIN BOITE
92C 3.......ZMIN BOITE
93C 4.......XMAX BOITE
94C 5.......YMAX BOITE
95C 6.......ZMAX BOITE
96C IRECT(4,*) ARRAY OF CONEC FACETTES E
97C X(3,*) COORDONNEES NODALES E
98C NB_NC NUMBER OF CANDIDATE NODES => Local
99C NB_EC NUMBER OF CANDIDATE ELEMENTS => Local
100C I_ADD position in the table of i/o addresses
101C NSV NOS SYSTEMES DES NODES E
102C Xmax larger abcisse existing e
103C XMAX largest order.existing E
104C Xmax larger existing side E
105C MAXSIZ TAILLE MEMOIRE MAX POSSIBLE E
106C I_STOK storage level for pairs
107C CANDIDATES impact E/S
108C ADNSTK current address in the node box
109C CAND_N boites resultats nodes C ADESTK current address in the element box
110C CAND_E adresses des boites resultat elements
111C MULNSN = MULTIMP*NSN maximum size now allowed for the
112C COUPLES NODES,ELT CANDIDATES
113C NOINT INTERFACE USER NUMBER
114C TZINF TAILLE ZONE INFLUENCE
115C MAXBOX TAILLE MAX BUCKET
116C MINBOX TAILLE MIN BUCKET
117C
118C Prov_n Provisional Cand_n (static variable in i7tri)
119C PROV_E CAND_E provisoire (variable static in i7tri)
120C-----------------------------------------------
121C D u m m y A r g u m e n t s
122C-----------------------------------------------
123 INTEGER I_ADD,MAXSIZ,I_MEM,ESHIFT,NSN,ISZNSNR,NRTM,NSNROLD,
124 . MULNSN,NB_N_B,NOINT,I_ADD_MAX,INACTI,IFQ,NSNR,IGAP,NIN,
125 . ADD(2,*),IRECT(4,*),
126 . NSV(*),CAND_N(*),CAND_E(*),CAND_A(*),IFPEN(*),RENUM(*),
127 . INTTH,II_STOK,ITIED
128 INTEGER KREMNOD(*),REMNOD(*),FLAGREMNODE
129 INTEGER, INTENT(IN) :: INTHEAT
130 INTEGER, INTENT(IN) :: IDT_THERM
131 INTEGER, INTENT(IN) :: NODADT_THERM
132C REAL
133 my_real
134 . x(3,*),xyzm(6,*),cand_p(*),stf(*),stfn(*),gap_s(*),gap_m(*),
135 . tzinf,maxbox,minbox,marge,gap,gapmin,gapmax,
136 . curv_max(*),gap_s_l(*),gap_m_l(*),cand_f(*)
137 my_real , INTENT(IN) :: drad,dgapload
138C-----------------------------------------------
139C L o c a l V a r i a b l e s
140C-----------------------------------------------
141 INTEGER NB_NCN,NB_NCN1,NB_ECN,ADDNN,ADDNE,I,J,DIR,NB_NC,NB_EC,
142 . N1,N2,N3,N4,NN,NE,K,L,NCAND_PROV,J_STOK,II,JJ,
143 . PROV_N(2*MVSIZ),PROV_E(2*MVSIZ),
144 .
145C BPE: Use on nrtm and not nrtm + 100 rigorously (here maxsiz = nrtm + 100)
146 . BPE(MAXSIZ/3),PE(MAXSIZ),BPN(NSN+NSNR),PN(NSN+NSNR),
147 . OLDNUM(ISZNSNR),IADD
148C REAL
149 my_real
150 . aaa,
151 . dx,dy,dz,dsup,trhreshold, xx1, xx2, xx3, xx4,
152 . xmin, xmax,ymin, ymax,zmin, zmax, tz, gapsmx, bgapsmx,
153 .
154 .
155 . smoins,splus,xx
156C REMNODE
157 INTEGER, DIMENSION(:), ALLOCATABLE :: TAGREMNODE
158 INTEGER DELNOD,M
159
160C-----------------------------------------------
161C initial construction phase of BPE and BPN moved from I7BUCE => I7TRI
162C
163
164 IF(flagremnode == 2) ALLOCATE(tagremnode(numnod+numfakenodigeo))
165
166 xmin = xyzm(1,i_add)
167 ymin = xyzm(2,i_add)
168 zmin = xyzm(3,i_add)
169 xmax = xyzm(4,i_add)
170 ymax = xyzm(5,i_add)
171 zmax = xyzm(6,i_add)
172C
173C DRAD = ZERO
174C
175C copy of segment and node numbers in BPE and BPN
176C
177 nb_ec = 0
178 DO i=1,nrtm
179C We do not retain the Destruit facets
180 IF(stf(i)/=zero)THEN
181 nb_ec = nb_ec + 1
182 bpe(nb_ec) = i
183 ENDIF
184 ENDDO
185
186 IF(igap==3) THEN
187 iadd = 10
188 ENDIF
189C
190C optimization // search for nodes within xmin xmax of the
191C processor elements
192C
193 nb_nc = 0
194 DO i=1,nsn
195 j=nsv(i)
196 IF(stfn(i)/=zero) THEN
197
198 IF(x(1,j)>=xmin.AND.x(1,j)<=xmax.AND.
199 . x(2,j)>=ymin.AND.x(2,j)<=ymax.AND.
200 . x(3,j)>=zmin.AND.x(3,j)<=zmax)THEN
201
202 nb_nc=nb_nc+1
203 bpn(nb_nc) = i
204 ENDIF
205 ENDIF
206 ENDDO
207C
208C Non -local candidate account in SPMD
209C
210 DO i = nsn+1, nsn+nsnr
211 IF( xrem(1,i-nsn)<xmin) cycle
212 IF( xrem(1,i-nsn)>xmax) cycle
213 IF( xrem(2,i-nsn)<ymin) cycle
214 IF( xrem(2,i-nsn)>ymax) cycle
215 IF( xrem(3,i-nsn)<zmin) cycle
216 IF( xrem(3,i-nsn)>zmax) cycle
217 nb_nc = nb_nc + 1
218 bpn(nb_nc) = i
219 ENDDO
220C
221C in SPMD, for inacti or IFQ, recovers old numbering of non-local candidates
222C
223 IF(nspmd>1.AND.
224 + (inacti==5.OR.inacti==6.OR.inacti==7.OR.ifq>0.OR.
225 + itied/=0)) THEN
226 CALL spmd_oldnumcd(renum,oldnum,isznsnr,nsnrold,intheat,idt_therm,nodadt_therm)
227 END IF
228C
229 j_stok = 0
230 GOTO 200
231C=======================================================================
232 100 CONTINUE
233C=======================================================================
234C-----------------------------------------------------------
235C
236C
237C 1- sorting phase on the median according to the largest direction
238C
239C
240C-----------------------------------------------------------
241C
242C 1- DETERMINER LA DIRECTION A DIVISER X,Y OU Z
243C
244 dir = 1
245 IF(dy==dsup) THEN
246 dir = 2
247 ELSE IF(dz==dsup) THEN
248 dir = 3
249 ENDIF
250 smoins = xyzm(dir,i_add)
251 splus = xyzm(dir+3,i_add)
252 trhreshold =(smoins+splus)*half
253C
254C 2- DIVISER LES NODES EN TWO ZONES
255C
256 nb_ncn= 0
257 nb_ncn1= 0
258 addnn= add(1,i_add)
259C
260 gapsmx = zero
261 DO i=1,nb_nc
262 j = bpn(i)
263 IF(j <= nsn) THEN
264 xx = x(dir,nsv(j))
265 IF(xx < trhreshold) THEN
266C we store at the bottom of the BP stack
267 nb_ncn1 = nb_ncn1 + 1
268 addnn = addnn + 1
269 pn(addnn) = j
270 IF(igap /=0) gapsmx = max(gapsmx,gap_s(j))
271 smoins = max(smoins,xx)
272 ENDIF
273 ENDIF
274 ENDDO
275 DO i=1,nb_nc
276 j = bpn(i)
277 IF(j > nsn) THEN
278 xx = xrem(dir,j-nsn)
279 IF(xx < trhreshold) THEN
280C we store at the bottom of the BP stack
281 nb_ncn1 = nb_ncn1 + 1
282 addnn = addnn + 1
283 pn(addnn) = j
284 IF(igap/=0) gapsmx = max(gapsmx,xrem(9,j-nsn))
285 smoins = max(smoins,xx)
286 ENDIF
287 ENDIF
288 ENDDO
289 bgapsmx = zero
290 DO i=1,nb_nc
291 j = bpn(i)
292 IF(j <= nsn) THEN
293 xx = x(dir,nsv(j))
294 IF(xx >= trhreshold) THEN
295C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPN
296 nb_ncn = nb_ncn + 1
297 bpn(nb_ncn) = j
298 IF(igap/=0) bgapsmx = max(bgapsmx,gap_s(j))
299 splus = min(splus,xx)
300 ENDIF
301 ENDIF
302 ENDDO
303 DO i=1,nb_nc
304 j = bpn(i)
305 IF(j > nsn) THEN
306 xx = xrem(dir,j-nsn)
307 IF(xx >= trhreshold) THEN
308C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPN
309 nb_ncn = nb_ncn + 1
310 bpn(nb_ncn) = j
311 IF(igap /= 0) bgapsmx = max(bgapsmx,xrem(9,j-nsn))
312 splus = min(splus,xx)
313 ENDIF
314 ENDIF
315 ENDDO
316C
317C 3- divide the elements
318C
319 nb_ecn= 0
320 addne= add(2,i_add)
321 IF(nb_ncn1==0) THEN
322 DO i=1,nb_ec
323 ne = bpe(i)
324 xx1=x(dir, irect(1,ne))
325 xx2=x(dir, irect(2,ne))
326 xx3=x(dir, irect(3,ne))
327 xx4=x(dir, irect(4,ne))
328 IF(igap == 0) THEN
329 aaa = tzinf+curv_max(ne)
330 ELSEIF(igap == 3) THEN
331 aaa = max(drad,dgapload+min(max(bgapsmx+max(gap_m(ne),gap_m_l(ne)),gapmin),gapmax))
332 + +marge+curv_max(ne)
333 ELSE
334 aaa = max(drad,dgapload+min(max(bgapsmx+gap_m(ne),gapmin),gapmax))
335 + +marge+curv_max(ne)
336 ENDIF
337 xmax = max(xx1,xx2,xx3,xx4) + aaa
338 IF(xmax >= splus) THEN
339C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPE
340 nb_ecn = nb_ecn + 1
341 bpe(nb_ecn) = ne
342 ENDIF
343 ENDDO
344 ELSEIF(nb_ncn == 0) THEN
345#include "vectorize.inc"
346 DO i=1,nb_ec
347 ne = bpe(i)
348 xx1=x(dir, irect(1,ne))
349 xx2=x(dir, irect(2,ne))
350 xx3=x(dir, irect(3,ne))
351 xx4=x(dir, irect(4,ne))
352 IF( igap == 0 ) THEN
353 aaa = -tzinf-curv_max(ne)
354 ELSEIF(igap == 3) THEN
355 aaa = -max(drad,dgapload+min(max(gapsmx+max(gap_m(ne),gap_m_l(ne)),gapmin),gapmax))
356 + -marge-curv_max(ne)
357 ELSE
358 aaa = -max(drad,dgapload+min(max(gapsmx+gap_m(ne),gapmin),gapmax))
359 - -marge-curv_max(ne)
360 ENDIF
361 xmin = min(xx1,xx2,xx3,xx4) + aaa
362
363 IF(xmin < smoins) THEN
364C we store at the bottom of the BP stack
365 addne = addne + 1
366 pe(addne) = ne
367 ENDIF
368 ENDDO
369 ELSE
370 DO i=1,nb_ec
371 ne = bpe(i)
372 xx1=x(dir, irect(1,ne))
373 xx2=x(dir, irect(2,ne))
374 xx3=x(dir, irect(3,ne))
375 xx4=x(dir, irect(4,ne))
376 IF( igap == 0 ) THEN
377 aaa=-tzinf-curv_max(ne)
378 ELSEIF(igap == 3) THEN
379 aaa= - max(drad,dgapload+min(max(gapsmx+max(gap_m(ne),gap_m_l(ne)),gapmin),gapmax))
380 + -marge-curv_max(ne)
381 ELSE
382 aaa= -max(drad,dgapload+min(max(gapsmx+gap_m(ne),gapmin),gapmax))
383 - -marge-curv_max(ne)
384 ENDIF
385 xmin = min(xx1,xx2,xx3,xx4) + aaa
386 IF(xmin < smoins) THEN
387C we store at the bottom of the BP stack
388 addne = addne + 1
389 pe(addne) = ne
390 ENDIF
391 ENDDO
392C
393 DO i=1,nb_ec
394 ne = bpe(i)
395 xx1=x(dir, irect(1,ne))
396 xx2=x(dir, irect(2,ne))
397 xx3=x(dir, irect(3,ne))
398 xx4=x(dir, irect(4,ne))
399 IF( igap == 0) THEN
400 aaa =tzinf+curv_max(ne)
401 ELSEIF( igap==3 ) THEN
402 aaa= max(drad,dgapload+min(max(bgapsmx+max(gap_m(ne),gap_m_l(ne)),gapmin),gapmax))
403 + +marge+curv_max(ne)
404 ELSE
405 aaa = max(drad,dgapload+min(max(bgapsmx+gap_m(ne),gapmin),gapmax))
406 + +marge+curv_max(ne)
407 ENDIF
408 xmax = max(xx1,xx2,xx3,xx4) + aaa
409
410 IF(xmax >= splus) THEN
411C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPE
412 nb_ecn = nb_ecn + 1
413 bpe(nb_ecn) = ne
414 ENDIF
415 ENDDO
416 ENDIF !NBNC1
417C
418C 4- REMPLIR LES TABLEAUX D'ADRESSES
419C
420 add(1,i_add+1) = addnn
421 add(2,i_add+1) = addne
422C-----we fill the min of the next box and the max of the current one
423C (i.e. THRESHOLD is a max for the current one)
424C We're going to go down and so we define a new box
425C we fill the max of the new box
426C initialises in i7buc1 a 1.E30 comme ca on recupere
427C either XMAX or the max of the box
428 xyzm(1,i_add+1) = xyzm(1,i_add)
429 xyzm(2,i_add+1) = xyzm(2,i_add)
430 xyzm(3,i_add+1) = xyzm(3,i_add)
431 xyzm(4,i_add+1) = xyzm(4,i_add)
432 xyzm(5,i_add+1) = xyzm(5,i_add)
433 xyzm(6,i_add+1) = xyzm(6,i_add)
434 xyzm(dir,i_add+1) = splus
435 xyzm(dir+3,i_add) = smoins
436C
437 nb_nc = nb_ncn
438 nb_ec = nb_ecn
439C we increment the descent level before exiting
440 i_add = i_add + 1
441 IF(i_add+1>=i_add_max) THEN
442c WRITE(6,*) __LINE__,__LINE__
443 i_mem = 3
444 RETURN
445 ENDIF
446C=======================================================================
447 200 CONTINUE
448C=======================================================================
449C-----------------------------------------------------------
450C
451C
452C 2- TEST ARRET = BOITE VIDE
453C BOITE TROP PETITE
454C BOITE NE CONTENANT QU'ONE NODE C No More Memory Available
455C
456C-------------------test for memory exceeded------------
457C
458 IF(add(2,i_add)+nb_ec>maxsiz) THEN
459C no more room in the stack of elements boxes too small
460 WRITE(6,*) __line__,__line__
461
462 i_mem = 1
463 RETURN
464 ENDIF
465C
466C--------------------test for empty boxes--------------
467C
468 IF(nb_ec/=0.AND.nb_nc/=0) THEN
469C
470 dx = xyzm(4,i_add) - xyzm(1,i_add)
471 dy = xyzm(5,i_add) - xyzm(2,i_add)
472 dz = xyzm(6,i_add) - xyzm(3,i_add)
473 dsup= max(dx,dy,dz)
474C
475C-------------------test for end of branch ------------
476C 1- storage of candidate node(s) and corresponding elements
477C remove the useless ones
478C
479C NCAND_PROV=NB_EC*NB_NC
480C NCAND_PROV negatif qd NB_EC*NB_NC > 2e31
481C
482 IF(nb_ec+nb_nc<=nvecsz) THEN
483 ncand_prov = nb_ec*nb_nc
484 ELSE
485 ncand_prov = nvecsz+1
486 ENDIF
487 IF(dsup<minbox.OR.(nb_nc<=nb_n_b)
488 & .OR.(ncand_prov<=nvecsz)) THEN
489 ncand_prov = nb_ec*nb_nc
490
491 IF(flagremnode==2) THEN
492 DO i=1,numnod+numfakenodigeo
493 tagremnode(i) = 0
494 ENDDO
495 ENDIF
496
497 DO k=1,ncand_prov,nvsiz
498 DO l=k,min(k-1+nvsiz,ncand_prov)
499 i = 1+(l-1)/nb_nc
500 j = l-(i-1)*nb_nc
501 ne = bpe(i)
502 n1=irect(1,ne)
503 n2=irect(2,ne)
504 n3=irect(3,ne)
505 n4=irect(4,ne)
506
507 IF(flagremnode==2) THEN
508 DO m= kremnod(2*(ne-1)+1)+1, kremnod(2*(ne-1)+2)
509 tagremnode(remnod(m)) = 1
510 ENDDO
511 ENDIF
512 jj = bpn(j)
513 IF( jj<=nsn ) THEN
514 IF( igap == 0 ) THEN
515 tz = tzinf+curv_max(ne)
516 ELSEIF( igap == 3 ) THEN
517 tz = max(drad,dgapload+max(min(gap_s_l(jj)+gap_m_l(ne),gapmax),gapmin)
518 . +marge+curv_max(ne))
519 ELSE
520 tz=max(drad,dgapload+max(min(gap_s(jj)+gap_m(ne),gapmax),gapmin)
521 + +marge+curv_max(ne))
522 ENDIF
523 ELSE
524 ii = jj-nsn
525 IF( igap == 0 ) THEN
526 tz = tzinf+curv_max(ne)
527 ELSEIF( igap == 3 ) THEN
528 tz = max(drad,dgapload+max(min(xrem(iadd,ii)+gap_m_l(ne)
529 . ,gapmax),gapmin))+marge+curv_max(ne)
530 ELSE
531 tz = max(drad,dgapload+max(min(xrem(9,ii)+gap_m(ne),gapmax),gapmin))
532 + +marge+curv_max(ne)
533 ENDIF
534 ENDIF
535 xx1=x(1, n1)
536 xx2=x(1, n2)
537 xx3=x(1, n3)
538 xx4=x(1, n4)
539 xmax=max(xx1,xx2,xx3,xx4)+tz
540 xmin=min(xx1,xx2,xx3,xx4)-tz
541 xx1=x(2, n1)
542 xx2=x(2, n2)
543 xx3=x(2, n3)
544 xx4=x(2, n4)
545 ymax=max(xx1,xx2,xx3,xx4)+tz
546 ymin=min(xx1,xx2,xx3,xx4)-tz
547 xx1=x(3, n1)
548 xx2=x(3, n2)
549 xx3=x(3, n3)
550 xx4=x(3, n4)
551 zmax=max(xx1,xx2,xx3,xx4)+tz
552 zmin=min(xx1,xx2,xx3,xx4)-tz
553 IF(jj<=nsn) THEN
554
555 IF(flagremnode==2) THEN
556 IF(tagremnode(nsv(jj)) == 1) cycle
557 ENDIF
558 nn=nsv(jj)
559 IF(nn/=n1.AND.nn/=n2.AND.nn/=n3.AND.nn/=n4.AND.
560 & x(1,nn)>xmin.AND.x(1,nn)<xmax.AND.
561 & x(2,nn)>ymin.AND.x(2,nn)<ymax.AND.
562 & x(3,nn)>zmin.AND.x(3,nn)<zmax ) THEN
563 j_stok = j_stok + 1
564 prov_n(j_stok) = jj
565 prov_e(j_stok) = ne
566 ENDIF
567 ELSE
568 ii = jj-nsn
569 IF(flagremnode==2) THEN
570 DO m= kremnod(2*(ne-1)+2) + 1, kremnod(2*(ne-1)+3)
571 IF(remnod(m) == -irem(2,ii) ) THEN
572 delnod = delnod + 1
573 EXIT
574 ENDIF
575 ENDDO
576 IF(delnod /= 0) cycle
577 ENDIF
578 IF(xrem(1,ii)>xmin.AND.
579 & xrem(1,ii)<xmax.AND.
580 & xrem(2,ii)>ymin.AND.
581 & xrem(2,ii)<ymax.AND.
582 & xrem(3,ii)>zmin.AND.
583 & xrem(3,ii)<zmax ) THEN
584 j_stok = j_stok + 1
585 prov_n(j_stok) = jj
586 prov_e(j_stok) = ne
587 ENDIF
588 ENDIF
589 ENDDO ! L=K,MIN(K-1+NVSIZ,NCAND_PROV)
590
591 IF(j_stok>=nvsiz)THEN
592 CALL i7sto(
593 1 nvsiz,irect ,x ,nsv ,ii_stok,
594 2 cand_n,cand_e ,mulnsn,noint ,marge ,
595 3 i_mem ,prov_n ,prov_e,eshift,inacti ,
596 4 ifq ,cand_a ,cand_p,ifpen ,nsn ,
597 5 oldnum,nsnrold,igap ,gap ,gap_s ,
598 6 gap_m ,gapmin ,gapmax,curv_max,nin ,
599 7 gap_s_l,gap_m_l,intth,drad,itied ,
600 8 cand_f,dgapload)
601 IF(i_mem==2) THEN
602 RETURN
603 ENDIF
604 j_stok = j_stok-nvsiz
605#include "vectorize.inc"
606 DO j=1,j_stok
607 prov_n(j) = prov_n(j+nvsiz)
608 prov_e(j) = prov_e(j+nvsiz)
609 ENDDO
610 ENDIF ! J_STOK >= NVSIZ
611
612 ENDDO ! LOOP OVER POSSIBLE CANDIDATES
613 ELSE
614C=======================================================================
615 GOTO 100
616C=======================================================================
617 ENDIF
618 ENDIF
619C-------------------------------------------------------------------------
620C empty box or
621C end of branch
622C we decrement the descent level before restarting
623C-------------------------------------------------------------------------
624 i_add = i_add - 1
625 IF (i_add/=0) THEN
626C-------------------------------------------------------------------------
627C the bottom of the stacks must be copied into corresponding stack bottoms
628C before going back down into the adjacent branch
629C-------------------------------------------------------------------------
630 CALL i7dstk(nb_nc,nb_ec,add(1,i_add),bpn,pn,bpe,pe)
631C=======================================================================
632 GOTO 200
633C=======================================================================
634 ENDIF
635C-------------------------------------------------------------------------
636C end of sorting
637C-------------------------------------------------------------------------
638 IF(j_stok/=0)CALL i7sto(
639 1 j_stok,irect ,x ,nsv ,ii_stok,
640 2 cand_n,cand_e ,mulnsn,noint ,marge ,
641 3 i_mem ,prov_n ,prov_e,eshift,inacti ,
642 4 ifq ,cand_a ,cand_p,ifpen ,nsn ,
643 5 oldnum,nsnrold,igap ,gap ,gap_s ,
644 6 gap_m ,gapmin ,gapmax,curv_max,nin ,
645 7 gap_s_l,gap_m_l,intth,drad,itied ,
646 8 cand_f ,dgapload)
647C-------------------------------------------------------------------------
648 IF(flagremnode==2) THEN
649 DEALLOCATE(tagremnode)
650 ENDIF
651 RETURN
#define my_real
Definition cppsort.cpp:32
subroutine i7sto(j_stok, irect, x, nsv, ii_stok, cand_n, cand_e, mulnsn, noint, marge, i_mem, prov_n, prov_e, eshift, inacti, ifq, cand_a, cand_p, ifpen, nsn, oldnum, nsnrold, igap, gap, gap_s, gap_m, gapmin, gapmax, curv_max, nin, gap_s_l, gap_m_l, intth, drad, itied, cand_f, dgapload)
Definition i7sto.F:41
subroutine ymax(idn, fac, npc, pld, stiffmin, stiffmax, stiffini, stiffavg)
Definition law100_upd.F:274
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
integer, dimension(:,:), allocatable irem
Definition tri7box.F:339
subroutine spmd_oldnumcd(renum, oldnum, nsnr, nsnrold, intheat, idt_therm, nodadt_therm)
Definition spmd_i7tool.F:38
subroutine i7dstk(i_add, nb_nc, nb_ec, add, bpn, pn, bpe, pe)
Definition i7dstk.F:33