xref: /petsc/src/mat/tests/ex51.c (revision 4b5966dda8e0e9d52f06e77e04a44e922c8134a2)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for MatBAIJ format.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8*4b5966ddSBarry Smith   Mat          A, B, E, *submatA, *submatB;
9c4762a1bSJed Brown   PetscInt     bs = 1, m = 43, ov = 1, i, j, k, *rows, *cols, M, nd = 5, *idx, mm, nn, lsize;
10c4762a1bSJed Brown   PetscScalar *vals, rval;
11c4762a1bSJed Brown   IS          *is1, *is2;
12c4762a1bSJed Brown   PetscRandom  rdm;
13c4762a1bSJed Brown   Vec          xx, s1, s2;
14c4762a1bSJed Brown   PetscReal    s1norm, s2norm, rnorm, tol = PETSC_SQRT_MACHINE_EPSILON;
15c4762a1bSJed Brown   PetscBool    flg;
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL));
23c4762a1bSJed Brown   M = m * bs;
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, bs, M, M, 1, NULL, &A));
269566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
279566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, M, M, 15, NULL, &B));
28*4b5966ddSBarry Smith   PetscCall(MatSetBlockSize(B, bs));
299566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
309566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
319566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &rows));
349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &cols));
359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * bs, &vals));
369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M, &idx));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Now set blocks of values */
39c4762a1bSJed Brown   for (i = 0; i < 20 * bs; i++) {
40*4b5966ddSBarry Smith     PetscInt nr = 1, nc = 1;
419566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
42c4762a1bSJed Brown     cols[0] = bs * (int)(PetscRealPart(rval) * m);
439566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
44c4762a1bSJed Brown     rows[0] = bs * (int)(PetscRealPart(rval) * m);
45c4762a1bSJed Brown     for (j = 1; j < bs; j++) {
46*4b5966ddSBarry Smith       PetscCall(PetscRandomGetValue(rdm, &rval));
47*4b5966ddSBarry Smith       if (PetscRealPart(rval) > .5) rows[nr++] = rows[0] + j - 1;
48*4b5966ddSBarry Smith     }
49*4b5966ddSBarry Smith     for (j = 1; j < bs; j++) {
50*4b5966ddSBarry Smith       PetscCall(PetscRandomGetValue(rdm, &rval));
51*4b5966ddSBarry Smith       if (PetscRealPart(rval) > .5) cols[nc++] = cols[0] + j - 1;
52c4762a1bSJed Brown     }
53c4762a1bSJed Brown 
54*4b5966ddSBarry Smith     for (j = 0; j < nr * nc; j++) {
559566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm, &rval));
56c4762a1bSJed Brown       vals[j] = rval;
57c4762a1bSJed Brown     }
58*4b5966ddSBarry Smith     PetscCall(MatSetValues(A, nr, rows, nc, cols, vals, ADD_VALUES));
59*4b5966ddSBarry Smith     PetscCall(MatSetValues(B, nr, rows, nc, cols, vals, ADD_VALUES));
60c4762a1bSJed Brown   }
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
66c4762a1bSJed Brown 
67*4b5966ddSBarry Smith   /* Test MatConvert_SeqAIJ_SeqBAIJ handles incompletely filled blocks */
68*4b5966ddSBarry Smith   PetscCall(MatConvert(B, MATBAIJ, MAT_INITIAL_MATRIX, &E));
69*4b5966ddSBarry Smith   PetscCall(MatDestroy(&E));
70*4b5966ddSBarry Smith 
71c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is1));
739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is2));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
769566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
77c4762a1bSJed Brown     lsize = (int)(PetscRealPart(rval) * m);
78c4762a1bSJed Brown     for (j = 0; j < lsize; j++) {
799566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm, &rval));
80c4762a1bSJed Brown       idx[j * bs] = bs * (int)(PetscRealPart(rval) * m);
81c4762a1bSJed Brown       for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k;
82c4762a1bSJed Brown     }
839566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is1 + i));
849566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is2 + i));
85c4762a1bSJed Brown   }
869566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(A, nd, is1, ov));
879566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(B, nd, is2, ov));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
909566063dSJacob Faibussowitsch     PetscCall(ISEqual(is1[i], is2[i], &flg));
9128b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "i=%" PetscInt_FMT ", flg =%d", i, (int)flg);
92c4762a1bSJed Brown   }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
959566063dSJacob Faibussowitsch     PetscCall(ISSort(is1[i]));
969566063dSJacob Faibussowitsch     PetscCall(ISSort(is2[i]));
97c4762a1bSJed Brown   }
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA));
1009566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_INITIAL_MATRIX, &submatB));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* Test MatMult() */
103c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1049566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i], &mm, &nn));
1059566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
1069566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s1));
1079566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s2));
108c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
1099566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx, rdm));
1109566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i], xx, s1));
1119566063dSJacob Faibussowitsch       PetscCall(MatMult(submatB[i], xx, s2));
1129566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_2, &s1norm));
1139566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &s2norm));
114c4762a1bSJed Brown       rnorm = s2norm - s1norm;
11548a46eb9SPierre Jolivet       if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm));
116c4762a1bSJed Brown     }
1179566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1189566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1199566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
120c4762a1bSJed Brown   }
121c4762a1bSJed Brown   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
1229566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA));
1239566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_REUSE_MATRIX, &submatB));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Test MatMult() */
126c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1279566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i], &mm, &nn));
1289566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
1299566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s1));
1309566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s2));
131c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
1329566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx, rdm));
1339566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i], xx, s1));
1349566063dSJacob Faibussowitsch       PetscCall(MatMult(submatB[i], xx, s2));
1359566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_2, &s1norm));
1369566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &s2norm));
137c4762a1bSJed Brown       rnorm = s2norm - s1norm;
13848a46eb9SPierre Jolivet       if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm));
139c4762a1bSJed Brown     }
1409566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1419566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1429566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
143c4762a1bSJed Brown   }
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* Free allocated memory */
146c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1479566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1[i]));
1489566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2[i]));
149c4762a1bSJed Brown   }
1509566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd, &submatA));
1519566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd, &submatB));
1529566063dSJacob Faibussowitsch   PetscCall(PetscFree(is1));
1539566063dSJacob Faibussowitsch   PetscCall(PetscFree(is2));
1549566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
1569566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals));
1589566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1609566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
1619566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
162b122ec5aSJacob Faibussowitsch   return 0;
163c4762a1bSJed Brown }
164c4762a1bSJed Brown 
165c4762a1bSJed Brown /*TEST
166c4762a1bSJed Brown 
167c4762a1bSJed Brown    test:
168c4762a1bSJed Brown       args: -mat_block_size {{1 2  5 7 8}} -ov {{1 3}} -mat_size {{11 13}} -nd {{7}}
169c4762a1bSJed Brown 
170c4762a1bSJed Brown TEST*/
171