xref: /petsc/src/mat/impls/baij/seq/baijsolvtrann.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h>
32c733ed4SBarry Smith 
42c733ed4SBarry Smith /* ----------------------------------------------------------- */
5d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_N_inplace(Mat A, Vec bb, Vec xx)
6d71ae5a4SJacob Faibussowitsch {
72c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
82c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
92c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout, *ai = a->i, *aj = a->j, *vi;
102c733ed4SBarry Smith   PetscInt           i, nz, j;
112c733ed4SBarry Smith   const PetscInt     n = a->mbs, bs = A->rmap->bs, bs2 = a->bs2;
122c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
132c733ed4SBarry Smith   PetscScalar       *x, *t, *ls;
142c733ed4SBarry Smith   const PetscScalar *b;
152c733ed4SBarry Smith 
162c733ed4SBarry Smith   PetscFunctionBegin;
179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
189566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
192c733ed4SBarry Smith   t = a->solve_work;
202c733ed4SBarry Smith 
219371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
229371c9d4SSatish Balay   r = rout;
239371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
249371c9d4SSatish Balay   c = cout;
252c733ed4SBarry Smith 
262c733ed4SBarry Smith   /* copy the b into temp work space according to permutation */
272c733ed4SBarry Smith   for (i = 0; i < n; i++) {
28ad540459SPierre Jolivet     for (j = 0; j < bs; j++) t[i * bs + j] = b[c[i] * bs + j];
292c733ed4SBarry Smith   }
302c733ed4SBarry Smith 
312c733ed4SBarry Smith   /* forward solve the upper triangular transpose */
322c733ed4SBarry Smith   ls = a->solve_work + A->cmap->n;
332c733ed4SBarry Smith   for (i = 0; i < n; i++) {
349566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
352c733ed4SBarry Smith     PetscKernel_w_gets_transA_times_v(bs, ls, aa + bs2 * a->diag[i], t + i * bs);
362c733ed4SBarry Smith     v  = aa + bs2 * (a->diag[i] + 1);
372c733ed4SBarry Smith     vi = aj + a->diag[i] + 1;
382c733ed4SBarry Smith     nz = ai[i + 1] - a->diag[i] - 1;
392c733ed4SBarry Smith     while (nz--) {
402c733ed4SBarry Smith       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (*vi++), v, t + i * bs);
412c733ed4SBarry Smith       v += bs2;
422c733ed4SBarry Smith     }
432c733ed4SBarry Smith   }
442c733ed4SBarry Smith 
452c733ed4SBarry Smith   /* backward solve the lower triangular transpose */
462c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
472c733ed4SBarry Smith     v  = aa + bs2 * ai[i];
482c733ed4SBarry Smith     vi = aj + ai[i];
492c733ed4SBarry Smith     nz = a->diag[i] - ai[i];
502c733ed4SBarry Smith     while (nz--) {
512c733ed4SBarry Smith       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (*vi++), v, t + i * bs);
522c733ed4SBarry Smith       v += bs2;
532c733ed4SBarry Smith     }
542c733ed4SBarry Smith   }
552c733ed4SBarry Smith 
562c733ed4SBarry Smith   /* copy t into x according to permutation */
572c733ed4SBarry Smith   for (i = 0; i < n; i++) {
58ad540459SPierre Jolivet     for (j = 0; j < bs; j++) x[bs * r[i] + j] = t[bs * i + j];
592c733ed4SBarry Smith   }
602c733ed4SBarry Smith 
619566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
629566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
659566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
66*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
672c733ed4SBarry Smith }
682c733ed4SBarry Smith 
69d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_N(Mat A, Vec bb, Vec xx)
70d71ae5a4SJacob Faibussowitsch {
712c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
722c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
732c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout;
742c733ed4SBarry Smith   const PetscInt     n = a->mbs, *ai = a->i, *aj = a->j, *vi, *diag = a->diag;
752c733ed4SBarry Smith   PetscInt           i, j, nz;
762c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
772c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
782c733ed4SBarry Smith   PetscScalar       *x, *t, *ls;
792c733ed4SBarry Smith   const PetscScalar *b;
802c733ed4SBarry Smith 
812c733ed4SBarry Smith   PetscFunctionBegin;
829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
839566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
842c733ed4SBarry Smith   t = a->solve_work;
852c733ed4SBarry Smith 
869371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
879371c9d4SSatish Balay   r = rout;
889371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
899371c9d4SSatish Balay   c = cout;
902c733ed4SBarry Smith 
912c733ed4SBarry Smith   /* copy the b into temp work space according to permutation */
922c733ed4SBarry Smith   for (i = 0; i < n; i++) {
93ad540459SPierre Jolivet     for (j = 0; j < bs; j++) t[i * bs + j] = b[c[i] * bs + j];
942c733ed4SBarry Smith   }
952c733ed4SBarry Smith 
962c733ed4SBarry Smith   /* forward solve the upper triangular transpose */
972c733ed4SBarry Smith   ls = a->solve_work + A->cmap->n;
982c733ed4SBarry Smith   for (i = 0; i < n; i++) {
999566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
1002c733ed4SBarry Smith     PetscKernel_w_gets_transA_times_v(bs, ls, aa + bs2 * diag[i], t + i * bs);
1012c733ed4SBarry Smith     v  = aa + bs2 * (diag[i] - 1);
1022c733ed4SBarry Smith     vi = aj + diag[i] - 1;
1032c733ed4SBarry Smith     nz = diag[i] - diag[i + 1] - 1;
1042c733ed4SBarry Smith     for (j = 0; j > -nz; j--) {
1052c733ed4SBarry Smith       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (vi[j]), v, t + i * bs);
1062c733ed4SBarry Smith       v -= bs2;
1072c733ed4SBarry Smith     }
1082c733ed4SBarry Smith   }
1092c733ed4SBarry Smith 
1102c733ed4SBarry Smith   /* backward solve the lower triangular transpose */
1112c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1122c733ed4SBarry Smith     v  = aa + bs2 * ai[i];
1132c733ed4SBarry Smith     vi = aj + ai[i];
1142c733ed4SBarry Smith     nz = ai[i + 1] - ai[i];
1152c733ed4SBarry Smith     for (j = 0; j < nz; j++) {
1162c733ed4SBarry Smith       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (vi[j]), v, t + i * bs);
1172c733ed4SBarry Smith       v += bs2;
1182c733ed4SBarry Smith     }
1192c733ed4SBarry Smith   }
1202c733ed4SBarry Smith 
1212c733ed4SBarry Smith   /* copy t into x according to permutation */
1222c733ed4SBarry Smith   for (i = 0; i < n; i++) {
123ad540459SPierre Jolivet     for (j = 0; j < bs; j++) x[bs * r[i] + j] = t[bs * i + j];
1242c733ed4SBarry Smith   }
1252c733ed4SBarry Smith 
1269566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
1279566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
1289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1309566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
131*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1322c733ed4SBarry Smith }
133