1c4762a1bSJed Brown static char help[] = "Tests converting a parallel AIJ formatted matrix to the parallel Row format.\n\ 2c4762a1bSJed Brown This also tests MatGetRow() and MatRestoreRow() for the parallel case.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 7d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown Mat C, A; 9c4762a1bSJed Brown PetscInt i, j, m = 3, n = 2, Ii, J, rstart, rend, nz; 10c4762a1bSJed Brown PetscMPIInt rank, size; 11c4762a1bSJed Brown const PetscInt *idx; 12c4762a1bSJed Brown PetscScalar v; 13c4762a1bSJed Brown const PetscScalar *values; 14c4762a1bSJed Brown 15327415f7SBarry Smith PetscFunctionBeginUser; 16*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 19c4762a1bSJed Brown n = 2 * size; 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN*/ 229566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 249566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 259566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(C, 5, NULL, 5, NULL)); 269566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(C, 5, NULL)); 27c4762a1bSJed Brown for (i = 0; i < m; i++) { 28c4762a1bSJed Brown for (j = 2 * rank; j < 2 * rank + 2; j++) { 299371c9d4SSatish Balay v = -1.0; 309371c9d4SSatish Balay Ii = j + n * i; 319371c9d4SSatish Balay if (i > 0) { 329371c9d4SSatish Balay J = Ii - n; 339371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 349371c9d4SSatish Balay } 359371c9d4SSatish Balay if (i < m - 1) { 369371c9d4SSatish Balay J = Ii + n; 379371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 389371c9d4SSatish Balay } 399371c9d4SSatish Balay if (j > 0) { 409371c9d4SSatish Balay J = Ii - 1; 419371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 429371c9d4SSatish Balay } 439371c9d4SSatish Balay if (j < n - 1) { 449371c9d4SSatish Balay J = Ii + 1; 459371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 469371c9d4SSatish Balay } 479371c9d4SSatish Balay v = 4.0; 489371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown } 519566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 529566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 539566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO)); 549566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 559566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 569566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 57c4762a1bSJed Brown 589566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend)); 59c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 609566063dSJacob Faibussowitsch PetscCall(MatGetRow(C, i, &nz, &idx, &values)); 619566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD, stdout, "[%d] get row %" PetscInt_FMT ": ", rank, i)); 6248a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD, stdout, "%" PetscInt_FMT " %g ", idx[j], (double)PetscRealPart(values[j]))); 639566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD, stdout, "\n")); 649566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(C, i, &nz, &idx, &values)); 65c4762a1bSJed Brown } 669566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, stdout)); 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATSAME, MAT_INITIAL_MATRIX, &A)); 699566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO)); 709566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 729566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 73c4762a1bSJed Brown 749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 769566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 77b122ec5aSJacob Faibussowitsch return 0; 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80c4762a1bSJed Brown /*TEST 81c4762a1bSJed Brown 82c4762a1bSJed Brown test: 83c4762a1bSJed Brown args: -mat_type seqaij 84c4762a1bSJed Brown 85c4762a1bSJed Brown TEST*/ 86