OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i7trivox1.F File Reference
#include "implicit_f.inc"
#include "mvsiz_p.inc"
#include "com04_c.inc"
#include "param_c.inc"

Go to the source code of this file.

Functions/Subroutines

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)

Function/Subroutine Documentation

◆ i7trivox1()

subroutine i7trivox1 ( integer nsn,
integer i_mem,
integer, dimension(4,nrtm) irect,
x,
stf,
stfn,
xyzm,
integer, dimension(nsn) nsv,
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, intent(in) nrtm_l,
integer igap,
gap,
gap_s,
gap_m,
gapmin,
gapmax,
marge,
curv_max,
bgapsmx,
integer istf,
integer i_stok,
integer, intent(in) nin,
integer id,
character(len=nchartitle) titr,
drad,
integer, dimension(*) index,
integer iremnode,
integer flagremnode,
integer, dimension(*) kremnode,
integer, dimension(*) remnode,
intent(in) dgapload,
integer, dimension(npari), intent(inout) ipari,
type(intbuf_struct_), intent(inout) intbuf_tab,
integer, dimension(nsn), intent(inout) iix,
integer, dimension(nsn), intent(inout) iiy,
integer, dimension(nsn), intent(inout) iiz,
integer, dimension(nsn), intent(inout) local_next_nod,
integer, intent(in) nrtm,
logical, intent(in) is_used_with_law151 )
Parameters
[in]nrtmnumber of segment
[in]nrtm_lnumber of segment on the current proc
[in]nininterface index
[in,out]ipariinterface data
[in,out]intbuf_tabinterface data

Definition at line 36 of file i7trivox1.F.

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 OF THE ROUTINE:
75C ===================
76C CLASSIFIES NODES INTO BOXES
77C SEARCHES FOR EACH FACETTE OF CONCERNED BOXES
78C SEARCH FOR CANDIDATES
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,*) ARRAY OF CONEC FACETTES E
85C X(3,*) COORDONNEES NODALES E
86C NSV NOS SYSTEMES DES NODES E
87C Xmax larger abcisse existing e
88C XMAX largest order.existing E
89C Xmax larger existing side E
90C I_STOK storage level of couples
91C CANDIDATES impact E/S
92C CAND_N boites resultats nodes C CAND_E adresses des boites resultat elements
93C MULNSN = MULTIMP*NSN MAX SIZE NOW ALLOWED FOR
94C COUPLES NODES,ELT CANDIDATES
95C NOINT INTERFACE USER NUMBER
96C TZINF TAILLE ZONE INFLUENCE
97C
98C Prov_n Provisional Cand_n (static variable in i7tri)
99C PROV_E CAND_E provisoire (variable static in i7tri)
100
101C VOXEL(ix,iy,iz) contains the local number of the first node of
102C the box
103C Local_next_nod (I) Next node in the same box (if /= 0)
104C Last_nod (i) Last node in the same box (if /= 0)
105C used only to go directly from the first
106C node at the last
107C INDEX index of IRECTM active for this domain
108C NRTM_L number of active segments for this domain
109C-----------------------------------------------
110C D u m m y A r g u m e n t s
111C-----------------------------------------------
112 integer, intent(in) :: nrtm !< number of segment
113 integer, intent(in) :: nrtm_l !< number of segment on the current proc
114 INTEGER I_MEM,NSN,
115 . MULNSN,NOINT,IGAP,NBX,NBY,NBZ,IREMNODE,FLAGREMNODE,
116 . NSV(NSN),
117 . IRECT(4,NRTM), VOXEL(NBX+2,NBY+2,NBZ+2),ISTF,
118 . I_STOK ,J_STOK,
119 . INDEX(*),KREMNODE(*),REMNODE(*)
120 my_real
121 . x(3,*),xyzm(6,2),stf(*),stfn(*),gap_s(*),gap_m(*),
122 . tzinf,marge,gap,gapmin,gapmax,bgapsmx,drad,
123 . curv_max(*),gap_s_l(*),gap_m_l(*)
124 INTEGER ID
125 CHARACTER(LEN=NCHARTITLE)::TITR
126 integer, intent(in) :: nin !< interface index
127 integer, dimension(npari), intent(inout) :: ipari !< interface data
128 type(intbuf_struct_), intent(inout) :: intbuf_tab !< interface data
129 INTEGER, dimension(nsn), intent(inout) :: iix,iiy,iiz,local_next_nod
130 LOGICAL,INTENT(IN) :: IS_USED_WITH_LAW151
131C-----------------------------------------------
132C L o c a l V a r i a b l e s
133C-----------------------------------------------
134 INTEGER NB_NCN,NB_NCN1,NB_ECN,I,J,DIR,NB_NC,NB_EC,
135 . N1,N2,N3,N4,NN,NE,K,L,NCAND_PROV,II,JJ,KK,
136 . NSNF, NSNL, I_BID,DELNOD
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,gapv(mvsiz)
142 INTEGER LAST_NOD(NSN)
143 INTEGER*8 meanjj8
144 INTEGER IX,IY,IZ,NEXT,M1,M2,M3,M4,
145 . IX1,IY1,IZ1,IX2,IY2,IZ2
146 integer, dimension(mvsiz) :: prov_n,prov_e
147 my_real
148 . xminb,yminb,zminb,xmaxb,ymaxb,zmaxb,
149 . xmine,ymine,zmine,xmaxe,ymaxe,zmaxe,aaa,tstart,tstop
150 my_real , INTENT(IN) :: dgapload
151 INTEGER FIRST,NEW,LAST,M
152 INTEGER, DIMENSION(MVSIZ) :: IX11,IX12,IX13,IX14,NSVG
153 my_real, DIMENSION(MVSIZ) :: x1,x2,x3,x4
154 my_real, DIMENSION(MVSIZ) :: y1,y2,y3,y4
155 my_real, DIMENSION(MVSIZ) :: z1,z2,z3,z4
156 my_real, DIMENSION(MVSIZ) :: xi,yi,zi
157 my_real, DIMENSION(MVSIZ) :: x0,y0,z0
158 my_real, DIMENSION(MVSIZ) :: nx1,ny1,nz1
159 my_real, DIMENSION(MVSIZ) :: nx2,ny2,nz2
160 my_real, DIMENSION(MVSIZ) :: nx3,ny3,nz3
161 my_real, DIMENSION(MVSIZ) :: nx4,ny4,nz4
162 my_real, DIMENSION(MVSIZ) :: p1,p2,p3,p4
163 my_real, DIMENSION(MVSIZ) :: lb1,lb2,lb3,lb4
164 my_real, DIMENSION(MVSIZ) :: lc1,lc2,lc3,lc4
165 my_real, DIMENSION(MVSIZ) :: n11,n21,n31
166 my_real, DIMENSION(MVSIZ) :: stif,pene
167
168 INTEGER, DIMENSION(:), ALLOCATABLE :: TAGNOD
169
170 integer , external :: omp_get_thread_num,omp_get_num_threads
171 integer :: itask,nthreads
172 integer :: my_old_size,my_address
173 integer :: local_i_stok,multimp
174 integer :: local_cand_n_size,local_cand_e_size
175
176 integer, save :: i_stok_old
177 integer, dimension(:), allocatable, save :: cand_n_size,cand_e_size
178 integer, dimension(:), allocatable, save :: address_cand_n,address_cand_e
179 type(array_type_int_1d) :: local_cand_n
180 type(array_type_int_1d) :: local_cand_e
181
182 integer :: my_size,mode
183 integer, dimension(:), allocatable :: my_index
184 integer, dimension(:,:), allocatable :: sort_array,save_array
185 integer, dimension(70000) :: work
186 integer stagnod
187C-----------------------------------------------
188 ! allocation of local arrays
189 itask = omp_get_thread_num()
190 nthreads = omp_get_num_threads()
191 local_cand_n_size = size(intbuf_tab%cand_n) / nthreads
192 local_cand_e_size = size(intbuf_tab%cand_e) / nthreads
193 local_i_stok = 0
194 local_cand_n%size_int_array_1d = local_cand_n_size
195 local_cand_e%size_int_array_1d = local_cand_e_size
196 call alloc_1d_array(local_cand_n)
197 call alloc_1d_array(local_cand_e)
198
199 xmin = xyzm(1,1)
200 ymin = xyzm(2,1)
201 zmin = xyzm(3,1)
202 xmax = xyzm(4,1)
203 ymax = xyzm(5,1)
204 zmax = xyzm(6,1)
205
206 xminb = xyzm(1,2)
207 yminb = xyzm(2,2)
208 zminb = xyzm(3,2)
209 xmaxb = xyzm(4,2)
210 ymaxb = xyzm(5,2)
211 zmaxb = xyzm(6,2)
212
213C=======================================================================
214C 1 placing nodes in boxes
215C=======================================================================
216
217!$OMP MASTER
218
219 allocate( cand_n_size(nthreads+1),cand_e_size(nthreads+1) )
220 allocate( address_cand_n(nthreads+1),address_cand_e(nthreads+1) )
221 cand_n_size(1:nthreads+1) = 0
222 cand_e_size(1:nthreads+1) = 0
223 address_cand_n(1:nthreads+1) = 0
224 address_cand_e(1:nthreads+1) = 0
225
226 DO i=1,nsn
227 iix(i)=0
228 iiy(i)=0
229 iiz(i)=0
230 IF(stfn(i) == zero)cycle
231 j=nsv(i)
232
233C First test vs CRVOXEL
234 ix=int(lrvoxel*(x(1,j)-xmin)/(xmax-xmin))
235 IF(ix < 0 .OR. ix > lrvoxel) cycle
236 iy=int(lrvoxel*(x(2,j)-ymin)/(ymax-ymin))
237 IF(iy < 0 .OR. iy > lrvoxel) cycle
238 iz=int(lrvoxel*(x(3,j)-zmin)/(zmax-zmin))
239 IF(iz < 0 .OR. iz > lrvoxel) cycle
240 IF(.NOT.(btest(crvoxel(iy,iz),ix))) cycle
241
242C Second test vs reduced Box
243 IF( (x(1,j)-xminb)/(xmaxb-xminb) > one ) THEN
244 iix(i) = nbx
245 ELSE
246 iix(i)=int(max(nbx*(x(1,j)-xminb)/(xmaxb-xminb),-one))
247 ENDIF
248 IF( (x(2,j)-yminb)/(ymaxb-yminb) > one ) THEN
249 iiy(i) = nby
250 ELSE
251 iiy(i)=int(max(nby*(x(2,j)-yminb)/(ymaxb-yminb),-one))
252 ENDIF
253 IF( (x(3,j)-zminb)/(zmaxb-zminb) > one ) THEN
254 iiz(i) = nbz
255 ELSE
256 iiz(i)=int(max(nbz*(x(3,j)-zminb)/(zmaxb-zminb),-one))
257 ENDIF
258
259
260 iix(i)=max(1,2+iix(i))
261 iiy(i)=max(1,2+iiy(i))
262 iiz(i)=max(1,2+iiz(i))
263
264 first = voxel(iix(i),iiy(i),iiz(i))
265 IF(first == 0)THEN
266c Empty Cell
267 voxel(iix(i),iiy(i),iiz(i)) = i ! first
268 local_next_nod(i) = 0 ! last one
269 last_nod(i) = 0 ! no last
270 ELSEIF(last_nod(first) == 0)THEN
271c cell containing one node
272c add as next node
273 local_next_nod(first) = i ! next
274 last_nod(first) = i ! last
275 local_next_nod(i) = 0 ! last one
276 ELSE
277c
278c jump to the last node of the cell
279 last = last_nod(first) ! last node in this voxel
280 local_next_nod(last) = i ! next
281 last_nod(first) = i ! last
282 local_next_nod(i) = 0 ! last one
283 ENDIF
284 ENDDO
285!$OMP END MASTER
286
287!$OMP BARRIER
288C=======================================================================
289C 3 Searching boxes concerning each facet
290C and creation of candidates
291C=======================================================================
292
293 stagnod = numnod+numfakenodigeo
294 IF(is_used_with_law151) stagnod = stagnod + numels
295 ALLOCATE(tagnod(stagnod)) ; tagnod(1:stagnod) = 0
296
297 local_i_stok = 0
298 j_stok = 0
299!$omp DO schedule(guided)
300 DO kk=1,nrtm_l
301 ne=index(kk)
302 IF(stf(ne) == zero)cycle
303
304 IF(flagremnode==2.AND.iremnode==2)THEN
305 k = kremnode(ne)+1
306 l = kremnode(ne+1)
307 DO m=k,l
308 tagnod(remnode(m)) = 1
309 ENDDO
310 ENDIF
311
312 IF(igap == 0)THEN
313 aaa = tzinf+curv_max(ne)
314 ELSE
315 aaa = marge+curv_max(ne)+max(min(gapmax,max(gapmin,bgapsmx+gap_m(ne)))+dgapload,drad)
316 ENDIF
317
318 m1 = irect(1,ne)
319 m2 = irect(2,ne)
320 m3 = irect(3,ne)
321 m4 = irect(4,ne)
322
323 xx1=x(1,m1)
324 xx2=x(1,m2)
325 xx3=x(1,m3)
326 xx4=x(1,m4)
327 xmaxe=max(xx1,xx2,xx3,xx4)
328 xmine=min(xx1,xx2,xx3,xx4)
329
330 yy1=x(2,m1)
331 yy2=x(2,m2)
332 yy3=x(2,m3)
333 yy4=x(2,m4)
334 ymaxe=max(yy1,yy2,yy3,yy4)
335 ymine=min(yy1,yy2,yy3,yy4)
336
337 zz1=x(3,m1)
338 zz2=x(3,m2)
339 zz3=x(3,m3)
340 zz4=x(3,m4)
341 zmaxe=max(zz1,zz2,zz3,zz4)
342 zmine=min(zz1,zz2,zz3,zz4)
343
344
345c calculation of the surface (for future elimination of candidates)
346
347 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
348 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
349 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
350 s2 = sx*sx + sy*sy + sz*sz
351
352c index of voxels occupied by the facet
353 IF( (xmine - aaa - xminb)/(xmaxb-xminb) > one ) THEN
354 ix1 = nbx
355 ELSE
356 ix1=int(max(nbx*(xmine-aaa-xminb)/(xmaxb-xminb),-one))
357 ENDIF
358 IF( (ymine - aaa - yminb)/(ymaxb-yminb) > one ) THEN
359 iy1 = nby
360 ELSE
361 iy1=int(max(nby*(ymine-aaa-yminb)/(ymaxb-yminb),-one))
362 ENDIF
363 IF( (zmine - aaa - zminb)/(zmaxb-zminb) > one ) THEN
364 iz1 = nbz
365 ELSE
366 iz1=int(max(nbz*(zmine-aaa-zminb)/(zmaxb-zminb),-one))
367 ENDIF
368
369 ix1=max(1,2+ix1)
370 iy1=max(1,2+iy1)
371 iz1=max(1,2+iz1)
372
373 IF( (xmaxe + aaa - xminb)/(xmaxb-xminb) > one ) THEN
374 ix2 = nbx
375 ELSE
376 ix2=int(max(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb),-one))
377 ENDIF
378 IF( (ymaxe + aaa - yminb)/(ymaxb-yminb) > one ) THEN
379 iy2 = nby
380 ELSE
381 iy2=int(max(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb),-one))
382 ENDIF
383 IF( (zmaxe + aaa - zminb)/(zmaxb-zminb) > one ) THEN
384 iz2 = nbz
385 ELSE
386 iz2=int(max(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb),-one))
387 ENDIF
388
389 ix2=max(1,2+ix2)
390 iy2=max(1,2+iy2)
391 iz2=max(1,2+iz2)
392
393 DO iz = iz1,iz2
394 DO iy = iy1,iy2
395 DO ix = ix1,ix2
396
397 jj = voxel(ix,iy,iz)
398
399 DO WHILE(jj /= 0)
400
401 nn=nsv(jj)
402 IF(tagnod(nn) == 1 ) GOTO 300
403 IF(nn == m1)GOTO 300
404 IF(nn == m2)GOTO 300
405 IF(nn == m3)GOTO 300
406 IF(nn == m4)GOTO 300
407
408 xs = x(1,nn)
409 ys = x(2,nn)
410 zs = x(3,nn)
411 IF(igap /= 0)THEN
412 aaa = marge+curv_max(ne)+max(min(gapmax,max(gapmin,gap_s(jj)+gap_m(ne)))+dgapload,drad)
413 ENDIF
414 IF(xs<=xmine-aaa)GOTO 300
415 IF(xs>=xmaxe+aaa)GOTO 300
416 IF(ys<=ymine-aaa)GOTO 300
417 IF(ys>=ymaxe+aaa)GOTO 300
418 IF(zs<=zmine-aaa)GOTO 300
419 IF(zs>=zmaxe+aaa)GOTO 300
420
421c underestimation of the distance**2 for candidate elimination
422
423 d1x = xs - xx1
424 d1y = ys - yy1
425 d1z = zs - zz1
426 d2x = xs - xx2
427 d2y = ys - yy2
428 d2z = zs - zz2
429 dd1 = d1x*sx+d1y*sy+d1z*sz
430 dd2 = d2x*sx+d2y*sy+d2z*sz
431 IF(dd1*dd2 > zero)THEN
432 d2 = min(dd1*dd1,dd2*dd2)
433 a2 = aaa*aaa*s2
434 IF(d2 > a2)GOTO 300
435 ENDIF
436 j_stok = j_stok + 1
437 prov_n(j_stok) = jj
438 prov_e(j_stok) = ne
439 IF(j_stok == nvsiz) THEN
440 CALL i7cor3(x ,irect,nsv ,prov_e ,prov_n,
441 . stf ,stfn ,gapv ,igap ,gap ,
442 . gap_s,gap_m,istf ,gapmin ,gapmax,
443 . gap_s_l,gap_m_l ,drad ,ix11,ix12,
444 5 ix13 ,ix14 ,nsvg,x1 ,x2 ,
445 6 x3 ,x4 ,y1 ,y2 ,y3 ,
446 7 y4 ,z1 ,z2 ,z3 ,z4 ,
447 8 xi ,yi ,zi ,stif ,dgapload,
448 9 j_stok)
449 CALL i7dst3(ix13,ix14,x1 ,x2 ,x3 ,
450 1 x4 ,y1 ,y2 ,y3 ,y4 ,
451 2 z1 ,z2 ,z3 ,z4 ,xi ,
452 3 yi ,zi ,x0 ,y0 ,z0 ,
453 4 nx1,ny1,nz1,nx2,ny2,
454 5 nz2,nx3,ny3,nz3,nx4,
455 6 ny4,nz4,p1 ,p2 ,p3 ,
456 7 p4 ,lb1,lb2,lb3,lb4,
457 8 lc1,lc2,lc3,lc4,j_stok)
458 CALL i7pen3(marge,gapv,n11,n21,n31,
459 1 pene ,nx1 ,ny1,nz1,nx2,
460 2 ny2 ,nz2 ,nx3,ny3,nz3,
461 3 nx4 ,ny4 ,nz4,p1 ,p2 ,
462 4 p3 ,p4,j_stok)
463 call inter_save_candidate( local_i_stok,j_stok,prov_n,prov_e,pene,local_cand_n,local_cand_e )
464 j_stok = 0
465 ENDIF
466
467 300 CONTINUE
468 jj = local_next_nod(jj)
469 ENDDO ! WHILE(JJ /= 0)
470 ENDDO
471 ENDDO
472 ENDDO
473
474 IF(flagremnode==2.AND.iremnode==2)THEN
475 k = kremnode(ne)+1
476 l = kremnode(ne+1)
477 DO m=k,l
478 tagnod(remnode(m)) = 0
479 ENDDO
480 ENDIF
481
482 ENDDO
483!$OMP END DO
484 DEALLOCATE(tagnod)
485
486 IF(j_stok/=0) THEN
487 CALL i7cor3(x ,irect,nsv ,prov_e ,prov_n,
488 . stf ,stfn ,gapv ,igap ,gap ,
489 . gap_s,gap_m,istf ,gapmin ,gapmax,
490 . gap_s_l,gap_m_l ,drad ,ix11,ix12,
491 5 ix13 ,ix14 ,nsvg,x1 ,x2 ,
492 6 x3 ,x4 ,y1 ,y2 ,y3 ,
493 7 y4 ,z1 ,z2 ,z3 ,z4 ,
494 8 xi ,yi ,zi ,stif ,dgapload,
495 9 j_stok)
496 CALL i7dst3(ix13,ix14,x1 ,x2 ,x3 ,
497 1 x4 ,y1 ,y2 ,y3 ,y4 ,
498 2 z1 ,z2 ,z3 ,z4 ,xi ,
499 3 yi ,zi ,x0 ,y0 ,z0 ,
500 4 nx1,ny1,nz1,nx2,ny2,
501 5 nz2,nx3,ny3,nz3,nx4,
502 6 ny4,nz4,p1 ,p2 ,p3 ,
503 7 p4 ,lb1,lb2,lb3,lb4,
504 8 lc1,lc2,lc3,lc4,j_stok)
505 CALL i7pen3(marge,gapv,n11,n21,n31,
506 1 pene ,nx1 ,ny1,nz1,nx2,
507 2 ny2 ,nz2 ,nx3,ny3,nz3,
508 3 nx4 ,ny4 ,nz4,p1 ,p2 ,
509 4 p3 ,p4,j_stok)
510 call inter_save_candidate( local_i_stok,j_stok,prov_n,prov_e,pene,local_cand_n,local_cand_e )
511 j_stok = 0
512 ENDIF
513
514 ! save the local number of candidates
515 cand_n_size(itask+1) = local_i_stok
516 cand_e_size(itask+1) = local_i_stok
517!$OMP BARRIER
518
519
520!$OMP SINGLE
521 ! compute the address for each threads & the total number of candidates
522 address_cand_n(1) = 0
523 address_cand_e(1) = 0
524 ! ------------
525 do i=1,nthreads
526 address_cand_n(i+1) = cand_n_size(i) + address_cand_n(i) ! address for cand_n
527 address_cand_e(i+1) = cand_e_size(i) + address_cand_e(i) ! address for cand_e
528
529 cand_n_size(nthreads+1) = cand_n_size(nthreads+1) + cand_n_size(i) ! total number of cand_n for this proc
530 cand_e_size(nthreads+1) = cand_e_size(nthreads+1) + cand_e_size(i) ! total number of cand_e for this proc
531 enddo
532 ! ------------
533
534 ! ------------
535 ! check the size of global cand_n & cand_e
536 my_old_size = ipari(18)*ipari(23)
537 i_stok_old = i_stok
538 i_stok = i_stok + cand_n_size(nthreads+1) ! total number of cand_n/cand_e (= sum( cand_n/e per proc))
539 if(i_stok > my_old_size) then ! total number of candidates > size of cand_n/e --> need to extend the 2 arrays
540 multimp = i_stok/ipari(18) + 1
541 call upgrade_multimp(nin,multimp,intbuf_tab)
542 endif
543 ! ------------
544!$OMP END SINGLE
545
546!$OMP BARRIER
547
548 ! ------------
549 my_address = i_stok_old + address_cand_n(itask+1) ! get the address for each thread
550 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
551 my_address = i_stok_old + address_cand_e(itask+1) ! get the address for each thread
552 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
553 ! ------------
554
555
556 call dealloc_1d_array(local_cand_n)
557 call dealloc_1d_array(local_cand_e)
558
559!$OMP BARRIER
560
561 ! Sort the candidates to ensure the same domain decomposition
562!$OMP SINGLE
563 ! ------------
564 my_size = cand_n_size(nthreads+1)
565 allocate(my_index(2*my_size))
566 allocate(sort_array(2,my_size))
567 allocate(save_array(2,my_size))
568
569 my_address = i_stok_old + address_cand_n(1) ! get the address of the first pair of candidate
570 sort_array(1,1:my_size) = intbuf_tab%cand_n(my_address+1:my_address+my_size)
571 my_address = i_stok_old + address_cand_e(1) ! get the address of the first pair of candidate
572 sort_array(2,1:my_size) = intbuf_tab%cand_e(my_address+1:my_address+my_size)
573 do i=1,my_size
574 my_index(i) = i
575 enddo
576 save_array(1:2,1:my_size) = sort_array(1:2,1:my_size)
577 mode = 0
578
579 call my_orders( mode,work,sort_array,my_index,my_size,2)
580 my_address = i_stok_old + address_cand_n(1) ! get the address of the first pair of candidate
581 do i=1,my_size
582 intbuf_tab%cand_n(my_address+i) = save_array(1,my_index(i))
583 enddo
584 my_address = i_stok_old + address_cand_e(1) ! get the address of the first pair of candidate
585 do i=1,my_size
586 intbuf_tab%cand_e(my_address+i) = save_array(2,my_index(i))
587 enddo
588 deallocate(my_index)
589 deallocate(sort_array)
590 deallocate(save_array)
591 ! ------------
592!$OMP END SINGLE
593
594!$OMP DO SCHEDULE(guided)
595 DO i=1,nsn
596 IF(iix(i)/=0)THEN
597 voxel(iix(i),iiy(i),iiz(i))=0
598 ENDIF
599 ENDDO
600!$OMP END DO
601
602!$OMP SINGLE
603 deallocate( cand_n_size,cand_e_size )
604 deallocate( address_cand_n,address_cand_e )
605!$OMP END SINGLE
606
607 RETURN
#define my_real
Definition cppsort.cpp:32
end diagonal values have been computed in the(sparse) matrix id.SOL
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
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
subroutine tagnod(ix, nix, nix1, nix2, numel, iparte, tagbuf, npart)
Definition tagnod.F:29
subroutine upgrade_multimp(ni, multimp_parameter, intbuf_tab)