OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
i25trivox1.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 i25trivox1 (nsn, irect, x, stfn, xyzm, nsv, ii_stok, eshift, bgapsmx, voxel, nbx, nby, nbz, nrtm, gap_s, gap_m, marge, nbinflg, mbinflg, ilev, msegtyp, igap, gap_s_l, gap_m_l, edge_l2, ledgmax, kremnode, remnode, iparts, npartns, lpartns, ielem, icode, iskew, drad, is_large_node, large_node, nb_large_nodes, dgapload, nrtmt, flag_removed_node, ielem_m, local_next_nod, iix, iiy, iiz, intbuf_tab, ipari, nin)

Function/Subroutine Documentation

◆ i25trivox1()

subroutine i25trivox1 ( integer nsn,
integer, dimension(4,*) irect,
x,
stfn,
xyzm,
integer, dimension(*) nsv,
integer ii_stok,
integer eshift,
bgapsmx,
integer, dimension(nbx+2,nby+2,nbz+2) voxel,
integer nbx,
integer nby,
integer nbz,
integer nrtm,
gap_s,
gap_m,
marge,
integer, dimension(*) nbinflg,
integer, dimension(*) mbinflg,
integer ilev,
integer, dimension(*) msegtyp,
integer igap,
gap_s_l,
gap_m_l,
edge_l2,
ledgmax,
integer, dimension(*) kremnode,
integer, dimension(*) remnode,
integer, dimension(*) iparts,
integer, dimension(*) npartns,
integer, dimension(*) lpartns,
integer, dimension(nrtm), intent(in) ielem,
integer, dimension(*) icode,
integer, dimension(*) iskew,
intent(in) drad,
integer, dimension(nsn) is_large_node,
integer, dimension(nsn) large_node,
integer nb_large_nodes,
intent(in) dgapload,
integer, intent(in) nrtmt,
logical, intent(in) flag_removed_node,
integer, dimension(2,nrtm), intent(in) ielem_m,
integer, dimension(nsn), intent(inout) local_next_nod,
integer, dimension(nsn), intent(inout) iix,
integer, dimension(nsn), intent(inout) iiy,
integer, dimension(nsn), intent(inout) iiz,
type(intbuf_struct_), intent(inout) intbuf_tab,
integer, dimension(npari), intent(inout) ipari,
integer, intent(in) nin )
Parameters
[in]flag_removed_nodeflag to remove some S node from the list of candidates
[in]nininterface index
[in,out]ipariinterface data
[in,out]intbuf_tabinterface data

Definition at line 32 of file i25trivox1.F.

46C============================================================================
47C M o d u l e s
48C-----------------------------------------------
49 USE tri7box
50 use array_mod
51 use intbufdef_mod
52C-----------------------------------------------
53C I m p l i c i t T y p e s
54C-----------------------------------------------
55#include "implicit_f.inc"
56C-----------------------------------------------
57C G l o b a l P a r a m e t e r s
58C-----------------------------------------------
59#include "mvsiz_p.inc"
60c parameter setting the size for the vector (orig version is 128)
61 INTEGER NVECSZ
62 parameter(nvecsz = mvsiz)
63C-----------------------------------------------
64C C o m m o n B l o c k s
65C-----------------------------------------------
66#include "com04_c.inc"
67#include "param_c.inc"
68C-----------------------------------------------
69C ROLE OF THE ROUTINE:
70C ===================
71C CLASSIFIES THE NODES INTO BOXES
72C SEARCHES FOR EACH FACET THE RELEVANT BOXES
73C SEARCH FOR CANDIDATES
74C-----------------------------------------------
75C D u m m y A r g u m e n t s
76C
77C NOM DESCRIPTION E/S
78C
79C IRECT(4,*) ARRAY OF CONEC FACETTES E
80C X(3,*) COORDONNEES NODALES E
81C NSV NOS SYSTEMES DES NODES E
82C Xmax larger abcisse existing e
83C XMAX largest order.existing E
84C Xmax larger existing side E
85C I_STOK storage level of the pairs
86C CANDIDATES impact E/S
87C CAND_N boites resultats nodes C CAND_E adresses des boites resultat elements
88C
89C Prov_n Provisional Cand_n (static variable in i7tri)
90C PROV_E CAND_E provisoire (variable static in i7tri)
91
92C VOXEL(ix,iy,iz) contains the local number of the first node of
93C the box
94C Local_next_nod (I) Next node in the same box (if /= 0)
95C Last_nod (i) Last node in the same box (if /= 0)
96C used only to go directly from the first
97C node at the last
98C-----------------------------------------------
99C D u m m y A r g u m e n t s
100C-----------------------------------------------
101 INTEGER ESHIFT,NSN,NRTM,IGAP,
102 . NBX,NBY,NBZ,
103 . NSV(*),
104 . IRECT(4,*), VOXEL(NBX+2,NBY+2,NBZ+2),II_STOK,
105 . NBINFLG(*),MBINFLG(*),ILEV,MSEGTYP(*),
106 . KREMNODE(*),REMNODE(*),
107 . IPARTS(*), NPARTNS(*), LPARTNS(*), ICODE(*), ISKEW(*)
108 LOGICAL, INTENT(in) :: FLAG_REMOVED_NODE !< flag to remove some S node from the list of candidates
109C REAL
110 my_real
111 . x(3,*),xyzm(6),stfn(*),gap_s(*),gap_m(*),
112 . gap_s_l(*),gap_m_l(*), edge_l2(*)
113 my_real
114 . ledgmax, marge, bgapsmx
115 my_real , INTENT(IN) :: drad, dgapload
116 INTEGER :: LARGE_NODE(NSN) ! list of nodes that have large research distance wrt solids
117 INTEGER :: IS_LARGE_NODE(NSN) ! tag nodes that have large research distance wrt solids
118 INTEGER :: NB_LARGE_NODES
119 INTEGER , INTENT(IN) :: NRTMT
120 INTEGER , INTENT(IN) :: IELEM_M(2,NRTM), IELEM(NRTM)
121 integer, intent(in) :: nin !< interface index
122 INTEGER, dimension(nsn), intent(inout) :: IIX,IIY,IIZ,LOCAL_NEXT_NOD
123 integer, dimension(npari), intent(inout) :: ipari !< interface data
124 type(intbuf_struct_), intent(inout) :: intbuf_tab !< interface data
125C-----------------------------------------------
126C L o c a l V a r i a b l e s
127C-----------------------------------------------
128 INTEGER I,J,
129 . NN,NE,K,L,J_STOK,JJ,
130 . PROV_N(MVSIZ),PROV_E(MVSIZ),
131 . M,IEM,IPM,IPS
132C REAL
133 my_real
134 . xs,ys,zs,sx,sy,sz,s2,
135 . xmin, xmax,ymin, ymax,zmin, zmax,
136 . xx1,xx2,xx3,xx4,yy1,yy2,yy3,yy4,zz1,zz2,zz3,zz4,
137 . d1x,d1y,d1z,d2x,d2y,d2z,dd1,dd2,d2,a2
138C
139 INTEGER, DIMENSION(:), ALLOCATABLE :: LAST_NOD
140 INTEGER IX,IY,IZ,M1,M2,M3,M4,
141 . IX1,IY1,IZ1,IX2,IY2,IZ2
142 my_real
143 . xminb,yminb,zminb,xmaxb,ymaxb,zmaxb,
144 . xmine,ymine,zmine,xmaxe,ymaxe,zmaxe,aaa
145 INTEGER FIRST,LAST
146 INTEGER, DIMENSION(:), ALLOCATABLE :: TAGNOD
147
148 integer , external :: omp_get_thread_num,omp_get_num_threads
149 integer :: itask,nthreads
150 integer :: my_old_size,my_address
151 integer :: local_i_stok,multimp
152 integer :: local_cand_n_size,local_cand_e_size
153
154 integer, dimension(:), allocatable, save :: cand_n_size,cand_e_size
155 integer, dimension(:), allocatable, save :: address_cand_n,address_cand_e
156 type(array_type_int_1d) :: local_cand_n
157 type(array_type_int_1d) :: local_cand_e
158
159 integer :: my_size,mode
160 integer, dimension(:), allocatable :: my_index
161 integer, dimension(:,:), allocatable :: sort_array,save_array
162 integer, dimension(70000) :: work
163C-----------------------------------------------
164 ! allocation of local arrays
165 itask = omp_get_thread_num()
166 nthreads = omp_get_num_threads()
167 local_cand_n_size = size(intbuf_tab%cand_n) / nthreads + 1
168 local_cand_e_size = size(intbuf_tab%cand_e) / nthreads + 1
169 local_i_stok = 0
170 local_cand_n%size_int_array_1d = local_cand_n_size
171 local_cand_e%size_int_array_1d = local_cand_e_size
172 call alloc_1d_array(local_cand_n)
173 call alloc_1d_array(local_cand_e)
174
175 xmin = xyzm(1)
176 ymin = xyzm(2)
177 zmin = xyzm(3)
178 xmax = xyzm(4)
179 ymax = xyzm(5)
180 zmax = xyzm(6)
181
182 xminb = xmin
183 yminb = ymin
184 zminb = zmin
185 xmaxb = xmax
186 ymaxb = ymax
187 zmaxb = zmax
188
189!$OMP MASTER
190 allocate( cand_n_size(nthreads+1),cand_e_size(nthreads+1) )
191 allocate( address_cand_n(nthreads+1),address_cand_e(nthreads+1) )
192 cand_n_size(1:nthreads+1) = 0
193 cand_e_size(1:nthreads+1) = 0
194 address_cand_n(1:nthreads+1) = 0
195 address_cand_e(1:nthreads+1) = 0
196 ALLOCATE(last_nod(nsn))
197C
198C Initial phase of construction of BPE and BPN moved from I7BUCE => I7TRI
199C
200
201C=======================================================================
202C 1 placing the nodes into the boxes
203C=======================================================================
204 DO i=1,nsn
205 iix(i)=0
206 iiy(i)=0
207 iiz(i)=0
208 IF(stfn(i) == zero)cycle
209 j=nsv(i)
210C Optimization // searches for nodes within xmin xmax of the
211C processor elements
212 IF(x(1,j) < xmin) cycle
213 IF(x(1,j) > xmax) cycle
214 IF(x(2,j) < ymin) cycle
215 IF(x(2,j) > ymax) cycle
216 IF(x(3,j) < zmin) cycle
217 IF(x(3,j) > zmax) cycle
218
219 iix(i)=int(nbx*(x(1,j)-xminb)/(xmaxb-xminb))
220 iiy(i)=int(nby*(x(2,j)-yminb)/(ymaxb-yminb))
221 iiz(i)=int(nbz*(x(3,j)-zminb)/(zmaxb-zminb))
222
223 iix(i)=max(1,2+min(nbx,iix(i)))
224 iiy(i)=max(1,2+min(nby,iiy(i)))
225 iiz(i)=max(1,2+min(nbz,iiz(i)))
226
227 first = voxel(iix(i),iiy(i),iiz(i))
228 IF(first == 0)THEN
229c Empty Cell
230 voxel(iix(i),iiy(i),iiz(i)) = i ! first
231 local_next_nod(i) = 0 ! last one
232 last_nod(i) = 0 ! no last
233 ELSEIF(last_nod(first) == 0)THEN
234c cell containing one node
235c add as next node
236 local_next_nod(first) = i ! next
237 last_nod(first) = i ! last
238 local_next_nod(i) = 0 ! last one
239 ELSE
240c
241c jump to the last node of the cell
242 last = last_nod(first) ! last node in this voxel
243 local_next_nod(last) = i ! next
244 last_nod(first) = i ! last
245 local_next_nod(i) = 0 ! last one
246 ENDIF
247 ENDDO
248!$OMP END MASTER
249
250!$OMP BARRIER
251
252C=======================================================================
253C 2 Searching boxes concerning each facet
254C and creation of candidates
255C=======================================================================
256 ALLOCATE( tagnod(numnod) )
257 tagnod(1:numnod) = 0
258C-----------------------------------------------
259 j_stok = 0
260!$OMP DO SCHEDULE(guided)
261 DO ne=1,nrtm
262 IF(ielem_m(2,ne) /=0) cycle
263
264c
265c It is possible to improve the algo by cutting the facet
266c in 2 (4,3,6,9 ...) if the facet is large in front of AAA and inclinee
267
268 m1 = irect(1,ne)
269 m2 = irect(2,ne)
270 m3 = irect(3,ne)
271 m4 = irect(4,ne)
272
273 IF(flag_removed_node)THEN
274 k = kremnode(ne)+1
275 l = kremnode(ne+1)
276 DO m=k,l
277 tagnod(remnode(m)) = 1
278 ENDDO
279 ENDIF
280
281 IF (msegtyp(ne)==0 .OR. msegtyp(ne)>nrtmt)THEN
282 ! ledgmax /=0 <=> inacti=5 or -1 and iddlevel=1 !
283 aaa = max(marge+max(bgapsmx+gap_m(ne)+dgapload,drad),ledgmax+bgapsmx+gap_m(ne)+dgapload)
284 ELSE
285 aaa = marge+max(bgapsmx+gap_m(ne)+dgapload,drad)
286 END IF
287
288
289 xx1=x(1,m1)
290 xx2=x(1,m2)
291 xx3=x(1,m3)
292 xx4=x(1,m4)
293 xmaxe=max(xx1,xx2,xx3,xx4)
294 xmine=min(xx1,xx2,xx3,xx4)
295
296 yy1=x(2,m1)
297 yy2=x(2,m2)
298 yy3=x(2,m3)
299 yy4=x(2,m4)
300 ymaxe=max(yy1,yy2,yy3,yy4)
301 ymine=min(yy1,yy2,yy3,yy4)
302
303 zz1=x(3,m1)
304 zz2=x(3,m2)
305 zz3=x(3,m3)
306 zz4=x(3,m4)
307 zmaxe=max(zz1,zz2,zz3,zz4)
308 zmine=min(zz1,zz2,zz3,zz4)
309
310
311c calculation of the surface (for future elimination of candidates)
312
313 sx = (yy3-yy1)*(zz4-zz2) - (zz3-zz1)*(yy4-yy2)
314 sy = (zz3-zz1)*(xx4-xx2) - (xx3-xx1)*(zz4-zz2)
315 sz = (xx3-xx1)*(yy4-yy2) - (yy3-yy1)*(xx4-xx2)
316 s2 = sx*sx + sy*sy + sz*sz
317
318c index of the voxels occupied by the facet
319
320 ix1=int(nbx*(xmine-aaa-xminb)/(xmaxb-xminb))
321 iy1=int(nby*(ymine-aaa-yminb)/(ymaxb-yminb))
322 iz1=int(nbz*(zmine-aaa-zminb)/(zmaxb-zminb))
323
324 ix1=max(1,2+min(nbx,ix1))
325 iy1=max(1,2+min(nby,iy1))
326 iz1=max(1,2+min(nbz,iz1))
327
328 ix2=int(nbx*(xmaxe+aaa-xminb)/(xmaxb-xminb))
329 iy2=int(nby*(ymaxe+aaa-yminb)/(ymaxb-yminb))
330 iz2=int(nbz*(zmaxe+aaa-zminb)/(zmaxb-zminb))
331
332 ix2=max(1,2+min(nbx,ix2))
333 iy2=max(1,2+min(nby,iy2))
334 iz2=max(1,2+min(nbz,iz2))
335
336 IF (msegtyp(ne)==0 .OR. msegtyp(ne)>nrtmt)THEN
337C Check for "large" nodes separately
338C if the current segment belongs to a solid or a coating shell
339 DO i = 1, nb_large_nodes
340 jj = large_node(i)
341 nn=nsv(jj)
342 IF(nn == m1)GOTO 400
343 IF(nn == m2)GOTO 400
344 IF(nn == m3)GOTO 400
345 IF(nn == m4)GOTO 400
346 IF(tagnod(nn) == 1)GOTO 400
347
348 xs = x(1,nn)
349 ys = x(2,nn)
350 zs = x(3,nn)
351c PMAX_GAP is a global overestimate penetration
352c NEED to communicate in SPMD
353c VMAXDT is a local overestimate of relative incremental displacement
354c NO need to communicate in SPMD
355
356 aaa = max(marge+max(gap_s(jj)+gap_m(ne)+dgapload,drad)+dgapload,edge_l2(jj)+gap_s(jj)+gap_m(ne)+dgapload)
357 IF(xs<=xmine-aaa)GOTO 400
358 IF(xs>=xmaxe+aaa)GOTO 400
359 IF(ys<=ymine-aaa)GOTO 400
360 IF(ys>=ymaxe+aaa)GOTO 400
361 IF(zs<=zmine-aaa)GOTO 400
362 IF(zs>=zmaxe+aaa)GOTO 400
363
364 iem=ielem(ne)
365 IF(iem/=0)THEN
366 ipm=iparts(iem)
367 ips=0
368 DO j=npartns(jj)+1,npartns(jj+1)
369 IF(lpartns(j)==ipm)THEN
370 ips=ipm
371 END IF
372 END DO
373 IF(ipm==ips) GOTO 400
374 END IF
375c END IF
376 d1x = xs - xx1
377 d1y = ys - yy1
378 d1z = zs - zz1
379 d2x = xs - xx2
380 d2y = ys - yy2
381 d2z = zs - zz2
382 dd1 = d1x*sx+d1y*sy+d1z*sz
383 dd2 = d2x*sx+d2y*sy+d2z*sz
384 IF(dd1*dd2 > zero)THEN
385 d2 = min(dd1*dd1,dd2*dd2)
386 a2 = aaa*aaa*s2
387 IF(d2 > a2)GOTO 400
388 ENDIF
389 j_stok = j_stok + 1
390 prov_n(j_stok) = jj
391 prov_e(j_stok) = ne
392 IF(j_stok == nvsiz)THEN
393 CALL i25sto(
394 1 nvsiz ,irect ,x ,nsv ,local_i_stok,
395 2 local_cand_n,local_cand_e ,marge ,
396 3 prov_n ,prov_e,eshift,nsn ,
397 4 nrtm ,gap_s ,gap_m ,nbinflg ,mbinflg,
398 5 ilev,msegtyp,igap ,gap_s_l,
399 6 gap_m_l,edge_l2,icode,iskew,drad ,
400 7 dgapload,nrtmt)
401 j_stok = 0
402 ENDIF
403 400 CONTINUE
404 ENDDO ! WHILE(JJ /= 0)
405 ENDIF ! solid or coating shell
406
407 DO iz = iz1,iz2
408 DO iy = iy1,iy2
409 DO ix = ix1,ix2
410
411cc nbpelem = nbpelem + 1
412
413 jj = voxel(ix,iy,iz)
414
415 DO WHILE(jj /= 0)
416
417cc nnpelem = nnpelem + 1
418
419 nn=nsv(jj)
420 IF(nn == m1)GOTO 300
421 IF(nn == m2)GOTO 300
422 IF(nn == m3)GOTO 300
423 IF(nn == m4)GOTO 300
424 IF(tagnod(nn) == 1)GOTO 300
425
426 xs = x(1,nn)
427 ys = x(2,nn)
428 zs = x(3,nn)
429c PMAX_GAP is a global overestimate penetration
430c NEED to communicate in SPMD
431c VMAXDT is a local overestimate of relative incremental displacement
432c NO need to communicate in SPMD
433
434 IF (msegtyp(ne)==0 .OR. msegtyp(ne)>nrtmt)THEN
435 IF(is_large_node(jj)==1) GOTO 300 ! node already checked before
436 ! LEDGMAX /=0 <=> INACTI=5 or -1 and IDDLEVEL=1 !
437 aaa = max(marge+max(gap_s(jj)+gap_m(ne)+dgapload,drad),edge_l2(jj)+gap_s(jj)+gap_m(ne)+dgapload)
438 ELSE
439 aaa = marge+max(gap_s(jj)+gap_m(ne)+dgapload,drad)
440 END IF
441
442 IF(xs<=xmine-aaa)GOTO 300
443 IF(xs>=xmaxe+aaa)GOTO 300
444 IF(ys<=ymine-aaa)GOTO 300
445 IF(ys>=ymaxe+aaa)GOTO 300
446 IF(zs<=zmine-aaa)GOTO 300
447 IF(zs>=zmaxe+aaa)GOTO 300
448
449 iem=ielem(ne)
450 IF(iem/=0)THEN
451 ipm=iparts(iem)
452 ips=0
453 DO j=npartns(jj)+1,npartns(jj+1)
454 IF(lpartns(j)==ipm)THEN
455 ips=ipm
456 END IF
457 END DO
458
459 IF(ipm==ips) GOTO 300
460 END IF
461c END IF
462
463c underestimation of the distance**2 for candidate elimination
464
465cc'nnr0pelem = nnr0pelem + 1
466
467 d1x = xs - xx1
468 d1y = ys - yy1
469 d1z = zs - zz1
470 d2x = xs - xx2
471 d2y = ys - yy2
472 d2z = zs - zz2
473 dd1 = d1x*sx+d1y*sy+d1z*sz
474 dd2 = d2x*sx+d2y*sy+d2z*sz
475 IF(dd1*dd2 > zero)THEN
476 d2 = min(dd1*dd1,dd2*dd2)
477 a2 = aaa*aaa*s2
478 IF(d2 > a2)GOTO 300
479 ENDIF
480
481 j_stok = j_stok + 1
482 prov_n(j_stok) = jj
483 prov_e(j_stok) = ne
484 IF(j_stok == nvsiz)THEN
485
486 CALL i25sto(
487 1 nvsiz ,irect ,x ,nsv ,local_i_stok,
488 2 local_cand_n,local_cand_e ,marge ,
489 3 prov_n ,prov_e,eshift,nsn ,
490 4 nrtm ,gap_s ,gap_m ,nbinflg ,mbinflg,
491 5 ilev,msegtyp,igap ,gap_s_l,
492 6 gap_m_l,edge_l2,icode,iskew,drad ,
493 7 dgapload,nrtmt)
494 j_stok = 0
495 ENDIF
496
497 300 CONTINUE
498
499 jj = local_next_nod(jj)
500
501 ENDDO ! WHILE(JJ /= 0)
502
503 ENDDO
504 ENDDO
505 ENDDO
506
507 IF(flag_removed_node)THEN
508 k = kremnode(ne)+1
509 l = kremnode(ne+1)
510 DO m=k,l
511 tagnod(remnode(m)) = 0
512 ENDDO
513 ENDIF
514
515 ENDDO
516!$OMP END DO
517
518!$OMP BARRIER
519C-------------------------------------------------------------------------
520C END OF SORTING
521C-------------------------------------------------------------------------
522 IF(j_stok/=0)CALL i25sto(
523 1 j_stok,irect ,x ,nsv ,local_i_stok,
524 2 local_cand_n,local_cand_e ,marge ,
525 3 prov_n ,prov_e,eshift,nsn ,
526 4 nrtm ,gap_s ,gap_m ,nbinflg ,mbinflg,
527 5 ilev,msegtyp,igap ,gap_s_l,
528 6 gap_m_l,edge_l2,icode,iskew,drad ,
529 7 dgapload,nrtmt)
530C 4 resetting the nodes in the boxes
531C=======================================================================
532 ! save the local number of candidates
533 cand_n_size(itask+1) = local_i_stok
534 cand_e_size(itask+1) = local_i_stok
535!$OMP BARRIER
536
537!$OMP SINGLE
538 ! compute the address for each threads & the total number of candidates
539 address_cand_n(1) = 0
540 address_cand_e(1) = 0
541 ! ------------
542 do i=1,nthreads
543 address_cand_n(i+1) = cand_n_size(i) + address_cand_n(i)
544 address_cand_e(i+1) = cand_e_size(i) + address_cand_e(i)
545
546 cand_n_size(nthreads+1) = cand_n_size(nthreads+1) + cand_n_size(i)
547 cand_e_size(nthreads+1) = cand_e_size(nthreads+1) + cand_e_size(i)
548 enddo
549 ! ------------
550
551 ! ------------
552 ! check the size of global cand_n & cand_e
553 my_old_size = ipari(18)*ipari(23)
554 if(cand_n_size(nthreads+1) > my_old_size) then ! total number of candidates > size of cand_n/e --> need to extend the 2 arrays
555 multimp = cand_n_size(nthreads+1)/ipari(18) + 1
556 call upgrade_multimp(nin,multimp,intbuf_tab)
557 endif
558 ii_stok = cand_n_size(nthreads+1) ! total number of cand_n/cand_e
559 ! ------------
560!$OMP END SINGLE
561
562 ! ------------
563 my_address = address_cand_n(itask+1) ! get the address for each thread
564 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
565 my_address = address_cand_e(itask+1) ! get the address for each thread
566 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
567 ! ------------
568
569 call dealloc_1d_array(local_cand_n)
570 call dealloc_1d_array(local_cand_e)
571
572!$OMP BARRIER
573
574
575
576
577 ! Sort the candidates to ensure the same domain decomposition
578!$OMP SINGLE
579 ! ------------
580 my_size = cand_n_size(nthreads+1)
581 allocate(my_index(2*my_size))
582 allocate(sort_array(2,my_size))
583 allocate(save_array(2,my_size))
584
585 my_address = address_cand_n(1) ! get the address of the first pair of candidate
586 sort_array(1,1:my_size) = intbuf_tab%cand_n(my_address+1:my_address+my_size)
587 my_address = address_cand_e(1) ! get the address of the first pair of candidate
588 sort_array(2,1:my_size) = intbuf_tab%cand_e(my_address+1:my_address+my_size)
589 do i=1,my_size
590 my_index(i) = i
591 enddo
592 save_array(1:2,1:my_size) = sort_array(1:2,1:my_size)
593 mode = 0
594
595 call my_orders( mode,work,sort_array,my_index,my_size,2)
596 my_address = address_cand_n(1) ! get the address of the first pair of candidate
597 do i=1,my_size
598 intbuf_tab%cand_n(my_address+i) = save_array(1,my_index(i))
599 enddo
600 my_address = address_cand_e(1) ! get the address of the first pair of candidate
601 do i=1,my_size
602 intbuf_tab%cand_e(my_address+i) = save_array(2,my_index(i))
603 enddo
604 deallocate(my_index)
605 deallocate(sort_array)
606 deallocate(save_array)
607 ! ------------
608!$OMP END SINGLE
609
610
611!$OMP DO SCHEDULE(guided)
612 DO i=1,nsn
613 IF(iix(i)/=0)THEN
614 voxel(iix(i),iiy(i),iiz(i))=0
615 ENDIF
616 ENDDO
617!$OMP END DO
618C
619!$OMP MASTER
620 DEALLOCATE(last_nod)
621 deallocate( cand_n_size,cand_e_size )
622 deallocate( address_cand_n,address_cand_e )
623!$OMP END MASTER
624 DEALLOCATE(tagnod)
625
626
627 RETURN
#define my_real
Definition cppsort.cpp:32
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
subroutine i25sto(j_stok, irect, x, nsv, local_i_stok, local_cand_n, local_cand_e, marge, prov_n, prov_e, eshift, nsn, nrtm, gap_s, gap_m, nbinflg, mbinflg, ilev, msegtyp, igap, gap_s_l, gap_m_l, edge_l2, icode, iskew, drad, dgapload, nrtmt)
Definition i25sto.F:42
subroutine tagnod(ix, nix, nix1, nix2, numel, iparte, tagbuf, npart)
Definition tagnod.F:29
subroutine upgrade_multimp(ni, multimp_parameter, intbuf_tab)