1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests the parallel case for MatIncreaseOverlap(). Input arguments are:\n\ 3*c4762a1bSJed Brown -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\ 4*c4762a1bSJed Brown -nd <size> : > 0 number of domains per processor \n\ 5*c4762a1bSJed Brown -ov <overlap> : >=0 amount of overlap between domains\n\n"; 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown #include <petscmat.h> 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown PetscErrorCode ISAllGatherDisjoint(IS iis, IS** ois) 10*c4762a1bSJed Brown { 11*c4762a1bSJed Brown IS *is2,is; 12*c4762a1bSJed Brown const PetscInt *idxs; 13*c4762a1bSJed Brown PetscInt i, ls,*sizes; 14*c4762a1bSJed Brown PetscErrorCode ierr; 15*c4762a1bSJed Brown PetscMPIInt size; 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown PetscFunctionBeginUser; 18*c4762a1bSJed Brown ierr = MPI_Comm_size(PetscObjectComm((PetscObject)iis),&size);CHKERRQ(ierr); 19*c4762a1bSJed Brown ierr = PetscMalloc1(size,&is2);CHKERRQ(ierr); 20*c4762a1bSJed Brown ierr = PetscMalloc1(size,&sizes);CHKERRQ(ierr); 21*c4762a1bSJed Brown ierr = ISGetLocalSize(iis,&ls);CHKERRQ(ierr); 22*c4762a1bSJed Brown /* we don't have a public ISGetLayout */ 23*c4762a1bSJed Brown ierr = MPI_Allgather(&ls,1,MPIU_INT,sizes,1,MPIU_INT,PetscObjectComm((PetscObject)iis));CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = ISAllGather(iis,&is);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr); 26*c4762a1bSJed Brown for (i = 0, ls = 0; i < size; i++) { 27*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,sizes[i],idxs+ls,PETSC_COPY_VALUES,&is2[i]);CHKERRQ(ierr); 28*c4762a1bSJed Brown ls += sizes[i]; 29*c4762a1bSJed Brown } 30*c4762a1bSJed Brown ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = PetscFree(sizes);CHKERRQ(ierr); 33*c4762a1bSJed Brown *ois = is2; 34*c4762a1bSJed Brown PetscFunctionReturn(0); 35*c4762a1bSJed Brown } 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown int main(int argc,char **args) 38*c4762a1bSJed Brown { 39*c4762a1bSJed Brown PetscErrorCode ierr; 40*c4762a1bSJed Brown PetscInt nd = 2,ov = 1,ndpar,i,start,m,n,end,lsize; 41*c4762a1bSJed Brown PetscMPIInt rank; 42*c4762a1bSJed Brown PetscBool flg, useND = PETSC_FALSE; 43*c4762a1bSJed Brown Mat A,B; 44*c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 45*c4762a1bSJed Brown PetscViewer fd; 46*c4762a1bSJed Brown IS *is1,*is2; 47*c4762a1bSJed Brown PetscRandom r; 48*c4762a1bSJed Brown PetscScalar rand; 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 51*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 53*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must use -f filename to indicate a file containing a PETSc binary matrix"); 54*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-nested_dissection",&useND,NULL);CHKERRQ(ierr); 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown /* Read matrix */ 59*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown /* Read the matrix again as a sequential matrix */ 67*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = MatLoad(B,fd);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = MatSetFromOptions(B);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown /* Create the IS corresponding to subdomains */ 75*c4762a1bSJed Brown if (useND) { 76*c4762a1bSJed Brown MatPartitioning part; 77*c4762a1bSJed Brown IS ndmap; 78*c4762a1bSJed Brown PetscMPIInt size; 79*c4762a1bSJed Brown 80*c4762a1bSJed Brown ndpar = 1; 81*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 82*c4762a1bSJed Brown nd = (PetscInt)size; 83*c4762a1bSJed Brown ierr = PetscMalloc1(ndpar,&is1);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = MatPartitioningCreate(PETSC_COMM_WORLD,&part);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = MatPartitioningSetAdjacency(part,A);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = MatPartitioningApplyND(part,&ndmap);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = ISBuildTwoSided(ndmap,NULL,&is1[0]);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = ISDestroy(&ndmap);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = ISAllGatherDisjoint(is1[0],&is2);CHKERRQ(ierr); 92*c4762a1bSJed Brown } else { 93*c4762a1bSJed Brown /* Create the random Index Sets */ 94*c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is1);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is2);CHKERRQ(ierr); 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&r);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 100*c4762a1bSJed Brown for (i=0; i<nd; i++) { 101*c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 102*c4762a1bSJed Brown start = (PetscInt)(rand*m); 103*c4762a1bSJed Brown ierr = PetscRandomGetValue(r,&rand);CHKERRQ(ierr); 104*c4762a1bSJed Brown end = (PetscInt)(rand*m); 105*c4762a1bSJed Brown lsize = end - start; 106*c4762a1bSJed Brown if (start > end) { start = end; lsize = -lsize;} 107*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,lsize,start,1,is1+i);CHKERRQ(ierr); 108*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,lsize,start,1,is2+i);CHKERRQ(ierr); 109*c4762a1bSJed Brown } 110*c4762a1bSJed Brown ndpar = nd; 111*c4762a1bSJed Brown ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 112*c4762a1bSJed Brown } 113*c4762a1bSJed Brown ierr = MatIncreaseOverlap(A,ndpar,is1,ov);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = MatIncreaseOverlap(B,nd,is2,ov);CHKERRQ(ierr); 115*c4762a1bSJed Brown if (useND) { 116*c4762a1bSJed Brown IS *is; 117*c4762a1bSJed Brown 118*c4762a1bSJed Brown ierr = ISAllGatherDisjoint(is1[0],&is);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = ISDestroy(&is1[0]);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = PetscFree(is1);CHKERRQ(ierr); 121*c4762a1bSJed Brown is1 = is; 122*c4762a1bSJed Brown } 123*c4762a1bSJed Brown /* Now see if the serial and parallel case have the same answers */ 124*c4762a1bSJed Brown for (i=0; i<nd; ++i) { 125*c4762a1bSJed Brown ierr = ISEqual(is1[i],is2[i],&flg);CHKERRQ(ierr); 126*c4762a1bSJed Brown if (!flg) { 127*c4762a1bSJed Brown ierr = ISViewFromOptions(is1[i],NULL,"-err_view");CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = ISViewFromOptions(is2[i],NULL,"-err_view");CHKERRQ(ierr); 129*c4762a1bSJed Brown SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"proc:[%d], i=%D, flg =%d\n",rank,i,(int)flg); 130*c4762a1bSJed Brown } 131*c4762a1bSJed Brown } 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown /* Free allocated memory */ 134*c4762a1bSJed Brown for (i=0; i<nd; ++i) { 135*c4762a1bSJed Brown ierr = ISDestroy(&is1[i]);CHKERRQ(ierr); 136*c4762a1bSJed Brown ierr = ISDestroy(&is2[i]);CHKERRQ(ierr); 137*c4762a1bSJed Brown } 138*c4762a1bSJed Brown ierr = PetscFree(is1);CHKERRQ(ierr); 139*c4762a1bSJed Brown ierr = PetscFree(is2);CHKERRQ(ierr); 140*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 141*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 142*c4762a1bSJed Brown ierr = PetscFinalize(); 143*c4762a1bSJed Brown return ierr; 144*c4762a1bSJed Brown } 145*c4762a1bSJed Brown 146*c4762a1bSJed Brown 147*c4762a1bSJed Brown 148*c4762a1bSJed Brown /*TEST 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown build: 151*c4762a1bSJed Brown requires: !complex 152*c4762a1bSJed Brown 153*c4762a1bSJed Brown testset: 154*c4762a1bSJed Brown nsize: 5 155*c4762a1bSJed Brown requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) !complex 156*c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info -ov 2 157*c4762a1bSJed Brown output_file: output/ex40_1.out 158*c4762a1bSJed Brown test: 159*c4762a1bSJed Brown suffix: 1 160*c4762a1bSJed Brown args: -nd 7 161*c4762a1bSJed Brown test: 162*c4762a1bSJed Brown requires: parmetis 163*c4762a1bSJed Brown suffix: 1_nd 164*c4762a1bSJed Brown args: -nested_dissection -mat_partitioning_type parmetis 165*c4762a1bSJed Brown 166*c4762a1bSJed Brown testset: 167*c4762a1bSJed Brown nsize: 3 168*c4762a1bSJed Brown requires: double !define(PETSC_USE_64BIT_INDICES) !complex 169*c4762a1bSJed Brown args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 -mat_increase_overlap_scalable 1 -ov 2 170*c4762a1bSJed Brown output_file: output/ex40_1.out 171*c4762a1bSJed Brown test: 172*c4762a1bSJed Brown suffix: 2 173*c4762a1bSJed Brown args: -nd 7 174*c4762a1bSJed Brown test: 175*c4762a1bSJed Brown requires: parmetis 176*c4762a1bSJed Brown suffix: 2_nd 177*c4762a1bSJed Brown args: -nested_dissection -mat_partitioning_type parmetis 178*c4762a1bSJed Brown 179*c4762a1bSJed Brown TEST*/ 180