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; 24ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 25589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL);CHKERRQ(ierr); 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* Read matrix and RHS */ 30c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 33c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* Read the matrix again as a seq matrix */ 37c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = MatLoad(B,fd);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Create the Random no generator */ 44c4762a1bSJed Brown ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&r);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Create the IS corresponding to subdomains */ 49c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is1);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is2);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = PetscMalloc1(m ,&idx);CHKERRQ(ierr); 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++) { 56c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 59c4762a1bSJed Brown lsize = (PetscInt)(rand*m); 60c4762a1bSJed Brown for (j=0; j<lsize; j++) { 61c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 62c4762a1bSJed Brown idx[j] = (PetscInt)(rand*m); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,lsize,idx,PETSC_COPY_VALUES,is1+i);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,lsize,idx,PETSC_COPY_VALUES,is2+i);CHKERRQ(ierr); 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown ierr = MatIncreaseOverlap(A,nd,is1,ov);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = MatIncreaseOverlap(B,nd,is2,ov);CHKERRQ(ierr); 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; 74c4762a1bSJed Brown ierr = ISEqual(is1[i],is2[i],&flg);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = ISGetSize(is1[i],&sz1);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = ISGetSize(is2[i],&sz2);CHKERRQ(ierr); 772758c9b9SBarry Smith 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); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Free Allocated Memory */ 81c4762a1bSJed Brown for (i=0; i<nd; ++i) { 82c4762a1bSJed Brown ierr = ISDestroy(&is1[i]);CHKERRQ(ierr); 83c4762a1bSJed Brown ierr = ISDestroy(&is2[i]);CHKERRQ(ierr); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = PetscFree(is1);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = PetscFree(is2);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = PetscFree(idx);CHKERRQ(ierr); 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 102*dfd57a17SPierre 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