1c4762a1bSJed Brown static char help[] = "Tests binary I/O of matrices and illustrates user-defined event logging.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown /* Note: Most applications would not read and write the same matrix within 6c4762a1bSJed Brown the same program. This example is intended only to demonstrate 7c4762a1bSJed Brown both input and output. */ 8c4762a1bSJed Brown 9d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 10d71ae5a4SJacob Faibussowitsch { 11c4762a1bSJed Brown Mat C; 12c4762a1bSJed Brown PetscScalar v; 13c4762a1bSJed Brown PetscInt i, j, Ii, J, Istart, Iend, N, m = 4, n = 4; 14c4762a1bSJed Brown PetscMPIInt rank, size; 15c4762a1bSJed Brown PetscViewer viewer; 16c4762a1bSJed Brown PetscLogEvent MATRIX_GENERATE, MATRIX_READ; 17c4762a1bSJed Brown 18327415f7SBarry Smith PetscFunctionBeginUser; 19*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 24c4762a1bSJed Brown N = m * n; 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* PART 1: Generate matrix, then write it in binary format */ 27c4762a1bSJed Brown 289566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("Generate Matrix", 0, &MATRIX_GENERATE)); 299566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MATRIX_GENERATE, 0, 0, 0, 0)); 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* Generate matrix */ 329566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 339566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, N, N)); 349566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 359566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 369566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &Istart, &Iend)); 37c4762a1bSJed Brown for (Ii = Istart; Ii < Iend; Ii++) { 389371c9d4SSatish Balay v = -1.0; 399371c9d4SSatish Balay i = Ii / n; 409371c9d4SSatish Balay j = Ii - i * n; 419371c9d4SSatish Balay if (i > 0) { 429371c9d4SSatish Balay J = Ii - n; 439371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, ADD_VALUES)); 449371c9d4SSatish Balay } 459371c9d4SSatish Balay if (i < m - 1) { 469371c9d4SSatish Balay J = Ii + n; 479371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, ADD_VALUES)); 489371c9d4SSatish Balay } 499371c9d4SSatish Balay if (j > 0) { 509371c9d4SSatish Balay J = Ii - 1; 519371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, ADD_VALUES)); 529371c9d4SSatish Balay } 539371c9d4SSatish Balay if (j < n - 1) { 549371c9d4SSatish Balay J = Ii + 1; 559371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, ADD_VALUES)); 569371c9d4SSatish Balay } 579371c9d4SSatish Balay v = 4.0; 589371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, ADD_VALUES)); 59c4762a1bSJed Brown } 609566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 619566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 629566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 63c4762a1bSJed Brown 649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "writing matrix in binary to matrix.dat ...\n")); 659566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &viewer)); 669566063dSJacob Faibussowitsch PetscCall(MatView(C, viewer)); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 699566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MATRIX_GENERATE, 0, 0, 0, 0)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* PART 2: Read in matrix in binary format */ 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* All processors wait until test matrix has been dumped */ 749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("Read Matrix", 0, &MATRIX_READ)); 779566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MATRIX_READ, 0, 0, 0, 0)); 789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "reading matrix in binary from matrix.dat ...\n")); 799566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &viewer)); 809566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 819566063dSJacob Faibussowitsch PetscCall(MatLoad(C, viewer)); 829566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 839566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MATRIX_READ, 0, 0, 0, 0)); 849566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Free data structures */ 879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 90b122ec5aSJacob Faibussowitsch return 0; 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93c4762a1bSJed Brown /*TEST 94c4762a1bSJed Brown 95c4762a1bSJed Brown test: 968cc725e6SPierre Jolivet filter: grep -v " MPI process" 97c4762a1bSJed Brown 98c4762a1bSJed Brown TEST*/ 99