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*327415f7SBarry Smith PetscFunctionBeginUser; 239566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 26c4762a1bSJed Brown 2708401ef6SPierre Jolivet PetscCheck(2 == size,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"Relevant with 2 processes only"); 289566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 29c4762a1bSJed Brown 30c4762a1bSJed Brown #ifdef SET_2nd_PROC_TO_HAVE_NO_LOCAL_LINES 31c4762a1bSJed Brown if (0 == rank) { 32c4762a1bSJed Brown locallines = m; 33c4762a1bSJed Brown d_nnz[0] = 1; 34c4762a1bSJed Brown d_nnz[1] = 1; 35c4762a1bSJed Brown d_nnz[2] = 1; 36c4762a1bSJed Brown } else { 37c4762a1bSJed Brown locallines = 0; 38c4762a1bSJed Brown } 39c4762a1bSJed Brown #else 40c4762a1bSJed Brown if (0 == rank) { 41c4762a1bSJed Brown locallines = m-1; 42c4762a1bSJed Brown d_nnz[0] = 1; 43c4762a1bSJed Brown d_nnz[1] = 1; 44c4762a1bSJed Brown } else { 45c4762a1bSJed Brown locallines = 1; 46c4762a1bSJed Brown d_nnz[0] = 1; 47c4762a1bSJed Brown } 48c4762a1bSJed Brown #endif 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,locallines,locallines,m,m)); 519566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 529566063dSJacob Faibussowitsch PetscCall(MatXAIJSetPreallocation(C,1,d_nnz,o_nnz,NULL,NULL)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown v = 2; 55c4762a1bSJed Brown /* Assembly on the diagonal: */ 56c4762a1bSJed Brown for (i=0; i<m; i++) { 579566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,1,&i,&v,ADD_VALUES)); 58c4762a1bSJed Brown } 599566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 609566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 619566063dSJacob Faibussowitsch PetscCall(MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 629566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 639566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 649566063dSJacob Faibussowitsch PetscCall(MatConvert(C,MATSAME, MAT_INITIAL_MATRIX, &lMatA)); 659566063dSJacob Faibussowitsch PetscCall(MatView(lMatA,PETSC_VIEWER_STDOUT_WORLD)); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(MatShift(lMatA,-1.0)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lMatA)); 709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 719566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 72b122ec5aSJacob Faibussowitsch return 0; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown /*TEST 76c4762a1bSJed Brown 77c4762a1bSJed Brown test: 78c4762a1bSJed Brown nsize: 2 79c4762a1bSJed Brown 80c4762a1bSJed Brown TEST*/ 81