12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h>
32c733ed4SBarry Smith
MatSolveTranspose_SeqBAIJ_7_inplace(Mat A,Vec bb,Vec xx)4d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_7_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, s7, x1, x2, x3, x4, x5, x6, x7, *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 = 7 * 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 t[ii + 6] = b[ic + 6];
362c733ed4SBarry Smith ii += 7;
372c733ed4SBarry Smith }
382c733ed4SBarry Smith
392c733ed4SBarry Smith /* forward solve the U^T */
402c733ed4SBarry Smith idx = 0;
412c733ed4SBarry Smith for (i = 0; i < n; i++) {
422c733ed4SBarry Smith v = aa + 49 * diag[i];
432c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */
449371c9d4SSatish Balay x1 = t[idx];
459371c9d4SSatish Balay x2 = t[1 + idx];
469371c9d4SSatish Balay x3 = t[2 + idx];
479371c9d4SSatish Balay x4 = t[3 + idx];
489371c9d4SSatish Balay x5 = t[4 + idx];
499371c9d4SSatish Balay x6 = t[5 + idx];
509371c9d4SSatish Balay x7 = t[6 + idx];
512c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6 + v[6] * x7;
522c733ed4SBarry Smith s2 = v[7] * x1 + v[8] * x2 + v[9] * x3 + v[10] * x4 + v[11] * x5 + v[12] * x6 + v[13] * x7;
532c733ed4SBarry Smith s3 = v[14] * x1 + v[15] * x2 + v[16] * x3 + v[17] * x4 + v[18] * x5 + v[19] * x6 + v[20] * x7;
542c733ed4SBarry Smith s4 = v[21] * x1 + v[22] * x2 + v[23] * x3 + v[24] * x4 + v[25] * x5 + v[26] * x6 + v[27] * x7;
552c733ed4SBarry Smith s5 = v[28] * x1 + v[29] * x2 + v[30] * x3 + v[31] * x4 + v[32] * x5 + v[33] * x6 + v[34] * x7;
562c733ed4SBarry Smith s6 = v[35] * x1 + v[36] * x2 + v[37] * x3 + v[38] * x4 + v[39] * x5 + v[40] * x6 + v[41] * x7;
572c733ed4SBarry Smith s7 = v[42] * x1 + v[43] * x2 + v[44] * x3 + v[45] * x4 + v[46] * x5 + v[47] * x6 + v[48] * x7;
582c733ed4SBarry Smith v += 49;
592c733ed4SBarry Smith
602c733ed4SBarry Smith vi = aj + diag[i] + 1;
612c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1;
622c733ed4SBarry Smith while (nz--) {
632c733ed4SBarry Smith oidx = 7 * (*vi++);
642c733ed4SBarry Smith t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7;
652c733ed4SBarry Smith t[oidx + 1] -= v[7] * s1 + v[8] * s2 + v[9] * s3 + v[10] * s4 + v[11] * s5 + v[12] * s6 + v[13] * s7;
662c733ed4SBarry Smith t[oidx + 2] -= v[14] * s1 + v[15] * s2 + v[16] * s3 + v[17] * s4 + v[18] * s5 + v[19] * s6 + v[20] * s7;
672c733ed4SBarry Smith t[oidx + 3] -= v[21] * s1 + v[22] * s2 + v[23] * s3 + v[24] * s4 + v[25] * s5 + v[26] * s6 + v[27] * s7;
682c733ed4SBarry Smith t[oidx + 4] -= v[28] * s1 + v[29] * s2 + v[30] * s3 + v[31] * s4 + v[32] * s5 + v[33] * s6 + v[34] * s7;
692c733ed4SBarry Smith t[oidx + 5] -= v[35] * s1 + v[36] * s2 + v[37] * s3 + v[38] * s4 + v[39] * s5 + v[40] * s6 + v[41] * s7;
702c733ed4SBarry Smith t[oidx + 6] -= v[42] * s1 + v[43] * s2 + v[44] * s3 + v[45] * s4 + v[46] * s5 + v[47] * s6 + v[48] * s7;
712c733ed4SBarry Smith v += 49;
722c733ed4SBarry Smith }
739371c9d4SSatish Balay t[idx] = s1;
749371c9d4SSatish Balay t[1 + idx] = s2;
759371c9d4SSatish Balay t[2 + idx] = s3;
769371c9d4SSatish Balay t[3 + idx] = s4;
779371c9d4SSatish Balay t[4 + idx] = s5;
789371c9d4SSatish Balay t[5 + idx] = s6;
799371c9d4SSatish Balay t[6 + idx] = s7;
802c733ed4SBarry Smith idx += 7;
812c733ed4SBarry Smith }
822c733ed4SBarry Smith /* backward solve the L^T */
832c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) {
842c733ed4SBarry Smith v = aa + 49 * diag[i] - 49;
852c733ed4SBarry Smith vi = aj + diag[i] - 1;
862c733ed4SBarry Smith nz = diag[i] - ai[i];
872c733ed4SBarry Smith idt = 7 * i;
889371c9d4SSatish Balay s1 = t[idt];
899371c9d4SSatish Balay s2 = t[1 + idt];
909371c9d4SSatish Balay s3 = t[2 + idt];
919371c9d4SSatish Balay s4 = t[3 + idt];
929371c9d4SSatish Balay s5 = t[4 + idt];
939371c9d4SSatish Balay s6 = t[5 + idt];
949371c9d4SSatish Balay s7 = t[6 + idt];
952c733ed4SBarry Smith while (nz--) {
962c733ed4SBarry Smith idx = 7 * (*vi--);
972c733ed4SBarry Smith t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7;
982c733ed4SBarry Smith t[idx + 1] -= v[7] * s1 + v[8] * s2 + v[9] * s3 + v[10] * s4 + v[11] * s5 + v[12] * s6 + v[13] * s7;
992c733ed4SBarry Smith t[idx + 2] -= v[14] * s1 + v[15] * s2 + v[16] * s3 + v[17] * s4 + v[18] * s5 + v[19] * s6 + v[20] * s7;
1002c733ed4SBarry Smith t[idx + 3] -= v[21] * s1 + v[22] * s2 + v[23] * s3 + v[24] * s4 + v[25] * s5 + v[26] * s6 + v[27] * s7;
1012c733ed4SBarry Smith t[idx + 4] -= v[28] * s1 + v[29] * s2 + v[30] * s3 + v[31] * s4 + v[32] * s5 + v[33] * s6 + v[34] * s7;
1022c733ed4SBarry Smith t[idx + 5] -= v[35] * s1 + v[36] * s2 + v[37] * s3 + v[38] * s4 + v[39] * s5 + v[40] * s6 + v[41] * s7;
1032c733ed4SBarry Smith t[idx + 6] -= v[42] * s1 + v[43] * s2 + v[44] * s3 + v[45] * s4 + v[46] * s5 + v[47] * s6 + v[48] * s7;
1042c733ed4SBarry Smith v -= 49;
1052c733ed4SBarry Smith }
1062c733ed4SBarry Smith }
1072c733ed4SBarry Smith
1082c733ed4SBarry Smith /* copy t into x according to permutation */
1092c733ed4SBarry Smith ii = 0;
1102c733ed4SBarry Smith for (i = 0; i < n; i++) {
1112c733ed4SBarry Smith ir = 7 * r[i];
1122c733ed4SBarry Smith x[ir] = t[ii];
1132c733ed4SBarry Smith x[ir + 1] = t[ii + 1];
1142c733ed4SBarry Smith x[ir + 2] = t[ii + 2];
1152c733ed4SBarry Smith x[ir + 3] = t[ii + 3];
1162c733ed4SBarry Smith x[ir + 4] = t[ii + 4];
1172c733ed4SBarry Smith x[ir + 5] = t[ii + 5];
1182c733ed4SBarry Smith x[ir + 6] = t[ii + 6];
1192c733ed4SBarry Smith ii += 7;
1202c733ed4SBarry Smith }
1212c733ed4SBarry Smith
1229566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout));
1239566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout));
1249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b));
1259566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
1269566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 49 * (a->nz) - 7.0 * A->cmap->n));
127*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1282c733ed4SBarry Smith }
MatSolveTranspose_SeqBAIJ_7(Mat A,Vec bb,Vec xx)129d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_7(Mat A, Vec bb, Vec xx)
130d71ae5a4SJacob Faibussowitsch {
1312c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
1322c733ed4SBarry Smith IS iscol = a->col, isrow = a->row;
1332c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
1342c733ed4SBarry Smith const PetscInt *r, *c, *rout, *cout;
1352c733ed4SBarry Smith PetscInt nz, idx, idt, j, i, oidx, ii, ic, ir;
1362c733ed4SBarry Smith const PetscInt bs = A->rmap->bs, bs2 = a->bs2;
1372c733ed4SBarry Smith const MatScalar *aa = a->a, *v;
1382c733ed4SBarry Smith PetscScalar s1, s2, s3, s4, s5, s6, s7, x1, x2, x3, x4, x5, x6, x7, *x, *t;
1392c733ed4SBarry Smith const PetscScalar *b;
1402c733ed4SBarry Smith
1412c733ed4SBarry Smith PetscFunctionBegin;
1429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b));
1439566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
1442c733ed4SBarry Smith t = a->solve_work;
1452c733ed4SBarry Smith
1469371c9d4SSatish Balay PetscCall(ISGetIndices(isrow, &rout));
1479371c9d4SSatish Balay r = rout;
1489371c9d4SSatish Balay PetscCall(ISGetIndices(iscol, &cout));
1499371c9d4SSatish Balay c = cout;
1502c733ed4SBarry Smith
1512c733ed4SBarry Smith /* copy b into temp work space according to permutation */
1522c733ed4SBarry Smith for (i = 0; i < n; i++) {
1539371c9d4SSatish Balay ii = bs * i;
1549371c9d4SSatish Balay ic = bs * c[i];
1559371c9d4SSatish Balay t[ii] = b[ic];
1569371c9d4SSatish Balay t[ii + 1] = b[ic + 1];
1579371c9d4SSatish Balay t[ii + 2] = b[ic + 2];
1589371c9d4SSatish Balay t[ii + 3] = b[ic + 3];
1599371c9d4SSatish Balay t[ii + 4] = b[ic + 4];
1609371c9d4SSatish Balay t[ii + 5] = b[ic + 5];
1619371c9d4SSatish Balay t[ii + 6] = b[ic + 6];
1622c733ed4SBarry Smith }
1632c733ed4SBarry Smith
1642c733ed4SBarry Smith /* forward solve the U^T */
1652c733ed4SBarry Smith idx = 0;
1662c733ed4SBarry Smith for (i = 0; i < n; i++) {
1672c733ed4SBarry Smith v = aa + bs2 * diag[i];
1682c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */
1699371c9d4SSatish Balay x1 = t[idx];
1709371c9d4SSatish Balay x2 = t[1 + idx];
1719371c9d4SSatish Balay x3 = t[2 + idx];
1729371c9d4SSatish Balay x4 = t[3 + idx];
1739371c9d4SSatish Balay x5 = t[4 + idx];
1749371c9d4SSatish Balay x6 = t[5 + idx];
1759371c9d4SSatish Balay x7 = t[6 + idx];
1762c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6 + v[6] * x7;
1772c733ed4SBarry Smith s2 = v[7] * x1 + v[8] * x2 + v[9] * x3 + v[10] * x4 + v[11] * x5 + v[12] * x6 + v[13] * x7;
1782c733ed4SBarry Smith s3 = v[14] * x1 + v[15] * x2 + v[16] * x3 + v[17] * x4 + v[18] * x5 + v[19] * x6 + v[20] * x7;
1792c733ed4SBarry Smith s4 = v[21] * x1 + v[22] * x2 + v[23] * x3 + v[24] * x4 + v[25] * x5 + v[26] * x6 + v[27] * x7;
1802c733ed4SBarry Smith s5 = v[28] * x1 + v[29] * x2 + v[30] * x3 + v[31] * x4 + v[32] * x5 + v[33] * x6 + v[34] * x7;
1812c733ed4SBarry Smith s6 = v[35] * x1 + v[36] * x2 + v[37] * x3 + v[38] * x4 + v[39] * x5 + v[40] * x6 + v[41] * x7;
1822c733ed4SBarry Smith s7 = v[42] * x1 + v[43] * x2 + v[44] * x3 + v[45] * x4 + v[46] * x5 + v[47] * x6 + v[48] * x7;
1832c733ed4SBarry Smith v -= bs2;
1842c733ed4SBarry Smith
1852c733ed4SBarry Smith vi = aj + diag[i] - 1;
1862c733ed4SBarry Smith nz = diag[i] - diag[i + 1] - 1;
1872c733ed4SBarry Smith for (j = 0; j > -nz; j--) {
1882c733ed4SBarry Smith oidx = bs * vi[j];
1892c733ed4SBarry Smith t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7;
1902c733ed4SBarry Smith t[oidx + 1] -= v[7] * s1 + v[8] * s2 + v[9] * s3 + v[10] * s4 + v[11] * s5 + v[12] * s6 + v[13] * s7;
1912c733ed4SBarry Smith t[oidx + 2] -= v[14] * s1 + v[15] * s2 + v[16] * s3 + v[17] * s4 + v[18] * s5 + v[19] * s6 + v[20] * s7;
1922c733ed4SBarry Smith t[oidx + 3] -= v[21] * s1 + v[22] * s2 + v[23] * s3 + v[24] * s4 + v[25] * s5 + v[26] * s6 + v[27] * s7;
1932c733ed4SBarry Smith t[oidx + 4] -= v[28] * s1 + v[29] * s2 + v[30] * s3 + v[31] * s4 + v[32] * s5 + v[33] * s6 + v[34] * s7;
1942c733ed4SBarry Smith t[oidx + 5] -= v[35] * s1 + v[36] * s2 + v[37] * s3 + v[38] * s4 + v[39] * s5 + v[40] * s6 + v[41] * s7;
1952c733ed4SBarry Smith t[oidx + 6] -= v[42] * s1 + v[43] * s2 + v[44] * s3 + v[45] * s4 + v[46] * s5 + v[47] * s6 + v[48] * s7;
1962c733ed4SBarry Smith v -= bs2;
1972c733ed4SBarry Smith }
1989371c9d4SSatish Balay t[idx] = s1;
1999371c9d4SSatish Balay t[1 + idx] = s2;
2009371c9d4SSatish Balay t[2 + idx] = s3;
2019371c9d4SSatish Balay t[3 + idx] = s4;
2029371c9d4SSatish Balay t[4 + idx] = s5;
2039371c9d4SSatish Balay t[5 + idx] = s6;
2049371c9d4SSatish Balay t[6 + idx] = s7;
2052c733ed4SBarry Smith idx += bs;
2062c733ed4SBarry Smith }
2072c733ed4SBarry Smith /* backward solve the L^T */
2082c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) {
2092c733ed4SBarry Smith v = aa + bs2 * ai[i];
2102c733ed4SBarry Smith vi = aj + ai[i];
2112c733ed4SBarry Smith nz = ai[i + 1] - ai[i];
2122c733ed4SBarry Smith idt = bs * i;
2139371c9d4SSatish Balay s1 = t[idt];
2149371c9d4SSatish Balay s2 = t[1 + idt];
2159371c9d4SSatish Balay s3 = t[2 + idt];
2169371c9d4SSatish Balay s4 = t[3 + idt];
2179371c9d4SSatish Balay s5 = t[4 + idt];
2189371c9d4SSatish Balay s6 = t[5 + idt];
2199371c9d4SSatish Balay s7 = t[6 + idt];
2202c733ed4SBarry Smith for (j = 0; j < nz; j++) {
2212c733ed4SBarry Smith idx = bs * vi[j];
2222c733ed4SBarry Smith t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7;
2232c733ed4SBarry Smith t[idx + 1] -= v[7] * s1 + v[8] * s2 + v[9] * s3 + v[10] * s4 + v[11] * s5 + v[12] * s6 + v[13] * s7;
2242c733ed4SBarry Smith t[idx + 2] -= v[14] * s1 + v[15] * s2 + v[16] * s3 + v[17] * s4 + v[18] * s5 + v[19] * s6 + v[20] * s7;
2252c733ed4SBarry Smith t[idx + 3] -= v[21] * s1 + v[22] * s2 + v[23] * s3 + v[24] * s4 + v[25] * s5 + v[26] * s6 + v[27] * s7;
2262c733ed4SBarry Smith t[idx + 4] -= v[28] * s1 + v[29] * s2 + v[30] * s3 + v[31] * s4 + v[32] * s5 + v[33] * s6 + v[34] * s7;
2272c733ed4SBarry Smith t[idx + 5] -= v[35] * s1 + v[36] * s2 + v[37] * s3 + v[38] * s4 + v[39] * s5 + v[40] * s6 + v[41] * s7;
2282c733ed4SBarry Smith t[idx + 6] -= v[42] * s1 + v[43] * s2 + v[44] * s3 + v[45] * s4 + v[46] * s5 + v[47] * s6 + v[48] * s7;
2292c733ed4SBarry Smith v += bs2;
2302c733ed4SBarry Smith }
2312c733ed4SBarry Smith }
2322c733ed4SBarry Smith
2332c733ed4SBarry Smith /* copy t into x according to permutation */
2342c733ed4SBarry Smith for (i = 0; i < n; i++) {
2359371c9d4SSatish Balay ii = bs * i;
2369371c9d4SSatish Balay ir = bs * r[i];
2379371c9d4SSatish Balay x[ir] = t[ii];
2389371c9d4SSatish Balay x[ir + 1] = t[ii + 1];
2399371c9d4SSatish Balay x[ir + 2] = t[ii + 2];
2409371c9d4SSatish Balay x[ir + 3] = t[ii + 3];
2419371c9d4SSatish Balay x[ir + 4] = t[ii + 4];
2429371c9d4SSatish Balay x[ir + 5] = t[ii + 5];
2439371c9d4SSatish Balay x[ir + 6] = t[ii + 6];
2442c733ed4SBarry Smith }
2452c733ed4SBarry Smith
2469566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout));
2479566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout));
2489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b));
2499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
2509566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
251*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2522c733ed4SBarry Smith }
253