1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests binary I/O of matrices and illustrates user-defined event logging.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown /* Note: Most applications would not read and write the same matrix within 7c4762a1bSJed Brown the same program. This example is intended only to demonstrate 8c4762a1bSJed Brown both input and output. */ 9c4762a1bSJed Brown 10*9371c9d4SSatish Balay int main(int argc, char **args) { 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 #if defined(PETSC_USE_LOG) 17c4762a1bSJed Brown PetscLogEvent MATRIX_GENERATE, MATRIX_READ; 18c4762a1bSJed Brown #endif 19c4762a1bSJed Brown 20327415f7SBarry Smith PetscFunctionBeginUser; 219566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 26c4762a1bSJed Brown N = m * n; 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* PART 1: Generate matrix, then write it in binary format */ 29c4762a1bSJed Brown 309566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("Generate Matrix", 0, &MATRIX_GENERATE)); 319566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MATRIX_GENERATE, 0, 0, 0, 0)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Generate matrix */ 349566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 359566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, N, N)); 369566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 379566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 389566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &Istart, &Iend)); 39c4762a1bSJed Brown for (Ii = Istart; Ii < Iend; Ii++) { 40*9371c9d4SSatish Balay v = -1.0; 41*9371c9d4SSatish Balay i = Ii / n; 42*9371c9d4SSatish Balay j = Ii - i * n; 43*9371c9d4SSatish Balay if (i > 0) { 44*9371c9d4SSatish Balay J = Ii - n; 45*9371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, ADD_VALUES)); 46*9371c9d4SSatish Balay } 47*9371c9d4SSatish Balay if (i < m - 1) { 48*9371c9d4SSatish Balay J = Ii + n; 49*9371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, ADD_VALUES)); 50*9371c9d4SSatish Balay } 51*9371c9d4SSatish Balay if (j > 0) { 52*9371c9d4SSatish Balay J = Ii - 1; 53*9371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, ADD_VALUES)); 54*9371c9d4SSatish Balay } 55*9371c9d4SSatish Balay if (j < n - 1) { 56*9371c9d4SSatish Balay J = Ii + 1; 57*9371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, ADD_VALUES)); 58*9371c9d4SSatish Balay } 59*9371c9d4SSatish Balay v = 4.0; 60*9371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, ADD_VALUES)); 61c4762a1bSJed Brown } 629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 649566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "writing matrix in binary to matrix.dat ...\n")); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &viewer)); 689566063dSJacob Faibussowitsch PetscCall(MatView(C, viewer)); 699566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 719566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MATRIX_GENERATE, 0, 0, 0, 0)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* PART 2: Read in matrix in binary format */ 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* All processors wait until test matrix has been dumped */ 769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 77c4762a1bSJed Brown 789566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("Read Matrix", 0, &MATRIX_READ)); 799566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MATRIX_READ, 0, 0, 0, 0)); 809566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "reading matrix in binary from matrix.dat ...\n")); 819566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &viewer)); 829566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 839566063dSJacob Faibussowitsch PetscCall(MatLoad(C, viewer)); 849566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 859566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MATRIX_READ, 0, 0, 0, 0)); 869566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Free data structures */ 899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 90c4762a1bSJed Brown 919566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 92b122ec5aSJacob Faibussowitsch return 0; 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown /*TEST 96c4762a1bSJed Brown 97c4762a1bSJed Brown test: 988cc725e6SPierre Jolivet filter: grep -v " MPI process" 99c4762a1bSJed Brown 100c4762a1bSJed Brown TEST*/ 101