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> 7d71ae5a4SJacob Faibussowitsch static PetscReal MakeValue(PetscInt i, PetscInt j, PetscInt M) 8d71ae5a4SJacob Faibussowitsch { 9c4762a1bSJed Brown PetscHash_t h = PetscHashCombine(PetscHashInt(i), PetscHashInt(j)); 10c4762a1bSJed Brown return (PetscReal)((h % 5 == 0) ? (1 + i + j * M) : 0); 11c4762a1bSJed Brown } 12c4762a1bSJed Brown 13d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckValuesAIJ(Mat A) 14d71ae5a4SJacob Faibussowitsch { 15c4762a1bSJed Brown PetscInt M, N, rstart, rend, i, j; 16c4762a1bSJed Brown PetscReal v, w; 17c4762a1bSJed Brown PetscScalar val; 18c4762a1bSJed Brown PetscBool seqsbaij, mpisbaij, sbaij; 19c4762a1bSJed Brown 20c4762a1bSJed Brown PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 229566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &seqsbaij)); 249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISBAIJ, &mpisbaij)); 25c4762a1bSJed Brown sbaij = (seqsbaij || mpisbaij) ? PETSC_TRUE : PETSC_FALSE; 26c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 27c4762a1bSJed Brown for (j = (sbaij ? i : 0); j < N; j++) { 289566063dSJacob Faibussowitsch PetscCall(MatGetValue(A, i, j, &val)); 299371c9d4SSatish Balay v = MakeValue(i, j, M); 309371c9d4SSatish Balay w = PetscRealPart(val); 3108401ef6SPierre 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); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown } 343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown 37d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 38d71ae5a4SJacob Faibussowitsch { 39c4762a1bSJed Brown Mat A; 40c4762a1bSJed Brown PetscInt M = 24, N = 24, bs = 3; 41c4762a1bSJed Brown PetscInt rstart, rend, i, j; 42c4762a1bSJed Brown PetscViewer view; 43c4762a1bSJed Brown 44327415f7SBarry Smith PetscFunctionBeginUser; 459566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, NULL, help)); 46c4762a1bSJed Brown /* 47c4762a1bSJed Brown Create a parallel SBAIJ matrix shared by all processors 48c4762a1bSJed Brown */ 499566063dSJacob Faibussowitsch PetscCall(MatCreateSBAIJ(PETSC_COMM_WORLD, bs, PETSC_DECIDE, PETSC_DECIDE, M, N, PETSC_DECIDE, NULL, PETSC_DECIDE, NULL, &A)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* 52c4762a1bSJed Brown Set values into the matrix 53c4762a1bSJed Brown */ 549566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 559566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 56c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 57c4762a1bSJed Brown for (j = 0; j < N; j++) { 58c4762a1bSJed Brown PetscReal v = MakeValue(i, j, M); 5948a46eb9SPierre Jolivet if (PetscAbsReal(v) > 0) PetscCall(MatSetValue(A, i, j, v, INSERT_VALUES)); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown } 629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 649566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-mat_base_view")); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* 67c4762a1bSJed Brown Store the binary matrix to a file 68c4762a1bSJed Brown */ 699566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view)); 7048a46eb9SPierre Jolivet for (i = 0; i < 3; i++) PetscCall(MatView(A, view)); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* 75c4762a1bSJed Brown Now reload the matrix and check its values 76c4762a1bSJed Brown */ 779566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view)); 789566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 799566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSBAIJ)); 80c4762a1bSJed Brown for (i = 0; i < 3; i++) { 819566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 829566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 839566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 84c4762a1bSJed Brown } 859566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 869566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-mat_load_view")); 879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* 90c4762a1bSJed Brown Reload in SEQSBAIJ matrix and check its values 91c4762a1bSJed Brown */ 929566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, "matrix.dat", FILE_MODE_READ, &view)); 939566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A)); 949566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQSBAIJ)); 95c4762a1bSJed Brown for (i = 0; i < 3; i++) { 969566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 979566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 989566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 99c4762a1bSJed Brown } 1009566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* 104c4762a1bSJed Brown Reload in MPISBAIJ matrix and check its values 105c4762a1bSJed Brown */ 1069566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view)); 1079566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1089566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPISBAIJ)); 109c4762a1bSJed Brown for (i = 0; i < 3; i++) { 1109566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A)); 1119566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 1129566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A)); 113c4762a1bSJed Brown } 1149566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 116c4762a1bSJed Brown 1179566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 118b122ec5aSJacob Faibussowitsch return 0; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown /*TEST 122c4762a1bSJed Brown 123c4762a1bSJed Brown testset: 124c4762a1bSJed Brown args: -viewer_binary_mpiio 0 125*3886731fSPierre Jolivet output_file: output/empty.out 126c4762a1bSJed Brown test: 127c4762a1bSJed Brown suffix: stdio_1 128c4762a1bSJed Brown nsize: 1 129c4762a1bSJed Brown test: 130c4762a1bSJed Brown suffix: stdio_2 131c4762a1bSJed Brown nsize: 2 132c4762a1bSJed Brown test: 133c4762a1bSJed Brown suffix: stdio_3 134c4762a1bSJed Brown nsize: 3 135c4762a1bSJed Brown test: 136c4762a1bSJed Brown suffix: stdio_4 137c4762a1bSJed Brown nsize: 4 138c4762a1bSJed Brown test: 139c4762a1bSJed Brown suffix: stdio_5 140c4762a1bSJed Brown nsize: 4 141c4762a1bSJed Brown 142c4762a1bSJed Brown testset: 143c4762a1bSJed Brown requires: mpiio 144c4762a1bSJed Brown args: -viewer_binary_mpiio 1 145*3886731fSPierre Jolivet output_file: output/empty.out 146c4762a1bSJed Brown test: 147c4762a1bSJed Brown suffix: mpiio_1 148c4762a1bSJed Brown nsize: 1 149c4762a1bSJed Brown test: 150c4762a1bSJed Brown suffix: mpiio_2 151c4762a1bSJed Brown nsize: 2 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: mpiio_3 154c4762a1bSJed Brown nsize: 3 155c4762a1bSJed Brown test: 156c4762a1bSJed Brown suffix: mpiio_4 157c4762a1bSJed Brown nsize: 4 158c4762a1bSJed Brown test: 159c4762a1bSJed Brown suffix: mpiio_5 160c4762a1bSJed Brown nsize: 5 161c4762a1bSJed Brown 162c4762a1bSJed Brown TEST*/ 163