OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i23trivox.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!|| i23trivox ../engine/source/interfaces/intsort/i23trivox.F
25!||--- called by ------------------------------------------------------
26!|| i23buce ../engine/source/interfaces/intsort/i23buce.F
27!||--- calls -----------------------------------------------------
28!|| i23sto ../engine/source/interfaces/int23/i23sto.F
29!|| my_barrier ../engine/source/system/machine.F
30!|| spmd_oldnumcd ../engine/source/mpi/interfaces/spmd_i7tool.F
31!||--- uses -----------------------------------------------------
32!|| tri7box ../engine/share/modules/tri7box.F
33!||====================================================================
34 SUBROUTINE i23trivox(
35 1 NSN ,RENUM ,NSNR ,ISZNSNR ,I_MEM ,
36 2 IRECT ,X ,STF ,STFN ,XYZM ,
37 3 NSV ,II_STOK,CAND_N ,ESHIFT ,CAND_E ,
38 4 MULNSN ,NOINT ,TZINF ,MSR ,
39 5 VOXEL ,NBX ,NBY ,NBZ ,
40 6 INACTI ,CAND_A ,CAND_P ,IFPEN ,
41 7 NRTM ,NSNROLD,IGAP ,GAP ,GAP_S ,
42 8 GAP_M ,GAPMIN ,GAPMAX ,MARGE ,CURV_MAX ,
43 9 NIN ,ITASK ,BGAPSMX ,INTHEAT,IDT_THERM,NODADT_THERM)
44C============================================================================
45C M o d u l e s
46C-----------------------------------------------
47 USE tri7box
48C-----------------------------------------------
49C I m p l i c i t T y p e s
50C-----------------------------------------------
51#include "implicit_f.inc"
52C-----------------------------------------------
53C G l o b a l P a r a m e t e r s
54C-----------------------------------------------
55#include "mvsiz_p.inc"
56c parameter setting the size for the vector (orig version is 128)
57 INTEGER NVECSZ
58 PARAMETER (NVECSZ = mvsiz)
59C-----------------------------------------------
60C C o m m o n B l o c k s
61C-----------------------------------------------
62#include "com01_c.inc"
63#include "param_c.inc"
64#include "task_c.inc"
65C-----------------------------------------------
66C role of the routine:
67C ===================
68C classifies the nodes in boxes
69C searches for each facet the concerned boxes
70C searches for candidates
71C-----------------------------------------------
72C D u m m y A r g u m e n t s
73C
74C NOM DESCRIPTION E/S
75C
76C IRECT(4,*) ARRAY OF CONEC FACETTES E
77C X(3,*) COORDONNEES NODALES E
78C NSV NOS SYSTEMES DES NODES E
79C Xmax larger abcisse existing e
80C XMAX largest order.existing E
81C Xmax larger existing side E
82C I_STOK storage level of pairs
83C CANDIDATES impact E/S
84C CAND_N boites resultats nodes C CAND_E adresses des boites resultat elements
85C MULNSN = MULTIMP*NSN max size admitted now for the
86C COUPLES NODES,ELT CANDIDATES
87C NOINT INTERFACE USER NUMBER
88C TZINF TAILLE ZONE INFLUENCE
89C
90C Prov_n Provisional Cand_n (static variable in i7tri)
91C PROV_E CAND_E provisoire (variable static in i7tri)
92
93C VOXEL(ix,iy,iz) contains the local number of the first node of
94C the box
95C Next_nod (i) Next node in the same box (if /= 0)
96C Last_nod (i) Last node in the same box (if /= 0)
97C used only to go directly from the first
98C node at the last
99C-----------------------------------------------
100C D u m m y A r g u m e n t s
101C-----------------------------------------------
102 INTEGER I_MEM,ESHIFT,NSN,ISZNSNR,NRTM,NSNROLD,NIN,ITASK,
103 . MULNSN,NOINT,INACTI,NSNR,IGAP,NBX,NBY,NBZ,
104 . NSV(*),CAND_N(*),CAND_E(*),CAND_A(*),IFPEN(*),RENUM(*),
105 . IRECT(4,*), VOXEL(NBX+2,NBY+2,NBZ+2), MSR(*),II_STOK
106 INTEGER, INTENT(IN) :: INTHEAT
107 INTEGER, INTENT(IN) :: IDT_THERM
108 INTEGER, INTENT(IN) :: NODADT_THERM
109C REAL
110 my_real
111 . x(3,*),xyzm(6),cand_p(*),stf(*),stfn(*),gap_s(*),gap_m(*),
112 . tzinf,marge,gap,gapmin,gapmax,bgapsmx,
113 . curv_max(*)
114C-----------------------------------------------
115C L o c a l V a r i a b l e s
116C-----------------------------------------------
117 INTEGER I,J,
118 . NN,NE,J_STOK,JJ,
119 . prov_n(mvsiz),prov_e(mvsiz),
120 . oldnum(isznsnr), nsnf, nsnl
121C REAL
122 my_real
123 . xs,ys,zs,sx,sy,sz,s2,
124 . xmin, xmax,ymin, ymax,zmin, zmax,
125 . xx1,xx2,xx3,xx4,yy1,yy2,yy3,yy4,zz1,zz2,zz3,zz4,
126 . d1x,d1y,d1z,d2x,d2y,d2z,dd1,dd2,d2,a2
127c provisional
128 INTEGER LAST_NOD(NSN+NSNR)
129 INTEGER IX,IY,IZ,M1,M2,M3,M4,
130 . IX1,IY1,IZ1,IX2,IY2,IZ2
131 INTEGER, DIMENSION(:),ALLOCATABLE :: IIX,IIY,IIZ
132 my_real
133 . XMINB,YMINB,ZMINB,XMAXB,YMAXB,ZMAXB,
134 . XMINE,YMINE,ZMINE,XMAXE,YMAXE,ZMAXE,AAA
135 INTEGER FIRST,LAST
136 SAVE IIX,IIY,IIZ
137C-----------------------------------------------
138 IF(itask == 0)THEN
139 ALLOCATE(next_nod(nsn+nsnr))
140 ALLOCATE(iix(nsn+nsnr))
141 ALLOCATE(iiy(nsn+nsnr))
142 ALLOCATE(iiz(nsn+nsnr))
143 END IF
144C Barrier to wait init voxel and allocation NEX_NOD
145 CALL my_barrier
146C initial phase of BPE and BPN construction moved from I7BUCE => I7TRI
147C
148 xmin = xyzm(1)
149 ymin = xyzm(2)
150 zmin = xyzm(3)
151 xmax = xyzm(4)
152 ymax = xyzm(5)
153 zmax = xyzm(6)
154
155c Dev Future: xminb larger than Xmin ...
156 xminb = xmin
157 yminb = ymin
158 zminb = zmin
159 xmaxb = xmax
160 ymaxb = ymax
161 zmaxb = zmax
162
163
164C En SPMD, retrouve ancienne numerotation des CANDIDATES non local
165 IF(nspmd>1) THEN
166 CALL spmd_oldnumcd(renum,oldnum,isznsnr,nsnrold,intheat,idt_therm,nodadt_therm)
167 END IF
168
169C=======================================================================
170C 1 placing nodes in boxes
171C=======================================================================
172 IF(itask == 0)THEN
173cc BGAPSMX = ZERO
174 DO i=1,nsn
175 iix(i)=0
176 iiy(i)=0
177 iiz(i)=0
178 IF(stfn(i) == zero)cycle
179 j=nsv(i)
180C optimization // searches for nodes included in xmin xmax of
181C processor elements
182 IF(x(1,j) < xmin) cycle
183 IF(x(1,j) > xmax) cycle
184 IF(x(2,j) < ymin) cycle
185 IF(x(2,j) > ymax) cycle
186 IF(x(3,j) < zmin) cycle
187 IF(x(3,j) > zmax) cycle
188
189 iix(i)=int(nbx*(x(1,j)-xminb)/(xmaxb-xminb))
190 iiy(i)=int(nby*(x(2,j)-yminb)/(ymaxb-yminb))
191 iiz(i)=int(nbz*(x(3,j)-zminb)/(zmaxb-zminb))
192
193 iix(i)=max(1,2+min(nbx,iix(i)))
194 iiy(i)=max(1,2+min(nby,iiy(i)))
195 iiz(i)=max(1,2+min(nbz,iiz(i)))
196
197 first = voxel(iix(i),iiy(i),iiz(i))
198 IF(first == 0)THEN
199c Empty Cell
200 voxel(iix(i),iiy(i),iiz(i)) = i ! first
201 next_nod(i) = 0 ! last one
202 last_nod(i) = 0 ! no last
203 ELSEIF(last_nod(first) == 0)THEN
204c cell containing one node
205c add as next node
206 next_nod(first) = i ! next
207 last_nod(first) = i ! last
208 next_nod(i) = 0 ! last one
209 ELSE
210c
211c jump to the last node of the cell
212 last = last_nod(first) ! last node in this voxel
213 next_nod(last) = i ! next
214 last_nod(first) = i ! last
215 next_nod(i) = 0 ! last one
216 ENDIF
217 ENDDO
218C=======================================================================
219C 2 placing nodes in boxes
220C non -local candidates in SPMD
221C=======================================================================
222 DO j = 1, nsnr
223 iix(nsn+j)=int(nbx*(xrem(1,j)-xminb)/(xmaxb-xminb))
224 iiy(nsn+j)=int(nby*(xrem(2,j)-yminb)/(ymaxb-yminb))
225 iiz(nsn+j)=int(nbz*(xrem(3,j)-zminb)/(zmaxb-zminb))
226cc IF(IGAP/=0) BGAPSMX = MAX(BGAPSMX,XREM(12,J))
227 iix(nsn+j)=max(1,2+min(nbx,iix(nsn+j)))
228 iiy(nsn+j)=max(1,2+min(nby,iiy(nsn+j)))
229 iiz(nsn+j)=max(1,2+min(nbz,iiz(nsn+j)))
230
231 first = voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))
232 IF(first == 0)THEN
233c Empty Cell
234 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j)) = nsn+j ! first
235 next_nod(nsn+j) = 0 ! last one
236 last_nod(nsn+j) = 0 ! no last
237 ELSEIF(last_nod(first) == 0)THEN
238c cell containing one node
239c add as next node
240 next_nod(first) = nsn+j ! next
241 last_nod(first) = nsn+j ! last
242 next_nod(nsn+j) = 0 ! last one
243 ELSE
244c
245c jump to the last node of the cell
246 last = last_nod(first) ! last node in this voxel
247 next_nod(last) = nsn+j ! next
248 last_nod(first) = nsn+j ! last
249 next_nod(nsn+j) = 0 ! last one
250 ENDIF
251 ENDDO
252cc BGAPSMXG=BGAPSMX
253 END IF
254C Barrier to wait task0 treatment
255 CALL my_barrier
256cc BGAPSMX=BGAPSMXG
257C=======================================================================
258C 3 Searching boxes concerning each facet
259C and creation of candidates
260C=======================================================================
261 j_stok = 0
262
263 DO ne=1,nrtm
264C We do not retain the Destruit facets
265 IF(stf(ne) == zero)cycle
266
267 IF(igap == 0)THEN
268 aaa = tzinf+sqrt(three)*curv_max(ne)
269 ELSE
270 aaa = marge+sqrt(three)*(curv_max(ne)+
271 . min(gapmax,max(gapmin,bgapsmx+gap_m(ne))))
272 ENDIF
273
274c It is possible to improve the algo by cutting the facet
275c in 2 (4,3,6,9 ...) if the facet is large in front of AAA and inclinee
276
277 m1 = irect(1,ne)
278 m2 = irect(2,ne)
279 m3 = irect(3,ne)
280 m4 = irect(4,ne)
281
282 xx1=x(1,m1)
283 xx2=x(1,m2)
284 xx3=x(1,m3)
285 xx4=x(1,m4)
286 xmaxe=max(xx1,xx2,xx3,xx4)
287 xmine=min(xx1,xx2,xx3,xx4)
288
289 yy1=x(2,m1)
290 yy2=x(2,m2)
291 yy3=x(2,m3)
292 yy4=x(2,m4)
293 ymaxe=max(yy1,yy2,yy3,yy4)
294 ymine=min(yy1,yy2,yy3,yy4)
295
296 zz1=x(3,m1)
297 zz2=x(3,m2)
298 zz3=x(3,m3)
299 zz4=x(3,m4)
300 zmaxe=max(zz1,zz2,zz3,zz4)
301 zmine=min(zz1,zz2,zz3,zz4)
302
303
304c surface calculation (for future elimination of candidates)
305
306 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
307 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
308 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
309 s2 = sx*sx + sy*sy + sz*sz
310
311c index of voxels occupied by the facet
312
313 ix1=int(nbx*(xmine-aaa-xminb)/(xmaxb-xminb))
314 iy1=int(nby*(ymine-aaa-yminb)/(ymaxb-yminb))
315 iz1=int(nbz*(zmine-aaa-zminb)/(zmaxb-zminb))
316
317 ix1=max(1,2+min(nbx,ix1))
318 iy1=max(1,2+min(nby,iy1))
319 iz1=max(1,2+min(nbz,iz1))
320
321 ix2=int(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb))
322 iy2=int(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb))
323 iz2=int(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb))
324
325 ix2=max(1,2+min(nbx,ix2))
326 iy2=max(1,2+min(nby,iy2))
327 iz2=max(1,2+min(nbz,iz2))
328
329cc nbpelem = 0
330cc'nnpelem = 0
331cc'nnr0pelem = 0
332cc'nnrpelem = 0
333
334 DO iz = iz1,iz2
335 DO iy = iy1,iy2
336 DO ix = ix1,ix2
337
338cc nbpelem = nbpelem + 1
339
340 jj = voxel(ix,iy,iz)
341
342 DO WHILE(jj /= 0)
343
344cc nnpelem = nnpelem + 1
345
346 IF(jj<=nsn)THEN
347 nn=nsv(jj)
348 IF(nn == m1)GOTO 200
349 IF(nn == m2)GOTO 200
350 IF(nn == m3)GOTO 200
351 IF(nn == m4)GOTO 200
352 xs = x(1,nn)
353 ys = x(2,nn)
354 zs = x(3,nn)
355 IF(igap /= 0)THEN
356 aaa = marge+
357 . sqrt(three)*(curv_max(ne)+min(gapmax,max(gapmin,
358 . gap_s(jj)+gap_m(ne))))
359 ENDIF
360 ELSE
361 j=jj-nsn
362 xs = xrem(1,j)
363 ys = xrem(2,j)
364 zs = xrem(3,j)
365 IF(igap /= 0)THEN
366 aaa = marge+
367 . sqrt(three)*(curv_max(ne)+min(gapmax,max(gapmin,
368 . xrem(9,j)+gap_m(ne))))
369 ENDIF
370 ENDIF
371
372
373 IF(xs<=xmine-aaa)GOTO 200
374 IF(xs>=xmaxe+aaa)GOTO 200
375 IF(ys<=ymine-aaa)GOTO 200
376 IF(ys>=ymaxe+aaa)GOTO 200
377 IF(zs<=zmine-aaa)GOTO 200
378 IF(zs>=zmaxe+aaa)GOTO 200
379
380c underestimation of distance**2 for elimination of candidates
381
382cc'nnr0pelem = nnr0pelem + 1
383
384 d1x = xs - xx1
385 d1y = ys - yy1
386 d1z = zs - zz1
387 d2x = xs - xx2
388 d2y = ys - yy2
389 d2z = zs - zz2
390 dd1 = d1x*sx+d1y*sy+d1z*sz
391 dd2 = d2x*sx+d2y*sy+d2z*sz
392 IF(dd1*dd2 > zero)THEN
393 d2 = min(dd1*dd1,dd2*dd2)
394 a2 = aaa*aaa*s2
395 IF(d2 > a2)GOTO 200
396 ENDIF
397
398cc nnrpelem = nnrpelem + 1
399
400 j_stok = j_stok + 1
401 prov_n(j_stok) = jj
402 prov_e(j_stok) = ne
403 IF(j_stok == nvsiz)THEN
404c CALL I7STOVOX(
405 CALL i23sto(
406 1 nvsiz ,irect ,x ,nsv ,ii_stok,
407 2 cand_n ,cand_e ,mulnsn ,noint ,marge ,
408 3 i_mem ,prov_n ,prov_e ,eshift ,inacti ,
409 4 igap ,gap ,gap_s ,gap_m ,gapmin ,
410 5 gapmax ,curv_max ,msr ,nsn ,oldnum ,
411 6 nsnrold,cand_a ,ifpen ,cand_p )
412 IF(i_mem==2)GOTO 100
413 j_stok = 0
414 ENDIF
415
416 200 CONTINUE
417
418 jj = next_nod(jj)
419
420 ENDDO ! WHILE(JJ /= 0)
421
422 ENDDO
423 ENDDO
424 ENDDO
425cc nbpelg = nbpelg + nbpelem
426cc nnpelg = nnpelg + nnpelem
427cc nnrpelg = nnrpelg + nnrpelem
428cc nnr0pelg = nnr0pelg + nnr0pelem
429 ENDDO
430
431C-------------------------------------------------------------------------
432C end of sorting
433C-------------------------------------------------------------------------
434c IF(J_STOK/=0)CALL I7STOVOX(
435 IF(j_stok/=0)CALL i23sto(
436 1 j_stok ,irect ,x ,nsv ,ii_stok,
437 2 cand_n ,cand_e ,mulnsn ,noint ,marge ,
438 3 i_mem ,prov_n ,prov_e ,eshift ,inacti ,
439 4 igap ,gap ,gap_s ,gap_m ,gapmin ,
440 5 gapmax ,curv_max ,msr ,nsn ,oldnum ,
441 6 nsnrold,cand_a ,ifpen ,cand_p )
442
443C=======================================================================
444C 4 reset nodes in boxes to zero
445C=======================================================================
446 100 CONTINUE
447
448C Barrier to avoid reinitialization before end of sorting
449 CALL my_barrier
450 nsnf = 1 + itask*nsn / nthread
451 nsnl = (itask+1)*nsn / nthread
452
453 DO i=nsnf,nsnl
454 IF(iix(i)/=0)THEN
455 voxel(iix(i),iiy(i),iiz(i))=0
456 ENDIF
457 ENDDO
458C=======================================================================
459C 5 reset nodes in boxes to zero
460C non -local candidates in SPMD
461C=======================================================================
462 nsnf = 1 + itask*nsnr / nthread
463 nsnl = (itask+1)*nsnr / nthread
464 DO j = nsnf, nsnl
465 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))=0
466 ENDDO
467
468C
469 CALL my_barrier()
470 IF(itask == 0)THEN
471 DEALLOCATE(next_nod)
472 DEALLOCATE(iix)
473 DEALLOCATE(iiy)
474 DEALLOCATE(iiz)
475 ENDIF
476
477 RETURN
478 END
479
subroutine i23sto(j_stok, irect, x, nsv, ii_stok, cand_n, cand_e, mulnsn, noint, marge, i_mem, prov_n, prov_e, eshift, inacti, igap, gap, gap_s, gap_m, gapmin, gapmax, curv_max, msr, nsn, oldnum, nsnrold, cand_a, ifpen, cand_p)
Definition i23sto.F:38
subroutine i23trivox(nsn, renum, nsnr, isznsnr, i_mem, irect, x, stf, stfn, xyzm, nsv, ii_stok, cand_n, eshift, cand_e, mulnsn, noint, tzinf, msr, voxel, nbx, nby, nbz, inacti, cand_a, cand_p, ifpen, nrtm, nsnrold, igap, gap, gap_s, gap_m, gapmin, gapmax, marge, curv_max, nin, itask, bgapsmx, intheat, idt_therm, nodadt_therm)
Definition i23trivox.F:44
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, dimension(:), allocatable next_nod
Definition tri7box.F:48
subroutine spmd_oldnumcd(renum, oldnum, nsnr, nsnrold, intheat, idt_therm, nodadt_therm)
Definition spmd_i7tool.F:38
subroutine my_barrier
Definition machine.F:31