12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith
MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)3d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx)
4d71ae5a4SJacob Faibussowitsch {
52c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
62c733ed4SBarry Smith PetscInt i, nz, idx, idt, oidx;
72c733ed4SBarry Smith const PetscInt *diag = a->diag, *vi, n = a->mbs, *ai = a->i, *aj = a->j;
82c733ed4SBarry Smith const MatScalar *aa = a->a, *v;
92c733ed4SBarry Smith PetscScalar s1, s2, x1, x2, *x;
102c733ed4SBarry Smith
112c733ed4SBarry Smith PetscFunctionBegin;
129566063dSJacob Faibussowitsch PetscCall(VecCopy(bb, xx));
139566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
142c733ed4SBarry Smith
152c733ed4SBarry Smith /* forward solve the U^T */
162c733ed4SBarry Smith idx = 0;
172c733ed4SBarry Smith for (i = 0; i < n; i++) {
182c733ed4SBarry Smith v = aa + 4 * diag[i];
192c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */
209371c9d4SSatish Balay x1 = x[idx];
219371c9d4SSatish Balay x2 = x[1 + idx];
222c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2;
232c733ed4SBarry Smith s2 = v[2] * x1 + v[3] * x2;
242c733ed4SBarry Smith v += 4;
252c733ed4SBarry Smith
262c733ed4SBarry Smith vi = aj + diag[i] + 1;
272c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1;
282c733ed4SBarry Smith while (nz--) {
292c733ed4SBarry Smith oidx = 2 * (*vi++);
302c733ed4SBarry Smith x[oidx] -= v[0] * s1 + v[1] * s2;
312c733ed4SBarry Smith x[oidx + 1] -= v[2] * s1 + v[3] * s2;
322c733ed4SBarry Smith v += 4;
332c733ed4SBarry Smith }
349371c9d4SSatish Balay x[idx] = s1;
359371c9d4SSatish Balay x[1 + idx] = s2;
362c733ed4SBarry Smith idx += 2;
372c733ed4SBarry Smith }
382c733ed4SBarry Smith /* backward solve the L^T */
392c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) {
402c733ed4SBarry Smith v = aa + 4 * diag[i] - 4;
412c733ed4SBarry Smith vi = aj + diag[i] - 1;
422c733ed4SBarry Smith nz = diag[i] - ai[i];
432c733ed4SBarry Smith idt = 2 * i;
449371c9d4SSatish Balay s1 = x[idt];
459371c9d4SSatish Balay s2 = x[1 + idt];
462c733ed4SBarry Smith while (nz--) {
472c733ed4SBarry Smith idx = 2 * (*vi--);
482c733ed4SBarry Smith x[idx] -= v[0] * s1 + v[1] * s2;
492c733ed4SBarry Smith x[idx + 1] -= v[2] * s1 + v[3] * s2;
502c733ed4SBarry Smith v -= 4;
512c733ed4SBarry Smith }
522c733ed4SBarry Smith }
539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
549566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 4.0 * (a->nz) - 2.0 * A->cmap->n));
55*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
562c733ed4SBarry Smith }
572c733ed4SBarry Smith
MatSolveTranspose_SeqBAIJ_2_NaturalOrdering(Mat A,Vec bb,Vec xx)58d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_2_NaturalOrdering(Mat A, Vec bb, Vec xx)
59d71ae5a4SJacob Faibussowitsch {
602c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
612c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
622c733ed4SBarry Smith PetscInt nz, idx, idt, j, i, oidx;
632c733ed4SBarry Smith const PetscInt bs = A->rmap->bs, bs2 = a->bs2;
642c733ed4SBarry Smith const MatScalar *aa = a->a, *v;
652c733ed4SBarry Smith PetscScalar s1, s2, x1, x2, *x;
662c733ed4SBarry Smith
672c733ed4SBarry Smith PetscFunctionBegin;
689566063dSJacob Faibussowitsch PetscCall(VecCopy(bb, xx));
699566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
702c733ed4SBarry Smith
712c733ed4SBarry Smith /* forward solve the U^T */
722c733ed4SBarry Smith idx = 0;
732c733ed4SBarry Smith for (i = 0; i < n; i++) {
742c733ed4SBarry Smith v = aa + bs2 * diag[i];
752c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */
769371c9d4SSatish Balay x1 = x[idx];
779371c9d4SSatish Balay x2 = x[1 + idx];
782c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2;
792c733ed4SBarry Smith s2 = v[2] * x1 + v[3] * x2;
802c733ed4SBarry Smith v -= bs2;
812c733ed4SBarry Smith
822c733ed4SBarry Smith vi = aj + diag[i] - 1;
832c733ed4SBarry Smith nz = diag[i] - diag[i + 1] - 1;
842c733ed4SBarry Smith for (j = 0; j > -nz; j--) {
852c733ed4SBarry Smith oidx = bs * vi[j];
862c733ed4SBarry Smith x[oidx] -= v[0] * s1 + v[1] * s2;
872c733ed4SBarry Smith x[oidx + 1] -= v[2] * s1 + v[3] * s2;
882c733ed4SBarry Smith v -= bs2;
892c733ed4SBarry Smith }
909371c9d4SSatish Balay x[idx] = s1;
919371c9d4SSatish Balay x[1 + idx] = s2;
922c733ed4SBarry Smith idx += bs;
932c733ed4SBarry Smith }
942c733ed4SBarry Smith /* backward solve the L^T */
952c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) {
962c733ed4SBarry Smith v = aa + bs2 * ai[i];
972c733ed4SBarry Smith vi = aj + ai[i];
982c733ed4SBarry Smith nz = ai[i + 1] - ai[i];
992c733ed4SBarry Smith idt = bs * i;
1009371c9d4SSatish Balay s1 = x[idt];
1019371c9d4SSatish Balay s2 = x[1 + idt];
1022c733ed4SBarry Smith for (j = 0; j < nz; j++) {
1032c733ed4SBarry Smith idx = bs * vi[j];
1042c733ed4SBarry Smith x[idx] -= v[0] * s1 + v[1] * s2;
1052c733ed4SBarry Smith x[idx + 1] -= v[2] * s1 + v[3] * s2;
1062c733ed4SBarry Smith v += bs2;
1072c733ed4SBarry Smith }
1082c733ed4SBarry Smith }
1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
1109566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
111*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1122c733ed4SBarry Smith }
113