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_4_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, x1, x2, x3, x4, *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 = 4 * 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 ii += 4; 332c733ed4SBarry Smith } 342c733ed4SBarry Smith 352c733ed4SBarry Smith /* forward solve the U^T */ 362c733ed4SBarry Smith idx = 0; 372c733ed4SBarry Smith for (i = 0; i < n; i++) { 382c733ed4SBarry Smith v = aa + 16 * diag[i]; 392c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 40*9371c9d4SSatish Balay x1 = t[idx]; 41*9371c9d4SSatish Balay x2 = t[1 + idx]; 42*9371c9d4SSatish Balay x3 = t[2 + idx]; 43*9371c9d4SSatish Balay x4 = t[3 + idx]; 442c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4; 452c733ed4SBarry Smith s2 = v[4] * x1 + v[5] * x2 + v[6] * x3 + v[7] * x4; 462c733ed4SBarry Smith s3 = v[8] * x1 + v[9] * x2 + v[10] * x3 + v[11] * x4; 472c733ed4SBarry Smith s4 = v[12] * x1 + v[13] * x2 + v[14] * x3 + v[15] * x4; 482c733ed4SBarry Smith v += 16; 492c733ed4SBarry Smith 502c733ed4SBarry Smith vi = aj + diag[i] + 1; 512c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1; 522c733ed4SBarry Smith while (nz--) { 532c733ed4SBarry Smith oidx = 4 * (*vi++); 542c733ed4SBarry Smith t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4; 552c733ed4SBarry Smith t[oidx + 1] -= v[4] * s1 + v[5] * s2 + v[6] * s3 + v[7] * s4; 562c733ed4SBarry Smith t[oidx + 2] -= v[8] * s1 + v[9] * s2 + v[10] * s3 + v[11] * s4; 572c733ed4SBarry Smith t[oidx + 3] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4; 582c733ed4SBarry Smith v += 16; 592c733ed4SBarry Smith } 60*9371c9d4SSatish Balay t[idx] = s1; 61*9371c9d4SSatish Balay t[1 + idx] = s2; 62*9371c9d4SSatish Balay t[2 + idx] = s3; 63*9371c9d4SSatish Balay t[3 + idx] = s4; 642c733ed4SBarry Smith idx += 4; 652c733ed4SBarry Smith } 662c733ed4SBarry Smith /* backward solve the L^T */ 672c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 682c733ed4SBarry Smith v = aa + 16 * diag[i] - 16; 692c733ed4SBarry Smith vi = aj + diag[i] - 1; 702c733ed4SBarry Smith nz = diag[i] - ai[i]; 712c733ed4SBarry Smith idt = 4 * i; 72*9371c9d4SSatish Balay s1 = t[idt]; 73*9371c9d4SSatish Balay s2 = t[1 + idt]; 74*9371c9d4SSatish Balay s3 = t[2 + idt]; 75*9371c9d4SSatish Balay s4 = t[3 + idt]; 762c733ed4SBarry Smith while (nz--) { 772c733ed4SBarry Smith idx = 4 * (*vi--); 782c733ed4SBarry Smith t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4; 792c733ed4SBarry Smith t[idx + 1] -= v[4] * s1 + v[5] * s2 + v[6] * s3 + v[7] * s4; 802c733ed4SBarry Smith t[idx + 2] -= v[8] * s1 + v[9] * s2 + v[10] * s3 + v[11] * s4; 812c733ed4SBarry Smith t[idx + 3] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4; 822c733ed4SBarry Smith v -= 16; 832c733ed4SBarry Smith } 842c733ed4SBarry Smith } 852c733ed4SBarry Smith 862c733ed4SBarry Smith /* copy t into x according to permutation */ 872c733ed4SBarry Smith ii = 0; 882c733ed4SBarry Smith for (i = 0; i < n; i++) { 892c733ed4SBarry Smith ir = 4 * r[i]; 902c733ed4SBarry Smith x[ir] = t[ii]; 912c733ed4SBarry Smith x[ir + 1] = t[ii + 1]; 922c733ed4SBarry Smith x[ir + 2] = t[ii + 2]; 932c733ed4SBarry Smith x[ir + 3] = t[ii + 3]; 942c733ed4SBarry Smith ii += 4; 952c733ed4SBarry Smith } 962c733ed4SBarry Smith 979566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout)); 989566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout)); 999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 1009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 1019566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 16 * (a->nz) - 4.0 * A->cmap->n)); 1022c733ed4SBarry Smith PetscFunctionReturn(0); 1032c733ed4SBarry Smith } 1042c733ed4SBarry Smith 105*9371c9d4SSatish Balay PetscErrorCode MatSolveTranspose_SeqBAIJ_4(Mat A, Vec bb, Vec xx) { 1062c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 1072c733ed4SBarry Smith IS iscol = a->col, isrow = a->row; 1082c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag; 1092c733ed4SBarry Smith const PetscInt *r, *c, *rout, *cout; 1102c733ed4SBarry Smith PetscInt nz, idx, idt, j, i, oidx, ii, ic, ir; 1112c733ed4SBarry Smith const PetscInt bs = A->rmap->bs, bs2 = a->bs2; 1122c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 1132c733ed4SBarry Smith PetscScalar s1, s2, s3, s4, x1, x2, x3, x4, *x, *t; 1142c733ed4SBarry Smith const PetscScalar *b; 1152c733ed4SBarry Smith 1162c733ed4SBarry Smith PetscFunctionBegin; 1179566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 1189566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 1192c733ed4SBarry Smith t = a->solve_work; 1202c733ed4SBarry Smith 121*9371c9d4SSatish Balay PetscCall(ISGetIndices(isrow, &rout)); 122*9371c9d4SSatish Balay r = rout; 123*9371c9d4SSatish Balay PetscCall(ISGetIndices(iscol, &cout)); 124*9371c9d4SSatish Balay c = cout; 1252c733ed4SBarry Smith 1262c733ed4SBarry Smith /* copy b into temp work space according to permutation */ 1272c733ed4SBarry Smith for (i = 0; i < n; i++) { 128*9371c9d4SSatish Balay ii = bs * i; 129*9371c9d4SSatish Balay ic = bs * c[i]; 130*9371c9d4SSatish Balay t[ii] = b[ic]; 131*9371c9d4SSatish Balay t[ii + 1] = b[ic + 1]; 132*9371c9d4SSatish Balay t[ii + 2] = b[ic + 2]; 133*9371c9d4SSatish Balay t[ii + 3] = b[ic + 3]; 1342c733ed4SBarry Smith } 1352c733ed4SBarry Smith 1362c733ed4SBarry Smith /* forward solve the U^T */ 1372c733ed4SBarry Smith idx = 0; 1382c733ed4SBarry Smith for (i = 0; i < n; i++) { 1392c733ed4SBarry Smith v = aa + bs2 * diag[i]; 1402c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 141*9371c9d4SSatish Balay x1 = t[idx]; 142*9371c9d4SSatish Balay x2 = t[1 + idx]; 143*9371c9d4SSatish Balay x3 = t[2 + idx]; 144*9371c9d4SSatish Balay x4 = t[3 + idx]; 1452c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4; 1462c733ed4SBarry Smith s2 = v[4] * x1 + v[5] * x2 + v[6] * x3 + v[7] * x4; 1472c733ed4SBarry Smith s3 = v[8] * x1 + v[9] * x2 + v[10] * x3 + v[11] * x4; 1482c733ed4SBarry Smith s4 = v[12] * x1 + v[13] * x2 + v[14] * x3 + v[15] * x4; 1492c733ed4SBarry Smith v -= bs2; 1502c733ed4SBarry Smith 1512c733ed4SBarry Smith vi = aj + diag[i] - 1; 1522c733ed4SBarry Smith nz = diag[i] - diag[i + 1] - 1; 1532c733ed4SBarry Smith for (j = 0; j > -nz; j--) { 1542c733ed4SBarry Smith oidx = bs * vi[j]; 1552c733ed4SBarry Smith t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4; 1562c733ed4SBarry Smith t[oidx + 1] -= v[4] * s1 + v[5] * s2 + v[6] * s3 + v[7] * s4; 1572c733ed4SBarry Smith t[oidx + 2] -= v[8] * s1 + v[9] * s2 + v[10] * s3 + v[11] * s4; 1582c733ed4SBarry Smith t[oidx + 3] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4; 1592c733ed4SBarry Smith v -= bs2; 1602c733ed4SBarry Smith } 161*9371c9d4SSatish Balay t[idx] = s1; 162*9371c9d4SSatish Balay t[1 + idx] = s2; 163*9371c9d4SSatish Balay t[2 + idx] = s3; 164*9371c9d4SSatish Balay t[3 + idx] = s4; 1652c733ed4SBarry Smith idx += bs; 1662c733ed4SBarry Smith } 1672c733ed4SBarry Smith /* backward solve the L^T */ 1682c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 1692c733ed4SBarry Smith v = aa + bs2 * ai[i]; 1702c733ed4SBarry Smith vi = aj + ai[i]; 1712c733ed4SBarry Smith nz = ai[i + 1] - ai[i]; 1722c733ed4SBarry Smith idt = bs * i; 173*9371c9d4SSatish Balay s1 = t[idt]; 174*9371c9d4SSatish Balay s2 = t[1 + idt]; 175*9371c9d4SSatish Balay s3 = t[2 + idt]; 176*9371c9d4SSatish Balay s4 = t[3 + idt]; 1772c733ed4SBarry Smith for (j = 0; j < nz; j++) { 1782c733ed4SBarry Smith idx = bs * vi[j]; 1792c733ed4SBarry Smith t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4; 1802c733ed4SBarry Smith t[idx + 1] -= v[4] * s1 + v[5] * s2 + v[6] * s3 + v[7] * s4; 1812c733ed4SBarry Smith t[idx + 2] -= v[8] * s1 + v[9] * s2 + v[10] * s3 + v[11] * s4; 1822c733ed4SBarry Smith t[idx + 3] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4; 1832c733ed4SBarry Smith v += bs2; 1842c733ed4SBarry Smith } 1852c733ed4SBarry Smith } 1862c733ed4SBarry Smith 1872c733ed4SBarry Smith /* copy t into x according to permutation */ 1882c733ed4SBarry Smith for (i = 0; i < n; i++) { 189*9371c9d4SSatish Balay ii = bs * i; 190*9371c9d4SSatish Balay ir = bs * r[i]; 191*9371c9d4SSatish Balay x[ir] = t[ii]; 192*9371c9d4SSatish Balay x[ir + 1] = t[ii + 1]; 193*9371c9d4SSatish Balay x[ir + 2] = t[ii + 2]; 194*9371c9d4SSatish Balay x[ir + 3] = t[ii + 3]; 1952c733ed4SBarry Smith } 1962c733ed4SBarry Smith 1979566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout)); 1989566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout)); 1999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 2009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 2019566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n)); 2022c733ed4SBarry Smith PetscFunctionReturn(0); 2032c733ed4SBarry Smith } 204