OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i7trivox1.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!|| i7trivox1 ../starter/source/interfaces/inter3d1/i7trivox1.F
25!||--- called by ------------------------------------------------------
26!|| i7buc_vox1 ../starter/source/interfaces/inter3d1/i7buc_vox1.F
27!||--- calls -----------------------------------------------------
28!|| i7cor3 ../starter/source/interfaces/inter3d1/i7cor3.F
29!|| i7dst3 ../starter/source/interfaces/inter3d1/i7dst3.F
30!|| i7pen3 ../starter/source/interfaces/inter3d1/i7pen3.F
31!|| inter_save_candidate ../starter/source/interfaces/inter3d1/inter_save_candidate.F90
32!||--- uses -----------------------------------------------------
33!|| inter_save_candidate_mod ../starter/source/interfaces/inter3d1/inter_save_candidate.F90
34!|| tri7box ../starter/share/modules1/tri7box.f
35!||====================================================================
36 SUBROUTINE i7trivox1(
37 1 NSN ,I_MEM ,IRECT ,X ,STF ,
38 2 STFN ,XYZM ,NSV ,
39 3 MULNSN ,NOINT ,TZINF ,GAP_S_L ,GAP_M_L ,
40 4 VOXEL ,NBX ,NBY ,NBZ ,NRTM_L ,
41 5 IGAP ,GAP ,GAP_S ,GAP_M ,GAPMIN ,
42 6 GAPMAX ,MARGE ,CURV_MAX ,BGAPSMX ,ISTF ,
43 7 I_STOK ,nin,
44 8 ID ,TITR ,DRAD ,INDEX ,
45 9 IREMNODE,FLAGREMNODE,KREMNODE,REMNODE,
46 1 DGAPLOAD,ipari,intbuf_tab,
47 2 iix,iiy,iiz,local_next_nod,nrtm,IS_USED_WITH_LAW151 )
48C============================================================================
49C M o d u l e s
50C-----------------------------------------------
51 USE my_alloc_mod
52 USE tri7box
54 use inter_save_candidate_mod , only : inter_save_candidate
55 use array_mod
56 use intbufdef_mod
57C-----------------------------------------------
58C I m p l i c i t T y p e s
59C-----------------------------------------------
60#include "implicit_f.inc"
61C-----------------------------------------------
62C G l o b a l P a r a m e t e r s
63C-----------------------------------------------
64#include "mvsiz_p.inc"
65c parameter setting the size for the vector (orig version is 128)
66 INTEGER NVECSZ
67 parameter(nvecsz = mvsiz)
68C-----------------------------------------------
69C C o m m o n B l o c k s
70C-----------------------------------------------
71#include "com04_c.inc"
72#include "param_c.inc"
73C-----------------------------------------------
74C ROLE DE LA ROUTINE:
75C ===================
76C CLASSE LES NOEUDS DANS DES BOITES
77C RECHERCHE POUR CHAQUE FACETTE DES BOITES CONCERNES
78C RECHERCHE DES CANDIDATS
79C-----------------------------------------------
80C D u m m y A r g u m e n t s
81C
82C NOM DESCRIPTION E/S
83C
84C IRECT(4,*) TABLEAU DES CONEC FACETTES E
85C X(3,*) COORDONNEES NODALES E
86C NSV NOS SYSTEMES DES NOEUDS E
87C XMAX plus grande abcisse existante E
88C XMAX plus grande ordonn. existante E
89C XMAX plus grande cote existante E
90C I_STOK niveau de stockage des couples
91C candidats impact E/S
92C CAND_N boites resultats noeuds
93C CAND_E adresses des boites resultat elements
94C MULNSN = MULTIMP*NSN TAILLE MAX ADMISE MAINTENANT POUR LES
95C COUPLES NOEUDS,ELT CANDIDATS
96C NOINT NUMERO USER DE L'INTERFACE
97C TZINF TAILLE ZONE INFLUENCE
98C
99C PROV_N CAND_N provisoire (variable static dans i7tri)
100C PROV_E CAND_E provisoire (variable static dans i7tri)
101
102C VOXEL(ix,iy,iz) contient le numero local du premier noeud de
103C la boite
104C LOCAL_NEXT_NOD(i) noeud suivant dans la meme boite (si /= 0)
105C LAST_NOD(i) dernier noeud dans la meme boite (si /= 0)
106C utilise uniquement pour aller directement du premier
107C noeud au dernier
108C INDEX index of IRECTM active for this domain
109C NRTM_L number of active segments for this domain
110C-----------------------------------------------
111C D u m m y A r g u m e n t s
112C-----------------------------------------------
113 integer, intent(in) :: nrtm !< number of segment
114 integer, intent(in) :: nrtm_l !< number of segment on the current proc
115 INTEGER I_MEM,NSN,
116 . MULNSN,NOINT,IGAP,NBX,NBY,NBZ,IREMNODE,FLAGREMNODE,
117 . NSV(NSN),
118 . IRECT(4,NRTM), VOXEL(NBX+2,NBY+2,NBZ+2),ISTF,
119 . I_STOK ,J_STOK,
120 . INDEX(*),KREMNODE(*),REMNODE(*)
121 my_real
122 . X(3,*),XYZM(6,2),STF(*),STFN(*),GAP_S(*),GAP_M(*),
123 . TZINF,MARGE,GAP,GAPMIN,GAPMAX,BGAPSMX,DRAD,
124 . curv_max(*),gap_s_l(*),gap_m_l(*)
125 INTEGER ID
126 CHARACTER(LEN=NCHARTITLE)::TITR
127 integer, intent(in) :: nin !< interface index
128 integer, dimension(npari), intent(inout) :: ipari !< interface data
129 type(intbuf_struct_), intent(inout) :: intbuf_tab !< interface data
130 INTEGER, dimension(nsn), intent(inout) :: iix,iiy,iiz,local_next_nod
131 LOGICAL,INTENT(IN) :: IS_USED_WITH_LAW151
132C-----------------------------------------------
133C L o c a l V a r i a b l e s
134C-----------------------------------------------
135 INTEGER NB_NCN,NB_NCN1,NB_ECN,I,J,DIR,NB_NC,NB_EC,
136 . n1,n2,n3,n4,nn,ne,k,l,ncand_prov,ii,jj,kk,
137 . nsnf, nsnl, i_bid,delnod
138 my_real
139 . dx,dy,dz,xs,ys,zs,xx,sx,sy,sz,s2,
140 . xmin, xmax,ymin, ymax,zmin, zmax, tz, gapsmx, gapl,
141 . xx1,xx2,xx3,xx4,yy1,yy2,yy3,yy4,zz1,zz2,zz3,zz4,
142 . d1x,d1y,d1z,d2x,d2y,d2z,dd1,dd2,d2,a2,gs,gapv(mvsiz)
143 INTEGER LAST_NOD(NSN)
144 INTEGER*8 meanjj8
145 INTEGER IX,IY,IZ,NEXT,M1,M2,M3,M4,
146 . ix1,iy1,iz1,ix2,iy2,iz2
147 integer, dimension(mvsiz) :: prov_n,prov_e
148 my_real
149 . xminb,yminb,zminb,xmaxb,ymaxb,zmaxb,
150 . xmine,ymine,zmine,xmaxe,ymaxe,zmaxe,aaa,tstart,tstop
151 my_real , INTENT(IN) :: dgapload
152 INTEGER FIRST,NEW,LAST,M
153 INTEGER, DIMENSION(MVSIZ) :: IX11,IX12,IX13,IX14,NSVG
154 my_real, DIMENSION(MVSIZ) :: x1,x2,x3,x4
155 my_real, DIMENSION(MVSIZ) :: y1,y2,y3,y4
156 my_real, DIMENSION(MVSIZ) :: z1,z2,z3,z4
157 my_real, DIMENSION(MVSIZ) :: xi,yi,zi
158 my_real, DIMENSION(MVSIZ) :: x0,y0,z0
159 my_real, DIMENSION(MVSIZ) :: nx1,ny1,nz1
160 my_real, DIMENSION(MVSIZ) :: nx2,ny2,nz2
161 my_real, DIMENSION(MVSIZ) :: nx3,ny3,nz3
162 my_real, DIMENSION(MVSIZ) :: nx4,ny4,nz4
163 my_real, DIMENSION(MVSIZ) :: p1,p2,p3,p4
164 my_real, DIMENSION(MVSIZ) :: lb1,lb2,lb3,lb4
165 my_real, DIMENSION(MVSIZ) :: lc1,lc2,lc3,lc4
166 my_real, DIMENSION(MVSIZ) :: n11,n21,n31
167 my_real, DIMENSION(MVSIZ) :: stif,pene
168
169 INTEGER, DIMENSION(:), ALLOCATABLE :: TAGNOD
170
171 integer , external :: omp_get_thread_num,omp_get_num_threads
172 integer :: itask,nthreads
173 integer :: my_old_size,my_address
174 integer :: local_i_stok,multimp
175 integer :: local_cand_n_size,local_cand_e_size
176
177 integer, save :: i_stok_old
178 integer, dimension(:), allocatable, save :: cand_n_size,cand_e_size
179 integer, dimension(:), allocatable, save :: address_cand_n,address_cand_e
180 type(array_type_int_1d) :: local_cand_n
181 type(array_type_int_1d) :: local_cand_e
182
183 integer :: my_size,mode
184 integer, dimension(:), allocatable :: my_index
185 integer, dimension(:,:), allocatable :: sort_array,save_array
186 integer, dimension(70000) :: work
187 integer stagnod
188C-----------------------------------------------
189 ! allocation of local arrays
190 itask = omp_get_thread_num()
191 nthreads = omp_get_num_threads()
192 local_cand_n_size = size(intbuf_tab%cand_n) / nthreads
193 local_cand_e_size = size(intbuf_tab%cand_e) / nthreads
194 local_i_stok = 0
195 local_cand_n%size_int_array_1d = local_cand_n_size
196 local_cand_e%size_int_array_1d = local_cand_e_size
197 call alloc_1d_array(local_cand_n)
198 call alloc_1d_array(local_cand_e)
199
200 xmin = xyzm(1,1)
201 ymin = xyzm(2,1)
202 zmin = xyzm(3,1)
203 xmax = xyzm(4,1)
204 ymax = xyzm(5,1)
205 zmax = xyzm(6,1)
206
207 xminb = xyzm(1,2)
208 yminb = xyzm(2,2)
209 zminb = xyzm(3,2)
210 xmaxb = xyzm(4,2)
211 ymaxb = xyzm(5,2)
212 zmaxb = xyzm(6,2)
213
214C=======================================================================
215C 1 mise des noeuds dans les boites
216C=======================================================================
217
218!$OMP MASTER
219
220 allocate( cand_n_size(nthreads+1),cand_e_size(nthreads+1) )
221 allocate( address_cand_n(nthreads+1),address_cand_e(nthreads+1) )
222 cand_n_size(1:nthreads+1) = 0
223 cand_e_size(1:nthreads+1) = 0
224 address_cand_n(1:nthreads+1) = 0
225 address_cand_e(1:nthreads+1) = 0
226
227 DO i=1,nsn
228 iix(i)=0
229 iiy(i)=0
230 iiz(i)=0
231 IF(stfn(i) == zero)cycle
232 j=nsv(i)
233
234C First test vs CRVOXEL
235 ix=int(lrvoxel*(x(1,j)-xmin)/(xmax-xmin))
236 IF(ix < 0 .OR. ix > lrvoxel) cycle
237 iy=int(lrvoxel*(x(2,j)-ymin)/(ymax-ymin))
238 IF(iy < 0 .OR. iy > lrvoxel) cycle
239 iz=int(lrvoxel*(x(3,j)-zmin)/(zmax-zmin))
240 IF(iz < 0 .OR. iz > lrvoxel) cycle
241 IF(.NOT.(btest(crvoxel(iy,iz),ix))) cycle
242
243C Second test vs reduced Box
244 IF( (x(1,j)-xminb)/(xmaxb-xminb) > one ) THEN
245 iix(i) = nbx
246 ELSE
247 iix(i)=int(max(nbx*(x(1,j)-xminb)/(xmaxb-xminb),-one))
248 ENDIF
249 IF( (x(2,j)-yminb)/(ymaxb-yminb) > one ) THEN
250 iiy(i) = nby
251 ELSE
252 iiy(i)=int(max(nby*(x(2,j)-yminb)/(ymaxb-yminb),-one))
253 ENDIF
254 IF( (x(3,j)-zminb)/(zmaxb-zminb) > one ) THEN
255 iiz(i) = nbz
256 ELSE
257 iiz(i)=int(max(nbz*(x(3,j)-zminb)/(zmaxb-zminb),-one))
258 ENDIF
259
260
261 iix(i)=max(1,2+iix(i))
262 iiy(i)=max(1,2+iiy(i))
263 iiz(i)=max(1,2+iiz(i))
264
265 first = voxel(iix(i),iiy(i),iiz(i))
266 IF(first == 0)THEN
267c empty cell
268 voxel(iix(i),iiy(i),iiz(i)) = i ! first
269 local_next_nod(i) = 0 ! last one
270 last_nod(i) = 0 ! no last
271 ELSEIF(last_nod(first) == 0)THEN
272c cell containing one node
273c add as next node
274 local_next_nod(first) = i ! next
275 last_nod(first) = i ! last
276 local_next_nod(i) = 0 ! last one
277 ELSE
278c
279c jump to the last node of the cell
280 last = last_nod(first) ! last node in this voxel
281 local_next_nod(last) = i ! next
282 last_nod(first) = i ! last
283 local_next_nod(i) = 0 ! last one
284 ENDIF
285 ENDDO
286!$OMP END MASTER
287
288!$OMP BARRIER
289C=======================================================================
290C 3 recherche des boites concernant chaque facette
291C et creation des candidats
292C=======================================================================
293
294 stagnod = numnod+numfakenodigeo
295 IF(is_used_with_law151) stagnod = stagnod + numels
296 ALLOCATE(tagnod(stagnod)) ; tagnod(1:stagnod) = 0
297
298 local_i_stok = 0
299 j_stok = 0
300!$OMP DO SCHEDULE(guided)
301 DO kk=1,nrtm_l
302 ne=index(kk)
303 IF(stf(ne) == zero)cycle
304
305 IF(flagremnode==2.AND.iremnode==2)THEN
306 k = kremnode(ne)+1
307 l = kremnode(ne+1)
308 DO m=k,l
309 tagnod(remnode(m)) = 1
310 ENDDO
311 ENDIF
312
313 IF(igap == 0)THEN
314 aaa = tzinf+curv_max(ne)
315 ELSE
316 aaa = marge+curv_max(ne)+max(min(gapmax,max(gapmin,bgapsmx+gap_m(ne)))+dgapload,drad)
317 ENDIF
318
319 m1 = irect(1,ne)
320 m2 = irect(2,ne)
321 m3 = irect(3,ne)
322 m4 = irect(4,ne)
323
324 xx1=x(1,m1)
325 xx2=x(1,m2)
326 xx3=x(1,m3)
327 xx4=x(1,m4)
328 xmaxe=max(xx1,xx2,xx3,xx4)
329 xmine=min(xx1,xx2,xx3,xx4)
330
331 yy1=x(2,m1)
332 yy2=x(2,m2)
333 yy3=x(2,m3)
334 yy4=x(2,m4)
335 ymaxe=max(yy1,yy2,yy3,yy4)
336 ymine=min(yy1,yy2,yy3,yy4)
337
338 zz1=x(3,m1)
339 zz2=x(3,m2)
340 zz3=x(3,m3)
341 zz4=x(3,m4)
342 zmaxe=max(zz1,zz2,zz3,zz4)
343 zmine=min(zz1,zz2,zz3,zz4)
344
345
346c calcul de la surface (pour elimination future de candidats)
347
348 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
349 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
350 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
351 s2 = sx*sx + sy*sy + sz*sz
352
353c indice des voxels occupes par la facette
354 IF( (xmine - aaa - xminb)/(xmaxb-xminb) > one ) THEN
355 ix1 = nbx
356 ELSE
357 ix1=int(max(nbx*(xmine-aaa-xminb)/(xmaxb-xminb),-one))
358 ENDIF
359 IF( (ymine - aaa - yminb)/(ymaxb-yminb) > one ) THEN
360 iy1 = nby
361 ELSE
362 iy1=int(max(nby*(ymine-aaa-yminb)/(ymaxb-yminb),-one))
363 ENDIF
364 IF( (zmine - aaa - zminb)/(zmaxb-zminb) > one ) THEN
365 iz1 = nbz
366 ELSE
367 iz1=int(max(nbz*(zmine-aaa-zminb)/(zmaxb-zminb),-one))
368 ENDIF
369
370 ix1=max(1,2+ix1)
371 iy1=max(1,2+iy1)
372 iz1=max(1,2+iz1)
373
374 IF( (xmaxe + aaa - xminb)/(xmaxb-xminb) > one ) THEN
375 ix2 = nbx
376 ELSE
377 ix2=int(max(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb),-one))
378 ENDIF
379 IF( (ymaxe + aaa - yminb)/(ymaxb-yminb) > one ) THEN
380 iy2 = nby
381 ELSE
382 iy2=int(max(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb),-one))
383 ENDIF
384 IF( (zmaxe + aaa - zminb)/(zmaxb-zminb) > one ) THEN
385 iz2 = nbz
386 ELSE
387 iz2=int(max(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb),-one))
388 ENDIF
389
390 ix2=max(1,2+ix2)
391 iy2=max(1,2+iy2)
392 iz2=max(1,2+iz2)
393
394 DO iz = iz1,iz2
395 DO iy = iy1,iy2
396 DO ix = ix1,ix2
397
398 jj = voxel(ix,iy,iz)
399
400 DO WHILE(jj /= 0)
401
402 nn=nsv(jj)
403 IF(tagnod(nn) == 1 ) GOTO 300
404 IF(nn == m1)GOTO 300
405 IF(nn == m2)GOTO 300
406 IF(nn == m3)GOTO 300
407 IF(nn == m4)GOTO 300
408
409 xs = x(1,nn)
410 ys = x(2,nn)
411 zs = x(3,nn)
412 IF(igap /= 0)THEN
413 aaa = marge+curv_max(ne)+max(min(gapmax,max(gapmin,gap_s(jj)+gap_m(ne)))+dgapload,drad)
414 ENDIF
415 IF(xs<=xmine-aaa)GOTO 300
416 IF(xs>=xmaxe+aaa)GOTO 300
417 IF(ys<=ymine-aaa)GOTO 300
418 IF(ys>=ymaxe+aaa)GOTO 300
419 IF(zs<=zmine-aaa)GOTO 300
420 IF(zs>=zmaxe+aaa)GOTO 300
421
422c sousestimation de la distance**2 pour elimination de candidats
423
424 d1x = xs - xx1
425 d1y = ys - yy1
426 d1z = zs - zz1
427 d2x = xs - xx2
428 d2y = ys - yy2
429 d2z = zs - zz2
430 dd1 = d1x*sx+d1y*sy+d1z*sz
431 dd2 = d2x*sx+d2y*sy+d2z*sz
432 IF(dd1*dd2 > zero)THEN
433 d2 = min(dd1*dd1,dd2*dd2)
434 a2 = aaa*aaa*s2
435 IF(d2 > a2)GOTO 300
436 ENDIF
437 j_stok = j_stok + 1
438 prov_n(j_stok) = jj
439 prov_e(j_stok) = ne
440 IF(j_stok == nvsiz) THEN
441 CALL i7cor3(x ,irect,nsv ,prov_e ,prov_n,
442 . stf ,stfn ,gapv ,igap ,gap ,
443 . gap_s,gap_m,istf ,gapmin ,gapmax,
444 . gap_s_l,gap_m_l ,drad ,ix11,ix12,
445 5 ix13 ,ix14 ,nsvg,x1 ,x2 ,
446 6 x3 ,x4 ,y1 ,y2 ,y3 ,
447 7 y4 ,z1 ,z2 ,z3 ,z4 ,
448 8 xi ,yi ,zi ,stif ,dgapload,
449 9 j_stok)
450 CALL i7dst3(ix13,ix14,x1 ,x2 ,x3 ,
451 1 x4 ,y1 ,y2 ,y3 ,y4 ,
452 2 z1 ,z2 ,z3 ,z4 ,xi ,
453 3 yi ,zi ,x0 ,y0 ,z0 ,
454 4 nx1,ny1,nz1,nx2,ny2,
455 5 nz2,nx3,ny3,nz3,nx4,
456 6 ny4,nz4,p1 ,p2 ,p3 ,
457 7 p4 ,lb1,lb2,lb3,lb4,
458 8 lc1,lc2,lc3,lc4,j_stok)
459 CALL i7pen3(marge,gapv,n11,n21,n31,
460 1 pene ,nx1 ,ny1,nz1,nx2,
461 2 ny2 ,nz2 ,nx3,ny3,nz3,
462 3 nx4 ,ny4 ,nz4,p1 ,p2 ,
463 4 p3 ,p4,j_stok)
464 call inter_save_candidate( local_i_stok,j_stok,prov_n,prov_e,pene,local_cand_n,local_cand_e )
465 j_stok = 0
466 ENDIF
467
468 300 CONTINUE
469 jj = local_next_nod(jj)
470 ENDDO ! WHILE(JJ /= 0)
471 ENDDO
472 ENDDO
473 ENDDO
474
475 IF(flagremnode==2.AND.iremnode==2)THEN
476 k = kremnode(ne)+1
477 l = kremnode(ne+1)
478 DO m=k,l
479 tagnod(remnode(m)) = 0
480 ENDDO
481 ENDIF
482
483 ENDDO
484!$OMP END DO
485 DEALLOCATE(tagnod)
486
487 IF(j_stok/=0) THEN
488 CALL i7cor3(x ,irect,nsv ,prov_e ,prov_n,
489 . stf ,stfn ,gapv ,igap ,gap ,
490 . gap_s,gap_m,istf ,gapmin ,gapmax,
491 . gap_s_l,gap_m_l ,drad ,ix11,ix12,
492 5 ix13 ,ix14 ,nsvg,x1 ,x2 ,
493 6 x3 ,x4 ,y1 ,y2 ,y3 ,
494 7 y4 ,z1 ,z2 ,z3 ,z4 ,
495 8 xi ,yi ,zi ,stif ,dgapload,
496 9 j_stok)
497 CALL i7dst3(ix13,ix14,x1 ,x2 ,x3 ,
498 1 x4 ,y1 ,y2 ,y3 ,y4 ,
499 2 z1 ,z2 ,z3 ,z4 ,xi ,
500 3 yi ,zi ,x0 ,y0 ,z0 ,
501 4 nx1,ny1,nz1,nx2,ny2,
502 5 nz2,nx3,ny3,nz3,nx4,
503 6 ny4,nz4,p1 ,p2 ,p3 ,
504 7 p4 ,lb1,lb2,lb3,lb4,
505 8 lc1,lc2,lc3,lc4,j_stok)
506 CALL i7pen3(marge,gapv,n11,n21,n31,
507 1 pene ,nx1 ,ny1,nz1,nx2,
508 2 ny2 ,nz2 ,nx3,ny3,nz3,
509 3 nx4 ,ny4 ,nz4,p1 ,p2 ,
510 4 p3 ,p4,j_stok)
511 call inter_save_candidate( local_i_stok,j_stok,prov_n,prov_e,pene,local_cand_n,local_cand_e )
512 j_stok = 0
513 ENDIF
514
515 ! save the local number of candidates
516 cand_n_size(itask+1) = local_i_stok
517 cand_e_size(itask+1) = local_i_stok
518!$OMP BARRIER
519
520
521!$OMP SINGLE
522 ! compute the address for each threads & the total number of candidates
523 address_cand_n(1) = 0
524 address_cand_e(1) = 0
525 ! ------------
526 do i=1,nthreads
527 address_cand_n(i+1) = cand_n_size(i) + address_cand_n(i) ! address for cand_n
528 address_cand_e(i+1) = cand_e_size(i) + address_cand_e(i) ! address for cand_e
529
530 cand_n_size(nthreads+1) = cand_n_size(nthreads+1) + cand_n_size(i) ! total number of cand_n for this proc
531 cand_e_size(nthreads+1) = cand_e_size(nthreads+1) + cand_e_size(i) ! total number of cand_e for this proc
532 enddo
533 ! ------------
534
535 ! ------------
536 ! check the size of global cand_n & cand_e
537 my_old_size = ipari(18)*ipari(23)
538 i_stok_old = i_stok
539 i_stok = i_stok + cand_n_size(nthreads+1) ! total number of cand_n/cand_e (= sum( cand_n/e per proc))
540 if(i_stok > my_old_size) then ! total number of candidates > size of cand_n/e --> need to extend the 2 arrays
541 multimp = i_stok/ipari(18) + 1
542 call upgrade_multimp(nin,multimp,intbuf_tab)
543 endif
544 ! ------------
545!$OMP END SINGLE
546
547!$OMP BARRIER
548
549 ! ------------
550 my_address = i_stok_old + address_cand_n(itask+1) ! get the address for each thread
551 intbuf_tab%cand_n(my_address+1:my_address+local_i_stok) = local_cand_n%int_array_1d(1:local_i_stok) ! save the cand_n into the global array
552 my_address = i_stok_old + address_cand_e(itask+1) ! get the address for each thread
553 intbuf_tab%cand_e(my_address+1:my_address+local_i_stok) = local_cand_e%int_array_1d(1:local_i_stok) ! save the cand_e into the global array
554 ! ------------
555
556
557 call dealloc_1d_array(local_cand_n)
558 call dealloc_1d_array(local_cand_e)
559
560!$OMP BARRIER
561
562 ! Sort the candidates to ensure the same domain decomposition
563!$OMP SINGLE
564 ! ------------
565 my_size = cand_n_size(nthreads+1)
566 allocate(my_index(2*my_size))
567 allocate(sort_array(2,my_size))
568 allocate(save_array(2,my_size))
569
570 my_address = i_stok_old + address_cand_n(1) ! get the address of the first pair of candidate
571 sort_array(1,1:my_size) = intbuf_tab%cand_n(my_address+1:my_address+my_size)
572 my_address = i_stok_old + address_cand_e(1) ! get the address of the first pair of candidate
573 sort_array(2,1:my_size) = intbuf_tab%cand_e(my_address+1:my_address+my_size)
574 do i=1,my_size
575 my_index(i) = i
576 enddo
577 save_array(1:2,1:my_size) = sort_array(1:2,1:my_size)
578 mode = 0
579
580 call my_orders( mode,work,sort_array,my_index,my_size,2)
581 my_address = i_stok_old + address_cand_n(1) ! get the address of the first pair of candidate
582 do i=1,my_size
583 intbuf_tab%cand_n(my_address+i) = save_array(1,my_index(i))
584 enddo
585 my_address = i_stok_old + address_cand_e(1) ! get the address of the first pair of candidate
586 do i=1,my_size
587 intbuf_tab%cand_e(my_address+i) = save_array(2,my_index(i))
588 enddo
589 deallocate(my_index)
590 deallocate(sort_array)
591 deallocate(save_array)
592 ! ------------
593!$OMP END SINGLE
594
595!$OMP DO SCHEDULE(guided)
596 DO i=1,nsn
597 IF(iix(i)/=0)THEN
598 voxel(iix(i),iiy(i),iiz(i))=0
599 ENDIF
600 ENDDO
601!$OMP END DO
602
603!$OMP SINGLE
604 deallocate( cand_n_size,cand_e_size )
605 deallocate( address_cand_n,address_cand_e )
606!$OMP END SINGLE
607
608 RETURN
609 END
subroutine i7trivox1(nsn, i_mem, irect, x, stf, stfn, xyzm, nsv, mulnsn, noint, tzinf, gap_s_l, gap_m_l, voxel, nbx, nby, nbz, nrtm_l, igap, gap, gap_s, gap_m, gapmin, gapmax, marge, curv_max, bgapsmx, istf, i_stok, nin, id, titr, drad, index, iremnode, flagremnode, kremnode, remnode, dgapload, ipari, intbuf_tab, iix, iiy, iiz, local_next_nod, nrtm, is_used_with_law151)
Definition i7trivox1.F:48
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
void my_orders(int *mode, int *iwork, int *data, int *index, int *n, int *irecl)
Definition my_orders.c:82
integer, parameter nchartitle
integer, dimension(0:lrvoxel, 0:lrvoxel) crvoxel
Definition tri7box.F:56
integer lrvoxel
Definition tri7box.F:54
subroutine i7cor3(x, irect, nsv, cand_e, cand_n, stf, stfn, gapv, igap, gap, gap_s, gap_m, istf, gapmin, gapmax, gap_s_l, gap_m_l, drad, ix1, ix2, ix3, ix4, nsvg, x1, x2, x3, x4, y1, y2, y3, y4, z1, z2, z3, z4, xi, yi, zi, stif, dgapload, last)
Definition i7cor3.F:43
subroutine i7dst3(ix3, ix4, x1, x2, x3, x4, y1, y2, y3, y4, z1, z2, z3, z4, xi, yi, zi, x0, y0, z0, nx1, ny1, nz1, nx2, ny2, nz2, nx3, ny3, nz3, nx4, ny4, nz4, p1, p2, p3, p4, lb1, lb2, lb3, lb4, lc1, lc2, lc3, lc4, last)
Definition i7dst3.F:46
subroutine i7pen3(marge, gapv, n1, n2, n3, pene, nx1, ny1, nz1, nx2, ny2, nz2, nx3, ny3, nz3, nx4, ny4, nz4, p1, p2, p3, p4, last)
Definition i7pen3.F:43
program starter
Definition starter.F:39
subroutine upgrade_multimp(ni, multimp_parameter, intbuf_tab)