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>
MakeValue(PetscInt i,PetscInt j,PetscInt M)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
CheckValuesAIJ(Mat A)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
19c4762a1bSJed Brown PetscFunctionBegin;
209566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N));
219566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
22c4762a1bSJed Brown for (i = rstart; i < rend; i++) {
23c4762a1bSJed Brown for (j = 0; j < N; j++) {
249566063dSJacob Faibussowitsch PetscCall(MatGetValue(A, i, j, &val));
259371c9d4SSatish Balay v = MakeValue(i, j, M);
269371c9d4SSatish Balay w = PetscRealPart(val);
2708401ef6SPierre 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);
28c4762a1bSJed Brown }
29c4762a1bSJed Brown }
303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
31c4762a1bSJed Brown }
32c4762a1bSJed Brown
main(int argc,char ** args)33d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
34d71ae5a4SJacob Faibussowitsch {
35c4762a1bSJed Brown Mat A;
36c4762a1bSJed Brown PetscInt M = 24, N = 48, bs = 2;
37c4762a1bSJed Brown PetscInt rstart, rend, i, j;
38c4762a1bSJed Brown PetscViewer view;
39c4762a1bSJed Brown
40327415f7SBarry Smith PetscFunctionBeginUser;
419566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, NULL, help));
42c4762a1bSJed Brown /*
43c4762a1bSJed Brown Create a parallel BAIJ matrix shared by all processors
44c4762a1bSJed Brown */
45d0609cedSBarry Smith PetscCall(MatCreateBAIJ(PETSC_COMM_WORLD, bs, PETSC_DECIDE, PETSC_DECIDE, M, N, PETSC_DECIDE, NULL, PETSC_DECIDE, NULL, &A));
46c4762a1bSJed Brown
47c4762a1bSJed Brown /*
48c4762a1bSJed Brown Set values into the matrix
49c4762a1bSJed Brown */
509566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N));
519566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
52c4762a1bSJed Brown for (i = rstart; i < rend; i++) {
53c4762a1bSJed Brown for (j = 0; j < N; j++) {
54c4762a1bSJed Brown PetscReal v = MakeValue(i, j, M);
5548a46eb9SPierre Jolivet if (PetscAbsReal(v) > 0) PetscCall(MatSetValue(A, i, j, v, INSERT_VALUES));
56c4762a1bSJed Brown }
57c4762a1bSJed Brown }
589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
609566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-mat_base_view"));
61c4762a1bSJed Brown
62c4762a1bSJed Brown /*
63c4762a1bSJed Brown Store the binary matrix to a file
64c4762a1bSJed Brown */
659566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view));
6648a46eb9SPierre Jolivet for (i = 0; i < 3; i++) PetscCall(MatView(A, view));
679566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view));
689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
69c4762a1bSJed Brown
70c4762a1bSJed Brown /*
71c4762a1bSJed Brown Now reload the matrix and check its values
72c4762a1bSJed Brown */
739566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view));
749566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
759566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATBAIJ));
76c4762a1bSJed Brown for (i = 0; i < 3; i++) {
779566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A));
789566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view));
799566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A));
80c4762a1bSJed Brown }
819566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view));
829566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-mat_load_view"));
839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
84c4762a1bSJed Brown
85c4762a1bSJed Brown /*
86c4762a1bSJed Brown Reload in SEQBAIJ matrix and check its values
87c4762a1bSJed Brown */
889566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, "matrix.dat", FILE_MODE_READ, &view));
899566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A));
909566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQBAIJ));
91c4762a1bSJed Brown for (i = 0; i < 3; i++) {
929566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A));
939566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view));
949566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A));
95c4762a1bSJed Brown }
969566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view));
979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
98c4762a1bSJed Brown
99c4762a1bSJed Brown /*
100c4762a1bSJed Brown Reload in MPIBAIJ matrix and check its values
101c4762a1bSJed Brown */
1029566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_READ, &view));
1039566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1049566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPIBAIJ));
105c4762a1bSJed Brown for (i = 0; i < 3; i++) {
1069566063dSJacob Faibussowitsch if (i > 0) PetscCall(MatZeroEntries(A));
1079566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view));
1089566063dSJacob Faibussowitsch PetscCall(CheckValuesAIJ(A));
109c4762a1bSJed Brown }
1109566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view));
1119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
112c4762a1bSJed Brown
1139566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
114b122ec5aSJacob Faibussowitsch return 0;
115c4762a1bSJed Brown }
116c4762a1bSJed Brown
117c4762a1bSJed Brown /*TEST
118c4762a1bSJed Brown
119c4762a1bSJed Brown testset:
120c4762a1bSJed Brown args: -viewer_binary_mpiio 0
121*3886731fSPierre Jolivet output_file: output/empty.out
122c4762a1bSJed Brown test:
123c4762a1bSJed Brown suffix: stdio_1
124c4762a1bSJed Brown nsize: 1
125c4762a1bSJed Brown test:
126c4762a1bSJed Brown suffix: stdio_2
127c4762a1bSJed Brown nsize: 2
128c4762a1bSJed Brown test:
129c4762a1bSJed Brown suffix: stdio_3
130c4762a1bSJed Brown nsize: 3
131c4762a1bSJed Brown test:
132c4762a1bSJed Brown suffix: stdio_4
133c4762a1bSJed Brown nsize: 4
134c4762a1bSJed Brown test:
135c4762a1bSJed Brown suffix: stdio_5
136c4762a1bSJed Brown nsize: 4
137c4762a1bSJed Brown
138c4762a1bSJed Brown testset:
139c4762a1bSJed Brown requires: mpiio
140c4762a1bSJed Brown args: -viewer_binary_mpiio 1
141*3886731fSPierre Jolivet output_file: output/empty.out
142c4762a1bSJed Brown test:
143c4762a1bSJed Brown suffix: mpiio_1
144c4762a1bSJed Brown nsize: 1
145c4762a1bSJed Brown test:
146c4762a1bSJed Brown suffix: mpiio_2
147c4762a1bSJed Brown nsize: 2
148c4762a1bSJed Brown test:
149c4762a1bSJed Brown suffix: mpiio_3
150c4762a1bSJed Brown nsize: 3
151c4762a1bSJed Brown test:
152c4762a1bSJed Brown suffix: mpiio_4
153c4762a1bSJed Brown nsize: 4
154c4762a1bSJed Brown test:
155c4762a1bSJed Brown suffix: mpiio_5
156c4762a1bSJed Brown nsize: 5
157c4762a1bSJed Brown
158c4762a1bSJed Brown TEST*/
159