xref: /petsc/src/mat/tests/ex40.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
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   PetscErrorCode ierr;
15c4762a1bSJed Brown   PetscMPIInt    size;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   PetscFunctionBeginUser;
18ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)iis),&size);CHKERRMPI(ierr);
19c4762a1bSJed Brown   ierr = PetscMalloc1(size,&is2);CHKERRQ(ierr);
20c4762a1bSJed Brown   ierr = PetscMalloc1(size,&sizes);CHKERRQ(ierr);
21c4762a1bSJed Brown   ierr = ISGetLocalSize(iis,&ls);CHKERRQ(ierr);
22c4762a1bSJed Brown   /* we don't have a public ISGetLayout */
23ffc4695bSBarry Smith   ierr = MPI_Allgather(&ls,1,MPIU_INT,sizes,1,MPIU_INT,PetscObjectComm((PetscObject)iis));CHKERRMPI(ierr);
24c4762a1bSJed Brown   ierr = ISAllGather(iis,&is);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
26c4762a1bSJed Brown   for (i = 0, ls = 0; i < size; i++) {
27c4762a1bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_SELF,sizes[i],idxs+ls,PETSC_COPY_VALUES,&is2[i]);CHKERRQ(ierr);
28c4762a1bSJed Brown     ls += sizes[i];
29c4762a1bSJed Brown   }
30c4762a1bSJed Brown   ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = ISDestroy(&is);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = PetscFree(sizes);CHKERRQ(ierr);
33c4762a1bSJed Brown   *ois = is2;
34c4762a1bSJed Brown   PetscFunctionReturn(0);
35c4762a1bSJed Brown }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown int main(int argc,char **args)
38c4762a1bSJed Brown {
39c4762a1bSJed Brown   PetscErrorCode ierr;
40c4762a1bSJed Brown   PetscInt       nd = 2,ov = 1,ndpar,i,start,m,n,end,lsize;
41c4762a1bSJed Brown   PetscMPIInt    rank;
42c4762a1bSJed Brown   PetscBool      flg, useND = PETSC_FALSE;
43c4762a1bSJed Brown   Mat            A,B;
44c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
45c4762a1bSJed Brown   PetscViewer    fd;
46c4762a1bSJed Brown   IS             *is1,*is2;
47c4762a1bSJed Brown   PetscRandom    r;
48c4762a1bSJed Brown   PetscScalar    rand;
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
51ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
52589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr);
53c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must use -f filename to indicate a file containing a PETSc binary matrix");
54c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL);CHKERRQ(ierr);
55c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL);CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-nested_dissection",&useND,NULL);CHKERRQ(ierr);
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* Read matrix */
59c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
61c4762a1bSJed Brown   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr = MatLoad(A,fd);CHKERRQ(ierr);
63c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
64c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* Read the matrix again as a sequential matrix */
67c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr);
70c4762a1bSJed Brown   ierr = MatLoad(B,fd);CHKERRQ(ierr);
71c4762a1bSJed Brown   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
72c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* Create the IS corresponding to subdomains */
75c4762a1bSJed Brown   if (useND) {
76c4762a1bSJed Brown     MatPartitioning part;
77c4762a1bSJed Brown     IS              ndmap;
78c4762a1bSJed Brown     PetscMPIInt     size;
79c4762a1bSJed Brown 
80c4762a1bSJed Brown     ndpar = 1;
81ffc4695bSBarry Smith     ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
82c4762a1bSJed Brown     nd   = (PetscInt)size;
83c4762a1bSJed Brown     ierr = PetscMalloc1(ndpar,&is1);CHKERRQ(ierr);
84c4762a1bSJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_WORLD,&part);CHKERRQ(ierr);
85c4762a1bSJed Brown     ierr = MatPartitioningSetAdjacency(part,A);CHKERRQ(ierr);
86c4762a1bSJed Brown     ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr);
87c4762a1bSJed Brown     ierr = MatPartitioningApplyND(part,&ndmap);CHKERRQ(ierr);
88c4762a1bSJed Brown     ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr);
89c4762a1bSJed Brown     ierr = ISBuildTwoSided(ndmap,NULL,&is1[0]);CHKERRQ(ierr);
90c4762a1bSJed Brown     ierr = ISDestroy(&ndmap);CHKERRQ(ierr);
91c4762a1bSJed Brown     ierr = ISAllGatherDisjoint(is1[0],&is2);CHKERRQ(ierr);
92c4762a1bSJed Brown   } else {
93c4762a1bSJed Brown     /* Create the random Index Sets */
94c4762a1bSJed Brown     ierr = PetscMalloc1(nd,&is1);CHKERRQ(ierr);
95c4762a1bSJed Brown     ierr = PetscMalloc1(nd,&is2);CHKERRQ(ierr);
96c4762a1bSJed Brown 
97c4762a1bSJed Brown     ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
98c4762a1bSJed Brown     ierr = PetscRandomCreate(PETSC_COMM_SELF,&r);CHKERRQ(ierr);
99c4762a1bSJed Brown     ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr);
100c4762a1bSJed Brown     for (i=0; i<nd; i++) {
101c4762a1bSJed Brown       ierr  = PetscRandomGetValue(r,&rand);CHKERRQ(ierr);
102c4762a1bSJed Brown       start = (PetscInt)(rand*m);
103c4762a1bSJed Brown       ierr  = PetscRandomGetValue(r,&rand);CHKERRQ(ierr);
104c4762a1bSJed Brown       end   = (PetscInt)(rand*m);
105c4762a1bSJed Brown       lsize =  end - start;
106c4762a1bSJed Brown       if (start > end) { start = end; lsize = -lsize;}
107c4762a1bSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,lsize,start,1,is1+i);CHKERRQ(ierr);
108c4762a1bSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,lsize,start,1,is2+i);CHKERRQ(ierr);
109c4762a1bSJed Brown     }
110c4762a1bSJed Brown     ndpar = nd;
111c4762a1bSJed Brown     ierr = PetscRandomDestroy(&r);CHKERRQ(ierr);
112c4762a1bSJed Brown   }
113c4762a1bSJed Brown   ierr = MatIncreaseOverlap(A,ndpar,is1,ov);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = MatIncreaseOverlap(B,nd,is2,ov);CHKERRQ(ierr);
115c4762a1bSJed Brown   if (useND) {
116c4762a1bSJed Brown     IS *is;
117c4762a1bSJed Brown 
118c4762a1bSJed Brown     ierr = ISAllGatherDisjoint(is1[0],&is);CHKERRQ(ierr);
119c4762a1bSJed Brown     ierr = ISDestroy(&is1[0]);CHKERRQ(ierr);
120c4762a1bSJed Brown     ierr = PetscFree(is1);CHKERRQ(ierr);
121c4762a1bSJed Brown     is1 = is;
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown   /* Now see if the serial and parallel case have the same answers */
124c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
125c4762a1bSJed Brown     ierr = ISEqual(is1[i],is2[i],&flg);CHKERRQ(ierr);
126c4762a1bSJed Brown     if (!flg) {
127c4762a1bSJed Brown       ierr = ISViewFromOptions(is1[i],NULL,"-err_view");CHKERRQ(ierr);
128c4762a1bSJed Brown       ierr = ISViewFromOptions(is2[i],NULL,"-err_view");CHKERRQ(ierr);
129c4762a1bSJed Brown       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"proc:[%d], i=%D, flg =%d\n",rank,i,(int)flg);
130c4762a1bSJed Brown     }
131c4762a1bSJed Brown   }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* Free allocated memory */
134c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
135c4762a1bSJed Brown     ierr = ISDestroy(&is1[i]);CHKERRQ(ierr);
136c4762a1bSJed Brown     ierr = ISDestroy(&is2[i]);CHKERRQ(ierr);
137c4762a1bSJed Brown   }
138c4762a1bSJed Brown   ierr = PetscFree(is1);CHKERRQ(ierr);
139c4762a1bSJed Brown   ierr = PetscFree(is2);CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
142c4762a1bSJed Brown   ierr = PetscFinalize();
143c4762a1bSJed Brown   return ierr;
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown /*TEST
147c4762a1bSJed Brown 
148c4762a1bSJed Brown    build:
149c4762a1bSJed Brown       requires: !complex
150c4762a1bSJed Brown 
151c4762a1bSJed Brown    testset:
152c4762a1bSJed Brown       nsize: 5
153*dfd57a17SPierre Jolivet       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
154c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info -ov 2
155c4762a1bSJed Brown       output_file: output/ex40_1.out
156c4762a1bSJed Brown       test:
157c4762a1bSJed Brown         suffix: 1
158c4762a1bSJed Brown         args: -nd 7
159c4762a1bSJed Brown       test:
160c4762a1bSJed Brown         requires: parmetis
161c4762a1bSJed Brown         suffix: 1_nd
162c4762a1bSJed Brown         args: -nested_dissection -mat_partitioning_type parmetis
163c4762a1bSJed Brown 
164c4762a1bSJed Brown    testset:
165c4762a1bSJed Brown       nsize: 3
166*dfd57a17SPierre Jolivet       requires: double !defined(PETSC_USE_64BIT_INDICES) !complex
167c4762a1bSJed Brown       args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 -mat_increase_overlap_scalable 1 -ov 2
168c4762a1bSJed Brown       output_file: output/ex40_1.out
169c4762a1bSJed Brown       test:
170c4762a1bSJed Brown         suffix: 2
171c4762a1bSJed Brown         args: -nd 7
172c4762a1bSJed Brown       test:
173c4762a1bSJed Brown         requires: parmetis
174c4762a1bSJed Brown         suffix: 2_nd
175c4762a1bSJed Brown         args: -nested_dissection -mat_partitioning_type parmetis
176c4762a1bSJed Brown 
177c4762a1bSJed Brown TEST*/
178