1 2 static char help[] = "Show MatShift BUG happening after copying a matrix with no rows on a process"; 3 /* 4 Contributed by: Eric Chamberland 5 */ 6 #include <petscmat.h> 7 8 /* DEFINE this to turn on/off the bug: */ 9 #define SET_2nd_PROC_TO_HAVE_NO_LOCAL_LINES 10 11 int main(int argc,char **args) 12 { 13 Mat C; 14 PetscInt i,m = 3; 15 PetscMPIInt rank,size; 16 PetscScalar v; 17 Mat lMatA; 18 PetscInt locallines; 19 PetscInt d_nnz[3] = {0,0,0}; 20 PetscInt o_nnz[3] = {0,0,0}; 21 22 PetscFunctionBeginUser; 23 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 24 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 25 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 26 27 PetscCheck(2 == size,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"Relevant with 2 processes only"); 28 PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 29 30 #ifdef SET_2nd_PROC_TO_HAVE_NO_LOCAL_LINES 31 if (0 == rank) { 32 locallines = m; 33 d_nnz[0] = 1; 34 d_nnz[1] = 1; 35 d_nnz[2] = 1; 36 } else { 37 locallines = 0; 38 } 39 #else 40 if (0 == rank) { 41 locallines = m-1; 42 d_nnz[0] = 1; 43 d_nnz[1] = 1; 44 } else { 45 locallines = 1; 46 d_nnz[0] = 1; 47 } 48 #endif 49 50 PetscCall(MatSetSizes(C,locallines,locallines,m,m)); 51 PetscCall(MatSetFromOptions(C)); 52 PetscCall(MatXAIJSetPreallocation(C,1,d_nnz,o_nnz,NULL,NULL)); 53 54 v = 2; 55 /* Assembly on the diagonal: */ 56 for (i=0; i<m; i++) { 57 PetscCall(MatSetValues(C,1,&i,1,&i,&v,ADD_VALUES)); 58 } 59 PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 60 PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 61 PetscCall(MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 62 PetscCall(MatSetOption(C, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 63 PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 64 PetscCall(MatConvert(C,MATSAME, MAT_INITIAL_MATRIX, &lMatA)); 65 PetscCall(MatView(lMatA,PETSC_VIEWER_STDOUT_WORLD)); 66 67 PetscCall(MatShift(lMatA,-1.0)); 68 69 PetscCall(MatDestroy(&lMatA)); 70 PetscCall(MatDestroy(&C)); 71 PetscCall(PetscFinalize()); 72 return 0; 73 } 74 75 /*TEST 76 77 test: 78 nsize: 2 79 80 TEST*/ 81