xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact10.c (revision e8e8640d1cb9a3a2f50c0c0d7b26e5c4d521e2e4)
1 #include <../src/mat/impls/sbaij/seq/sbaij.h>
2 #include <petsc/private/kernels/blockinvert.h>
3 
4 /*
5       Version for when blocks are 6 by 6 Using natural ordering
6 */
MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat C,Mat A,const MatFactorInfo * info)7 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering(Mat C, Mat A, const MatFactorInfo *info)
8 {
9   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data;
10   PetscInt      i, j, mbs = a->mbs, *bi = b->i, *bj = b->j;
11   PetscInt     *ai, *aj, k, k1, jmin, jmax, *jl, *il, vj, nexti, ili;
12   MatScalar    *ba = b->a, *aa, *ap, *dk, *uik;
13   MatScalar    *u, *d, *w, *wp, u0, u1, u2, u3, u4, u5, u6, u7, u8, u9, u10, u11, u12;
14   MatScalar     u13, u14, u15, u16, u17, u18, u19, u20, u21, u22, u23, u24, u25, u26, u27;
15   MatScalar     u28, u29, u30, u31, u32, u33, u34, u35;
16   MatScalar     d0, d1, d2, d3, d4, d5, d6, d7, d8, d9, d10, d11, d12;
17   MatScalar     d13, d14, d15, d16, d17, d18, d19, d20, d21, d22, d23, d24, d25, d26, d27;
18   MatScalar     d28, d29, d30, d31, d32, d33, d34, d35;
19   MatScalar     m0, m1, m2, m3, m4, m5, m6, m7, m8, m9, m10, m11, m12;
20   MatScalar     m13, m14, m15, m16, m17, m18, m19, m20, m21, m22, m23, m24, m25, m26, m27;
21   MatScalar     m28, m29, m30, m31, m32, m33, m34, m35;
22   PetscReal     shift = info->shiftamount;
23   PetscBool     allowzeropivot, zeropivotdetected;
24 
25   PetscFunctionBegin;
26   /* initialization */
27   allowzeropivot = PetscNot(A->erroriffailure);
28   PetscCall(PetscCalloc1(36 * mbs, &w));
29   PetscCall(PetscMalloc2(mbs, &il, mbs, &jl));
30   il[0] = 0;
31   for (i = 0; i < mbs; i++) jl[i] = mbs;
32 
33   PetscCall(PetscMalloc2(36, &dk, 36, &uik));
34   ai = a->i;
35   aj = a->j;
36   aa = a->a;
37 
38   /* for each row k */
39   for (k = 0; k < mbs; k++) {
40     /*initialize k-th row with elements nonzero in row k of A */
41     jmin = ai[k];
42     jmax = ai[k + 1];
43     if (jmin < jmax) {
44       ap = aa + jmin * 36;
45       for (j = jmin; j < jmax; j++) {
46         vj = aj[j]; /* block col. index */
47         wp = w + vj * 36;
48         for (i = 0; i < 36; i++) *wp++ = *ap++;
49       }
50     }
51 
52     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
53     PetscCall(PetscArraycpy(dk, w + k * 36, 36));
54     i = jl[k]; /* first row to be added to k_th row  */
55 
56     while (i < mbs) {
57       nexti = jl[i]; /* next row to be added to k_th row */
58 
59       /* compute multiplier */
60       ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
61 
62       /* uik = -inv(Di)*U_bar(i,k) */
63       d = ba + i * 36;
64       u = ba + ili * 36;
65 
66       u0  = u[0];
67       u1  = u[1];
68       u2  = u[2];
69       u3  = u[3];
70       u4  = u[4];
71       u5  = u[5];
72       u6  = u[6];
73       u7  = u[7];
74       u8  = u[8];
75       u9  = u[9];
76       u10 = u[10];
77       u11 = u[11];
78       u12 = u[12];
79       u13 = u[13];
80       u14 = u[14];
81       u15 = u[15];
82       u16 = u[16];
83       u17 = u[17];
84       u18 = u[18];
85       u19 = u[19];
86       u20 = u[20];
87       u21 = u[21];
88       u22 = u[22];
89       u23 = u[23];
90       u24 = u[24];
91       u25 = u[25];
92       u26 = u[26];
93       u27 = u[27];
94       u28 = u[28];
95       u29 = u[29];
96       u30 = u[30];
97       u31 = u[31];
98       u32 = u[32];
99       u33 = u[33];
100       u34 = u[34];
101       u35 = u[35];
102 
103       d0  = d[0];
104       d1  = d[1];
105       d2  = d[2];
106       d3  = d[3];
107       d4  = d[4];
108       d5  = d[5];
109       d6  = d[6];
110       d7  = d[7];
111       d8  = d[8];
112       d9  = d[9];
113       d10 = d[10];
114       d11 = d[11];
115       d12 = d[12];
116       d13 = d[13];
117       d14 = d[14];
118       d15 = d[15];
119       d16 = d[16];
120       d17 = d[17];
121       d18 = d[18];
122       d19 = d[19];
123       d20 = d[20];
124       d21 = d[21];
125       d22 = d[22];
126       d23 = d[23];
127       d24 = d[24];
128       d25 = d[25];
129       d26 = d[26];
130       d27 = d[27];
131       d28 = d[28];
132       d29 = d[29];
133       d30 = d[30];
134       d31 = d[31];
135       d32 = d[32];
136       d33 = d[33];
137       d34 = d[34];
138       d35 = d[35];
139 
140       m0 = uik[0] = -(d0 * u0 + d6 * u1 + d12 * u2 + d18 * u3 + d24 * u4 + d30 * u5);
141       m1 = uik[1] = -(d1 * u0 + d7 * u1 + d13 * u2 + d19 * u3 + d25 * u4 + d31 * u5);
142       m2 = uik[2] = -(d2 * u0 + d8 * u1 + d14 * u2 + d20 * u3 + d26 * u4 + d32 * u5);
143       m3 = uik[3] = -(d3 * u0 + d9 * u1 + d15 * u2 + d21 * u3 + d27 * u4 + d33 * u5);
144       m4 = uik[4] = -(d4 * u0 + d10 * u1 + d16 * u2 + d22 * u3 + d28 * u4 + d34 * u5);
145       m5 = uik[5] = -(d5 * u0 + d11 * u1 + d17 * u2 + d23 * u3 + d29 * u4 + d35 * u5);
146 
147       m6 = uik[6] = -(d0 * u6 + d6 * u7 + d12 * u8 + d18 * u9 + d24 * u10 + d30 * u11);
148       m7 = uik[7] = -(d1 * u6 + d7 * u7 + d13 * u8 + d19 * u9 + d25 * u10 + d31 * u11);
149       m8 = uik[8] = -(d2 * u6 + d8 * u7 + d14 * u8 + d20 * u9 + d26 * u10 + d32 * u11);
150       m9 = uik[9] = -(d3 * u6 + d9 * u7 + d15 * u8 + d21 * u9 + d27 * u10 + d33 * u11);
151       m10 = uik[10] = -(d4 * u6 + d10 * u7 + d16 * u8 + d22 * u9 + d28 * u10 + d34 * u11);
152       m11 = uik[11] = -(d5 * u6 + d11 * u7 + d17 * u8 + d23 * u9 + d29 * u10 + d35 * u11);
153 
154       m12 = uik[12] = -(d0 * u12 + d6 * u13 + d12 * u14 + d18 * u15 + d24 * u16 + d30 * u17);
155       m13 = uik[13] = -(d1 * u12 + d7 * u13 + d13 * u14 + d19 * u15 + d25 * u16 + d31 * u17);
156       m14 = uik[14] = -(d2 * u12 + d8 * u13 + d14 * u14 + d20 * u15 + d26 * u16 + d32 * u17);
157       m15 = uik[15] = -(d3 * u12 + d9 * u13 + d15 * u14 + d21 * u15 + d27 * u16 + d33 * u17);
158       m16 = uik[16] = -(d4 * u12 + d10 * u13 + d16 * u14 + d22 * u15 + d28 * u16 + d34 * u17);
159       m17 = uik[17] = -(d5 * u12 + d11 * u13 + d17 * u14 + d23 * u15 + d29 * u16 + d35 * u17);
160 
161       m18 = uik[18] = -(d0 * u18 + d6 * u19 + d12 * u20 + d18 * u21 + d24 * u22 + d30 * u23);
162       m19 = uik[19] = -(d1 * u18 + d7 * u19 + d13 * u20 + d19 * u21 + d25 * u22 + d31 * u23);
163       m20 = uik[20] = -(d2 * u18 + d8 * u19 + d14 * u20 + d20 * u21 + d26 * u22 + d32 * u23);
164       m21 = uik[21] = -(d3 * u18 + d9 * u19 + d15 * u20 + d21 * u21 + d27 * u22 + d33 * u23);
165       m22 = uik[22] = -(d4 * u18 + d10 * u19 + d16 * u20 + d22 * u21 + d28 * u22 + d34 * u23);
166       m23 = uik[23] = -(d5 * u18 + d11 * u19 + d17 * u20 + d23 * u21 + d29 * u22 + d35 * u23);
167 
168       m24 = uik[24] = -(d0 * u24 + d6 * u25 + d12 * u26 + d18 * u27 + d24 * u28 + d30 * u29);
169       m25 = uik[25] = -(d1 * u24 + d7 * u25 + d13 * u26 + d19 * u27 + d25 * u28 + d31 * u29);
170       m26 = uik[26] = -(d2 * u24 + d8 * u25 + d14 * u26 + d20 * u27 + d26 * u28 + d32 * u29);
171       m27 = uik[27] = -(d3 * u24 + d9 * u25 + d15 * u26 + d21 * u27 + d27 * u28 + d33 * u29);
172       m28 = uik[28] = -(d4 * u24 + d10 * u25 + d16 * u26 + d22 * u27 + d28 * u28 + d34 * u29);
173       m29 = uik[29] = -(d5 * u24 + d11 * u25 + d17 * u26 + d23 * u27 + d29 * u28 + d35 * u29);
174 
175       m30 = uik[30] = -(d0 * u30 + d6 * u31 + d12 * u32 + d18 * u33 + d24 * u34 + d30 * u35);
176       m31 = uik[31] = -(d1 * u30 + d7 * u31 + d13 * u32 + d19 * u33 + d25 * u34 + d31 * u35);
177       m32 = uik[32] = -(d2 * u30 + d8 * u31 + d14 * u32 + d20 * u33 + d26 * u34 + d32 * u35);
178       m33 = uik[33] = -(d3 * u30 + d9 * u31 + d15 * u32 + d21 * u33 + d27 * u34 + d33 * u35);
179       m34 = uik[34] = -(d4 * u30 + d10 * u31 + d16 * u32 + d22 * u33 + d28 * u34 + d34 * u35);
180       m35 = uik[35] = -(d5 * u30 + d11 * u31 + d17 * u32 + d23 * u33 + d29 * u34 + d35 * u35);
181 
182       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
183       dk[0] += m0 * u0 + m1 * u1 + m2 * u2 + m3 * u3 + m4 * u4 + m5 * u5;
184       dk[1] += m6 * u0 + m7 * u1 + m8 * u2 + m9 * u3 + m10 * u4 + m11 * u5;
185       dk[2] += m12 * u0 + m13 * u1 + m14 * u2 + m15 * u3 + m16 * u4 + m17 * u5;
186       dk[3] += m18 * u0 + m19 * u1 + m20 * u2 + m21 * u3 + m22 * u4 + m23 * u5;
187       dk[4] += m24 * u0 + m25 * u1 + m26 * u2 + m27 * u3 + m28 * u4 + m29 * u5;
188       dk[5] += m30 * u0 + m31 * u1 + m32 * u2 + m33 * u3 + m34 * u4 + m35 * u5;
189 
190       dk[6] += m0 * u6 + m1 * u7 + m2 * u8 + m3 * u9 + m4 * u10 + m5 * u11;
191       dk[7] += m6 * u6 + m7 * u7 + m8 * u8 + m9 * u9 + m10 * u10 + m11 * u11;
192       dk[8] += m12 * u6 + m13 * u7 + m14 * u8 + m15 * u9 + m16 * u10 + m17 * u11;
193       dk[9] += m18 * u6 + m19 * u7 + m20 * u8 + m21 * u9 + m22 * u10 + m23 * u11;
194       dk[10] += m24 * u6 + m25 * u7 + m26 * u8 + m27 * u9 + m28 * u10 + m29 * u11;
195       dk[11] += m30 * u6 + m31 * u7 + m32 * u8 + m33 * u9 + m34 * u10 + m35 * u11;
196 
197       dk[12] += m0 * u12 + m1 * u13 + m2 * u14 + m3 * u15 + m4 * u16 + m5 * u17;
198       dk[13] += m6 * u12 + m7 * u13 + m8 * u14 + m9 * u15 + m10 * u16 + m11 * u17;
199       dk[14] += m12 * u12 + m13 * u13 + m14 * u14 + m15 * u15 + m16 * u16 + m17 * u17;
200       dk[15] += m18 * u12 + m19 * u13 + m20 * u14 + m21 * u15 + m22 * u16 + m23 * u17;
201       dk[16] += m24 * u12 + m25 * u13 + m26 * u14 + m27 * u15 + m28 * u16 + m29 * u17;
202       dk[17] += m30 * u12 + m31 * u13 + m32 * u14 + m33 * u15 + m34 * u16 + m35 * u17;
203 
204       dk[18] += m0 * u18 + m1 * u19 + m2 * u20 + m3 * u21 + m4 * u22 + m5 * u23;
205       dk[19] += m6 * u18 + m7 * u19 + m8 * u20 + m9 * u21 + m10 * u22 + m11 * u23;
206       dk[20] += m12 * u18 + m13 * u19 + m14 * u20 + m15 * u21 + m16 * u22 + m17 * u23;
207       dk[21] += m18 * u18 + m19 * u19 + m20 * u20 + m21 * u21 + m22 * u22 + m23 * u23;
208       dk[22] += m24 * u18 + m25 * u19 + m26 * u20 + m27 * u21 + m28 * u22 + m29 * u23;
209       dk[23] += m30 * u18 + m31 * u19 + m32 * u20 + m33 * u21 + m34 * u22 + m35 * u23;
210 
211       dk[24] += m0 * u24 + m1 * u25 + m2 * u26 + m3 * u27 + m4 * u28 + m5 * u29;
212       dk[25] += m6 * u24 + m7 * u25 + m8 * u26 + m9 * u27 + m10 * u28 + m11 * u29;
213       dk[26] += m12 * u24 + m13 * u25 + m14 * u26 + m15 * u27 + m16 * u28 + m17 * u29;
214       dk[27] += m18 * u24 + m19 * u25 + m20 * u26 + m21 * u27 + m22 * u28 + m23 * u29;
215       dk[28] += m24 * u24 + m25 * u25 + m26 * u26 + m27 * u27 + m28 * u28 + m29 * u29;
216       dk[29] += m30 * u24 + m31 * u25 + m32 * u26 + m33 * u27 + m34 * u28 + m35 * u29;
217 
218       dk[30] += m0 * u30 + m1 * u31 + m2 * u32 + m3 * u33 + m4 * u34 + m5 * u35;
219       dk[31] += m6 * u30 + m7 * u31 + m8 * u32 + m9 * u33 + m10 * u34 + m11 * u35;
220       dk[32] += m12 * u30 + m13 * u31 + m14 * u32 + m15 * u33 + m16 * u34 + m17 * u35;
221       dk[33] += m18 * u30 + m19 * u31 + m20 * u32 + m21 * u33 + m22 * u34 + m23 * u35;
222       dk[34] += m24 * u30 + m25 * u31 + m26 * u32 + m27 * u33 + m28 * u34 + m29 * u35;
223       dk[35] += m30 * u30 + m31 * u31 + m32 * u32 + m33 * u33 + m34 * u34 + m35 * u35;
224 
225       PetscCall(PetscLogFlops(216.0 * 4.0));
226 
227       /* update -U(i,k) */
228       PetscCall(PetscArraycpy(ba + ili * 36, uik, 36));
229 
230       /* add multiple of row i to k-th row ... */
231       jmin = ili + 1;
232       jmax = bi[i + 1];
233       if (jmin < jmax) {
234         for (j = jmin; j < jmax; j++) {
235           /* w += -U(i,k)^T * U_bar(i,j) */
236           wp = w + bj[j] * 36;
237           u  = ba + j * 36;
238 
239           u0  = u[0];
240           u1  = u[1];
241           u2  = u[2];
242           u3  = u[3];
243           u4  = u[4];
244           u5  = u[5];
245           u6  = u[6];
246           u7  = u[7];
247           u8  = u[8];
248           u9  = u[9];
249           u10 = u[10];
250           u11 = u[11];
251           u12 = u[12];
252           u13 = u[13];
253           u14 = u[14];
254           u15 = u[15];
255           u16 = u[16];
256           u17 = u[17];
257           u18 = u[18];
258           u19 = u[19];
259           u20 = u[20];
260           u21 = u[21];
261           u22 = u[22];
262           u23 = u[23];
263           u24 = u[24];
264           u25 = u[25];
265           u26 = u[26];
266           u27 = u[27];
267           u28 = u[28];
268           u29 = u[29];
269           u30 = u[30];
270           u31 = u[31];
271           u32 = u[32];
272           u33 = u[33];
273           u34 = u[34];
274           u35 = u[35];
275 
276           wp[0] += m0 * u0 + m1 * u1 + m2 * u2 + m3 * u3 + m4 * u4 + m5 * u5;
277           wp[1] += m6 * u0 + m7 * u1 + m8 * u2 + m9 * u3 + m10 * u4 + m11 * u5;
278           wp[2] += m12 * u0 + m13 * u1 + m14 * u2 + m15 * u3 + m16 * u4 + m17 * u5;
279           wp[3] += m18 * u0 + m19 * u1 + m20 * u2 + m21 * u3 + m22 * u4 + m23 * u5;
280           wp[4] += m24 * u0 + m25 * u1 + m26 * u2 + m27 * u3 + m28 * u4 + m29 * u5;
281           wp[5] += m30 * u0 + m31 * u1 + m32 * u2 + m33 * u3 + m34 * u4 + m35 * u5;
282 
283           wp[6] += m0 * u6 + m1 * u7 + m2 * u8 + m3 * u9 + m4 * u10 + m5 * u11;
284           wp[7] += m6 * u6 + m7 * u7 + m8 * u8 + m9 * u9 + m10 * u10 + m11 * u11;
285           wp[8] += m12 * u6 + m13 * u7 + m14 * u8 + m15 * u9 + m16 * u10 + m17 * u11;
286           wp[9] += m18 * u6 + m19 * u7 + m20 * u8 + m21 * u9 + m22 * u10 + m23 * u11;
287           wp[10] += m24 * u6 + m25 * u7 + m26 * u8 + m27 * u9 + m28 * u10 + m29 * u11;
288           wp[11] += m30 * u6 + m31 * u7 + m32 * u8 + m33 * u9 + m34 * u10 + m35 * u11;
289 
290           wp[12] += m0 * u12 + m1 * u13 + m2 * u14 + m3 * u15 + m4 * u16 + m5 * u17;
291           wp[13] += m6 * u12 + m7 * u13 + m8 * u14 + m9 * u15 + m10 * u16 + m11 * u17;
292           wp[14] += m12 * u12 + m13 * u13 + m14 * u14 + m15 * u15 + m16 * u16 + m17 * u17;
293           wp[15] += m18 * u12 + m19 * u13 + m20 * u14 + m21 * u15 + m22 * u16 + m23 * u17;
294           wp[16] += m24 * u12 + m25 * u13 + m26 * u14 + m27 * u15 + m28 * u16 + m29 * u17;
295           wp[17] += m30 * u12 + m31 * u13 + m32 * u14 + m33 * u15 + m34 * u16 + m35 * u17;
296 
297           wp[18] += m0 * u18 + m1 * u19 + m2 * u20 + m3 * u21 + m4 * u22 + m5 * u23;
298           wp[19] += m6 * u18 + m7 * u19 + m8 * u20 + m9 * u21 + m10 * u22 + m11 * u23;
299           wp[20] += m12 * u18 + m13 * u19 + m14 * u20 + m15 * u21 + m16 * u22 + m17 * u23;
300           wp[21] += m18 * u18 + m19 * u19 + m20 * u20 + m21 * u21 + m22 * u22 + m23 * u23;
301           wp[22] += m24 * u18 + m25 * u19 + m26 * u20 + m27 * u21 + m28 * u22 + m29 * u23;
302           wp[23] += m30 * u18 + m31 * u19 + m32 * u20 + m33 * u21 + m34 * u22 + m35 * u23;
303 
304           wp[24] += m0 * u24 + m1 * u25 + m2 * u26 + m3 * u27 + m4 * u28 + m5 * u29;
305           wp[25] += m6 * u24 + m7 * u25 + m8 * u26 + m9 * u27 + m10 * u28 + m11 * u29;
306           wp[26] += m12 * u24 + m13 * u25 + m14 * u26 + m15 * u27 + m16 * u28 + m17 * u29;
307           wp[27] += m18 * u24 + m19 * u25 + m20 * u26 + m21 * u27 + m22 * u28 + m23 * u29;
308           wp[28] += m24 * u24 + m25 * u25 + m26 * u26 + m27 * u27 + m28 * u28 + m29 * u29;
309           wp[29] += m30 * u24 + m31 * u25 + m32 * u26 + m33 * u27 + m34 * u28 + m35 * u29;
310 
311           wp[30] += m0 * u30 + m1 * u31 + m2 * u32 + m3 * u33 + m4 * u34 + m5 * u35;
312           wp[31] += m6 * u30 + m7 * u31 + m8 * u32 + m9 * u33 + m10 * u34 + m11 * u35;
313           wp[32] += m12 * u30 + m13 * u31 + m14 * u32 + m15 * u33 + m16 * u34 + m17 * u35;
314           wp[33] += m18 * u30 + m19 * u31 + m20 * u32 + m21 * u33 + m22 * u34 + m23 * u35;
315           wp[34] += m24 * u30 + m25 * u31 + m26 * u32 + m27 * u33 + m28 * u34 + m29 * u35;
316           wp[35] += m30 * u30 + m31 * u31 + m32 * u32 + m33 * u33 + m34 * u34 + m35 * u35;
317         }
318         PetscCall(PetscLogFlops(2.0 * 216.0 * (jmax - jmin)));
319 
320         /* ... add i to row list for next nonzero entry */
321         il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
322         j     = bj[jmin];
323         jl[i] = jl[j];
324         jl[j] = i; /* update jl */
325       }
326       i = nexti;
327     }
328 
329     /* save nonzero entries in k-th row of U ... */
330 
331     /* invert diagonal block */
332     d = ba + k * 36;
333     PetscCall(PetscArraycpy(d, dk, 36));
334     PetscCall(PetscKernel_A_gets_inverse_A_6(d, shift, allowzeropivot, &zeropivotdetected));
335     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
336 
337     jmin = bi[k];
338     jmax = bi[k + 1];
339     if (jmin < jmax) {
340       for (j = jmin; j < jmax; j++) {
341         vj = bj[j]; /* block col. index of U */
342         u  = ba + j * 36;
343         wp = w + vj * 36;
344         for (k1 = 0; k1 < 36; k1++) {
345           *u++  = *wp;
346           *wp++ = 0.0;
347         }
348       }
349 
350       /* ... add k to row list for first nonzero entry in k-th row */
351       il[k] = jmin;
352       i     = bj[jmin];
353       jl[k] = jl[i];
354       jl[i] = k;
355     }
356   }
357 
358   PetscCall(PetscFree(w));
359   PetscCall(PetscFree2(il, jl));
360   PetscCall(PetscFree2(dk, uik));
361 
362   C->ops->solve          = MatSolve_SeqSBAIJ_6_NaturalOrdering_inplace;
363   C->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering_inplace;
364   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace;
365   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering_inplace;
366   C->assembled           = PETSC_TRUE;
367   C->preallocated        = PETSC_TRUE;
368 
369   PetscCall(PetscLogFlops(1.3333 * 216 * b->mbs)); /* from inverting diagonal blocks */
370   PetscFunctionReturn(PETSC_SUCCESS);
371 }
372