xref: /petsc/src/vec/is/tests/ex2.c (revision be096a4675ce3b20dd7f9c195e6046e69d9955cf)
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 
179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20*be096a46SBarry Smith   PetscCheck(size < 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Example only works with up to three processes");
21c4762a1bSJed Brown 
229566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(size*n,&izero));
23c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
249566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n,ix[i][rank],PETSC_COPY_VALUES,&isx[i]));
25c4762a1bSJed Brown   }
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   for (j = 0; j < 3; j++) {
289566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_WRITE,&vx));
299566063dSJacob Faibussowitsch     PetscCall(ISView(isx[0],vx));
309566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vx));
31c4762a1bSJed Brown 
329566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl));
339566063dSJacob Faibussowitsch     PetscCall(ISCreate(PETSC_COMM_WORLD,&il));
349566063dSJacob Faibussowitsch     PetscCall(ISLoad(il,vl));
359566063dSJacob Faibussowitsch     PetscCall(ISEqual(il,isx[0],&equal));
3628b400f6SJacob Faibussowitsch     PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set loaded from file does not match",j);
379566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&il));
389566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
39c4762a1bSJed Brown 
409566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_APPEND,&vx));
419566063dSJacob Faibussowitsch     PetscCall(ISView(isx[1],vx));
429566063dSJacob Faibussowitsch     PetscCall(ISView(isx[2],vx));
439566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vx));
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl));
46c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
479566063dSJacob Faibussowitsch       PetscCall(ISCreate(PETSC_COMM_WORLD,&il));
489566063dSJacob Faibussowitsch       PetscCall(ISLoad(il,vl));
499566063dSJacob Faibussowitsch       PetscCall(ISEqual(il,isx[i],&equal));
5028b400f6SJacob Faibussowitsch       PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
519566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&il));
52c4762a1bSJed Brown     }
539566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
54c4762a1bSJed Brown 
559566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl));
56c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
579566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n,izero,PETSC_COPY_VALUES,&il));
589566063dSJacob Faibussowitsch       PetscCall(ISLoad(il,vl));
599566063dSJacob Faibussowitsch       PetscCall(ISEqual(il,isx[i],&equal));
6028b400f6SJacob Faibussowitsch       PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
619566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&il));
62c4762a1bSJed Brown     }
639566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
64c4762a1bSJed Brown   }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   for (j = 0; j < 3; j++) {
679566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile_noheader",FILE_MODE_WRITE,&vx));
689566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinarySetSkipHeader(vx,PETSC_TRUE));
69c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
709566063dSJacob Faibussowitsch       PetscCall(ISView(isx[i],vx));
71c4762a1bSJed Brown     }
729566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vx));
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile_noheader",FILE_MODE_READ,&vl));
759566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinarySetSkipHeader(vl,PETSC_TRUE));
76c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
779566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n,izero,PETSC_COPY_VALUES,&il));
789566063dSJacob Faibussowitsch       PetscCall(ISLoad(il,vl));
799566063dSJacob Faibussowitsch       PetscCall(ISEqual(il,isx[i],&equal));
8028b400f6SJacob Faibussowitsch       PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
819566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&il));
82c4762a1bSJed Brown     }
839566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
84c4762a1bSJed Brown   }
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
879566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isx[i]));
88c4762a1bSJed Brown   }
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
91c4762a1bSJed Brown     const char *filename  = (j == 0) ? "testfile_isstride" : "testfile_isblock";
92c4762a1bSJed Brown     PetscInt    blocksize = (j == 0) ? 1 : size;
939566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&vx));
94c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
95c4762a1bSJed Brown       if (j == 0) {
969566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,rank,rank+1,&isx[i]));
97c4762a1bSJed Brown       } else {
989566063dSJacob Faibussowitsch         PetscCall(ISCreateBlock(PETSC_COMM_WORLD,blocksize,n,ix[i][rank],PETSC_COPY_VALUES,&isx[i]));
99c4762a1bSJed Brown       }
1009566063dSJacob Faibussowitsch       PetscCall(ISView(isx[i],vx));
1019566063dSJacob Faibussowitsch       PetscCall(ISToGeneral(isx[i]));
102c4762a1bSJed Brown     }
1039566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vx));
1049566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&vl));
105c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
1069566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,blocksize*n,izero,PETSC_COPY_VALUES,&il));
1079566063dSJacob Faibussowitsch       PetscCall(ISLoad(il,vl));
1089566063dSJacob Faibussowitsch       PetscCall(ISEqual(il,isx[i],&equal));
10928b400f6SJacob Faibussowitsch       PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
1109566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&il));
111c4762a1bSJed Brown     }
1129566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vl));
113c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
1149566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&isx[i]));
115c4762a1bSJed Brown     }
116c4762a1bSJed Brown   }
1179566063dSJacob Faibussowitsch   PetscCall(PetscFree(izero));
118c4762a1bSJed Brown 
1199566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
120b122ec5aSJacob Faibussowitsch   return 0;
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown /*TEST
124c4762a1bSJed Brown 
125c4762a1bSJed Brown    testset:
126c4762a1bSJed Brown       args: -viewer_binary_mpiio 0
127c4762a1bSJed Brown       output_file: output/ex2_1.out
128c4762a1bSJed Brown       test:
129c4762a1bSJed Brown         suffix: stdio_1
130c4762a1bSJed Brown         nsize: 1
131c4762a1bSJed Brown       test:
132c4762a1bSJed Brown         suffix: stdio_2
133c4762a1bSJed Brown         nsize: 2
134c4762a1bSJed Brown       test:
135c4762a1bSJed Brown         suffix: stdio_3
136c4762a1bSJed Brown         nsize: 3
137c4762a1bSJed Brown 
138c4762a1bSJed Brown    testset:
139c4762a1bSJed Brown       requires: mpiio
140c4762a1bSJed Brown       args: -viewer_binary_mpiio 1
141c4762a1bSJed Brown       output_file: output/ex2_1.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 
152c4762a1bSJed Brown TEST*/
153