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_2_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, x1, x2, *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 = 2 * c[i]; 282c733ed4SBarry Smith t[ii] = b[ic]; 292c733ed4SBarry Smith t[ii + 1] = b[ic + 1]; 302c733ed4SBarry Smith ii += 2; 312c733ed4SBarry Smith } 322c733ed4SBarry Smith 332c733ed4SBarry Smith /* forward solve the U^T */ 342c733ed4SBarry Smith idx = 0; 352c733ed4SBarry Smith for (i = 0; i < n; i++) { 362c733ed4SBarry Smith v = aa + 4 * diag[i]; 372c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 38*9371c9d4SSatish Balay x1 = t[idx]; 39*9371c9d4SSatish Balay x2 = t[1 + idx]; 402c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2; 412c733ed4SBarry Smith s2 = v[2] * x1 + v[3] * x2; 422c733ed4SBarry Smith v += 4; 432c733ed4SBarry Smith 442c733ed4SBarry Smith vi = aj + diag[i] + 1; 452c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1; 462c733ed4SBarry Smith while (nz--) { 472c733ed4SBarry Smith oidx = 2 * (*vi++); 482c733ed4SBarry Smith t[oidx] -= v[0] * s1 + v[1] * s2; 492c733ed4SBarry Smith t[oidx + 1] -= v[2] * s1 + v[3] * s2; 502c733ed4SBarry Smith v += 4; 512c733ed4SBarry Smith } 52*9371c9d4SSatish Balay t[idx] = s1; 53*9371c9d4SSatish Balay t[1 + idx] = s2; 542c733ed4SBarry Smith idx += 2; 552c733ed4SBarry Smith } 562c733ed4SBarry Smith /* backward solve the L^T */ 572c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 582c733ed4SBarry Smith v = aa + 4 * diag[i] - 4; 592c733ed4SBarry Smith vi = aj + diag[i] - 1; 602c733ed4SBarry Smith nz = diag[i] - ai[i]; 612c733ed4SBarry Smith idt = 2 * i; 62*9371c9d4SSatish Balay s1 = t[idt]; 63*9371c9d4SSatish Balay s2 = t[1 + idt]; 642c733ed4SBarry Smith while (nz--) { 652c733ed4SBarry Smith idx = 2 * (*vi--); 662c733ed4SBarry Smith t[idx] -= v[0] * s1 + v[1] * s2; 672c733ed4SBarry Smith t[idx + 1] -= v[2] * s1 + v[3] * s2; 682c733ed4SBarry Smith v -= 4; 692c733ed4SBarry Smith } 702c733ed4SBarry Smith } 712c733ed4SBarry Smith 722c733ed4SBarry Smith /* copy t into x according to permutation */ 732c733ed4SBarry Smith ii = 0; 742c733ed4SBarry Smith for (i = 0; i < n; i++) { 752c733ed4SBarry Smith ir = 2 * r[i]; 762c733ed4SBarry Smith x[ir] = t[ii]; 772c733ed4SBarry Smith x[ir + 1] = t[ii + 1]; 782c733ed4SBarry Smith ii += 2; 792c733ed4SBarry Smith } 802c733ed4SBarry Smith 819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout)); 829566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout)); 839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 859566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 4 * (a->nz) - 2.0 * A->cmap->n)); 862c733ed4SBarry Smith PetscFunctionReturn(0); 872c733ed4SBarry Smith } 882c733ed4SBarry Smith 89*9371c9d4SSatish Balay PetscErrorCode MatSolveTranspose_SeqBAIJ_2(Mat A, Vec bb, Vec xx) { 902c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 912c733ed4SBarry Smith IS iscol = a->col, isrow = a->row; 922c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag; 932c733ed4SBarry Smith const PetscInt *r, *c, *rout, *cout; 942c733ed4SBarry Smith PetscInt nz, idx, idt, j, i, oidx, ii, ic, ir; 952c733ed4SBarry Smith const PetscInt bs = A->rmap->bs, bs2 = a->bs2; 962c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 972c733ed4SBarry Smith PetscScalar s1, s2, x1, x2, *x, *t; 982c733ed4SBarry Smith const PetscScalar *b; 992c733ed4SBarry Smith 1002c733ed4SBarry Smith PetscFunctionBegin; 1019566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 1029566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 1032c733ed4SBarry Smith t = a->solve_work; 1042c733ed4SBarry Smith 105*9371c9d4SSatish Balay PetscCall(ISGetIndices(isrow, &rout)); 106*9371c9d4SSatish Balay r = rout; 107*9371c9d4SSatish Balay PetscCall(ISGetIndices(iscol, &cout)); 108*9371c9d4SSatish Balay c = cout; 1092c733ed4SBarry Smith 1102c733ed4SBarry Smith /* copy b into temp work space according to permutation */ 1112c733ed4SBarry Smith for (i = 0; i < n; i++) { 112*9371c9d4SSatish Balay ii = bs * i; 113*9371c9d4SSatish Balay ic = bs * c[i]; 114*9371c9d4SSatish Balay t[ii] = b[ic]; 115*9371c9d4SSatish Balay t[ii + 1] = b[ic + 1]; 1162c733ed4SBarry Smith } 1172c733ed4SBarry Smith 1182c733ed4SBarry Smith /* forward solve the U^T */ 1192c733ed4SBarry Smith idx = 0; 1202c733ed4SBarry Smith for (i = 0; i < n; i++) { 1212c733ed4SBarry Smith v = aa + bs2 * diag[i]; 1222c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 123*9371c9d4SSatish Balay x1 = t[idx]; 124*9371c9d4SSatish Balay x2 = t[1 + idx]; 1252c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2; 1262c733ed4SBarry Smith s2 = v[2] * x1 + v[3] * x2; 1272c733ed4SBarry Smith v -= bs2; 1282c733ed4SBarry Smith 1292c733ed4SBarry Smith vi = aj + diag[i] - 1; 1302c733ed4SBarry Smith nz = diag[i] - diag[i + 1] - 1; 1312c733ed4SBarry Smith for (j = 0; j > -nz; j--) { 1322c733ed4SBarry Smith oidx = bs * vi[j]; 1332c733ed4SBarry Smith t[oidx] -= v[0] * s1 + v[1] * s2; 1342c733ed4SBarry Smith t[oidx + 1] -= v[2] * s1 + v[3] * s2; 1352c733ed4SBarry Smith v -= bs2; 1362c733ed4SBarry Smith } 137*9371c9d4SSatish Balay t[idx] = s1; 138*9371c9d4SSatish Balay t[1 + idx] = s2; 1392c733ed4SBarry Smith idx += bs; 1402c733ed4SBarry Smith } 1412c733ed4SBarry Smith /* backward solve the L^T */ 1422c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 1432c733ed4SBarry Smith v = aa + bs2 * ai[i]; 1442c733ed4SBarry Smith vi = aj + ai[i]; 1452c733ed4SBarry Smith nz = ai[i + 1] - ai[i]; 1462c733ed4SBarry Smith idt = bs * i; 147*9371c9d4SSatish Balay s1 = t[idt]; 148*9371c9d4SSatish Balay s2 = t[1 + idt]; 1492c733ed4SBarry Smith for (j = 0; j < nz; j++) { 1502c733ed4SBarry Smith idx = bs * vi[j]; 1512c733ed4SBarry Smith t[idx] -= v[0] * s1 + v[1] * s2; 1522c733ed4SBarry Smith t[idx + 1] -= v[2] * s1 + v[3] * s2; 1532c733ed4SBarry Smith v += bs2; 1542c733ed4SBarry Smith } 1552c733ed4SBarry Smith } 1562c733ed4SBarry Smith 1572c733ed4SBarry Smith /* copy t into x according to permutation */ 1582c733ed4SBarry Smith for (i = 0; i < n; i++) { 159*9371c9d4SSatish Balay ii = bs * i; 160*9371c9d4SSatish Balay ir = bs * r[i]; 161*9371c9d4SSatish Balay x[ir] = t[ii]; 162*9371c9d4SSatish Balay x[ir + 1] = t[ii + 1]; 1632c733ed4SBarry Smith } 1642c733ed4SBarry Smith 1659566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout)); 1669566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout)); 1679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 1689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 1699566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n)); 1702c733ed4SBarry Smith PetscFunctionReturn(0); 1712c733ed4SBarry Smith } 172