1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Show MatShift BUG happening after copying a matrix with no rows on a process"; 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Contributed by: Eric Chamberland 5c4762a1bSJed Brown */ 6c4762a1bSJed Brown #include <petscmat.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown /* DEFINE this to turn on/off the bug: */ 9c4762a1bSJed Brown #define SET_2nd_PROC_TO_HAVE_NO_LOCAL_LINES 10c4762a1bSJed Brown 11c4762a1bSJed Brown int main(int argc,char **args) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown Mat C; 14c4762a1bSJed Brown PetscInt i,m = 3; 15c4762a1bSJed Brown PetscMPIInt rank,size; 16c4762a1bSJed Brown PetscScalar v; 17c4762a1bSJed Brown Mat lMatA; 18c4762a1bSJed Brown PetscInt locallines; 19c4762a1bSJed Brown PetscInt d_nnz[3] = {0,0,0}; 20c4762a1bSJed Brown PetscInt o_nnz[3] = {0,0,0}; 21c4762a1bSJed Brown 22*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 235f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 245f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 25c4762a1bSJed Brown 262c71b3e2SJacob Faibussowitsch PetscCheckFalse(2 != size,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"Relevant with 2 processes only"); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown #ifdef SET_2nd_PROC_TO_HAVE_NO_LOCAL_LINES 30c4762a1bSJed Brown if (0 == rank) { 31c4762a1bSJed Brown locallines = m; 32c4762a1bSJed Brown d_nnz[0] = 1; 33c4762a1bSJed Brown d_nnz[1] = 1; 34c4762a1bSJed Brown d_nnz[2] = 1; 35c4762a1bSJed Brown } else { 36c4762a1bSJed Brown locallines = 0; 37c4762a1bSJed Brown } 38c4762a1bSJed Brown #else 39c4762a1bSJed Brown if (0 == rank) { 40c4762a1bSJed Brown locallines = m-1; 41c4762a1bSJed Brown d_nnz[0] = 1; 42c4762a1bSJed Brown d_nnz[1] = 1; 43c4762a1bSJed Brown } else { 44c4762a1bSJed Brown locallines = 1; 45c4762a1bSJed Brown d_nnz[0] = 1; 46c4762a1bSJed Brown } 47c4762a1bSJed Brown #endif 48c4762a1bSJed Brown 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,locallines,locallines,m,m)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatXAIJSetPreallocation(C,1,d_nnz,o_nnz,NULL,NULL)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown v = 2; 54c4762a1bSJed Brown /* Assembly on the diagonal: */ 55c4762a1bSJed Brown for (i=0; i<m; i++) { 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,1,&i,&v,ADD_VALUES)); 57c4762a1bSJed Brown } 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(C, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATSAME, MAT_INITIAL_MATRIX, &lMatA)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(lMatA,PETSC_VIEWER_STDOUT_WORLD)); 65c4762a1bSJed Brown 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(lMatA,-1.0)); 67c4762a1bSJed Brown 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&lMatA)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 70*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 71*b122ec5aSJacob Faibussowitsch return 0; 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown /*TEST 75c4762a1bSJed Brown 76c4762a1bSJed Brown test: 77c4762a1bSJed Brown nsize: 2 78c4762a1bSJed Brown 79c4762a1bSJed Brown TEST*/ 80