12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h> 22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h> 32c733ed4SBarry Smith 4d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_1(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 *rout, *cout, *r, *c, *adiag = a->diag, *ai = a->i, *aj = a->j, *vi; 92c733ed4SBarry Smith PetscInt i, n = a->mbs, j; 102c733ed4SBarry Smith PetscInt nz; 112c733ed4SBarry Smith PetscScalar *x, *tmp, s1; 122c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 132c733ed4SBarry Smith const PetscScalar *b; 142c733ed4SBarry Smith 152c733ed4SBarry Smith PetscFunctionBegin; 169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 179566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 182c733ed4SBarry Smith tmp = 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 for (i = 0; i < n; i++) tmp[i] = b[c[i]]; 272c733ed4SBarry Smith 282c733ed4SBarry Smith /* forward solve the U^T */ 292c733ed4SBarry Smith for (i = 0; i < n; i++) { 302c733ed4SBarry Smith v = aa + adiag[i + 1] + 1; 312c733ed4SBarry Smith vi = aj + adiag[i + 1] + 1; 322c733ed4SBarry Smith nz = adiag[i] - adiag[i + 1] - 1; 332c733ed4SBarry Smith s1 = tmp[i]; 342c733ed4SBarry Smith s1 *= v[nz]; /* multiply by inverse of diagonal entry */ 352c733ed4SBarry Smith for (j = 0; j < nz; j++) tmp[vi[j]] -= s1 * v[j]; 362c733ed4SBarry Smith tmp[i] = s1; 372c733ed4SBarry Smith } 382c733ed4SBarry Smith 392c733ed4SBarry Smith /* backward solve the L^T */ 402c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 412c733ed4SBarry Smith v = aa + ai[i]; 422c733ed4SBarry Smith vi = aj + ai[i]; 432c733ed4SBarry Smith nz = ai[i + 1] - ai[i]; 442c733ed4SBarry Smith s1 = tmp[i]; 452c733ed4SBarry Smith for (j = 0; j < nz; j++) tmp[vi[j]] -= s1 * v[j]; 462c733ed4SBarry Smith } 472c733ed4SBarry Smith 482c733ed4SBarry Smith /* copy tmp into x according to permutation */ 492c733ed4SBarry Smith for (i = 0; i < n; i++) x[r[i]] = tmp[i]; 502c733ed4SBarry Smith 519566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout)); 529566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout)); 539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 552c733ed4SBarry Smith 569566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * a->nz - A->cmap->n)); 57*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 582c733ed4SBarry Smith } 592c733ed4SBarry Smith 60d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_1_inplace(Mat A, Vec bb, Vec xx) 61d71ae5a4SJacob Faibussowitsch { 622c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 632c733ed4SBarry Smith IS iscol = a->col, isrow = a->row; 642c733ed4SBarry Smith const PetscInt *r, *c, *rout, *cout; 652c733ed4SBarry Smith const PetscInt *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j; 662c733ed4SBarry Smith PetscInt i, nz; 672c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 682c733ed4SBarry Smith PetscScalar s1, *x, *t; 692c733ed4SBarry Smith const PetscScalar *b; 702c733ed4SBarry Smith 712c733ed4SBarry Smith PetscFunctionBegin; 729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 739566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 742c733ed4SBarry Smith t = a->solve_work; 752c733ed4SBarry Smith 769371c9d4SSatish Balay PetscCall(ISGetIndices(isrow, &rout)); 779371c9d4SSatish Balay r = rout; 789371c9d4SSatish Balay PetscCall(ISGetIndices(iscol, &cout)); 799371c9d4SSatish Balay c = cout; 802c733ed4SBarry Smith 812c733ed4SBarry Smith /* copy the b into temp work space according to permutation */ 822c733ed4SBarry Smith for (i = 0; i < n; i++) t[i] = b[c[i]]; 832c733ed4SBarry Smith 842c733ed4SBarry Smith /* forward solve the U^T */ 852c733ed4SBarry Smith for (i = 0; i < n; i++) { 862c733ed4SBarry Smith v = aa + diag[i]; 872c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */ 882c733ed4SBarry Smith s1 = (*v++) * t[i]; 892c733ed4SBarry Smith vi = aj + diag[i] + 1; 902c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1; 91ad540459SPierre Jolivet while (nz--) t[*vi++] -= (*v++) * s1; 922c733ed4SBarry Smith t[i] = s1; 932c733ed4SBarry Smith } 942c733ed4SBarry Smith /* backward solve the L^T */ 952c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 962c733ed4SBarry Smith v = aa + diag[i] - 1; 972c733ed4SBarry Smith vi = aj + diag[i] - 1; 982c733ed4SBarry Smith nz = diag[i] - ai[i]; 992c733ed4SBarry Smith s1 = t[i]; 100ad540459SPierre Jolivet while (nz--) t[*vi--] -= (*v--) * s1; 1012c733ed4SBarry Smith } 1022c733ed4SBarry Smith 1032c733ed4SBarry Smith /* copy t into x according to permutation */ 1042c733ed4SBarry Smith for (i = 0; i < n; i++) x[r[i]] = t[i]; 1052c733ed4SBarry Smith 1069566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &rout)); 1079566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &cout)); 1089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 1109566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * (a->nz) - A->cmap->n)); 111*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1122c733ed4SBarry Smith } 113