xref: /petsc/src/vec/is/tests/ex2.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
6c4762a1bSJed Brown int main(int argc,char **argv)
7c4762a1bSJed Brown {
82cf5aabcSBarry Smith   PetscInt               n = 3, *izero, j, i;
9c4762a1bSJed Brown   PetscInt               ix[3][3][3] = {{{3,5,4},{1,7,9},{0,2,8}},
10c4762a1bSJed Brown                                         {{0,2,8},{3,5,4},{1,7,9}},
11c4762a1bSJed Brown                                         {{1,7,9},{0,2,8},{3,5,4}}};
12c4762a1bSJed Brown   IS                     isx[3],il;
13c4762a1bSJed Brown   PetscMPIInt            size,rank;
14c4762a1bSJed Brown   PetscViewer            vx,vl;
15c4762a1bSJed Brown   PetscBool              equal;
16c4762a1bSJed Brown 
17*327415f7SBarry Smith   PetscFunctionBeginUser;
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
21be096a46SBarry Smith   PetscCheck(size < 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Example only works with up to three processes");
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(size*n,&izero));
24c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
259566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n,ix[i][rank],PETSC_COPY_VALUES,&isx[i]));
26c4762a1bSJed Brown   }
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));
70c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
719566063dSJacob Faibussowitsch       PetscCall(ISView(isx[i],vx));
72c4762a1bSJed Brown     }
739566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vx));
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile_noheader",FILE_MODE_READ,&vl));
769566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinarySetSkipHeader(vl,PETSC_TRUE));
77c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
789566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n,izero,PETSC_COPY_VALUES,&il));
799566063dSJacob Faibussowitsch       PetscCall(ISLoad(il,vl));
809566063dSJacob Faibussowitsch       PetscCall(ISEqual(il,isx[i],&equal));
8128b400f6SJacob Faibussowitsch       PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
829566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&il));
83c4762a1bSJed Brown     }
849566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
85c4762a1bSJed Brown   }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
889566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isx[i]));
89c4762a1bSJed Brown   }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
92c4762a1bSJed Brown     const char *filename  = (j == 0) ? "testfile_isstride" : "testfile_isblock";
93c4762a1bSJed Brown     PetscInt    blocksize = (j == 0) ? 1 : size;
949566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&vx));
95c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
96c4762a1bSJed Brown       if (j == 0) {
979566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,rank,rank+1,&isx[i]));
98c4762a1bSJed Brown       } else {
999566063dSJacob Faibussowitsch         PetscCall(ISCreateBlock(PETSC_COMM_WORLD,blocksize,n,ix[i][rank],PETSC_COPY_VALUES,&isx[i]));
100c4762a1bSJed Brown       }
1019566063dSJacob Faibussowitsch       PetscCall(ISView(isx[i],vx));
1029566063dSJacob Faibussowitsch       PetscCall(ISToGeneral(isx[i]));
103c4762a1bSJed Brown     }
1049566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vx));
1059566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&vl));
106c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
1079566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,blocksize*n,izero,PETSC_COPY_VALUES,&il));
1089566063dSJacob Faibussowitsch       PetscCall(ISLoad(il,vl));
1099566063dSJacob Faibussowitsch       PetscCall(ISEqual(il,isx[i],&equal));
11028b400f6SJacob Faibussowitsch       PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
1119566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&il));
112c4762a1bSJed Brown     }
1139566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
114c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
1159566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&isx[i]));
116c4762a1bSJed Brown     }
117c4762a1bSJed Brown   }
1189566063dSJacob Faibussowitsch   PetscCall(PetscFree(izero));
119c4762a1bSJed Brown 
1209566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
121b122ec5aSJacob Faibussowitsch   return 0;
122c4762a1bSJed Brown }
123c4762a1bSJed Brown 
124c4762a1bSJed Brown /*TEST
125c4762a1bSJed Brown 
126c4762a1bSJed Brown    testset:
127c4762a1bSJed Brown       args: -viewer_binary_mpiio 0
128c4762a1bSJed Brown       output_file: output/ex2_1.out
129c4762a1bSJed Brown       test:
130c4762a1bSJed Brown         suffix: stdio_1
131c4762a1bSJed Brown         nsize: 1
132c4762a1bSJed Brown       test:
133c4762a1bSJed Brown         suffix: stdio_2
134c4762a1bSJed Brown         nsize: 2
135c4762a1bSJed Brown       test:
136c4762a1bSJed Brown         suffix: stdio_3
137c4762a1bSJed Brown         nsize: 3
138c4762a1bSJed Brown 
139c4762a1bSJed Brown    testset:
140c4762a1bSJed Brown       requires: mpiio
141c4762a1bSJed Brown       args: -viewer_binary_mpiio 1
142c4762a1bSJed Brown       output_file: output/ex2_1.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 
153c4762a1bSJed Brown TEST*/
154