xref: /petsc/src/vec/is/tests/ex2.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Tests ISView() and ISLoad() \n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscis.h>
4c4762a1bSJed Brown #include <petscviewer.h>
5c4762a1bSJed Brown 
6d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
7d71ae5a4SJacob Faibussowitsch {
82cf5aabcSBarry Smith   PetscInt n           = 3, *izero, j, i;
99371c9d4SSatish Balay   PetscInt ix[3][3][3] = {
109371c9d4SSatish Balay     {{3, 5, 4}, {1, 7, 9}, {0, 2, 8}},
11c4762a1bSJed Brown     {{0, 2, 8}, {3, 5, 4}, {1, 7, 9}},
129371c9d4SSatish Balay     {{1, 7, 9}, {0, 2, 8}, {3, 5, 4}}
139371c9d4SSatish Balay   };
14c4762a1bSJed Brown   IS          isx[3], il;
15c4762a1bSJed Brown   PetscMPIInt size, rank;
16c4762a1bSJed Brown   PetscViewer vx, vl;
17c4762a1bSJed Brown   PetscBool   equal;
18c4762a1bSJed Brown 
19327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
23be096a46SBarry Smith   PetscCheck(size < 4, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example only works with up to three processes");
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(size * n, &izero));
2648a46eb9SPierre Jolivet   for (i = 0; i < 3; i++) PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, ix[i][rank], PETSC_COPY_VALUES, &isx[i]));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   for (j = 0; j < 3; j++) {
299566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "testfile", FILE_MODE_WRITE, &vx));
309566063dSJacob Faibussowitsch     PetscCall(ISView(isx[0], vx));
319566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vx));
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "testfile", FILE_MODE_READ, &vl));
349566063dSJacob Faibussowitsch     PetscCall(ISCreate(PETSC_COMM_WORLD, &il));
359566063dSJacob Faibussowitsch     PetscCall(ISLoad(il, vl));
369566063dSJacob Faibussowitsch     PetscCall(ISEqual(il, isx[0], &equal));
3728b400f6SJacob Faibussowitsch     PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Iteration %" PetscInt_FMT " - Index set loaded from file does not match", j);
389566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&il));
399566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "testfile", FILE_MODE_APPEND, &vx));
429566063dSJacob Faibussowitsch     PetscCall(ISView(isx[1], vx));
439566063dSJacob Faibussowitsch     PetscCall(ISView(isx[2], vx));
449566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vx));
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "testfile", FILE_MODE_READ, &vl));
47c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
489566063dSJacob Faibussowitsch       PetscCall(ISCreate(PETSC_COMM_WORLD, &il));
499566063dSJacob Faibussowitsch       PetscCall(ISLoad(il, vl));
509566063dSJacob Faibussowitsch       PetscCall(ISEqual(il, isx[i], &equal));
5128b400f6SJacob Faibussowitsch       PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match", j, i);
529566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&il));
53c4762a1bSJed Brown     }
549566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "testfile", FILE_MODE_READ, &vl));
57c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
589566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, izero, PETSC_COPY_VALUES, &il));
599566063dSJacob Faibussowitsch       PetscCall(ISLoad(il, vl));
609566063dSJacob Faibussowitsch       PetscCall(ISEqual(il, isx[i], &equal));
6128b400f6SJacob Faibussowitsch       PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match", j, i);
629566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&il));
63c4762a1bSJed Brown     }
649566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
65c4762a1bSJed Brown   }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   for (j = 0; j < 3; j++) {
689566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "testfile_noheader", FILE_MODE_WRITE, &vx));
699566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinarySetSkipHeader(vx, PETSC_TRUE));
7048a46eb9SPierre Jolivet     for (i = 0; i < 3; i++) PetscCall(ISView(isx[i], vx));
719566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vx));
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "testfile_noheader", FILE_MODE_READ, &vl));
749566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinarySetSkipHeader(vl, PETSC_TRUE));
75c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
769566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, izero, PETSC_COPY_VALUES, &il));
779566063dSJacob Faibussowitsch       PetscCall(ISLoad(il, vl));
789566063dSJacob Faibussowitsch       PetscCall(ISEqual(il, isx[i], &equal));
7928b400f6SJacob Faibussowitsch       PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match", j, i);
809566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&il));
81c4762a1bSJed Brown     }
829566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown 
8548a46eb9SPierre Jolivet   for (i = 0; i < 3; i++) PetscCall(ISDestroy(&isx[i]));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
88c4762a1bSJed Brown     const char *filename  = (j == 0) ? "testfile_isstride" : "testfile_isblock";
89c4762a1bSJed Brown     PetscInt    blocksize = (j == 0) ? 1 : size;
909566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_WRITE, &vx));
91c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
92c4762a1bSJed Brown       if (j == 0) {
939566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, rank, rank + 1, &isx[i]));
94c4762a1bSJed Brown       } else {
959566063dSJacob Faibussowitsch         PetscCall(ISCreateBlock(PETSC_COMM_WORLD, blocksize, n, ix[i][rank], PETSC_COPY_VALUES, &isx[i]));
96c4762a1bSJed Brown       }
979566063dSJacob Faibussowitsch       PetscCall(ISView(isx[i], vx));
989566063dSJacob Faibussowitsch       PetscCall(ISToGeneral(isx[i]));
99c4762a1bSJed Brown     }
1009566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vx));
1019566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &vl));
102c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
1039566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, blocksize * n, izero, PETSC_COPY_VALUES, &il));
1049566063dSJacob Faibussowitsch       PetscCall(ISLoad(il, vl));
1059566063dSJacob Faibussowitsch       PetscCall(ISEqual(il, isx[i], &equal));
10628b400f6SJacob Faibussowitsch       PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match", j, i);
1079566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&il));
108c4762a1bSJed Brown     }
1099566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
11048a46eb9SPierre Jolivet     for (i = 0; i < 3; i++) PetscCall(ISDestroy(&isx[i]));
111c4762a1bSJed Brown   }
1129566063dSJacob Faibussowitsch   PetscCall(PetscFree(izero));
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
122*3886731fSPierre Jolivet       output_file: output/empty.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 
133c4762a1bSJed Brown    testset:
134c4762a1bSJed Brown       requires: mpiio
135c4762a1bSJed Brown       args: -viewer_binary_mpiio 1
136*3886731fSPierre Jolivet       output_file: output/empty.out
137c4762a1bSJed Brown       test:
138c4762a1bSJed Brown         suffix: mpiio_1
139c4762a1bSJed Brown         nsize: 1
140c4762a1bSJed Brown       test:
141c4762a1bSJed Brown         suffix: mpiio_2
142c4762a1bSJed Brown         nsize: 2
143c4762a1bSJed Brown       test:
144c4762a1bSJed Brown         suffix: mpiio_3
145c4762a1bSJed Brown         nsize: 3
146c4762a1bSJed Brown 
147c4762a1bSJed Brown TEST*/
148