1c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap() - the parallel case. This example\n\
2c4762a1bSJed Brown is similar to ex40.c; here the index sets used are random. Input arguments are:\n\
3c4762a1bSJed Brown -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\
4c4762a1bSJed Brown -nd <size> : > 0 no of domains per processor \n\
5c4762a1bSJed Brown -ov <overlap> : >=0 amount of overlap between domains\n\n";
6c4762a1bSJed Brown
7c4762a1bSJed Brown #include <petscmat.h>
8c4762a1bSJed Brown
main(int argc,char ** args)9d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
10d71ae5a4SJacob Faibussowitsch {
11c4762a1bSJed Brown PetscInt nd = 2, ov = 1, i, j, m, n, *idx, lsize;
12c4762a1bSJed Brown PetscMPIInt rank;
13c4762a1bSJed Brown PetscBool flg;
14c4762a1bSJed Brown Mat A, B;
15c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN];
16c4762a1bSJed Brown PetscViewer fd;
17c4762a1bSJed Brown IS *is1, *is2;
18c4762a1bSJed Brown PetscRandom r;
19c4762a1bSJed Brown PetscScalar rand;
20c4762a1bSJed Brown
21327415f7SBarry Smith PetscFunctionBeginUser;
22c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), NULL));
259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL));
269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL));
27c4762a1bSJed Brown
28c4762a1bSJed Brown /* Read matrix and RHS */
299566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
309566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
319566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPIAIJ));
329566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd));
339566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd));
34c4762a1bSJed Brown
35c4762a1bSJed Brown /* Read the matrix again as a seq matrix */
369566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, file, FILE_MODE_READ, &fd));
379566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &B));
389566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQAIJ));
399566063dSJacob Faibussowitsch PetscCall(MatLoad(B, fd));
409566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd));
41c4762a1bSJed Brown
42c4762a1bSJed Brown /* Create the Random no generator */
439566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n));
449566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
459566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r));
46c4762a1bSJed Brown
47c4762a1bSJed Brown /* Create the IS corresponding to subdomains */
489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is1));
499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is2));
509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idx));
51c4762a1bSJed Brown
52c4762a1bSJed Brown /* Create the random Index Sets */
53c4762a1bSJed Brown for (i = 0; i < nd; i++) {
5448a46eb9SPierre Jolivet for (j = 0; j < rank; j++) PetscCall(PetscRandomGetValue(r, &rand));
559566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rand));
56c4762a1bSJed Brown lsize = (PetscInt)(rand * m);
57c4762a1bSJed Brown for (j = 0; j < lsize; j++) {
589566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &rand));
59c4762a1bSJed Brown idx[j] = (PetscInt)(rand * m);
60c4762a1bSJed Brown }
619566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize, idx, PETSC_COPY_VALUES, is1 + i));
629566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize, idx, PETSC_COPY_VALUES, is2 + i));
63c4762a1bSJed Brown }
64c4762a1bSJed Brown
659566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, nd, is1, ov));
669566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(B, nd, is2, ov));
67c4762a1bSJed Brown
68c4762a1bSJed Brown /* Now see if the serial and parallel case have the same answers */
69c4762a1bSJed Brown for (i = 0; i < nd; ++i) {
70c4762a1bSJed Brown PetscInt sz1, sz2;
719566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i], is2[i], &flg));
729566063dSJacob Faibussowitsch PetscCall(ISGetSize(is1[i], &sz1));
739566063dSJacob Faibussowitsch PetscCall(ISGetSize(is2[i], &sz2));
7428b400f6SJacob 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);
75c4762a1bSJed Brown }
76c4762a1bSJed Brown
77c4762a1bSJed Brown /* Free Allocated Memory */
78c4762a1bSJed Brown for (i = 0; i < nd; ++i) {
799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i]));
809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i]));
81c4762a1bSJed Brown }
829566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r));
839566063dSJacob Faibussowitsch PetscCall(PetscFree(is1));
849566063dSJacob Faibussowitsch PetscCall(PetscFree(is2));
859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
879566063dSJacob Faibussowitsch PetscCall(PetscFree(idx));
889566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
89b122ec5aSJacob Faibussowitsch return 0;
90c4762a1bSJed Brown }
91c4762a1bSJed Brown
92c4762a1bSJed Brown /*TEST
93c4762a1bSJed Brown
94c4762a1bSJed Brown build:
95c4762a1bSJed Brown requires: !complex
96c4762a1bSJed Brown
97c4762a1bSJed Brown test:
98c4762a1bSJed Brown nsize: 3
99dfd57a17SPierre Jolivet requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
100c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/arco1 -nd 3 -ov 1
101*3886731fSPierre Jolivet output_file: output/empty.out
102c4762a1bSJed Brown
103c4762a1bSJed Brown TEST*/
104