xref: /petsc/src/vec/is/tests/ex2.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 {
8c4762a1bSJed Brown   PetscErrorCode         ierr;
92cf5aabcSBarry Smith   PetscInt               n = 3, *izero, j, i;
10c4762a1bSJed Brown   PetscInt               ix[3][3][3] = {{{3,5,4},{1,7,9},{0,2,8}},
11c4762a1bSJed Brown                                         {{0,2,8},{3,5,4},{1,7,9}},
12c4762a1bSJed Brown                                         {{1,7,9},{0,2,8},{3,5,4}}};
13c4762a1bSJed Brown   IS                     isx[3],il;
14c4762a1bSJed Brown   PetscMPIInt            size,rank;
15c4762a1bSJed Brown   PetscViewer            vx,vl;
16c4762a1bSJed Brown   PetscBool              equal;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
19*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
20*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
212c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 3,PETSC_COMM_WORLD,PETSC_ERR_ARG_SIZ,"Example only works with up to three processes");
22c4762a1bSJed Brown 
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(size*n,&izero));
24c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
25*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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++) {
29*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_WRITE,&vx));
30*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(isx[0],vx));
31*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&vx));
32c4762a1bSJed Brown 
33*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl));
34*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreate(PETSC_COMM_WORLD,&il));
35*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLoad(il,vl));
36*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISEqual(il,isx[0],&equal));
372c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set loaded from file does not match",j);
38*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&il));
39*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&vl));
40c4762a1bSJed Brown 
41*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_APPEND,&vx));
42*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(isx[1],vx));
43*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(isx[2],vx));
44*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&vx));
45c4762a1bSJed Brown 
46*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl));
47c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
48*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreate(PETSC_COMM_WORLD,&il));
49*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISLoad(il,vl));
50*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISEqual(il,isx[i],&equal));
512c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
52*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&il));
53c4762a1bSJed Brown     }
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&vl));
55c4762a1bSJed Brown 
56*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl));
57c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
58*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,n,izero,PETSC_COPY_VALUES,&il));
59*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISLoad(il,vl));
60*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISEqual(il,isx[i],&equal));
612c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
62*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&il));
63c4762a1bSJed Brown     }
64*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&vl));
65c4762a1bSJed Brown   }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   for (j = 0; j < 3; j++) {
68*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile_noheader",FILE_MODE_WRITE,&vx));
69*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinarySetSkipHeader(vx,PETSC_TRUE));
70c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
71*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISView(isx[i],vx));
72c4762a1bSJed Brown     }
73*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&vx));
74c4762a1bSJed Brown 
75*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile_noheader",FILE_MODE_READ,&vl));
76*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinarySetSkipHeader(vl,PETSC_TRUE));
77c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
78*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,n,izero,PETSC_COPY_VALUES,&il));
79*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISLoad(il,vl));
80*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISEqual(il,isx[i],&equal));
812c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
82*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&il));
83c4762a1bSJed Brown     }
84*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&vl));
85c4762a1bSJed Brown   }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
88*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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;
94*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&vx));
95c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
96c4762a1bSJed Brown       if (j == 0) {
97*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,n,rank,rank+1,&isx[i]));
98c4762a1bSJed Brown       } else {
99*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISCreateBlock(PETSC_COMM_WORLD,blocksize,n,ix[i][rank],PETSC_COPY_VALUES,&isx[i]));
100c4762a1bSJed Brown       }
101*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISView(isx[i],vx));
102*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISToGeneral(isx[i]));
103c4762a1bSJed Brown     }
104*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&vx));
105*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&vl));
106c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
107*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,blocksize*n,izero,PETSC_COPY_VALUES,&il));
108*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISLoad(il,vl));
109*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISEqual(il,isx[i],&equal));
1102c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
111*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&il));
112c4762a1bSJed Brown     }
113*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&vl));
114c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
115*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&isx[i]));
116c4762a1bSJed Brown     }
117c4762a1bSJed Brown   }
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(izero));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   ierr = PetscFinalize();
121c4762a1bSJed Brown   return ierr;
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