OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
sfac_distrib_distentry.F
Go to the documentation of this file.
1C
2C This file is part of MUMPS 5.5.1, released
3C on Tue Jul 12 13:17:24 UTC 2022
4C
5C
6C Copyright 1991-2022 CERFACS, CNRS, ENS Lyon, INP Toulouse, Inria,
7C Mumps Technologies, University of Bordeaux.
8C
9C This version of MUMPS is provided to you free of charge. It is
10C released under the CeCILL-C license
11C (see doc/CeCILL-C_V1-en.txt, doc/CeCILL-C_V1-fr.txt, and
12C https://cecill.info/licences/Licence_CeCILL-C_V1-en.html)
13C
15 & ( n, mapping, nnz, irn, jcn, procnode, step,
16 & slavef, perm, fils,
17 & rg2l, keep,keep8, mblock, nblock, nprow, npcol )
19 IMPLICIT NONE
20 INTEGER N, SLAVEF, MBLOCK, NBLOCK, NPROW, NPCOL
21 iNTEGER(8) :: NNZ
22 INTEGER KEEP(500)
23 INTEGER(8) KEEP8(150)
24 INTEGER IRN( NNZ ), JCN( NNZ )
25 INTEGER MAPPING( NNZ ), STEP( N )
26 INTEGER PROCNODE( KEEP(28) ), PERM( N ), FILS( N ), RG2L( N )
27 INTEGER MUMPS_PROCNODE, MUMPS_TYPENODE
28 EXTERNAL mumps_procnode, mumps_typenode
29 INTEGER K4, IOLD, JOLD, INEW, JNEW, ISEND, JSEND, IARR, INODE
30 INTEGER(8) :: K8
31 INTEGER TYPE_NODE, DEST
32 INTEGER IPOSROOT, JPOSROOT, IROW_GRID, JCOL_GRID
33 inode = keep(38)
34 k4 = 1
35 DO WHILE ( inode .GT. 0 )
36 rg2l( inode ) = k4
37 inode = fils( inode )
38 k4 = k4 + 1
39 END DO
40 DO k8 = 1_8, nnz
41 iold = irn( k8 )
42 jold = jcn( k8 )
43 IF ( iold .GT. n .OR. iold .LT. 1 .OR.
44 & jold .GT. n .OR. jold .LT. 1 ) THEN
45 mapping( k8 ) = -1
46 cycle
47 END IF
48 IF ( iold .eq. jold ) THEN
49 isend = iold
50 jsend = jold
51 ELSE
52 inew = perm( iold )
53 jnew = perm( jold )
54 IF ( inew .LT. jnew ) THEN
55 isend = iold
56 IF ( keep(50) .ne. 0 ) isend = -iold
57 jsend = jold
58 ELSE
59 isend = -jold
60 jsend = iold
61 END IF
62 END IF
63 iarr = abs( isend )
64 type_node = mumps_typenode( procnode(abs(step(iarr))),
65 & keep(199) )
66 IF ( type_node .eq. 1 .or. type_node .eq. 2 ) THEN
67 IF ( keep(46) .eq. 0 ) THEN
68 dest = mumps_procnode( procnode(abs(step(iarr))),
69 & keep(199) ) + 1
70 ELSE
71 dest = mumps_procnode( procnode(abs(step(iarr))),
72 & keep(199) )
73 END IF
74 ELSE
75 IF ( isend .LT. 0 ) THEN
76 iposroot = rg2l( jsend )
77 jposroot = rg2l( iarr )
78 ELSE
79 iposroot = rg2l( iarr )
80 jposroot = rg2l( jsend )
81 END IF
82 irow_grid = mod( ( iposroot - 1 )/mblock, nprow )
83 jcol_grid = mod( ( jposroot - 1 )/nblock, npcol )
84 IF ( keep( 46 ) .eq. 0 ) THEN
85 dest = irow_grid * npcol + jcol_grid + 1
86 ELSE
87 dest = irow_grid * npcol + jcol_grid
88 END IF
89 END IF
90 mapping( k8 ) = dest
91 END DO
92 RETURN
93 END SUBROUTINE smumps_build_mapping
95 & N, NZ_loc8, id,
96 & DBLARR, LDBLARR, INTARR, LINTARR,
97 & PTRAIW, PTRARW, KEEP,KEEP8, MYID, COMM, NBRECORDS,
98 &
99 & A, LA, root, PROCNODE_STEPS, SLAVEF, PERM, STEP,
100 & ICNTL, INFO, NSEND8, NLOCAL8,
101 & ISTEP_TO_INIV2, CANDIDATES
102 & )
103!$ USE OMP_LIB
105 IMPLICIT NONE
106 INTEGER N
107 INTEGER(8) :: NZ_loc8
108 TYPE (SMUMPS_STRUC) :: id
109 INTEGER(8) :: LDBLARR, LINTARR
110 REAL DBLARR( LDBLARR )
111 INTEGER INTARR( LINTARR )
112 INTEGER(8), INTENT(IN) :: PTRAIW( N ), PTRARW( N )
113 INTEGER KEEP(500)
114 INTEGER(8) KEEP8(150)
115 INTEGER MYID, COMM, NBRECORDS
116 INTEGER(8) :: LA
117 INTEGER SLAVEF
118 INTEGER ISTEP_TO_INIV2(KEEP(71))
119 INTEGER CANDIDATES(SLAVEF+1, max(1,KEEP(56)))
120 REAL A( LA )
121 TYPE (SMUMPS_ROOT_STRUC) :: root
122 INTEGER PROCNODE_STEPS(KEEP(28)), PERM( N ), STEP( N )
123 INTEGER INFO( 80 ), ICNTL(60)
124 INTEGER MUMPS_PROCNODE, MUMPS_TYPENODE, numroc,
126 EXTERNAL mumps_procnode, mumps_typenode, numroc,
128 include 'mumps_tags.h'
129 include 'mpif.h'
130 INTEGER :: IERR, MSGSOU
131 INTEGER :: STATUS(MPI_STATUS_SIZE)
132 REAL ZERO
133 PARAMETER( ZERO = 0.0e0 )
134 INTEGER, ALLOCATABLE, DIMENSION(:,:) :: IW4
135 INTEGER END_MSG_2_RECV
136 INTEGER I
137 INTEGER(8) :: I18, IA8
138 INTEGER(8) :: K8
139 INTEGER TYPE_NODE, DEST, DEST_SHR
140 INTEGER IOLD, JOLD, IARR, ISEND, JSEND
141 INTEGER ISEND_SHR, JSEND_SHR
142 INTEGER allocok, TYPESPLIT, T4MASTER, INIV2, NCAND
143 LOGICAL T4_MASTER_CONCERNED, EARLYT3ROOTINS
144 REAL VAL, VAL_SHR
145 INTEGER(8) :: PTR_ROOT
146 INTEGER LOCAL_M, LOCAL_N, ARROW_ROOT
147 INTEGER IROW_GRID, JCOL_GRID, IPOSROOT, JPOSROOT,
148 & ilocroot, jlocroot
149 INTEGER MP,LP
150 INTEGER KPROBE, FREQPROBE
151 INTEGER(8) :: IS18, IIW8, IS8, IAS8
152 INTEGER ISHIFT
153 INTEGER, ALLOCATABLE, DIMENSION(:,:,:) :: BUFI
154 REAL, ALLOCATABLE, DIMENSION(:,:,:) :: BUFR
155 INTEGER, ALLOCATABLE, DIMENSION(:) :: BUFRECI
156 REAL, ALLOCATABLE, DIMENSION(:) :: BUFRECR
157 INTEGER, ALLOCATABLE, DIMENSION(:) :: IACT, IREQI, IREQR
158 LOGICAL, ALLOCATABLE, DIMENSION(:) :: SEND_ACTIVE
159 LOGICAL :: FLAG
160 INTEGER(8), INTENT(OUT) :: NSEND8, NLOCAL8
161 INTEGER MASTER_NODE, ISTEP
162 LOGICAL :: DOIT, OMP_FLAG, OMP_FLAG_P
163 INTEGER NOMP, NOMP_P, IOMP, P2
164 nsend8 = 0_8
165 nlocal8 = 0_8
166 lp = icntl(1)
167 mp = icntl(2)
168 end_msg_2_recv = slavef
169 ALLOCATE( iact(slavef), stat=allocok)
170 IF ( allocok .GT. 0 ) THEN
171 IF ( lp > 0 ) THEN
172 WRITE(lp,*)
173 & '** Error allocating IACT in matrix distribution'
174 END IF
175 info(1) = -13
176 info(2) = slavef
177 GOTO 20
178 END IF
179 ALLOCATE( ireqi(slavef), stat=allocok)
180 IF ( allocok .GT. 0 ) THEN
181 IF ( lp > 0 ) THEN
182 WRITE(lp,*)
183 & '** Error allocating IREQI in matrix distribution'
184 END IF
185 info(1) = -13
186 info(2) = slavef
187 GOTO 20
188 END IF
189 ALLOCATE( ireqr(slavef), stat=allocok)
190 IF ( allocok .GT. 0 ) THEN
191 IF ( lp > 0 ) THEN
192 WRITE(lp,*)
193 & '** Error allocating IREQR in matrix distribution'
194 END IF
195 info(1) = -13
196 info(2) = slavef
197 GOTO 20
198 END IF
199 ALLOCATE( send_active(slavef), stat=allocok)
200 IF ( allocok .GT. 0 ) THEN
201 IF ( lp > 0 ) THEN
202 WRITE(lp,*)
203 & '** Error allocating SEND_ACTIVE in matrix distribution'
204 END IF
205 info(1) = -13
206 info(2) = slavef
207 GOTO 20
208 END IF
209 ALLOCATE( bufi( nbrecords * 2 + 1, 2, slavef ), stat=allocok)
210 IF ( allocok .GT. 0 ) THEN
211 IF ( lp > 0 ) THEN
212 WRITE(lp,*)
213 & '** Error allocating int buffer for matrix distribution'
214 END IF
215 info(1) = -13
216 info(2) = ( nbrecords * 2 + 1 ) * slavef * 2
217 GOTO 20
218 END IF
219 ALLOCATE( bufr( nbrecords, 2, slavef), stat = allocok)
220 IF ( allocok .GT. 0 ) THEN
221 IF ( lp > 0 ) THEN
222 WRITE(lp,*)
223 & '** Error allocating real buffer for matrix distribution'
224 END IF
225 info(1) = -13
226 info(2) = nbrecords * slavef * 2
227 GOTO 20
228 END IF
229 ALLOCATE( bufreci( nbrecords * 2 + 1 ), stat = allocok )
230 IF ( allocok .GT. 0 ) THEN
231 IF ( lp > 0 ) THEN
232 WRITE(lp,*)
233 & '** Error allocating int recv buffer for matrix distribution'
234 END IF
235 info(1) = -13
236 info(2) = nbrecords * 2 + 1
237 GOTO 20
238 END IF
239 ALLOCATE( bufrecr( nbrecords ), stat = allocok )
240 IF ( allocok .GT. 0 ) THEN
241 IF ( lp > 0 ) THEN
242 WRITE(lp,*)
243 & '** Error allocating int recv buffer for matrix distribution'
244 END IF
245 info(1) = -13
246 info(2) = nbrecords
247 GOTO 20
248 END IF
249 ALLOCATE( iw4( n, 2 ), stat = allocok )
250 IF ( allocok .GT. 0 ) THEN
251 WRITE(lp,*) '** Error allocating IW4 for matrix distribution'
252 info(1) = -13
253 info(2) = n * 2
254 END IF
255 20 CONTINUE
256 CALL mumps_propinfo( icntl, info, comm, myid )
257 IF ( info(1) .LT. 0 ) GOTO 100
258 arrow_root = 0
259 DO i = 1, n
260 i18 = ptraiw( i )
261 ia8 = ptrarw( i )
262 IF ( ia8 .GT. 0_8 ) THEN
263 dblarr( ia8 ) = zero
264 iw4( i, 1 ) = intarr( i18 )
265 iw4( i, 2 ) = -intarr( i18 + 1_8 )
266 intarr( i18 + 2_8 ) = i
267 END IF
268 END DO
269 earlyt3rootins = keep(200) .EQ.0
270 & .OR. ( keep(200) .LT. 0 .AND. keep(400) .EQ. 0 )
271 IF ( keep(38) .NE. 0 .AND. earlyt3rootins ) THEN
272 CALL smumps_get_root_info(root,local_m, local_n, ptr_root, la)
273 CALL smumps_set_root_to_zero(root, keep, a, la)
274 ELSE
275 local_m = -19999; local_n = -29999; ptr_root = -99999_8
276 END IF
277 DO i = 1, slavef
278 bufi( 1, 1, i ) = 0
279 END DO
280 DO i = 1, slavef
281 bufi( 1, 2, i ) = 0
282 END DO
283 DO i = 1, slavef
284 send_active( i ) = .false.
285 iact( i ) = 1
286 END DO
287 kprobe = 0
288 freqprobe = max(1,nbrecords/10)
289 IF (slavef .EQ. 1) freqprobe = huge(freqprobe)
290 nomp = 1
291!$ NOMP=omp_get_max_threads()
292 omp_flag = keep(399).EQ.1 .AND. nomp .GE.2 .AND. slavef.EQ.1
293!$OMP PARALLEL PRIVATE( K8, I, DEST,
294!$omp& t4master, t4_master_concerned,
295!$OMP& INIV2, NCAND, IROW_GRID, JCOL_GRID, IPOSROOT, JPOSROOT,
296!$OMP& ILOCROOT, JLOCROOT,
297!$OMP& TYPE_NODE, TYPESPLIT, MASTER_NODE,
298!$OMP& IA8, ISHIFT, IIW8, IS18, IS8, IAS8, VAL,
299!$OMP& IARR, ISTEP, ISEND, JSEND,
300!$OMP& IOLD, JOLD, IOMP, DOIT, P2, NOMP_P, OMP_FLAG_P )
301!$OMP& REDUCTION(+:NSEND8, NLOCAL8, ARROW_ROOT) IF (OMP_FLAG)
302 iomp=0
303!$ IOMP=omp_get_thread_num()
304 nomp_p=1
305!$ NOMP_P=omp_get_num_threads()
306 omp_flag_p = .false.
307!$ OMP_FLAG_P = OMP_FLAG .AND. NOMP_P .GT. 1
308 IF (omp_flag_p) THEN
309 IF ( nomp_p .GE. 16 ) THEN
310 nomp_p=16
311 p2 = 4
312 ELSE IF (nomp_p.GE.8) THEN
313 nomp_p=8
314 p2 = 3
315 ELSE IF (nomp_p.GE.4) THEN
316 nomp_p=4
317 p2 = 2
318 ELSE IF (nomp_p.GE.2) THEN
319 nomp_p=2
320 p2 = 1
321 ENDIF
322 ELSE
323 nomp_p = 1
324 p2 = 0
325 ENDIF
326 IF ( iomp .LT. nomp_p ) THEN
327 DO k8 = 1_8, nz_loc8
328 IF ( slavef .GT. 1 ) THEN
329!$OMP MASTER
330 kprobe = kprobe + 1
331 IF ( kprobe .eq. freqprobe ) THEN
332 kprobe = 0
333 CALL mpi_iprobe( mpi_any_source, arr_int, comm,
334 & flag, status, ierr )
335 IF ( flag ) THEN
336 msgsou = status( mpi_source )
337 CALL mpi_recv( bufreci(1), nbrecords * 2 + 1,
338 & mpi_integer,
339 & msgsou, arr_int, comm, status, ierr )
340 CALL mpi_recv( bufrecr(1), nbrecords,
341 & mpi_real,
342 & msgsou, arr_real, comm, status, ierr )
344 & bufreci, bufrecr, nbrecords, n, iw4(1,1),
345 & keep,keep8, local_m, local_n, root, ptr_root,
346 & a, la,
347 & end_msg_2_recv, myid, procnode_steps, slavef,
348 & ptraiw, ptrarw, perm, step,
349 & intarr, lintarr, dblarr, ldblarr
350 & )
351 END IF
352 END IF
353!$OMP END MASTER
354 ENDIF
355 iold = id%IRN_loc(k8)
356 jold = id%JCN_loc(k8)
357 IF ( (iold.GT.n).OR.(jold.GT.n).OR.(iold.LT.1)
358 & .OR.(jold.LT.1) ) THEN
359 cycle
360 ENDIF
361 IF (omp_flag_p) THEN
362 IF (iold.EQ.jold) THEN
363 iarr = iold
364 ELSE IF (perm(iold).LT.perm(jold)) THEN
365 iarr = iold
366 ELSE
367 iarr = jold
368 ENDIF
369 doit = ( iomp .EQ. ibits(iarr, p2-1, p2))
370 ELSE
371 doit = .true.
372 ENDIF
373 IF (doit) THEN
374 IF (iold.EQ.jold) THEN
375 isend = iold
376 jsend = iold
377 iarr = iold
378 ELSE IF (perm(iold).LT.perm(jold)) THEN
379 iarr = iold
380 IF ( keep(50) .NE. 0 ) THEN
381 isend = -iold
382 ELSE
383 isend = iold
384 ENDIF
385 jsend = jold
386 ELSE
387 iarr = jold
388 isend = -jold
389 jsend = iold
390 ENDIF
391 istep = abs(step(iarr))
392 CALL mumps_typeandprocnode( type_node, master_node,
393 & procnode_steps(istep), keep(199) )
394 t4_master_concerned = .false.
395 t4master = -9999
396 val = id%A_loc(k8)
397 IF ((keep(52).EQ.7).OR.(keep(52).EQ.8)) THEN
398 val = val * id%ROWSCA(iold)*id%COLSCA(jold)
399 ENDIF
400 IF ( type_node .eq. 1 ) THEN
401 dest = master_node
402 IF (dest.EQ.myid) THEN
403 nlocal8 = nlocal8 + 1_8
404 IF (isend.EQ.jsend) THEN
405 ia8 = ptrarw(isend)
406 dblarr(ia8) = dblarr(ia8) + val
407 ELSE IF (isend.GE.0) THEN
408 is18 = ptraiw(iarr)
409 ishift = intarr(is18) + iw4(iarr,2)
410 intarr(is18+ishift+2) = jsend
411 dblarr(ptrarw(iarr)+ishift) = val
412 iw4(iarr,2) = iw4(iarr,2) - 1
413 ELSE
414 ishift = iw4(iarr,1)
415 intarr(ptraiw(iarr)+ishift+2) = jsend
416 dblarr(ptrarw(iarr)+ishift) = val
417 iw4(iarr,1) = iw4(iarr,1) - 1
418 IF ( iw4(iarr,1) .EQ. 0
419 & .AND. step(iarr) > 0 ) THEN
420 CALL smumps_quick_sort_arrowheads( n, perm,
421 & intarr( ptraiw(iarr) + 3 ),
422 & dblarr( ptrarw(iarr) + 1 ),
423 & intarr( ptraiw(iarr) ), 1,
424 & intarr( ptraiw(iarr) ) )
425 END IF
426 ENDIF
427 cycle
428 ENDIF
429 ELSE IF ( type_node .eq. 2 ) THEN
430 IF ( isend .LT. 0 ) THEN
431 dest = -1
432 ELSE
433 dest = master_node
434 END IF
435 iniv2 = istep_to_iniv2(istep)
436 IF ( keep(79) .GT. 0) THEN
437 typesplit = mumps_typesplit( procnode_steps(istep),
438 & keep(199) )
439 IF ( (typesplit.EQ.5).OR.(typesplit.EQ.6)) THEN
440 t4_master_concerned = .true.
441 t4master=candidates(candidates(slavef+1,iniv2)+1,iniv2)
442 ENDIF
443 ENDIF
444 ELSE
445 arrow_root = arrow_root + 1
446 IF (earlyt3rootins) THEN
447 IF ( isend < 0 ) THEN
448 iposroot = root%RG2L_ROW(jsend)
449 jposroot = root%RG2L_ROW(iarr )
450 ELSE
451 iposroot = root%RG2L_ROW(iarr )
452 jposroot = root%RG2L_ROW(jsend)
453 END IF
454 irow_grid = mod( ( iposroot-1 )/root%MBLOCK, root%NPROW )
455 jcol_grid = mod( ( jposroot-1 )/root%NBLOCK, root%NPCOL )
456 dest = irow_grid * root%NPCOL + jcol_grid
457 ELSE
458 dest = -2
459 ENDIF
460 IF ( omp_flag_p ) THEN
461 IF ( earlyt3rootins ) THEN
462 ilocroot = root%MBLOCK * ( ( iposroot - 1 ) /
463 & ( root%MBLOCK * root%NPROW ) )
464 & + mod( iposroot - 1, root%MBLOCK ) + 1
465 jlocroot = root%NBLOCK * ( ( jposroot - 1 ) /
466 & ( root%NBLOCK * root%NPCOL ) )
467 & + mod( jposroot - 1, root%NBLOCK ) + 1
468 IF (keep(60)==0) THEN
469 a( ptr_root + int(jlocroot-1,8) * int(local_m,8)
470 & + int(ilocroot-1,8)) = a( ptr_root
471 & + int(jlocroot - 1,8) * int(local_m,8)
472 & + int(ilocroot - 1,8) )
473 & + val
474 ELSE
475 root%SCHUR_POINTER( int(jlocroot-1,8)
476 & * int(root%SCHUR_LLD,8)
477 & + int(ilocroot,8) )
478 & = root%SCHUR_POINTER( int(jlocroot - 1,8)
479 & * int(root%SCHUR_LLD,8)
480 & + int(ilocroot,8))
481 & + val
482 ENDIF
483 ELSE
484 IF (isend.EQ.jsend) THEN
485 ia8 = ptrarw(isend)
486 dblarr(ia8) = dblarr(ia8) + val
487 ELSE IF (isend.GE.0) THEN
488 is18 = ptraiw(iarr)
489 ishift = intarr(is18) + iw4(iarr,2)
490 iw4(iarr,2) = iw4(iarr,2) - 1
491 iiw8 = is18 + ishift + 2
492 intarr(iiw8) = jsend
493 is8 = ptrarw(iarr)
494 ias8 = is8 + ishift
495 dblarr(ias8) = val
496 ELSE
497 is8 = ptraiw(iarr)+iw4(iarr,1)+2
498 intarr(is8) = jsend
499 ias8 = ptrarw(iarr)+iw4(iarr,1)
500 iw4(iarr,1) = iw4(iarr,1) - 1
501 dblarr(ias8) = val
502 IF ( iw4(iarr,1) .EQ. 0
503 & .AND. step(iarr) > 0 ) THEN
504 CALL smumps_quick_sort_arrowheads( n, perm,
505 & intarr( ptraiw(iarr) + 3 ),
506 & dblarr( ptrarw(iarr) + 1 ),
507 & intarr( ptraiw(iarr) ), 1,
508 & intarr( ptraiw(iarr) ) )
509 END IF
510 ENDIF
511 ENDIF
512 cycle
513 ENDIF
514 END IF
515 IF (dest .eq. -1) THEN
516 nlocal8 = nlocal8 + 1_8
517 nsend8 = nsend8 + int(slavef -1,8)
518 ELSE IF (dest .EQ. -2) THEN
519 nlocal8 = nlocal8 + 1_8
520 nsend8 = nsend8 + int(slavef -1,8)
521 ELSE
522 IF (dest .eq.myid ) THEN
523 nlocal8 = nlocal8 + 1_8
524 ELSE
525 nsend8 = nsend8 + 1_8
526 ENDIF
527 ENDIF
528 IF ( dest.EQ.-1) THEN
529 iniv2 = istep_to_iniv2(istep)
530 ncand = candidates(slavef+1,iniv2)
531 IF (keep(79) .GT. 0) THEN
532 DO i=1, slavef
533 dest=candidates(i,iniv2)
534 IF (dest.LT.0) EXIT
535 IF (i.EQ.ncand+1) cycle
536 dest_shr=dest;isend_shr=isend
537 jsend_shr=jsend;val_shr=val
539 ENDDO
540 ELSE
541 DO i=1, ncand
542 dest=candidates(i,iniv2)
543 dest_shr=dest;isend_shr=isend
544 jsend_shr=jsend;val_shr=val
546 ENDDO
547 ENDIF
548 dest=master_node
549 dest_shr=dest;isend_shr=isend
550 jsend_shr=jsend;val_shr=val
552 IF (t4_master_concerned) THEN
553 dest = t4master
554 dest_shr=dest;isend_shr=isend
555 jsend_shr=jsend;val_shr=val
557 ENDIF
558 ELSE IF (dest .GE. 0) THEN
559 dest_shr=dest;isend_shr=isend
560 jsend_shr=jsend;val_shr=val
562 IF (t4_master_concerned) THEN
563 dest = t4master
564 dest_shr=dest;isend_shr=isend
565 jsend_shr=jsend;val_shr=val
567 ENDIF
568 ELSE IF (dest .EQ. -2) THEN
569 DO i = 0, slavef-1
570 dest=i
571 dest_shr=dest;isend_shr=isend
572 jsend_shr=jsend;val_shr=val
574 ENDDO
575 ENDIF
576 ENDIF
577 END DO
578 ENDIF
579!$OMP END PARALLEL
580 dest_shr = -3
582 DO WHILE ( end_msg_2_recv .NE. 0 )
583 CALL mpi_recv( bufreci(1), nbrecords * 2 + 1, mpi_integer,
584 & mpi_any_source, arr_int, comm, status, ierr )
585 msgsou = status( mpi_source )
586 CALL mpi_recv( bufrecr(1), nbrecords, mpi_real,
587 & msgsou, arr_real, comm, status, ierr )
589 & bufreci, bufrecr, nbrecords, n, iw4(1,1),
590 & keep,keep8, local_m, local_n, root, ptr_root,
591 & a, la,
592 & end_msg_2_recv, myid, procnode_steps, slavef,
593 & ptraiw, ptrarw, perm, step,
594 & intarr, lintarr, dblarr, ldblarr
595 & )
596 END DO
597 DO i = 1, slavef
598 IF ( send_active( i ) ) THEN
599 CALL mpi_wait( ireqi( i ), status, ierr )
600 CALL mpi_wait( ireqr( i ), status, ierr )
601 END IF
602 END DO
603 keep(49) = arrow_root
604 100 CONTINUE
605 IF (ALLOCATED(iw4)) DEALLOCATE( iw4 )
606 IF (ALLOCATED(bufi)) DEALLOCATE( bufi )
607 IF (ALLOCATED(bufr)) DEALLOCATE( bufr )
608 IF (ALLOCATED(bufreci)) DEALLOCATE( bufreci )
609 IF (ALLOCATED(bufrecr)) DEALLOCATE( bufrecr )
610 IF (ALLOCATED(iact)) DEALLOCATE( iact )
611 IF (ALLOCATED(ireqi)) DEALLOCATE( ireqi )
612 IF (ALLOCATED(ireqr)) DEALLOCATE( ireqr )
613 IF (ALLOCATED(send_active)) DEALLOCATE( send_active )
614 RETURN
615 CONTAINS
617 IMPLICIT NONE
618 INTEGER ISLAVE, IBEG, IEND, NBREC, IREQ
619 INTEGER TAILLE_SEND_I, TAILLE_SEND_R
620 LOGICAL SEND_LOCAL
621 IF ( dest_shr .eq. -3 ) THEN
622 ibeg = 1
623 iend = slavef
624 ELSE
625 ibeg = dest_shr + 1
626 iend = dest_shr + 1
627 END IF
628 send_local = .false.
629 DO islave = ibeg, iend
630 nbrec = bufi(1,iact(islave),islave)
631 IF ( dest_shr .eq. -3 ) THEN
632 bufi(1,iact(islave),islave) = - nbrec
633 END IF
634 IF ( dest_shr .eq. -3 .or. nbrec + 1 > nbrecords ) THEN
635 DO WHILE ( send_active( islave ) )
636 CALL mpi_test( ireqr( islave ), flag, status, ierr )
637 IF ( .NOT. flag ) THEN
638 CALL mpi_iprobe( mpi_any_source, arr_int, comm,
639 & flag, status, ierr )
640 IF ( flag ) THEN
641 msgsou = status(mpi_source)
642 CALL mpi_recv( bufreci(1), 2*nbrecords+1,
643 & mpi_integer, msgsou, arr_int, comm,
644 & status, ierr )
645 CALL mpi_recv( bufrecr(1), nbrecords,
646 & mpi_real, msgsou,
647 & arr_real, comm, status, ierr )
649 & bufreci, bufrecr, nbrecords, n, iw4(1,1),
650 & keep,keep8, local_m, local_n, root, ptr_root,
651 & a, la,
652 & end_msg_2_recv, myid, procnode_steps, slavef,
653 & ptraiw, ptrarw, perm, step,
654 & intarr, lintarr, dblarr, ldblarr
655 & )
656 END IF
657 ELSE
658 CALL mpi_wait( ireqi( islave ), status, ierr )
659 send_active( islave ) = .false.
660 END IF
661 END DO
662 IF ( islave - 1 .ne. myid ) THEN
663 taille_send_i = nbrec * 2 + 1
664 taille_send_r = nbrec
665 CALL mpi_isend( bufi(1, iact(islave), islave ),
666 & taille_send_i,
667 & mpi_integer, islave - 1, arr_int, comm,
668 & ireqi( islave ), ierr )
669 CALL mpi_isend( bufr(1, iact(islave), islave ),
670 & taille_send_r,
671 & mpi_real, islave - 1, arr_real, comm,
672 & ireqr( islave ), ierr )
673 send_active( islave ) = .true.
674 ELSE
675 send_local = .true.
676 END IF
677 iact( islave ) = 3 - iact( islave )
678 bufi( 1, iact( islave ), islave ) = 0
679 END IF
680 IF ( dest_shr .ne. -3 ) THEN
681 ireq = bufi(1,iact(islave),islave) + 1
682 bufi(1,iact(islave),islave) = ireq
683 bufi(ireq*2,iact(islave),islave) = isend_shr
684 bufi(ireq*2+1,iact(islave),islave) = jsend_shr
685 bufr(ireq,iact(islave),islave ) = val_shr
686 END IF
687 END DO
688 IF ( send_local ) THEN
689 islave = myid + 1
691 & bufi(1,3-iact(islave),islave),
692 & bufr(1,3-iact(islave),islave),
693 & nbrecords, n, iw4(1,1),
694 & keep,keep8, local_m, local_n, root, ptr_root,
695 & a, la,
696 & end_msg_2_recv, myid, procnode_steps, slavef,
697 & ptraiw, ptrarw, perm, step,
698 & intarr, lintarr, dblarr, ldblarr
699 & )
700 END IF
701 RETURN
702 END SUBROUTINE smumps_dist_fill_buffer
703 END SUBROUTINE smumps_redistribution
705 & ( bufi, bufr, nbrecords, n, iw4,
706 & keep,keep8, local_m, local_n, root, ptr_root, a, la,
707 & end_msg_2_recv, myid, procnode_steps,
708 & slavef,
709 & ptraiw, ptrarw, perm, step,
710 & intarr, lintarr, dblarr, ldblarr )
711 USE smumps_struc_def, ONLY : smumps_root_struc
712 IMPLICIT NONE
713 TYPE (SMUMPS_ROOT_STRUC) :: root
714 INTEGER NBRECORDS, N, MYID, SLAVEF
715 INTEGER BUFI( NBRECORDS * 2 + 1 )
716 REAL BUFR( NBRECORDS )
717 INTEGER IW4( N, 2 )
718 INTEGER KEEP(500)
719 INTEGER(8) KEEP8(150)
720 INTEGER END_MSG_2_RECV
721 INTEGER(8) :: PTRAIW( N ), PTRARW( N )
722 INTEGER :: PERM( N ), STEP( N )
723 INTEGER PROCNODE_STEPS( KEEP(28) )
724 INTEGER(8), INTENT(IN) :: LINTARR, LDBLARR
725 INTEGER INTARR( LINTARR )
726 INTEGER LOCAL_M, LOCAL_N
727 INTEGER(8) :: PTR_ROOT, LA
728 REAL A( LA ), DBLARR( LDBLARR )
729 INTEGER MUMPS_TYPENODE, MUMPS_PROCNODE
730 EXTERNAL mumps_typenode, mumps_procnode
731 INTEGER IREC, NB_REC, NODE_TYPE, IPROC
732 INTEGER IPOSROOT, JPOSROOT, ILOCROOT, JLOCROOT
733 INTEGER(8) :: IA8, IS18, IIW8, IS8, IAS8
734 INTEGER ISHIFT, IARR, JARR
735 INTEGER TAILLE
736 LOGICAL :: EARLYT3ROOTINS
737 REAL VAL
738 earlyt3rootins = keep(200) .EQ.0
739 & .OR. ( keep(200) .LT. 0 .AND. keep(400) .EQ. 0 )
740 nb_rec = bufi( 1 )
741 IF ( nb_rec .LE. 0 ) THEN
742 end_msg_2_recv = end_msg_2_recv - 1
743 nb_rec = - nb_rec
744 END IF
745 IF ( nb_rec .eq. 0 ) GOTO 100
746 DO irec = 1, nb_rec
747 iarr = bufi( irec * 2 )
748 jarr = bufi( irec * 2 + 1 )
749 val = bufr( irec )
750 node_type = mumps_typenode(
751 & procnode_steps(abs(step(abs( iarr )))),
752 & keep(199) )
753 IF ( node_type .eq. 3 .AND. earlyt3rootins ) THEN
754 IF ( iarr .GT. 0 ) THEN
755 iposroot = root%RG2L_ROW( iarr )
756 jposroot = root%RG2L_COL( jarr )
757 ELSE
758 iposroot = root%RG2L_ROW( jarr )
759 jposroot = root%RG2L_COL( -iarr )
760 END IF
761 ilocroot = root%MBLOCK * ( ( iposroot - 1 ) /
762 & ( root%MBLOCK * root%NPROW ) )
763 & + mod( iposroot - 1, root%MBLOCK ) + 1
764 jlocroot = root%NBLOCK * ( ( jposroot - 1 ) /
765 & ( root%NBLOCK * root%NPCOL ) )
766 & + mod( jposroot - 1, root%NBLOCK ) + 1
767 IF (keep(60)==0) THEN
768 a( ptr_root + int(jlocroot-1,8) * int(local_m,8)
769 & + int(ilocroot-1,8)) = a( ptr_root
770 & + int(jlocroot - 1,8) * int(local_m,8)
771 & + int(ilocroot - 1,8) )
772 & + val
773 ELSE
774 root%SCHUR_POINTER( int(jlocroot-1,8)
775 & * int(root%SCHUR_LLD,8)
776 & + int(ilocroot,8) )
777 & = root%SCHUR_POINTER( int(jlocroot - 1,8)
778 & * int(root%SCHUR_LLD,8)
779 & + int(ilocroot,8))
780 & + val
781 ENDIF
782 ELSE IF (iarr.GE.0) THEN
783 IF (iarr.EQ.jarr) THEN
784 ia8 = ptrarw(iarr)
785 dblarr(ia8) = dblarr(ia8) + val
786 ELSE
787 is18 = ptraiw(iarr)
788 ishift = intarr(is18) + iw4(iarr,2)
789 iw4(iarr,2) = iw4(iarr,2) - 1
790 iiw8 = is18 + ishift + 2
791 intarr(iiw8) = jarr
792 is8 = ptrarw(iarr)
793 ias8 = is8 + ishift
794 dblarr(ias8) = val
795 ENDIF
796 ELSE
797 iarr = -iarr
798 is8 = ptraiw(iarr)+iw4(iarr,1)+2
799 intarr(is8) = jarr
800 ias8 = ptrarw(iarr)+iw4(iarr,1)
801 iw4(iarr,1) = iw4(iarr,1) - 1
802 dblarr(ias8) = val
803 IF ( iw4(iarr,1) .EQ. 0
804 & .AND. step(iarr) > 0 ) THEN
805 iproc = mumps_procnode( procnode_steps(step(iarr)),
806 & keep(199) )
807 IF ( iproc .EQ. myid ) THEN
808 taille = intarr( ptraiw(iarr) )
809 CALL smumps_quick_sort_arrowheads( n, perm,
810 & intarr( ptraiw(iarr) + 3 ),
811 & dblarr( ptrarw(iarr) + 1 ),
812 & taille, 1, taille )
813 ENDIF
814 END IF
815 ENDIF
816 ENDDO
817 100 CONTINUE
818 RETURN
819 END SUBROUTINE smumps_dist_treat_recv_buf
subroutine mumps_propinfo(icntl, info, comm, id)
subroutine mpi_recv(buf, cnt, datatype, source, tag, comm, status, ierr)
Definition mpi.f:461
subroutine mpi_isend(buf, cnt, datatype, dest, tag, comm, ireq, ierr)
Definition mpi.f:382
subroutine mpi_test(ireq, flag, status, ierr)
Definition mpi.f:502
subroutine mpi_iprobe(source, tag, comm, flag, status, ierr)
Definition mpi.f:360
subroutine mpi_wait(ireq, status, ierr)
Definition mpi.f:525
recursive subroutine smumps_quick_sort_arrowheads(n, perm, intlist, dbllist, taille, lo, hi)
subroutine smumps_get_root_info(root, local_m, local_n, ptr_root, la)
subroutine smumps_set_root_to_zero(root, keep, a, la)
subroutine smumps_redistribution(n, nz_loc8, id, dblarr, ldblarr, intarr, lintarr, ptraiw, ptrarw, keep, keep8, myid, comm, nbrecords a, la, root, procnode_steps, slavef, perm, step, icntl, info, nsend8, nlocal8, istep_to_iniv2, candidates)
subroutine smumps_dist_fill_buffer()
subroutine smumps_dist_treat_recv_buf(bufi, bufr, nbrecords, n, iw4, keep, keep8, local_m, local_n, root, ptr_root, a, la, end_msg_2_recv, myid, procnode_steps, slavef, ptraiw, ptrarw, perm, step, intarr, lintarr, dblarr, ldblarr)
subroutine smumps_build_mapping(n, mapping, nnz, irn, jcn, procnode, step, slavef, perm, fils, rg2l, keep, keep8, mblock, nblock, nprow, npcol)
integer function mumps_typesplit(procinfo_inode, k199)
subroutine mumps_typeandprocnode(tpn, mumps_procnode, procinfo_inode, k199)