OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i25trivox.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!|| i25trivox ../engine/source/interfaces/intsort/i25trivox.F
25!||--- called by ------------------------------------------------------
26!|| i25buce ../engine/source/interfaces/intsort/i25buce.F
27!||--- calls -----------------------------------------------------
28!|| i25sto ../engine/source/interfaces/intsort/i25sto.F
29!|| my_barrier ../engine/source/system/machine.F
30!||--- uses -----------------------------------------------------
31!|| my_alloc_mod ../common_source/tools/memory/my_alloc.F90
32!|| tri7box ../engine/share/modules/tri7box.F
33!||====================================================================
34 SUBROUTINE i25trivox(
35 1 NSN ,NSNR ,ISZNSNR ,I_MEM ,VMAXDT ,
36 2 IRECT ,X ,STF ,STFN ,XYZM ,
37 3 NSV ,II_STOK ,CAND_N ,ESHIFT ,CAND_E ,
38 4 MULNSN ,NOINT ,V ,BGAPSMX ,
39 5 VOXEL ,NBX ,NBY ,NBZ ,PMAX_GAP ,
40 6 NRTM ,GAP_S ,GAP_M ,MARGE ,CURV_MAX ,
41 7 NIN ,ITASK ,PENE_OLD,ITAB ,NBINFLG ,
42 8 MBINFLG,ILEV ,MSEGTYP ,
43 9 FLAGREMNODE,KREMNOD,REMNOD ,
44 A IGAP ,GAP_S_L,GAP_M_L ,ICODT ,ISKEW ,
45 B DRAD ,DGAPLOAD )
46C============================================================================
47C M o d u l e s
48C-----------------------------------------------
49 USE tri7box
50 use my_alloc_mod
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 "com04_c.inc"
66#include "param_c.inc"
67#include "task_c.inc"
68C-----------------------------------------------
69C ROLE DE LA ROUTINE:
70C ===================
71C CLASSE LES NOEUDS DANS DES BOITES
72C RECHERCHE POUR CHAQUE FACETTE DES BOITES CONCERNES
73C RECHERCHE DES CANDIDATS
74C-----------------------------------------------
75C D u m m y A r g u m e n t s
76C
77C NOM DESCRIPTION E/S
78C
79C IRECT(4,*) TABLEAU DES CONEC FACETTES E
80C X(3,*) COORDONNEES NODALES E
81C NSV NOS SYSTEMES DES NOEUDS E
82C XMAX plus grande abcisse existante E
83C XMAX plus grande ordonn. existante E
84C XMAX plus grande cote existante E
85C I_STOK niveau de stockage des couples
86C candidats impact E/S
87C CAND_N boites resultats noeuds
88C CAND_E adresses des boites resultat elements
89C MULNSN = MULTIMP*NSN TAILLE MAX ADMISE MAINTENANT POUR LES
90C COUPLES NOEUDS,ELT CANDIDATS
91C NOINT NUMERO USER DE L'INTERFACE
92C
93C PROV_N CAND_N provisoire (variable static dans i7tri)
94C PROV_E CAND_E provisoire (variable static dans i7tri)
95
96C VOXEL(ix,iy,iz) contient le numero local du premier noeud de
97C la boite
98C NEXT_NOD(i) noeud suivant dans la meme boite (si /= 0)
99C LAST_NOD(i) dernier noeud dans la meme boite (si /= 0)
100C utilise uniquement pour aller directement du premier
101C noeud au dernier
102C-----------------------------------------------
103C D u m m y A r g u m e n t s
104C-----------------------------------------------
105 INTEGER I_MEM,ESHIFT,NSN,ISZNSNR,NRTM,NIN,ITASK,IGAP,
106 . MULNSN,NOINT,NSNR,NBX,NBY,NBZ,
107 . NSV(*),CAND_N(*),CAND_E(*),
108 . IRECT(4,*), VOXEL(NBX+2,NBY+2,NBZ+2),II_STOK,ITAB(*),
109 . NBINFLG(*),MBINFLG(*),ILEV,MSEGTYP(*),
110 . FLAGREMNODE,KREMNOD(*),REMNOD(*), ICODT(*), ISKEW(*)
111C REAL
112 my_real
113 . X(3,*),V(3,*),XYZM(6),STF(*),STFN(*),GAP_S(*),GAP_M(*),
114 . CURV_MAX(*),PENE_OLD(5,NSN),GAP_S_L(*),GAP_M_L(*),
115 . MARGE,BGAPSMX,PMAX_GAP,VMAXDT
116 my_real , INTENT(IN) :: dgapload ,drad
117C-----------------------------------------------
118C L o c a l V a r i a b l e s
119C-----------------------------------------------
120 INTEGER NB_NCN,NB_NCN1,NB_ECN,I,J,DIR,NB_NC,NB_EC,
121 . N1,N2,N3,N4,NN,NE,K,L,J_STOK,II,JJ,
122 . PROV_N(MVSIZ),PROV_E(MVSIZ),
123 . nsnf, nsnl,m,ns1,ns2,nse,ns,ip,delnod
124C REAL
125 my_real
126 . dx,dy,dz,xs,ys,zs,xx,sx,sy,sz,s2,
127 . xmin, xmax,ymin, ymax,zmin, zmax, tz, gapsmx, gapl,
128 . xx1,xx2,xx3,xx4,yy1,yy2,yy3,yy4,zz1,zz2,zz3,zz4,
129 . d1x,d1y,d1z,d2x,d2y,d2z,dd1,dd2,d2,a2,gs
130c
131 INTEGER IX,IY,IZ,NEXT,M1,M2,M3,M4,
132 . IX1,IY1,IZ1,IX2,IY2,IZ2
133 INTEGER, DIMENSION(:),ALLOCATABLE :: LAST_NOD
134 INTEGER, DIMENSION(:),ALLOCATABLE :: IIX,IIY,IIZ
135 my_real
136 . xminb,yminb,zminb,xmaxb,ymaxb,zmaxb,
137 . xmine,ymine,zmine,xmaxe,ymaxe,zmaxe,aaa
138 INTEGER FIRST,NEW,LAST
139 SAVE IIX,IIY,IIZ
140 INTEGER, DIMENSION(:),ALLOCATABLE :: TAG
141C-----------------------------------------------
142 CALL my_alloc(tag,numnod)
143 CALL my_alloc(last_nod,nsn+nsnr)
144 IF(itask == 0)THEN
145 ALLOCATE(next_nod(nsn+nsnr))
146 ALLOCATE(iix(nsn+nsnr))
147 ALLOCATE(iiy(nsn+nsnr))
148 ALLOCATE(iiz(nsn+nsnr))
149 END IF
150C Barrier to wait init voxel and allocation NEX_NOD
151 CALL my_barrier
152C Phase initiale de construction de BPE et BPN deplacee de I7BUCE => I7TRI
153C
154 xmin = xyzm(1)
155 ymin = xyzm(2)
156 zmin = xyzm(3)
157 xmax = xyzm(4)
158 ymax = xyzm(5)
159 zmax = xyzm(6)
160
161c dev future: xminb plus grand que xmin...
162 xminb = xmin
163 yminb = ymin
164 zminb = zmin
165 xmaxb = xmax
166 ymaxb = ymax
167 zmaxb = zmax
168
169C=======================================================================
170C 1 mise des noeuds dans les boites
171C=======================================================================
172 IF(itask == 0)THEN
173 DO i=1,nsn
174 iix(i)=0
175 iiy(i)=0
176 iiz(i)=0
177 IF(stfn(i) <= zero)cycle
178 j=nsv(i)
179C Optimisation // recherche les noeuds compris dans xmin xmax des
180C elements du processeur
181 IF(x(1,j) < xmin) cycle
182 IF(x(1,j) > xmax) cycle
183 IF(x(2,j) < ymin) cycle
184 IF(x(2,j) > ymax) cycle
185 IF(x(3,j) < zmin) cycle
186 IF(x(3,j) > zmax) cycle
187
188 iix(i)=int(nbx*(x(1,j)-xminb)/(xmaxb-xminb))
189 iiy(i)=int(nby*(x(2,j)-yminb)/(ymaxb-yminb))
190 iiz(i)=int(nbz*(x(3,j)-zminb)/(zmaxb-zminb))
191
192 iix(i)=max(1,2+min(nbx,iix(i)))
193 iiy(i)=max(1,2+min(nby,iiy(i)))
194 iiz(i)=max(1,2+min(nbz,iiz(i)))
195
196 first = voxel(iix(i),iiy(i),iiz(i))
197 IF(first == 0)THEN
198c empty cell
199 voxel(iix(i),iiy(i),iiz(i)) = i ! first
200 next_nod(i) = 0 ! last one
201 last_nod(i) = 0 ! no last
202 ELSEIF(last_nod(first) == 0)THEN
203c cell containing one node
204c add as next node
205 next_nod(first) = i ! next
206 last_nod(first) = i ! last
207 next_nod(i) = 0 ! last one
208 ELSE
209c
210c jump to the last node of the cell
211 last = last_nod(first) ! last node in this voxel
212 next_nod(last) = i ! next
213 last_nod(first) = i ! last
214 next_nod(i) = 0 ! last one
215 ENDIF
216 ENDDO
217C=======================================================================
218C 2 mise des noeuds dans les boites
219C candidats non locaux en SPMD
220C=======================================================================
221 DO j = 1, nsnr
222 iix(nsn+j)=int(nbx*(xrem(1,j)-xminb)/(xmaxb-xminb))
223 iiy(nsn+j)=int(nby*(xrem(2,j)-yminb)/(ymaxb-yminb))
224 iiz(nsn+j)=int(nbz*(xrem(3,j)-zminb)/(zmaxb-zminb))
225 iix(nsn+j)=max(1,2+min(nbx,iix(nsn+j)))
226 iiy(nsn+j)=max(1,2+min(nby,iiy(nsn+j)))
227 iiz(nsn+j)=max(1,2+min(nbz,iiz(nsn+j)))
228
229 first = voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))
230 IF(first == 0)THEN
231c empty cell
232 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j)) = nsn+j ! first
233 next_nod(nsn+j) = 0 ! last one
234 last_nod(nsn+j) = 0 ! no last
235 ELSEIF(last_nod(first) == 0)THEN
236c cell containing one node
237c add as next node
238 next_nod(first) = nsn+j ! next
239 last_nod(first) = nsn+j ! last
240 next_nod(nsn+j) = 0 ! last one
241 ELSE
242c
243c jump to the last node of the cell
244 last = last_nod(first) ! last node in this voxel
245 next_nod(last) = nsn+j ! next
246 last_nod(first) = nsn+j ! last
247 next_nod(nsn+j) = 0 ! last one
248 ENDIF
249 ENDDO
250 END IF
251C Barrier to wait task0 treatment
252 CALL my_barrier
253C=======================================================================
254C 3 recherche des boites concernant chaque facette
255C et creation des candidats
256C=======================================================================
257 j_stok = 0
258C
259 IF(flagremnode == 2)THEN
260 DO i=1,numnod
261 tag(i) = 0
262 ENDDO
263
264 DO ne=1,nrtm
265C don't keep deleted segments
266 IF(stf(ne) <= zero)cycle
267 k = kremnod(2*(ne-1)+1)+1
268 l = kremnod(2*(ne-1)+2)
269 DO i=k,l
270 tag(remnod(i)) = 1
271 ENDDO
272c
273 aaa = marge+curv_max(ne)+max(max(bgapsmx+gap_m(ne),pmax_gap)+dgapload,drad)+vmaxdt
274 +
275
276
277c il est possible d'ameliorer l'algo en decoupant la facette
278c en 2(4,3,6,9...) si la facette est grande devant AAA et inclinee
279
280 m1 = irect(1,ne)
281 m2 = irect(2,ne)
282 m3 = irect(3,ne)
283 m4 = irect(4,ne)
284
285 xx1=x(1,m1)
286 xx2=x(1,m2)
287 xx3=x(1,m3)
288 xx4=x(1,m4)
289 xmaxe=max(xx1,xx2,xx3,xx4)
290 xmine=min(xx1,xx2,xx3,xx4)
291
292 yy1=x(2,m1)
293 yy2=x(2,m2)
294 yy3=x(2,m3)
295 yy4=x(2,m4)
296 ymaxe=max(yy1,yy2,yy3,yy4)
297 ymine=min(yy1,yy2,yy3,yy4)
298
299 zz1=x(3,m1)
300 zz2=x(3,m2)
301 zz3=x(3,m3)
302 zz4=x(3,m4)
303 zmaxe=max(zz1,zz2,zz3,zz4)
304 zmine=min(zz1,zz2,zz3,zz4)
305
306
307c calcul de la surface (pour elimination future de candidats)
308
309 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
310 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
311 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
312 s2 = sx*sx + sy*sy + sz*sz
313
314c indice des voxels occupes par la facette
315
316 ix1=int(nbx*(xmine-aaa-xminb)/(xmaxb-xminb))
317 iy1=int(nby*(ymine-aaa-yminb)/(ymaxb-yminb))
318 iz1=int(nbz*(zmine-aaa-zminb)/(zmaxb-zminb))
319
320 ix1=max(1,2+min(nbx,ix1))
321 iy1=max(1,2+min(nby,iy1))
322 iz1=max(1,2+min(nbz,iz1))
323
324 ix2=int(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb))
325 iy2=int(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb))
326 iz2=int(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb))
327
328 ix2=max(1,2+min(nbx,ix2))
329 iy2=max(1,2+min(nby,iy2))
330 iz2=max(1,2+min(nbz,iz2))
331
332cc nbpelem = 0
333cc nnpelem = 0
334cc nnr0pelem = 0
335cc nnrpelem = 0
336
337 DO iz = iz1,iz2
338 DO iy = iy1,iy2
339 DO ix = ix1,ix2
340
341cc nbpelem = nbpelem + 1
342
343 jj = voxel(ix,iy,iz)
344
345 DO WHILE(jj /= 0)
346 delnod = 0
347cc nnpelem = nnpelem + 1
348
349 IF(jj<=nsn)THEN
350 nn=nsv(jj)
351 IF(nn == m1)GOTO 200
352 IF(nn == m2)GOTO 200
353 IF(nn == m3)GOTO 200
354 IF(nn == m4)GOTO 200
355 IF(tag(nn) == 1)GOTO 200
356 xs = x(1,nn)
357 ys = x(2,nn)
358 zs = x(3,nn)
359c PMAX_GAP is a global overestimate penetration
360c NEED to communicate in SPMD
361c VMAXDT is a local overestimate of relative incremental displacement
362c NO need to communicate in SPMD
363
364 aaa = marge + curv_max(ne)
365 + + max(gap_s(jj)+gap_m(ne)+dgapload,drad)
366 + +vmaxdt
367 ELSE
368 j=jj-nsn
369
370 k = kremnod(2*(ne-1)+2) + 1
371 l = kremnod(2*(ne-1)+3)
372
373 DO m=k,l
374 IF(remnod(m) == -irem(2,j) ) THEN
375 delnod = delnod + 1
376 EXIT
377 ENDIF
378 ENDDO
379
380 IF(delnod /= 0)GOTO 200
381
382 xs = xrem(1,j)
383 ys = xrem(2,j)
384 zs = xrem(3,j)
385 aaa = marge+curv_max(ne)
386c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
387c +EDGE_L2(JJ) remote
388c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
389 + + max(xrem(igapxremp,j)+gap_m(ne)+dgapload,drad)
390 + + vmaxdt
391 ENDIF
392
393 IF(xs<=xmine-aaa)GOTO 200
394 IF(xs>=xmaxe+aaa)GOTO 200
395 IF(ys<=ymine-aaa)GOTO 200
396 IF(ys>=ymaxe+aaa)GOTO 200
397 IF(zs<=zmine-aaa)GOTO 200
398 IF(zs>=zmaxe+aaa)GOTO 200
399
400c sousestimation de la distance**2 pour elimination de candidats
401
402cc nnr0pelem = nnr0pelem + 1
403
404 d1x = xs - xx1
405 d1y = ys - yy1
406 d1z = zs - zz1
407 d2x = xs - xx2
408 d2y = ys - yy2
409 d2z = zs - zz2
410 dd1 = d1x*sx+d1y*sy+d1z*sz
411 dd2 = d2x*sx+d2y*sy+d2z*sz
412 IF(dd1*dd2 > zero)THEN
413 d2 = min(dd1*dd1,dd2*dd2)
414 a2 = aaa*aaa*s2
415 IF(d2 > a2)GOTO 200
416 ENDIF
417
418cc nnrpelem = nnrpelem + 1
419
420 j_stok = j_stok + 1
421 prov_n(j_stok) = jj
422 prov_e(j_stok) = ne
423 IF(j_stok == nvsiz)THEN
424
425 CALL i25sto(
426 1 nvsiz ,irect ,x ,nsv ,ii_stok,
427 2 cand_n,cand_e ,mulnsn,noint ,marge ,
428 3 i_mem ,prov_n ,prov_e,eshift,v ,
429 4 nsn ,nrtm ,gap_s ,gap_m ,curv_max,nin ,
430 5 pene_old,nbinflg ,mbinflg,ilev,msegtyp,
431 6 itab ,igap ,gap_s_l,gap_m_l,icodt,iskew,
432 7 drad ,dgapload)
433 IF(i_mem==2)GOTO 100
434 j_stok = 0
435 ENDIF
436
437 200 CONTINUE
438
439 jj = next_nod(jj)
440
441 ENDDO ! WHILE(JJ /= 0)
442
443 ENDDO
444 ENDDO
445 ENDDO
446 k = kremnod(2*(ne-1)+1)+1
447 l = kremnod(2*(ne-1)+2)
448 DO i=k,l
449 tag(remnod(i)) = 0
450 ENDDO
451
452cc nbpelg = nbpelg + nbpelem
453cc nnpelg = nnpelg + nnpelem
454cc nnrpelg = nnrpelg + nnrpelem
455cc nnr0pelg = nnr0pelg + nnr0pelem
456 ENDDO
457
458 ELSE !FLAGREMNODE
459 DO ne=1,nrtm
460C don't keep deleted segments
461 IF(stf(ne) <= zero)cycle
462c
463 aaa = marge+curv_max(ne)+max(max(bgapsmx+gap_m(ne),pmax_gap)+dgapload,drad)+vmaxdt
464 +
465
466
467c il est possible d'ameliorer l'algo en decoupant la facette
468c en 2(4,3,6,9...) si la facette est grande devant AAA et inclinee
469
470 m1 = irect(1,ne)
471 m2 = irect(2,ne)
472 m3 = irect(3,ne)
473 m4 = irect(4,ne)
474
475 xx1=x(1,m1)
476 xx2=x(1,m2)
477 xx3=x(1,m3)
478 xx4=x(1,m4)
479 xmaxe=max(xx1,xx2,xx3,xx4)
480 xmine=min(xx1,xx2,xx3,xx4)
481
482 yy1=x(2,m1)
483 yy2=x(2,m2)
484 yy3=x(2,m3)
485 yy4=x(2,m4)
486 ymaxe=max(yy1,yy2,yy3,yy4)
487 ymine=min(yy1,yy2,yy3,yy4)
488
489 zz1=x(3,m1)
490 zz2=x(3,m2)
491 zz3=x(3,m3)
492 zz4=x(3,m4)
493 zmaxe=max(zz1,zz2,zz3,zz4)
494 zmine=min(zz1,zz2,zz3,zz4)
495
496
497c calcul de la surface (pour elimination future de candidats)
498
499 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
500 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
501 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
502 s2 = sx*sx + sy*sy + sz*sz
503
504c indice des voxels occupes par la facette
505
506 ix1=int(nbx*(xmine-aaa-xminb)/(xmaxb-xminb))
507 iy1=int(nby*(ymine-aaa-yminb)/(ymaxb-yminb))
508 iz1=int(nbz*(zmine-aaa-zminb)/(zmaxb-zminb))
509
510 ix1=max(1,2+min(nbx,ix1))
511 iy1=max(1,2+min(nby,iy1))
512 iz1=max(1,2+min(nbz,iz1))
513
514 ix2=int(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb))
515 iy2=int(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb))
516 iz2=int(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb))
517
518 ix2=max(1,2+min(nbx,ix2))
519 iy2=max(1,2+min(nby,iy2))
520 iz2=max(1,2+min(nbz,iz2))
521
522cc nbpelem = 0
523cc nnpelem = 0
524cc nnr0pelem = 0
525cc nnrpelem = 0
526
527 DO iz = iz1,iz2
528 DO iy = iy1,iy2
529 DO ix = ix1,ix2
530
531cc nbpelem = nbpelem + 1
532
533 jj = voxel(ix,iy,iz)
534
535 DO WHILE(jj /= 0)
536
537cc nnpelem = nnpelem + 1
538
539 IF(jj<=nsn)THEN
540 nn=nsv(jj)
541 IF(nn == m1)GOTO 300
542 IF(nn == m2)GOTO 300
543 IF(nn == m3)GOTO 300
544 IF(nn == m4)GOTO 300
545 xs = x(1,nn)
546 ys = x(2,nn)
547 zs = x(3,nn)
548c PMAX_GAP is a global overestimate penetration
549c NEED to communicate in SPMD
550c VMAXDT is a local overestimate of relative incremental displacement
551c NO need to communicate in SPMD
552
553 aaa = marge + curv_max(ne)
554 + + max(gap_s(jj)+gap_m(ne)+dgapload,drad)
555 + +vmaxdt
556 ELSE
557 j=jj-nsn
558
559 xs = xrem(1,j)
560 ys = xrem(2,j)
561 zs = xrem(3,j)
562 aaa = marge+curv_max(ne)
563c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
564c +EDGE_L2(JJ) remote
565c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
566 + + max(xrem(igapxremp,j)+gap_m(ne)+dgapload,drad)
567 + + vmaxdt
568 ENDIF
569
570 IF(xs<=xmine-aaa)GOTO 300
571 IF(xs>=xmaxe+aaa)GOTO 300
572 IF(ys<=ymine-aaa)GOTO 300
573 IF(ys>=ymaxe+aaa)GOTO 300
574 IF(zs<=zmine-aaa)GOTO 300
575 IF(zs>=zmaxe+aaa)GOTO 300
576
577c sousestimation de la distance**2 pour elimination de candidats
578
579cc nnr0pelem = nnr0pelem + 1
580
581 d1x = xs - xx1
582 d1y = ys - yy1
583 d1z = zs - zz1
584 d2x = xs - xx2
585 d2y = ys - yy2
586 d2z = zs - zz2
587 dd1 = d1x*sx+d1y*sy+d1z*sz
588 dd2 = d2x*sx+d2y*sy+d2z*sz
589 IF(dd1*dd2 > zero)THEN
590 d2 = min(dd1*dd1,dd2*dd2)
591 a2 = aaa*aaa*s2
592 IF(d2 > a2)GOTO 300
593 ENDIF
594
595cc nnrpelem = nnrpelem + 1
596
597 j_stok = j_stok + 1
598 prov_n(j_stok) = jj
599 prov_e(j_stok) = ne
600 IF(j_stok == nvsiz)THEN
601
602 CALL i25sto(
603 1 nvsiz ,irect ,x ,nsv ,ii_stok,
604 2 cand_n,cand_e ,mulnsn,noint ,marge ,
605 3 i_mem ,prov_n ,prov_e,eshift,v ,
606 4 nsn ,nrtm ,gap_s ,gap_m ,curv_max,nin ,
607 5 pene_old,nbinflg ,mbinflg,ilev,msegtyp,
608 6 itab ,igap ,gap_s_l,gap_m_l,icodt,iskew,
609 7 drad ,dgapload)
610 IF(i_mem==2)GOTO 100
611 j_stok = 0
612 ENDIF
613
614 300 CONTINUE
615
616 jj = next_nod(jj)
617
618 ENDDO ! WHILE(JJ /= 0)
619
620 ENDDO
621 ENDDO
622 ENDDO
623
624cc nbpelg = nbpelg + nbpelem
625cc nnpelg = nnpelg + nnpelem
626cc nnrpelg = nnrpelg + nnrpelem
627cc nnr0pelg = nnr0pelg + nnr0pelem
628 ENDDO
629 END IF !FLAGREMNODE
630C-------------------------------------------------------------------------
631C FIN DU TRI
632C-------------------------------------------------------------------------
633 IF(j_stok/=0)CALL i25sto(
634 1 j_stok,irect ,x ,nsv ,ii_stok,
635 2 cand_n,cand_e ,mulnsn,noint ,marge ,
636 3 i_mem ,prov_n ,prov_e,eshift,v ,
637 4 nsn ,nrtm ,gap_s ,gap_m ,curv_max,nin ,
638 5 pene_old,nbinflg,mbinflg,ilev ,msegtyp,
639 6 itab ,igap ,gap_s_l,gap_m_l,icodt,iskew,
640 7 drad ,dgapload)
641
642C=======================================================================
643C 4 remise a zero des noeuds dans les boites
644C=======================================================================
645 100 CONTINUE
646
647C Barrier to avoid reinitialization before end of sorting
648 CALL my_barrier
649 nsnf = 1 + itask*nsn / nthread
650 nsnl = (itask+1)*nsn / nthread
651
652 DO i=nsnf,nsnl
653 IF(iix(i)/=0)THEN
654 voxel(iix(i),iiy(i),iiz(i))=0
655 ENDIF
656 ENDDO
657C=======================================================================
658C 5 remise a zero des noeuds dans les boites
659C candidats non locaux en SPMD
660C=======================================================================
661 nsnf = 1 + itask*nsnr / nthread
662 nsnl = (itask+1)*nsnr / nthread
663 DO j = nsnf, nsnl
664 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))=0
665 ENDDO
666
667C
668 CALL my_barrier()
669 IF(itask == 0)THEN
670 DEALLOCATE(next_nod)
671 DEALLOCATE(iix)
672 DEALLOCATE(iiy)
673 DEALLOCATE(iiz)
674 ENDIF
675
676 DEALLOCATE(tag)
677 DEALLOCATE(last_nod)
678 RETURN
679 END
680
subroutine i25trivox(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, flagremnode, kremnod, remnod, igap, gap_s_l, gap_m_l, icodt, iskew, drad, dgapload)
Definition i25trivox.F:46
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 igapxremp
Definition tri7box.F:423
integer, dimension(:), allocatable next_nod
Definition tri7box.F:48
integer, dimension(:,:), allocatable irem
Definition tri7box.F:339
subroutine i25sto(j_stok, irect, x, nsv, local_i_stok, local_cand_n, local_cand_e, marge, prov_n, prov_e, eshift, nsn, nrtm, gap_s, gap_m, nbinflg, mbinflg, ilev, msegtyp, igap, gap_s_l, gap_m_l, edge_l2, icode, iskew, drad, dgapload, nrtmt)
Definition i25sto.F:42
subroutine my_barrier
Definition machine.F:31