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 of the routine:
70C ===================
71C classify nodes into boxes
72C search for concerned boxes for each facet
73C search for candidates
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,*) ARRAY OF CONEC FACETTES E
80C X(3,*) COORDONNEES NODALES E
81C NSV NOS SYSTEMES DES NODES E
82C Xmax larger abcisse existing e
83C XMAX largest order.existing E
84C Xmax larger existing side E
85C i_stok storage level of pairs
86C CANDIDATES impact E/S
87C CAND_N boites resultats nodes C CAND_E adresses des boites resultat elements
88C mulnsn = multimp*nsn max size now allowed for the
89C COUPLES NODES,ELT CANDIDATES
90C NOINT INTERFACE USER NUMBER
91C
92C Prov_n Provisional Cand_n (static variable in i7tri)
93C PROV_E CAND_E provisoire (variable static in i7tri)
94
95C voxel(ix,iy,iz) contains the local number of the first node of
96C the box
97C Next_nod (i) Next node in the same box (if /= 0)
98C Last_nod (i) Last node in the same box (if /= 0)
99C used only to go directly from the first
100C node at the last
101C-----------------------------------------------
102C D u m m y A r g u m e n t s
103C-----------------------------------------------
104 INTEGER I_MEM,ESHIFT,NSN,ISZNSNR,NRTM,NIN,ITASK,IGAP,
105 . MULNSN,NOINT,NSNR,NBX,NBY,NBZ,
106 . NSV(*),CAND_N(*),CAND_E(*),
107 . IRECT(4,*), VOXEL(NBX+2,NBY+2,NBZ+2),II_STOK,ITAB(*),
108 . NBINFLG(*),MBINFLG(*),ILEV,MSEGTYP(*),
109 . FLAGREMNODE,KREMNOD(*),REMNOD(*), ICODT(*), ISKEW(*)
110C REAL
111 my_real
112 . X(3,*),V(3,*),XYZM(6),STF(*),STFN(*),GAP_S(*),GAP_M(*),
113 . CURV_MAX(*),PENE_OLD(5,NSN),GAP_S_L(*),GAP_M_L(*),
114 . MARGE,BGAPSMX,PMAX_GAP,VMAXDT
115 my_real , INTENT(IN) :: dgapload ,drad
116C-----------------------------------------------
117C L o c a l V a r i a b l e s
118C-----------------------------------------------
119 INTEGER I,J,
120 . NN,NE,K,L,J_STOK,JJ,
121 . PROV_N(MVSIZ),PROV_E(MVSIZ),
122 . nsnf, nsnl,m,delnod
123C REAL
124 my_real
125 . xs,ys,zs,sx,sy,sz,s2,
126 . xmin, xmax,ymin, ymax,zmin, zmax,
127 . xx1,xx2,xx3,xx4,yy1,yy2,yy3,yy4,zz1,zz2,zz3,zz4,
128 . d1x,d1y,d1z,d2x,d2y,d2z,dd1,dd2,d2,a2
129c
130 INTEGER IX,IY,IZ,M1,M2,M3,M4,
131 . IX1,IY1,IZ1,IX2,IY2,IZ2
132 INTEGER, DIMENSION(:),ALLOCATABLE :: LAST_NOD
133 INTEGER, DIMENSION(:),ALLOCATABLE :: IIX,IIY,IIZ
134 my_real
135 . xminb,yminb,zminb,xmaxb,ymaxb,zmaxb,
136 . xmine,ymine,zmine,xmaxe,ymaxe,zmaxe,aaa
137 INTEGER FIRST,LAST
138 SAVE IIX,IIY,IIZ
139 INTEGER, DIMENSION(:),ALLOCATABLE :: TAG
140C-----------------------------------------------
141 CALL my_alloc(tag,numnod)
142 CALL my_alloc(last_nod,nsn+nsnr)
143 IF(itask == 0)THEN
144 ALLOCATE(next_nod(nsn+nsnr))
145 ALLOCATE(iix(nsn+nsnr))
146 ALLOCATE(iiy(nsn+nsnr))
147 ALLOCATE(iiz(nsn+nsnr))
148 END IF
149C Barrier to wait init voxel and allocation NEX_NOD
150 CALL my_barrier
151C initial construction phase of bpe and bpn moved from i7buce => i7tri
152C
153 xmin = xyzm(1)
154 ymin = xyzm(2)
155 zmin = xyzm(3)
156 xmax = xyzm(4)
157 ymax = xyzm(5)
158 zmax = xyzm(6)
159
160c Dev Future: xminb larger than Xmin ...
161 xminb = xmin
162 yminb = ymin
163 zminb = zmin
164 xmaxb = xmax
165 ymaxb = ymax
166 zmaxb = zmax
167
168C=======================================================================
169C 1 placing nodes into boxes
170C=======================================================================
171 IF(itask == 0)THEN
172 DO i=1,nsn
173 iix(i)=0
174 iiy(i)=0
175 iiz(i)=0
176 IF(stfn(i) <= zero)cycle
177 j=nsv(i)
178C optimization // search for nodes within xmin xmax of the
179C processor elements
180 IF(x(1,j) < xmin) cycle
181 IF(x(1,j) > xmax) cycle
182 IF(x(2,j) < ymin) cycle
183 IF(x(2,j) > ymax) cycle
184 IF(x(3,j) < zmin) cycle
185 IF(x(3,j) > zmax) cycle
186
187 iix(i)=int(nbx*(x(1,j)-xminb)/(xmaxb-xminb))
188 iiy(i)=int(nby*(x(2,j)-yminb)/(ymaxb-yminb))
189 iiz(i)=int(nbz*(x(3,j)-zminb)/(zmaxb-zminb))
190
191 iix(i)=max(1,2+min(nbx,iix(i)))
192 iiy(i)=max(1,2+min(nby,iiy(i)))
193 iiz(i)=max(1,2+min(nbz,iiz(i)))
194
195 first = voxel(iix(i),iiy(i),iiz(i))
196 IF(first == 0)THEN
197c Empty Cell
198 voxel(iix(i),iiy(i),iiz(i)) = i ! first
199 next_nod(i) = 0 ! last one
200 last_nod(i) = 0 ! no last
201 ELSEIF(last_nod(first) == 0)THEN
202c cell containing one node
203c add as next node
204 next_nod(first) = i ! next
205 last_nod(first) = i ! last
206 next_nod(i) = 0 ! last one
207 ELSE
208c
209c jump to the last node of the cell
210 last = last_nod(first) ! last node in this voxel
211 next_nod(last) = i ! next
212 last_nod(first) = i ! last
213 next_nod(i) = 0 ! last one
214 ENDIF
215 ENDDO
216C=======================================================================
217C 2 placing nodes into boxes
218C non -local candidates in SPMD
219C=======================================================================
220 DO j = 1, nsnr
221 iix(nsn+j)=int(nbx*(xrem(1,j)-xminb)/(xmaxb-xminb))
222 iiy(nsn+j)=int(nby*(xrem(2,j)-yminb)/(ymaxb-yminb))
223 iiz(nsn+j)=int(nbz*(xrem(3,j)-zminb)/(zmaxb-zminb))
224 iix(nsn+j)=max(1,2+min(nbx,iix(nsn+j)))
225 iiy(nsn+j)=max(1,2+min(nby,iiy(nsn+j)))
226 iiz(nsn+j)=max(1,2+min(nbz,iiz(nsn+j)))
227
228 first = voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))
229 IF(first == 0)THEN
230c Empty Cell
231 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j)) = nsn+j ! first
232 next_nod(nsn+j) = 0 ! last one
233 last_nod(nsn+j) = 0 ! no last
234 ELSEIF(last_nod(first) == 0)THEN
235c cell containing one node
236c add as next node
237 next_nod(first) = nsn+j ! next
238 last_nod(first) = nsn+j ! last
239 next_nod(nsn+j) = 0 ! last one
240 ELSE
241c
242c jump to the last node of the cell
243 last = last_nod(first) ! last node in this voxel
244 next_nod(last) = nsn+j ! next
245 last_nod(first) = nsn+j ! last
246 next_nod(nsn+j) = 0 ! last one
247 ENDIF
248 ENDDO
249 END IF
250C Barrier to wait task0 treatment
251 CALL my_barrier
252C=======================================================================
253C 3 Searching boxes concerning each facet
254C and creation of candidates
255C=======================================================================
256 j_stok = 0
257C
258 IF(flagremnode == 2)THEN
259 DO i=1,numnod
260 tag(i) = 0
261 ENDDO
262
263 DO ne=1,nrtm
264C don't keep deleted segments
265 IF(stf(ne) <= zero)cycle
266 k = kremnod(2*(ne-1)+1)+1
267 l = kremnod(2*(ne-1)+2)
268 DO i=k,l
269 tag(remnod(i)) = 1
270 ENDDO
271c
272 aaa = marge+curv_max(ne)+max(max(bgapsmx+gap_m(ne),pmax_gap)+dgapload,drad)+vmaxdt
273 +
274
275
276c It is possible to improve the algo by cutting the facet
277c in 2 (4,3,6,9 ...) if the facet is large in front of AAA and inclinee
278
279 m1 = irect(1,ne)
280 m2 = irect(2,ne)
281 m3 = irect(3,ne)
282 m4 = irect(4,ne)
283
284 xx1=x(1,m1)
285 xx2=x(1,m2)
286 xx3=x(1,m3)
287 xx4=x(1,m4)
288 xmaxe=max(xx1,xx2,xx3,xx4)
289 xmine=min(xx1,xx2,xx3,xx4)
290
291 yy1=x(2,m1)
292 yy2=x(2,m2)
293 yy3=x(2,m3)
294 yy4=x(2,m4)
295 ymaxe=max(yy1,yy2,yy3,yy4)
296 ymine=min(yy1,yy2,yy3,yy4)
297
298 zz1=x(3,m1)
299 zz2=x(3,m2)
300 zz3=x(3,m3)
301 zz4=x(3,m4)
302 zmaxe=max(zz1,zz2,zz3,zz4)
303 zmine=min(zz1,zz2,zz3,zz4)
304
305
306c surface calculation (for future elimination of candidates)
307
308 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
309 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
310 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
311 s2 = sx*sx + sy*sy + sz*sz
312
313c index of voxels occupied by the facet
314
315 ix1=int(nbx*(xmine-aaa-xminb)/(xmaxb-xminb))
316 iy1=int(nby*(ymine-aaa-yminb)/(ymaxb-yminb))
317 iz1=int(nbz*(zmine-aaa-zminb)/(zmaxb-zminb))
318
319 ix1=max(1,2+min(nbx,ix1))
320 iy1=max(1,2+min(nby,iy1))
321 iz1=max(1,2+min(nbz,iz1))
322
323 ix2=int(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb))
324 iy2=int(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb))
325 iz2=int(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb))
326
327 ix2=max(1,2+min(nbx,ix2))
328 iy2=max(1,2+min(nby,iy2))
329 iz2=max(1,2+min(nbz,iz2))
330
331cc nbpelem = 0
332cc'nnpelem = 0
333cc'nnr0pelem = 0
334cc'nnrpelem = 0
335
336 DO iz = iz1,iz2
337 DO iy = iy1,iy2
338 DO ix = ix1,ix2
339
340cc nbpelem = nbpelem + 1
341
342 jj = voxel(ix,iy,iz)
343
344 DO WHILE(jj /= 0)
345 delnod = 0
346cc nnpelem = nnpelem + 1
347
348 IF(jj<=nsn)THEN
349 nn=nsv(jj)
350 IF(nn == m1)GOTO 200
351 IF(nn == m2)GOTO 200
352 IF(nn == m3)GOTO 200
353 IF(nn == m4)GOTO 200
354 IF(tag(nn) == 1)GOTO 200
355 xs = x(1,nn)
356 ys = x(2,nn)
357 zs = x(3,nn)
358c PMAX_GAP is a global overestimate penetration
359c NEED to communicate in SPMD
360c VMAXDT is a local overestimate of relative incremental displacement
361c NO need to communicate in SPMD
362
363 aaa = marge + curv_max(ne)
364 + + max(gap_s(jj)+gap_m(ne)+dgapload,drad)
365 + +vmaxdt
366 ELSE
367 j=jj-nsn
368
369 k = kremnod(2*(ne-1)+2) + 1
370 l = kremnod(2*(ne-1)+3)
371
372 DO m=k,l
373 IF(remnod(m) == -irem(2,j) ) THEN
374 delnod = delnod + 1
375 EXIT
376 ENDIF
377 ENDDO
378
379 IF(delnod /= 0)GOTO 200
380
381 xs = xrem(1,j)
382 ys = xrem(2,j)
383 zs = xrem(3,j)
384 aaa = marge+curv_max(ne)
385c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
386c +EDGE_L2(JJ) remote
387c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
388 + + max(xrem(igapxremp,j)+gap_m(ne)+dgapload,drad)
389 + + vmaxdt
390 ENDIF
391
392 IF(xs<=xmine-aaa)GOTO 200
393 IF(xs>=xmaxe+aaa)GOTO 200
394 IF(ys<=ymine-aaa)GOTO 200
395 IF(ys>=ymaxe+aaa)GOTO 200
396 IF(zs<=zmine-aaa)GOTO 200
397 IF(zs>=zmaxe+aaa)GOTO 200
398
399c underestimation of distance**2 for candidate elimination
400
401cc'nnr0pelem = nnr0pelem + 1
402
403 d1x = xs - xx1
404 d1y = ys - yy1
405 d1z = zs - zz1
406 d2x = xs - xx2
407 d2y = ys - yy2
408 d2z = zs - zz2
409 dd1 = d1x*sx+d1y*sy+d1z*sz
410 dd2 = d2x*sx+d2y*sy+d2z*sz
411 IF(dd1*dd2 > zero)THEN
412 d2 = min(dd1*dd1,dd2*dd2)
413 a2 = aaa*aaa*s2
414 IF(d2 > a2)GOTO 200
415 ENDIF
416
417cc nnrpelem = nnrpelem + 1
418
419 j_stok = j_stok + 1
420 prov_n(j_stok) = jj
421 prov_e(j_stok) = ne
422 IF(j_stok == nvsiz)THEN
423
424 CALL i25sto(
425 1 nvsiz ,irect ,x ,nsv ,ii_stok,
426 2 cand_n,cand_e ,mulnsn,noint ,marge ,
427 3 i_mem ,prov_n ,prov_e,eshift,v ,
428 4 nsn ,nrtm ,gap_s ,gap_m ,curv_max,nin ,
429 5 pene_old,nbinflg ,mbinflg,ilev,msegtyp,
430 6 itab ,igap ,gap_s_l,gap_m_l,icodt,iskew,
431 7 drad ,dgapload)
432 IF(i_mem==2)GOTO 100
433 j_stok = 0
434 ENDIF
435
436 200 CONTINUE
437
438 jj = next_nod(jj)
439
440 ENDDO ! WHILE(JJ /= 0)
441
442 ENDDO
443 ENDDO
444 ENDDO
445 k = kremnod(2*(ne-1)+1)+1
446 l = kremnod(2*(ne-1)+2)
447 DO i=k,l
448 tag(remnod(i)) = 0
449 ENDDO
450
451cc nbpelg = nbpelg + nbpelem
452cc nnpelg = nnpelg + nnpelem
453cc nnrpelg = nnrpelg + nnrpelem
454cc nnr0pelg = nnr0pelg + nnr0pelem
455 ENDDO
456
457 ELSE !FLAGREMNODE
458 DO ne=1,nrtm
459C don't keep deleted segments
460 IF(stf(ne) <= zero)cycle
461c
462 aaa = marge+curv_max(ne)+max(max(bgapsmx+gap_m(ne),pmax_gap)+dgapload,drad)+vmaxdt
463 +
464
465
466c It is possible to improve the algo by cutting the facet
467c in 2 (4,3,6,9 ...) if the facet is large in front of AAA and inclinee
468
469 m1 = irect(1,ne)
470 m2 = irect(2,ne)
471 m3 = irect(3,ne)
472 m4 = irect(4,ne)
473
474 xx1=x(1,m1)
475 xx2=x(1,m2)
476 xx3=x(1,m3)
477 xx4=x(1,m4)
478 xmaxe=max(xx1,xx2,xx3,xx4)
479 xmine=min(xx1,xx2,xx3,xx4)
480
481 yy1=x(2,m1)
482 yy2=x(2,m2)
483 yy3=x(2,m3)
484 yy4=x(2,m4)
485 ymaxe=max(yy1,yy2,yy3,yy4)
486 ymine=min(yy1,yy2,yy3,yy4)
487
488 zz1=x(3,m1)
489 zz2=x(3,m2)
490 zz3=x(3,m3)
491 zz4=x(3,m4)
492 zmaxe=max(zz1,zz2,zz3,zz4)
493 zmine=min(zz1,zz2,zz3,zz4)
494
495
496c surface calculation (for future elimination of candidates)
497
498 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
499 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
500 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
501 s2 = sx*sx + sy*sy + sz*sz
502
503c index of voxels occupied by the facet
504
505 ix1=int(nbx*(xmine-aaa-xminb)/(xmaxb-xminb))
506 iy1=int(nby*(ymine-aaa-yminb)/(ymaxb-yminb))
507 iz1=int(nbz*(zmine-aaa-zminb)/(zmaxb-zminb))
508
509 ix1=max(1,2+min(nbx,ix1))
510 iy1=max(1,2+min(nby,iy1))
511 iz1=max(1,2+min(nbz,iz1))
512
513 ix2=int(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb))
514 iy2=int(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb))
515 iz2=int(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb))
516
517 ix2=max(1,2+min(nbx,ix2))
518 iy2=max(1,2+min(nby,iy2))
519 iz2=max(1,2+min(nbz,iz2))
520
521cc nbpelem = 0
522cc'nnpelem = 0
523cc'nnr0pelem = 0
524cc'nnrpelem = 0
525
526 DO iz = iz1,iz2
527 DO iy = iy1,iy2
528 DO ix = ix1,ix2
529
530cc nbpelem = nbpelem + 1
531
532 jj = voxel(ix,iy,iz)
533
534 DO WHILE(jj /= 0)
535
536cc nnpelem = nnpelem + 1
537
538 IF(jj<=nsn)THEN
539 nn=nsv(jj)
540 IF(nn == m1)GOTO 300
541 IF(nn == m2)GOTO 300
542 IF(nn == m3)GOTO 300
543 IF(nn == m4)GOTO 300
544 xs = x(1,nn)
545 ys = x(2,nn)
546 zs = x(3,nn)
547c PMAX_GAP is a global overestimate penetration
548c NEED to communicate in SPMD
549c VMAXDT is a local overestimate of relative incremental displacement
550c NO need to communicate in SPMD
551
552 aaa = marge + curv_max(ne)
553 + + max(gap_s(jj)+gap_m(ne)+dgapload,drad)
554 + +vmaxdt
555 ELSE
556 j=jj-nsn
557
558 xs = xrem(1,j)
559 ys = xrem(2,j)
560 zs = xrem(3,j)
561 aaa = marge+curv_max(ne)
562c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
563c +EDGE_L2(JJ) remote
564c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
565 + + max(xrem(igapxremp,j)+gap_m(ne)+dgapload,drad)
566 + + vmaxdt
567 ENDIF
568
569 IF(xs<=xmine-aaa)GOTO 300
570 IF(xs>=xmaxe+aaa)GOTO 300
571 IF(ys<=ymine-aaa)GOTO 300
572 IF(ys>=ymaxe+aaa)GOTO 300
573 IF(zs<=zmine-aaa)GOTO 300
574 IF(zs>=zmaxe+aaa)GOTO 300
575
576c underestimation of distance**2 for candidate elimination
577
578cc'nnr0pelem = nnr0pelem + 1
579
580 d1x = xs - xx1
581 d1y = ys - yy1
582 d1z = zs - zz1
583 d2x = xs - xx2
584 d2y = ys - yy2
585 d2z = zs - zz2
586 dd1 = d1x*sx+d1y*sy+d1z*sz
587 dd2 = d2x*sx+d2y*sy+d2z*sz
588 IF(dd1*dd2 > zero)THEN
589 d2 = min(dd1*dd1,dd2*dd2)
590 a2 = aaa*aaa*s2
591 IF(d2 > a2)GOTO 300
592 ENDIF
593
594cc nnrpelem = nnrpelem + 1
595
596 j_stok = j_stok + 1
597 prov_n(j_stok) = jj
598 prov_e(j_stok) = ne
599 IF(j_stok == nvsiz)THEN
600
601 CALL i25sto(
602 1 nvsiz ,irect ,x ,nsv ,ii_stok,
603 2 cand_n,cand_e ,mulnsn,noint ,marge ,
604 3 i_mem ,prov_n ,prov_e,eshift,v ,
605 4 nsn ,nrtm ,gap_s ,gap_m ,curv_max,nin ,
606 5 pene_old,nbinflg ,mbinflg,ilev,msegtyp,
607 6 itab ,igap ,gap_s_l,gap_m_l,icodt,iskew,
608 7 drad ,dgapload)
609 IF(i_mem==2)GOTO 100
610 j_stok = 0
611 ENDIF
612
613 300 CONTINUE
614
615 jj = next_nod(jj)
616
617 ENDDO ! WHILE(JJ /= 0)
618
619 ENDDO
620 ENDDO
621 ENDDO
622
623cc nbpelg = nbpelg + nbpelem
624cc nnpelg = nnpelg + nnpelem
625cc nnrpelg = nnrpelg + nnrpelem
626cc nnr0pelg = nnr0pelg + nnr0pelem
627 ENDDO
628 END IF !FLAGREMNODE
629C-------------------------------------------------------------------------
630C end of sorting
631C-------------------------------------------------------------------------
632 IF(j_stok/=0)CALL i25sto(
633 1 j_stok,irect ,x ,nsv ,ii_stok,
634 2 cand_n,cand_e ,mulnsn,noint ,marge ,
635 3 i_mem ,prov_n ,prov_e,eshift,v ,
636 4 nsn ,nrtm ,gap_s ,gap_m ,curv_max,nin ,
637 5 pene_old,nbinflg,mbinflg,ilev ,msegtyp,
638 6 itab ,igap ,gap_s_l,gap_m_l,icodt,iskew,
639 7 drad ,dgapload)
640
641C=======================================================================
642C 4 reset nodes to zero in boxes
643C=======================================================================
644 100 CONTINUE
645
646C Barrier to avoid reinitialization before end of sorting
647 CALL my_barrier
648 nsnf = 1 + itask*nsn / nthread
649 nsnl = (itask+1)*nsn / nthread
650
651 DO i=nsnf,nsnl
652 IF(iix(i)/=0)THEN
653 voxel(iix(i),iiy(i),iiz(i))=0
654 ENDIF
655 ENDDO
656C=======================================================================
657C 5 reset nodes to zero in boxes
658C non -local candidates in SPMD
659C=======================================================================
660 nsnf = 1 + itask*nsnr / nthread
661 nsnl = (itask+1)*nsnr / nthread
662 DO j = nsnf, nsnl
663 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))=0
664 ENDDO
665
666C
667 CALL my_barrier()
668 IF(itask == 0)THEN
669 DEALLOCATE(next_nod)
670 DEALLOCATE(iix)
671 DEALLOCATE(iiy)
672 DEALLOCATE(iiz)
673 ENDIF
674
675 DEALLOCATE(tag)
676 DEALLOCATE(last_nod)
677 RETURN
678 END
679
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:274
#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