1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap() - the parallel case. This example\n\ 3*c4762a1bSJed Brown is similar to ex40.c; here the index sets used are random. Input arguments are:\n\ 4*c4762a1bSJed Brown -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\ 5*c4762a1bSJed Brown -nd <size> : > 0 no of domains per processor \n\ 6*c4762a1bSJed Brown -ov <overlap> : >=0 amount of overlap between domains\n\n"; 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown #include <petscmat.h> 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown int main(int argc,char **args) 11*c4762a1bSJed Brown { 12*c4762a1bSJed Brown PetscInt nd = 2,ov=1,i,j,m,n,*idx,lsize; 13*c4762a1bSJed Brown PetscErrorCode ierr; 14*c4762a1bSJed Brown PetscMPIInt rank; 15*c4762a1bSJed Brown PetscBool flg; 16*c4762a1bSJed Brown Mat A,B; 17*c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 18*c4762a1bSJed Brown PetscViewer fd; 19*c4762a1bSJed Brown IS *is1,*is2; 20*c4762a1bSJed Brown PetscRandom r; 21*c4762a1bSJed Brown PetscScalar rand; 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 24*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL);CHKERRQ(ierr); 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown /* Read matrix and RHS */ 30*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown /* Read the matrix again as a seq matrix */ 37*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = MatLoad(B,fd);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown /* Create the Random no generator */ 44*c4762a1bSJed Brown ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&r);CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown /* Create the IS corresponding to subdomains */ 49*c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is1);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is2);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = PetscMalloc1(m ,&idx);CHKERRQ(ierr); 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown /* Create the random Index Sets */ 54*c4762a1bSJed Brown for (i=0; i<nd; i++) { 55*c4762a1bSJed Brown for (j=0; j<rank; j++) { 56*c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 57*c4762a1bSJed Brown } 58*c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 59*c4762a1bSJed Brown lsize = (PetscInt)(rand*m); 60*c4762a1bSJed Brown for (j=0; j<lsize; j++) { 61*c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 62*c4762a1bSJed Brown idx[j] = (PetscInt)(rand*m); 63*c4762a1bSJed Brown } 64*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,lsize,idx,PETSC_COPY_VALUES,is1+i);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,lsize,idx,PETSC_COPY_VALUES,is2+i);CHKERRQ(ierr); 66*c4762a1bSJed Brown } 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown ierr = MatIncreaseOverlap(A,nd,is1,ov);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = MatIncreaseOverlap(B,nd,is2,ov);CHKERRQ(ierr); 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown /* Now see if the serial and parallel case have the same answers */ 72*c4762a1bSJed Brown for (i=0; i<nd; ++i) { 73*c4762a1bSJed Brown PetscInt sz1,sz2; 74*c4762a1bSJed Brown ierr = ISEqual(is1[i],is2[i],&flg);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = ISGetSize(is1[i],&sz1);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = ISGetSize(is2[i],&sz2);CHKERRQ(ierr); 77*c4762a1bSJed Brown if (!flg) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"proc:[%d], i=%D, flg =%d sz1 = %D sz2 = %D\n",rank,i,(int)flg,sz1,sz2);CHKERRQ(ierr); 78*c4762a1bSJed Brown } 79*c4762a1bSJed Brown 80*c4762a1bSJed Brown /* Free Allocated Memory */ 81*c4762a1bSJed Brown for (i=0; i<nd; ++i) { 82*c4762a1bSJed Brown ierr = ISDestroy(&is1[i]);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = ISDestroy(&is2[i]);CHKERRQ(ierr); 84*c4762a1bSJed Brown } 85*c4762a1bSJed Brown ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = PetscFree(is1);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = PetscFree(is2);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = PetscFree(idx);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = PetscFinalize(); 92*c4762a1bSJed Brown return ierr; 93*c4762a1bSJed Brown } 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown /*TEST 98*c4762a1bSJed Brown 99*c4762a1bSJed Brown build: 100*c4762a1bSJed Brown requires: !complex 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown test: 103*c4762a1bSJed Brown nsize: 3 104*c4762a1bSJed Brown requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) !complex 105*c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 3 -ov 1 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown TEST*/ 109