OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i24tri.F File Reference
#include "implicit_f.inc"
#include "mvsiz_p.inc"
#include "param_c.inc"
#include "com04_c.inc"
#include "vect07_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine i24tri (bpe, pe, bpn, pn, add, irect, x, nb_nc, nb_ec, xyzm, i_add, nsv, i_amax, xmax, ymax, zmax, maxsiz, i_stok, i_mem, nb_n_b, cand_n, cand_e, nsn, noint, tzinf, maxbox, minbox, stf, stfn, j_stok, multimp, istf, itab, gap, gap_s, gap_m, igap, gapmin, gapmax, marge, gap_s_l, gap_m_l, id, titr, ilev, nbinflg, mbinflg, mvoisn, ixs, ixs10, ixs16, ixs20, ipartns, ipen0, inacti, msegtyp, marge_sh, nrtm, irtse, is2se, ix1, ix2, ix3, ix4, nsvg, x1, x2, x3, x4, y1, y2, y3, y4, z1, z2, z3, z4, xi, yi, zi, x0, y0, z0, stif, nx1, ny1, nz1, nx2, ny2, nz2, nx3, ny3, nz3, nx4, ny4, nz4, p1, p2, p3, p4, lb1, lb2, lb3, lb4, lc1, lc2, lc3, lc4, pene, prov_n, prov_e, n11, n21, n31, dgapload, s_kremnode, s_remnode, kremnode, remnode, tag_removed_node, flag_removed_node)

Function/Subroutine Documentation

◆ i24tri()

subroutine i24tri ( integer, dimension(*) bpe,
integer, dimension(*) pe,
integer, dimension(*) bpn,
integer, dimension(*) pn,
integer, dimension(2,0:*) add,
integer, dimension(4,*) irect,
x,
integer nb_nc,
integer nb_ec,
xyzm,
integer i_add,
integer, dimension(*) nsv,
integer i_amax,
xmax,
ymax,
zmax,
integer maxsiz,
integer i_stok,
integer i_mem,
integer nb_n_b,
integer, dimension(*) cand_n,
integer, dimension(*) cand_e,
integer nsn,
integer noint,
tzinf,
maxbox,
minbox,
stf,
stfn,
integer j_stok,
integer multimp,
integer istf,
integer, dimension(*) itab,
gap,
gap_s,
gap_m,
integer igap,
gapmin,
gapmax,
marge,
gap_s_l,
gap_m_l,
integer id,
character(len=nchartitle) titr,
integer ilev,
integer, dimension(*) nbinflg,
integer, dimension(*) mbinflg,
integer, dimension(4,*) mvoisn,
integer, dimension(nixs,*) ixs,
integer, dimension(6,*) ixs10,
integer, dimension(8,*) ixs16,
integer, dimension(12,*) ixs20,
integer, dimension(*) ipartns,
integer ipen0,
integer inacti,
integer, dimension(*) msegtyp,
marge_sh,
integer nrtm,
integer, dimension(*) irtse,
integer, dimension(*) is2se,
integer, dimension(mvsiz), intent(inout) ix1,
integer, dimension(mvsiz), intent(inout) ix2,
integer, dimension(mvsiz), intent(inout) ix3,
integer, dimension(mvsiz), intent(inout) ix4,
integer, dimension(mvsiz), intent(inout) nsvg,
intent(inout) x1,
intent(inout) x2,
intent(inout) x3,
intent(inout) x4,
intent(inout) y1,
intent(inout) y2,
intent(inout) y3,
intent(inout) y4,
intent(inout) z1,
intent(inout) z2,
intent(inout) z3,
intent(inout) z4,
intent(inout) xi,
intent(inout) yi,
intent(inout) zi,
intent(inout) x0,
intent(inout) y0,
intent(inout) z0,
intent(inout) stif,
intent(inout) nx1,
intent(inout) ny1,
intent(inout) nz1,
intent(inout) nx2,
intent(inout) ny2,
intent(inout) nz2,
intent(inout) nx3,
intent(inout) ny3,
intent(inout) nz3,
intent(inout) nx4,
intent(inout) ny4,
intent(inout) nz4,
intent(inout) p1,
intent(inout) p2,
intent(inout) p3,
intent(inout) p4,
intent(inout) lb1,
intent(inout) lb2,
intent(inout) lb3,
intent(inout) lb4,
intent(inout) lc1,
intent(inout) lc2,
intent(inout) lc3,
intent(inout) lc4,
intent(inout) pene,
integer, dimension(mvsiz), intent(inout) prov_n,
integer, dimension(mvsiz), intent(inout) prov_e,
intent(inout) n11,
intent(inout) n21,
intent(inout) n31,
intent(in) dgapload,
integer, intent(in) s_kremnode,
integer, intent(in) s_remnode,
integer, dimension(s_kremnode), intent(in) kremnode,
integer, dimension(s_remnode), intent(in) remnode,
integer, dimension(numnod), intent(inout) tag_removed_node,
logical, intent(in) flag_removed_node )
Parameters
[in]flag_removed_nodeflag to remove some S node from the list of candidates
[in]s_kremnodesize of KREMNODE
[in]s_remnodesize of REMNODE

Definition at line 39 of file i24tri.F.

65#ifndef HYPERMESH_LIB
66 USE message_mod
67#endif
69 use element_mod , only :nixs
70C-----------------------------------------------
71C I m p l i c i t T y p e s
72C-----------------------------------------------
73#include "implicit_f.inc"
74C-----------------------------------------------
75C G l o b a l P a r a m e t e r s
76C-----------------------------------------------
77#include "mvsiz_p.inc"
78#include "param_c.inc"
79C-----------------------------------------------
80C C o m m o n B l o c k s
81C-----------------------------------------------
82#include "com04_c.inc"
83#include "vect07_c.inc"
84C-----------------------------------------------
85C D u m m y A r g u m e n t s
86C-----------------------------------------------
87 INTEGER NB_NC,NB_EC,I_ADD,MAXSIZ,I_STOK,J_STOK,I_MEM,ISTF
88 INTEGER I_BID, I_AMAX,NB_N_B, NOINT, NSN,MULTIMP, IGAP
89 INTEGER ADD(2,0:*),IRECT(4,*),BPE(*),PE(*),BPN(*),PN(*)
90 INTEGER NSV(*),CAND_N(*),CAND_E(*), ITAB(*),NBINFLG(*),MBINFLG(*),
91 * ILEV,MVOISN(4,*),IPARTNS(*),IPEN0,INACTI,NRTM
92 INTEGER IXS(NIXS,*), IXS10(6,*), IXS16(8,*), IXS20(12,*),IRTSE(*),IS2SE(*)
93C REAL
95 . x(3,*),xyzm(6,*),tzinf,stf(*),stfn(*),
96 . maxbox,minbox, xmax, ymax, zmax,
97 . gap, gap_s(*), gap_m(*),
98 . gapmin, gapmax, marge, gapsmx, bgapsmx,
99 . gap_s_l(*),gap_m_l(*),marge_sh
100 my_real , INTENT(IN) :: dgapload
101 INTEGER ID,MSEGTYP(*)
102 LOGICAL, INTENT(in) :: FLAG_REMOVED_NODE !< flag to remove some S node from the list of candidates
103 INTEGER, INTENT(in) :: S_KREMNODE !< size of KREMNODE
104 INTEGER, INTENT(in) :: S_REMNODE !< size of REMNODE
105 INTEGER, DIMENSION(S_KREMNODE), INTENT(in) :: KREMNODE
106 INTEGER, DIMENSION(S_REMNODE), INTENT(in) :: REMNODE
107 INTEGER, DIMENSION(NUMNOD), INTENT(inout) :: TAG_REMOVED_NODE
108 CHARACTER(LEN=NCHARTITLE) :: TITR
109 INTEGER, DIMENSION(MVSIZ), INTENT(INOUT) ::PROV_N,PROV_E
110 INTEGER, DIMENSION(MVSIZ), INTENT(INOUT) :: IX1,IX2,IX3,IX4,NSVG
111 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: x1,x2,x3,x4
112 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: y1,y2,y3,y4
113 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: z1,z2,z3,z4
114 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: xi,yi,zi
115 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: x0,y0,z0,stif
116 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: n11,n21,n31,pene
117 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: nx1,ny1,nz1
118 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: nx2,ny2,nz2
119 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: nx3,ny3,nz3
120 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: nx4,ny4,nz4
121 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: p1,p2,p3,p4
122 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: lb1,lb2,lb3,lb4
123 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: lc1,lc2,lc3,lc4
124C-----------------------------------------------
125C L o c a l V a r i a b l e s
126C-----------------------------------------------
127 INTEGER NB_NCN,NB_ECN,ADDNN,ADDNE,I,J
128 INTEGER INF,SUP,DIR,N1,N2,N3,N4,NN,NE
129 INTEGER SKIP,NS,NS1,NS2,NSE
130 INTEGER :: FIRST,LAST
131 INTEGER :: IJK
132C REAL
133 my_real
134 . dx,dy,dz,dsup,seuil,xmx,xmn,
135 . gapv(mvsiz),marge_e
136C-----------------------------------------------
137C ROLE OF THE ROUTINE:
138C ===================
139C CLASSIFIES THE ELEMENTS OF BPE AND THE NODES OF BPN INTO TWO ZONES
140C > OR < TO A BORDER HERE DETERMINED AND OUTPUTS EVERYTHING
141C INTO bpe,hpe, and bpn,hpn
142C-----------------------------------------------
143C D u m m y A r g u m e n t s
144C
145C NOM DESCRIPTION E/S
146C
147C BPE ARRAY OF FACETTES TO SORT E/S
148C AND OF THE RESULT MAX SIDE
149C PE ARRAY OF FACETTES S
150C RESULTAT COTE MIN
151C BPN SORTED NODES ARRAY E/S
152C AND OF THE RESULT MAX SIDE
153C PN NODES ARRAY S
154C RESULTAT COTE MIN
155C ADD(2,*) ARRAY OF ADRESSES E/S
156C 1.......ADRESSES NODES C 2.......ADRESSES ELEMENTS
157C ZYZM(6,*) ARRAY OF XYZMIN E/S
158C 1.......XMIN BOITE
159C 2.......YMIN BOITE
160C 3.......ZMIN BOITE
161C 4.......XMAX BOITE
162C 5.......YMAX BOITE
163C 6.......ZMAX BOITE
164C IRECT(4,*) ARRAY OF CONEC FACETTES E
165C X(3,*) COORDONNEES NODALES E
166C NB_NC NUMBER OF CANDIDATE NODES E/S
167C NB_EC NUMBER OF CANDIDATE ELEMENTS E/S
168C I_ADD POSITION IN THE TAB OF E/S ADDRESSES
169C NSV NOS SYSTEMES DES NODES E
170C Xmax larger abcisse existing e
171C XMAX largest order.existing E
172C Xmax larger existing side E
173C MAXSIZ TAILLE MEMOIRE MAX POSSIBLE E
174C I_STOK storage level of couples
175C CANDIDATES impact E/S
176C CAND_N boites resultats nodes C CAND_E adresses des boites resultat elements
177C MULTIMP*NSN MAX SIZE NOW ALLOWED FOR THE
178C COUPLES NODES,ELT CANDIDATES
179C NOINT INTERFACE USER NUMBER
180C TZINF TAILLE ZONE INFLUENCE
181C MAXBOX TAILLE MAX BUCKET
182C MINBOX TAILLE MIN BUCKET
183C=======================================================================
184C
185C
186C 1- TEST ARRET = BOITE VIDE
187C BOITE TROP PETITE
188C BOITE NE CONTENANT QU'ONE NODE C NO MORE MEMORY AVAILABLE
189C
190C-----------------------------------------------------------
191C
192C IF(MEMX>ADD(2,1)+NB_EC)THEN
193C WRITE(ISTDO,*)' *******MEM MAX=',MEMX
194C MEMX=-1
195C ELSEIF(MEMX/=-1)THEN
196C MEMX=ADD(2,1)+NB_EC
197C ENDIF
198C--------------------TEST ON EMPTY BOXES--------------
199C
200 IF(nb_ec==0.OR.nb_nc==0) THEN
201C write(6,*)" BOITE VIDE"
202C MUST COPY THE BOTTOMS OF THE STACKS INTO CORRESPONDING BAS_DE_PILE
203C BEFORE DESCENDING INTO THE ADJACENT BRANCH
204 CALL i7dstk(i_add,nb_nc,nb_ec,add,bpn,pn,bpe,pe)
205 RETURN
206 ENDIF
207C
208C-------------------TEST ON END OF BRANCH / OVERFLOW------------
209C
210 dx = xyzm(4,i_add) - xyzm(1,i_add)
211 dy = xyzm(5,i_add) - xyzm(2,i_add)
212 dz = xyzm(6,i_add) - xyzm(3,i_add)
213 dsup= max(dx,dy,dz)
214C
215 IF(add(2,1)+nb_ec>=maxsiz) THEN
216C NO MORE SPACE IN THE STACK OF ELEMENTS BOXES TOO SMALL
217 IF ( nb_n_b == numnod) THEN
218#ifndef HYPERMESH_LIB
219 CALL ancmsg(msgid=83,
220 . msgtype=msgerror,
221 . anmode=aninfo,
222 . i1=id,
223 . c1=titr)
224#endif
225 ENDIF
226 i_mem = 1
227 RETURN
228 ENDIF
229 IF(dsup<minbox.OR.
230 . nb_nc<=nb_n_b.AND.dsup<maxbox.OR.
231 . nb_nc<=nb_n_b.AND.nb_ec==1) THEN
232C
233C write(6,*)" NOUVELLE BOITE "
234C 1- STORAGE OF THE CANDIDATE NODE AND CORRESP. ELEMENTS.
235C REMOVE THE USELESS ONES
236 DO i=1,nb_ec
237 ne = bpe(i)
238 ! --------------
239 ! do not take into account node_id = REMNODE(i)
240 IF(flag_removed_node) THEN
241 first = kremnode(ne)+1
242 last = kremnode(ne+1)
243 DO ijk=first,last
244 IF(remnode(ijk)<=numnod) tag_removed_node(remnode(ijk)) = 1
245 ENDDO
246 ENDIF
247 ! --------------
248 n1=irect(1,ne)
249 n2=irect(2,ne)
250 n3=irect(3,ne)
251 n4=irect(4,ne)
252 IF (msegtyp(ne)==0.OR.msegtyp(ne)>nrtm) THEN
253 marge_e = marge
254 ELSE
255 marge_e = marge_sh
256 END IF
257 DO j=1,nb_nc
258 nn=nsv(bpn(j))
259C-----------------
260C-------- not sure if use Inacti----
261C INS = 1
262C IF (INACTI > 0)
263C . CALL INSOLBOX(X ,MVOISN(1,NE),MVOISN(2,NE),NOINT ,IXS,
264C . IXS10 ,IXS16 ,IXS20 ,NN ,GAP ,
265C . MVOISN(3,NE),IPARTNS(BPN(J)),IPEN0,INS )
266 IF(nn/=n1.AND.nn/=n2.AND.nn/=n3.AND.nn/=n4) THEN
267C Fictive Nodes to filter
268 skip=0
269 IF(nn > numnod)THEN
270 ns=nn-numnod
271 CALL i24fic_getn(ns ,irtse ,is2se ,nse ,
272 + ns1 ,ns2 )
273 IF(ns1 == n1 .OR. ns2 == n1) skip=1
274 IF(ns1 == n2 .OR. ns2 == n2) skip=1
275 IF(ns1 == n3 .OR. ns2 == n3) skip=1
276 IF(ns1 == n4 .OR. ns2 == n4) skip=1
277
278 IF(skip==0) THEN
279 IF(flag_removed_node) THEN
280 first = kremnode(ne)+1
281 last = kremnode(ne+1)
282 DO ijk=first,last
283 IF(remnode(ijk)==nn) skip = 1
284 ENDDO
285 ENDIF
286 ENDIF
287 ELSE
288 ! --------------
289 ! S node must be remove from the list of candidates if TAG_REMOVED_NODE = 1
290 IF(flag_removed_node) THEN
291 IF(tag_removed_node(nn)==1) skip = 1
292 ENDIF
293 ! --------------
294 ENDIF
295
296 IF (skip==0)THEN
297 j_stok = j_stok + 1
298 prov_n(j_stok) = bpn(j)
299 prov_e(j_stok) = ne
300 IF(j_stok==nvsiz) THEN
301 lft = 1
302 llt = nvsiz
303 nft = 0
304 j_stok = 0
305 CALL i7cor3(x ,irect,nsv ,prov_e ,prov_n,
306 . stf ,stfn ,gapv ,igap ,gap ,
307 . gap_s,gap_m,istf ,gapmin ,gapmax,
308 . gap_s_l,gap_m_l ,zero ,ix1 ,ix2 ,
309 5 ix3 ,ix4 ,nsvg,x1 ,x2 ,
310 6 x3 ,x4 ,y1 ,y2 ,y3 ,
311 7 y4 ,z1 ,z2 ,z3 ,z4 ,
312 8 xi ,yi ,zi ,stif ,dgapload,
313 9 llt)
314 CALL i7dst3(ix3,ix4,x1 ,x2 ,x3 ,
315 1 x4 ,y1 ,y2 ,y3 ,y4 ,
316 2 z1 ,z2 ,z3 ,z4 ,xi ,
317 3 yi ,zi ,x0 ,y0 ,z0 ,
318 4 nx1,ny1,nz1,nx2,ny2,
319 5 nz2,nx3,ny3,nz3,nx4,
320 6 ny4,nz4,p1 ,p2 ,p3 ,
321 7 p4 ,lb1,lb2,lb3,lb4,
322 8 lc1,lc2,lc3,lc4,llt)
323
324 CALL i7pen3(marge_e,gapv,n11,n21,n31,
325 1 pene ,nx1 ,ny1,nz1,nx2,
326 2 ny2 ,nz2 ,nx3,ny3,nz3,
327 3 nx4 ,ny4 ,nz4,p1 ,p2 ,
328 4 p3 ,p4,llt)
329
330 IF (ilev==2) CALL i24s1s2(prov_n,prov_e,nbinflg,mbinflg,pene)
331 IF(i_stok+nvsiz<multimp*nsn) THEN
332 CALL i7cmp3(i_stok,cand_e ,cand_n,1,pene,
333 1 prov_n,prov_e)
334 ELSE
335 i_bid = 0
336 CALL i7cmp3(i_bid,cand_e,cand_n,0,pene,
337 1 prov_n,prov_e)
338 IF(i_stok+i_bid<multimp*nsn) THEN
339 CALL i7cmp3(i_stok,cand_e,cand_n,1,pene,
340 1 prov_n,prov_e)
341 ELSE
342 i_mem = 2
343c CALL ANSTCKI(NOINT)
344c CALL ANCWARN(103,ANINFO_BLIND_2)
345 RETURN
346 ENDIF
347 ENDIF
348 ENDIF
349C write(6,*)"Node candidat",BPN(J)
350C write(6,*)"Element candidat",NE
351 ENDIF
352 ENDIF
353 ENDDO
354
355 ! --------------
356 ! flush to 0
357 IF(flag_removed_node) THEN
358 first = kremnode(ne)+1
359 last = kremnode(ne+1)
360 DO ijk=first,last
361 IF(remnode(ijk)<=numnod) tag_removed_node(remnode(ijk)) = 0
362 ENDDO
363 ENDIF
364 ! --------------
365 ENDDO
366C MUST COPY THE BOTTOMS OF THE STACKS INTO CORRESPONDING BAS_DE_PILE
367C BEFORE DESCENDING INTO THE ADJACENT BRANCH
368 CALL i7dstk(i_add,nb_nc,nb_ec,add,bpn,pn,bpe,pe)
369 RETURN
370 ENDIF
371C
372C-----------------------------------------------------------
373C
374C
375C 2- SORTING PHASE ON THE MEDIAN ACCORDING TO THE LARGEST DIRECTION
376C
377C
378C-----------------------------------------------------------
379C
380C
381C 1- DETERMINER LA DIRECTION A DIVISER X,Y OU Z
382C
383 dir = 1
384 IF(dy==dsup) THEN
385 dir = 2
386 ELSE IF(dz==dsup) THEN
387 dir = 3
388 ENDIF
389 seuil =(xyzm(dir+3,i_add)+xyzm(dir,i_add))/2
390C
391C 2- DIVISER LES NODES EN TWO ZONES
392C
393 nb_ncn= 0
394 addnn= add(1,1)
395 inf = 0
396 sup = 0
397 IF(igap==0)THEN
398 DO i=1,nb_nc
399 IF(x(dir,nsv(bpn(i)))<seuil) THEN
400C STORED IN THE BOTTOM OF THE BP STACK
401 addnn = addnn + 1
402 pn(addnn) = bpn(i)
403 inf = 1
404 ELSE
405 nb_ncn = nb_ncn + 1
406 bpn(nb_ncn) = bpn(i)
407C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPN
408 sup = 1
409 ENDIF
410 END DO
411 ELSE
412 gapsmx = zero
413 bgapsmx = zero
414 DO i=1,nb_nc
415 IF(x(dir,nsv(bpn(i)))<seuil) THEN
416C STORED IN THE BOTTOM OF THE BP STACK
417 addnn = addnn + 1
418 pn(addnn) = bpn(i)
419 gapsmx = max(gapsmx,gap_s(bpn(i)))
420 inf = 1
421 ELSE
422C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPN
423 nb_ncn = nb_ncn + 1
424 bpn(nb_ncn) = bpn(i)
425 bgapsmx = max(bgapsmx,gap_s(bpn(i)))
426 sup = 1
427 ENDIF
428 END DO
429 END IF
430
431CC
432CC 3- DIVISER LES ELEMENTS
433CC
434C NB_ECN= 0
435C ADDNE= ADD(2,1)
436C SEUILI = SEUIL-TZINF
437C SEUILS = SEUIL+TZINF
438C DO 85 I=1,NB_EC
439C INF = 0
440C SUP = 0
441C DO 80 J=1,4
442C IP = IRECT(J,BPE(I))
443C IF(X(DIR,IP)<SEUILS) THEN
444C INF = 1
445C IF(SUP==1) GOTO 81
446C ENDIF
447C IF(X(DIR,IP)>=SEUILI) THEN
448C SUP = 1
449C IF(INF==1) GOTO 81
450C ENDIF
451C 80 CONTINUE
452C
453C 81 CONTINUE
454C IF(INF==1) THEN
455C STORED IN THE BOTTOM OF THE BP STACK
456C ADDNE = ADDNE + 1
457C PE(ADDNE) = BPE(I)
458C ENDIF
459C IF(SUP==1) THEN
460C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPE
461C NB_ECN = NB_ECN + 1
462C BPE(NB_ECN) = BPE(I)
463C ENDIF
464C 85 CONTINUE
465C
466C 3- DIVIDE THE ELEMENTS
467C
468C 2 LIGNES PROV FOR TEST
469C INF = 1
470C SUP = 1
471C
472 nb_ecn= 0
473 addne= add(2,1)
474 IF(igap==0)THEN
475 DO i=1,nb_ec
476 xmx = max(x(dir,irect(1,bpe(i))),x(dir,irect(2,bpe(i))),
477 . x(dir,irect(3,bpe(i))),x(dir,irect(4,bpe(i))))
478 . + tzinf
479 xmn = min(x(dir,irect(1,bpe(i))),x(dir,irect(2,bpe(i))),
480 . x(dir,irect(3,bpe(i))),x(dir,irect(4,bpe(i))))
481 . - tzinf
482 IF(xmn<seuil.AND.inf==1) THEN
483C STORED IN THE BOTTOM OF THE BP STACK
484 addne = addne + 1
485 pe(addne) = bpe(i)
486 ENDIF
487 IF(xmx>=seuil.AND.sup==1) THEN
488C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPE
489 nb_ecn = nb_ecn + 1
490 bpe(nb_ecn) = bpe(i)
491 ENDIF
492 ENDDO
493 ELSE
494 DO i=1,nb_ec
495 ne = bpe(i)
496 IF (msegtyp(ne)==0.OR.msegtyp(ne)>nrtm) THEN
497 marge_e = marge
498 ELSE
499 marge_e = marge_sh
500 END IF
501 xmn = min(x(dir,irect(1,ne)),x(dir,irect(2,ne)),
502 . x(dir,irect(3,ne)),x(dir,irect(4,ne)))
503 . - max(min(gapsmx+gap_m(ne),gapmax),gapmin)+dgapload-marge_e
504 IF(xmn<seuil.AND.inf==1) THEN
505C WE STORE IN THE BOTTOM OF THE STACK BP
506 addne = addne + 1
507 pe(addne) = bpe(i)
508 ENDIF
509 xmx = max(x(dir,irect(1,ne)),x(dir,irect(2,ne)),
510 . x(dir,irect(3,ne)),x(dir,irect(4,ne)))
511 . + max(min(bgapsmx+gap_m(ne),gapmax),gapmin)+dgapload+marge_e
512 IF(xmx>=seuil.AND.sup==1) THEN
513C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPE
514 nb_ecn = nb_ecn + 1
515 bpe(nb_ecn) = bpe(i)
516 ENDIF
517 ENDDO
518 END IF
519C
520C 4- REMPLIR LES TABLEAUX D'ADRESSES
521C
522 add(1,2) = addnn
523 add(2,2) = addne
524C-----we fill the mins of the next box and the maxs of the current one
525C (i.e. threshold is a max for the current one)
526C We're going to go down and so we define a new box
527C we fill the maxs of the new box
528C initialises in i7buc1 a 1.E30 comme ca on recupere
529C either XMAX or the max of the box
530 xyzm(1,i_add+1) = xyzm(1,i_add)
531 xyzm(2,i_add+1) = xyzm(2,i_add)
532 xyzm(3,i_add+1) = xyzm(3,i_add)
533 xyzm(4,i_add+1) = xyzm(4,i_add)
534 xyzm(5,i_add+1) = xyzm(5,i_add)
535 xyzm(6,i_add+1) = xyzm(6,i_add)
536 xyzm(dir,i_add+1) = seuil
537 xyzm(dir+3,i_add) = seuil
538C
539 nb_nc = nb_ncn
540 nb_ec = nb_ecn
541C we increment the descent level before exiting
542 i_add = i_add + 1
543 IF(i_add>=1000) THEN
544C TROP NIVEAUX PILE ON VAS LES PRENDRE PLUS GRANDES...
545 IF ( nb_n_b == numnod) THEN
546#ifndef HYPERMESH_LIB
547 CALL ancmsg(msgid=83,
548 . msgtype=msgerror,
549 . anmode=aninfo,
550 . i1=id,
551 . c1=titr)
552#endif
553 ENDIF
554 i_mem = 1
555 RETURN
556 ENDIF
557C
558C this return is only reached when ok = 0
559 RETURN
#define my_real
Definition cppsort.cpp:32
subroutine i24s1s2(prov_n, prov_e, nbinflg, mbinflg, pene)
Definition i24buc1.F:599
subroutine i24fic_getn(ns, irtse, is2se, ie, ns1, ns2)
Definition i24surfi.F:1923
subroutine i7cmp3(i_stok, cand_e, cand_n, iflag, pene, prov_n, prov_e)
Definition i7cmp3.F:82
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
initmumps id
integer, parameter nchartitle
subroutine i7cor3(x, irect, nsv, cand_e, cand_n, stf, stfn, gapv, igap, gap, gap_s, gap_m, istf, gapmin, gapmax, gap_s_l, gap_m_l, drad, ix1, ix2, ix3, ix4, nsvg, x1, x2, x3, x4, y1, y2, y3, y4, z1, z2, z3, z4, xi, yi, zi, stif, dgapload, last)
Definition i7cor3.F:43
subroutine i7dst3(ix3, ix4, x1, x2, x3, x4, y1, y2, y3, y4, z1, z2, z3, z4, xi, yi, zi, x0, y0, z0, nx1, ny1, nz1, nx2, ny2, nz2, nx3, ny3, nz3, nx4, ny4, nz4, p1, p2, p3, p4, lb1, lb2, lb3, lb4, lc1, lc2, lc3, lc4, last)
Definition i7dst3.F:46
subroutine i7dstk(i_add, nb_nc, nb_ec, add, bpn, pn, bpe, pe)
Definition i7dstk.F:33
subroutine i7pen3(marge, gapv, n1, n2, n3, pene, nx1, ny1, nz1, nx2, ny2, nz2, nx3, ny3, nz3, nx4, ny4, nz4, p1, p2, p3, p4, last)
Definition i7pen3.F:43
subroutine ancmsg(msgid, msgtype, anmode, i1, i2, i3, i4, i5, i6, i7, i8, i9, i10, i11, i12, i13, i14, i15, i16, i17, i18, i19, i20, r1, r2, r3, r4, r5, r6, r7, r8, r9, c1, c2, c3, c4, c5, c6, c7, c8, c9, prmode)
Definition message.F:895