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