xref: /petsc/src/mat/impls/baij/seq/baijsolvtran5.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 MatSolveTranspose_SeqBAIJ_5_inplace(Mat A, Vec bb, Vec xx) {
52c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
62c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
72c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout;
82c733ed4SBarry Smith   const PetscInt    *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j;
92c733ed4SBarry Smith   PetscInt           i, nz, idx, idt, ii, ic, ir, oidx;
102c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
112c733ed4SBarry Smith   PetscScalar        s1, s2, s3, s4, s5, x1, x2, x3, x4, x5, *x, *t;
122c733ed4SBarry Smith   const PetscScalar *b;
132c733ed4SBarry Smith 
142c733ed4SBarry Smith   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
169566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
172c733ed4SBarry Smith   t = a->solve_work;
182c733ed4SBarry Smith 
19*9371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
20*9371c9d4SSatish Balay   r = rout;
21*9371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
22*9371c9d4SSatish Balay   c = cout;
232c733ed4SBarry Smith 
242c733ed4SBarry Smith   /* copy the b into temp work space according to permutation */
252c733ed4SBarry Smith   ii = 0;
262c733ed4SBarry Smith   for (i = 0; i < n; i++) {
272c733ed4SBarry Smith     ic        = 5 * c[i];
282c733ed4SBarry Smith     t[ii]     = b[ic];
292c733ed4SBarry Smith     t[ii + 1] = b[ic + 1];
302c733ed4SBarry Smith     t[ii + 2] = b[ic + 2];
312c733ed4SBarry Smith     t[ii + 3] = b[ic + 3];
322c733ed4SBarry Smith     t[ii + 4] = b[ic + 4];
332c733ed4SBarry Smith     ii += 5;
342c733ed4SBarry Smith   }
352c733ed4SBarry Smith 
362c733ed4SBarry Smith   /* forward solve the U^T */
372c733ed4SBarry Smith   idx = 0;
382c733ed4SBarry Smith   for (i = 0; i < n; i++) {
392c733ed4SBarry Smith     v  = aa + 25 * diag[i];
402c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
41*9371c9d4SSatish Balay     x1 = t[idx];
42*9371c9d4SSatish Balay     x2 = t[1 + idx];
43*9371c9d4SSatish Balay     x3 = t[2 + idx];
44*9371c9d4SSatish Balay     x4 = t[3 + idx];
45*9371c9d4SSatish Balay     x5 = t[4 + idx];
462c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5;
472c733ed4SBarry Smith     s2 = v[5] * x1 + v[6] * x2 + v[7] * x3 + v[8] * x4 + v[9] * x5;
482c733ed4SBarry Smith     s3 = v[10] * x1 + v[11] * x2 + v[12] * x3 + v[13] * x4 + v[14] * x5;
492c733ed4SBarry Smith     s4 = v[15] * x1 + v[16] * x2 + v[17] * x3 + v[18] * x4 + v[19] * x5;
502c733ed4SBarry Smith     s5 = v[20] * x1 + v[21] * x2 + v[22] * x3 + v[23] * x4 + v[24] * x5;
512c733ed4SBarry Smith     v += 25;
522c733ed4SBarry Smith 
532c733ed4SBarry Smith     vi = aj + diag[i] + 1;
542c733ed4SBarry Smith     nz = ai[i + 1] - diag[i] - 1;
552c733ed4SBarry Smith     while (nz--) {
562c733ed4SBarry Smith       oidx = 5 * (*vi++);
572c733ed4SBarry Smith       t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
582c733ed4SBarry Smith       t[oidx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
592c733ed4SBarry Smith       t[oidx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
602c733ed4SBarry Smith       t[oidx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
612c733ed4SBarry Smith       t[oidx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
622c733ed4SBarry Smith       v += 25;
632c733ed4SBarry Smith     }
64*9371c9d4SSatish Balay     t[idx]     = s1;
65*9371c9d4SSatish Balay     t[1 + idx] = s2;
66*9371c9d4SSatish Balay     t[2 + idx] = s3;
67*9371c9d4SSatish Balay     t[3 + idx] = s4;
68*9371c9d4SSatish Balay     t[4 + idx] = s5;
692c733ed4SBarry Smith     idx += 5;
702c733ed4SBarry Smith   }
712c733ed4SBarry Smith   /* backward solve the L^T */
722c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
732c733ed4SBarry Smith     v   = aa + 25 * diag[i] - 25;
742c733ed4SBarry Smith     vi  = aj + diag[i] - 1;
752c733ed4SBarry Smith     nz  = diag[i] - ai[i];
762c733ed4SBarry Smith     idt = 5 * i;
77*9371c9d4SSatish Balay     s1  = t[idt];
78*9371c9d4SSatish Balay     s2  = t[1 + idt];
79*9371c9d4SSatish Balay     s3  = t[2 + idt];
80*9371c9d4SSatish Balay     s4  = t[3 + idt];
81*9371c9d4SSatish Balay     s5  = t[4 + idt];
822c733ed4SBarry Smith     while (nz--) {
832c733ed4SBarry Smith       idx = 5 * (*vi--);
842c733ed4SBarry Smith       t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
852c733ed4SBarry Smith       t[idx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
862c733ed4SBarry Smith       t[idx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
872c733ed4SBarry Smith       t[idx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
882c733ed4SBarry Smith       t[idx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
892c733ed4SBarry Smith       v -= 25;
902c733ed4SBarry Smith     }
912c733ed4SBarry Smith   }
922c733ed4SBarry Smith 
932c733ed4SBarry Smith   /* copy t into x according to permutation */
942c733ed4SBarry Smith   ii = 0;
952c733ed4SBarry Smith   for (i = 0; i < n; i++) {
962c733ed4SBarry Smith     ir        = 5 * r[i];
972c733ed4SBarry Smith     x[ir]     = t[ii];
982c733ed4SBarry Smith     x[ir + 1] = t[ii + 1];
992c733ed4SBarry Smith     x[ir + 2] = t[ii + 2];
1002c733ed4SBarry Smith     x[ir + 3] = t[ii + 3];
1012c733ed4SBarry Smith     x[ir + 4] = t[ii + 4];
1022c733ed4SBarry Smith     ii += 5;
1032c733ed4SBarry Smith   }
1042c733ed4SBarry Smith 
1059566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
1069566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
1079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1099566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 25 * (a->nz) - 5.0 * A->cmap->n));
1102c733ed4SBarry Smith   PetscFunctionReturn(0);
1112c733ed4SBarry Smith }
1122c733ed4SBarry Smith 
113*9371c9d4SSatish Balay PetscErrorCode MatSolveTranspose_SeqBAIJ_5(Mat A, Vec bb, Vec xx) {
1142c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
1152c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
1162c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
1172c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout;
1182c733ed4SBarry Smith   PetscInt           nz, idx, idt, j, i, oidx, ii, ic, ir;
1192c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
1202c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
1212c733ed4SBarry Smith   PetscScalar        s1, s2, s3, s4, s5, x1, x2, x3, x4, x5, *x, *t;
1222c733ed4SBarry Smith   const PetscScalar *b;
1232c733ed4SBarry Smith 
1242c733ed4SBarry Smith   PetscFunctionBegin;
1259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1272c733ed4SBarry Smith   t = a->solve_work;
1282c733ed4SBarry Smith 
129*9371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
130*9371c9d4SSatish Balay   r = rout;
131*9371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
132*9371c9d4SSatish Balay   c = cout;
1332c733ed4SBarry Smith 
1342c733ed4SBarry Smith   /* copy b into temp work space according to permutation */
1352c733ed4SBarry Smith   for (i = 0; i < n; i++) {
136*9371c9d4SSatish Balay     ii        = bs * i;
137*9371c9d4SSatish Balay     ic        = bs * c[i];
138*9371c9d4SSatish Balay     t[ii]     = b[ic];
139*9371c9d4SSatish Balay     t[ii + 1] = b[ic + 1];
140*9371c9d4SSatish Balay     t[ii + 2] = b[ic + 2];
141*9371c9d4SSatish Balay     t[ii + 3] = b[ic + 3];
1422c733ed4SBarry Smith     t[ii + 4] = b[ic + 4];
1432c733ed4SBarry Smith   }
1442c733ed4SBarry Smith 
1452c733ed4SBarry Smith   /* forward solve the U^T */
1462c733ed4SBarry Smith   idx = 0;
1472c733ed4SBarry Smith   for (i = 0; i < n; i++) {
1482c733ed4SBarry Smith     v  = aa + bs2 * diag[i];
1492c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
150*9371c9d4SSatish Balay     x1 = t[idx];
151*9371c9d4SSatish Balay     x2 = t[1 + idx];
152*9371c9d4SSatish Balay     x3 = t[2 + idx];
153*9371c9d4SSatish Balay     x4 = t[3 + idx];
154*9371c9d4SSatish Balay     x5 = t[4 + idx];
1552c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5;
1562c733ed4SBarry Smith     s2 = v[5] * x1 + v[6] * x2 + v[7] * x3 + v[8] * x4 + v[9] * x5;
1572c733ed4SBarry Smith     s3 = v[10] * x1 + v[11] * x2 + v[12] * x3 + v[13] * x4 + v[14] * x5;
1582c733ed4SBarry Smith     s4 = v[15] * x1 + v[16] * x2 + v[17] * x3 + v[18] * x4 + v[19] * x5;
1592c733ed4SBarry Smith     s5 = v[20] * x1 + v[21] * x2 + v[22] * x3 + v[23] * x4 + v[24] * x5;
1602c733ed4SBarry Smith     v -= bs2;
1612c733ed4SBarry Smith 
1622c733ed4SBarry Smith     vi = aj + diag[i] - 1;
1632c733ed4SBarry Smith     nz = diag[i] - diag[i + 1] - 1;
1642c733ed4SBarry Smith     for (j = 0; j > -nz; j--) {
1652c733ed4SBarry Smith       oidx = bs * vi[j];
1662c733ed4SBarry Smith       t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
1672c733ed4SBarry Smith       t[oidx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
1682c733ed4SBarry Smith       t[oidx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
1692c733ed4SBarry Smith       t[oidx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
1702c733ed4SBarry Smith       t[oidx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
1712c733ed4SBarry Smith       v -= bs2;
1722c733ed4SBarry Smith     }
173*9371c9d4SSatish Balay     t[idx]     = s1;
174*9371c9d4SSatish Balay     t[1 + idx] = s2;
175*9371c9d4SSatish Balay     t[2 + idx] = s3;
176*9371c9d4SSatish Balay     t[3 + idx] = s4;
177*9371c9d4SSatish Balay     t[4 + idx] = s5;
1782c733ed4SBarry Smith     idx += bs;
1792c733ed4SBarry Smith   }
1802c733ed4SBarry Smith   /* backward solve the L^T */
1812c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1822c733ed4SBarry Smith     v   = aa + bs2 * ai[i];
1832c733ed4SBarry Smith     vi  = aj + ai[i];
1842c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1852c733ed4SBarry Smith     idt = bs * i;
186*9371c9d4SSatish Balay     s1  = t[idt];
187*9371c9d4SSatish Balay     s2  = t[1 + idt];
188*9371c9d4SSatish Balay     s3  = t[2 + idt];
189*9371c9d4SSatish Balay     s4  = t[3 + idt];
190*9371c9d4SSatish Balay     s5  = t[4 + idt];
1912c733ed4SBarry Smith     for (j = 0; j < nz; j++) {
1922c733ed4SBarry Smith       idx = bs * vi[j];
1932c733ed4SBarry Smith       t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
1942c733ed4SBarry Smith       t[idx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
1952c733ed4SBarry Smith       t[idx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
1962c733ed4SBarry Smith       t[idx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
1972c733ed4SBarry Smith       t[idx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
1982c733ed4SBarry Smith       v += bs2;
1992c733ed4SBarry Smith     }
2002c733ed4SBarry Smith   }
2012c733ed4SBarry Smith 
2022c733ed4SBarry Smith   /* copy t into x according to permutation */
2032c733ed4SBarry Smith   for (i = 0; i < n; i++) {
204*9371c9d4SSatish Balay     ii        = bs * i;
205*9371c9d4SSatish Balay     ir        = bs * r[i];
206*9371c9d4SSatish Balay     x[ir]     = t[ii];
207*9371c9d4SSatish Balay     x[ir + 1] = t[ii + 1];
208*9371c9d4SSatish Balay     x[ir + 2] = t[ii + 2];
209*9371c9d4SSatish Balay     x[ir + 3] = t[ii + 3];
2102c733ed4SBarry Smith     x[ir + 4] = t[ii + 4];
2112c733ed4SBarry Smith   }
2122c733ed4SBarry Smith 
2139566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
2149566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
2159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
2169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
2179566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
2182c733ed4SBarry Smith   PetscFunctionReturn(0);
2192c733ed4SBarry Smith }
220