xref: /petsc/src/mat/impls/baij/seq/baijsolvtran3.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_3_inplace(Mat A,Vec bb,Vec xx)4d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_3_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, s3, x1, x2, x3, *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        = 3 * c[i];
292c733ed4SBarry Smith     t[ii]     = b[ic];
302c733ed4SBarry Smith     t[ii + 1] = b[ic + 1];
312c733ed4SBarry Smith     t[ii + 2] = b[ic + 2];
322c733ed4SBarry Smith     ii += 3;
332c733ed4SBarry Smith   }
342c733ed4SBarry Smith 
352c733ed4SBarry Smith   /* forward solve the U^T */
362c733ed4SBarry Smith   idx = 0;
372c733ed4SBarry Smith   for (i = 0; i < n; i++) {
382c733ed4SBarry Smith     v = aa + 9 * diag[i];
392c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
409371c9d4SSatish Balay     x1 = t[idx];
419371c9d4SSatish Balay     x2 = t[1 + idx];
429371c9d4SSatish Balay     x3 = t[2 + idx];
432c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3;
442c733ed4SBarry Smith     s2 = v[3] * x1 + v[4] * x2 + v[5] * x3;
452c733ed4SBarry Smith     s3 = v[6] * x1 + v[7] * x2 + v[8] * x3;
462c733ed4SBarry Smith     v += 9;
472c733ed4SBarry Smith 
482c733ed4SBarry Smith     vi = aj + diag[i] + 1;
492c733ed4SBarry Smith     nz = ai[i + 1] - diag[i] - 1;
502c733ed4SBarry Smith     while (nz--) {
512c733ed4SBarry Smith       oidx = 3 * (*vi++);
522c733ed4SBarry Smith       t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3;
532c733ed4SBarry Smith       t[oidx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3;
542c733ed4SBarry Smith       t[oidx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3;
552c733ed4SBarry Smith       v += 9;
562c733ed4SBarry Smith     }
579371c9d4SSatish Balay     t[idx]     = s1;
589371c9d4SSatish Balay     t[1 + idx] = s2;
599371c9d4SSatish Balay     t[2 + idx] = s3;
602c733ed4SBarry Smith     idx += 3;
612c733ed4SBarry Smith   }
622c733ed4SBarry Smith   /* backward solve the L^T */
632c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
642c733ed4SBarry Smith     v   = aa + 9 * diag[i] - 9;
652c733ed4SBarry Smith     vi  = aj + diag[i] - 1;
662c733ed4SBarry Smith     nz  = diag[i] - ai[i];
672c733ed4SBarry Smith     idt = 3 * i;
689371c9d4SSatish Balay     s1  = t[idt];
699371c9d4SSatish Balay     s2  = t[1 + idt];
709371c9d4SSatish Balay     s3  = t[2 + idt];
712c733ed4SBarry Smith     while (nz--) {
722c733ed4SBarry Smith       idx = 3 * (*vi--);
732c733ed4SBarry Smith       t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3;
742c733ed4SBarry Smith       t[idx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3;
752c733ed4SBarry Smith       t[idx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3;
762c733ed4SBarry Smith       v -= 9;
772c733ed4SBarry Smith     }
782c733ed4SBarry Smith   }
792c733ed4SBarry Smith 
802c733ed4SBarry Smith   /* copy t into x according to permutation */
812c733ed4SBarry Smith   ii = 0;
822c733ed4SBarry Smith   for (i = 0; i < n; i++) {
832c733ed4SBarry Smith     ir        = 3 * r[i];
842c733ed4SBarry Smith     x[ir]     = t[ii];
852c733ed4SBarry Smith     x[ir + 1] = t[ii + 1];
862c733ed4SBarry Smith     x[ir + 2] = t[ii + 2];
872c733ed4SBarry Smith     ii += 3;
882c733ed4SBarry Smith   }
892c733ed4SBarry Smith 
909566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
919566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
949566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 9 * (a->nz) - 3.0 * A->cmap->n));
95*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
962c733ed4SBarry Smith }
972c733ed4SBarry Smith 
MatSolveTranspose_SeqBAIJ_3(Mat A,Vec bb,Vec xx)98d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_3(Mat A, Vec bb, Vec xx)
99d71ae5a4SJacob Faibussowitsch {
1002c733ed4SBarry Smith   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
1012c733ed4SBarry Smith   IS                 iscol = a->col, isrow = a->row;
1022c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
1032c733ed4SBarry Smith   const PetscInt    *r, *c, *rout, *cout;
1042c733ed4SBarry Smith   PetscInt           nz, idx, idt, j, i, oidx, ii, ic, ir;
1052c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
1062c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
1072c733ed4SBarry Smith   PetscScalar        s1, s2, s3, x1, x2, x3, *x, *t;
1082c733ed4SBarry Smith   const PetscScalar *b;
1092c733ed4SBarry Smith 
1102c733ed4SBarry Smith   PetscFunctionBegin;
1119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1129566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1132c733ed4SBarry Smith   t = a->solve_work;
1142c733ed4SBarry Smith 
1159371c9d4SSatish Balay   PetscCall(ISGetIndices(isrow, &rout));
1169371c9d4SSatish Balay   r = rout;
1179371c9d4SSatish Balay   PetscCall(ISGetIndices(iscol, &cout));
1189371c9d4SSatish Balay   c = cout;
1192c733ed4SBarry Smith 
1202c733ed4SBarry Smith   /* copy b into temp work space according to permutation */
1212c733ed4SBarry Smith   for (i = 0; i < n; i++) {
1229371c9d4SSatish Balay     ii        = bs * i;
1239371c9d4SSatish Balay     ic        = bs * c[i];
1249371c9d4SSatish Balay     t[ii]     = b[ic];
1259371c9d4SSatish Balay     t[ii + 1] = b[ic + 1];
1269371c9d4SSatish Balay     t[ii + 2] = b[ic + 2];
1272c733ed4SBarry Smith   }
1282c733ed4SBarry Smith 
1292c733ed4SBarry Smith   /* forward solve the U^T */
1302c733ed4SBarry Smith   idx = 0;
1312c733ed4SBarry Smith   for (i = 0; i < n; i++) {
1322c733ed4SBarry Smith     v = aa + bs2 * diag[i];
1332c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
1349371c9d4SSatish Balay     x1 = t[idx];
1359371c9d4SSatish Balay     x2 = t[1 + idx];
1369371c9d4SSatish Balay     x3 = t[2 + idx];
1372c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3;
1382c733ed4SBarry Smith     s2 = v[3] * x1 + v[4] * x2 + v[5] * x3;
1392c733ed4SBarry Smith     s3 = v[6] * x1 + v[7] * x2 + v[8] * x3;
1402c733ed4SBarry Smith     v -= bs2;
1412c733ed4SBarry Smith 
1422c733ed4SBarry Smith     vi = aj + diag[i] - 1;
1432c733ed4SBarry Smith     nz = diag[i] - diag[i + 1] - 1;
1442c733ed4SBarry Smith     for (j = 0; j > -nz; j--) {
1452c733ed4SBarry Smith       oidx = bs * vi[j];
1462c733ed4SBarry Smith       t[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3;
1472c733ed4SBarry Smith       t[oidx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3;
1482c733ed4SBarry Smith       t[oidx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3;
1492c733ed4SBarry Smith       v -= bs2;
1502c733ed4SBarry Smith     }
1519371c9d4SSatish Balay     t[idx]     = s1;
1529371c9d4SSatish Balay     t[1 + idx] = s2;
1539371c9d4SSatish Balay     t[2 + idx] = s3;
1542c733ed4SBarry Smith     idx += bs;
1552c733ed4SBarry Smith   }
1562c733ed4SBarry Smith   /* backward solve the L^T */
1572c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1582c733ed4SBarry Smith     v   = aa + bs2 * ai[i];
1592c733ed4SBarry Smith     vi  = aj + ai[i];
1602c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1612c733ed4SBarry Smith     idt = bs * i;
1629371c9d4SSatish Balay     s1  = t[idt];
1639371c9d4SSatish Balay     s2  = t[1 + idt];
1649371c9d4SSatish Balay     s3  = t[2 + idt];
1652c733ed4SBarry Smith     for (j = 0; j < nz; j++) {
1662c733ed4SBarry Smith       idx = bs * vi[j];
1672c733ed4SBarry Smith       t[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3;
1682c733ed4SBarry Smith       t[idx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3;
1692c733ed4SBarry Smith       t[idx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3;
1702c733ed4SBarry Smith       v += bs2;
1712c733ed4SBarry Smith     }
1722c733ed4SBarry Smith   }
1732c733ed4SBarry Smith 
1742c733ed4SBarry Smith   /* copy t into x according to permutation */
1752c733ed4SBarry Smith   for (i = 0; i < n; i++) {
1769371c9d4SSatish Balay     ii        = bs * i;
1779371c9d4SSatish Balay     ir        = bs * r[i];
1789371c9d4SSatish Balay     x[ir]     = t[ii];
1799371c9d4SSatish Balay     x[ir + 1] = t[ii + 1];
1809371c9d4SSatish Balay     x[ir + 2] = t[ii + 2];
1812c733ed4SBarry Smith   }
1822c733ed4SBarry Smith 
1839566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &rout));
1849566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol, &cout));
1859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1879566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
188*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1892c733ed4SBarry Smith }
190