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
69C-----------------------------------------------
70C I m p l i c i t T y p e s
71C-----------------------------------------------
72#include "implicit_f.inc"
73C-----------------------------------------------
74C G l o b a l P a r a m e t e r s
75C-----------------------------------------------
76#include "mvsiz_p.inc"
77#include "param_c.inc"
78C-----------------------------------------------
79C C o m m o n B l o c k s
80C-----------------------------------------------
81#include "com04_c.inc"
82#include "vect07_c.inc"
83C-----------------------------------------------
84C D u m m y A r g u m e n t s
85C-----------------------------------------------
86 INTEGER NB_NC,NB_EC,I_ADD,MAXSIZ,I_STOK,J_STOK,I_MEM,ISTF
87 INTEGER I_BID, I_AMAX,NB_N_B, NOINT, NSN,MULTIMP, IGAP
88 INTEGER ADD(2,0:*),IRECT(4,*),BPE(*),PE(*),BPN(*),PN(*)
89 INTEGER NSV(*),CAND_N(*),CAND_E(*), ITAB(*),NBINFLG(*),MBINFLG(*),
90 * ILEV,MVOISN(4,*),IPARTNS(*),IPEN0,INACTI,NRTM
91 INTEGER IXS(NIXS,*), IXS10(6,*), IXS16(8,*), IXS20(12,*),IRTSE(*),IS2SE(*)
92C REAL
94 . x(3,*),xyzm(6,*),tzinf,dbuc,stf(*),stfn(*),
95 . maxbox,minbox, xmax, ymax, zmax,
96 . gap, gap_s(*), gap_m(*),
97 . gapmin, gapmax, marge, gapsmx, bgapsmx,
98 . gap_s_l(*),gap_m_l(*),marge_sh
99 my_real , INTENT(IN) :: dgapload
100 INTEGER ID,MSEGTYP(*)
101 LOGICAL, INTENT(in) :: FLAG_REMOVED_NODE !< flag to remove some S node from the list of candidates
102 INTEGER, INTENT(in) :: S_KREMNODE !< size of KREMNODE
103 INTEGER, INTENT(in) :: S_REMNODE !< size of REMNODE
104 INTEGER, DIMENSION(S_KREMNODE), INTENT(in) :: KREMNODE
105 INTEGER, DIMENSION(S_REMNODE), INTENT(in) :: REMNODE
106 INTEGER, DIMENSION(NUMNOD), INTENT(inout) :: TAG_REMOVED_NODE
107 CHARACTER(LEN=NCHARTITLE) :: TITR
108 INTEGER, DIMENSION(MVSIZ), INTENT(INOUT) ::PROV_N,PROV_E
109 INTEGER, DIMENSION(MVSIZ), INTENT(INOUT) :: IX1,IX2,IX3,IX4,NSVG
110 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: x1,x2,x3,x4
111 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: y1,y2,y3,y4
112 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: z1,z2,z3,z4
113 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: xi,yi,zi
114 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: x0,y0,z0,stif
115 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: n11,n21,n31,pene
116 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: nx1,ny1,nz1
117 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: nx2,ny2,nz2
118 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: nx3,ny3,nz3
119 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: nx4,ny4,nz4
120 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: p1,p2,p3,p4
121 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: lb1,lb2,lb3,lb4
122 my_real, DIMENSION(MVSIZ), INTENT(INOUT) :: lc1,lc2,lc3,lc4
123C-----------------------------------------------
124C L o c a l V a r i a b l e s
125C-----------------------------------------------
126 INTEGER NB_NCN,NB_ECN,ADDNN,ADDNE,IPOS,I,IP,J
127 INTEGER INF,SUP,DIR,N1,N2,N3,N4,NN,NE,INS
128 INTEGER SKIP,NS,NS1,NS2,NSE
129 INTEGER :: FIRST,LAST
130 INTEGER :: IJK
131C REAL
132 my_real
133 . bid,dx,dy,dz,dsup,seuil,xmx,xmn,gapsmax,
134 . gapv(mvsiz),marge_e
135C-----------------------------------------------
136C ROLE DE LA ROUTINE:
137C ===================
138C CLASSE LES ELETS DE BPE ET LES NOEUDS DE BPN EN TWO ZONES
139C > OU < A UNE FRONTIERE ICI DETERMINEE ET SORT LE TOUT
140C DANS bpe,hpe, et bpn,hpn
141C-----------------------------------------------
142C D u m m y A r g u m e n t s
143C
144C NOM DESCRIPTION E/S
145C
146C BPE TABLEAU DES FACETTES A TRIER E/S
147C ET DU RESULTAT COTE MAX
148C PE TABLEAU DES FACETTES S
149C RESULTAT COTE MIN
150C BPN TABLEAU DES NOEUDS A TRIER E/S
151C ET DU RESULTAT COTE MAX
152C PN TABLEAU DES NOEUDS S
153C RESULTAT COTE MIN
154C ADD(2,*) TABLEAU DES ADRESSES E/S
155C 1.......ADRESSES NOEUDS
156C 2.......ADRESSES ELEMENTS
157C ZYZM(6,*) TABLEAU DES 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,*) TABLEAU DES CONEC FACETTES E
165C X(3,*) COORDONNEES NODALES E
166C NB_NC NOMBRE DE NOEUDS CANDIDATS E/S
167C NB_EC NOMBRE D'ELTS CANDIDATS E/S
168C I_ADD POSITION DANS LE TAB DES ADRESSES E/S
169C NSV NOS SYSTEMES DES NOEUDS E
170C XMAX plus grande abcisse existante E
171C XMAX plus grande ordonn. existante E
172C XMAX plus grande cote existante E
173C MAXSIZ TAILLE MEMOIRE MAX POSSIBLE E
174C I_STOK niveau de stockage des couples
175C candidats impact E/S
176C CAND_N boites resultats noeuds
177C CAND_E adresses des boites resultat elements
178C MULTIMP*NSN TAILLE MAX ADMISE MAINTENANT POUR LES
179C COUPLES NOEUDS,ELT CANDIDATS
180C NOINT NUMERO USER DE L'INTERFACE
181C TZINF TAILLE ZONE INFLUENCE
182C MAXBOX TAILLE MAX BUCKET
183C MINBOX TAILLE MIN BUCKET
184C=======================================================================
185C
186C
187C 1- TEST ARRET = BOITE VIDE
188C BOITE TROP PETITE
189C BOITE NE CONTENANT QU'ONE NOEUD
190C PLUS DE MEMOIRE DISPONIBLE
191C
192C-----------------------------------------------------------
193C
194C IF(MEMX>ADD(2,1)+NB_EC)THEN
195C WRITE(ISTDO,*)' *******MEM MAX=',MEMX
196C MEMX=-1
197C ELSEIF(MEMX/=-1)THEN
198C MEMX=ADD(2,1)+NB_EC
199C ENDIF
200C--------------------TEST SUR BOITE VIDES--------------
201C
202 IF(nb_ec==0.OR.nb_nc==0) THEN
203C write(6,*)" BOITE VIDE"
204C IL FAUT COPIER LES BAS DES PILES DANS BAS_DE_PILE CORRESPONDANTS
205C AVANT DE REDESCENDRE DANS LA BRANCHE MITOYENNE
206 CALL i7dstk(i_add,nb_nc,nb_ec,add,bpn,pn,bpe,pe)
207 RETURN
208 ENDIF
209C
210C-------------------TEST SUR FIN DE BRANCHE / MEMOIRE DEPASSEE------------
211C
212 dx = xyzm(4,i_add) - xyzm(1,i_add)
213 dy = xyzm(5,i_add) - xyzm(2,i_add)
214 dz = xyzm(6,i_add) - xyzm(3,i_add)
215 dsup= max(dx,dy,dz)
216C
217 IF(add(2,1)+nb_ec>=maxsiz) THEN
218C PLUS DE PLACE DANS LA PILE DES ELEMENTS BOITES TROP PETITES
219 IF ( nb_n_b == numnod) THEN
220#ifndef HYPERMESH_LIB
221 CALL ancmsg(msgid=83,
222 . msgtype=msgerror,
223 . anmode=aninfo,
224 . i1=id,
225 . c1=titr)
226#endif
227 ENDIF
228 i_mem = 1
229 RETURN
230 ENDIF
231 IF(dsup<minbox.OR.
232 . nb_nc<=nb_n_b.AND.dsup<maxbox.OR.
233 . nb_nc<=nb_n_b.AND.nb_ec==1) THEN
234C
235C write(6,*)" NOUVELLE BOITE "
236C 1- STOCKAGE DU OU DES NOEUD CANDIDAT ET DES ELTS CORRESP.
237C VIRER LES INUTILES
238 DO i=1,nb_ec
239 ne = bpe(i)
240 ! --------------
241 ! do not take into account node_id = remnode(i)
242 IF(flag_removed_node) THEN
243 first = kremnode(ne)+1
244 last = kremnode(ne+1)
245 DO ijk=first,last
246 IF(remnode(ijk)<=numnod) tag_removed_node(remnode(ijk)) = 1
247 ENDDO
248 ENDIF
249 ! --------------
250 n1=irect(1,ne)
251 n2=irect(2,ne)
252 n3=irect(3,ne)
253 n4=irect(4,ne)
254 IF (msegtyp(ne)==0.OR.msegtyp(ne)>nrtm) THEN
255 marge_e = marge
256 ELSE
257 marge_e = marge_sh
258 END IF
259 DO j=1,nb_nc
260 nn=nsv(bpn(j))
261C-----------------
262C-------- not sure if use Inacti----
263C INS = 1
264C IF (INACTI > 0)
265C . CALL INSOLBOX(X ,MVOISN(1,NE),MVOISN(2,NE),NOINT ,IXS,
266C . IXS10 ,IXS16 ,IXS20 ,NN ,GAP ,
267C . MVOISN(3,NE),IPARTNS(BPN(J)),IPEN0,INS )
268 IF(nn/=n1.AND.nn/=n2.AND.nn/=n3.AND.nn/=n4) THEN
269C Fictive Nodes to filter
270 skip=0
271 IF(nn > numnod)THEN
272 ns=nn-numnod
273 CALL i24fic_getn(ns ,irtse ,is2se ,nse ,
274 + ns1 ,ns2 )
275 IF(ns1 == n1 .OR. ns2 == n1) skip=1
276 IF(ns1 == n2 .OR. ns2 == n2) skip=1
277 IF(ns1 == n3 .OR. ns2 == n3) skip=1
278 IF(ns1 == n4 .OR. ns2 == n4) skip=1
279
280 IF(skip==0) THEN
281 IF(flag_removed_node) THEN
282 first = kremnode(ne)+1
283 last = kremnode(ne+1)
284 DO ijk=first,last
285 IF(remnode(ijk)==nn) skip = 1
286 ENDDO
287 ENDIF
288 ENDIF
289 ELSE
290 ! --------------
291 ! S node must be remove from the list of candidates if TAG_REMOVED_NODE = 1
292 IF(flag_removed_node) THEN
293 IF(tag_removed_node(nn)==1) skip = 1
294 ENDIF
295 ! --------------
296 ENDIF
297
298 IF (skip==0)THEN
299 j_stok = j_stok + 1
300 prov_n(j_stok) = bpn(j)
301 prov_e(j_stok) = ne
302 IF(j_stok==nvsiz) THEN
303 lft = 1
304 llt = nvsiz
305 nft = 0
306 j_stok = 0
307 CALL i7cor3(x ,irect,nsv ,prov_e ,prov_n,
308 . stf ,stfn ,gapv ,igap ,gap ,
309 . gap_s,gap_m,istf ,gapmin ,gapmax,
310 . gap_s_l,gap_m_l ,zero ,ix1 ,ix2 ,
311 5 ix3 ,ix4 ,nsvg,x1 ,x2 ,
312 6 x3 ,x4 ,y1 ,y2 ,y3 ,
313 7 y4 ,z1 ,z2 ,z3 ,z4 ,
314 8 xi ,yi ,zi ,stif ,dgapload,
315 9 llt)
316 CALL i7dst3(ix3,ix4,x1 ,x2 ,x3 ,
317 1 x4 ,y1 ,y2 ,y3 ,y4 ,
318 2 z1 ,z2 ,z3 ,z4 ,xi ,
319 3 yi ,zi ,x0 ,y0 ,z0 ,
320 4 nx1,ny1,nz1,nx2,ny2,
321 5 nz2,nx3,ny3,nz3,nx4,
322 6 ny4,nz4,p1 ,p2 ,p3 ,
323 7 p4 ,lb1,lb2,lb3,lb4,
324 8 lc1,lc2,lc3,lc4,llt)
325
326 CALL i7pen3(marge_e,gapv,n11,n21,n31,
327 1 pene ,nx1 ,ny1,nz1,nx2,
328 2 ny2 ,nz2 ,nx3,ny3,nz3,
329 3 nx4 ,ny4 ,nz4,p1 ,p2 ,
330 4 p3 ,p4,llt)
331
332 IF (ilev==2) CALL i24s1s2(prov_n,prov_e,nbinflg,mbinflg,pene)
333 IF(i_stok+nvsiz<multimp*nsn) THEN
334 CALL i7cmp3(i_stok,cand_e ,cand_n,1,pene,
335 1 prov_n,prov_e)
336 ELSE
337 i_bid = 0
338 CALL i7cmp3(i_bid,cand_e,cand_n,0,pene,
339 1 prov_n,prov_e)
340 IF(i_stok+i_bid<multimp*nsn) THEN
341 CALL i7cmp3(i_stok,cand_e,cand_n,1,pene,
342 1 prov_n,prov_e)
343 ELSE
344 i_mem = 2
345c CALL ANSTCKI(NOINT)
346c CALL ANCWARN(103,ANINFO_BLIND_2)
347 RETURN
348 ENDIF
349 ENDIF
350 ENDIF
351C write(6,*)"Noeud candidat",BPN(J)
352C write(6,*)"Element candidat",NE
353 ENDIF
354 ENDIF
355 ENDDO
356
357 ! --------------
358 ! flush to 0
359 IF(flag_removed_node) THEN
360 first = kremnode(ne)+1
361 last = kremnode(ne+1)
362 DO ijk=first,last
363 IF(remnode(ijk)<=numnod) tag_removed_node(remnode(ijk)) = 0
364 ENDDO
365 ENDIF
366 ! --------------
367 ENDDO
368C IL FAUT COPIER LES BAS DES PILES DANS BAS_DE_PILE CORRESPONDANTS
369C AVANT DE REDESCENDRE DANS LA BRANCHE MITOYENNE
370 CALL i7dstk(i_add,nb_nc,nb_ec,add,bpn,pn,bpe,pe)
371 RETURN
372 ENDIF
373C
374C-----------------------------------------------------------
375C
376C
377C 2- PHASE DE TRI SUR LA MEDIANE SELON LA + GDE DIRECTION
378C
379C
380C-----------------------------------------------------------
381C
382C
383C 1- DETERMINER LA DIRECTION A DIVISER X,Y OU Z
384C
385 dir = 1
386 IF(dy==dsup) THEN
387 dir = 2
388 ELSE IF(dz==dsup) THEN
389 dir = 3
390 ENDIF
391 seuil =(xyzm(dir+3,i_add)+xyzm(dir,i_add))/2
392C
393C 2- DIVISER LES NOEUDS EN TWO ZONES
394C
395 nb_ncn= 0
396 addnn= add(1,1)
397 inf = 0
398 sup = 0
399 IF(igap==0)THEN
400 DO i=1,nb_nc
401 IF(x(dir,nsv(bpn(i)))<seuil) THEN
402C ON STOCKE DANS LE BAS DE LA PILE BP
403 addnn = addnn + 1
404 pn(addnn) = bpn(i)
405 inf = 1
406 ELSE
407 nb_ncn = nb_ncn + 1
408 bpn(nb_ncn) = bpn(i)
409C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPN
410 sup = 1
411 ENDIF
412 END DO
413 ELSE
414 gapsmx = zero
415 bgapsmx = zero
416 DO i=1,nb_nc
417 IF(x(dir,nsv(bpn(i)))<seuil) THEN
418C ON STOCKE DANS LE BAS DE LA PILE BP
419 addnn = addnn + 1
420 pn(addnn) = bpn(i)
421 gapsmx = max(gapsmx,gap_s(bpn(i)))
422 inf = 1
423 ELSE
424C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPN
425 nb_ncn = nb_ncn + 1
426 bpn(nb_ncn) = bpn(i)
427 bgapsmx = max(bgapsmx,gap_s(bpn(i)))
428 sup = 1
429 ENDIF
430 END DO
431 END IF
432
433CC
434CC 3- DIVISER LES ELEMENTS
435CC
436C NB_ECN= 0
437C ADDNE= ADD(2,1)
438C SEUILI = SEUIL-TZINF
439C SEUILS = SEUIL+TZINF
440C DO 85 I=1,NB_EC
441C INF = 0
442C SUP = 0
443C DO 80 J=1,4
444C IP = IRECT(J,BPE(I))
445C IF(X(DIR,IP)<SEUILS) THEN
446C INF = 1
447C IF(SUP==1) GOTO 81
448C ENDIF
449C IF(X(DIR,IP)>=SEUILI) THEN
450C SUP = 1
451C IF(INF==1) GOTO 81
452C ENDIF
453C 80 CONTINUE
454C
455C 81 CONTINUE
456C IF(INF==1) THEN
457C ON STOCKE DANS LE BAS DE LA PILE BP
458C ADDNE = ADDNE + 1
459C PE(ADDNE) = BPE(I)
460C ENDIF
461C IF(SUP==1) THEN
462C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPE
463C NB_ECN = NB_ECN + 1
464C BPE(NB_ECN) = BPE(I)
465C ENDIF
466C 85 CONTINUE
467C
468C 3- DIVISER LES ELEMENTS
469C
470C 2 LIGNES PROV POUR TEST
471C INF = 1
472C SUP = 1
473C
474 nb_ecn= 0
475 addne= add(2,1)
476 IF(igap==0)THEN
477 DO i=1,nb_ec
478 xmx = max(x(dir,irect(1,bpe(i))),x(dir,irect(2,bpe(i))),
479 . x(dir,irect(3,bpe(i))),x(dir,irect(4,bpe(i))))
480 . + tzinf
481 xmn = min(x(dir,irect(1,bpe(i))),x(dir,irect(2,bpe(i))),
482 . x(dir,irect(3,bpe(i))),x(dir,irect(4,bpe(i))))
483 . - tzinf
484 IF(xmn<seuil.AND.inf==1) THEN
485C ON STOCKE DANS LE BAS DE LA PILE BP
486 addne = addne + 1
487 pe(addne) = bpe(i)
488 ENDIF
489 IF(xmx>=seuil.AND.sup==1) THEN
490C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPE
491 nb_ecn = nb_ecn + 1
492 bpe(nb_ecn) = bpe(i)
493 ENDIF
494 ENDDO
495 ELSE
496 DO i=1,nb_ec
497 ne = bpe(i)
498 IF (msegtyp(ne)==0.OR.msegtyp(ne)>nrtm) THEN
499 marge_e = marge
500 ELSE
501 marge_e = marge_sh
502 END IF
503 xmn = min(x(dir,irect(1,ne)),x(dir,irect(2,ne)),
504 . x(dir,irect(3,ne)),x(dir,irect(4,ne)))
505 . - max(min(gapsmx+gap_m(ne),gapmax),gapmin)+dgapload-marge_e
506 IF(xmn<seuil.AND.inf==1) THEN
507C ON STOCKE DANS LE BAS DE LA PILE BP
508 addne = addne + 1
509 pe(addne) = bpe(i)
510 ENDIF
511 xmx = max(x(dir,irect(1,ne)),x(dir,irect(2,ne)),
512 . x(dir,irect(3,ne)),x(dir,irect(4,ne)))
513 . + max(min(bgapsmx+gap_m(ne),gapmax),gapmin)+dgapload+marge_e
514 IF(xmx>=seuil.AND.sup==1) THEN
515C ON STOCKE EN ECRASANT PROGRESSIVEMENT BPE
516 nb_ecn = nb_ecn + 1
517 bpe(nb_ecn) = bpe(i)
518 ENDIF
519 ENDDO
520 END IF
521C
522C 4- REMPLIR LES TABLEAUX D'ADRESSES
523C
524 add(1,2) = addnn
525 add(2,2) = addne
526C-----on remplit les min de la boite suivante et les max de la courante
527C (i.e. seuil est un max pour la courante)
528C on va redescendre et donc on definit une nouvelle boite
529C on remplit les max de la nouvelle boite
530C initialises dans i7buc1 a 1.E30 comme ca on recupere
531C soit XMAX soit le max de la boite
532 xyzm(1,i_add+1) = xyzm(1,i_add)
533 xyzm(2,i_add+1) = xyzm(2,i_add)
534 xyzm(3,i_add+1) = xyzm(3,i_add)
535 xyzm(4,i_add+1) = xyzm(4,i_add)
536 xyzm(5,i_add+1) = xyzm(5,i_add)
537 xyzm(6,i_add+1) = xyzm(6,i_add)
538 xyzm(dir,i_add+1) = seuil
539 xyzm(dir+3,i_add) = seuil
540C
541 nb_nc = nb_ncn
542 nb_ec = nb_ecn
543C on incremente le niveau de descente avant de sortir
544 i_add = i_add + 1
545 IF(i_add>=1000) THEN
546C TROP NIVEAUX PILE ON VAS LES PRENDRE PLUS GRANDES...
547 IF ( nb_n_b == numnod) THEN
548#ifndef HYPERMESH_LIB
549 CALL ancmsg(msgid=83,
550 . msgtype=msgerror,
551 . anmode=aninfo,
552 . i1=id,
553 . c1=titr)
554#endif
555 ENDIF
556 i_mem = 1
557 RETURN
558 ENDIF
559C
560C ce return n'est atteint que dans le cas ok = 0
561 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:1921
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:272
#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:34
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:889