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