xref: /petsc/src/mat/impls/baij/seq/baijsolvnat6.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h>
32c733ed4SBarry Smith 
4*9371c9d4SSatish Balay PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx) {
52c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
62c733ed4SBarry Smith   PetscInt           i, nz, idx, idt, jdx;
72c733ed4SBarry Smith   const PetscInt    *diag = a->diag, *vi, n = a->mbs, *ai = a->i, *aj = a->j;
82c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
92c733ed4SBarry Smith   PetscScalar       *x, s1, s2, s3, s4, s5, s6, x1, x2, x3, x4, x5, x6;
102c733ed4SBarry Smith   const PetscScalar *b;
112c733ed4SBarry Smith 
122c733ed4SBarry Smith   PetscFunctionBegin;
139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
149566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
152c733ed4SBarry Smith   /* forward solve the lower triangular */
162c733ed4SBarry Smith   idx  = 0;
17*9371c9d4SSatish Balay   x[0] = b[idx];
18*9371c9d4SSatish Balay   x[1] = b[1 + idx];
19*9371c9d4SSatish Balay   x[2] = b[2 + idx];
20*9371c9d4SSatish Balay   x[3] = b[3 + idx];
21*9371c9d4SSatish Balay   x[4] = b[4 + idx];
22*9371c9d4SSatish Balay   x[5] = b[5 + idx];
232c733ed4SBarry Smith   for (i = 1; i < n; i++) {
242c733ed4SBarry Smith     v   = aa + 36 * ai[i];
252c733ed4SBarry Smith     vi  = aj + ai[i];
262c733ed4SBarry Smith     nz  = diag[i] - ai[i];
272c733ed4SBarry Smith     idx = 6 * i;
28*9371c9d4SSatish Balay     s1  = b[idx];
29*9371c9d4SSatish Balay     s2  = b[1 + idx];
30*9371c9d4SSatish Balay     s3  = b[2 + idx];
31*9371c9d4SSatish Balay     s4  = b[3 + idx];
32*9371c9d4SSatish Balay     s5  = b[4 + idx];
33*9371c9d4SSatish Balay     s6  = b[5 + idx];
342c733ed4SBarry Smith     while (nz--) {
352c733ed4SBarry Smith       jdx = 6 * (*vi++);
36*9371c9d4SSatish Balay       x1  = x[jdx];
37*9371c9d4SSatish Balay       x2  = x[1 + jdx];
38*9371c9d4SSatish Balay       x3  = x[2 + jdx];
39*9371c9d4SSatish Balay       x4  = x[3 + jdx];
40*9371c9d4SSatish Balay       x5  = x[4 + jdx];
41*9371c9d4SSatish Balay       x6  = x[5 + jdx];
422c733ed4SBarry Smith       s1 -= v[0] * x1 + v[6] * x2 + v[12] * x3 + v[18] * x4 + v[24] * x5 + v[30] * x6;
432c733ed4SBarry Smith       s2 -= v[1] * x1 + v[7] * x2 + v[13] * x3 + v[19] * x4 + v[25] * x5 + v[31] * x6;
442c733ed4SBarry Smith       s3 -= v[2] * x1 + v[8] * x2 + v[14] * x3 + v[20] * x4 + v[26] * x5 + v[32] * x6;
452c733ed4SBarry Smith       s4 -= v[3] * x1 + v[9] * x2 + v[15] * x3 + v[21] * x4 + v[27] * x5 + v[33] * x6;
462c733ed4SBarry Smith       s5 -= v[4] * x1 + v[10] * x2 + v[16] * x3 + v[22] * x4 + v[28] * x5 + v[34] * x6;
472c733ed4SBarry Smith       s6 -= v[5] * x1 + v[11] * x2 + v[17] * x3 + v[23] * x4 + v[29] * x5 + v[35] * x6;
482c733ed4SBarry Smith       v += 36;
492c733ed4SBarry Smith     }
502c733ed4SBarry Smith     x[idx]     = s1;
512c733ed4SBarry Smith     x[1 + idx] = s2;
522c733ed4SBarry Smith     x[2 + idx] = s3;
532c733ed4SBarry Smith     x[3 + idx] = s4;
542c733ed4SBarry Smith     x[4 + idx] = s5;
552c733ed4SBarry Smith     x[5 + idx] = s6;
562c733ed4SBarry Smith   }
572c733ed4SBarry Smith   /* backward solve the upper triangular */
582c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
592c733ed4SBarry Smith     v   = aa + 36 * diag[i] + 36;
602c733ed4SBarry Smith     vi  = aj + diag[i] + 1;
612c733ed4SBarry Smith     nz  = ai[i + 1] - diag[i] - 1;
622c733ed4SBarry Smith     idt = 6 * i;
63*9371c9d4SSatish Balay     s1  = x[idt];
64*9371c9d4SSatish Balay     s2  = x[1 + idt];
65*9371c9d4SSatish Balay     s3  = x[2 + idt];
66*9371c9d4SSatish Balay     s4  = x[3 + idt];
67*9371c9d4SSatish Balay     s5  = x[4 + idt];
68*9371c9d4SSatish Balay     s6  = x[5 + idt];
692c733ed4SBarry Smith     while (nz--) {
702c733ed4SBarry Smith       idx = 6 * (*vi++);
71*9371c9d4SSatish Balay       x1  = x[idx];
72*9371c9d4SSatish Balay       x2  = x[1 + idx];
73*9371c9d4SSatish Balay       x3  = x[2 + idx];
74*9371c9d4SSatish Balay       x4  = x[3 + idx];
75*9371c9d4SSatish Balay       x5  = x[4 + idx];
76*9371c9d4SSatish Balay       x6  = x[5 + idx];
772c733ed4SBarry Smith       s1 -= v[0] * x1 + v[6] * x2 + v[12] * x3 + v[18] * x4 + v[24] * x5 + v[30] * x6;
782c733ed4SBarry Smith       s2 -= v[1] * x1 + v[7] * x2 + v[13] * x3 + v[19] * x4 + v[25] * x5 + v[31] * x6;
792c733ed4SBarry Smith       s3 -= v[2] * x1 + v[8] * x2 + v[14] * x3 + v[20] * x4 + v[26] * x5 + v[32] * x6;
802c733ed4SBarry Smith       s4 -= v[3] * x1 + v[9] * x2 + v[15] * x3 + v[21] * x4 + v[27] * x5 + v[33] * x6;
812c733ed4SBarry Smith       s5 -= v[4] * x1 + v[10] * x2 + v[16] * x3 + v[22] * x4 + v[28] * x5 + v[34] * x6;
822c733ed4SBarry Smith       s6 -= v[5] * x1 + v[11] * x2 + v[17] * x3 + v[23] * x4 + v[29] * x5 + v[35] * x6;
832c733ed4SBarry Smith       v += 36;
842c733ed4SBarry Smith     }
852c733ed4SBarry Smith     v          = aa + 36 * diag[i];
862c733ed4SBarry Smith     x[idt]     = v[0] * s1 + v[6] * s2 + v[12] * s3 + v[18] * s4 + v[24] * s5 + v[30] * s6;
872c733ed4SBarry Smith     x[1 + idt] = v[1] * s1 + v[7] * s2 + v[13] * s3 + v[19] * s4 + v[25] * s5 + v[31] * s6;
882c733ed4SBarry Smith     x[2 + idt] = v[2] * s1 + v[8] * s2 + v[14] * s3 + v[20] * s4 + v[26] * s5 + v[32] * s6;
892c733ed4SBarry Smith     x[3 + idt] = v[3] * s1 + v[9] * s2 + v[15] * s3 + v[21] * s4 + v[27] * s5 + v[33] * s6;
902c733ed4SBarry Smith     x[4 + idt] = v[4] * s1 + v[10] * s2 + v[16] * s3 + v[22] * s4 + v[28] * s5 + v[34] * s6;
912c733ed4SBarry Smith     x[5 + idt] = v[5] * s1 + v[11] * s2 + v[17] * s3 + v[23] * s4 + v[29] * s5 + v[35] * s6;
922c733ed4SBarry Smith   }
932c733ed4SBarry Smith 
949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
969566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 36 * (a->nz) - 6.0 * A->cmap->n));
972c733ed4SBarry Smith   PetscFunctionReturn(0);
982c733ed4SBarry Smith }
992c733ed4SBarry Smith 
100*9371c9d4SSatish Balay PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering(Mat A, Vec bb, Vec xx) {
1012c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
1022c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *ai = a->i, *aj = a->j, *adiag = a->diag;
1032c733ed4SBarry Smith   PetscInt           i, k, nz, idx, jdx, idt;
1042c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
1052c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
1062c733ed4SBarry Smith   PetscScalar       *x;
1072c733ed4SBarry Smith   const PetscScalar *b;
1082c733ed4SBarry Smith   PetscScalar        s1, s2, s3, s4, s5, s6, x1, x2, x3, x4, x5, x6;
1092c733ed4SBarry Smith 
1102c733ed4SBarry Smith   PetscFunctionBegin;
1119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1129566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1132c733ed4SBarry Smith   /* forward solve the lower triangular */
1142c733ed4SBarry Smith   idx  = 0;
115*9371c9d4SSatish Balay   x[0] = b[idx];
116*9371c9d4SSatish Balay   x[1] = b[1 + idx];
117*9371c9d4SSatish Balay   x[2] = b[2 + idx];
118*9371c9d4SSatish Balay   x[3] = b[3 + idx];
119*9371c9d4SSatish Balay   x[4] = b[4 + idx];
120*9371c9d4SSatish Balay   x[5] = b[5 + idx];
1212c733ed4SBarry Smith   for (i = 1; i < n; i++) {
1222c733ed4SBarry Smith     v   = aa + bs2 * ai[i];
1232c733ed4SBarry Smith     vi  = aj + ai[i];
1242c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1252c733ed4SBarry Smith     idx = bs * i;
126*9371c9d4SSatish Balay     s1  = b[idx];
127*9371c9d4SSatish Balay     s2  = b[1 + idx];
128*9371c9d4SSatish Balay     s3  = b[2 + idx];
129*9371c9d4SSatish Balay     s4  = b[3 + idx];
130*9371c9d4SSatish Balay     s5  = b[4 + idx];
131*9371c9d4SSatish Balay     s6  = b[5 + idx];
1322c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
1332c733ed4SBarry Smith       jdx = bs * vi[k];
134*9371c9d4SSatish Balay       x1  = x[jdx];
135*9371c9d4SSatish Balay       x2  = x[1 + jdx];
136*9371c9d4SSatish Balay       x3  = x[2 + jdx];
137*9371c9d4SSatish Balay       x4  = x[3 + jdx];
138*9371c9d4SSatish Balay       x5  = x[4 + jdx];
139*9371c9d4SSatish Balay       x6  = x[5 + jdx];
1402c733ed4SBarry Smith       s1 -= v[0] * x1 + v[6] * x2 + v[12] * x3 + v[18] * x4 + v[24] * x5 + v[30] * x6;
141feb237baSPierre Jolivet       s2 -= v[1] * x1 + v[7] * x2 + v[13] * x3 + v[19] * x4 + v[25] * x5 + v[31] * x6;
1422c733ed4SBarry Smith       s3 -= v[2] * x1 + v[8] * x2 + v[14] * x3 + v[20] * x4 + v[26] * x5 + v[32] * x6;
1432c733ed4SBarry Smith       s4 -= v[3] * x1 + v[9] * x2 + v[15] * x3 + v[21] * x4 + v[27] * x5 + v[33] * x6;
1442c733ed4SBarry Smith       s5 -= v[4] * x1 + v[10] * x2 + v[16] * x3 + v[22] * x4 + v[28] * x5 + v[34] * x6;
1452c733ed4SBarry Smith       s6 -= v[5] * x1 + v[11] * x2 + v[17] * x3 + v[23] * x4 + v[29] * x5 + v[35] * x6;
1462c733ed4SBarry Smith       v += bs2;
1472c733ed4SBarry Smith     }
1482c733ed4SBarry Smith 
1492c733ed4SBarry Smith     x[idx]     = s1;
1502c733ed4SBarry Smith     x[1 + idx] = s2;
1512c733ed4SBarry Smith     x[2 + idx] = s3;
1522c733ed4SBarry Smith     x[3 + idx] = s4;
1532c733ed4SBarry Smith     x[4 + idx] = s5;
1542c733ed4SBarry Smith     x[5 + idx] = s6;
1552c733ed4SBarry Smith   }
1562c733ed4SBarry Smith 
1572c733ed4SBarry Smith   /* backward solve the upper triangular */
1582c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1592c733ed4SBarry Smith     v   = aa + bs2 * (adiag[i + 1] + 1);
1602c733ed4SBarry Smith     vi  = aj + adiag[i + 1] + 1;
1612c733ed4SBarry Smith     nz  = adiag[i] - adiag[i + 1] - 1;
1622c733ed4SBarry Smith     idt = bs * i;
163*9371c9d4SSatish Balay     s1  = x[idt];
164*9371c9d4SSatish Balay     s2  = x[1 + idt];
165*9371c9d4SSatish Balay     s3  = x[2 + idt];
166*9371c9d4SSatish Balay     s4  = x[3 + idt];
167*9371c9d4SSatish Balay     s5  = x[4 + idt];
168*9371c9d4SSatish Balay     s6  = x[5 + idt];
1692c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
1702c733ed4SBarry Smith       idx = bs * vi[k];
171*9371c9d4SSatish Balay       x1  = x[idx];
172*9371c9d4SSatish Balay       x2  = x[1 + idx];
173*9371c9d4SSatish Balay       x3  = x[2 + idx];
174*9371c9d4SSatish Balay       x4  = x[3 + idx];
175*9371c9d4SSatish Balay       x5  = x[4 + idx];
176*9371c9d4SSatish Balay       x6  = x[5 + idx];
1772c733ed4SBarry Smith       s1 -= v[0] * x1 + v[6] * x2 + v[12] * x3 + v[18] * x4 + v[24] * x5 + v[30] * x6;
178feb237baSPierre Jolivet       s2 -= v[1] * x1 + v[7] * x2 + v[13] * x3 + v[19] * x4 + v[25] * x5 + v[31] * x6;
1792c733ed4SBarry Smith       s3 -= v[2] * x1 + v[8] * x2 + v[14] * x3 + v[20] * x4 + v[26] * x5 + v[32] * x6;
1802c733ed4SBarry Smith       s4 -= v[3] * x1 + v[9] * x2 + v[15] * x3 + v[21] * x4 + v[27] * x5 + v[33] * x6;
1812c733ed4SBarry Smith       s5 -= v[4] * x1 + v[10] * x2 + v[16] * x3 + v[22] * x4 + v[28] * x5 + v[34] * x6;
1822c733ed4SBarry Smith       s6 -= v[5] * x1 + v[11] * x2 + v[17] * x3 + v[23] * x4 + v[29] * x5 + v[35] * x6;
1832c733ed4SBarry Smith       v += bs2;
1842c733ed4SBarry Smith     }
1852c733ed4SBarry Smith     /* x = inv_diagonal*x */
1862c733ed4SBarry Smith     x[idt]     = v[0] * s1 + v[6] * s2 + v[12] * s3 + v[18] * s4 + v[24] * s5 + v[30] * s6;
1872c733ed4SBarry Smith     x[1 + idt] = v[1] * s1 + v[7] * s2 + v[13] * s3 + v[19] * s4 + v[25] * s5 + v[31] * s6;
1882c733ed4SBarry Smith     x[2 + idt] = v[2] * s1 + v[8] * s2 + v[14] * s3 + v[20] * s4 + v[26] * s5 + v[32] * s6;
1892c733ed4SBarry Smith     x[3 + idt] = v[3] * s1 + v[9] * s2 + v[15] * s3 + v[21] * s4 + v[27] * s5 + v[33] * s6;
1902c733ed4SBarry Smith     x[4 + idt] = v[4] * s1 + v[10] * s2 + v[16] * s3 + v[22] * s4 + v[28] * s5 + v[34] * s6;
1912c733ed4SBarry Smith     x[5 + idt] = v[5] * s1 + v[11] * s2 + v[17] * s3 + v[23] * s4 + v[29] * s5 + v[35] * s6;
1922c733ed4SBarry Smith   }
1932c733ed4SBarry Smith 
1949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1969566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
1972c733ed4SBarry Smith   PetscFunctionReturn(0);
1982c733ed4SBarry Smith }
199