xref: /petsc/src/mat/tests/ex41.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap() - the parallel case. This example\n\
3c4762a1bSJed Brown is similar to ex40.c; here the index sets used are random. Input arguments are:\n\
4c4762a1bSJed Brown   -f <input_file> : file to load.  For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\
5c4762a1bSJed Brown   -nd <size>      : > 0  no of domains per processor \n\
6c4762a1bSJed Brown   -ov <overlap>   : >=0  amount of overlap between domains\n\n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petscmat.h>
9c4762a1bSJed Brown 
10c4762a1bSJed Brown int main(int argc,char **args)
11c4762a1bSJed Brown {
12c4762a1bSJed Brown   PetscInt       nd = 2,ov=1,i,j,m,n,*idx,lsize;
13c4762a1bSJed Brown   PetscErrorCode ierr;
14c4762a1bSJed Brown   PetscMPIInt    rank;
15c4762a1bSJed Brown   PetscBool      flg;
16c4762a1bSJed Brown   Mat            A,B;
17c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
18c4762a1bSJed Brown   PetscViewer    fd;
19c4762a1bSJed Brown   IS             *is1,*is2;
20c4762a1bSJed Brown   PetscRandom    r;
21c4762a1bSJed Brown   PetscScalar    rand;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
245f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* Read matrix and RHS */
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATMPIAIJ));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /* Read the matrix again as a seq matrix */
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&B));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATSEQAIJ));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(B,fd));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Create the Random no generator */
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,&m,&n));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&r));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(r));
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   /* Create the IS corresponding to subdomains */
495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nd,&is1));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nd,&is2));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(m ,&idx));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Create the random Index Sets */
54c4762a1bSJed Brown   for (i=0; i<nd; i++) {
55c4762a1bSJed Brown     for (j=0; j<rank; j++) {
565f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(r,&rand));
57c4762a1bSJed Brown     }
585f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(r,&rand));
59c4762a1bSJed Brown     lsize = (PetscInt)(rand*m);
60c4762a1bSJed Brown     for (j=0; j<lsize; j++) {
615f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(r,&rand));
62c4762a1bSJed Brown       idx[j] = (PetscInt)(rand*m);
63c4762a1bSJed Brown     }
645f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,lsize,idx,PETSC_COPY_VALUES,is1+i));
655f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,lsize,idx,PETSC_COPY_VALUES,is2+i));
66c4762a1bSJed Brown   }
67c4762a1bSJed Brown 
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(A,nd,is1,ov));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIncreaseOverlap(B,nd,is2,ov));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* Now see if the serial and parallel case have the same answers */
72c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
73c4762a1bSJed Brown     PetscInt sz1,sz2;
745f80ce2aSJacob Faibussowitsch     CHKERRQ(ISEqual(is1[i],is2[i],&flg));
755f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetSize(is1[i],&sz1));
765f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetSize(is2[i],&sz2));
77*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"proc:[%d], i=%" PetscInt_FMT ", flg =%d  sz1 = %" PetscInt_FMT " sz2 = %" PetscInt_FMT,rank,i,(int)flg,sz1,sz2);
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Free Allocated Memory */
81c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
825f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is1[i]));
835f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is2[i]));
84c4762a1bSJed Brown   }
855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&r));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is1));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is2));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(idx));
91c4762a1bSJed Brown   ierr = PetscFinalize();
92c4762a1bSJed Brown   return ierr;
93c4762a1bSJed Brown }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown /*TEST
96c4762a1bSJed Brown 
97c4762a1bSJed Brown    build:
98c4762a1bSJed Brown       requires: !complex
99c4762a1bSJed Brown 
100c4762a1bSJed Brown    test:
101c4762a1bSJed Brown       nsize: 3
102dfd57a17SPierre Jolivet       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
103c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1 -nd 3 -ov 1
104c4762a1bSJed Brown 
105c4762a1bSJed Brown TEST*/
106