1c4762a1bSJed Brown static char help[] = "Show MatShift BUG happening after copying a matrix with no rows on a process"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown Contributed by: Eric Chamberland 4c4762a1bSJed Brown */ 5c4762a1bSJed Brown #include <petscmat.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown /* DEFINE this to turn on/off the bug: */ 8c4762a1bSJed Brown #define SET_2nd_PROC_TO_HAVE_NO_LOCAL_LINES 9c4762a1bSJed Brown 10d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 11d71ae5a4SJacob Faibussowitsch { 12c4762a1bSJed Brown Mat C; 13c4762a1bSJed Brown PetscInt i, m = 3; 14c4762a1bSJed Brown PetscMPIInt rank, size; 15c4762a1bSJed Brown PetscScalar v; 16c4762a1bSJed Brown Mat lMatA; 17c4762a1bSJed Brown PetscInt locallines; 18c4762a1bSJed Brown PetscInt d_nnz[3] = {0, 0, 0}; 19c4762a1bSJed Brown PetscInt o_nnz[3] = {0, 0, 0}; 20c4762a1bSJed Brown 21327415f7SBarry Smith PetscFunctionBeginUser; 22*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 25c4762a1bSJed Brown 2608401ef6SPierre Jolivet PetscCheck(2 == size, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Relevant with 2 processes only"); 279566063dSJacob Faibussowitsch PetscCall(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 499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, locallines, locallines, m, m)); 509566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 519566063dSJacob Faibussowitsch PetscCall(MatXAIJSetPreallocation(C, 1, d_nnz, o_nnz, NULL, NULL)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown v = 2; 54c4762a1bSJed Brown /* Assembly on the diagonal: */ 5548a46eb9SPierre Jolivet for (i = 0; i < m; i++) PetscCall(MatSetValues(C, 1, &i, 1, &i, &v, ADD_VALUES)); 569566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 579566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 589566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 599566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 609566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 619566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATSAME, MAT_INITIAL_MATRIX, &lMatA)); 629566063dSJacob Faibussowitsch PetscCall(MatView(lMatA, PETSC_VIEWER_STDOUT_WORLD)); 63c4762a1bSJed Brown 649566063dSJacob Faibussowitsch PetscCall(MatShift(lMatA, -1.0)); 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&lMatA)); 679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 689566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 69b122ec5aSJacob Faibussowitsch return 0; 70c4762a1bSJed Brown } 71c4762a1bSJed Brown 72c4762a1bSJed Brown /*TEST 73c4762a1bSJed Brown 74c4762a1bSJed Brown test: 75c4762a1bSJed Brown nsize: 2 76c4762a1bSJed Brown 77c4762a1bSJed Brown TEST*/ 78