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