1 #include <../src/mat/impls/sbaij/seq/sbaij.h>
2 #include <petsc/private/kernels/blockinvert.h>
3
4 /* Version for when blocks are 3 by 3 */
MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat C,Mat A,const MatFactorInfo * info)5 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat C, Mat A, const MatFactorInfo *info)
6 {
7 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data;
8 IS perm = b->row;
9 const PetscInt *ai, *aj, *perm_ptr, mbs = a->mbs, *bi = b->i, *bj = b->j;
10 PetscInt *a2anew, i, j, k, k1, jmin, jmax, *jl, *il, vj, nexti, ili;
11 MatScalar *ba = b->a, *aa, *ap, *dk, *uik;
12 MatScalar *u, *diag, *rtmp, *rtmp_ptr;
13 PetscReal shift = info->shiftamount;
14 PetscBool allowzeropivot, zeropivotdetected;
15
16 PetscFunctionBegin;
17 /* initialization */
18 allowzeropivot = PetscNot(A->erroriffailure);
19 PetscCall(PetscCalloc1(9 * mbs, &rtmp));
20 PetscCall(PetscMalloc2(mbs, &il, mbs, &jl));
21 il[0] = 0;
22 for (i = 0; i < mbs; i++) jl[i] = mbs;
23
24 PetscCall(PetscMalloc2(9, &dk, 9, &uik));
25 PetscCall(ISGetIndices(perm, &perm_ptr));
26
27 /* check permutation */
28 if (!a->permute) {
29 ai = a->i;
30 aj = a->j;
31 aa = a->a;
32 } else {
33 ai = a->inew;
34 aj = a->jnew;
35 PetscCall(PetscMalloc1(9 * ai[mbs], &aa));
36 PetscCall(PetscArraycpy(aa, a->a, 9 * ai[mbs]));
37 PetscCall(PetscMalloc1(ai[mbs], &a2anew));
38 PetscCall(PetscArraycpy(a2anew, a->a2anew, ai[mbs]));
39
40 for (i = 0; i < mbs; i++) {
41 jmin = ai[i];
42 jmax = ai[i + 1];
43 for (j = jmin; j < jmax; j++) {
44 while (a2anew[j] != j) {
45 k = a2anew[j];
46 a2anew[j] = a2anew[k];
47 a2anew[k] = k;
48 for (k1 = 0; k1 < 9; k1++) {
49 dk[k1] = aa[k * 9 + k1];
50 aa[k * 9 + k1] = aa[j * 9 + k1];
51 aa[j * 9 + k1] = dk[k1];
52 }
53 }
54 /* transform column-oriented blocks that lie in the lower triangle to row-oriented blocks */
55 if (i > aj[j]) {
56 /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
57 ap = aa + j * 9; /* ptr to the beginning of j-th block of aa */
58 for (k = 0; k < 9; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
59 for (k = 0; k < 3; k++) { /* j-th block of aa <- dk^T */
60 for (k1 = 0; k1 < 3; k1++) *ap++ = dk[k + 3 * k1];
61 }
62 }
63 }
64 }
65 PetscCall(PetscFree(a2anew));
66 }
67
68 /* for each row k */
69 for (k = 0; k < mbs; k++) {
70 /*initialize k-th row with elements nonzero in row perm(k) of A */
71 jmin = ai[perm_ptr[k]];
72 jmax = ai[perm_ptr[k] + 1];
73 if (jmin < jmax) {
74 ap = aa + jmin * 9;
75 for (j = jmin; j < jmax; j++) {
76 vj = perm_ptr[aj[j]]; /* block col. index */
77 rtmp_ptr = rtmp + vj * 9;
78 for (i = 0; i < 9; i++) *rtmp_ptr++ = *ap++;
79 }
80 }
81
82 /* modify k-th row by adding in those rows i with U(i,k) != 0 */
83 PetscCall(PetscArraycpy(dk, rtmp + k * 9, 9));
84 i = jl[k]; /* first row to be added to k_th row */
85
86 while (i < mbs) {
87 nexti = jl[i]; /* next row to be added to k_th row */
88
89 /* compute multiplier */
90 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
91
92 /* uik = -inv(Di)*U_bar(i,k) */
93 diag = ba + i * 9;
94 u = ba + ili * 9;
95
96 uik[0] = -(diag[0] * u[0] + diag[3] * u[1] + diag[6] * u[2]);
97 uik[1] = -(diag[1] * u[0] + diag[4] * u[1] + diag[7] * u[2]);
98 uik[2] = -(diag[2] * u[0] + diag[5] * u[1] + diag[8] * u[2]);
99
100 uik[3] = -(diag[0] * u[3] + diag[3] * u[4] + diag[6] * u[5]);
101 uik[4] = -(diag[1] * u[3] + diag[4] * u[4] + diag[7] * u[5]);
102 uik[5] = -(diag[2] * u[3] + diag[5] * u[4] + diag[8] * u[5]);
103
104 uik[6] = -(diag[0] * u[6] + diag[3] * u[7] + diag[6] * u[8]);
105 uik[7] = -(diag[1] * u[6] + diag[4] * u[7] + diag[7] * u[8]);
106 uik[8] = -(diag[2] * u[6] + diag[5] * u[7] + diag[8] * u[8]);
107
108 /* update D(k) += -U(i,k)^T * U_bar(i,k) */
109 dk[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2];
110 dk[1] += uik[3] * u[0] + uik[4] * u[1] + uik[5] * u[2];
111 dk[2] += uik[6] * u[0] + uik[7] * u[1] + uik[8] * u[2];
112
113 dk[3] += uik[0] * u[3] + uik[1] * u[4] + uik[2] * u[5];
114 dk[4] += uik[3] * u[3] + uik[4] * u[4] + uik[5] * u[5];
115 dk[5] += uik[6] * u[3] + uik[7] * u[4] + uik[8] * u[5];
116
117 dk[6] += uik[0] * u[6] + uik[1] * u[7] + uik[2] * u[8];
118 dk[7] += uik[3] * u[6] + uik[4] * u[7] + uik[5] * u[8];
119 dk[8] += uik[6] * u[6] + uik[7] * u[7] + uik[8] * u[8];
120
121 PetscCall(PetscLogFlops(27.0 * 4.0));
122
123 /* update -U(i,k) */
124 PetscCall(PetscArraycpy(ba + ili * 9, uik, 9));
125
126 /* add multiple of row i to k-th row ... */
127 jmin = ili + 1;
128 jmax = bi[i + 1];
129 if (jmin < jmax) {
130 for (j = jmin; j < jmax; j++) {
131 /* rtmp += -U(i,k)^T * U_bar(i,j) */
132 rtmp_ptr = rtmp + bj[j] * 9;
133 u = ba + j * 9;
134 rtmp_ptr[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2];
135 rtmp_ptr[1] += uik[3] * u[0] + uik[4] * u[1] + uik[5] * u[2];
136 rtmp_ptr[2] += uik[6] * u[0] + uik[7] * u[1] + uik[8] * u[2];
137
138 rtmp_ptr[3] += uik[0] * u[3] + uik[1] * u[4] + uik[2] * u[5];
139 rtmp_ptr[4] += uik[3] * u[3] + uik[4] * u[4] + uik[5] * u[5];
140 rtmp_ptr[5] += uik[6] * u[3] + uik[7] * u[4] + uik[8] * u[5];
141
142 rtmp_ptr[6] += uik[0] * u[6] + uik[1] * u[7] + uik[2] * u[8];
143 rtmp_ptr[7] += uik[3] * u[6] + uik[4] * u[7] + uik[5] * u[8];
144 rtmp_ptr[8] += uik[6] * u[6] + uik[7] * u[7] + uik[8] * u[8];
145 }
146 PetscCall(PetscLogFlops(2.0 * 27.0 * (jmax - jmin)));
147
148 /* ... add i to row list for next nonzero entry */
149 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
150 j = bj[jmin];
151 jl[i] = jl[j];
152 jl[j] = i; /* update jl */
153 }
154 i = nexti;
155 }
156
157 /* save nonzero entries in k-th row of U ... */
158
159 /* invert diagonal block */
160 diag = ba + k * 9;
161 PetscCall(PetscArraycpy(diag, dk, 9));
162 PetscCall(PetscKernel_A_gets_inverse_A_3(diag, shift, allowzeropivot, &zeropivotdetected));
163 if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
164
165 jmin = bi[k];
166 jmax = bi[k + 1];
167 if (jmin < jmax) {
168 for (j = jmin; j < jmax; j++) {
169 vj = bj[j]; /* block col. index of U */
170 u = ba + j * 9;
171 rtmp_ptr = rtmp + vj * 9;
172 for (k1 = 0; k1 < 9; k1++) {
173 *u++ = *rtmp_ptr;
174 *rtmp_ptr++ = 0.0;
175 }
176 }
177
178 /* ... add k to row list for first nonzero entry in k-th row */
179 il[k] = jmin;
180 i = bj[jmin];
181 jl[k] = jl[i];
182 jl[i] = k;
183 }
184 }
185
186 PetscCall(PetscFree(rtmp));
187 PetscCall(PetscFree2(il, jl));
188 PetscCall(PetscFree2(dk, uik));
189 if (a->permute) PetscCall(PetscFree(aa));
190
191 PetscCall(ISRestoreIndices(perm, &perm_ptr));
192
193 C->ops->solve = MatSolve_SeqSBAIJ_3_inplace;
194 C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_inplace;
195 C->assembled = PETSC_TRUE;
196 C->preallocated = PETSC_TRUE;
197
198 PetscCall(PetscLogFlops(1.3333 * 27 * b->mbs)); /* from inverting diagonal blocks */
199 PetscFunctionReturn(PETSC_SUCCESS);
200 }
201