xref: /petsc/src/mat/tests/ex40.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests the parallel case for MatIncreaseOverlap(). Input arguments are:\n\
3c4762a1bSJed Brown   -f <input_file> : file to load.  For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\
4c4762a1bSJed Brown   -nd <size>      : > 0  number of domains per processor \n\
5c4762a1bSJed Brown   -ov <overlap>   : >=0  amount of overlap between domains\n\n";
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #include <petscmat.h>
8c4762a1bSJed Brown 
9c4762a1bSJed Brown PetscErrorCode ISAllGatherDisjoint(IS iis, IS** ois)
10c4762a1bSJed Brown {
11c4762a1bSJed Brown   IS             *is2,is;
12c4762a1bSJed Brown   const PetscInt *idxs;
13c4762a1bSJed Brown   PetscInt       i, ls,*sizes;
14c4762a1bSJed Brown   PetscMPIInt    size;
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   PetscFunctionBeginUser;
175f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)iis),&size));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size,&is2));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size,&sizes));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(iis,&ls));
21c4762a1bSJed Brown   /* we don't have a public ISGetLayout */
225f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allgather(&ls,1,MPIU_INT,sizes,1,MPIU_INT,PetscObjectComm((PetscObject)iis)));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(ISAllGather(iis,&is));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(is,&idxs));
25c4762a1bSJed Brown   for (i = 0, ls = 0; i < size; i++) {
265f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,sizes[i],idxs+ls,PETSC_COPY_VALUES,&is2[i]));
27c4762a1bSJed Brown     ls += sizes[i];
28c4762a1bSJed Brown   }
295f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(is,&idxs));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(sizes));
32c4762a1bSJed Brown   *ois = is2;
33c4762a1bSJed Brown   PetscFunctionReturn(0);
34c4762a1bSJed Brown }
35c4762a1bSJed Brown 
36c4762a1bSJed Brown int main(int argc,char **args)
37c4762a1bSJed Brown {
38c4762a1bSJed Brown   PetscInt       nd = 2,ov = 1,ndpar,i,start,m,n,end,lsize;
39c4762a1bSJed Brown   PetscMPIInt    rank;
40c4762a1bSJed Brown   PetscBool      flg, useND = PETSC_FALSE;
41c4762a1bSJed Brown   Mat            A,B;
42c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
43c4762a1bSJed Brown   PetscViewer    fd;
44c4762a1bSJed Brown   IS             *is1,*is2;
45c4762a1bSJed Brown   PetscRandom    r;
46c4762a1bSJed Brown   PetscScalar    rand;
47c4762a1bSJed Brown 
48*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
495f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
5128b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must use -f filename to indicate a file containing a PETSc binary matrix");
525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-nested_dissection",&useND,NULL));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Read matrix */
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATMPIAIJ));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Read the matrix again as a sequential matrix */
655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&B));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATSEQAIJ));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(B,fd));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(B));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Create the IS corresponding to subdomains */
73c4762a1bSJed Brown   if (useND) {
74c4762a1bSJed Brown     MatPartitioning part;
75c4762a1bSJed Brown     IS              ndmap;
76c4762a1bSJed Brown     PetscMPIInt     size;
77c4762a1bSJed Brown 
78c4762a1bSJed Brown     ndpar = 1;
795f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
80c4762a1bSJed Brown     nd   = (PetscInt)size;
815f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(ndpar,&is1));
825f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningCreate(PETSC_COMM_WORLD,&part));
835f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningSetAdjacency(part,A));
845f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningSetFromOptions(part));
855f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningApplyND(part,&ndmap));
865f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPartitioningDestroy(&part));
875f80ce2aSJacob Faibussowitsch     CHKERRQ(ISBuildTwoSided(ndmap,NULL,&is1[0]));
885f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&ndmap));
895f80ce2aSJacob Faibussowitsch     CHKERRQ(ISAllGatherDisjoint(is1[0],&is2));
90c4762a1bSJed Brown   } else {
91c4762a1bSJed Brown     /* Create the random Index Sets */
925f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nd,&is1));
935f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nd,&is2));
94c4762a1bSJed Brown 
955f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetSize(A,&m,&n));
965f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&r));
975f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomSetFromOptions(r));
98c4762a1bSJed Brown     for (i=0; i<nd; i++) {
995f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(r,&rand));
100c4762a1bSJed Brown       start = (PetscInt)(rand*m);
1015f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(r,&rand));
102c4762a1bSJed Brown       end   = (PetscInt)(rand*m);
103c4762a1bSJed Brown       lsize =  end - start;
104c4762a1bSJed Brown       if (start > end) { start = end; lsize = -lsize;}
1055f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateStride(PETSC_COMM_SELF,lsize,start,1,is1+i));
1065f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateStride(PETSC_COMM_SELF,lsize,start,1,is2+i));
107c4762a1bSJed Brown     }
108c4762a1bSJed Brown     ndpar = nd;
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomDestroy(&r));
110c4762a1bSJed Brown   }
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(A,ndpar,is1,ov));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(B,nd,is2,ov));
113c4762a1bSJed Brown   if (useND) {
114c4762a1bSJed Brown     IS *is;
115c4762a1bSJed Brown 
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(ISAllGatherDisjoint(is1[0],&is));
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is1[0]));
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(is1));
119c4762a1bSJed Brown     is1 = is;
120c4762a1bSJed Brown   }
121c4762a1bSJed Brown   /* Now see if the serial and parallel case have the same answers */
122c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1235f80ce2aSJacob Faibussowitsch     CHKERRQ(ISEqual(is1[i],is2[i],&flg));
124c4762a1bSJed Brown     if (!flg) {
1255f80ce2aSJacob Faibussowitsch       CHKERRQ(ISViewFromOptions(is1[i],NULL,"-err_view"));
1265f80ce2aSJacob Faibussowitsch       CHKERRQ(ISViewFromOptions(is2[i],NULL,"-err_view"));
12798921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"proc:[%d], i=%" PetscInt_FMT ", flg =%d",rank,i,(int)flg);
128c4762a1bSJed Brown     }
129c4762a1bSJed Brown   }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /* Free allocated memory */
132c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is1[i]));
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is2[i]));
135c4762a1bSJed Brown   }
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is1));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is2));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
140*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
141*b122ec5aSJacob Faibussowitsch   return 0;
142c4762a1bSJed Brown }
143c4762a1bSJed Brown 
144c4762a1bSJed Brown /*TEST
145c4762a1bSJed Brown 
146c4762a1bSJed Brown    build:
147c4762a1bSJed Brown       requires: !complex
148c4762a1bSJed Brown 
149c4762a1bSJed Brown    testset:
150c4762a1bSJed Brown       nsize: 5
151dfd57a17SPierre Jolivet       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
152c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info -ov 2
153c4762a1bSJed Brown       output_file: output/ex40_1.out
154c4762a1bSJed Brown       test:
155c4762a1bSJed Brown         suffix: 1
156c4762a1bSJed Brown         args: -nd 7
157c4762a1bSJed Brown       test:
158c4762a1bSJed Brown         requires: parmetis
159c4762a1bSJed Brown         suffix: 1_nd
160c4762a1bSJed Brown         args: -nested_dissection -mat_partitioning_type parmetis
161c4762a1bSJed Brown 
162c4762a1bSJed Brown    testset:
163c4762a1bSJed Brown       nsize: 3
164dfd57a17SPierre Jolivet       requires: double !defined(PETSC_USE_64BIT_INDICES) !complex
165c4762a1bSJed Brown       args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 -mat_increase_overlap_scalable 1 -ov 2
166c4762a1bSJed Brown       output_file: output/ex40_1.out
167c4762a1bSJed Brown       test:
168c4762a1bSJed Brown         suffix: 2
169c4762a1bSJed Brown         args: -nd 7
170c4762a1bSJed Brown       test:
171c4762a1bSJed Brown         requires: parmetis
172c4762a1bSJed Brown         suffix: 2_nd
173c4762a1bSJed Brown         args: -nested_dissection -mat_partitioning_type parmetis
174c4762a1bSJed Brown 
175c4762a1bSJed Brown TEST*/
176