1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap() - the parallel case. This example\n\ 3c4762a1bSJed Brown is similar to ex40.c; here the index sets used are random. Input arguments are:\n\ 4c4762a1bSJed Brown -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\ 5c4762a1bSJed Brown -nd <size> : > 0 no of domains per processor \n\ 6c4762a1bSJed Brown -ov <overlap> : >=0 amount of overlap between domains\n\n"; 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include <petscmat.h> 9c4762a1bSJed Brown 10*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 11*d71ae5a4SJacob Faibussowitsch { 12c4762a1bSJed Brown PetscInt nd = 2, ov = 1, i, j, m, n, *idx, lsize; 13c4762a1bSJed Brown PetscMPIInt rank; 14c4762a1bSJed Brown PetscBool flg; 15c4762a1bSJed Brown Mat A, B; 16c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 17c4762a1bSJed Brown PetscViewer fd; 18c4762a1bSJed Brown IS *is1, *is2; 19c4762a1bSJed Brown PetscRandom r; 20c4762a1bSJed Brown PetscScalar rand; 21c4762a1bSJed Brown 22327415f7SBarry Smith PetscFunctionBeginUser; 239566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), NULL)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* Read matrix and RHS */ 309566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 319566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 329566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPIAIJ)); 339566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* Read the matrix again as a seq matrix */ 379566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, file, FILE_MODE_READ, &fd)); 389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &B)); 399566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQAIJ)); 409566063dSJacob Faibussowitsch PetscCall(MatLoad(B, fd)); 419566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Create the Random no generator */ 449566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n)); 459566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r)); 469566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Create the IS corresponding to subdomains */ 499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is1)); 509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is2)); 519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idx)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Create the random Index Sets */ 54c4762a1bSJed Brown for (i = 0; i < nd; i++) { 5548a46eb9SPierre Jolivet for (j = 0; j < rank; j++) PetscCall(PetscRandomGetValue(r, &rand)); 569566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rand)); 57c4762a1bSJed Brown lsize = (PetscInt)(rand * m); 58c4762a1bSJed Brown for (j = 0; j < lsize; j++) { 599566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rand)); 60c4762a1bSJed Brown idx[j] = (PetscInt)(rand * m); 61c4762a1bSJed Brown } 629566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize, idx, PETSC_COPY_VALUES, is1 + i)); 639566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize, idx, PETSC_COPY_VALUES, is2 + i)); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, nd, is1, ov)); 679566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(B, nd, is2, ov)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Now see if the serial and parallel case have the same answers */ 70c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 71c4762a1bSJed Brown PetscInt sz1, sz2; 729566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i], is2[i], &flg)); 739566063dSJacob Faibussowitsch PetscCall(ISGetSize(is1[i], &sz1)); 749566063dSJacob Faibussowitsch PetscCall(ISGetSize(is2[i], &sz2)); 7528b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "proc:[%d], i=%" PetscInt_FMT ", flg =%d sz1 = %" PetscInt_FMT " sz2 = %" PetscInt_FMT, rank, i, (int)flg, sz1, sz2); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Free Allocated Memory */ 79c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 82c4762a1bSJed Brown } 839566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 849566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 859566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 889566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 90b122ec5aSJacob Faibussowitsch return 0; 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93c4762a1bSJed Brown /*TEST 94c4762a1bSJed Brown 95c4762a1bSJed Brown build: 96c4762a1bSJed Brown requires: !complex 97c4762a1bSJed Brown 98c4762a1bSJed Brown test: 99c4762a1bSJed Brown nsize: 3 100dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex 101c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 3 -ov 1 102c4762a1bSJed Brown 103c4762a1bSJed Brown TEST*/ 104