xref: /petsc/src/mat/impls/baij/seq/baijsolvtran6.c (revision 31d78bcd2b98084dc1368b20eb1129c8b9fb39fe)
12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h>
32c733ed4SBarry Smith 
MatSolveTranspose_SeqBAIJ_6_inplace(Mat A,Vec bb,Vec xx)4d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_6_inplace(Mat A, Vec bb, Vec xx)
5d71ae5a4SJacob Faibussowitsch {
62c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
72c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
82c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout;
92c733ed4SBarry Smith   const PetscInt    *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j;
102c733ed4SBarry Smith   PetscInt           i, nz, idx, idt, ii, ic, ir, oidx;
112c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
122c733ed4SBarry Smith   PetscScalar        s1, s2, s3, s4, s5, s6, x1, x2, x3, x4, x5, x6, *x, *t;
132c733ed4SBarry Smith   const PetscScalar *b;
142c733ed4SBarry Smith 
152c733ed4SBarry Smith   PetscFunctionBegin;
169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
182c733ed4SBarry Smith   t = a->solve_work;
192c733ed4SBarry Smith 
209371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
219371c9d4SSatish Balay   r = rout;
229371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
239371c9d4SSatish Balay   c = cout;
242c733ed4SBarry Smith 
252c733ed4SBarry Smith   /* copy the b into temp work space according to permutation */
262c733ed4SBarry Smith   ii = 0;
272c733ed4SBarry Smith   for (i = 0; i < n; i++) {
282c733ed4SBarry Smith     ic        = 6 * c[i];
292c733ed4SBarry Smith     t[ii]     = b[ic];
302c733ed4SBarry Smith     t[ii + 1] = b[ic + 1];
312c733ed4SBarry Smith     t[ii + 2] = b[ic + 2];
322c733ed4SBarry Smith     t[ii + 3] = b[ic + 3];
332c733ed4SBarry Smith     t[ii + 4] = b[ic + 4];
342c733ed4SBarry Smith     t[ii + 5] = b[ic + 5];
352c733ed4SBarry Smith     ii += 6;
362c733ed4SBarry Smith   }
372c733ed4SBarry Smith 
382c733ed4SBarry Smith   /* forward solve the U^T */
392c733ed4SBarry Smith   idx = 0;
402c733ed4SBarry Smith   for (i = 0; i < n; i++) {
412c733ed4SBarry Smith     v = aa + 36 * diag[i];
422c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
439371c9d4SSatish Balay     x1 = t[idx];
449371c9d4SSatish Balay     x2 = t[1 + idx];
459371c9d4SSatish Balay     x3 = t[2 + idx];
469371c9d4SSatish Balay     x4 = t[3 + idx];
479371c9d4SSatish Balay     x5 = t[4 + idx];
482c733ed4SBarry Smith     x6 = t[5 + idx];
492c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6;
502c733ed4SBarry Smith     s2 = v[6] * x1 + v[7] * x2 + v[8] * x3 + v[9] * x4 + v[10] * x5 + v[11] * x6;
512c733ed4SBarry Smith     s3 = v[12] * x1 + v[13] * x2 + v[14] * x3 + v[15] * x4 + v[16] * x5 + v[17] * x6;
522c733ed4SBarry Smith     s4 = v[18] * x1 + v[19] * x2 + v[20] * x3 + v[21] * x4 + v[22] * x5 + v[23] * x6;
532c733ed4SBarry Smith     s5 = v[24] * x1 + v[25] * x2 + v[26] * x3 + v[27] * x4 + v[28] * x5 + v[29] * x6;
542c733ed4SBarry Smith     s6 = v[30] * x1 + v[31] * x2 + v[32] * x3 + v[33] * x4 + v[34] * x5 + v[35] * x6;
552c733ed4SBarry Smith     v += 36;
562c733ed4SBarry Smith 
572c733ed4SBarry Smith     vi = aj + diag[i] + 1;
582c733ed4SBarry Smith     nz = ai[i + 1] - diag[i] - 1;
592c733ed4SBarry Smith     while (nz--) {
602c733ed4SBarry Smith       oidx = 6 * (*vi++);
612c733ed4SBarry Smith       t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
622c733ed4SBarry Smith       t[oidx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
632c733ed4SBarry Smith       t[oidx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
642c733ed4SBarry Smith       t[oidx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
652c733ed4SBarry Smith       t[oidx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
662c733ed4SBarry Smith       t[oidx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
672c733ed4SBarry Smith       v += 36;
682c733ed4SBarry Smith     }
699371c9d4SSatish Balay     t[idx]     = s1;
709371c9d4SSatish Balay     t[1 + idx] = s2;
719371c9d4SSatish Balay     t[2 + idx] = s3;
729371c9d4SSatish Balay     t[3 + idx] = s4;
739371c9d4SSatish Balay     t[4 + idx] = s5;
742c733ed4SBarry Smith     t[5 + idx] = s6;
752c733ed4SBarry Smith     idx += 6;
762c733ed4SBarry Smith   }
772c733ed4SBarry Smith   /* backward solve the L^T */
782c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
792c733ed4SBarry Smith     v   = aa + 36 * diag[i] - 36;
802c733ed4SBarry Smith     vi  = aj + diag[i] - 1;
812c733ed4SBarry Smith     nz  = diag[i] - ai[i];
822c733ed4SBarry Smith     idt = 6 * i;
839371c9d4SSatish Balay     s1  = t[idt];
849371c9d4SSatish Balay     s2  = t[1 + idt];
859371c9d4SSatish Balay     s3  = t[2 + idt];
869371c9d4SSatish Balay     s4  = t[3 + idt];
879371c9d4SSatish Balay     s5  = t[4 + idt];
882c733ed4SBarry Smith     s6  = t[5 + idt];
892c733ed4SBarry Smith     while (nz--) {
902c733ed4SBarry Smith       idx = 6 * (*vi--);
912c733ed4SBarry Smith       t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
922c733ed4SBarry Smith       t[idx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
932c733ed4SBarry Smith       t[idx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
942c733ed4SBarry Smith       t[idx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
952c733ed4SBarry Smith       t[idx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
962c733ed4SBarry Smith       t[idx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
972c733ed4SBarry Smith       v -= 36;
982c733ed4SBarry Smith     }
992c733ed4SBarry Smith   }
1002c733ed4SBarry Smith 
1012c733ed4SBarry Smith   /* copy t into x according to permutation */
1022c733ed4SBarry Smith   ii = 0;
1032c733ed4SBarry Smith   for (i = 0; i < n; i++) {
1042c733ed4SBarry Smith     ir        = 6 * r[i];
1052c733ed4SBarry Smith     x[ir]     = t[ii];
1062c733ed4SBarry Smith     x[ir + 1] = t[ii + 1];
1072c733ed4SBarry Smith     x[ir + 2] = t[ii + 2];
1082c733ed4SBarry Smith     x[ir + 3] = t[ii + 3];
1092c733ed4SBarry Smith     x[ir + 4] = t[ii + 4];
1102c733ed4SBarry Smith     x[ir + 5] = t[ii + 5];
1112c733ed4SBarry Smith     ii += 6;
1122c733ed4SBarry Smith   }
1132c733ed4SBarry Smith 
1149566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
1159566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
1169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1189566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 36 * (a->nz) - 6.0 * A->cmap->n));
119*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1202c733ed4SBarry Smith }
1212c733ed4SBarry Smith 
MatSolveTranspose_SeqBAIJ_6(Mat A,Vec bb,Vec xx)122d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_6(Mat A, Vec bb, Vec xx)
123d71ae5a4SJacob Faibussowitsch {
1242c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
1252c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
1262c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
1272c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout;
1282c733ed4SBarry Smith   PetscInt           nz, idx, idt, j, i, oidx, ii, ic, ir;
1292c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
1302c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
1312c733ed4SBarry Smith   PetscScalar        s1, s2, s3, s4, s5, s6, x1, x2, x3, x4, x5, x6, *x, *t;
1322c733ed4SBarry Smith   const PetscScalar *b;
1332c733ed4SBarry Smith 
1342c733ed4SBarry Smith   PetscFunctionBegin;
1359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1372c733ed4SBarry Smith   t = a->solve_work;
1382c733ed4SBarry Smith 
1399371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
1409371c9d4SSatish Balay   r = rout;
1419371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
1429371c9d4SSatish Balay   c = cout;
1432c733ed4SBarry Smith 
1442c733ed4SBarry Smith   /* copy b into temp work space according to permutation */
1452c733ed4SBarry Smith   for (i = 0; i < n; i++) {
1469371c9d4SSatish Balay     ii        = bs * i;
1479371c9d4SSatish Balay     ic        = bs * c[i];
1489371c9d4SSatish Balay     t[ii]     = b[ic];
1499371c9d4SSatish Balay     t[ii + 1] = b[ic + 1];
1509371c9d4SSatish Balay     t[ii + 2] = b[ic + 2];
1519371c9d4SSatish Balay     t[ii + 3] = b[ic + 3];
1529371c9d4SSatish Balay     t[ii + 4] = b[ic + 4];
1539371c9d4SSatish Balay     t[ii + 5] = b[ic + 5];
1542c733ed4SBarry Smith   }
1552c733ed4SBarry Smith 
1562c733ed4SBarry Smith   /* forward solve the U^T */
1572c733ed4SBarry Smith   idx = 0;
1582c733ed4SBarry Smith   for (i = 0; i < n; i++) {
1592c733ed4SBarry Smith     v = aa + bs2 * diag[i];
1602c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
1619371c9d4SSatish Balay     x1 = t[idx];
1629371c9d4SSatish Balay     x2 = t[1 + idx];
1639371c9d4SSatish Balay     x3 = t[2 + idx];
1649371c9d4SSatish Balay     x4 = t[3 + idx];
1659371c9d4SSatish Balay     x5 = t[4 + idx];
1662c733ed4SBarry Smith     x6 = t[5 + idx];
1672c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6;
1682c733ed4SBarry Smith     s2 = v[6] * x1 + v[7] * x2 + v[8] * x3 + v[9] * x4 + v[10] * x5 + v[11] * x6;
1692c733ed4SBarry Smith     s3 = v[12] * x1 + v[13] * x2 + v[14] * x3 + v[15] * x4 + v[16] * x5 + v[17] * x6;
1702c733ed4SBarry Smith     s4 = v[18] * x1 + v[19] * x2 + v[20] * x3 + v[21] * x4 + v[22] * x5 + v[23] * x6;
1712c733ed4SBarry Smith     s5 = v[24] * x1 + v[25] * x2 + v[26] * x3 + v[27] * x4 + v[28] * x5 + v[29] * x6;
1722c733ed4SBarry Smith     s6 = v[30] * x1 + v[31] * x2 + v[32] * x3 + v[33] * x4 + v[34] * x5 + v[35] * x6;
1732c733ed4SBarry Smith     v -= bs2;
1742c733ed4SBarry Smith 
1752c733ed4SBarry Smith     vi = aj + diag[i] - 1;
1762c733ed4SBarry Smith     nz = diag[i] - diag[i + 1] - 1;
1772c733ed4SBarry Smith     for (j = 0; j > -nz; j--) {
1782c733ed4SBarry Smith       oidx = bs * vi[j];
1792c733ed4SBarry Smith       t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
1802c733ed4SBarry Smith       t[oidx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
1812c733ed4SBarry Smith       t[oidx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
1822c733ed4SBarry Smith       t[oidx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
1832c733ed4SBarry Smith       t[oidx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
1842c733ed4SBarry Smith       t[oidx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
1852c733ed4SBarry Smith       v -= bs2;
1862c733ed4SBarry Smith     }
1879371c9d4SSatish Balay     t[idx]     = s1;
1889371c9d4SSatish Balay     t[1 + idx] = s2;
1899371c9d4SSatish Balay     t[2 + idx] = s3;
1909371c9d4SSatish Balay     t[3 + idx] = s4;
1919371c9d4SSatish Balay     t[4 + idx] = s5;
1922c733ed4SBarry Smith     t[5 + idx] = s6;
1932c733ed4SBarry Smith     idx += bs;
1942c733ed4SBarry Smith   }
1952c733ed4SBarry Smith   /* backward solve the L^T */
1962c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1972c733ed4SBarry Smith     v   = aa + bs2 * ai[i];
1982c733ed4SBarry Smith     vi  = aj + ai[i];
1992c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
2002c733ed4SBarry Smith     idt = bs * i;
2019371c9d4SSatish Balay     s1  = t[idt];
2029371c9d4SSatish Balay     s2  = t[1 + idt];
2039371c9d4SSatish Balay     s3  = t[2 + idt];
2049371c9d4SSatish Balay     s4  = t[3 + idt];
2059371c9d4SSatish Balay     s5  = t[4 + idt];
2062c733ed4SBarry Smith     s6  = t[5 + idt];
2072c733ed4SBarry Smith     for (j = 0; j < nz; j++) {
2082c733ed4SBarry Smith       idx = bs * vi[j];
2092c733ed4SBarry Smith       t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
2102c733ed4SBarry Smith       t[idx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
2112c733ed4SBarry Smith       t[idx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
2122c733ed4SBarry Smith       t[idx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
2132c733ed4SBarry Smith       t[idx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
2142c733ed4SBarry Smith       t[idx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
2152c733ed4SBarry Smith       v += bs2;
2162c733ed4SBarry Smith     }
2172c733ed4SBarry Smith   }
2182c733ed4SBarry Smith 
2192c733ed4SBarry Smith   /* copy t into x according to permutation */
2202c733ed4SBarry Smith   for (i = 0; i < n; i++) {
2219371c9d4SSatish Balay     ii        = bs * i;
2229371c9d4SSatish Balay     ir        = bs * r[i];
2239371c9d4SSatish Balay     x[ir]     = t[ii];
2249371c9d4SSatish Balay     x[ir + 1] = t[ii + 1];
2259371c9d4SSatish Balay     x[ir + 2] = t[ii + 2];
2269371c9d4SSatish Balay     x[ir + 3] = t[ii + 3];
2279371c9d4SSatish Balay     x[ir + 4] = t[ii + 4];
2289371c9d4SSatish Balay     x[ir + 5] = t[ii + 5];
2292c733ed4SBarry Smith   }
2302c733ed4SBarry Smith 
2319566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
2329566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
2339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
2349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
2359566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
236*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2372c733ed4SBarry Smith }
238