xref: /petsc/src/mat/impls/baij/seq/baijsolvnat3.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h>
32c733ed4SBarry Smith 
42c733ed4SBarry Smith /*
52c733ed4SBarry Smith       Special case where the matrix was ILU(0) factored in the natural
62c733ed4SBarry Smith    ordering. This eliminates the need for the column and row permutation.
72c733ed4SBarry Smith */
8*9371c9d4SSatish Balay PetscErrorCode MatSolve_SeqBAIJ_3_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx) {
92c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
102c733ed4SBarry Smith   const PetscInt     n = a->mbs, *ai = a->i, *aj = a->j;
112c733ed4SBarry Smith   const PetscInt    *diag = a->diag, *vi;
122c733ed4SBarry Smith   const MatScalar   *aa   = a->a, *v;
132c733ed4SBarry Smith   PetscScalar       *x, s1, s2, s3, x1, x2, x3;
142c733ed4SBarry Smith   const PetscScalar *b;
152c733ed4SBarry Smith   PetscInt           jdx, idt, idx, nz, i;
162c733ed4SBarry Smith 
172c733ed4SBarry Smith   PetscFunctionBegin;
189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
199566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
202c733ed4SBarry Smith 
212c733ed4SBarry Smith   /* forward solve the lower triangular */
222c733ed4SBarry Smith   idx  = 0;
23*9371c9d4SSatish Balay   x[0] = b[0];
24*9371c9d4SSatish Balay   x[1] = b[1];
25*9371c9d4SSatish Balay   x[2] = b[2];
262c733ed4SBarry Smith   for (i = 1; i < n; i++) {
272c733ed4SBarry Smith     v  = aa + 9 * ai[i];
282c733ed4SBarry Smith     vi = aj + ai[i];
292c733ed4SBarry Smith     nz = diag[i] - ai[i];
302c733ed4SBarry Smith     idx += 3;
31*9371c9d4SSatish Balay     s1 = b[idx];
32*9371c9d4SSatish Balay     s2 = b[1 + idx];
33*9371c9d4SSatish Balay     s3 = b[2 + idx];
342c733ed4SBarry Smith     while (nz--) {
352c733ed4SBarry Smith       jdx = 3 * (*vi++);
36*9371c9d4SSatish Balay       x1  = x[jdx];
37*9371c9d4SSatish Balay       x2  = x[1 + jdx];
38*9371c9d4SSatish Balay       x3  = x[2 + jdx];
392c733ed4SBarry Smith       s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
402c733ed4SBarry Smith       s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
412c733ed4SBarry Smith       s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
422c733ed4SBarry Smith       v += 9;
432c733ed4SBarry Smith     }
442c733ed4SBarry Smith     x[idx]     = s1;
452c733ed4SBarry Smith     x[1 + idx] = s2;
462c733ed4SBarry Smith     x[2 + idx] = s3;
472c733ed4SBarry Smith   }
482c733ed4SBarry Smith   /* backward solve the upper triangular */
492c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
502c733ed4SBarry Smith     v   = aa + 9 * diag[i] + 9;
512c733ed4SBarry Smith     vi  = aj + diag[i] + 1;
522c733ed4SBarry Smith     nz  = ai[i + 1] - diag[i] - 1;
532c733ed4SBarry Smith     idt = 3 * i;
54*9371c9d4SSatish Balay     s1  = x[idt];
55*9371c9d4SSatish Balay     s2  = x[1 + idt];
562c733ed4SBarry Smith     s3  = x[2 + idt];
572c733ed4SBarry Smith     while (nz--) {
582c733ed4SBarry Smith       idx = 3 * (*vi++);
59*9371c9d4SSatish Balay       x1  = x[idx];
60*9371c9d4SSatish Balay       x2  = x[1 + idx];
61*9371c9d4SSatish Balay       x3  = x[2 + idx];
622c733ed4SBarry Smith       s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
632c733ed4SBarry Smith       s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
642c733ed4SBarry Smith       s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
652c733ed4SBarry Smith       v += 9;
662c733ed4SBarry Smith     }
672c733ed4SBarry Smith     v          = aa + 9 * diag[i];
682c733ed4SBarry Smith     x[idt]     = v[0] * s1 + v[3] * s2 + v[6] * s3;
692c733ed4SBarry Smith     x[1 + idt] = v[1] * s1 + v[4] * s2 + v[7] * s3;
702c733ed4SBarry Smith     x[2 + idt] = v[2] * s1 + v[5] * s2 + v[8] * s3;
712c733ed4SBarry Smith   }
722c733ed4SBarry Smith 
739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
759566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 9 * (a->nz) - 3.0 * A->cmap->n));
762c733ed4SBarry Smith   PetscFunctionReturn(0);
772c733ed4SBarry Smith }
782c733ed4SBarry Smith 
79*9371c9d4SSatish Balay PetscErrorCode MatSolve_SeqBAIJ_3_NaturalOrdering(Mat A, Vec bb, Vec xx) {
802c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
812c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *ai = a->i, *aj = a->j, *adiag = a->diag;
822c733ed4SBarry Smith   PetscInt           i, k, nz, idx, jdx, idt;
832c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
842c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
852c733ed4SBarry Smith   PetscScalar       *x;
862c733ed4SBarry Smith   const PetscScalar *b;
872c733ed4SBarry Smith   PetscScalar        s1, s2, s3, x1, x2, x3;
882c733ed4SBarry Smith 
892c733ed4SBarry Smith   PetscFunctionBegin;
909566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
919566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
922c733ed4SBarry Smith   /* forward solve the lower triangular */
932c733ed4SBarry Smith   idx  = 0;
94*9371c9d4SSatish Balay   x[0] = b[idx];
95*9371c9d4SSatish Balay   x[1] = b[1 + idx];
96*9371c9d4SSatish Balay   x[2] = b[2 + idx];
972c733ed4SBarry Smith   for (i = 1; i < n; i++) {
982c733ed4SBarry Smith     v   = aa + bs2 * ai[i];
992c733ed4SBarry Smith     vi  = aj + ai[i];
1002c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1012c733ed4SBarry Smith     idx = bs * i;
102*9371c9d4SSatish Balay     s1  = b[idx];
103*9371c9d4SSatish Balay     s2  = b[1 + idx];
104*9371c9d4SSatish Balay     s3  = b[2 + idx];
1052c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
1062c733ed4SBarry Smith       jdx = bs * vi[k];
107*9371c9d4SSatish Balay       x1  = x[jdx];
108*9371c9d4SSatish Balay       x2  = x[1 + jdx];
109*9371c9d4SSatish Balay       x3  = x[2 + jdx];
1102c733ed4SBarry Smith       s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
1112c733ed4SBarry Smith       s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
1122c733ed4SBarry Smith       s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
1132c733ed4SBarry Smith 
1142c733ed4SBarry Smith       v += bs2;
1152c733ed4SBarry Smith     }
1162c733ed4SBarry Smith 
1172c733ed4SBarry Smith     x[idx]     = s1;
1182c733ed4SBarry Smith     x[1 + idx] = s2;
1192c733ed4SBarry Smith     x[2 + idx] = s3;
1202c733ed4SBarry Smith   }
1212c733ed4SBarry Smith 
1222c733ed4SBarry Smith   /* backward solve the upper triangular */
1232c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1242c733ed4SBarry Smith     v   = aa + bs2 * (adiag[i + 1] + 1);
1252c733ed4SBarry Smith     vi  = aj + adiag[i + 1] + 1;
1262c733ed4SBarry Smith     nz  = adiag[i] - adiag[i + 1] - 1;
1272c733ed4SBarry Smith     idt = bs * i;
128*9371c9d4SSatish Balay     s1  = x[idt];
129*9371c9d4SSatish Balay     s2  = x[1 + idt];
130*9371c9d4SSatish Balay     s3  = x[2 + idt];
1312c733ed4SBarry Smith 
1322c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
1332c733ed4SBarry Smith       idx = bs * vi[k];
134*9371c9d4SSatish Balay       x1  = x[idx];
135*9371c9d4SSatish Balay       x2  = x[1 + idx];
136*9371c9d4SSatish Balay       x3  = x[2 + idx];
1372c733ed4SBarry Smith       s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
1382c733ed4SBarry Smith       s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
1392c733ed4SBarry Smith       s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
1402c733ed4SBarry Smith 
1412c733ed4SBarry Smith       v += bs2;
1422c733ed4SBarry Smith     }
1432c733ed4SBarry Smith     /* x = inv_diagonal*x */
1442c733ed4SBarry Smith     x[idt]     = v[0] * s1 + v[3] * s2 + v[6] * s3;
1452c733ed4SBarry Smith     x[1 + idt] = v[1] * s1 + v[4] * s2 + v[7] * s3;
1462c733ed4SBarry Smith     x[2 + idt] = v[2] * s1 + v[5] * s2 + v[8] * s3;
1472c733ed4SBarry Smith   }
1482c733ed4SBarry Smith 
1499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1519566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
1522c733ed4SBarry Smith   PetscFunctionReturn(0);
1532c733ed4SBarry Smith }
1542c733ed4SBarry Smith 
155*9371c9d4SSatish Balay PetscErrorCode MatForwardSolve_SeqBAIJ_3_NaturalOrdering(Mat A, Vec bb, Vec xx) {
1562c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
1572c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *ai = a->i, *aj = a->j;
1582c733ed4SBarry Smith   PetscInt           i, k, nz, idx, jdx;
1592c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
1602c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
1612c733ed4SBarry Smith   PetscScalar       *x;
1622c733ed4SBarry Smith   const PetscScalar *b;
1632c733ed4SBarry Smith   PetscScalar        s1, s2, s3, x1, x2, x3;
1642c733ed4SBarry Smith 
1652c733ed4SBarry Smith   PetscFunctionBegin;
1669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1682c733ed4SBarry Smith   /* forward solve the lower triangular */
1692c733ed4SBarry Smith   idx  = 0;
170*9371c9d4SSatish Balay   x[0] = b[idx];
171*9371c9d4SSatish Balay   x[1] = b[1 + idx];
172*9371c9d4SSatish Balay   x[2] = b[2 + idx];
1732c733ed4SBarry Smith   for (i = 1; i < n; i++) {
1742c733ed4SBarry Smith     v   = aa + bs2 * ai[i];
1752c733ed4SBarry Smith     vi  = aj + ai[i];
1762c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1772c733ed4SBarry Smith     idx = bs * i;
178*9371c9d4SSatish Balay     s1  = b[idx];
179*9371c9d4SSatish Balay     s2  = b[1 + idx];
180*9371c9d4SSatish Balay     s3  = b[2 + idx];
1812c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
1822c733ed4SBarry Smith       jdx = bs * vi[k];
183*9371c9d4SSatish Balay       x1  = x[jdx];
184*9371c9d4SSatish Balay       x2  = x[1 + jdx];
185*9371c9d4SSatish Balay       x3  = x[2 + jdx];
1862c733ed4SBarry Smith       s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
1872c733ed4SBarry Smith       s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
1882c733ed4SBarry Smith       s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
1892c733ed4SBarry Smith 
1902c733ed4SBarry Smith       v += bs2;
1912c733ed4SBarry Smith     }
1922c733ed4SBarry Smith 
1932c733ed4SBarry Smith     x[idx]     = s1;
1942c733ed4SBarry Smith     x[1 + idx] = s2;
1952c733ed4SBarry Smith     x[2 + idx] = s3;
1962c733ed4SBarry Smith   }
1972c733ed4SBarry Smith 
1989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
2009566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.0 * bs2 * (a->nz) - bs * A->cmap->n));
2012c733ed4SBarry Smith   PetscFunctionReturn(0);
2022c733ed4SBarry Smith }
2032c733ed4SBarry Smith 
204*9371c9d4SSatish Balay PetscErrorCode MatBackwardSolve_SeqBAIJ_3_NaturalOrdering(Mat A, Vec bb, Vec xx) {
2052c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
2062c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *aj = a->j, *adiag = a->diag;
2072c733ed4SBarry Smith   PetscInt           i, k, nz, idx, idt;
2082c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
2092c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
2102c733ed4SBarry Smith   PetscScalar       *x;
2112c733ed4SBarry Smith   const PetscScalar *b;
2122c733ed4SBarry Smith   PetscScalar        s1, s2, s3, x1, x2, x3;
2132c733ed4SBarry Smith 
2142c733ed4SBarry Smith   PetscFunctionBegin;
2159566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
2169566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
2172c733ed4SBarry Smith 
2182c733ed4SBarry Smith   /* backward solve the upper triangular */
2192c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
2202c733ed4SBarry Smith     v   = aa + bs2 * (adiag[i + 1] + 1);
2212c733ed4SBarry Smith     vi  = aj + adiag[i + 1] + 1;
2222c733ed4SBarry Smith     nz  = adiag[i] - adiag[i + 1] - 1;
2232c733ed4SBarry Smith     idt = bs * i;
224*9371c9d4SSatish Balay     s1  = b[idt];
225*9371c9d4SSatish Balay     s2  = b[1 + idt];
226*9371c9d4SSatish Balay     s3  = b[2 + idt];
2272c733ed4SBarry Smith 
2282c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
2292c733ed4SBarry Smith       idx = bs * vi[k];
230*9371c9d4SSatish Balay       x1  = x[idx];
231*9371c9d4SSatish Balay       x2  = x[1 + idx];
232*9371c9d4SSatish Balay       x3  = x[2 + idx];
2332c733ed4SBarry Smith       s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
2342c733ed4SBarry Smith       s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
2352c733ed4SBarry Smith       s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
2362c733ed4SBarry Smith 
2372c733ed4SBarry Smith       v += bs2;
2382c733ed4SBarry Smith     }
2392c733ed4SBarry Smith     /* x = inv_diagonal*x */
2402c733ed4SBarry Smith     x[idt]     = v[0] * s1 + v[3] * s2 + v[6] * s3;
2412c733ed4SBarry Smith     x[1 + idt] = v[1] * s1 + v[4] * s2 + v[7] * s3;
2422c733ed4SBarry Smith     x[2 + idt] = v[2] * s1 + v[5] * s2 + v[8] * s3;
2432c733ed4SBarry Smith   }
2442c733ed4SBarry Smith 
2459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
2469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
2479566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
2482c733ed4SBarry Smith   PetscFunctionReturn(0);
2492c733ed4SBarry Smith }
250