1c4762a1bSJed Brown static char help[] = "Tests MatView()/MatLoad() with binary viewers for BAIJ matrices.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown #include <petscviewer.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown #include <petsc/private/hashtable.h> 79371c9d4SSatish Balay static PetscReal MakeValue(PetscInt i, PetscInt j, PetscInt M) { 8c4762a1bSJed Brown PetscHash_t h = PetscHashCombine(PetscHashInt(i), PetscHashInt(j)); 9c4762a1bSJed Brown return (PetscReal)((h % 5 == 0) ? (1 + i + j * M) : 0); 10c4762a1bSJed Brown } 11c4762a1bSJed Brown 129371c9d4SSatish Balay static PetscErrorCode CheckValuesAIJ(Mat A) { 13c4762a1bSJed Brown PetscInt M, N, rstart, rend, i, j; 14c4762a1bSJed Brown PetscReal v, w; 15c4762a1bSJed Brown PetscScalar val; 16c4762a1bSJed Brown 17c4762a1bSJed Brown PetscFunctionBegin; 189566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 199566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 20c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 21c4762a1bSJed Brown for (j = 0; j < N; j++) { 229566063dSJacob Faibussowitsch PetscCall(MatGetValue(A, i, j, &val)); 239371c9d4SSatish Balay v = MakeValue(i, j, M); 249371c9d4SSatish Balay w = PetscRealPart(val); 2508401ef6SPierre Jolivet PetscCheck(PetscAbsReal(v - w) <= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Matrix entry (%" PetscInt_FMT ",%" PetscInt_FMT ") should be %g, got %g", i, j, (double)v, (double)w); 26c4762a1bSJed Brown } 27c4762a1bSJed Brown } 28c4762a1bSJed Brown PetscFunctionReturn(0); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 319371c9d4SSatish Balay int main(int argc, char **args) { 32c4762a1bSJed Brown Mat A; 33c4762a1bSJed Brown PetscInt M = 24, N = 48, bs = 2; 34c4762a1bSJed Brown PetscInt rstart, rend, i, j; 35c4762a1bSJed Brown PetscViewer view; 36c4762a1bSJed Brown 37327415f7SBarry Smith PetscFunctionBeginUser; 389566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, NULL, help)); 39c4762a1bSJed Brown /* 40c4762a1bSJed Brown Create a parallel BAIJ matrix shared by all processors 41c4762a1bSJed Brown */ 42d0609cedSBarry Smith PetscCall(MatCreateBAIJ(PETSC_COMM_WORLD, bs, PETSC_DECIDE, PETSC_DECIDE, M, N, PETSC_DECIDE, NULL, PETSC_DECIDE, NULL, &A)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* 45c4762a1bSJed Brown Set values into the matrix 46c4762a1bSJed Brown */ 479566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 489566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 49c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 50c4762a1bSJed Brown for (j = 0; j < N; j++) { 51c4762a1bSJed Brown PetscReal v = MakeValue(i, j, M); 52*48a46eb9SPierre Jolivet if (PetscAbsReal(v) > 0) PetscCall(MatSetValue(A, i, j, v, INSERT_VALUES)); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown } 559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 579566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-mat_base_view")); 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* 60c4762a1bSJed Brown Store the binary matrix to a file 61c4762a1bSJed Brown */ 629566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view)); 63*48a46eb9SPierre Jolivet for (i = 0; i < 3; i++) PetscCall(MatView(A, view)); 649566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* 68c4762a1bSJed Brown Now reload the matrix and check its values 69c4762a1bSJed Brown */ 709566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view)); 719566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 729566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATBAIJ)); 73c4762a1bSJed Brown for (i = 0; i < 3; i++) { 749566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 759566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 769566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 77c4762a1bSJed Brown } 789566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 799566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-mat_load_view")); 809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* 83c4762a1bSJed Brown Reload in SEQBAIJ matrix and check its values 84c4762a1bSJed Brown */ 859566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, "matrix.dat", FILE_MODE_READ, &view)); 869566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A)); 879566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQBAIJ)); 88c4762a1bSJed Brown for (i = 0; i < 3; i++) { 899566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 909566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 919566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 92c4762a1bSJed Brown } 939566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* 97c4762a1bSJed Brown Reload in MPIBAIJ matrix and check its values 98c4762a1bSJed Brown */ 999566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view)); 1009566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1019566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPIBAIJ)); 102c4762a1bSJed Brown for (i = 0; i < 3; i++) { 1039566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 1049566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 1059566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 106c4762a1bSJed Brown } 1079566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 111b122ec5aSJacob Faibussowitsch return 0; 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown /*TEST 115c4762a1bSJed Brown 116c4762a1bSJed Brown testset: 117c4762a1bSJed Brown args: -viewer_binary_mpiio 0 118c4762a1bSJed Brown output_file: output/ex45.out 119c4762a1bSJed Brown test: 120c4762a1bSJed Brown suffix: stdio_1 121c4762a1bSJed Brown nsize: 1 122c4762a1bSJed Brown test: 123c4762a1bSJed Brown suffix: stdio_2 124c4762a1bSJed Brown nsize: 2 125c4762a1bSJed Brown test: 126c4762a1bSJed Brown suffix: stdio_3 127c4762a1bSJed Brown nsize: 3 128c4762a1bSJed Brown test: 129c4762a1bSJed Brown suffix: stdio_4 130c4762a1bSJed Brown nsize: 4 131c4762a1bSJed Brown test: 132c4762a1bSJed Brown suffix: stdio_5 133c4762a1bSJed Brown nsize: 4 134c4762a1bSJed Brown 135c4762a1bSJed Brown testset: 136c4762a1bSJed Brown requires: mpiio 137c4762a1bSJed Brown args: -viewer_binary_mpiio 1 138c4762a1bSJed Brown output_file: output/ex45.out 139c4762a1bSJed Brown test: 140c4762a1bSJed Brown suffix: mpiio_1 141c4762a1bSJed Brown nsize: 1 142c4762a1bSJed Brown test: 143c4762a1bSJed Brown suffix: mpiio_2 144c4762a1bSJed Brown nsize: 2 145c4762a1bSJed Brown test: 146c4762a1bSJed Brown suffix: mpiio_3 147c4762a1bSJed Brown nsize: 3 148c4762a1bSJed Brown test: 149c4762a1bSJed Brown suffix: mpiio_4 150c4762a1bSJed Brown nsize: 4 151c4762a1bSJed Brown test: 152c4762a1bSJed Brown suffix: mpiio_5 153c4762a1bSJed Brown nsize: 5 154c4762a1bSJed Brown 155c4762a1bSJed Brown TEST*/ 156