OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i7trivox.F File Reference
#include "implicit_f.inc"
#include "mvsiz_p.inc"
#include "com01_c.inc"
#include "com04_c.inc"
#include "param_c.inc"
#include "task_c.inc"
#include "ige3d_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine i7trivox (nsn, renum, nsnr, isznsnr, i_mem, irect, x, stf, stfn, xyzm, nsv, ii_stok, cand_n, eshift, cand_e, mulnsn, noint, tzinf, gap_s_l, gap_m_l, voxel, nbx, nby, nbz, intth, inacti, ifq, cand_a, cand_p, ifpen, nrtm, nsnrold, igap, gap, gap_s, gap_m, gapmin, gapmax, marge, curv_max, nin, itask, bgapsmx, kremnod, remnod, itab, flagremnode, drad, itied, cand_f, dgapload, remote_s_node, list_remote_s_node, total_nb_nrtm, intheat, idt_therm, nodadt_therm)

Function/Subroutine Documentation

◆ i7trivox()

subroutine i7trivox ( integer nsn,
integer, dimension(*) renum,
integer nsnr,
integer isznsnr,
integer i_mem,
integer, dimension(4,*) irect,
x,
stf,
stfn,
xyzm,
integer, dimension(*) nsv,
integer ii_stok,
integer, dimension(*) cand_n,
integer eshift,
integer, dimension(*) cand_e,
integer mulnsn,
integer noint,
tzinf,
gap_s_l,
gap_m_l,
integer, dimension(nbx+2,nby+2,nbz+2) voxel,
integer nbx,
integer nby,
integer nbz,
integer intth,
integer inacti,
integer ifq,
integer, dimension(*) cand_a,
cand_p,
integer, dimension(*) ifpen,
integer, intent(in) nrtm,
integer nsnrold,
integer igap,
gap,
gap_s,
gap_m,
gapmin,
gapmax,
marge,
curv_max,
integer nin,
integer itask,
bgapsmx,
integer, dimension(*) kremnod,
integer, dimension(*) remnod,
integer, dimension(*) itab,
integer flagremnode,
intent(in) drad,
integer itied,
cand_f,
intent(in) dgapload,
integer, intent(inout) remote_s_node,
integer, dimension(nsnr), intent(inout) list_remote_s_node,
integer, intent(in) total_nb_nrtm,
integer, intent(in) intheat,
integer, intent(in) idt_therm,
integer, intent(in) nodadt_therm )
Parameters
[in]nrtmnumber of segments per threads
[in]total_nb_nrtmtotal number of segments

Definition at line 37 of file i7trivox.F.

50C============================================================================
51C M o d u l e s
52C-----------------------------------------------
53 USE message_mod
54 USE tri7box
55C-----------------------------------------------
56C I m p l i c i t T y p e s
57C-----------------------------------------------
58#include "implicit_f.inc"
59C-----------------------------------------------
60C G l o b a l P a r a m e t e r s
61C-----------------------------------------------
62#include "mvsiz_p.inc"
63c parameter setting the size for the vector (orig version is 128)
64 INTEGER NVECSZ
65 parameter(nvecsz = mvsiz)
66C-----------------------------------------------
67C C o m m o n B l o c k s
68C-----------------------------------------------
69#include "com01_c.inc"
70#include "com04_c.inc"
71#include "param_c.inc"
72#include "task_c.inc"
73#include "ige3d_c.inc"
74C-----------------------------------------------
75C ROUTINE PURPOSE:
76C ================
77C CLASSIFIES NODES INTO BOXES
78C SEARCHES FOR BOXES CONCERNED FOR EACH FACET
79C SEARCHES FOR CANDIDATES
80C-----------------------------------------------
81C D u m m y A r g u m e n t s
82C
83C NAME DESCRIPTION I/O
84C
85C IRECT(4,*) FACET CONNECTIVITY ARRAY I
86C X(3,*) NODAL COORDINATES I
87C NSV NODE NUMBERS I
88C XMAX largest existing abscissa I
89C XMAX largest existing ordinate I
90C XMAX largest existing coordinate I
91C I_STOK storage level of candidate
92C impact couples I/O
93C CAND_N node result boxes
94C CAND_E element result box addresses
95C MULNSN = MULTIMP*NSN MAX SIZE ALLOWED NOW FOR
96C NODE,ELEMENT CANDIDATE COUPLES
97C NOINT USER INTERFACE NUMBER
98C TZINF INFLUENCE ZONE SIZE
99C
100C PROV_N provisional CAND_N (static variable in i7tri)
101C PROV_E provisional CAND_E (static variable in i7tri)
102
103C VOXEL(ix,iy,iz) contains the local number of the first node in
104C the box
105C NEXT_NOD(i) next node in the same box (if /= 0)
106C LAST_NOD(i) last node in the same box (if /= 0)
107C used only to go directly from first
108C node to last
109C-----------------------------------------------
110C D u m m y A r g u m e n t s
111C-----------------------------------------------
112 INTEGER I_MEM,ESHIFT,NSN,ISZNSNR,NSNROLD,NIN,ITASK,
113 . MULNSN,NOINT,INACTI,IFQ,NSNR,IGAP,NBX,NBY,NBZ,
114 . NSV(*),CAND_N(*),CAND_E(*),CAND_A(*),IFPEN(*),RENUM(*),
115 . INTTH,IRECT(4,*), VOXEL(NBX+2,NBY+2,NBZ+2),II_STOK,
116 . KREMNOD(*),REMNOD(*),ITAB(*),FLAGREMNODE,ITIED
117 INTEGER, INTENT(in) :: NRTM !< number of segments per threads
118 INTEGER, INTENT(in) :: TOTAL_NB_NRTM !< total number of segments
119 INTEGER, INTENT(IN) :: INTHEAT
120 INTEGER, INTENT(IN) :: IDT_THERM
121 INTEGER, INTENT(IN) :: NODADT_THERM
122 my_real
123 . x(3,*),xyzm(12),cand_p(*),stf(*),stfn(*),gap_s(*),gap_m(*),
124 . tzinf,marge,gap,gapmin,gapmax,bgapsmx,
125 . curv_max(*),gap_s_l(*),gap_m_l(*),cand_f(*)
126 my_real , INTENT(IN) :: drad,dgapload
127 INTEGER, INTENT(inout) :: REMOTE_S_NODE
128 INTEGER, DIMENSION(NSNR), INTENT(inout) :: LIST_REMOTE_S_NODE
129C-----------------------------------------------
130C L o c a l V a r i a b l e s
131C-----------------------------------------------
132 INTEGER NB_NCN,NB_NCN1,NB_ECN,I,J,DIR,NB_NC,NB_EC,
133 . N1,N2,N3,N4,NN,NE,K,L,NCAND_PROV,J_STOK,II,JJ,
134 . PROV_N(MVSIZ),PROV_E(MVSIZ),
135 . OLDNUM(ISZNSNR), NSNF, NSNL,DELNOD,M
136 INTEGER, DIMENSION(:), ALLOCATABLE :: TAGREMNODE
137 my_real
138 . dx,dy,dz,xs,ys,zs,xx,sx,sy,sz,s2,
139 . xmin, xmax,ymin, ymax,zmin, zmax, tz, gapsmx, gapl,
140 . xx1,xx2,xx3,xx4,yy1,yy2,yy3,yy4,zz1,zz2,zz3,zz4,
141 . d1x,d1y,d1z,d2x,d2y,d2z,dd1,dd2,d2,a2,gs
142c temporary
143 INTEGER LAST_NOD(NSN+NSNR)
144 INTEGER IX,IY,IZ,NEXT,M1,M2,M3,M4,
145 . IX1,IY1,IZ1,IX2,IY2,IZ2
146 INTEGER, DIMENSION(:),ALLOCATABLE :: IIX,IIY,IIZ
147 my_real
148 . xminb,yminb,zminb,xmaxb,ymaxb,zmaxb,
149 . xmine,ymine,zmine,xmaxe,ymaxe,zmaxe,aaa
150 INTEGER CPT_VOX0
151 LOGICAL TEST_V0
152 INTEGER FIRST,NEW,LAST,IERROR
153 LOGICAL DBG_type18_fvm
154 SAVE iix,iiy,iiz
155C-----------------------------------------------
156 IF(itask == 0)THEN
157 remote_s_node = 0
158 cpt_vox0 = 0
159 ALLOCATE(next_nod(nsn+nsnr),stat=ierror)
160 IF(ierror/=0) THEN
161 CALL ancmsg(msgid=19,anmode=aninfo,
162 . c1='(/INTER/TYPE7)')
163 CALL arret(2)
164 ENDIF
165 ALLOCATE(iix(nsn+nsnr),iiy(nsn+nsnr),iiz(nsn+nsnr),stat=ierror)
166 IF(ierror/=0) THEN
167 CALL ancmsg(msgid=19,anmode=aninfo,
168 . c1='(/INTER/TYPE7)')
169 CALL arret(2)
170 ENDIF
171 END IF
172
173 CALL my_barrier !Barrier to wait init voxel and allocation NEXT_NOD
174
175! Initial phase of BPE and BPN construction moved from I7BUCE => I7TRI
176
177 xmin = xyzm(1)
178 ymin = xyzm(2)
179 zmin = xyzm(3)
180 xmax = xyzm(4)
181 ymax = xyzm(5)
182 zmax = xyzm(6)
183
184 !dev future: xminb larger than xmin...
185 IF( to_trim(nin) ) THEN
186 xminb = xyzm(7)
187 yminb = xyzm(8)
188 zminb = xyzm(9)
189 xmaxb = xyzm(10)
190 ymaxb = xyzm(11)
191 zmaxb = xyzm(12)
192 ELSE
193! we use the global bounding box
194! the box reduction seemed inefficient
195! at the first sorting
196! because too many SECONDARY nodes were outside
197! the reduced box
198 xminb = xmin
199 yminb = ymin
200 zminb = zmin
201 xmaxb = xmax
202 ymaxb = ymax
203 zmaxb = zmax
204 ENDIF
205
206 ! In SPMD, for inacti or IFQ, retrieves old numbering of non-local candidates
207 IF(nspmd>1.AND.(inacti==5.OR.inacti==6.OR.inacti==7.OR.ifq>0.OR.itied/=0)) THEN
208 CALL spmd_oldnumcd(renum,oldnum,isznsnr,nsnrold,intheat,idt_therm,nodadt_therm)
209 END IF
210
211C=======================================================================
212C 1 NODE CLUSTERING
213C=======================================================================
214 IF(itask==0.AND.total_nb_nrtm>0)THEN
215 !==========LAGRANGIAN and STAGGERED ALE==========!
216 test_v0 = first_test(nin)
217 DO i=1,nsn
218 iix(i)=0
219 iiy(i)=0
220 iiz(i)=0
221 IF(stfn(i) == zero)cycle
222 j=nsv(i)
223 !C Optimization // search for nodes within xmin xmax of
224 !C processor elements
225 IF(x(1,j) < xmin) cycle
226 IF(x(1,j) > xmax) cycle
227 IF(x(2,j) < ymin) cycle
228 IF(x(2,j) > ymax) cycle
229 IF(x(3,j) < zmin) cycle
230 IF(x(3,j) > zmax) cycle
231 iix(i)=int(nbx*(x(1,j)-xminb)/(xmaxb-xminb))
232 iiy(i)=int(nby*(x(2,j)-yminb)/(ymaxb-yminb))
233 iiz(i)=int(nbz*(x(3,j)-zminb)/(zmaxb-zminb))
234 iix(i)=max(1,2+min(nbx,iix(i)))
235 iiy(i)=max(1,2+min(nby,iiy(i)))
236 iiz(i)=max(1,2+min(nbz,iiz(i)))
237 first = voxel(iix(i),iiy(i),iiz(i))
238 IF(test_v0) THEN
239! count the number of SECONDARY nodes outiside the reduced box
240 IF(iix(i) == 1 .OR. iiy(i) == 1 .OR. iiz(i) == 1 .AND.
241 . iix(i) == nbx+2 .OR. iiy(i) == nby+2 .OR. iiz(i) == nbz+2) THEN
242 cpt_vox0 = cpt_vox0 +1
243 ENDIF
244 ENDIF
245 IF(first == 0)THEN
246 !empty cell
247 voxel(iix(i),iiy(i),iiz(i)) = i ! first
248 next_nod(i) = 0 ! last one
249 last_nod(i) = 0 ! no last
250 ELSEIF(last_nod(first) == 0)THEN
251 !cell containing one node
252 !add as next node
253 next_nod(first) = i ! next
254 last_nod(first) = i ! last
255 next_nod(i) = 0 ! last one
256 ELSE
257 !
258 !jump to the last node of the cell
259 last = last_nod(first) ! last node in this voxel
260 next_nod(last) = i ! next
261 last_nod(first) = i ! last
262 next_nod(i) = 0 ! last one
263 ENDIF
264 ENDDO
265
266C=======================================================================
267C 2 REMOTE NODE CLUSTERING
268C=======================================================================
269 DO j = 1, nsnr
270
271 IF(xrem(1,j) < xmin) cycle
272 IF(xrem(1,j) > xmax) cycle
273 IF(xrem(2,j) < ymin) cycle
274 IF(xrem(2,j) > ymax) cycle
275 IF(xrem(3,j) < zmin) cycle
276 IF(xrem(3,j) > zmax) cycle
277
278 remote_s_node = remote_s_node + 1
279 list_remote_s_node( remote_s_node ) = j
280 iix(nsn+j)=int(nbx*(xrem(1,j)-xminb)/(xmaxb-xminb))
281 iiy(nsn+j)=int(nby*(xrem(2,j)-yminb)/(ymaxb-yminb))
282 iiz(nsn+j)=int(nbz*(xrem(3,j)-zminb)/(zmaxb-zminb))
283 iix(nsn+j)=max(1,2+min(nbx,iix(nsn+j)))
284 iiy(nsn+j)=max(1,2+min(nby,iiy(nsn+j)))
285 iiz(nsn+j)=max(1,2+min(nbz,iiz(nsn+j)))
286
287 first = voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))
288 IF(test_v0) THEN
289! count the number of SECONDARY nodes outiside the reduced box
290 IF(iix(j+nsn) == 1 .OR. iiy(j+nsn) == 1 .OR. iiz(j+nsn) == 1 .AND.
291 . iix(j+nsn) == nbx+2 .OR. iiy(j+nsn) == nby+2 .OR. iiz(j+nsn) == nbz+2) THEN
292 cpt_vox0 = cpt_vox0 +1
293 ENDIF
294 ENDIF
295
296 IF(first == 0)THEN
297 ! empty cell
298 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j)) = nsn+j ! first
299 next_nod(nsn+j) = 0 ! last one
300 last_nod(nsn+j) = 0 ! no last
301 ELSEIF(last_nod(first) == 0)THEN
302 ! cell containing one node, add it as next node
303 next_nod(first) = nsn+j ! next
304 last_nod(first) = nsn+j ! last
305 next_nod(nsn+j) = 0 ! last one
306 ELSE
307 ! , jump to the last node of the cell
308 last = last_nod(first) ! last node in this voxel
309 next_nod(last) = nsn+j ! next
310 last_nod(first) = nsn+j ! last
311 next_nod(nsn+j) = 0 ! last one
312 ENDIF
313 ENDDO
314 END IF !ITASK == 0
315
316 CALL my_barrier ! Barrier to wait task0 treatment
317
318 IF(first_test(nin)) THEN
319!During the first sorting, we use the reduced box
320!If too many (>5 %) nodes are outside the reduced box
321!then the box is not reduced anymore
322 CALL my_barrier
323 IF(itask == 0) THEN
324 IF(cpt_vox0 > 5*(remote_s_node + nsn)/100) to_trim(nin) = .false.
325 first_test(nin) = .false.
326 ENDIF
327 CALL my_barrier
328 ENDIF
329
330
331C=======================================================================
332C 3 FACE RECOVERY AND ENUMERATION OF CANDIDATE COUPLES
333C=======================================================================
334 j_stok = 0
335 IF(flagremnode == 2) THEN
336 ALLOCATE(tagremnode(numnod+numfakenodigeo))
337 DO i=1,numnod+numfakenodigeo
338 tagremnode(i) = 0
339 ENDDO
340 ENDIF
341 DO ne=1,nrtm
342 IF(stf(ne) == zero)cycle ! do not keep destroyed facets
343 IF(flagremnode == 2) THEN
344 k = kremnod(2*(ne-1)+1)+1
345 l = kremnod(2*(ne-1)+2)
346 DO i=k,l
347 tagremnode(remnod(i)) = 1
348 ENDDO
349 ENDIF
350 IF(igap == 0)THEN
351 aaa = tzinf+curv_max(ne)
352 ELSE
353 aaa = marge+curv_max(ne)+max(min(gapmax,max(gapmin,bgapsmx+gap_m(ne)))+dgapload,drad)
354 ENDIF
355 !it is possible to improve the algorithm by subdividing the facet
356 !into 2(4,3,6,9...) if the facet is large compared to AAA and inclined
357 m1 = irect(1,ne)
358 m2 = irect(2,ne)
359 m3 = irect(3,ne)
360 m4 = irect(4,ne)
361
362 xx1=x(1,m1)
363 xx2=x(1,m2)
364 xx3=x(1,m3)
365 xx4=x(1,m4)
366 xmaxe=max(xx1,xx2,xx3,xx4)
367 xmine=min(xx1,xx2,xx3,xx4)
368
369 yy1=x(2,m1)
370 yy2=x(2,m2)
371 yy3=x(2,m3)
372 yy4=x(2,m4)
373 ymaxe=max(yy1,yy2,yy3,yy4)
374 ymine=min(yy1,yy2,yy3,yy4)
375
376 zz1=x(3,m1)
377 zz2=x(3,m2)
378 zz3=x(3,m3)
379 zz4=x(3,m4)
380 zmaxe=max(zz1,zz2,zz3,zz4)
381 zmine=min(zz1,zz2,zz3,zz4)
382
383 !surface calculation (for future candidate elimination)
384 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
385 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
386 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
387 s2 = sx*sx + sy*sy + sz*sz
388
389 !index of voxels occupied by the facet
390 IF(nbx>1) THEN
391 ix1=int(nbx*(xmine-aaa-xminb)/(xmaxb-xminb))
392 ix2=int(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb))
393 ELSE
394 ix1=-2
395 ix2=1
396 ENDIF
397
398 IF(nby>1) THEN
399 iy1=int(nby*(ymine-aaa-yminb)/(ymaxb-yminb))
400 iy2=int(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb))
401 ELSE
402 iy1=-2
403 iy2=1
404 ENDIF
405
406 IF(nbz>1) THEN
407 iz1=int(nbz*(zmine-aaa-zminb)/(zmaxb-zminb))
408 iz2=int(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb))
409 ELSE
410 iz1=-2
411 iz2=1
412 ENDIF
413
414 ix1=max(1,2+min(nbx,ix1))
415 iy1=max(1,2+min(nby,iy1))
416 iz1=max(1,2+min(nbz,iz1))
417
418 ix2=max(1,2+min(nbx,ix2))
419 iy2=max(1,2+min(nby,iy2))
420 iz2=max(1,2+min(nbz,iz2))
421
422 !nbpelem = 0
423 !nnpelem = 0
424 !nnr0pelem = 0
425 !nnrpelem = 0
426
427 DO iz = iz1,iz2
428 DO iy = iy1,iy2
429 DO ix = ix1,ix2
430
431 ! nbpelem = nbpelem + 1
432
433 jj = voxel(ix,iy,iz)
434 DO WHILE(jj /= 0)
435 delnod = 0
436
437 !nnpelem = nnpelem + 1
438 IF(jj<=nsn)THEN
439 nn=nsv(jj)
440
441 IF(nn == m1)GOTO 200
442 IF(nn == m2)GOTO 200
443 IF(nn == m3)GOTO 200
444 IF(nn == m4)GOTO 200
445
446 IF(flagremnode == 2) THEN
447 IF( tagremnode(nsv(jj)) == 1)GOTO 200
448 ENDIF
449 xs = x(1,nn)
450 ys = x(2,nn)
451 zs = x(3,nn)
452 IF(igap /= 0)THEN
453 aaa = marge+curv_max(ne)+max(min(gapmax,max(gapmin,gap_s(jj)+gap_m(ne)))+dgapload,drad)
454 ENDIF
455 ELSE
456 j=jj-nsn
457 IF(flagremnode == 2) THEN
458 nn = irem(1,j)
459 k = kremnod(2*(ne-1)+2) + 1
460 l = kremnod(2*(ne-1)+3)
461 DO m=k,l
462 IF(remnod(m) == -irem(2,j) ) THEN
463 delnod = delnod + 1
464 EXIT
465 ENDIF
466 ENDDO
467 IF(delnod /= 0)GOTO 200
468 ENDIF
469
470 xs = xrem(1,j)
471 ys = xrem(2,j)
472 zs = xrem(3,j)
473 IF(igap /= 0)THEN
474 aaa = marge+curv_max(ne)+max(min(gapmax,max(gapmin,xrem(9,j)+gap_m(ne)))+dgapload,drad)
475 ENDIF
476 ENDIF
477
478 IF(xs<=xmine-aaa)GOTO 200
479 IF(xs>=xmaxe+aaa)GOTO 200
480 IF(ys<=ymine-aaa)GOTO 200
481 IF(ys>=ymaxe+aaa)GOTO 200
482 IF(zs<=zmine-aaa)GOTO 200
483 IF(zs>=zmaxe+aaa)GOTO 200
484
485 !underestimation of distance**2 for candidate elimination
486 !nnr0pelem = nnr0pelem + 1
487
488 d1x = xs - xx1
489 d1y = ys - yy1
490 d1z = zs - zz1
491 d2x = xs - xx2
492 d2y = ys - yy2
493 d2z = zs - zz2
494 dd1 = d1x*sx+d1y*sy+d1z*sz
495 dd2 = d2x*sx+d2y*sy+d2z*sz
496 IF(dd1*dd2 > zero)THEN
497 d2 = min(dd1*dd1,dd2*dd2)
498 a2 = aaa*aaa*s2
499 IF(d2 > a2)GOTO 200
500 ENDIF
501
502 !nnrpelem = nnrpelem + 1
503
504 j_stok = j_stok + 1
505 prov_n(j_stok) = jj
506 prov_e(j_stok) = ne
507 IF(j_stok == nvsiz)THEN
508 CALL i7sto(
509 1 nvsiz ,irect ,x ,nsv ,ii_stok,
510 2 cand_n,cand_e ,mulnsn,noint ,marge ,
511 3 i_mem ,prov_n ,prov_e,eshift,inacti ,
512 4 ifq ,cand_a ,cand_p,ifpen ,nsn ,
513 5 oldnum,nsnrold,igap ,gap ,gap_s ,
514 6 gap_m ,gapmin ,gapmax,curv_max,nin ,
515 7 gap_s_l,gap_m_l,intth,drad,itied ,
516 8 cand_f ,dgapload)
517 IF(i_mem==2)GOTO 100
518 j_stok = 0
519 ENDIF
520 200 CONTINUE
521 jj = next_nod(jj)
522 ENDDO ! WHILE(JJ /= 0)
523 ENDDO ! X
524 ENDDO ! Y
525 ENDDO ! Z
526 IF(flagremnode == 2) THEN
527 k = kremnod(2*(ne-1)+1)+1
528 l = kremnod(2*(ne-1)+2)
529 DO i=k,l
530 tagremnode(remnod(i)) = 0
531 ENDDO
532 ENDIF
533 ENDDO
534C=======================================================================
535C 4 STORAGE
536C=======================================================================
537 IF(j_stok/=0)CALL i7sto(
538 1 j_stok,irect ,x ,nsv ,ii_stok,
539 2 cand_n,cand_e ,mulnsn,noint ,marge ,
540 3 i_mem ,prov_n ,prov_e,eshift,inacti ,
541 4 ifq ,cand_a ,cand_p,ifpen ,nsn ,
542 5 oldnum,nsnrold,igap ,gap ,gap_s ,
543 6 gap_m ,gapmin ,gapmax,curv_max,nin ,
544 7 gap_s_l,gap_m_l,intth,drad ,itied ,
545 8 cand_f ,dgapload)
546
547C=======================================================================
548C 5 VOXEL CLUSTERING RESET
549C=======================================================================
550 100 CONTINUE
551 ! Barrier to avoid reinitialization before end of sorting
552 CALL my_barrier
553 IF(total_nb_nrtm>0) THEN
554 nsnf = 1 + itask*nsn / nthread
555 nsnl = (itask+1)*nsn / nthread
556 DO i=nsnf,nsnl
557 IF(iix(i)/=0)THEN
558 voxel(iix(i),iiy(i),iiz(i))=0
559 ENDIF
560 ENDDO
561 nsnf = 1 + itask*remote_s_node / nthread
562 nsnl = (itask+1)*remote_s_node / nthread
563 IF(itask+1==nthread) nsnl=remote_s_node
564 DO jj = nsnf, nsnl
565 j = list_remote_s_node(jj)
566 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))=0
567 ENDDO
568 ENDIF
569 CALL my_barrier()
570
571C=======================================================================
572C 6 DEBUG OUTPUT
573C=======================================================================
574
575 dbg_type18_fvm=.false.
576 if(inacti==7 .AND. dbg_type18_fvm)then
577 write(*,fmt='(A)')"------------------------------------------"
578 write(*,*)"RESULT : Search Algorithm with VOXEL partitioning"
579 write(*,*)" Number of couples =", ii_stok
580 if(ii_stok>0)then
581 write(*,fmt='(A,(I10))')" --> SECONDARY Node ids: ", cand_n(1:ii_stok)
582 write(*,fmt='(A,(I10))')" --> Local Face ids: ", cand_e(1:ii_stok)
583 endif
584 write(*,*)" Structure domain :"
585 write(*,fmt='(A,F30.16,A,F30.16)')" Xmin=",xmin," Xmax=",xmax
586 write(*,fmt='(A,F30.16,A,F30.16)')" Ymin=",ymin," Ymax=",ymax
587 write(*,fmt='(A,F30.16,A,F30.16)')" Zmin=",zmin," Zmax=",zmax
588 write(*,*)" Partitioning domain :"
589 write(*,*)" TZINF,AAA=",tzinf,aaa
590 write(*,fmt='(A,F30.16,A,F30.16)')" Xmin=",xmin-aaa," Xmax=",xmax+aaa
591 write(*,fmt='(A,F30.16,A,F30.16)')" Ymin=",ymin-aaa," Ymax=",ymax+aaa
592 write(*,fmt='(A,F30.16,A,F30.16)')" Zmin=",zmin-aaa," Zmax=",zmax+aaa
593 write(*,fmt='(A)')"------------------------------------------"
594 endif
595
596C=======================================================================
597C 7 DEALLOCATE
598C=======================================================================
599 IF(itask == 0)THEN
600 DEALLOCATE(next_nod)
601 DEALLOCATE(iix)
602 DEALLOCATE(iiy)
603 DEALLOCATE(iiz)
604 ENDIF
605 IF(flagremnode == 2) THEN
606 IF(ALLOCATED(tagremnode)) DEALLOCATE(tagremnode)
607 ENDIF
608
609
610 RETURN
#define my_real
Definition cppsort.cpp:32
subroutine i7sto(j_stok, irect, x, nsv, ii_stok, cand_n, cand_e, mulnsn, noint, marge, i_mem, prov_n, prov_e, eshift, inacti, ifq, cand_a, cand_p, ifpen, nsn, oldnum, nsnrold, igap, gap, gap_s, gap_m, gapmin, gapmax, curv_max, nin, gap_s_l, gap_m_l, intth, drad, itied, cand_f, dgapload)
Definition i7sto.F:41
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
logical, dimension(:), allocatable first_test
Definition tri7box.F:420
integer, dimension(:), allocatable next_nod
Definition tri7box.F:48
logical, dimension(:), allocatable to_trim
Definition tri7box.F:420
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 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
subroutine arret(nn)
Definition arret.F:87
subroutine my_barrier
Definition machine.F:31