<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=iso-8859-1">
<style type="text/css" style="display:none;"> P {margin-top:0;margin-bottom:0;} </style>
</head>
<body dir="ltr">
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
Dear All,</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
             As part of preparing a code to call the SLEPC eigenvalue solving library,</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
I am constructing a matrix in sparse CSR format row-by-row. Just for debugging <br>
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
purposes I write out the column values for a given row, which are stored in a <br>
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
PetscInt allocatable vector, using PetscIntView.</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<br>
</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
Everything works fine when the number of MPI processes exactly divide the</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
number of rows of the matrix, and so each process owns the same number of rows.</div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
However, when the <span style="font-size:12pt">number of MPI processes does not exactly divide the</span><br>
<span style="font-size:12pt">number of rows of the matrix, and so each process owns a different number of rows,</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt">the code hangs when it reaches the line that calls PetscIntView.</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt">To be precise the code hangs on the final row that a process, other than root,<b><i>
</i></b>owns<b><i>.</i></b></span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<b><i><span style="font-size:12pt">If I however comment out the call to <span style="font-size:12pt">
PetscIntView the code completes without error,</span></span></i></b></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<b><i><span style="font-size:12pt"><span style="font-size:12pt"> and produces the correct eigenvalues (hence we are not missing a row / miswriting a row).</span></span></i></b></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt"><span style="font-size:12pt">  </span> Note also that a simple direct writeout of this same array using a plain fortran command</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt">will write out the array without problem.<br>
</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt"><br>
</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt">I have attached below a small code that reproduces the problem.</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt">For this code we have nominally assigned 200 rows to our matrix. The code runs without</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt">problem using 1,2,4,5,8 or 10 MPI processes, all of which precisely divide 200,</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt"> but will hang for 3 MPI processes for example.</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt">For the case of 3 MPI processes the subroutine WHOSE_ROW_IS_IT allocates the rows
<br>
</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt">to each process as :</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt">  process no       first row           last row       no. of rows
<div>   0                             1                     66               66</div>
<div>   1                            67                   133             67</div>
<div>   2                          134                   200             67</div>
   <br>
</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt">The code will hang when process 1 calls <span style="font-size:12pt">
<span>PetscIntView</span></span> for its last row, row 133 for example.<br>
</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<span style="font-size:12pt"><br>
</span></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<b><span style="font-size:12pt"><i>One piece of additional information that may be relevant</i> is that the code
<i>does</i> run to completion <br>
</span></b></div>
<div style="font-family: Calibri, Helvetica, sans-serif; font-size: 12pt; color: rgb(0, 0, 0);">
<b><span style="font-size:12pt"> without hanging if I comment out the final slepc/MPI finalisation command</span></b></div>
<div><span style="font-size:12pt"> CALL SlepcFinalize(ierr_pets)</span> <br>
</div>
<div>(I of course I get ' bad termination' errors, but the otherwise the run is successful.)<br>
</div>
<div> </div>
<div> I would appreciate it if anyone has any ideas on what is going wrong!</div>
<div>  Many thanks,</div>
<div>                       Dan.<br>
</div>
<div><br>
</div>
<div><br>
</div>
<div>code:</div>
<div><br>
</div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      MODULE ALL_STAB_ROUTINES</span>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      IMPLICIT NONE</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      CONTAINS</span></div>
<div><br>
</div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      SUBROUTINE WHOSE_ROW_IS_IT(ROW_NO, TOTAL_NO_ROWS, NO_PROCESSES,     &</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">     &      OWNER)</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">!     THIS ROUTINE ALLOCATES ROWS EVENLY BETWEEN mpi PROCESSES</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">#include <slepc/finclude/slepceps.h></span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      use slepceps</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      IMPLICIT NONE</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      PetscInt, INTENT(IN) :: ROW_NO, TOTAL_NO_ROWS, NO_PROCESSES</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      PetscInt, INTENT(OUT) :: OWNER</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      PetscInt :: P, REM</span></div>
<div><br>
</div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      P = TOTAL_NO_ROWS / NO_PROCESSES ! NOTE INTEGER DIVISION</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      REM = TOTAL_NO_ROWS - P*NO_PROCESSES</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      IF (ROW_NO < (NO_PROCESSES - REM)*P + 1 ) THEN</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">        OWNER = (ROW_NO - 1)/P ! NOTE INTEGER DIVISION</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      ELSE</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">        OWNER = (  ROW_NO  +   NO_PROCESSES - REM -1 )/(P+1) ! NOTE INTEGER DIVISION</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      ENDIF    </span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      END SUBROUTINE WHOSE_ROW_IS_IT</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      END MODULE ALL_STAB_ROUTINES</span></div>
<div><br>
</div>
<div><br>
</div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      PROGRAM trialer</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      USE MPI</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">#include <slepc/finclude/slepceps.h></span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      use slepceps</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      USE ALL_STAB_ROUTINES</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      IMPLICIT NONE</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      PetscMPIInt    rank3, total_mpi_size</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      PetscInt nl3, code,  PROC_ROW, ISTATUS, jm, N_rows,NO_A_ENTRIES</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      PetscInt, ALLOCATABLE, DIMENSION(:) :: JALOC</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      PetscInt, PARAMETER  ::  ZERO = 0 , ONE = 1, TWO = 2, THREE = 3  </span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      PetscErrorCode ierr_pets</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">     
</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">! Initialise sleps/mpi</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      call SlepcInitialize(PETSC_NULL_CHARACTER,ierr_pets) ! note that this initialises MPI</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      call MPI_COMM_SIZE(MPI_COMM_WORLD, total_mpi_size, ierr_pets) !! find total no of MPI processes  
</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      nL3= total_mpi_size</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      call MPI_COMM_RANK(MPI_COMM_WORLD,rank3,ierr_pets) !! find my overall rank -> rank3</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      write(*,*)'Welcome: PROCESS NO , TOTAL NO. OF PROCESSES =  ',rank3, nl3</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;"> </span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      N_rows = 200 ! NUMBER OF ROWS OF A NOTIONAL MATRIX</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      NO_A_ENTRIES = 12 ! NUMBER OF ENTRIES FOR JALOC</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">     
</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">!     LOOP OVER ROWS      </span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      do jm = 1, N_rows</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;"> 
</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      CALL whose_row_is_it(JM,  N_rows , NL3, PROC_ROW) ! FIND OUT WHICH PROCESS OWNS ROW</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      if (rank3 == PROC_ROW) then ! IF mpi PROCESS OWNS THIS ROW THEN ..</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">!       ALLOCATE jaloc ARRAY AND INITIALISE</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">       
</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">        allocate(jaloc(NO_A_ENTRIES), STAT=ISTATUS )</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">        jaloc = three</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">       
</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">       
</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">        WRITE(*,*)'JALOC',JALOC ! THIS SIMPLE PLOT ALWAYS WORKS</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">        write(*,*)'calling PetscIntView: PROCESS NO. ROW NO.',rank3, jm</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">        ! THIS CALL TO PetscIntView CAUSES CODE TO HANG WHEN E.G. total_mpi_size=3, JM=133</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">        call PetscIntView(NO_A_ENTRIES,JALOC(1:NO_A_ENTRIES),              &</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">     &       PETSC_VIEWER_STDOUT_WORLD, ierr_pets)</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">        CHKERRA(ierr_pets)</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">        deallocate(jaloc)</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      endif</span></div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      enddo</span></div>
<div><br>
</div>
<div><span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      CALL SlepcFinalize(ierr_pets)</span></div>
<span style="color: rgb(23, 78, 134); font-family: Calibri, Helvetica, sans-serif;">      end program trialer</span><br>
</div>
<div><br>
</div>
</body>
</html>