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_3_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, x1, x2, x3, *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 = 3 * 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 ii += 3; 322c733ed4SBarry Smith } 332c733ed4SBarry Smith 342c733ed4SBarry Smith /* forward solve the U^T */ 352c733ed4SBarry Smith idx = 0; 362c733ed4SBarry Smith for (i = 0; i < n; i++) { 372c733ed4SBarry Smith v = aa + 9 * diag[i]; 382c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 39*9371c9d4SSatish Balay x1 = t[idx]; 40*9371c9d4SSatish Balay x2 = t[1 + idx]; 41*9371c9d4SSatish Balay x3 = t[2 + idx]; 422c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2 + v[2] * x3; 432c733ed4SBarry Smith s2 = v[3] * x1 + v[4] * x2 + v[5] * x3; 442c733ed4SBarry Smith s3 = v[6] * x1 + v[7] * x2 + v[8] * x3; 452c733ed4SBarry Smith v += 9; 462c733ed4SBarry Smith 472c733ed4SBarry Smith vi = aj + diag[i] + 1; 482c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1; 492c733ed4SBarry Smith while (nz--) { 502c733ed4SBarry Smith oidx = 3 * (*vi++); 512c733ed4SBarry Smith t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3; 522c733ed4SBarry Smith t[oidx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3; 532c733ed4SBarry Smith t[oidx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3; 542c733ed4SBarry Smith v += 9; 552c733ed4SBarry Smith } 56*9371c9d4SSatish Balay t[idx] = s1; 57*9371c9d4SSatish Balay t[1 + idx] = s2; 58*9371c9d4SSatish Balay t[2 + idx] = s3; 592c733ed4SBarry Smith idx += 3; 602c733ed4SBarry Smith } 612c733ed4SBarry Smith /* backward solve the L^T */ 622c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 632c733ed4SBarry Smith v = aa + 9 * diag[i] - 9; 642c733ed4SBarry Smith vi = aj + diag[i] - 1; 652c733ed4SBarry Smith nz = diag[i] - ai[i]; 662c733ed4SBarry Smith idt = 3 * i; 67*9371c9d4SSatish Balay s1 = t[idt]; 68*9371c9d4SSatish Balay s2 = t[1 + idt]; 69*9371c9d4SSatish Balay s3 = t[2 + idt]; 702c733ed4SBarry Smith while (nz--) { 712c733ed4SBarry Smith idx = 3 * (*vi--); 722c733ed4SBarry Smith t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3; 732c733ed4SBarry Smith t[idx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3; 742c733ed4SBarry Smith t[idx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3; 752c733ed4SBarry Smith v -= 9; 762c733ed4SBarry Smith } 772c733ed4SBarry Smith } 782c733ed4SBarry Smith 792c733ed4SBarry Smith /* copy t into x according to permutation */ 802c733ed4SBarry Smith ii = 0; 812c733ed4SBarry Smith for (i = 0; i < n; i++) { 822c733ed4SBarry Smith ir = 3 * r[i]; 832c733ed4SBarry Smith x[ir] = t[ii]; 842c733ed4SBarry Smith x[ir + 1] = t[ii + 1]; 852c733ed4SBarry Smith x[ir + 2] = t[ii + 2]; 862c733ed4SBarry Smith ii += 3; 872c733ed4SBarry Smith } 882c733ed4SBarry Smith 899566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout)); 909566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout)); 919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 929566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 939566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 9 * (a->nz) - 3.0 * A->cmap->n)); 942c733ed4SBarry Smith PetscFunctionReturn(0); 952c733ed4SBarry Smith } 962c733ed4SBarry Smith 97*9371c9d4SSatish Balay PetscErrorCode MatSolveTranspose_SeqBAIJ_3(Mat A, Vec bb, Vec xx) { 982c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 992c733ed4SBarry Smith IS iscol = a->col, isrow = a->row; 1002c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag; 1012c733ed4SBarry Smith const PetscInt *r, *c, *rout, *cout; 1022c733ed4SBarry Smith PetscInt nz, idx, idt, j, i, oidx, ii, ic, ir; 1032c733ed4SBarry Smith const PetscInt bs = A->rmap->bs, bs2 = a->bs2; 1042c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 1052c733ed4SBarry Smith PetscScalar s1, s2, s3, x1, x2, x3, *x, *t; 1062c733ed4SBarry Smith const PetscScalar *b; 1072c733ed4SBarry Smith 1082c733ed4SBarry Smith PetscFunctionBegin; 1099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 1109566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 1112c733ed4SBarry Smith t = a->solve_work; 1122c733ed4SBarry Smith 113*9371c9d4SSatish Balay PetscCall(ISGetIndices(isrow, &rout)); 114*9371c9d4SSatish Balay r = rout; 115*9371c9d4SSatish Balay PetscCall(ISGetIndices(iscol, &cout)); 116*9371c9d4SSatish Balay c = cout; 1172c733ed4SBarry Smith 1182c733ed4SBarry Smith /* copy b into temp work space according to permutation */ 1192c733ed4SBarry Smith for (i = 0; i < n; i++) { 120*9371c9d4SSatish Balay ii = bs * i; 121*9371c9d4SSatish Balay ic = bs * c[i]; 122*9371c9d4SSatish Balay t[ii] = b[ic]; 123*9371c9d4SSatish Balay t[ii + 1] = b[ic + 1]; 124*9371c9d4SSatish Balay t[ii + 2] = b[ic + 2]; 1252c733ed4SBarry Smith } 1262c733ed4SBarry Smith 1272c733ed4SBarry Smith /* forward solve the U^T */ 1282c733ed4SBarry Smith idx = 0; 1292c733ed4SBarry Smith for (i = 0; i < n; i++) { 1302c733ed4SBarry Smith v = aa + bs2 * diag[i]; 1312c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 132*9371c9d4SSatish Balay x1 = t[idx]; 133*9371c9d4SSatish Balay x2 = t[1 + idx]; 134*9371c9d4SSatish Balay x3 = t[2 + idx]; 1352c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2 + v[2] * x3; 1362c733ed4SBarry Smith s2 = v[3] * x1 + v[4] * x2 + v[5] * x3; 1372c733ed4SBarry Smith s3 = v[6] * x1 + v[7] * x2 + v[8] * x3; 1382c733ed4SBarry Smith v -= bs2; 1392c733ed4SBarry Smith 1402c733ed4SBarry Smith vi = aj + diag[i] - 1; 1412c733ed4SBarry Smith nz = diag[i] - diag[i + 1] - 1; 1422c733ed4SBarry Smith for (j = 0; j > -nz; j--) { 1432c733ed4SBarry Smith oidx = bs * vi[j]; 1442c733ed4SBarry Smith t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3; 1452c733ed4SBarry Smith t[oidx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3; 1462c733ed4SBarry Smith t[oidx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3; 1472c733ed4SBarry Smith v -= bs2; 1482c733ed4SBarry Smith } 149*9371c9d4SSatish Balay t[idx] = s1; 150*9371c9d4SSatish Balay t[1 + idx] = s2; 151*9371c9d4SSatish Balay t[2 + idx] = s3; 1522c733ed4SBarry Smith idx += bs; 1532c733ed4SBarry Smith } 1542c733ed4SBarry Smith /* backward solve the L^T */ 1552c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 1562c733ed4SBarry Smith v = aa + bs2 * ai[i]; 1572c733ed4SBarry Smith vi = aj + ai[i]; 1582c733ed4SBarry Smith nz = ai[i + 1] - ai[i]; 1592c733ed4SBarry Smith idt = bs * i; 160*9371c9d4SSatish Balay s1 = t[idt]; 161*9371c9d4SSatish Balay s2 = t[1 + idt]; 162*9371c9d4SSatish Balay s3 = t[2 + idt]; 1632c733ed4SBarry Smith for (j = 0; j < nz; j++) { 1642c733ed4SBarry Smith idx = bs * vi[j]; 1652c733ed4SBarry Smith t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3; 1662c733ed4SBarry Smith t[idx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3; 1672c733ed4SBarry Smith t[idx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3; 1682c733ed4SBarry Smith v += bs2; 1692c733ed4SBarry Smith } 1702c733ed4SBarry Smith } 1712c733ed4SBarry Smith 1722c733ed4SBarry Smith /* copy t into x according to permutation */ 1732c733ed4SBarry Smith for (i = 0; i < n; i++) { 174*9371c9d4SSatish Balay ii = bs * i; 175*9371c9d4SSatish Balay ir = bs * r[i]; 176*9371c9d4SSatish Balay x[ir] = t[ii]; 177*9371c9d4SSatish Balay x[ir + 1] = t[ii + 1]; 178*9371c9d4SSatish Balay x[ir + 2] = t[ii + 2]; 1792c733ed4SBarry Smith } 1802c733ed4SBarry Smith 1819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout)); 1829566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout)); 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 1849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 1859566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n)); 1862c733ed4SBarry Smith PetscFunctionReturn(0); 1872c733ed4SBarry Smith } 188