xref: /petsc/src/vec/is/tests/ex2.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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;
19ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
20ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
21*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 3,PETSC_COMM_WORLD,PETSC_ERR_ARG_SIZ,"Example only works with up to three processes");
22c4762a1bSJed Brown 
232cf5aabcSBarry Smith   ierr = PetscCalloc1(size*n,&izero);CHKERRQ(ierr);
24c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
25c4762a1bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_WORLD,n,ix[i][rank],PETSC_COPY_VALUES,&isx[i]);CHKERRQ(ierr);
26c4762a1bSJed Brown   }
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   for (j = 0; j < 3; j++) {
29c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_WRITE,&vx);CHKERRQ(ierr);
30c4762a1bSJed Brown     ierr = ISView(isx[0],vx);CHKERRQ(ierr);
31c4762a1bSJed Brown     ierr = PetscViewerDestroy(&vx);CHKERRQ(ierr);
32c4762a1bSJed Brown 
33c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl);CHKERRQ(ierr);
34c4762a1bSJed Brown     ierr = ISCreate(PETSC_COMM_WORLD,&il);CHKERRQ(ierr);
35c4762a1bSJed Brown     ierr = ISLoad(il,vl);CHKERRQ(ierr);
36c4762a1bSJed Brown     ierr = ISEqual(il,isx[0],&equal);CHKERRQ(ierr);
37*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set loaded from file does not match",j);
38c4762a1bSJed Brown     ierr = ISDestroy(&il);CHKERRQ(ierr);
39c4762a1bSJed Brown     ierr = PetscViewerDestroy(&vl);CHKERRQ(ierr);
40c4762a1bSJed Brown 
41c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_APPEND,&vx);CHKERRQ(ierr);
42c4762a1bSJed Brown     ierr = ISView(isx[1],vx);CHKERRQ(ierr);
43c4762a1bSJed Brown     ierr = ISView(isx[2],vx);CHKERRQ(ierr);
44c4762a1bSJed Brown     ierr = PetscViewerDestroy(&vx);CHKERRQ(ierr);
45c4762a1bSJed Brown 
46c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl);CHKERRQ(ierr);
47c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
48c4762a1bSJed Brown       ierr = ISCreate(PETSC_COMM_WORLD,&il);CHKERRQ(ierr);
49c4762a1bSJed Brown       ierr = ISLoad(il,vl);CHKERRQ(ierr);
50c4762a1bSJed Brown       ierr = ISEqual(il,isx[i],&equal);CHKERRQ(ierr);
51*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
52c4762a1bSJed Brown       ierr = ISDestroy(&il);CHKERRQ(ierr);
53c4762a1bSJed Brown     }
54c4762a1bSJed Brown     ierr = PetscViewerDestroy(&vl);CHKERRQ(ierr);
55c4762a1bSJed Brown 
56c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile",FILE_MODE_READ,&vl);CHKERRQ(ierr);
57c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
58c4762a1bSJed Brown       ierr = ISCreateGeneral(PETSC_COMM_WORLD,n,izero,PETSC_COPY_VALUES,&il);CHKERRQ(ierr);
59c4762a1bSJed Brown       ierr = ISLoad(il,vl);CHKERRQ(ierr);
60c4762a1bSJed Brown       ierr = ISEqual(il,isx[i],&equal);CHKERRQ(ierr);
61*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
62c4762a1bSJed Brown       ierr = ISDestroy(&il);CHKERRQ(ierr);
63c4762a1bSJed Brown     }
64c4762a1bSJed Brown     ierr = PetscViewerDestroy(&vl);CHKERRQ(ierr);
65c4762a1bSJed Brown   }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   for (j = 0; j < 3; j++) {
68c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile_noheader",FILE_MODE_WRITE,&vx);CHKERRQ(ierr);
69c4762a1bSJed Brown     ierr = PetscViewerBinarySetSkipHeader(vx,PETSC_TRUE);CHKERRQ(ierr);
70c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
71c4762a1bSJed Brown       ierr = ISView(isx[i],vx);CHKERRQ(ierr);
72c4762a1bSJed Brown     }
73c4762a1bSJed Brown     ierr = PetscViewerDestroy(&vx);CHKERRQ(ierr);
74c4762a1bSJed Brown 
75c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"testfile_noheader",FILE_MODE_READ,&vl);CHKERRQ(ierr);
76c4762a1bSJed Brown     ierr = PetscViewerBinarySetSkipHeader(vl,PETSC_TRUE);CHKERRQ(ierr);
77c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
78c4762a1bSJed Brown       ierr = ISCreateGeneral(PETSC_COMM_WORLD,n,izero,PETSC_COPY_VALUES,&il);CHKERRQ(ierr);
79c4762a1bSJed Brown       ierr = ISLoad(il,vl);CHKERRQ(ierr);
80c4762a1bSJed Brown       ierr = ISEqual(il,isx[i],&equal);CHKERRQ(ierr);
81*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
82c4762a1bSJed Brown       ierr = ISDestroy(&il);CHKERRQ(ierr);
83c4762a1bSJed Brown     }
84c4762a1bSJed Brown     ierr = PetscViewerDestroy(&vl);CHKERRQ(ierr);
85c4762a1bSJed Brown   }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
88c4762a1bSJed Brown     ierr = ISDestroy(&isx[i]);CHKERRQ(ierr);
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;
94c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&vx);CHKERRQ(ierr);
95c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
96c4762a1bSJed Brown       if (j == 0) {
97c4762a1bSJed Brown         ierr = ISCreateStride(PETSC_COMM_WORLD,n,rank,rank+1,&isx[i]);CHKERRQ(ierr);
98c4762a1bSJed Brown       } else {
99c4762a1bSJed Brown         ierr = ISCreateBlock(PETSC_COMM_WORLD,blocksize,n,ix[i][rank],PETSC_COPY_VALUES,&isx[i]);CHKERRQ(ierr);
100c4762a1bSJed Brown       }
101c4762a1bSJed Brown       ierr = ISView(isx[i],vx);CHKERRQ(ierr);
102c4762a1bSJed Brown       ierr = ISToGeneral(isx[i]);CHKERRQ(ierr);
103c4762a1bSJed Brown     }
104c4762a1bSJed Brown     ierr = PetscViewerDestroy(&vx);CHKERRQ(ierr);
105c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&vl);CHKERRQ(ierr);
106c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
107c4762a1bSJed Brown       ierr = ISCreateGeneral(PETSC_COMM_WORLD,blocksize*n,izero,PETSC_COPY_VALUES,&il);CHKERRQ(ierr);
108c4762a1bSJed Brown       ierr = ISLoad(il,vl);CHKERRQ(ierr);
109c4762a1bSJed Brown       ierr = ISEqual(il,isx[i],&equal);CHKERRQ(ierr);
110*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Iteration %" PetscInt_FMT " - Index set %" PetscInt_FMT " loaded from file does not match",j,i);
111c4762a1bSJed Brown       ierr = ISDestroy(&il);CHKERRQ(ierr);
112c4762a1bSJed Brown     }
113c4762a1bSJed Brown     ierr = PetscViewerDestroy(&vl);CHKERRQ(ierr);
114c4762a1bSJed Brown     for (i = 0; i < 3; i++) {
115c4762a1bSJed Brown       ierr = ISDestroy(&isx[i]);CHKERRQ(ierr);
116c4762a1bSJed Brown     }
117c4762a1bSJed Brown   }
1182cf5aabcSBarry Smith   ierr = PetscFree(izero);CHKERRQ(ierr);
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