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