xref: /petsc/src/mat/impls/baij/seq/baijsolvtran2.c (revision 31d78bcd2b98084dc1368b20eb1129c8b9fb39fe)
12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h>
32c733ed4SBarry Smith 
MatSolveTranspose_SeqBAIJ_2_inplace(Mat A,Vec bb,Vec xx)4d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_2_inplace(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    *r, *c, *rout, *cout;
92c733ed4SBarry Smith   const PetscInt    *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j;
102c733ed4SBarry Smith   PetscInt           i, nz, idx, idt, ii, ic, ir, oidx;
112c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
122c733ed4SBarry Smith   PetscScalar        s1, s2, x1, x2, *x, *t;
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 
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   ii = 0;
272c733ed4SBarry Smith   for (i = 0; i < n; i++) {
282c733ed4SBarry Smith     ic        = 2 * c[i];
292c733ed4SBarry Smith     t[ii]     = b[ic];
302c733ed4SBarry Smith     t[ii + 1] = b[ic + 1];
312c733ed4SBarry Smith     ii += 2;
322c733ed4SBarry Smith   }
332c733ed4SBarry Smith 
342c733ed4SBarry Smith   /* forward solve the U^T */
352c733ed4SBarry Smith   idx = 0;
362c733ed4SBarry Smith   for (i = 0; i < n; i++) {
372c733ed4SBarry Smith     v = aa + 4 * diag[i];
382c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
399371c9d4SSatish Balay     x1 = t[idx];
409371c9d4SSatish Balay     x2 = t[1 + idx];
412c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2;
422c733ed4SBarry Smith     s2 = v[2] * x1 + v[3] * x2;
432c733ed4SBarry Smith     v += 4;
442c733ed4SBarry Smith 
452c733ed4SBarry Smith     vi = aj + diag[i] + 1;
462c733ed4SBarry Smith     nz = ai[i + 1] - diag[i] - 1;
472c733ed4SBarry Smith     while (nz--) {
482c733ed4SBarry Smith       oidx = 2 * (*vi++);
492c733ed4SBarry Smith       t[oidx] -= v[0] * s1 + v[1] * s2;
502c733ed4SBarry Smith       t[oidx + 1] -= v[2] * s1 + v[3] * s2;
512c733ed4SBarry Smith       v += 4;
522c733ed4SBarry Smith     }
539371c9d4SSatish Balay     t[idx]     = s1;
549371c9d4SSatish Balay     t[1 + idx] = s2;
552c733ed4SBarry Smith     idx += 2;
562c733ed4SBarry Smith   }
572c733ed4SBarry Smith   /* backward solve the L^T */
582c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
592c733ed4SBarry Smith     v   = aa + 4 * diag[i] - 4;
602c733ed4SBarry Smith     vi  = aj + diag[i] - 1;
612c733ed4SBarry Smith     nz  = diag[i] - ai[i];
622c733ed4SBarry Smith     idt = 2 * i;
639371c9d4SSatish Balay     s1  = t[idt];
649371c9d4SSatish Balay     s2  = t[1 + idt];
652c733ed4SBarry Smith     while (nz--) {
662c733ed4SBarry Smith       idx = 2 * (*vi--);
672c733ed4SBarry Smith       t[idx] -= v[0] * s1 + v[1] * s2;
682c733ed4SBarry Smith       t[idx + 1] -= v[2] * s1 + v[3] * s2;
692c733ed4SBarry Smith       v -= 4;
702c733ed4SBarry Smith     }
712c733ed4SBarry Smith   }
722c733ed4SBarry Smith 
732c733ed4SBarry Smith   /* copy t into x according to permutation */
742c733ed4SBarry Smith   ii = 0;
752c733ed4SBarry Smith   for (i = 0; i < n; i++) {
762c733ed4SBarry Smith     ir        = 2 * r[i];
772c733ed4SBarry Smith     x[ir]     = t[ii];
782c733ed4SBarry Smith     x[ir + 1] = t[ii + 1];
792c733ed4SBarry Smith     ii += 2;
802c733ed4SBarry Smith   }
812c733ed4SBarry Smith 
829566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
839566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
869566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 4 * (a->nz) - 2.0 * A->cmap->n));
87*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
882c733ed4SBarry Smith }
892c733ed4SBarry Smith 
MatSolveTranspose_SeqBAIJ_2(Mat A,Vec bb,Vec xx)90d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_2(Mat A, Vec bb, Vec xx)
91d71ae5a4SJacob Faibussowitsch {
922c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
932c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
942c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
952c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout;
962c733ed4SBarry Smith   PetscInt           nz, idx, idt, j, i, oidx, ii, ic, ir;
972c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
982c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
992c733ed4SBarry Smith   PetscScalar        s1, s2, x1, x2, *x, *t;
1002c733ed4SBarry Smith   const PetscScalar *b;
1012c733ed4SBarry Smith 
1022c733ed4SBarry Smith   PetscFunctionBegin;
1039566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1052c733ed4SBarry Smith   t = a->solve_work;
1062c733ed4SBarry Smith 
1079371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
1089371c9d4SSatish Balay   r = rout;
1099371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
1109371c9d4SSatish Balay   c = cout;
1112c733ed4SBarry Smith 
1122c733ed4SBarry Smith   /* copy b into temp work space according to permutation */
1132c733ed4SBarry Smith   for (i = 0; i < n; i++) {
1149371c9d4SSatish Balay     ii        = bs * i;
1159371c9d4SSatish Balay     ic        = bs * c[i];
1169371c9d4SSatish Balay     t[ii]     = b[ic];
1179371c9d4SSatish Balay     t[ii + 1] = b[ic + 1];
1182c733ed4SBarry Smith   }
1192c733ed4SBarry Smith 
1202c733ed4SBarry Smith   /* forward solve the U^T */
1212c733ed4SBarry Smith   idx = 0;
1222c733ed4SBarry Smith   for (i = 0; i < n; i++) {
1232c733ed4SBarry Smith     v = aa + bs2 * diag[i];
1242c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
1259371c9d4SSatish Balay     x1 = t[idx];
1269371c9d4SSatish Balay     x2 = t[1 + idx];
1272c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2;
1282c733ed4SBarry Smith     s2 = v[2] * x1 + v[3] * x2;
1292c733ed4SBarry Smith     v -= bs2;
1302c733ed4SBarry Smith 
1312c733ed4SBarry Smith     vi = aj + diag[i] - 1;
1322c733ed4SBarry Smith     nz = diag[i] - diag[i + 1] - 1;
1332c733ed4SBarry Smith     for (j = 0; j > -nz; j--) {
1342c733ed4SBarry Smith       oidx = bs * vi[j];
1352c733ed4SBarry Smith       t[oidx] -= v[0] * s1 + v[1] * s2;
1362c733ed4SBarry Smith       t[oidx + 1] -= v[2] * s1 + v[3] * s2;
1372c733ed4SBarry Smith       v -= bs2;
1382c733ed4SBarry Smith     }
1399371c9d4SSatish Balay     t[idx]     = s1;
1409371c9d4SSatish Balay     t[1 + idx] = s2;
1412c733ed4SBarry Smith     idx += bs;
1422c733ed4SBarry Smith   }
1432c733ed4SBarry Smith   /* backward solve the L^T */
1442c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1452c733ed4SBarry Smith     v   = aa + bs2 * ai[i];
1462c733ed4SBarry Smith     vi  = aj + ai[i];
1472c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1482c733ed4SBarry Smith     idt = bs * i;
1499371c9d4SSatish Balay     s1  = t[idt];
1509371c9d4SSatish Balay     s2  = t[1 + idt];
1512c733ed4SBarry Smith     for (j = 0; j < nz; j++) {
1522c733ed4SBarry Smith       idx = bs * vi[j];
1532c733ed4SBarry Smith       t[idx] -= v[0] * s1 + v[1] * s2;
1542c733ed4SBarry Smith       t[idx + 1] -= v[2] * s1 + v[3] * s2;
1552c733ed4SBarry Smith       v += bs2;
1562c733ed4SBarry Smith     }
1572c733ed4SBarry Smith   }
1582c733ed4SBarry Smith 
1592c733ed4SBarry Smith   /* copy t into x according to permutation */
1602c733ed4SBarry Smith   for (i = 0; i < n; i++) {
1619371c9d4SSatish Balay     ii        = bs * i;
1629371c9d4SSatish Balay     ir        = bs * r[i];
1639371c9d4SSatish Balay     x[ir]     = t[ii];
1649371c9d4SSatish Balay     x[ir + 1] = t[ii + 1];
1652c733ed4SBarry Smith   }
1662c733ed4SBarry Smith 
1679566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
1689566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
1699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1719566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
172*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1732c733ed4SBarry Smith }
174