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