xref: /petsc/src/mat/tests/ex51.c (revision 4b5966dda8e0e9d52f06e77e04a44e922c8134a2)
1 
2 static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for MatBAIJ format.\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc, char **args)
7 {
8   Mat          A, B, E, *submatA, *submatB;
9   PetscInt     bs = 1, m = 43, ov = 1, i, j, k, *rows, *cols, M, nd = 5, *idx, mm, nn, lsize;
10   PetscScalar *vals, rval;
11   IS          *is1, *is2;
12   PetscRandom  rdm;
13   Vec          xx, s1, s2;
14   PetscReal    s1norm, s2norm, rnorm, tol = PETSC_SQRT_MACHINE_EPSILON;
15   PetscBool    flg;
16 
17   PetscFunctionBeginUser;
18   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
19   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
20   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL));
21   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL));
22   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL));
23   M = m * bs;
24 
25   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, bs, M, M, 1, NULL, &A));
26   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
27   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, M, M, 15, NULL, &B));
28   PetscCall(MatSetBlockSize(B, bs));
29   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
30   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
31   PetscCall(PetscRandomSetFromOptions(rdm));
32 
33   PetscCall(PetscMalloc1(bs, &rows));
34   PetscCall(PetscMalloc1(bs, &cols));
35   PetscCall(PetscMalloc1(bs * bs, &vals));
36   PetscCall(PetscMalloc1(M, &idx));
37 
38   /* Now set blocks of values */
39   for (i = 0; i < 20 * bs; i++) {
40     PetscInt nr = 1, nc = 1;
41     PetscCall(PetscRandomGetValue(rdm, &rval));
42     cols[0] = bs * (int)(PetscRealPart(rval) * m);
43     PetscCall(PetscRandomGetValue(rdm, &rval));
44     rows[0] = bs * (int)(PetscRealPart(rval) * m);
45     for (j = 1; j < bs; j++) {
46       PetscCall(PetscRandomGetValue(rdm, &rval));
47       if (PetscRealPart(rval) > .5) rows[nr++] = rows[0] + j - 1;
48     }
49     for (j = 1; j < bs; j++) {
50       PetscCall(PetscRandomGetValue(rdm, &rval));
51       if (PetscRealPart(rval) > .5) cols[nc++] = cols[0] + j - 1;
52     }
53 
54     for (j = 0; j < nr * nc; j++) {
55       PetscCall(PetscRandomGetValue(rdm, &rval));
56       vals[j] = rval;
57     }
58     PetscCall(MatSetValues(A, nr, rows, nc, cols, vals, ADD_VALUES));
59     PetscCall(MatSetValues(B, nr, rows, nc, cols, vals, ADD_VALUES));
60   }
61 
62   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
63   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
64   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
65   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
66 
67   /* Test MatConvert_SeqAIJ_SeqBAIJ handles incompletely filled blocks */
68   PetscCall(MatConvert(B, MATBAIJ, MAT_INITIAL_MATRIX, &E));
69   PetscCall(MatDestroy(&E));
70 
71   /* Test MatIncreaseOverlap() */
72   PetscCall(PetscMalloc1(nd, &is1));
73   PetscCall(PetscMalloc1(nd, &is2));
74 
75   for (i = 0; i < nd; i++) {
76     PetscCall(PetscRandomGetValue(rdm, &rval));
77     lsize = (int)(PetscRealPart(rval) * m);
78     for (j = 0; j < lsize; j++) {
79       PetscCall(PetscRandomGetValue(rdm, &rval));
80       idx[j * bs] = bs * (int)(PetscRealPart(rval) * m);
81       for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k;
82     }
83     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is1 + i));
84     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is2 + i));
85   }
86   PetscCall(MatIncreaseOverlap(A, nd, is1, ov));
87   PetscCall(MatIncreaseOverlap(B, nd, is2, ov));
88 
89   for (i = 0; i < nd; ++i) {
90     PetscCall(ISEqual(is1[i], is2[i], &flg));
91     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "i=%" PetscInt_FMT ", flg =%d", i, (int)flg);
92   }
93 
94   for (i = 0; i < nd; ++i) {
95     PetscCall(ISSort(is1[i]));
96     PetscCall(ISSort(is2[i]));
97   }
98 
99   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA));
100   PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_INITIAL_MATRIX, &submatB));
101 
102   /* Test MatMult() */
103   for (i = 0; i < nd; i++) {
104     PetscCall(MatGetSize(submatA[i], &mm, &nn));
105     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
106     PetscCall(VecDuplicate(xx, &s1));
107     PetscCall(VecDuplicate(xx, &s2));
108     for (j = 0; j < 3; j++) {
109       PetscCall(VecSetRandom(xx, rdm));
110       PetscCall(MatMult(submatA[i], xx, s1));
111       PetscCall(MatMult(submatB[i], xx, s2));
112       PetscCall(VecNorm(s1, NORM_2, &s1norm));
113       PetscCall(VecNorm(s2, NORM_2, &s2norm));
114       rnorm = s2norm - s1norm;
115       if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm));
116     }
117     PetscCall(VecDestroy(&xx));
118     PetscCall(VecDestroy(&s1));
119     PetscCall(VecDestroy(&s2));
120   }
121   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
122   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA));
123   PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_REUSE_MATRIX, &submatB));
124 
125   /* Test MatMult() */
126   for (i = 0; i < nd; i++) {
127     PetscCall(MatGetSize(submatA[i], &mm, &nn));
128     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
129     PetscCall(VecDuplicate(xx, &s1));
130     PetscCall(VecDuplicate(xx, &s2));
131     for (j = 0; j < 3; j++) {
132       PetscCall(VecSetRandom(xx, rdm));
133       PetscCall(MatMult(submatA[i], xx, s1));
134       PetscCall(MatMult(submatB[i], xx, s2));
135       PetscCall(VecNorm(s1, NORM_2, &s1norm));
136       PetscCall(VecNorm(s2, NORM_2, &s2norm));
137       rnorm = s2norm - s1norm;
138       if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm));
139     }
140     PetscCall(VecDestroy(&xx));
141     PetscCall(VecDestroy(&s1));
142     PetscCall(VecDestroy(&s2));
143   }
144 
145   /* Free allocated memory */
146   for (i = 0; i < nd; ++i) {
147     PetscCall(ISDestroy(&is1[i]));
148     PetscCall(ISDestroy(&is2[i]));
149   }
150   PetscCall(MatDestroySubMatrices(nd, &submatA));
151   PetscCall(MatDestroySubMatrices(nd, &submatB));
152   PetscCall(PetscFree(is1));
153   PetscCall(PetscFree(is2));
154   PetscCall(PetscFree(idx));
155   PetscCall(PetscFree(rows));
156   PetscCall(PetscFree(cols));
157   PetscCall(PetscFree(vals));
158   PetscCall(MatDestroy(&A));
159   PetscCall(MatDestroy(&B));
160   PetscCall(PetscRandomDestroy(&rdm));
161   PetscCall(PetscFinalize());
162   return 0;
163 }
164 
165 /*TEST
166 
167    test:
168       args: -mat_block_size {{1 2  5 7 8}} -ov {{1 3}} -mat_size {{11 13}} -nd {{7}}
169 
170 TEST*/
171