19 IF (
associated(lmat%COL))
THEN
21 IF (
associated(lmat%COL(j)%IRN))
THEN
22 DEALLOCATE(lmat%COL(j)%IRN)
23 NULLIFY(lmat%COL(j)%IRN)
35 IF (
associated(gcomp%IPE))
THEN
39 IF (
associated(gcomp%ADJ))
THEN
46 & NBLK, NDOF, BLKPTR, BLKVAR,
47 & SIZEOFBLOCKS, DOF2BLOCK )
49 INTEGER,
INTENT(IN) :: NBLK, NDOF
50 INTEGER,
INTENT(IN) :: BLKPTR(NBLK+1), BLKVAR(NDOF)
51 INTEGER,
INTENT(OUT):: SIZEOFBLOCKS(NBLK), DOF2BLOCK(NDOF)
52 INTEGER :: IB, I, IDOF
54 sizeofblocks(ib)= blkptr(ib+1)-blkptr(ib)
55 DO i=blkptr(ib), blkptr(ib+1)-1
63 & NBLK, NDOF, NNZ, IRN, JCN,
65 & IFLAG, IERROR, LP, LPOK,
69 INTEGER,
INTENT(IN) :: MYID, NBLK, NDOF
70 INTEGER(8),
INTENT(IN) :: NNZ
71 INTEGER,
INTENT(IN) :: IRN(max(1_8,NNZ)), JCN(max(1_8,NNZ))
72 INTEGER,
INTENT(IN) :: DOF2BLOCK(NDOF)
73 INTEGER :: LP, IFLAG, IERROR
74 LOGICAL,
INTENT(IN) :: LPOK
76 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: FLAG
78 INTEGER :: I, J, JJB, IIB, IB, JB, NB, PT
82 ALLOCATE(lmat%COL(nblk),flag(nblk), stat=allocok)
83 IF (allocok.NE.0)
THEN
87 WRITE(lp, *)
" ERROR allocate of LMAT%COL"
92 lmat%COL(ib)%NBINCOL = 0
99 IF ( (i.GT.ndof).OR.(j.GT.ndof).OR.(i.LT.1)
107 lmat%NZL = lmat%NZL+1_8
108 lmat%COL(jjb)%NBINCOL = lmat%COL(jjb)%NBINCOL + 1
112 IF (ierror.GE.1)
THEN
113 IF (mod(iflag,2) .EQ. 0) iflag = iflag+1
116 nb = lmat%COL(jb)%NBINCOL
118 ALLOCATE(lmat%COL(jb)%IRN(nb), stat=allocok)
119 IF (allocok.NE.0)
THEN
123 WRITE(lp, *)
" ERROR allocate of LMAT%COL"
132 IF ( (i.LE.ndof).AND.(j.LE.ndof).AND.(i.GE.1)
133 & .AND.(j.GE.1))
THEN
141 lmat%COL(jjb)%IRN(pt) = iib
146 & nblk, lmat, flag(1), iflag, ierror, lp, lpok
152 & NBLK, LMAT, FLAG, IFLAG, IERROR, LP, LPOK
156 INTEGER,
INTENT(IN) :: MYID, NBLK, LP
157 LOGICAL,
INTENT(IN) :: LPOK
158 INTEGER,
INTENT(OUT) :: FLAG(NBLK)
159 INTEGER,
INTENT(INOUT) :: IFLAG, IERROR
160 TYPE(LMATRIX_T),
INTENT(INOUT) :: LMAT
161 INTEGER,
POINTER,
DIMENSION(:) :: PTCLEAN
162 INTEGER :: allocok, IB, JB, NB
168 IF ( lmat%COL(jb)%NBINCOL.EQ.0) cycle
170 DO ib=1, lmat%COL(jb)%NBINCOL
171 IF (flag(lmat%COL(jb)%IRN(ib)).EQ.jb)
THEN
172 lmat%COL(jb)%IRN(ib)=0
175 lmat%NZL = lmat%NZL+1_8
176 flag(lmat%COL(jb)%IRN(ib)) = jb
180 ALLOCATE(ptclean(nb), stat=allocok)
181 IF (allocok.NE.0)
THEN
185 WRITE(lp, *)
" ERROR allocate PTCLEAN of size",
191 DO ib=1, lmat%COL(jb)%NBINCOL
192 IF (lmat%COL(jb)%IRN(ib).NE.0)
THEN
194 ptclean(nb)=lmat%COL(jb)%IRN(ib)
197 lmat%COL(jb)%NBINCOL = nb
198 deallocate(lmat%COL(jb)%IRN)
199 lmat%COL(jb)%IRN => ptclean
202 deallocate(lmat%COL(jb)%IRN)
203 NULLIFY(lmat%COL(jb)%IRN)
209 & LMAT, LUMAT, INFO, ICNTL )
213 INTEGER,
INTENT(IN) :: ICNTL(60)
214 INTEGER,
INTENT(INOUT) :: INFO(80)
215 INTEGER :: IB, IIB, JB, allocok, LP, MPG, NB, IERR
218 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
220 prokg = ( mpg .GT. 0 .and. (icntl(4).GE.2) )
221 lumat%NBCOL = lmat%NBCOL
222 lumat%NZL = 2_8*lmat%NZL
223 ALLOCATE( lumat%COL(lmat%NBCOL),stat=allocok)
224 IF (allocok.NE.0)
THEN
226 info( 2 ) = lmat%NBCOL
228 WRITE(lp, *)
" ERROR allocating LUMAT%COL "
233 lumat%COL(jb)%NBINCOL = lmat%COL(jb)%NBINCOL
236 DO ib=1, lmat%COL(jb)%NBINCOL
237 iib=lmat%COL(jb)%IRN(ib)
238 lumat%COL(iib)%NBINCOL = lumat%COL(iib)%NBINCOL + 1
242 nb = lumat%COL(jb)%NBINCOL
243 ALLOCATE(lumat%COL(jb)%IRN(nb), stat=ierr)
248 WRITE(lp, *)
" ERROR allocating columns of LUMAT"
254 lumat%COL(jb)%NBINCOL = 0
257 DO ib=1, lmat%COL(jb)%NBINCOL
258 iib=lmat%COL(jb)%IRN(ib)
259 nb = lumat%COL(jb)%NBINCOL+1
260 lumat%COL(jb)%NBINCOL = nb
261 lumat%COL(jb)%IRN(nb)= iib
262 nb = lumat%COL(iib)%NBINCOL+1
263 lumat%COL(iib)%NBINCOL = nb
264 lumat%COL(iib)%IRN(nb)= jb
273 INTEGER,
INTENT(IN) :: MYID, LP
275 write(lp,*) myid,
" ... LMATRIX %NBCOL, %NZL= ",
276 & lmat%NBCOL, lmat%NZL
277 IF (lmat%NBCOL.GE.0.AND.
associated(lmat%COL))
THEN
279 IF (lmat%COL(jb)%NBINCOL.GT.0)
THEN
280 WRITE(lp,*) myid,
" ... Column=", jb ,
" nb entries =",
281 & lmat%COL(jb)%NBINCOL,
" List of entries:",
282 & lmat%COL(jb)%IRN(1:lmat%COL(jb)%NBINCOL)
290 & LMAT, GCOMP, INFO, ICNTL )
293 INTEGER,
INTENT(IN) :: MYID
294 LOGICAL,
INTENT(IN) :: UNFOLD, READY_FOR_ANA_F
297 INTEGER,
INTENT(IN) :: ICNTL(60)
298 INTEGER,
INTENT(INOUT) :: INFO(80)
299 INTEGER :: IB, IIB, JJB, allocok, LP, MPG
300 INTEGER(8) :: JPOS, SIZEGCOMPALLOCATED
301 INTEGER(8),
ALLOCATABLE,
DIMENSION(:) :: IQ
302#if defined(DETERMINISTIC_PARALLEL_GRAPH)
303 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: WORK
304 INTEGER(8) :: IFIRST, ILAST
309 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
311 prokg = ( mpg .GT. 0 .and. (icntl(4).GE.2) )
312 gcomp%NG = lmat%NBCOL
314 gcomp%NZG = 2_8*lmat%NZL
315 sizegcompallocated = gcomp%NZG + int(gcomp%NG,8)+1_8
316 ELSE IF (ready_for_ana_f)
THEN
318 sizegcompallocated = gcomp%NZG + int(gcomp%NG,8)+1_8
321 sizegcompallocated = gcomp%NZG
323 gcomp%SIZEADJALLOCATED= sizegcompallocated
324 ALLOCATE( gcomp%ADJ(sizegcompallocated),
325 & gcomp%IPE(gcomp%NG+1),
326 & iq(gcomp%NG),stat=allocok)
327 IF (allocok.NE.0)
THEN
330 & gcomp%NZG + 3_8*int(gcomp%NG,8)+1_8, info(2))
332 WRITE(lp, *)
" ERROR allocating graph in",
333 &
" MUMPS_AB_LMAT_TO_CLEAN_G"
342 DO ib=1, lmat%COL(jjb)%NBINCOL
343 iib=lmat%COL(jjb)%IRN(ib)
350 iq(jjb) = lmat%COL(jjb)%NBINCOL
355 gcomp%IPE(jjb+1) = gcomp%IPE(jjb)+iq(jjb)
359 iq(jjb)= gcomp%IPE(jjb)
362 DO ib=1, lmat%COL(jjb)%NBINCOL
363 iib=lmat%COL(jjb)%IRN(ib)
364 gcomp%ADJ(iq(iib))= jjb
365 iq(iib) = iq(iib)+1_8
366 gcomp%ADJ(iq(jjb))= iib
367 iq(jjb) = iq(jjb)+1_8
372 jpos = gcomp%IPE(jjb)
373 DO ib=1, lmat%COL(jjb)%NBINCOL
374 iib=lmat%COL(jjb)%IRN(ib)
381#if defined(DETERMINISTIC_PARALLEL_GRAPH)
382 IF (.NOT.ready_for_ana_f)
THEN
383 ALLOCATE(work(0:gcomp%NG),stat=allocok)
384 IF (allocok.NE.0)
THEN
388 WRITE(lp, *)
" ERROR allocating WORK in",
389 &
" MUMPS_AB_LMAT_TO_CLEAN_G"
394 ifirst = gcomp%IPE(jjb)
395 ilast= gcomp%IPE(jjb+1)-1
396 l = int(ilast-ifirst+1)
398 IF (l .GE. gcomp%NG )
THEN
399 WRITE(*,*)
" Internal error in MUMPS_AB_LMAT_TO_CLEAN_G",
403 CALL mumps_mergesort( l,
404 & gcomp%ADJ(ifirst:ilast), work(0:l+1) )
405 CALL mumps_mergeswap1( l,
406 & work(0:l+1), gcomp%ADJ(ifirst:ilast) )
413#if defined(DETERMINISTIC_PARALLEL_GRAPH)
415 SUBROUTINE mumps_mergesort(N, K, L)
417 INTEGER :: K(:), L(0:)
418 INTEGER :: P, Q, S, T
423 IF (k(p) <= k(p+1))
THEN
432 IF (l(n+1) == 0)
THEN
435 l(n+1) = iabs(l(n+1))
444 IF(k(p) .GT. k(q))
GOTO 600
449 IF (p .GT. 0)
GOTO 300
463 IF (q .GT. 0)
GOTO 300
481 END SUBROUTINE mumps_mergesort
482 SUBROUTINE mumps_mergeswap1(N, L, A)
483 INTEGER :: I, LP, ISWAP, N
484 INTEGER :: L(0:), A(:)
488 IF ((lp==0).OR.(i>n))
EXIT
502 END SUBROUTINE mumps_mergeswap1
506 & INFO, ICNTL, COMM, NBLK, MYID, NPROCS,
511 include
'mumps_tags.h'
513 INTEGER,
INTENT(IN) :: OPTION, NBLK
514 INTEGER,
INTENT(IN) :: ICNTL(60), COMM, MYID, NPROCS
516 TYPE(LMATRIX_T) :: LMAT
517 INTEGER,
INTENT(OUT):: MAPCOL(NBLK)
518 INTEGER :: LP, SIZE_NZROW, I
520 INTEGER(8) :: NZL, NNZ
521 INTEGER,
DIMENSION(:),
ALLOCATABLE :: NZ_ROW
523 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
524 IF (option.EQ.1)
THEN
531 ALLOCATE(nz_row(nblk), stat=ierr)
537 &
" ERROR allocate in MUMPS_AB_COL_DISTRIBUTION ", info(2)
542 IF (info(1).LT.0)
GOTO 500
543 IF (option.NE.1)
THEN
545 mapcol(i) = lmat%COL(i)%NBINCOL
548 & mpi_integer, mpi_sum, comm, ierr)
550 & mpi_integer8, mpi_sum, comm, ierr)
553 & nnz, nz_row(1), size_nzrow, nblk, nprocs, mapcol(1))
555 IF (
allocated(nz_row))
DEALLOCATE(nz_row)
559 & MYID, NNZ, NZ_ROW, SIZE_NZROW, NBLK, NPROCS, MAPCOL )
560 INTEGER,
INTENT(IN) :: OPTION, MYID, SIZE_NZROW, NBLK
561 INTEGER,
INTENT(IN) :: ICNTL(60), NPROCS
564 INTEGER,
INTENT(IN) :: NZ_ROW(SIZE_NZROW)
565 INTEGER,
INTENT(OUT):: MAPCOL(NBLK)
566 INTEGER :: I, J, P, F, LP, IERR
568 INTEGER(8) :: SHARE, T
569 INTEGER,
DIMENSION(:),
ALLOCATABLE :: FIRST
571 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
572 ALLOCATE(first(nprocs+1), stat=ierr)
578 &
" ERROR allocate in MUMPS_AB_COL_DISTRIBUTION ", info(2)
585 IF (option.EQ.1)
THEN
586 share = int(nblk/nprocs,8)
588 first(i) = (i-1)*int(share)+1
590 first(nprocs+1)=nblk+1
592 share = (nnz-1_8)/int(nprocs,8) + 1_8
597 t = t+int(nz_row(i),8)
599 & (t .GE. share) .OR.
600 & ((nblk-i).EQ.(nprocs-p-1)) .OR.
613 IF ((i.EQ.nblk).AND.(p.NE.nprocs))
THEN
619 first(nprocs+1) = nblk+1
622 DO j=first(i), first(i+1)-1
626 IF (
allocated(first))
DEALLOCATE(first)
631 & MAPCOLonLUMAT, MAPCOL_IN_NSTEPS,
632 & INFO, ICNTL, KEEP, COMM, MYID, NBLK, NPROCS,
633 & LMAT, MAPCOL, SIZEMAPCOL,
639 include
'mumps_tags.h'
640 LOGICAL,
INTENT(IN) :: MAPCOLonLUMAT, MAPCOL_IN_NSTEPS
641 INTEGER,
INTENT(IN) :: MYID, NPROCS, NBLK, SIZEMAPCOL
642 INTEGER,
INTENT(IN) :: ICNTL(60), COMM, KEEP(500)
643 INTEGER,
INTENT(IN) :: SIZESTEP
644 INTEGER,
INTENT(IN) :: STEP(SIZESTEP)
645 INTEGER,
INTENT(INOUT) :: INFO(80)
647 INTEGER,
INTENT(INOUT) :: MAPCOL(SIZEMAPCOL)
649 INTEGER :: NBLKloc, IERR, JB, IB, LP, NB, I,
651 INTEGER(8) :: NNZ, NZ_locMAX8, NSEND8, NLOCAL8
653 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: WT, WNBINCOL
657 IF (nblkloc.NE.nblk)
THEN
658 write(6,*)
"Internal error in MUMPS_AB_BUILD_DCLEAN_LUMATRIX ",
659 &
"NBLKloc, NBLK=", nblkloc, nblk
662 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
663 ALLOCATE(wt(nblk), wnbincol(nblk), stat=ierr)
668 WRITE(lp, *)
" ERROR allocate of LUMAT%COL; WT"
674 IF ( info(1) .LT. 0 )
GOTO 500
676 wt(jb) = lmat%COL(jb)%NBINCOL
679 IF ( lmat%COL(jb)%NBINCOL.EQ.0) cycle
680 DO ib=1, lmat%COL(jb)%NBINCOL
681 i = lmat%COL(jb)%IRN(ib)
686 & mpi_integer, mpi_sum, comm, ierr)
687 IF (
allocated(wt))
DEALLOCATE(wt)
688 IF (mapcolonlumat)
THEN
691 nnz=nnz+int(wnbincol(i),8)
694 & myid, nnz, wnbincol(1), nblk,
695 & nblk, nprocs, mapcol(1))
698 IF ( info(1) .LT. 0 )
GOTO 500
702 ALLOCATE(lumat%COL(nblk), stat=ierr)
707 WRITE(lp, *)
" ERROR allocate of LUMAT%COL; WT"
710 IF ( info(1) .GE. 0 )
THEN
713 IF (mapcol_in_nsteps)
THEN
714 IF (mapcol(abs(step(jb))).EQ.myid)
THEN
715 lumat%NZL = lumat%NZL + int(nb,8)
720 IF (mapcol(jb).EQ.myid)
THEN
721 lumat%NZL = lumat%NZL + int(nb,8)
726 lumat%COL(jb)%NBINCOL = nb
728 ALLOCATE(lumat%COL(jb)%IRN(nb), stat=ierr)
733 WRITE(lp, *)
" ERROR allocate of LUMAT%COL"
742 IF ( info(1) .LT. 0 )
GOTO 500
743 IF (
allocated(wnbincol))
DEALLOCATE(wnbincol)
745 & mpi_max, comm, ierr)
747 IF (nz_locmax8 .LT. int(nbrecords,8))
THEN
748 nbrecords = int(nz_locmax8)
753 & info, icntl, comm, myid, nblk, nprocs,
754 & lmat, mapcol, sizemapcol, step, sizestep,
755 & lumat, nbrecords, nsend8, nlocal8
760 IF ( info(1) .LT. 0 )
GOTO 500
761 ALLOCATE(wt(nblk), stat=ierr)
766 WRITE(lp, *)
" ERROR allocate of LUMAT%COL; WT"
771 & nblk, lumat, wt(1), info(1), info(2), lp, lpok
775 IF ( info(1) .LT. 0 )
GOTO 500
779 IF (
allocated(wt))
DEALLOCATE(wt)
780 IF (
allocated(wnbincol))
DEALLOCATE(wnbincol)
785 & INFO, ICNTL, KEEP, COMM, MYID, NBLK,
786 & LUMAT, PROCNODE_STEPS, NSTEPS, MAPCOL,
787 & LUMAT_REMAP, NBRECORDS, STEP
792 include
'mumps_tags.h'
793 INTEGER :: IERR, MASTER
795 integer,
INTENT(IN) :: myid, nblk, nsteps, keep(500)
796 INTEGER,
INTENT(IN) :: ICNTL(60), COMM
798 INTEGER,
INTENT(IN) :: PROCNODE_STEPS(NSTEPS)
800 INTEGER,
INTENT(IN) :: STEP(NBLK)
801 TYPE(
lmatrix_t),
INTENT(INOUT) :: LUMAT_REMAP
802 INTEGER,
INTENT(OUT) :: NBRECORDS
803 INTEGER,
INTENT(OUT) :: MAPCOL(NSTEPS)
804 INTEGER :: LP, MP, ISTEP, JB, NB
806 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: WT, WNBINCOL
807 INTEGER MUMPS_PROCNODE
808 INTEGER(8) :: NZ_locMAX8
811 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
812 ALLOCATE(wt(nblk), wnbincol(nblk), stat=ierr)
817 WRITE(lp, *)
" ERROR allocate WT"
822 IF ( info(1) .LT. 0 )
GOTO 500
824 wt(jb) = lumat%COL(jb)%NBINCOL
827 & mpi_integer, mpi_sum, comm, ierr)
828 IF (
allocated(wt))
DEALLOCATE(wt)
829 IF (myid.EQ.master)
THEN
832 & mumps_procnode(procnode_steps(istep),keep(199))
835 CALL mpi_bcast( mapcol(1), nsteps, mpi_integer,
836 & master, comm, ierr )
837 CALL mpi_bcast( step(1), nblk, mpi_integer,
838 & master, comm, ierr )
839 lumat_remap%NBCOL = nblk
840 ALLOCATE(lumat_remap%COL(nblk), stat=ierr)
845 WRITE(lp, *)
" ERROR allocate of LUMAT_REMAP%COL"
848 IF ( info(1) .GE. 0 )
THEN
849 lumat_remap%NZL = 0_8
852 IF (mapcol(abs(step(jb))).EQ.myid)
THEN
853 lumat_remap%NZL = lumat_remap%NZL + int(nb,8)
857 lumat_remap%COL(jb)%NBINCOL = nb
859 ALLOCATE(lumat_remap%COL(jb)%IRN(nb), stat=ierr)
864 WRITE(lp, *)
" ERROR allocate of LUMAT_REMAP%COL"
873 IF ( info(1) .LT. 0 )
GOTO 500
874 IF (
allocated(wnbincol))
DEALLOCATE(wnbincol)
875 CALL mpi_allreduce(lumat_remap%NZL, nz_locmax8, 1, mpi_integer8,
876 & mpi_max, comm, ierr)
878 IF (nz_locmax8 .LT. int(nbrecords,8))
THEN
879 nbrecords = int(nz_locmax8)
883 IF (
allocated(wt))
DEALLOCATE(wt)
884 IF (
allocated(wnbincol))
DEALLOCATE(wnbincol)
889 & MYID, NPROCS, COMM,
891 & IRN, JCN, DOF2BLOCK,
893 & LUMAT, GCOMP, READY_FOR_ANA_F)
897 include
'mumps_tags.h'
899 parameter( master = 0 )
900 INTEGER,
INTENT(IN) :: MYID, NPROCS, NBLK, NDOF
901 INTEGER(8),
INTENT(IN) :: NNZ
902 INTEGER,
INTENT(IN) :: IRN(max(1_8,NNZ)), JCN(max(1_8,NNZ))
903 LOGICAL,
INTENT(IN) :: READY_FOR_ANA_F
904 INTEGER,
INTENT(INOUT) :: DOF2BLOCK(NDOF)
905 INTEGER,
INTENT(IN) :: ICNTL(60), COMM
906 INTEGER,
INTENT(INOUT) :: KEEP(500), INFO(80)
910 INTEGER :: IDUMMY_ARRAY(1)
911 INTEGER :: allocok, LP, MPG
912 LOGICAL :: LPOK, PROKG
913 INTEGER,
DIMENSION(:),
ALLOCATABLE :: MAPCOL
914 LOGICAL :: MAPCOLonLUMAT, MAPCOL_IN_NSTEPS
918 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
920 prokg = ( mpg .GT. 0 .and. myid .eq. master )
921 mapcolonlumat = .false.
922 mapcol_in_nsteps = .false.
923 IF (keep(14).EQ.1)
THEN
926 IF (keep(14).EQ.0)
THEN
927 CALL mpi_bcast( dof2block, ndof, mpi_integer, master,
931 & nblk, ndof, nnz, irn, jcn,
933 & info(1), info(2), lp, lpok,
937 IF ( info(1) .LT. 0 )
GOTO 500
938 ALLOCATE(mapcol(nblk), stat=allocok)
939 IF (allocok.NE.0)
THEN
943 WRITE(lp, *)
" ERROR allocate MAPCOL of size",
949 IF ( info(1) .LT. 0 )
GOTO 500
950 IF (.NOT.mapcolonlumat)
THEN
952 & info, icntl, comm, nblk, myid, nprocs,
956 IF ( info(1) .LT. 0 )
GOTO 500
959 & mapcolonlumat, mapcol_in_nsteps,
960 & info, icntl, keep, comm, myid, nblk, nprocs,
961 & lmat, mapcol, nblk,
966 IF ( info(1) .LT. 0 )
GOTO 500
967 IF (
allocated(mapcol))
DEALLOCATE(mapcol)
970 & lumat, gcomp, info, icntl
974 IF ( info(1) .LT. 0 )
GOTO 500
975 IF (keep(494).EQ.0)
THEN
980 IF (
allocated(mapcol))
DEALLOCATE(mapcol)
987 & MYID, NPROCS, COMM,
990 & PROCNODE_STEPS, NSTEPS, STEP,
996 include
'mumps_tags.h'
998 parameter( master = 0 )
999 INTEGER,
INTENT(IN) :: MYID, NPROCS, NBLK, NDOF, NSTEPS
1000 INTEGER(8),
INTENT(IN) :: NNZ
1001 INTEGER,
INTENT(IN) :: IRN(max(1_8,NNZ)), JCN(max(1_8,NNZ))
1002 INTEGER,
INTENT(IN) :: ICNTL(60), COMM
1003 INTEGER,
INTENT(IN) :: PROCNODE_STEPS(NSTEPS)
1004 INTEGER,
INTENT(IN) :: STEP(NBLK)
1005 INTEGER,
INTENT(INOUT) :: KEEP(500), INFO(80)
1006 INTEGER,
INTENT(OUT) :: MAPCOL(NSTEPS)
1007 TYPE(LMATRIX_T) :: LUMAT
1008 INTEGER,
DIMENSION(:),
allocatable:: DOF2BLOCK
1010 INTEGER :: allocok, LP
1012 INTEGER :: IDOF, ISTEP
1013 LOGICAL :: MAPCOL_IN_NSTEPS, MAPCOLonLUMAT
1016 INTEGER MUMPS_PROCNODE
1018 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
1019 mapcolonlumat = .false.
1020 mapcol_in_nsteps = .true.
1021 IF (keep(14).EQ.1)
THEN
1024 allocate(dof2block(ndof), stat=allocok)
1025 IF (allocok.NE.0)
THEN
1028 IF ( lpok )
WRITE(lp, 150)
' DOF2BLOCK'
1032 IF ( info(1) .LT. 0 )
GOTO 500
1034 dof2block(idof) = idof
1037 & nblk, ndof, nnz, irn, jcn,
1039 & info(1), info(2), lp, lpok,
1043 IF ( info(1) .LT. 0 )
GOTO 500
1044 IF (
allocated(dof2block))
DEALLOCATE(dof2block)
1045 IF (myid.EQ.master)
THEN
1048 & mumps_procnode(procnode_steps(istep),keep(199))
1051 CALL mpi_bcast( mapcol(1), nsteps, mpi_integer,
1052 & master, comm, ierr )
1053 CALL mpi_bcast( step(1), nblk, mpi_integer,
1054 & master, comm, ierr )
1056 & mapcolonlumat, mapcol_in_nsteps,
1057 & info, icntl, keep, comm, myid, nblk, nprocs,
1058 & lmat, mapcol, nsteps,
1059 & step, nblk, lumat)
1062 IF ( info(1) .LT. 0 )
GOTO 500
1065 IF (
allocated(dof2block))
DEALLOCATE(dof2block)
1071 & /
' ** FAILURE IN MUMPS_AB_DCOORD_TO_DTREE_LUMAT, ',
1072 &
' DYNAMIC ALLOCATION OF ',
1078 & INFO, ICNTL, COMM, MYID, NBLK, SLAVEF,
1079 & LMAT, MAPCOL, SIZEMAPCOL, STEP, SIZESTEP,
1080 & LUMAT, NBRECORDS, NSEND8, NLOCAL8
1085 include
'mumps_tags.h'
1086 INTEGER :: IERR, MASTER, MSGSOU
1088 INTEGER :: STATUS(MPI_STATUS_SIZE)
1089 LOGICAL,
INTENT(IN) :: UNFOLD, MAPCOL_IN_NSTEPS
1090 INTEGER,
INTENT(IN) :: MYID, SLAVEF, NBLK
1091 INTEGER,
INTENT(IN) :: SIZEMAPCOL, SIZESTEP
1092 INTEGER,
INTENT(IN) :: ICNTL(60), COMM, NBRECORDS
1094 TYPE(LMATRIX_T),
INTENT(IN) :: LMAT
1095 INTEGER,
INTENT(IN) :: MAPCOL(SIZEMAPCOL)
1096 INTEGER,
INTENT(IN) :: STEP(SIZESTEP)
1097 TYPE(LMATRIX_T),
INTENT(INOUT) :: LUMAT
1098 INTEGER(8),
INTENT(OUT) :: NSEND8, NLOCAL8
1099 INTEGER :: LP, MP, allocok
1100 INTEGER :: IB, JB, I, II, ISEND, JSEND, ITOSEND
1103 INTEGER END_MSG_2_RECV
1104 INTEGER KPROBE, FREQPROBE
1105 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: PTLOC
1106 INTEGER,
ALLOCATABLE,
DIMENSION(:,:,:) :: BUFI
1107 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: BUFRECI
1108 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: IACT, IREQI
1109 LOGICAL,
ALLOCATABLE,
DIMENSION(:) :: SEND_ACTIVE
1114 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
1122 end_msg_2_recv = slavef-1
1123 ALLOCATE( iact(slavef), stat=allocok)
1124 IF ( allocok .GT. 0 )
THEN
1127 &
'** Error allocating IACT in matrix distribution'
1133 ALLOCATE( ireqi(slavef), stat=allocok)
1134 IF ( allocok .GT. 0 )
THEN
1137 &
'** Error allocating IREQI in matrix distribution'
1143 ALLOCATE( send_active(slavef), stat=allocok)
1144 IF ( allocok .GT. 0 )
THEN
1147 &
'** Error allocating SEND_ACTIVE in matrix distribution'
1153 ALLOCATE( bufi( nbrecords * 2 + 1, 2, slavef ), stat=allocok)
1154 IF ( allocok .GT. 0 )
THEN
1157 &
'** Error allocating int buffer for matrix distribution'
1160 info(2) = ( nbrecords * 2 + 1 ) * slavef * 2
1163 ALLOCATE( bufreci( nbrecords * 2 + 1 ), stat = allocok )
1164 IF ( allocok .GT. 0 )
THEN
1167 &
'** Error allocating int recv buffer for matrix distribution'
1170 info(2) = nbrecords * 2 + 1
1173 ALLOCATE( ptloc( nblk ), stat = allocok )
1174 IF ( allocok .GT. 0 )
THEN
1177 &
'** Error allocating int recv buffer for matrix distribution'
1185 IF ( info(1) .LT. 0 )
GOTO 100
1193 send_active( i ) = .false.
1200 freqprobe =
max(1,nbrecords/10)
1201 IF (slavef .EQ. 1) freqprobe = huge(freqprobe)
1203 IF ( lmat%COL(jb)%NBINCOL.EQ.0) cycle
1204 DO ii=1, lmat%COL(jb)%NBINCOL
1206 IF ( kprobe .eq. freqprobe )
THEN
1208 CALL mpi_iprobe( mpi_any_source, lmatdist, comm,
1209 & flag, status, ierr )
1211 msgsou = status( mpi_source )
1212 CALL mpi_recv( bufreci(1), nbrecords * 2 + 1,
1214 & msgsou, lmatdist, comm, status, ierr )
1216 & myid, bufreci(1), nbrecords, lumat,
1217 & nblk, ptloc(1), end_msg_2_recv
1221 ib = lmat%COL(jb)%IRN(ii)
1222 DO itosend=1,nbtosend
1223 IF (itosend.EQ.1)
THEN
1224 IF (mapcol_in_nsteps)
THEN
1225 dest = mapcol(abs(step(jb)))
1232 IF (mapcol_in_nsteps)
THEN
1233 dest = mapcol(abs(step(ib)))
1240 IF (dest.EQ.myid)
THEN
1241 lumat%COL(jsend)%IRN(1+ptloc(jsend))= isend
1242 ptloc(jsend) = ptloc(jsend) + 1
1243 nlocal8 = nlocal8 + 1_8
1245 nsend8 = nsend8 + 1_8
1247 & dest, isend, jsend, nblk,
1248 & bufi, bufreci, ptloc,
1249 & nbrecords, slavef, comm, myid, iact, ireqi,
1250 & send_active, lmat, lumat, end_msg_2_recv
1258 & nblk, bufi, bufreci, ptloc,
1259 & nbrecords, slavef, comm, myid, iact, ireqi,
1260 & send_active, lmat, lumat, end_msg_2_recv
1262 DO WHILE ( end_msg_2_recv .NE. 0 )
1263 CALL mpi_recv( bufreci(1), nbrecords * 2 + 1, mpi_integer,
1264 & mpi_any_source, lmatdist, comm, status, ierr )
1266 & myid, bufreci(1), nbrecords, lumat,
1267 & nblk, ptloc(1), end_msg_2_recv
1271 IF ( send_active( i ) )
THEN
1272 CALL mpi_wait( ireqi( i ), status, ierr )
1276 IF (
ALLOCATED(ptloc))
DEALLOCATE( ptloc )
1277 IF (
ALLOCATED(bufi))
DEALLOCATE( bufi )
1278 IF (
ALLOCATED(bufreci))
DEALLOCATE( bufreci )
1279 IF (
ALLOCATED(iact))
DEALLOCATE( iact )
1280 IF (
ALLOCATED(ireqi))
DEALLOCATE( ireqi )
1281 IF (
ALLOCATED(send_active))
DEALLOCATE( send_active )
1285 & MYID, BUFI, NBRECORDS, LUMAT,
1286 & NBLK, PTLOC, END_MSG_2_RECV
1291 include
'mumps_tags.h'
1292 INTEGER,
INTENT(IN) :: NBLK, MYID, NBRECORDS
1293 INTEGER,
INTENT(IN) :: BUFI( NBRECORDS * 2 + 1 )
1294 INTEGER,
INTENT(INOUT):: END_MSG_2_RECV, PTLOC(NBLK)
1296 INTEGER :: IREC, NB_REC, IB, JB
1298 IF ( nb_rec .LE. 0 )
THEN
1299 end_msg_2_recv = end_msg_2_recv - 1
1302 IF ( nb_rec .eq. 0 )
RETURN
1304 ib = bufi( irec * 2 )
1305 jb = bufi( irec * 2 + 1 )
1306 lumat%COL(jb)%IRN(1+ptloc(jb))= ib
1307 ptloc(jb) = ptloc(jb) + 1
1312 & DEST, ISEND, JSEND, NBLK,
1313 & BUFI, BUFRECI, PTLOC,
1314 & NBRECORDS, SLAVEF, COMM, MYID, IACT, IREQI,
1315 & SEND_ACTIVE, LMAT, LUMAT, END_MSG_2_RECV
1320 include
'mumps_tags.h'
1321 INTEGER :: STATUS(MPI_STATUS_SIZE)
1322 INTEGER,
INTENT(IN) :: DEST, ISEND, JSEND, SLAVEF, COMM, MYID,
1324 INTEGER,
INTENT(INOUT) :: END_MSG_2_RECV, PTLOC(NBLK)
1327 LOGICAL,
INTENT(INOUT) :: SEND_ACTIVE(SLAVEF)
1328 INTEGER,
INTENT(INOUT) :: IREQI(SLAVEF), IACT(SLAVEF)
1329 INTEGER,
INTENT(INOUT) :: BUFI( NBRECORDS * 2 + 1, 2, SLAVEF )
1330 INTEGER,
INTENT(INOUT) :: BUFRECI( NBRECORDS * 2 + 1)
1331 INTEGER :: IBEG, IEND, ISLAVE, TAILLE_SEND_I, IREQ, MSGSOU,
1334 IF ( DEST .eq. -3 ) THEN
1341 DO islave = ibeg, iend
1342 nbrec = bufi(1,iact(islave),islave)
1343 IF ( dest .eq. -3 )
THEN
1344 bufi(1,iact(islave),islave) = - nbrec
1346 IF ( dest .eq. -3 .or. nbrec + 1 > nbrecords )
THEN
1347 DO WHILE ( send_active( islave ) )
1348 CALL mpi_test( ireqi( islave ), flag, status, ierr )
1349 IF ( .NOT. flag )
THEN
1350 CALL mpi_iprobe( mpi_any_source, lmatdist, comm,
1351 & flag, status, ierr )
1353 msgsou = status(mpi_source)
1354 CALL mpi_recv( bufreci(1), 2*nbrecords+1,
1355 & mpi_integer, msgsou, lmatdist, comm,
1358 & myid, bufreci, nbrecords, lumat,
1359 & nblk, ptloc(1), end_msg_2_recv
1363 send_active( islave ) = .false.
1366 IF ( islave - 1 .ne. myid )
THEN
1367 taille_send_i = nbrec * 2 + 1
1368 CALL mpi_isend( bufi(1, iact(islave), islave ),
1370 & mpi_integer, islave - 1, lmatdist, comm,
1371 & ireqi( islave ), ierr )
1372 send_active( islave ) = .true.
1374 IF (nbrec.NE.0)
THEN
1375 write(*,*)
" Internal error in ",
1376 &
" MUMPS_AB_LMAT_FILL_BUFFER "
1380 iact( islave ) = 3 - iact( islave )
1381 bufi( 1, iact( islave ), islave ) = 0
1383 IF ( dest .ne. -3 )
THEN
1384 ireq = bufi(1,iact(islave),islave) + 1
1385 bufi(1,iact(islave),islave) = ireq
1386 bufi(ireq*2,iact(islave),islave) = isend
1387 bufi(ireq*2+1,iact(islave),islave) = jsend
1393 & ICNTL, KEEP, COMM, MYID, NPROCS, INFO,
1394 & GCOMP_DIST, GCOMP)
1398 include
'mumps_tags.h'
1399 INTEGER IERR, MASTER
1400 parameter( master = 0 )
1401 INTEGER :: STATUS(MPI_STATUS_SIZE)
1403 INTEGER,
INTENT(IN) :: MYID, NPROCS, ICNTL(60), COMM,
1405 INTEGER,
INTENT(INOUT) :: INFO(80)
1406 TYPE(COMPACT_GRAPH_T) :: GCOMP
1407 INTEGER :: NG, allocok, LP, MPG, I, J, K
1408 INTEGER :: INDX, NB_BLOCK_SENT, MAX_NBBLOCK_loc, NRECV,
1409 & BLOCKSIZE, SIZE_SENT, NB_BLOCKS, NBNONEMPTY,
1410 & FIRSTNONEMPTY, LASTNONEMPTY, NBBLOCK_loc
1411 INTEGER(4) :: IOVFLO
1412 INTEGER(8) :: NZG, NZG_CENT, I8, IBEG8, IEND8,
1413 & sizegcompallocated
1414 LOGICAL :: LPOK, PROKG
1415 INTEGER(8),
ALLOCATABLE,
DIMENSION(:) ::
1416 INTEGER,
ALLOCATABLE :: REQPTR(:)
1417 INTEGER(8),
ALLOCATABLE :: GPTR(:), GPTR_cp(:)
1419 lpok = ((lp.GT.0).AND.(icntl(4).GE.1))
1421 prokg = ( mpg .GT. 0 .and. myid .eq. master )
1422 prokg = (prokg.AND.(icntl(4).GE.2))
1423 iovflo = huge(iovflo)
1424 blocksize = int(
max(100000_8,int(iovflo,8)/200_8))
1425 nzg = gcomp_dist%NZG
1427 CALL mpi_reduce( nzg, nzg_cent, 1, mpi_integer8,
1428 & mpi_sum, master, comm, ierr )
1429 IF (myid.EQ.master)
THEN
1430 gcomp%NZG = nzg_cent
1432 sizegcompallocated = nzg_cent+int(ng,8)+1_8
1433 gcomp%SIZEADJALLOCATED = sizegcompallocated
1434 ALLOCATE( gcomp%ADJ(sizegcompallocated),
1437 & gptr_cp( nprocs ),
1438 & reqptr( nprocs-1 ),
1439 & iq(ng+1),stat=allocok)
1440 IF (allocok.NE.0)
THEN
1443 & nzg_cent + 3_8*int(ng,8)+3_8+3_8*int(nprocs,8)-1_8,
1446 &
WRITE(lp, *)
" ERROR allocating graph in",
1447 &
" MUMPS_AB_GATHER_GRAPH"
1450 ALLOCATE( iq(ng+1), stat=allocok)
1451 IF (allocok.NE.0)
THEN
1455 &
WRITE(lp, *)
" ERROR allocating pointers",
1456 &
" MUMPS_AB_GATHER_GRAPH"
1461 IF (info(1).LT.0)
GOTO 500
1465 iq(i) = int(gcomp_dist%IPE(i+1)-gcomp_dist%IPE(i))
1466 IF (iq(i).NE.0)
THEN
1467 IF (firstnonempty.EQ.0) firstnonempty=i
1471 nbnonempty = lastnonempty-firstnonempty+1
1472 IF (myid.EQ.master)
THEN
1477 IF (nbnonempty.GT.0)
THEN
1478 DO i=firstnonempty, lastnonempty
1479 gcomp%IPE(j) = iq(i)
1483 DO i = 1, nprocs - 1
1486 & gatherg_nb, comm, status, ierr )
1487 IF (nbnonempty.GT.0)
THEN
1490 & gatherg_first, comm, status, ierr )
1491 CALL mpi_recv( gcomp%IPE(j), nbnonempty,
1493 & gatherg_ipe, comm, status, ierr )
1497 CALL mpi_send( nbnonempty, 1, mpi_integer, master,
1498 & gatherg_nb, comm, ierr )
1499 IF (nbnonempty.GT.0)
THEN
1500 CALL mpi_send( firstnonempty, 1, mpi_integer, master,
1501 & gatherg_first, comm, ierr )
1502 CALL mpi_send( iq(firstnonempty), nbnonempty,
1503 & mpi_integer8, master,
1504 & gatherg_ipe, comm, ierr )
1507 IF (myid.EQ.master)
THEN
1510 iq(i+1) = iq(i) + gcomp%IPE(i)
1511 gcomp%IPE(i) = iq(i)
1513 gcomp%IPE(ng+1) = iq(ng+1)
1518 IF (myid.EQ.master)
THEN
1521 DO i = 1, nprocs - 1
1524 & gatherg_nzg, comm, status, ierr )
1525 nbblock_loc = ceiling(dble(gptr(i+1))/dble(blocksize))
1526 max_nbblock_loc =
max(max_nbblock_loc, nbblock_loc)
1527 nb_block_sent = nb_block_sent + nbblock_loc
1529 gptr( 1 ) = nzg + 1_8
1531 gptr( i ) = gptr( i ) + gptr( i-1 )
1534 CALL mpi_send( nzg, 1, mpi_integer8, master,
1535 & gatherg_nzg, comm, ierr )
1537 IF (myid.EQ.master)
THEN
1539 gptr_cp(i) = gptr(i)
1541 IF (nzg.GT.0_8)
THEN
1543 gcomp%ADJ(i8) = gcomp_dist%ADJ(i8)
1547 DO k = 1, max_nbblock_loc
1549 DO i = 1, nprocs - 1
1550 ibeg8 = gptr_cp( i )
1551 IF ( ibeg8 .LT. gptr(i+1))
THEN
1553 iend8 =
min(ibeg8+int(blocksize,8)-1_8,
1555 gptr_cp( i ) = iend8 + 1_8
1556 size_sent = int(iend8 - ibeg8 + 1_8)
1557 nb_blocks = nb_blocks + 1
1558 CALL mpi_irecv( gcomp%ADJ(ibeg8), size_sent,
1560 & i, gatherg_adj, comm, reqptr(i), ierr )
1562 reqptr( i ) = mpi_request_null
1567 & ( nprocs-1, reqptr, indx,
1571 DEALLOCATE( reqptr )
1573 DEALLOCATE( gptr_cp )
1575 IF (nzg.EQ.0)
GOTO 600
1576 DO i8=1_8, nzg, int(blocksize,8)
1577 size_sent = blocksize
1578 IF (nzg-i8+1_8.LT.int(blocksize,8))
THEN
1579 size_sent = int(nzg-i8+1_8)
1582 & gcomp_dist%ADJ(i8), size_sent,
1583 & mpi_integer, master,
1584 & gatherg_adj, comm, ierr )
1589 IF (myid.EQ.master)
THEN
1590 IF (
associated(gcomp%ADJ))
THEN
1591 DEALLOCATE(gcomp%ADJ)
1594 IF (
associated(gcomp%IPE))
THEN
1595 DEALLOCATE(gcomp%IPE)
1600 IF (
allocated(iq))
DEALLOCATE(iq)
subroutine mumps_ab_dcoord_to_dtree_lumat(myid, nprocs, comm, nblk, ndof, nnz, irn, jcn, procnode_steps, nsteps, step, icntl, info, keep, mapcol, lumat)
subroutine mumps_ab_print_lmatrix(lmat, myid, lp)
subroutine mumps_ab_lmat_to_clean_g(myid, unfold, ready_for_ana_f, lmat, gcomp, info, icntl)
subroutine mumps_ab_col_distribution(option, info, icntl, comm, nblk, myid, nprocs, lmat, mapcol)
subroutine mumps_ab_build_dclean_lumatrix(mapcolonlumat, mapcol_in_nsteps, info, icntl, keep, comm, myid, nblk, nprocs, lmat, mapcol, sizemapcol, step, sizestep, lumat)
subroutine mumps_ab_compute_sizeofblock(nblk, ndof, blkptr, blkvar, sizeofblocks, dof2block)
subroutine mumps_inialize_redist_lumat(info, icntl, keep, comm, myid, nblk, lumat, procnode_steps, nsteps, mapcol, lumat_remap, nbrecords, step)
subroutine mumps_ab_lmat_to_lumat(lmat, lumat, info, icntl)
subroutine mumps_ab_gather_graph(icntl, keep, comm, myid, nprocs, info, gcomp_dist, gcomp)
subroutine mumps_ab_lmat_treat_recv_buf(myid, bufi, nbrecords, lumat, nblk, ptloc, end_msg_2_recv)
subroutine mumps_ab_free_gcomp(gcomp)
subroutine mumps_ab_dist_lmat_to_lumat(unfold, mapcol_in_nsteps, info, icntl, comm, myid, nblk, slavef, lmat, mapcol, sizemapcol, step, sizestep, lumat, nbrecords, nsend8, nlocal8)
subroutine mumps_ab_localclean_lmat(myid, nblk, lmat, flag, iflag, ierror, lp, lpok)
subroutine mumps_ab_lmat_fill_buffer(dest, isend, jsend, nblk, bufi, bufreci, ptloc, nbrecords, slavef, comm, myid, iact, ireqi, send_active, lmat, lumat, end_msg_2_recv)
subroutine mumps_ab_coord_to_lmat(myid, nblk, ndof, nnz, irn, jcn, dof2block, iflag, ierror, lp, lpok, lmat)
subroutine mumps_ab_dcoord_to_dcompg(myid, nprocs, comm, nblk, ndof, nnz, irn, jcn, dof2block, icntl, info, keep, lumat, gcomp, ready_for_ana_f)
subroutine mumps_ab_free_lmat(lmat)
subroutine mumps_ab_compute_mapcol(option, info, icntl, myid, nnz, nz_row, size_nzrow, nblk, nprocs, mapcol)
subroutine mumps_propinfo(icntl, info, comm, id)
subroutine mpi_recv(buf, cnt, datatype, source, tag, comm, status, ierr)
subroutine mpi_isend(buf, cnt, datatype, dest, tag, comm, ireq, ierr)
subroutine mpi_test(ireq, flag, status, ierr)
subroutine mpi_iprobe(source, tag, comm, flag, status, ierr)
subroutine mpi_wait(ireq, status, ierr)
subroutine mpi_reduce(sendbuf, recvbuf, cnt, datatype, op, root, comm, ierr)
subroutine mpi_send(buf, cnt, datatype, dest, tag, comm, ierr)
subroutine mpi_allreduce(sendbuf, recvbuf, cnt, datatype, operation, comm, ierr)
subroutine mpi_waitany(cnt, array_of_requests, index, status, ierr)
subroutine mpi_bcast(buffer, cnt, datatype, root, comm, ierr)
subroutine mpi_irecv(buf, cnt, datatype, source, tag, comm, ireq, ierr)