xref: /petsc/src/mat/impls/baij/seq/baijsolvtran1.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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_1(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    *rout, *cout, *r, *c, *adiag = a->diag, *ai = a->i, *aj = a->j, *vi;
82c733ed4SBarry Smith   PetscInt           i, n = a->mbs, j;
92c733ed4SBarry Smith   PetscInt           nz;
102c733ed4SBarry Smith   PetscScalar       *x, *tmp, s1;
112c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
122c733ed4SBarry Smith   const PetscScalar *b;
132c733ed4SBarry Smith 
142c733ed4SBarry Smith   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
169566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
172c733ed4SBarry Smith   tmp = 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   for (i = 0; i < n; i++) tmp[i] = b[c[i]];
262c733ed4SBarry Smith 
272c733ed4SBarry Smith   /* forward solve the U^T */
282c733ed4SBarry Smith   for (i = 0; i < n; i++) {
292c733ed4SBarry Smith     v  = aa + adiag[i + 1] + 1;
302c733ed4SBarry Smith     vi = aj + adiag[i + 1] + 1;
312c733ed4SBarry Smith     nz = adiag[i] - adiag[i + 1] - 1;
322c733ed4SBarry Smith     s1 = tmp[i];
332c733ed4SBarry Smith     s1 *= v[nz]; /* multiply by inverse of diagonal entry */
342c733ed4SBarry Smith     for (j = 0; j < nz; j++) tmp[vi[j]] -= s1 * v[j];
352c733ed4SBarry Smith     tmp[i] = s1;
362c733ed4SBarry Smith   }
372c733ed4SBarry Smith 
382c733ed4SBarry Smith   /* backward solve the L^T */
392c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
402c733ed4SBarry Smith     v  = aa + ai[i];
412c733ed4SBarry Smith     vi = aj + ai[i];
422c733ed4SBarry Smith     nz = ai[i + 1] - ai[i];
432c733ed4SBarry Smith     s1 = tmp[i];
442c733ed4SBarry Smith     for (j = 0; j < nz; j++) tmp[vi[j]] -= s1 * v[j];
452c733ed4SBarry Smith   }
462c733ed4SBarry Smith 
472c733ed4SBarry Smith   /* copy tmp into x according to permutation */
482c733ed4SBarry Smith   for (i = 0; i < n; i++) x[r[i]] = tmp[i];
492c733ed4SBarry Smith 
509566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
519566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
542c733ed4SBarry Smith 
559566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * a->nz - A->cmap->n));
562c733ed4SBarry Smith   PetscFunctionReturn(0);
572c733ed4SBarry Smith }
582c733ed4SBarry Smith 
59*9371c9d4SSatish Balay PetscErrorCode MatSolveTranspose_SeqBAIJ_1_inplace(Mat A, Vec bb, Vec xx) {
602c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
612c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
622c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout;
632c733ed4SBarry Smith   const PetscInt    *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j;
642c733ed4SBarry Smith   PetscInt           i, nz;
652c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
662c733ed4SBarry Smith   PetscScalar        s1, *x, *t;
672c733ed4SBarry Smith   const PetscScalar *b;
682c733ed4SBarry Smith 
692c733ed4SBarry Smith   PetscFunctionBegin;
709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
719566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
722c733ed4SBarry Smith   t = a->solve_work;
732c733ed4SBarry Smith 
74*9371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
75*9371c9d4SSatish Balay   r = rout;
76*9371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
77*9371c9d4SSatish Balay   c = cout;
782c733ed4SBarry Smith 
792c733ed4SBarry Smith   /* copy the b into temp work space according to permutation */
802c733ed4SBarry Smith   for (i = 0; i < n; i++) t[i] = b[c[i]];
812c733ed4SBarry Smith 
822c733ed4SBarry Smith   /* forward solve the U^T */
832c733ed4SBarry Smith   for (i = 0; i < n; i++) {
842c733ed4SBarry Smith     v  = aa + diag[i];
852c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
862c733ed4SBarry Smith     s1 = (*v++) * t[i];
872c733ed4SBarry Smith     vi = aj + diag[i] + 1;
882c733ed4SBarry Smith     nz = ai[i + 1] - diag[i] - 1;
89*9371c9d4SSatish Balay     while (nz--) { t[*vi++] -= (*v++) * s1; }
902c733ed4SBarry Smith     t[i] = s1;
912c733ed4SBarry Smith   }
922c733ed4SBarry Smith   /* backward solve the L^T */
932c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
942c733ed4SBarry Smith     v  = aa + diag[i] - 1;
952c733ed4SBarry Smith     vi = aj + diag[i] - 1;
962c733ed4SBarry Smith     nz = diag[i] - ai[i];
972c733ed4SBarry Smith     s1 = t[i];
98*9371c9d4SSatish Balay     while (nz--) { t[*vi--] -= (*v--) * s1; }
992c733ed4SBarry Smith   }
1002c733ed4SBarry Smith 
1012c733ed4SBarry Smith   /* copy t into x according to permutation */
1022c733ed4SBarry Smith   for (i = 0; i < n; i++) x[r[i]] = t[i];
1032c733ed4SBarry Smith 
1049566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
1059566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
1069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1089566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * (a->nz) - A->cmap->n));
1092c733ed4SBarry Smith   PetscFunctionReturn(0);
1102c733ed4SBarry Smith }
111