1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap() and MatCreateSubmatrices() for the parallel case.\n\ 3c4762a1bSJed Brown This example is similar to ex40.c; here the index sets used are random.\n\ 4c4762a1bSJed Brown Input arguments are:\n\ 5c4762a1bSJed Brown -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\ 6c4762a1bSJed Brown -nd <size> : > 0 no of domains per processor \n\ 7c4762a1bSJed Brown -ov <overlap> : >=0 amount of overlap between domains\n\n"; 8c4762a1bSJed Brown 9c4762a1bSJed Brown #include <petscmat.h> 10c4762a1bSJed Brown 11c4762a1bSJed Brown int main(int argc,char **args) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown PetscErrorCode ierr; 14c4762a1bSJed Brown PetscInt nd = 2,ov=1,i,j,lsize,m,n,*idx,bs; 15c4762a1bSJed Brown PetscMPIInt rank, size; 16c4762a1bSJed Brown PetscBool flg; 17c4762a1bSJed Brown Mat A,B,*submatA,*submatB; 18c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 19c4762a1bSJed Brown PetscViewer fd; 20c4762a1bSJed Brown IS *is1,*is2; 21c4762a1bSJed Brown PetscRandom r; 22c4762a1bSJed Brown PetscBool test_unsorted = PETSC_FALSE; 23c4762a1bSJed Brown PetscScalar rand; 24c4762a1bSJed Brown 25c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 26ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 27ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 28589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_unsorted",&test_unsorted,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Read matrix A and RHS */ 34c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Read the same matrix as a seq matrix B */ 42c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = MatSetFromOptions(B);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = MatLoad(B,fd);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 48c4762a1bSJed Brown 49c4762a1bSJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Create the Random no generator */ 52c4762a1bSJed Brown ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&r);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Create the IS corresponding to subdomains */ 57c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is1);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is2);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = PetscMalloc1(m ,&idx);CHKERRQ(ierr); 60c4762a1bSJed Brown for (i = 0; i < m; i++) {idx[i] = i;CHKERRQ(ierr);} 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* Create the random Index Sets */ 63c4762a1bSJed Brown for (i=0; i<nd; i++) { 64c4762a1bSJed Brown /* Skip a few,so that the IS on different procs are diffeent*/ 65c4762a1bSJed Brown for (j=0; j<rank; j++) { 66c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 67c4762a1bSJed Brown } 68c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 69c4762a1bSJed Brown lsize = (PetscInt)(rand*(m/bs)); 70c4762a1bSJed Brown /* shuffle */ 71c4762a1bSJed Brown for (j=0; j<lsize; j++) { 72c4762a1bSJed Brown PetscInt k, swap, l; 73c4762a1bSJed Brown 74c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 75c4762a1bSJed Brown k = j + (PetscInt)(rand*((m/bs)-j)); 76c4762a1bSJed Brown for (l = 0; l < bs; l++) { 77c4762a1bSJed Brown swap = idx[bs*j+l]; 78c4762a1bSJed Brown idx[bs*j+l] = idx[bs*k+l]; 79c4762a1bSJed Brown idx[bs*k+l] = swap; 80c4762a1bSJed Brown } 81c4762a1bSJed Brown } 82c4762a1bSJed Brown if (!test_unsorted) {ierr = PetscSortInt(lsize*bs,idx);CHKERRQ(ierr);} 83c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is1+i);CHKERRQ(ierr); 84c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is2+i);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = ISSetBlockSize(is1[i],bs);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = ISSetBlockSize(is2[i],bs);CHKERRQ(ierr); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89c4762a1bSJed Brown if (!test_unsorted) { 90c4762a1bSJed Brown ierr = MatIncreaseOverlap(A,nd,is1,ov);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = MatIncreaseOverlap(B,nd,is2,ov);CHKERRQ(ierr); 92c4762a1bSJed Brown 93c4762a1bSJed Brown for (i=0; i<nd; ++i) { 94c4762a1bSJed Brown ierr = ISSort(is1[i]);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = ISSort(is2[i]);CHKERRQ(ierr); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 99c4762a1bSJed Brown ierr = MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = MatCreateSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB);CHKERRQ(ierr); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* Now see if the serial and parallel case have the same answers */ 103c4762a1bSJed Brown for (i=0; i<nd; ++i) { 104c4762a1bSJed Brown ierr = MatEqual(submatA[i],submatB[i],&flg);CHKERRQ(ierr); 105c4762a1bSJed Brown if (!flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"%D-th paralle submatA != seq submatB",i); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Free Allocated Memory */ 109c4762a1bSJed Brown for (i=0; i<nd; ++i) { 110c4762a1bSJed Brown ierr = ISDestroy(&is1[i]);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = ISDestroy(&is2[i]);CHKERRQ(ierr); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown ierr = MatDestroySubMatrices(nd,&submatA);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = MatDestroySubMatrices(nd,&submatB);CHKERRQ(ierr); 115c4762a1bSJed Brown 116c4762a1bSJed Brown ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = PetscFree(is1);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = PetscFree(is2);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = PetscFree(idx);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = PetscFinalize(); 123c4762a1bSJed Brown return ierr; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown /*TEST 127c4762a1bSJed Brown 128c4762a1bSJed Brown build: 129c4762a1bSJed Brown requires: !complex 130c4762a1bSJed Brown 131c4762a1bSJed Brown test: 132c4762a1bSJed Brown nsize: 3 133*dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 134c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 5 -ov 2 135c4762a1bSJed Brown 136c4762a1bSJed Brown test: 137c4762a1bSJed Brown suffix: 2 138c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -ov 2 139*dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 140c4762a1bSJed Brown 141c4762a1bSJed Brown test: 142c4762a1bSJed Brown suffix: unsorted_baij_mpi 143c4762a1bSJed Brown nsize: 3 144*dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 145c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -nd 8 -mat_type baij -test_unsorted 146c4762a1bSJed Brown 147c4762a1bSJed Brown test: 148c4762a1bSJed Brown suffix: unsorted_baij_seq 149*dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 150c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -nd 8 -mat_type baij -test_unsorted 151c4762a1bSJed Brown 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: unsorted_mpi 154c4762a1bSJed Brown nsize: 3 155*dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 156c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -test_unsorted 157c4762a1bSJed Brown 158c4762a1bSJed Brown test: 159c4762a1bSJed Brown suffix: unsorted_seq 160*dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 161c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -test_unsorted 162c4762a1bSJed Brown 163c4762a1bSJed Brown TEST*/ 164