Actual source code: ex42.c


  2: static char help[] = "Tests MatIncreaseOverlap() and MatCreateSubmatrices() for the parallel case.\n\
  3: This example is similar to ex40.c; here the index sets used are random.\n\
  4: Input arguments are:\n\
  5:   -f <input_file> : file to load.  For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\
  6:   -nd <size>      : > 0  no of domains per processor \n\
  7:   -ov <overlap>   : >=0  amount of overlap between domains\n\n";

  9: #include <petscmat.h>

 11: int main(int argc,char **args)
 12: {
 13:   PetscInt       nd = 2,ov=1,i,j,lsize,m,n,*idx,bs;
 14:   PetscMPIInt    rank, size;
 15:   PetscBool      flg;
 16:   Mat            A,B,*submatA,*submatB;
 17:   char           file[PETSC_MAX_PATH_LEN];
 18:   PetscViewer    fd;
 19:   IS             *is1,*is2;
 20:   PetscRandom    r;
 21:   PetscBool      test_unsorted = PETSC_FALSE;
 22:   PetscScalar    rand;

 24:   PetscInitialize(&argc,&args,(char*)0,help);
 25:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 26:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 27:   PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL);
 28:   PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL);
 30:   PetscOptionsGetBool(NULL,NULL,"-test_unsorted",&test_unsorted,NULL);

 32:   /* Read matrix A and RHS */
 33:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);
 34:   MatCreate(PETSC_COMM_WORLD,&A);
 35:   MatSetType(A,MATAIJ);
 36:   MatSetFromOptions(A);
 37:   MatLoad(A,fd);
 38:   PetscViewerDestroy(&fd);

 40:   /* Read the same matrix as a seq matrix B */
 41:   PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd);
 42:   MatCreate(PETSC_COMM_SELF,&B);
 43:   MatSetType(B,MATSEQAIJ);
 44:   MatSetFromOptions(B);
 45:   MatLoad(B,fd);
 46:   PetscViewerDestroy(&fd);

 48:   MatGetBlockSize(A,&bs);

 50:   /* Create the Random no generator */
 51:   MatGetSize(A,&m,&n);
 52:   PetscRandomCreate(PETSC_COMM_SELF,&r);
 53:   PetscRandomSetFromOptions(r);

 55:   /* Create the IS corresponding to subdomains */
 56:   PetscMalloc1(nd,&is1);
 57:   PetscMalloc1(nd,&is2);
 58:   PetscMalloc1(m ,&idx);
 59:   for (i = 0; i < m; i++) {idx[i] = i;}

 61:   /* Create the random Index Sets */
 62:   for (i=0; i<nd; i++) {
 63:     /* Skip a few,so that the IS on different procs are diffeent*/
 64:     for (j=0; j<rank; j++) {
 65:       PetscRandomGetValue(r,&rand);
 66:     }
 67:     PetscRandomGetValue(r,&rand);
 68:     lsize = (PetscInt)(rand*(m/bs));
 69:     /* shuffle */
 70:     for (j=0; j<lsize; j++) {
 71:       PetscInt k, swap, l;

 73:       PetscRandomGetValue(r,&rand);
 74:       k      = j + (PetscInt)(rand*((m/bs)-j));
 75:       for (l = 0; l < bs; l++) {
 76:         swap        = idx[bs*j+l];
 77:         idx[bs*j+l] = idx[bs*k+l];
 78:         idx[bs*k+l] = swap;
 79:       }
 80:     }
 81:     if (!test_unsorted) PetscSortInt(lsize*bs,idx);
 82:     ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is1+i);
 83:     ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is2+i);
 84:     ISSetBlockSize(is1[i],bs);
 85:     ISSetBlockSize(is2[i],bs);
 86:   }

 88:   if (!test_unsorted) {
 89:     MatIncreaseOverlap(A,nd,is1,ov);
 90:     MatIncreaseOverlap(B,nd,is2,ov);

 92:     for (i=0; i<nd; ++i) {
 93:       ISSort(is1[i]);
 94:       ISSort(is2[i]);
 95:     }
 96:   }

 98:   MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA);
 99:   MatCreateSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB);

101:   /* Now see if the serial and parallel case have the same answers */
102:   for (i=0; i<nd; ++i) {
103:     MatEqual(submatA[i],submatB[i],&flg);
105:   }

107:   /* Free Allocated Memory */
108:   for (i=0; i<nd; ++i) {
109:     ISDestroy(&is1[i]);
110:     ISDestroy(&is2[i]);
111:   }
112:   MatDestroySubMatrices(nd,&submatA);
113:   MatDestroySubMatrices(nd,&submatB);

115:   PetscRandomDestroy(&r);
116:   PetscFree(is1);
117:   PetscFree(is2);
118:   MatDestroy(&A);
119:   MatDestroy(&B);
120:   PetscFree(idx);
121:   PetscFinalize();
122:   return 0;
123: }

125: /*TEST

127:    build:
128:       requires: !complex

130:    test:
131:       nsize: 3
132:       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
133:       args: -f ${DATAFILESPATH}/matrices/arco1 -nd 5 -ov 2

135:    test:
136:       suffix: 2
137:       args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -ov 2
138:       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex

140:    test:
141:       suffix: unsorted_baij_mpi
142:       nsize: 3
143:       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
144:       args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -nd 8 -mat_type baij -test_unsorted

146:    test:
147:       suffix: unsorted_baij_seq
148:       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
149:       args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -nd 8 -mat_type baij -test_unsorted

151:    test:
152:       suffix: unsorted_mpi
153:       nsize: 3
154:       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
155:       args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -test_unsorted

157:    test:
158:       suffix: unsorted_seq
159:       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
160:       args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -test_unsorted

162: TEST*/