xref: /petsc/src/mat/impls/baij/seq/baijsolvtrann.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h>
32c733ed4SBarry Smith 
42c733ed4SBarry Smith /* ----------------------------------------------------------- */
5*9371c9d4SSatish Balay PetscErrorCode MatSolveTranspose_SeqBAIJ_N_inplace(Mat A, Vec bb, Vec xx) {
62c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
72c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
82c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout, *ai = a->i, *aj = a->j, *vi;
92c733ed4SBarry Smith   PetscInt           i, nz, j;
102c733ed4SBarry Smith   const PetscInt     n = a->mbs, bs = A->rmap->bs, bs2 = a->bs2;
112c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
122c733ed4SBarry Smith   PetscScalar       *x, *t, *ls;
132c733ed4SBarry Smith   const PetscScalar *b;
142c733ed4SBarry Smith 
152c733ed4SBarry Smith   PetscFunctionBegin;
169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
182c733ed4SBarry Smith   t = a->solve_work;
192c733ed4SBarry Smith 
20*9371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
21*9371c9d4SSatish Balay   r = rout;
22*9371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
23*9371c9d4SSatish 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++) {
27*9371c9d4SSatish Balay     for (j = 0; j < bs; j++) { t[i * bs + j] = b[c[i] * bs + j]; }
282c733ed4SBarry Smith   }
292c733ed4SBarry Smith 
302c733ed4SBarry Smith   /* forward solve the upper triangular transpose */
312c733ed4SBarry Smith   ls = a->solve_work + A->cmap->n;
322c733ed4SBarry Smith   for (i = 0; i < n; i++) {
339566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
342c733ed4SBarry Smith     PetscKernel_w_gets_transA_times_v(bs, ls, aa + bs2 * a->diag[i], t + i * bs);
352c733ed4SBarry Smith     v  = aa + bs2 * (a->diag[i] + 1);
362c733ed4SBarry Smith     vi = aj + a->diag[i] + 1;
372c733ed4SBarry Smith     nz = ai[i + 1] - a->diag[i] - 1;
382c733ed4SBarry Smith     while (nz--) {
392c733ed4SBarry Smith       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (*vi++), v, t + i * bs);
402c733ed4SBarry Smith       v += bs2;
412c733ed4SBarry Smith     }
422c733ed4SBarry Smith   }
432c733ed4SBarry Smith 
442c733ed4SBarry Smith   /* backward solve the lower triangular transpose */
452c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
462c733ed4SBarry Smith     v  = aa + bs2 * ai[i];
472c733ed4SBarry Smith     vi = aj + ai[i];
482c733ed4SBarry Smith     nz = a->diag[i] - ai[i];
492c733ed4SBarry Smith     while (nz--) {
502c733ed4SBarry Smith       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (*vi++), v, t + i * bs);
512c733ed4SBarry Smith       v += bs2;
522c733ed4SBarry Smith     }
532c733ed4SBarry Smith   }
542c733ed4SBarry Smith 
552c733ed4SBarry Smith   /* copy t into x according to permutation */
562c733ed4SBarry Smith   for (i = 0; i < n; i++) {
57*9371c9d4SSatish Balay     for (j = 0; j < bs; j++) { x[bs * r[i] + j] = t[bs * i + j]; }
582c733ed4SBarry Smith   }
592c733ed4SBarry Smith 
609566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
619566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
649566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
652c733ed4SBarry Smith   PetscFunctionReturn(0);
662c733ed4SBarry Smith }
672c733ed4SBarry Smith 
68*9371c9d4SSatish Balay PetscErrorCode MatSolveTranspose_SeqBAIJ_N(Mat A, Vec bb, Vec xx) {
692c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
702c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
712c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout;
722c733ed4SBarry Smith   const PetscInt     n = a->mbs, *ai = a->i, *aj = a->j, *vi, *diag = a->diag;
732c733ed4SBarry Smith   PetscInt           i, j, nz;
742c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
752c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
762c733ed4SBarry Smith   PetscScalar       *x, *t, *ls;
772c733ed4SBarry Smith   const PetscScalar *b;
782c733ed4SBarry Smith 
792c733ed4SBarry Smith   PetscFunctionBegin;
809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
819566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
822c733ed4SBarry Smith   t = a->solve_work;
832c733ed4SBarry Smith 
84*9371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
85*9371c9d4SSatish Balay   r = rout;
86*9371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
87*9371c9d4SSatish Balay   c = cout;
882c733ed4SBarry Smith 
892c733ed4SBarry Smith   /* copy the b into temp work space according to permutation */
902c733ed4SBarry Smith   for (i = 0; i < n; i++) {
91*9371c9d4SSatish Balay     for (j = 0; j < bs; j++) { t[i * bs + j] = b[c[i] * bs + j]; }
922c733ed4SBarry Smith   }
932c733ed4SBarry Smith 
942c733ed4SBarry Smith   /* forward solve the upper triangular transpose */
952c733ed4SBarry Smith   ls = a->solve_work + A->cmap->n;
962c733ed4SBarry Smith   for (i = 0; i < n; i++) {
979566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
982c733ed4SBarry Smith     PetscKernel_w_gets_transA_times_v(bs, ls, aa + bs2 * diag[i], t + i * bs);
992c733ed4SBarry Smith     v  = aa + bs2 * (diag[i] - 1);
1002c733ed4SBarry Smith     vi = aj + diag[i] - 1;
1012c733ed4SBarry Smith     nz = diag[i] - diag[i + 1] - 1;
1022c733ed4SBarry Smith     for (j = 0; j > -nz; j--) {
1032c733ed4SBarry Smith       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (vi[j]), v, t + i * bs);
1042c733ed4SBarry Smith       v -= bs2;
1052c733ed4SBarry Smith     }
1062c733ed4SBarry Smith   }
1072c733ed4SBarry Smith 
1082c733ed4SBarry Smith   /* backward solve the lower triangular transpose */
1092c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1102c733ed4SBarry Smith     v  = aa + bs2 * ai[i];
1112c733ed4SBarry Smith     vi = aj + ai[i];
1122c733ed4SBarry Smith     nz = ai[i + 1] - ai[i];
1132c733ed4SBarry Smith     for (j = 0; j < nz; j++) {
1142c733ed4SBarry Smith       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (vi[j]), v, t + i * bs);
1152c733ed4SBarry Smith       v += bs2;
1162c733ed4SBarry Smith     }
1172c733ed4SBarry Smith   }
1182c733ed4SBarry Smith 
1192c733ed4SBarry Smith   /* copy t into x according to permutation */
1202c733ed4SBarry Smith   for (i = 0; i < n; i++) {
121*9371c9d4SSatish Balay     for (j = 0; j < bs; j++) { x[bs * r[i] + j] = t[bs * i + j]; }
1222c733ed4SBarry Smith   }
1232c733ed4SBarry Smith 
1249566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
1259566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
1269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1289566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
1292c733ed4SBarry Smith   PetscFunctionReturn(0);
1302c733ed4SBarry Smith }
131