1c4762a1bSJed Brown static char help[] = "Tests MatView()/MatLoad() with binary viewers for SBAIJ 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 PetscBool seqsbaij, mpisbaij, sbaij; 17c4762a1bSJed Brown 18c4762a1bSJed Brown PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 209566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &seqsbaij)); 229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISBAIJ, &mpisbaij)); 23c4762a1bSJed Brown sbaij = (seqsbaij || mpisbaij) ? PETSC_TRUE : PETSC_FALSE; 24c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 25c4762a1bSJed Brown for (j = (sbaij ? i : 0); j < N; j++) { 269566063dSJacob Faibussowitsch PetscCall(MatGetValue(A, i, j, &val)); 279371c9d4SSatish Balay v = MakeValue(i, j, M); 289371c9d4SSatish Balay w = PetscRealPart(val); 2908401ef6SPierre 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); 30c4762a1bSJed Brown } 31c4762a1bSJed Brown } 32c4762a1bSJed Brown PetscFunctionReturn(0); 33c4762a1bSJed Brown } 34c4762a1bSJed Brown 359371c9d4SSatish Balay int main(int argc, char **args) { 36c4762a1bSJed Brown Mat A; 37c4762a1bSJed Brown PetscInt M = 24, N = 24, bs = 3; 38c4762a1bSJed Brown PetscInt rstart, rend, i, j; 39c4762a1bSJed Brown PetscViewer view; 40c4762a1bSJed Brown 41327415f7SBarry Smith PetscFunctionBeginUser; 429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, NULL, help)); 43c4762a1bSJed Brown /* 44c4762a1bSJed Brown Create a parallel SBAIJ matrix shared by all processors 45c4762a1bSJed Brown */ 469566063dSJacob Faibussowitsch PetscCall(MatCreateSBAIJ(PETSC_COMM_WORLD, bs, PETSC_DECIDE, PETSC_DECIDE, M, N, PETSC_DECIDE, NULL, PETSC_DECIDE, NULL, &A)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* 49c4762a1bSJed Brown Set values into the matrix 50c4762a1bSJed Brown */ 519566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 529566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 53c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 54c4762a1bSJed Brown for (j = 0; j < N; j++) { 55c4762a1bSJed Brown PetscReal v = MakeValue(i, j, M); 56*48a46eb9SPierre Jolivet if (PetscAbsReal(v) > 0) PetscCall(MatSetValue(A, i, j, v, INSERT_VALUES)); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown } 599566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 609566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 619566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-mat_base_view")); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* 64c4762a1bSJed Brown Store the binary matrix to a file 65c4762a1bSJed Brown */ 669566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view)); 67*48a46eb9SPierre Jolivet for (i = 0; i < 3; i++) PetscCall(MatView(A, view)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* 72c4762a1bSJed Brown Now reload the matrix and check its values 73c4762a1bSJed Brown */ 749566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view)); 759566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 769566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSBAIJ)); 77c4762a1bSJed Brown for (i = 0; i < 3; i++) { 789566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 799566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 809566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 81c4762a1bSJed Brown } 829566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 839566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-mat_load_view")); 849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* 87c4762a1bSJed Brown Reload in SEQSBAIJ matrix and check its values 88c4762a1bSJed Brown */ 899566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, "matrix.dat", FILE_MODE_READ, &view)); 909566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A)); 919566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQSBAIJ)); 92c4762a1bSJed Brown for (i = 0; i < 3; i++) { 939566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 949566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 959566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 96c4762a1bSJed Brown } 979566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* 101c4762a1bSJed Brown Reload in MPISBAIJ matrix and check its values 102c4762a1bSJed Brown */ 1039566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view)); 1049566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1059566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPISBAIJ)); 106c4762a1bSJed Brown for (i = 0; i < 3; i++) { 1079566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 1089566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 1099566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 110c4762a1bSJed Brown } 1119566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 115b122ec5aSJacob Faibussowitsch return 0; 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118c4762a1bSJed Brown /*TEST 119c4762a1bSJed Brown 120c4762a1bSJed Brown testset: 121c4762a1bSJed Brown args: -viewer_binary_mpiio 0 122c4762a1bSJed Brown output_file: output/ex50.out 123c4762a1bSJed Brown test: 124c4762a1bSJed Brown suffix: stdio_1 125c4762a1bSJed Brown nsize: 1 126c4762a1bSJed Brown test: 127c4762a1bSJed Brown suffix: stdio_2 128c4762a1bSJed Brown nsize: 2 129c4762a1bSJed Brown test: 130c4762a1bSJed Brown suffix: stdio_3 131c4762a1bSJed Brown nsize: 3 132c4762a1bSJed Brown test: 133c4762a1bSJed Brown suffix: stdio_4 134c4762a1bSJed Brown nsize: 4 135c4762a1bSJed Brown test: 136c4762a1bSJed Brown suffix: stdio_5 137c4762a1bSJed Brown nsize: 4 138c4762a1bSJed Brown 139c4762a1bSJed Brown testset: 140c4762a1bSJed Brown requires: mpiio 141c4762a1bSJed Brown args: -viewer_binary_mpiio 1 142c4762a1bSJed Brown output_file: output/ex50.out 143c4762a1bSJed Brown test: 144c4762a1bSJed Brown suffix: mpiio_1 145c4762a1bSJed Brown nsize: 1 146c4762a1bSJed Brown test: 147c4762a1bSJed Brown suffix: mpiio_2 148c4762a1bSJed Brown nsize: 2 149c4762a1bSJed Brown test: 150c4762a1bSJed Brown suffix: mpiio_3 151c4762a1bSJed Brown nsize: 3 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: mpiio_4 154c4762a1bSJed Brown nsize: 4 155c4762a1bSJed Brown test: 156c4762a1bSJed Brown suffix: mpiio_5 157c4762a1bSJed Brown nsize: 5 158c4762a1bSJed Brown 159c4762a1bSJed Brown TEST*/ 160