OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i24trivox.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!|| i24trivox ../engine/source/interfaces/intsort/i24trivox.F
25!||--- called by ------------------------------------------------------
26!|| i24buce ../engine/source/interfaces/intsort/i24buce.F
27!||--- calls -----------------------------------------------------
28!|| i24fic_getn ../engine/source/interfaces/int24/i24for3e.F
29!|| i24sto ../engine/source/interfaces/intsort/i24sto.F
30!|| my_barrier ../engine/source/system/machine.F
31!|| spmd_oldnumcd ../engine/source/mpi/interfaces/spmd_i7tool.F
32!||--- uses -----------------------------------------------------
33!|| tri7box ../engine/share/modules/tri7box.F
34!||====================================================================
35 SUBROUTINE i24trivox(
36 1 NSN ,NSNR ,ISZNSNR ,I_MEM ,VMAXDT ,
37 2 IRECT ,X ,STF ,STFN ,XYZM ,
38 3 NSV ,II_STOK ,CAND_N ,ESHIFT ,CAND_E ,
39 4 MULNSN ,NOINT ,V ,BGAPSMX ,
40 5 VOXEL ,NBX ,NBY ,NBZ ,PMAX_GAP ,
41 6 NRTM ,GAP_S ,GAP_M ,MARGE ,CURV_MAX ,
42 7 NIN ,ITASK ,PENE_OLD,ITAB ,NBINFLG ,
43 8 MBINFLG,ILEV ,MSEGTYP ,EDGE_L2 ,IEDGE ,
44 9 ISEADD ,ISEDGE ,CAND_T ,FLAGREMNODE,KREMNOD,
45 A REMNOD ,CAND_A ,RENUM ,NSNROLD ,IRTSE ,
46 B IS2SE ,NSNE ,DGAPLOAD,INTHEAT,IDT_THERM,NODADT_THERM)
47C============================================================================
48C M o d u l e s
49C-----------------------------------------------
50 USE tri7box
51C-----------------------------------------------
52C I m p l i c i t T y p e s
53C-----------------------------------------------
54#include "implicit_f.inc"
55C-----------------------------------------------
56C G l o b a l P a r a m e t e r s
57C-----------------------------------------------
58#include "mvsiz_p.inc"
59c parameter setting the size for the vector (orig version is 128)
60 INTEGER NVECSZ
61 parameter(nvecsz = mvsiz)
62C-----------------------------------------------
63C C o m m o n B l o c k s
64C-----------------------------------------------
65#include "com01_c.inc"
66#include "com04_c.inc"
67#include "param_c.inc"
68#include "task_c.inc"
69C-----------------------------------------------
70C ROLE DE LA ROUTINE:
71C ===================
72C CLASSE LES NOEUDS DANS DES BOITES
73C RECHERCHE POUR CHAQUE FACETTE DES BOITES CONCERNES
74C RECHERCHE DES CANDIDATS
75C-----------------------------------------------
76C D u m m y A r g u m e n t s
77C
78C NOM DESCRIPTION E/S
79C
80C IRECT(4,*) TABLEAU DES CONEC FACETTES E
81C X(3,*) COORDONNEES NODALES E
82C NSV NOS SYSTEMES DES NOEUDS E
83C XMAX plus grande abcisse existante E
84C XMAX plus grande ordonn. existante E
85C XMAX plus grande cote existante E
86C I_STOK niveau de stockage des couples
87C candidats impact E/S
88C CAND_N boites resultats noeuds
89C CAND_E adresses des boites resultat elements
90C MULNSN = MULTIMP*NSN TAILLE MAX ADMISE MAINTENANT POUR LES
91C COUPLES NOEUDS,ELT CANDIDATS
92C NOINT NUMERO USER DE L'INTERFACE
93C
94C PROV_N CAND_N provisoire (variable static dans i7tri)
95C PROV_E CAND_E provisoire (variable static dans i7tri)
96
97C VOXEL(ix,iy,iz) contient le numero local du premier noeud de
98C la boite
99C NEXT_NOD(i) noeud suivant dans la meme boite (si /= 0)
100C LAST_NOD(i) dernier noeud dans la meme boite (si /= 0)
101C utilise uniquement pour aller directement du premier
102C noeud au dernier
103C-----------------------------------------------
104C D u m m y A r g u m e n t s
105C-----------------------------------------------
106 INTEGER I_MEM,ESHIFT,NSN,ISZNSNR,NRTM,NIN,ITASK,
107 . MULNSN,NOINT,NSNR,NBX,NBY,NBZ,IEDGE,NSNE,
108 . NSV(*),CAND_N(*),CAND_E(*),
109 . IRECT(4,*), VOXEL(NBX+2,NBY+2,NBZ+2),II_STOK,ITAB(*),
110 . NBINFLG(*),MBINFLG(*),ILEV,MSEGTYP(*),CAND_T(*),
111 . ISEADD(*) ,ISEDGE(*),FLAGREMNODE,KREMNOD(*),REMNOD(*),CAND_A(*),
112 . RENUM(*),NSNROLD,IRTSE(5,*),IS2SE(2,*)
113 INTEGER, INTENT(IN) :: INTHEAT
114 INTEGER, INTENT(IN) :: IDT_THERM
115 INTEGER, INTENT(IN) :: NODADT_THERM
116C REAL
117 my_real
118 . x(3,*),v(3,*),xyzm(6),stf(*),stfn(*),gap_s(*),
119 . gap_m(*),curv_max(*),pene_old(5,nsn),edge_l2(*),
120 . marge,bgapsmx,pmax_gap,vmaxdt
121 my_real , INTENT(IN) :: dgapload
122C-----------------------------------------------
123C L o c a l V a r i a b l e s
124C-----------------------------------------------
125 INTEGER NB_NCN,NB_NCN1,NB_ECN,I,J,DIR,NB_NC,NB_EC,
126 . n1,n2,n3,n4,nn,ne,k,l,j_stok,ii,jj,
127 . prov_n(mvsiz),prov_e(mvsiz),
128 . oldnum(isznsnr), nsnf, nsnl,m,nse,ns,ip
129C REAL
130 my_real
131 . dx,dy,dz,xs,ys,zs,xx,sx,sy,sz,s2,
132 . xmin, xmax,ymin, ymax,zmin, zmax, tz, gapsmx, gapl,
133 . xx1,xx2,xx3,xx4,yy1,yy2,yy3,yy4,zz1,zz2,zz3,zz4,
134 . d1x,d1y,d1z,d2x,d2y,d2z,dd1,dd2,d2,a2,gs
135c provisoire
136 INTEGER LAST_NOD(NSN+NSNR)
137 INTEGER IX,IY,IZ,NEXT,M1,M2,M3,M4,
138 . IX1,IY1,IZ1,IX2,IY2,IZ2
139 INTEGER, DIMENSION(:),ALLOCATABLE :: IIX,IIY,IIZ
140 my_real
141 . xminb,yminb,zminb,xmaxb,ymaxb,zmaxb,
142 . xmine,ymine,zmine,xmaxe,ymaxe,zmaxe,aaa
143 INTEGER FIRST,NEW,LAST
144 SAVE IIX,IIY,IIZ
145 INTEGER, DIMENSION(NUMNOD+NSNE) :: TAG
146C --------------------------------
147C TYPE24 E2E - I24FIC_GETN method
148C --------------------------------
149 INTEGER IK1(4),IK2(4),IE1,IE2,IED,NS1,NS2,NS1ID,NS2ID
150 DATA ik1 /1,2,3,4/
151 DATA ik2 /2,3,4,1/
152C-----------------------------------------------
153 IF(itask == 0)THEN
154 ALLOCATE(next_nod(nsn+nsnr))
155 ALLOCATE(iix(nsn+nsnr))
156 ALLOCATE(iiy(nsn+nsnr))
157 ALLOCATE(iiz(nsn+nsnr))
158 END IF
159C Barrier to wait init voxel and allocation NEX_NOD
160 CALL my_barrier
161C Phase initiale de construction de BPE et BPN deplacee de I7BUCE => I7TRI
162C
163 xmin = xyzm(1)
164 ymin = xyzm(2)
165 zmin = xyzm(3)
166 xmax = xyzm(4)
167 ymax = xyzm(5)
168 zmax = xyzm(6)
169
170c dev future: xminb plus grand que xmin...
171 xminb = xmin
172 yminb = ymin
173 zminb = zmin
174 xmaxb = xmax
175 ymaxb = ymax
176 zmaxb = zmax
177
178c!!!!!!!!!!!!!!! A VERIFIER !!!!!!!!!!!!!!!
179C En SPMD, pour IFQ, retrouve ancienne numerotation des candidats non locaux
180c INUTIL POUR INT 24 !!!!!!!!!!!!!!!!!
181 IF(nspmd>1) THEN
182 CALL spmd_oldnumcd(renum,oldnum,isznsnr,nsnrold,intheat,idt_therm,nodadt_therm)
183 END IF
184
185C=======================================================================
186C 1 mise des noeuds dans les boites
187C=======================================================================
188C Note for Edge2Edge : X is no more the Radioss X Array but an extension
189C NUMNOD+SNE
190C It is updated at any cycle
191 IF(itask == 0)THEN
192 DO i=1,nsn
193 iix(i)=0
194 iiy(i)=0
195 iiz(i)=0
196 IF(stfn(i) == zero)cycle
197 j=nsv(i)
198C Optimisation // recherche les noeuds compris dans xmin xmax des
199C elements du processeur
200 IF(x(1,j) < xmin) cycle
201 IF(x(1,j) > xmax) cycle
202 IF(x(2,j) < ymin) cycle
203 IF(x(2,j) > ymax) cycle
204 IF(x(3,j) < zmin) cycle
205 IF(x(3,j) > zmax) cycle
206
207 iix(i)=int(nbx*(x(1,j)-xminb)/(xmaxb-xminb))
208 iiy(i)=int(nby*(x(2,j)-yminb)/(ymaxb-yminb))
209 iiz(i)=int(nbz*(x(3,j)-zminb)/(zmaxb-zminb))
210
211 iix(i)=max(1,2+min(nbx,iix(i)))
212 iiy(i)=max(1,2+min(nby,iiy(i)))
213 iiz(i)=max(1,2+min(nbz,iiz(i)))
214
215 first = voxel(iix(i),iiy(i),iiz(i))
216 IF(first == 0)THEN
217c empty cell
218 voxel(iix(i),iiy(i),iiz(i)) = i ! first
219 next_nod(i) = 0 ! last one
220 last_nod(i) = 0 ! no last
221 ELSEIF(last_nod(first) == 0)THEN
222c cell containing one node
223c add as next node
224 next_nod(first) = i ! next
225 last_nod(first) = i ! last
226 next_nod(i) = 0 ! last one
227 ELSE
228c
229c jump to the last node of the cell
230 last = last_nod(first) ! last node in this voxel
231 next_nod(last) = i ! next
232 last_nod(first) = i ! last
233 next_nod(i) = 0 ! last one
234 ENDIF
235 ENDDO
236C=======================================================================
237C 2 mise des noeuds dans les boites
238C candidats non locaux en SPMD
239C=======================================================================
240 DO j = 1, nsnr
241
242 IF(irem(8,j)==-1) cycle ! case IREM / ISEDGE_FI==-1 : Node was added due to Fictive Remote Node only
243 ! Do not retain in sorting, otherwise node can be candidate twice
244
245 iix(nsn+j)=int(nbx*(xrem(1,j)-xminb)/(xmaxb-xminb))
246 iiy(nsn+j)=int(nby*(xrem(2,j)-yminb)/(ymaxb-yminb))
247 iiz(nsn+j)=int(nbz*(xrem(3,j)-zminb)/(zmaxb-zminb))
248 iix(nsn+j)=max(1,2+min(nbx,iix(nsn+j)))
249 iiy(nsn+j)=max(1,2+min(nby,iiy(nsn+j)))
250 iiz(nsn+j)=max(1,2+min(nbz,iiz(nsn+j)))
251
252 first = voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))
253 IF(first == 0)THEN
254c empty cell
255 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j)) = nsn+j ! first
256 next_nod(nsn+j) = 0 ! last one
257 last_nod(nsn+j) = 0 ! no last
258 ELSEIF(last_nod(first) == 0)THEN
259c cell containing one node
260c add as next node
261 next_nod(first) = nsn+j ! next
262 last_nod(first) = nsn+j ! last
263 next_nod(nsn+j) = 0 ! last one
264 ELSE
265c
266c jump to the last node of the cell
267 last = last_nod(first) ! last node in this voxel
268 next_nod(last) = nsn+j ! next
269 last_nod(first) = nsn+j ! last
270 next_nod(nsn+j) = 0 ! last one
271 ENDIF
272 ENDDO
273 END IF
274C Barrier to wait task0 treatment
275 CALL my_barrier
276C=======================================================================
277C 3 recherche des boites concernant chaque facette
278C et creation des candidats
279C=======================================================================
280 j_stok = 0
281 IF(flagremnode == 2)THEN
282 DO i=1,numnod+nsne
283 tag(i) = 0
284 ENDDO
285 END IF
286
287 DO ne=1,nrtm
288C on ne retient pas les facettes detruites
289 IF(stf(ne) == zero)cycle
290
291 aaa = marge+curv_max(ne)+bgapsmx+pmax_gap+vmaxdt
292 + + gap_m(ne)+dgapload
293
294
295c il est possible d'ameliorer l'algo en decoupant la facette
296c en 2(4,3,6,9...) si la facette est grande devant AAA et inclinee
297
298 m1 = irect(1,ne)
299 m2 = irect(2,ne)
300 m3 = irect(3,ne)
301 m4 = irect(4,ne)
302
303 xx1=x(1,m1)
304 xx2=x(1,m2)
305 xx3=x(1,m3)
306 xx4=x(1,m4)
307 xmaxe=max(xx1,xx2,xx3,xx4)
308 xmine=min(xx1,xx2,xx3,xx4)
309
310 yy1=x(2,m1)
311 yy2=x(2,m2)
312 yy3=x(2,m3)
313 yy4=x(2,m4)
314 ymaxe=max(yy1,yy2,yy3,yy4)
315 ymine=min(yy1,yy2,yy3,yy4)
316
317 zz1=x(3,m1)
318 zz2=x(3,m2)
319 zz3=x(3,m3)
320 zz4=x(3,m4)
321 zmaxe=max(zz1,zz2,zz3,zz4)
322 zmine=min(zz1,zz2,zz3,zz4)
323
324
325c calcul de la surface (pour elimination future de candidats)
326
327 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
328 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
329 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
330 s2 = sx*sx + sy*sy + sz*sz
331
332c indice des voxels occupes par la facette
333
334 ix1=int(nbx*(xmine-aaa-xminb)/(xmaxb-xminb))
335 iy1=int(nby*(ymine-aaa-yminb)/(ymaxb-yminb))
336 iz1=int(nbz*(zmine-aaa-zminb)/(zmaxb-zminb))
337
338 ix1=max(1,2+min(nbx,ix1))
339 iy1=max(1,2+min(nby,iy1))
340 iz1=max(1,2+min(nbz,iz1))
341
342 ix2=int(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb))
343 iy2=int(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb))
344 iz2=int(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb))
345
346 ix2=max(1,2+min(nbx,ix2))
347 iy2=max(1,2+min(nby,iy2))
348 iz2=max(1,2+min(nbz,iz2))
349
350 IF(flagremnode == 2)THEN
351 k = kremnod(2*(ne-1)+1)+1
352 l = kremnod(2*(ne-1)+2)
353 DO i=k,l
354 tag(remnod(i)) = 1
355 ENDDO
356 END if!(FLAGREMNODE == 2)THEN
357cc nbpelem = 0
358cc nnpelem = 0
359cc nnr0pelem = 0
360cc nnrpelem = 0
361
362 DO iz = iz1,iz2
363 DO iy = iy1,iy2
364 DO ix = ix1,ix2
365
366cc nbpelem = nbpelem + 1
367
368 jj = voxel(ix,iy,iz)
369
370 DO WHILE(jj /= 0)
371
372cc nnpelem = nnpelem + 1
373
374 IF(jj<=nsn)THEN
375 nn=nsv(jj)
376 IF(nn == m1)GOTO 200
377 IF(nn == m2)GOTO 200
378 IF(nn == m3)GOTO 200
379 IF(nn == m4)GOTO 200
380 IF(flagremnode == 2)THEN
381 IF(tag(nn) == 1)GOTO 200
382 END IF
383C-----Fictive nodes on edges-: exclude auto-impact------
384 IF (nn >numnod) THEN
385 ns = nn-numnod
386 CALL i24fic_getn(ns ,irtse ,is2se ,nse ,
387 + ns1 ,ns2 )
388 IF(ns1 == m1 .OR. ns2 == m1) GOTO 200
389 IF(ns1 == m2 .OR. ns2 == m2) GOTO 200
390 IF(ns1 == m3 .OR. ns2 == m3) GOTO 200
391 IF(ns1 == m4 .OR. ns2 == m4) GOTO 200
392 END IF
393 xs = x(1,nn)
394 ys = x(2,nn)
395 zs = x(3,nn)
396c PMAX_GAP is a global overestimate penetration
397c NEED to communicate in SPMD
398c VMAXDT is a local overestimate of relative incremental displacement
399c NO need to communicate in SPMD
400
401 IF (iedge > 0) THEN
402 aaa = marge + curv_max(ne)
403 + + max(gap_s(jj)+gap_m(ne)+edge_l2(jj)+dgapload
404 + ,pene_old(3,jj))+vmaxdt
405 ELSE
406 aaa = marge + curv_max(ne)
407 + + max(gap_s(jj)+gap_m(ne)+dgapload
408 + ,pene_old(3,jj))+vmaxdt
409 END IF
410 ELSE
411 j=jj-nsn
412 IF(flagremnode == 2)THEN
413 k = kremnod(2*(ne-1)+2) + 1
414 l = kremnod(2*(ne-1)+3)
415 IF(irem(8,j)==1) THEN
416 DO m=k,l
417 IF(remnod(m) == -irem(2,j) ) GOTO 200
418 ENDDO
419 ELSE
420 DO m=k,l
421 IF(remnod(m) == -irem(2,j) ) GOTO 200
422 ENDDO
423 ENDIF
424 END if!(FLAGREMNODE == 2)THEN
425C
426C Auto impact between main surface and secnd Edge
427C can happen with Remote nodes when Secnd Nodes are on border between 2 domains
428C and Fictive node is remote
429 IF(irem(8,j)==1) THEN
430 ! Same than in I24FIC_GETN but for Remote Node
431 i24irempnsne=irem(7,j) ! in IREM IRTSE is located in IREM(I24IREMPNSNE,J) to IREM(I24IREMPNSNE+4,J)
432 ied = irem(i24irempnsne+4,j) ! IED = IRTSE(5,xx)
433 ns1 = irem(i24irempnsne-1+ik1(ied),j) ! NS1 = IRTSE(IK1(IED))
434 ns2 = irem(i24irempnsne-1+ik2(ied),j) ! NS2 = IRTSE(IK2(IED))
435 ns1id = irem(2,ns1) ! ITAB Remote NS1
436 ns2id = irem(2,ns2) ! ITAB Remote NS2
437 IF (ns1id == itab(m1) .OR. ns2id == itab(m1)) GOTO 200
438 IF (ns1id == itab(m2) .OR. ns2id == itab(m2)) GOTO 200
439 IF (ns1id == itab(m3) .OR. ns2id == itab(m3)) GOTO 200
440 IF (ns1id == itab(m4) .OR. ns2id == itab(m4)) GOTO 200
441 ENDIF
442 xs = xrem(1,j)
443 ys = xrem(2,j)
444 zs = xrem(3,j)
445 aaa = marge+curv_max(ne)
446c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
447c +EDGE_L2(JJ) remote
448c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
449 + + max(xrem(igapxremp,j)+gap_m(ne)+dgapload,xrem(i24xremp+6,j))
450 + + vmaxdt
451 ENDIF
452
453 IF(xs<=xmine-aaa)GOTO 200
454 IF(xs>=xmaxe+aaa)GOTO 200
455 IF(ys<=ymine-aaa)GOTO 200
456 IF(ys>=ymaxe+aaa)GOTO 200
457 IF(zs<=zmine-aaa)GOTO 200
458 IF(zs>=zmaxe+aaa)GOTO 200
459
460c sousestimation de la distance**2 pour elimination de candidats
461
462cc nnr0pelem = nnr0pelem + 1
463
464 d1x = xs - xx1
465 d1y = ys - yy1
466 d1z = zs - zz1
467 d2x = xs - xx2
468 d2y = ys - yy2
469 d2z = zs - zz2
470 dd1 = d1x*sx+d1y*sy+d1z*sz
471 dd2 = d2x*sx+d2y*sy+d2z*sz
472 IF(dd1*dd2 > zero)THEN
473 d2 = min(dd1*dd1,dd2*dd2)
474 a2 = aaa*aaa*s2
475 IF(d2 > a2)GOTO 200
476 ENDIF
477
478cc nnrpelem = nnrpelem + 1
479
480 j_stok = j_stok + 1
481 prov_n(j_stok) = jj
482 prov_e(j_stok) = ne
483 IF(j_stok == nvsiz)THEN
484
485 CALL i24sto(
486 1 nvsiz ,irect ,x ,nsv ,ii_stok,
487 2 cand_n,cand_e ,mulnsn,noint ,marge ,
488 3 i_mem ,prov_n ,prov_e,eshift,v ,
489 4 nsn ,gap_s ,gap_m ,curv_max,nin ,
490 5 pene_old,nbinflg ,mbinflg,ilev,msegtyp,
491 6 edge_l2,iedge,iseadd ,isedge ,cand_t,itab,
492 7 cand_a,oldnum,nsnrold,dgapload)
493 IF(i_mem==2)GOTO 100
494 j_stok = 0
495 ENDIF
496
497 200 CONTINUE
498
499 jj = next_nod(jj)
500
501 ENDDO ! WHILE(JJ /= 0)
502
503 ENDDO
504 ENDDO
505 ENDDO
506cc nbpelg = nbpelg + nbpelem
507cc nnpelg = nnpelg + nnpelem
508cc nnrpelg = nnrpelg + nnrpelem
509cc nnr0pelg = nnr0pelg + nnr0pelem
510 IF(flagremnode == 2)THEN
511 k = kremnod(2*(ne-1)+1)+1
512 l = kremnod(2*(ne-1)+2)
513 DO i=k,l
514 tag(remnod(i)) = 0
515 ENDDO
516 END IF
517 ENDDO
518
519C-------------------------------------------------------------------------
520C FIN DU TRI
521C-------------------------------------------------------------------------
522 IF(j_stok/=0)CALL i24sto(
523 1 j_stok,irect ,x ,nsv ,ii_stok,
524 2 cand_n,cand_e ,mulnsn,noint ,marge ,
525 3 i_mem ,prov_n ,prov_e,eshift,v ,
526 4 nsn ,gap_s ,gap_m ,curv_max,nin ,
527 5 pene_old,nbinflg,mbinflg,ilev ,msegtyp,
528 6 edge_l2,iedge,iseadd ,isedge ,cand_t,itab,
529 7 cand_a,oldnum,nsnrold,dgapload)
530
531C=======================================================================
532C 4 remise a zero des noeuds dans les boites
533C=======================================================================
534 100 CONTINUE
535
536C Barrier to avoid reinitialization before end of sorting
537 CALL my_barrier
538 nsnf = 1 + itask*nsn / nthread
539 nsnl = (itask+1)*nsn / nthread
540
541 DO i=nsnf,nsnl
542 IF(iix(i)/=0)THEN
543 voxel(iix(i),iiy(i),iiz(i))=0
544 ENDIF
545 ENDDO
546C=======================================================================
547C 5 remise a zero des noeuds dans les boites
548C candidats non locaux en SPMD
549C=======================================================================
550 nsnf = 1 + itask*nsnr / nthread
551 nsnl = (itask+1)*nsnr / nthread
552 DO j = nsnf, nsnl
553 IF(irem(8,j)==-1)cycle
554 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))=0
555 ENDDO
556
557C
558 CALL my_barrier()
559 IF(itask == 0)THEN
560 DEALLOCATE(next_nod)
561 DEALLOCATE(iix)
562 DEALLOCATE(iiy)
563 DEALLOCATE(iiz)
564 ENDIF
565
566 RETURN
567 END
568
#define my_real
Definition cppsort.cpp:32
if(complex_arithmetic) id
subroutine i24sto(j_stok, irect, x, nsv, ii_stok, cand_n, cand_e, mulnsn, noint, marge, i_mem, prov_n, prov_e, eshift, v, nsn, gap_s, gap_m, curv_max, nin, pene_old, nbinflg, mbinflg, ilev, msegtyp, edge_l2, iedge, iseadd, isedge, cand_t, itab, cand_a, oldnum, nsnrold, dgapload)
Definition i24sto.F:43
subroutine i24fic_getn(ns, irtse, is2se, ie, ns1, ns2)
Definition i24surfi.F:1921
subroutine i24trivox(nsn, nsnr, isznsnr, i_mem, vmaxdt, irect, x, stf, stfn, xyzm, nsv, ii_stok, cand_n, eshift, cand_e, mulnsn, noint, v, bgapsmx, voxel, nbx, nby, nbz, pmax_gap, nrtm, gap_s, gap_m, marge, curv_max, nin, itask, pene_old, itab, nbinflg, mbinflg, ilev, msegtyp, edge_l2, iedge, iseadd, isedge, cand_t, flagremnode, kremnod, remnod, cand_a, renum, nsnrold, irtse, is2se, nsne, dgapload, intheat, idt_therm, nodadt_therm)
Definition i24trivox.F:47
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
integer, dimension(:), allocatable next_nod
Definition tri7box.F:48
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 my_barrier
Definition machine.F:31