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 DE LA ROUTINE:
67C ===================
68C CLASSE LES NOEUDS DANS DES BOITES
69C RECHERCHE POUR CHAQUE FACETTE DES BOITES CONCERNES
70C RECHERCHE DES CANDIDATS
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,*) TABLEAU DES CONEC FACETTES E
77C X(3,*) COORDONNEES NODALES E
78C NSV NOS SYSTEMES DES NOEUDS E
79C XMAX plus grande abcisse existante E
80C XMAX plus grande ordonn. existante E
81C XMAX plus grande cote existante E
82C I_STOK niveau de stockage des couples
83C candidats impact E/S
84C CAND_N boites resultats noeuds
85C CAND_E adresses des boites resultat elements
86C MULNSN = MULTIMP*NSN TAILLE MAX ADMISE MAINTENANT POUR LES
87C COUPLES NOEUDS,ELT CANDIDATS
88C NOINT NUMERO USER DE L'INTERFACE
89C TZINF TAILLE ZONE INFLUENCE
90C
91C PROV_N CAND_N provisoire (variable static dans i7tri)
92C PROV_E CAND_E provisoire (variable static dans i7tri)
93
94C VOXEL(ix,iy,iz) contient le numero local du premier noeud de
95C la boite
96C NEXT_NOD(i) noeud suivant dans la meme boite (si /= 0)
97C LAST_NOD(i) dernier noeud dans la meme boite (si /= 0)
98C utilise uniquement pour aller directement du premier
99C noeud au dernier
100C-----------------------------------------------
101C D u m m y A r g u m e n t s
102C-----------------------------------------------
103 INTEGER I_MEM,ESHIFT,NSN,ISZNSNR,NRTM,NSNROLD,NIN,ITASK,
104 . MULNSN,NOINT,INACTI,NSNR,IGAP,NBX,NBY,NBZ,
105 . NSV(*),CAND_N(*),CAND_E(*),CAND_A(*),IFPEN(*),RENUM(*),
106 . IRECT(4,*), VOXEL(NBX+2,NBY+2,NBZ+2), MSR(*),II_STOK
107 INTEGER, INTENT(IN) :: INTHEAT
108 INTEGER, INTENT(IN) :: IDT_THERM
109 INTEGER, INTENT(IN) :: NODADT_THERM
110C REAL
111 my_real
112 . x(3,*),xyzm(6),cand_p(*),stf(*),stfn(*),gap_s(*),gap_m(*),
113 . tzinf,marge,gap,gapmin,gapmax,bgapsmx,
114 . curv_max(*)
115C-----------------------------------------------
116C L o c a l V a r i a b l e s
117C-----------------------------------------------
118 INTEGER NB_NCN,NB_NCN1,NB_ECN,I,J,DIR,NB_NC,NB_EC,
119 . N1,N2,N3,N4,NN,NE,K,L,NCAND_PROV,J_STOK,II,JJ,
120 . prov_n(mvsiz),prov_e(mvsiz),
121 . oldnum(isznsnr), nsnf, nsnl
122C REAL
123 my_real
124 . dx,dy,dz,xs,ys,zs,xx,sx,sy,sz,s2,
125 . xmin, xmax,ymin, ymax,zmin, zmax, tz, gapsmx, gapl,
126 . xx1,xx2,xx3,xx4,yy1,yy2,yy3,yy4,zz1,zz2,zz3,zz4,
127 . d1x,d1y,d1z,d2x,d2y,d2z,dd1,dd2,d2,a2,gs
128c provisoire
129 INTEGER LAST_NOD(NSN+NSNR)
130 INTEGER IX,IY,IZ,NEXT,M1,M2,M3,M4,
131 . IX1,IY1,IZ1,IX2,IY2,IZ2
132 INTEGER, DIMENSION(:),ALLOCATABLE :: IIX,IIY,IIZ
133 my_real
134 . XMINB,YMINB,ZMINB,XMAXB,YMAXB,ZMAXB,
135 . XMINE,YMINE,ZMINE,XMAXE,YMAXE,ZMAXE,AAA
136 INTEGER FIRST,NEW,LAST
137 SAVE IIX,IIY,IIZ
138C-----------------------------------------------
139 IF(itask == 0)THEN
140 ALLOCATE(next_nod(nsn+nsnr))
141 ALLOCATE(iix(nsn+nsnr))
142 ALLOCATE(iiy(nsn+nsnr))
143 ALLOCATE(iiz(nsn+nsnr))
144 END IF
145C Barrier to wait init voxel and allocation NEX_NOD
146 CALL my_barrier
147C Phase initiale de construction de BPE et BPN deplacee de I7BUCE => I7TRI
148C
149 xmin = xyzm(1)
150 ymin = xyzm(2)
151 zmin = xyzm(3)
152 xmax = xyzm(4)
153 ymax = xyzm(5)
154 zmax = xyzm(6)
155
156c dev future: xminb plus grand que xmin...
157 xminb = xmin
158 yminb = ymin
159 zminb = zmin
160 xmaxb = xmax
161 ymaxb = ymax
162 zmaxb = zmax
163
164
165C En SPMD, retrouve ancienne numerotation des candidats non locaux
166 IF(nspmd>1) THEN
167 CALL spmd_oldnumcd(renum,oldnum,isznsnr,nsnrold,intheat,idt_therm,nodadt_therm)
168 END IF
169
170C=======================================================================
171C 1 mise des noeuds dans les boites
172C=======================================================================
173 IF(itask == 0)THEN
174cc BGAPSMX = ZERO
175 DO i=1,nsn
176 iix(i)=0
177 iiy(i)=0
178 iiz(i)=0
179 IF(stfn(i) == zero)cycle
180 j=nsv(i)
181C Optimisation // recherche les noeuds compris dans xmin xmax des
182C elements du processeur
183 IF(x(1,j) < xmin) cycle
184 IF(x(1,j) > xmax) cycle
185 IF(x(2,j) < ymin) cycle
186 IF(x(2,j) > ymax) cycle
187 IF(x(3,j) < zmin) cycle
188 IF(x(3,j) > zmax) cycle
189
190 iix(i)=int(nbx*(x(1,j)-xminb)/(xmaxb-xminb))
191 iiy(i)=int(nby*(x(2,j)-yminb)/(ymaxb-yminb))
192 iiz(i)=int(nbz*(x(3,j)-zminb)/(zmaxb-zminb))
193
194 iix(i)=max(1,2+min(nbx,iix(i)))
195 iiy(i)=max(1,2+min(nby,iiy(i)))
196 iiz(i)=max(1,2+min(nbz,iiz(i)))
197
198 first = voxel(iix(i),iiy(i),iiz(i))
199 IF(first == 0)THEN
200c empty cell
201 voxel(iix(i),iiy(i),iiz(i)) = i ! first
202 next_nod(i) = 0 ! last one
203 last_nod(i) = 0 ! no last
204 ELSEIF(last_nod(first) == 0)THEN
205c cell containing one node
206c add as next node
207 next_nod(first) = i ! next
208 last_nod(first) = i ! last
209 next_nod(i) = 0 ! last one
210 ELSE
211c
212c jump to the last node of the cell
213 last = last_nod(first) ! last node in this voxel
214 next_nod(last) = i ! next
215 last_nod(first) = i ! last
216 next_nod(i) = 0 ! last one
217 ENDIF
218 ENDDO
219C=======================================================================
220C 2 mise des noeuds dans les boites
221C candidats non locaux en SPMD
222C=======================================================================
223 DO j = 1, nsnr
224 iix(nsn+j)=int(nbx*(xrem(1,j)-xminb)/(xmaxb-xminb))
225 iiy(nsn+j)=int(nby*(xrem(2,j)-yminb)/(ymaxb-yminb))
226 iiz(nsn+j)=int(nbz*(xrem(3,j)-zminb)/(zmaxb-zminb))
227cc IF(IGAP/=0) BGAPSMX = MAX(BGAPSMX,XREM(12,J))
228 iix(nsn+j)=max(1,2+min(nbx,iix(nsn+j)))
229 iiy(nsn+j)=max(1,2+min(nby,iiy(nsn+j)))
230 iiz(nsn+j)=max(1,2+min(nbz,iiz(nsn+j)))
231
232 first = voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))
233 IF(first == 0)THEN
234c empty cell
235 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j)) = nsn+j ! first
236 next_nod(nsn+j) = 0 ! last one
237 last_nod(nsn+j) = 0 ! no last
238 ELSEIF(last_nod(first) == 0)THEN
239c cell containing one node
240c add as next node
241 next_nod(first) = nsn+j ! next
242 last_nod(first) = nsn+j ! last
243 next_nod(nsn+j) = 0 ! last one
244 ELSE
245c
246c jump to the last node of the cell
247 last = last_nod(first) ! last node in this voxel
248 next_nod(last) = nsn+j ! next
249 last_nod(first) = nsn+j ! last
250 next_nod(nsn+j) = 0 ! last one
251 ENDIF
252 ENDDO
253cc BGAPSMXG=BGAPSMX
254 END IF
255C Barrier to wait task0 treatment
256 CALL my_barrier
257cc BGAPSMX=BGAPSMXG
258C=======================================================================
259C 3 recherche des boites concernant chaque facette
260C et creation des candidats
261C=======================================================================
262 j_stok = 0
263
264 DO ne=1,nrtm
265C on ne retient pas les facettes detruites
266 IF(stf(ne) == zero)cycle
267
268 IF(igap == 0)THEN
269 aaa = tzinf+sqrt(three)*curv_max(ne)
270 ELSE
271 aaa = marge+sqrt(three)*(curv_max(ne)+
272 . min(gapmax,max(gapmin,bgapsmx+gap_m(ne))))
273 ENDIF
274
275c il est possible d'ameliorer l'algo en decoupant la facette
276c en 2(4,3,6,9...) si la facette est grande devant AAA et inclinee
277
278 m1 = irect(1,ne)
279 m2 = irect(2,ne)
280 m3 = irect(3,ne)
281 m4 = irect(4,ne)
282
283 xx1=x(1,m1)
284 xx2=x(1,m2)
285 xx3=x(1,m3)
286 xx4=x(1,m4)
287 xmaxe=max(xx1,xx2,xx3,xx4)
288 xmine=min(xx1,xx2,xx3,xx4)
289
290 yy1=x(2,m1)
291 yy2=x(2,m2)
292 yy3=x(2,m3)
293 yy4=x(2,m4)
294 ymaxe=max(yy1,yy2,yy3,yy4)
295 ymine=min(yy1,yy2,yy3,yy4)
296
297 zz1=x(3,m1)
298 zz2=x(3,m2)
299 zz3=x(3,m3)
300 zz4=x(3,m4)
301 zmaxe=max(zz1,zz2,zz3,zz4)
302 zmine=min(zz1,zz2,zz3,zz4)
303
304
305c calcul de la surface (pour elimination future de candidats)
306
307 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
308 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
309 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
310 s2 = sx*sx + sy*sy + sz*sz
311
312c indice des voxels occupes par la facette
313
314 ix1=int(nbx*(xmine-aaa-xminb)/(xmaxb-xminb))
315 iy1=int(nby*(ymine-aaa-yminb)/(ymaxb-yminb))
316 iz1=int(nbz*(zmine-aaa-zminb)/(zmaxb-zminb))
317
318 ix1=max(1,2+min(nbx,ix1))
319 iy1=max(1,2+min(nby,iy1))
320 iz1=max(1,2+min(nbz,iz1))
321
322 ix2=int(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb))
323 iy2=int(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb))
324 iz2=int(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb))
325
326 ix2=max(1,2+min(nbx,ix2))
327 iy2=max(1,2+min(nby,iy2))
328 iz2=max(1,2+min(nbz,iz2))
329
330cc nbpelem = 0
331cc nnpelem = 0
332cc nnr0pelem = 0
333cc nnrpelem = 0
334
335 DO iz = iz1,iz2
336 DO iy = iy1,iy2
337 DO ix = ix1,ix2
338
339cc nbpelem = nbpelem + 1
340
341 jj = voxel(ix,iy,iz)
342
343 DO WHILE(jj /= 0)
344
345cc nnpelem = nnpelem + 1
346
347 IF(jj<=nsn)THEN
348 nn=nsv(jj)
349 IF(nn == m1)GOTO 200
350 IF(nn == m2)GOTO 200
351 IF(nn == m3)GOTO 200
352 IF(nn == m4)GOTO 200
353 xs = x(1,nn)
354 ys = x(2,nn)
355 zs = x(3,nn)
356 IF(igap /= 0)THEN
357 aaa = marge+
358 . sqrt(three)*(curv_max(ne)+min(gapmax,max(gapmin,
359 . gap_s(jj)+gap_m(ne))))
360 ENDIF
361 ELSE
362 j=jj-nsn
363 xs = xrem(1,j)
364 ys = xrem(2,j)
365 zs = xrem(3,j)
366 IF(igap /= 0)THEN
367 aaa = marge+
368 . sqrt(three)*(curv_max(ne)+min(gapmax,max(gapmin,
369 . xrem(9,j)+gap_m(ne))))
370 ENDIF
371 ENDIF
372
373
374 IF(xs<=xmine-aaa)GOTO 200
375 IF(xs>=xmaxe+aaa)GOTO 200
376 IF(ys<=ymine-aaa)GOTO 200
377 IF(ys>=ymaxe+aaa)GOTO 200
378 IF(zs<=zmine-aaa)GOTO 200
379 IF(zs>=zmaxe+aaa)GOTO 200
380
381c sousestimation de la distance**2 pour elimination de candidats
382
383cc nnr0pelem = nnr0pelem + 1
384
385 d1x = xs - xx1
386 d1y = ys - yy1
387 d1z = zs - zz1
388 d2x = xs - xx2
389 d2y = ys - yy2
390 d2z = zs - zz2
391 dd1 = d1x*sx+d1y*sy+d1z*sz
392 dd2 = d2x*sx+d2y*sy+d2z*sz
393 IF(dd1*dd2 > zero)THEN
394 d2 = min(dd1*dd1,dd2*dd2)
395 a2 = aaa*aaa*s2
396 IF(d2 > a2)GOTO 200
397 ENDIF
398
399cc nnrpelem = nnrpelem + 1
400
401 j_stok = j_stok + 1
402 prov_n(j_stok) = jj
403 prov_e(j_stok) = ne
404 IF(j_stok == nvsiz)THEN
405c CALL I7STOVOX(
406 CALL i23sto(
407 1 nvsiz ,irect ,x ,nsv ,ii_stok,
408 2 cand_n ,cand_e ,mulnsn ,noint ,marge ,
409 3 i_mem ,prov_n ,prov_e ,eshift ,inacti ,
410 4 igap ,gap ,gap_s ,gap_m ,gapmin ,
411 5 gapmax ,curv_max ,msr ,nsn ,oldnum ,
412 6 nsnrold,cand_a ,ifpen ,cand_p )
413 IF(i_mem==2)GOTO 100
414 j_stok = 0
415 ENDIF
416
417 200 CONTINUE
418
419 jj = next_nod(jj)
420
421 ENDDO ! WHILE(JJ /= 0)
422
423 ENDDO
424 ENDDO
425 ENDDO
426cc nbpelg = nbpelg + nbpelem
427cc nnpelg = nnpelg + nnpelem
428cc nnrpelg = nnrpelg + nnrpelem
429cc nnr0pelg = nnr0pelg + nnr0pelem
430 ENDDO
431
432C-------------------------------------------------------------------------
433C FIN DU TRI
434C-------------------------------------------------------------------------
435c IF(J_STOK/=0)CALL I7STOVOX(
436 IF(j_stok/=0)CALL i23sto(
437 1 j_stok ,irect ,x ,nsv ,ii_stok,
438 2 cand_n ,cand_e ,mulnsn ,noint ,marge ,
439 3 i_mem ,prov_n ,prov_e ,eshift ,inacti ,
440 4 igap ,gap ,gap_s ,gap_m ,gapmin ,
441 5 gapmax ,curv_max ,msr ,nsn ,oldnum ,
442 6 nsnrold,cand_a ,ifpen ,cand_p )
443
444C=======================================================================
445C 4 remise a zero des noeuds dans les boites
446C=======================================================================
447 100 CONTINUE
448
449C Barrier to avoid reinitialization before end of sorting
450 CALL my_barrier
451 nsnf = 1 + itask*nsn / nthread
452 nsnl = (itask+1)*nsn / nthread
453
454 DO i=nsnf,nsnl
455 IF(iix(i)/=0)THEN
456 voxel(iix(i),iiy(i),iiz(i))=0
457 ENDIF
458 ENDDO
459C=======================================================================
460C 5 remise a zero des noeuds dans les boites
461C candidats non locaux en SPMD
462C=======================================================================
463 nsnf = 1 + itask*nsnr / nthread
464 nsnl = (itask+1)*nsnr / nthread
465 DO j = nsnf, nsnl
466 voxel(iix(nsn+j),iiy(nsn+j),iiz(nsn+j))=0
467 ENDDO
468
469C
470 CALL my_barrier()
471 IF(itask == 0)THEN
472 DEALLOCATE(next_nod)
473 DEALLOCATE(iix)
474 DEALLOCATE(iiy)
475 DEALLOCATE(iiz)
476 ENDIF
477
478 RETURN
479 END
480
subroutine i23buce(x, irect, nsv, inacti, nrtm, nsn, cand_e, cand_n, gap, noint, ii_stok, tzinf, maxbox, minbox, ncontact, nb_n_b, eshift, cand_p, ncont, ild, weight, stfn, nin, stf, igap, gap_s, gapmin, gapmax, icurv, num_imp, itask, i_mem, msr, gap_m, nsnr, curv_max, renum, nsnrold, ifpen, mwag, bminma, nmn, irectg, bgapsmx, intheat, idt_therm, nodadt_therm)
Definition i23buce.F:46
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: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
subroutine spmd_oldnumcd(renum, oldnum, nsnr, nsnrold, intheat, idt_therm, nodadt_therm)
Definition spmd_i7tool.F:38
subroutine my_barrier
Definition machine.F:31