xref: /petsc/src/mat/tests/ex42.c (revision da81f9329be15cc55f054c8a00978087195c9247)
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 
11d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
12d71ae5a4SJacob Faibussowitsch {
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 
24327415f7SBarry Smith   PetscFunctionBeginUser;
259566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), NULL));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL));
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_unsorted", &test_unsorted, NULL));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /* Read matrix A and RHS */
349566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
359566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
369566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
379566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
389566063dSJacob Faibussowitsch   PetscCall(MatLoad(A, fd));
399566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* Read the same matrix as a seq matrix B */
429566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, file, FILE_MODE_READ, &fd));
439566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &B));
449566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATSEQAIJ));
459566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
469566063dSJacob Faibussowitsch   PetscCall(MatLoad(B, fd));
479566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Create the Random no generator */
529566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &m, &n));
539566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
549566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(r));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Create the IS corresponding to subdomains */
579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is1));
589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is2));
599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &idx));
60ad540459SPierre Jolivet   for (i = 0; i < m; i++) idx[i] = i;
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* Create the random Index Sets */
63c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
64*da81f932SPierre Jolivet     /* Skip a few,so that the IS on different procs are different*/
6548a46eb9SPierre Jolivet     for (j = 0; j < rank; j++) PetscCall(PetscRandomGetValue(r, &rand));
669566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(r, &rand));
67c4762a1bSJed Brown     lsize = (PetscInt)(rand * (m / bs));
68c4762a1bSJed Brown     /* shuffle */
69c4762a1bSJed Brown     for (j = 0; j < lsize; j++) {
70c4762a1bSJed Brown       PetscInt k, swap, l;
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(r, &rand));
73c4762a1bSJed Brown       k = j + (PetscInt)(rand * ((m / bs) - j));
74c4762a1bSJed Brown       for (l = 0; l < bs; l++) {
75c4762a1bSJed Brown         swap            = idx[bs * j + l];
76c4762a1bSJed Brown         idx[bs * j + l] = idx[bs * k + l];
77c4762a1bSJed Brown         idx[bs * k + l] = swap;
78c4762a1bSJed Brown       }
79c4762a1bSJed Brown     }
809566063dSJacob Faibussowitsch     if (!test_unsorted) PetscCall(PetscSortInt(lsize * bs, idx));
819566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is1 + i));
829566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is2 + i));
839566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize(is1[i], bs));
849566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize(is2[i], bs));
85c4762a1bSJed Brown   }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   if (!test_unsorted) {
889566063dSJacob Faibussowitsch     PetscCall(MatIncreaseOverlap(A, nd, is1, ov));
899566063dSJacob Faibussowitsch     PetscCall(MatIncreaseOverlap(B, nd, is2, ov));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown     for (i = 0; i < nd; ++i) {
929566063dSJacob Faibussowitsch       PetscCall(ISSort(is1[i]));
939566063dSJacob Faibussowitsch       PetscCall(ISSort(is2[i]));
94c4762a1bSJed Brown     }
95c4762a1bSJed Brown   }
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA));
989566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_INITIAL_MATRIX, &submatB));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* Now see if the serial and parallel case have the same answers */
101c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1029566063dSJacob Faibussowitsch     PetscCall(MatEqual(submatA[i], submatB[i], &flg));
1036aad120cSJose E. Roman     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "%" PetscInt_FMT "-th parallel submatA != seq submatB", i);
104c4762a1bSJed Brown   }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* Free Allocated Memory */
107c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1089566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1[i]));
1099566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2[i]));
110c4762a1bSJed Brown   }
1119566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd, &submatA));
1129566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd, &submatB));
113c4762a1bSJed Brown 
1149566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&r));
1159566063dSJacob Faibussowitsch   PetscCall(PetscFree(is1));
1169566063dSJacob Faibussowitsch   PetscCall(PetscFree(is2));
1179566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1199566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
1209566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
121b122ec5aSJacob Faibussowitsch   return 0;
122c4762a1bSJed Brown }
123c4762a1bSJed Brown 
124c4762a1bSJed Brown /*TEST
125c4762a1bSJed Brown 
126c4762a1bSJed Brown    build:
127c4762a1bSJed Brown       requires: !complex
128c4762a1bSJed Brown 
129c4762a1bSJed Brown    test:
130c4762a1bSJed Brown       nsize: 3
131dfd57a17SPierre Jolivet       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
132c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1 -nd 5 -ov 2
133c4762a1bSJed Brown 
134c4762a1bSJed Brown    test:
135c4762a1bSJed Brown       suffix: 2
136c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -ov 2
137dfd57a17SPierre Jolivet       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
138c4762a1bSJed Brown 
139c4762a1bSJed Brown    test:
140c4762a1bSJed Brown       suffix: unsorted_baij_mpi
141c4762a1bSJed Brown       nsize: 3
142dfd57a17SPierre Jolivet       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
143c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -nd 8 -mat_type baij -test_unsorted
144c4762a1bSJed Brown 
145c4762a1bSJed Brown    test:
146c4762a1bSJed Brown       suffix: unsorted_baij_seq
147dfd57a17SPierre Jolivet       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
148c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -nd 8 -mat_type baij -test_unsorted
149c4762a1bSJed Brown 
150c4762a1bSJed Brown    test:
151c4762a1bSJed Brown       suffix: unsorted_mpi
152c4762a1bSJed Brown       nsize: 3
153dfd57a17SPierre Jolivet       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
154c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -test_unsorted
155c4762a1bSJed Brown 
156c4762a1bSJed Brown    test:
157c4762a1bSJed Brown       suffix: unsorted_seq
158dfd57a17SPierre Jolivet       requires: datafilespath double !defined(PETSC_USE_64BIT_INDICES) !complex
159c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/arco1 -nd 8 -test_unsorted
160c4762a1bSJed Brown 
161c4762a1bSJed Brown TEST*/
162