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_7_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, s6, s7, x1, x2, x3, x4, x5, x6, x7, *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 = 7 * 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 t[ii + 5] = b[ic + 5]; 342c733ed4SBarry Smith t[ii + 6] = b[ic + 6]; 352c733ed4SBarry Smith ii += 7; 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 + 49 * diag[i]; 422c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 43*9371c9d4SSatish Balay x1 = t[idx]; 44*9371c9d4SSatish Balay x2 = t[1 + idx]; 45*9371c9d4SSatish Balay x3 = t[2 + idx]; 46*9371c9d4SSatish Balay x4 = t[3 + idx]; 47*9371c9d4SSatish Balay x5 = t[4 + idx]; 48*9371c9d4SSatish Balay x6 = t[5 + idx]; 49*9371c9d4SSatish Balay x7 = t[6 + idx]; 502c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6 + v[6] * x7; 512c733ed4SBarry Smith s2 = v[7] * x1 + v[8] * x2 + v[9] * x3 + v[10] * x4 + v[11] * x5 + v[12] * x6 + v[13] * x7; 522c733ed4SBarry Smith s3 = v[14] * x1 + v[15] * x2 + v[16] * x3 + v[17] * x4 + v[18] * x5 + v[19] * x6 + v[20] * x7; 532c733ed4SBarry Smith s4 = v[21] * x1 + v[22] * x2 + v[23] * x3 + v[24] * x4 + v[25] * x5 + v[26] * x6 + v[27] * x7; 542c733ed4SBarry Smith s5 = v[28] * x1 + v[29] * x2 + v[30] * x3 + v[31] * x4 + v[32] * x5 + v[33] * x6 + v[34] * x7; 552c733ed4SBarry Smith s6 = v[35] * x1 + v[36] * x2 + v[37] * x3 + v[38] * x4 + v[39] * x5 + v[40] * x6 + v[41] * x7; 562c733ed4SBarry Smith s7 = v[42] * x1 + v[43] * x2 + v[44] * x3 + v[45] * x4 + v[46] * x5 + v[47] * x6 + v[48] * x7; 572c733ed4SBarry Smith v += 49; 582c733ed4SBarry Smith 592c733ed4SBarry Smith vi = aj + diag[i] + 1; 602c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1; 612c733ed4SBarry Smith while (nz--) { 622c733ed4SBarry Smith oidx = 7 * (*vi++); 632c733ed4SBarry Smith t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7; 642c733ed4SBarry 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; 652c733ed4SBarry 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; 662c733ed4SBarry 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; 672c733ed4SBarry 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; 682c733ed4SBarry 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; 692c733ed4SBarry 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; 702c733ed4SBarry Smith v += 49; 712c733ed4SBarry Smith } 72*9371c9d4SSatish Balay t[idx] = s1; 73*9371c9d4SSatish Balay t[1 + idx] = s2; 74*9371c9d4SSatish Balay t[2 + idx] = s3; 75*9371c9d4SSatish Balay t[3 + idx] = s4; 76*9371c9d4SSatish Balay t[4 + idx] = s5; 77*9371c9d4SSatish Balay t[5 + idx] = s6; 78*9371c9d4SSatish Balay t[6 + idx] = s7; 792c733ed4SBarry Smith idx += 7; 802c733ed4SBarry Smith } 812c733ed4SBarry Smith /* backward solve the L^T */ 822c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 832c733ed4SBarry Smith v = aa + 49 * diag[i] - 49; 842c733ed4SBarry Smith vi = aj + diag[i] - 1; 852c733ed4SBarry Smith nz = diag[i] - ai[i]; 862c733ed4SBarry Smith idt = 7 * i; 87*9371c9d4SSatish Balay s1 = t[idt]; 88*9371c9d4SSatish Balay s2 = t[1 + idt]; 89*9371c9d4SSatish Balay s3 = t[2 + idt]; 90*9371c9d4SSatish Balay s4 = t[3 + idt]; 91*9371c9d4SSatish Balay s5 = t[4 + idt]; 92*9371c9d4SSatish Balay s6 = t[5 + idt]; 93*9371c9d4SSatish Balay s7 = t[6 + idt]; 942c733ed4SBarry Smith while (nz--) { 952c733ed4SBarry Smith idx = 7 * (*vi--); 962c733ed4SBarry Smith t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7; 972c733ed4SBarry 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; 982c733ed4SBarry 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; 992c733ed4SBarry 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; 1002c733ed4SBarry 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; 1012c733ed4SBarry 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; 1022c733ed4SBarry 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; 1032c733ed4SBarry Smith v -= 49; 1042c733ed4SBarry Smith } 1052c733ed4SBarry Smith } 1062c733ed4SBarry Smith 1072c733ed4SBarry Smith /* copy t into x according to permutation */ 1082c733ed4SBarry Smith ii = 0; 1092c733ed4SBarry Smith for (i = 0; i < n; i++) { 1102c733ed4SBarry Smith ir = 7 * r[i]; 1112c733ed4SBarry Smith x[ir] = t[ii]; 1122c733ed4SBarry Smith x[ir + 1] = t[ii + 1]; 1132c733ed4SBarry Smith x[ir + 2] = t[ii + 2]; 1142c733ed4SBarry Smith x[ir + 3] = t[ii + 3]; 1152c733ed4SBarry Smith x[ir + 4] = t[ii + 4]; 1162c733ed4SBarry Smith x[ir + 5] = t[ii + 5]; 1172c733ed4SBarry Smith x[ir + 6] = t[ii + 6]; 1182c733ed4SBarry Smith ii += 7; 1192c733ed4SBarry Smith } 1202c733ed4SBarry Smith 1219566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout)); 1229566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout)); 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 1249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 1259566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 49 * (a->nz) - 7.0 * A->cmap->n)); 1262c733ed4SBarry Smith PetscFunctionReturn(0); 1272c733ed4SBarry Smith } 128*9371c9d4SSatish Balay PetscErrorCode MatSolveTranspose_SeqBAIJ_7(Mat A, Vec bb, Vec xx) { 1292c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 1302c733ed4SBarry Smith IS iscol = a->col, isrow = a->row; 1312c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag; 1322c733ed4SBarry Smith const PetscInt *r, *c, *rout, *cout; 1332c733ed4SBarry Smith PetscInt nz, idx, idt, j, i, oidx, ii, ic, ir; 1342c733ed4SBarry Smith const PetscInt bs = A->rmap->bs, bs2 = a->bs2; 1352c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 1362c733ed4SBarry Smith PetscScalar s1, s2, s3, s4, s5, s6, s7, x1, x2, x3, x4, x5, x6, x7, *x, *t; 1372c733ed4SBarry Smith const PetscScalar *b; 1382c733ed4SBarry Smith 1392c733ed4SBarry Smith PetscFunctionBegin; 1409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 1419566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 1422c733ed4SBarry Smith t = a->solve_work; 1432c733ed4SBarry Smith 144*9371c9d4SSatish Balay PetscCall(ISGetIndices(isrow, &rout)); 145*9371c9d4SSatish Balay r = rout; 146*9371c9d4SSatish Balay PetscCall(ISGetIndices(iscol, &cout)); 147*9371c9d4SSatish Balay c = cout; 1482c733ed4SBarry Smith 1492c733ed4SBarry Smith /* copy b into temp work space according to permutation */ 1502c733ed4SBarry Smith for (i = 0; i < n; i++) { 151*9371c9d4SSatish Balay ii = bs * i; 152*9371c9d4SSatish Balay ic = bs * c[i]; 153*9371c9d4SSatish Balay t[ii] = b[ic]; 154*9371c9d4SSatish Balay t[ii + 1] = b[ic + 1]; 155*9371c9d4SSatish Balay t[ii + 2] = b[ic + 2]; 156*9371c9d4SSatish Balay t[ii + 3] = b[ic + 3]; 157*9371c9d4SSatish Balay t[ii + 4] = b[ic + 4]; 158*9371c9d4SSatish Balay t[ii + 5] = b[ic + 5]; 159*9371c9d4SSatish Balay t[ii + 6] = b[ic + 6]; 1602c733ed4SBarry Smith } 1612c733ed4SBarry Smith 1622c733ed4SBarry Smith /* forward solve the U^T */ 1632c733ed4SBarry Smith idx = 0; 1642c733ed4SBarry Smith for (i = 0; i < n; i++) { 1652c733ed4SBarry Smith v = aa + bs2 * diag[i]; 1662c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 167*9371c9d4SSatish Balay x1 = t[idx]; 168*9371c9d4SSatish Balay x2 = t[1 + idx]; 169*9371c9d4SSatish Balay x3 = t[2 + idx]; 170*9371c9d4SSatish Balay x4 = t[3 + idx]; 171*9371c9d4SSatish Balay x5 = t[4 + idx]; 172*9371c9d4SSatish Balay x6 = t[5 + idx]; 173*9371c9d4SSatish Balay x7 = t[6 + idx]; 1742c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6 + v[6] * x7; 1752c733ed4SBarry Smith s2 = v[7] * x1 + v[8] * x2 + v[9] * x3 + v[10] * x4 + v[11] * x5 + v[12] * x6 + v[13] * x7; 1762c733ed4SBarry Smith s3 = v[14] * x1 + v[15] * x2 + v[16] * x3 + v[17] * x4 + v[18] * x5 + v[19] * x6 + v[20] * x7; 1772c733ed4SBarry Smith s4 = v[21] * x1 + v[22] * x2 + v[23] * x3 + v[24] * x4 + v[25] * x5 + v[26] * x6 + v[27] * x7; 1782c733ed4SBarry Smith s5 = v[28] * x1 + v[29] * x2 + v[30] * x3 + v[31] * x4 + v[32] * x5 + v[33] * x6 + v[34] * x7; 1792c733ed4SBarry Smith s6 = v[35] * x1 + v[36] * x2 + v[37] * x3 + v[38] * x4 + v[39] * x5 + v[40] * x6 + v[41] * x7; 1802c733ed4SBarry Smith s7 = v[42] * x1 + v[43] * x2 + v[44] * x3 + v[45] * x4 + v[46] * x5 + v[47] * x6 + v[48] * x7; 1812c733ed4SBarry Smith v -= bs2; 1822c733ed4SBarry Smith 1832c733ed4SBarry Smith vi = aj + diag[i] - 1; 1842c733ed4SBarry Smith nz = diag[i] - diag[i + 1] - 1; 1852c733ed4SBarry Smith for (j = 0; j > -nz; j--) { 1862c733ed4SBarry Smith oidx = bs * vi[j]; 1872c733ed4SBarry Smith t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7; 1882c733ed4SBarry 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; 1892c733ed4SBarry 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; 1902c733ed4SBarry 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; 1912c733ed4SBarry 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; 1922c733ed4SBarry 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; 1932c733ed4SBarry 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; 1942c733ed4SBarry Smith v -= bs2; 1952c733ed4SBarry Smith } 196*9371c9d4SSatish Balay t[idx] = s1; 197*9371c9d4SSatish Balay t[1 + idx] = s2; 198*9371c9d4SSatish Balay t[2 + idx] = s3; 199*9371c9d4SSatish Balay t[3 + idx] = s4; 200*9371c9d4SSatish Balay t[4 + idx] = s5; 201*9371c9d4SSatish Balay t[5 + idx] = s6; 202*9371c9d4SSatish Balay t[6 + idx] = s7; 2032c733ed4SBarry Smith idx += bs; 2042c733ed4SBarry Smith } 2052c733ed4SBarry Smith /* backward solve the L^T */ 2062c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 2072c733ed4SBarry Smith v = aa + bs2 * ai[i]; 2082c733ed4SBarry Smith vi = aj + ai[i]; 2092c733ed4SBarry Smith nz = ai[i + 1] - ai[i]; 2102c733ed4SBarry Smith idt = bs * i; 211*9371c9d4SSatish Balay s1 = t[idt]; 212*9371c9d4SSatish Balay s2 = t[1 + idt]; 213*9371c9d4SSatish Balay s3 = t[2 + idt]; 214*9371c9d4SSatish Balay s4 = t[3 + idt]; 215*9371c9d4SSatish Balay s5 = t[4 + idt]; 216*9371c9d4SSatish Balay s6 = t[5 + idt]; 217*9371c9d4SSatish Balay s7 = t[6 + idt]; 2182c733ed4SBarry Smith for (j = 0; j < nz; j++) { 2192c733ed4SBarry Smith idx = bs * vi[j]; 2202c733ed4SBarry Smith t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7; 2212c733ed4SBarry 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; 2222c733ed4SBarry 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; 2232c733ed4SBarry 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; 2242c733ed4SBarry 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; 2252c733ed4SBarry 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; 2262c733ed4SBarry 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; 2272c733ed4SBarry Smith v += bs2; 2282c733ed4SBarry Smith } 2292c733ed4SBarry Smith } 2302c733ed4SBarry Smith 2312c733ed4SBarry Smith /* copy t into x according to permutation */ 2322c733ed4SBarry Smith for (i = 0; i < n; i++) { 233*9371c9d4SSatish Balay ii = bs * i; 234*9371c9d4SSatish Balay ir = bs * r[i]; 235*9371c9d4SSatish Balay x[ir] = t[ii]; 236*9371c9d4SSatish Balay x[ir + 1] = t[ii + 1]; 237*9371c9d4SSatish Balay x[ir + 2] = t[ii + 2]; 238*9371c9d4SSatish Balay x[ir + 3] = t[ii + 3]; 239*9371c9d4SSatish Balay x[ir + 4] = t[ii + 4]; 240*9371c9d4SSatish Balay x[ir + 5] = t[ii + 5]; 241*9371c9d4SSatish Balay x[ir + 6] = t[ii + 6]; 2422c733ed4SBarry Smith } 2432c733ed4SBarry Smith 2449566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout)); 2459566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout)); 2469566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 2479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 2489566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n)); 2492c733ed4SBarry Smith PetscFunctionReturn(0); 2502c733ed4SBarry Smith } 251