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 PetscMPIInt rank; 14c4762a1bSJed Brown PetscBool flg; 15c4762a1bSJed Brown Mat A,B; 16c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 17c4762a1bSJed Brown PetscViewer fd; 18c4762a1bSJed Brown IS *is1,*is2; 19c4762a1bSJed Brown PetscRandom r; 20c4762a1bSJed Brown PetscScalar rand; 21c4762a1bSJed Brown 22*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 23*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 24*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL)); 25*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL)); 26*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* Read matrix and RHS */ 29*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 30*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 31*9566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATMPIAIJ)); 32*9566063dSJacob Faibussowitsch PetscCall(MatLoad(A,fd)); 33*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* Read the matrix again as a seq matrix */ 36*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd)); 37*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&B)); 38*9566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATSEQAIJ)); 39*9566063dSJacob Faibussowitsch PetscCall(MatLoad(B,fd)); 40*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* Create the Random no generator */ 43*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&m,&n)); 44*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&r)); 45*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* Create the IS corresponding to subdomains */ 48*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd,&is1)); 49*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd,&is2)); 50*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m ,&idx)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Create the random Index Sets */ 53c4762a1bSJed Brown for (i=0; i<nd; i++) { 54c4762a1bSJed Brown for (j=0; j<rank; j++) { 55*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r,&rand)); 56c4762a1bSJed Brown } 57*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r,&rand)); 58c4762a1bSJed Brown lsize = (PetscInt)(rand*m); 59c4762a1bSJed Brown for (j=0; j<lsize; j++) { 60*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r,&rand)); 61c4762a1bSJed Brown idx[j] = (PetscInt)(rand*m); 62c4762a1bSJed Brown } 63*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,lsize,idx,PETSC_COPY_VALUES,is1+i)); 64*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,lsize,idx,PETSC_COPY_VALUES,is2+i)); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67*9566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A,nd,is1,ov)); 68*9566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(B,nd,is2,ov)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* Now see if the serial and parallel case have the same answers */ 71c4762a1bSJed Brown for (i=0; i<nd; ++i) { 72c4762a1bSJed Brown PetscInt sz1,sz2; 73*9566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i],is2[i],&flg)); 74*9566063dSJacob Faibussowitsch PetscCall(ISGetSize(is1[i],&sz1)); 75*9566063dSJacob Faibussowitsch PetscCall(ISGetSize(is2[i],&sz2)); 7628b400f6SJacob 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); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Free Allocated Memory */ 80c4762a1bSJed Brown for (i=0; i<nd; ++i) { 81*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 82*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 83c4762a1bSJed Brown } 84*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 85*9566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 86*9566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 87*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 88*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 89*9566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 90*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 91b122ec5aSJacob Faibussowitsch return 0; 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 94c4762a1bSJed Brown /*TEST 95c4762a1bSJed Brown 96c4762a1bSJed Brown build: 97c4762a1bSJed Brown requires: !complex 98c4762a1bSJed Brown 99c4762a1bSJed Brown test: 100c4762a1bSJed Brown nsize: 3 101dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 102c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 3 -ov 1 103c4762a1bSJed Brown 104c4762a1bSJed Brown TEST*/ 105