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; 265f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 275f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_unsorted",&test_unsorted,NULL)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Read matrix A and RHS */ 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,fd)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Read the same matrix as a seq matrix B */ 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&B)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATSEQAIJ)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(B,fd)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 48c4762a1bSJed Brown 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(A,&bs)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Create the Random no generator */ 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&m,&n)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&r)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(r)); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Create the IS corresponding to subdomains */ 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nd,&is1)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nd,&is2)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m ,&idx)); 602f613bf5SBarry Smith for (i = 0; i < m; i++) {idx[i] = i;} 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++) { 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(r,&rand)); 67c4762a1bSJed Brown } 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(r,&rand)); 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 745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(r,&rand)); 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 } 825f80ce2aSJacob Faibussowitsch if (!test_unsorted) CHKERRQ(PetscSortInt(lsize*bs,idx)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is1+i)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is2+i)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(ISSetBlockSize(is1[i],bs)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(ISSetBlockSize(is2[i],bs)); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89c4762a1bSJed Brown if (!test_unsorted) { 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatIncreaseOverlap(A,nd,is1,ov)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatIncreaseOverlap(B,nd,is2,ov)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown for (i=0; i<nd; ++i) { 945f80ce2aSJacob Faibussowitsch CHKERRQ(ISSort(is1[i])); 955f80ce2aSJacob Faibussowitsch CHKERRQ(ISSort(is2[i])); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* Now see if the serial and parallel case have the same answers */ 103c4762a1bSJed Brown for (i=0; i<nd; ++i) { 1045f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(submatA[i],submatB[i],&flg)); 105*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"%" PetscInt_FMT "-th paralle submatA != seq submatB",i); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Free Allocated Memory */ 109c4762a1bSJed Brown for (i=0; i<nd; ++i) { 1105f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is1[i])); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is2[i])); 112c4762a1bSJed Brown } 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroySubMatrices(nd,&submatA)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroySubMatrices(nd,&submatB)); 115c4762a1bSJed Brown 1165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&r)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(is1)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(is2)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(idx)); 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 133dfd57a17SPierre 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 139dfd57a17SPierre 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 144dfd57a17SPierre 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 149dfd57a17SPierre 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 155dfd57a17SPierre 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 160dfd57a17SPierre 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