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