12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h>
32c733ed4SBarry Smith
42c733ed4SBarry Smith /*
52c733ed4SBarry Smith Special case where the matrix was ILU(0) factored in the natural
62c733ed4SBarry Smith ordering. This eliminates the need for the column and row permutation.
72c733ed4SBarry Smith */
MatSolve_SeqBAIJ_3_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)8d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_3_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx)
9d71ae5a4SJacob Faibussowitsch {
102c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
112c733ed4SBarry Smith const PetscInt n = a->mbs, *ai = a->i, *aj = a->j;
122c733ed4SBarry Smith const PetscInt *diag = a->diag, *vi;
132c733ed4SBarry Smith const MatScalar *aa = a->a, *v;
142c733ed4SBarry Smith PetscScalar *x, s1, s2, s3, x1, x2, x3;
152c733ed4SBarry Smith const PetscScalar *b;
162c733ed4SBarry Smith PetscInt jdx, idt, idx, nz, i;
172c733ed4SBarry Smith
182c733ed4SBarry Smith PetscFunctionBegin;
199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b));
209566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
212c733ed4SBarry Smith
222c733ed4SBarry Smith /* forward solve the lower triangular */
232c733ed4SBarry Smith idx = 0;
249371c9d4SSatish Balay x[0] = b[0];
259371c9d4SSatish Balay x[1] = b[1];
269371c9d4SSatish Balay x[2] = b[2];
272c733ed4SBarry Smith for (i = 1; i < n; i++) {
282c733ed4SBarry Smith v = aa + 9 * ai[i];
292c733ed4SBarry Smith vi = aj + ai[i];
302c733ed4SBarry Smith nz = diag[i] - ai[i];
312c733ed4SBarry Smith idx += 3;
329371c9d4SSatish Balay s1 = b[idx];
339371c9d4SSatish Balay s2 = b[1 + idx];
349371c9d4SSatish Balay s3 = b[2 + idx];
352c733ed4SBarry Smith while (nz--) {
362c733ed4SBarry Smith jdx = 3 * (*vi++);
379371c9d4SSatish Balay x1 = x[jdx];
389371c9d4SSatish Balay x2 = x[1 + jdx];
399371c9d4SSatish Balay x3 = x[2 + jdx];
402c733ed4SBarry Smith s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
412c733ed4SBarry Smith s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
422c733ed4SBarry Smith s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
432c733ed4SBarry Smith v += 9;
442c733ed4SBarry Smith }
452c733ed4SBarry Smith x[idx] = s1;
462c733ed4SBarry Smith x[1 + idx] = s2;
472c733ed4SBarry Smith x[2 + idx] = s3;
482c733ed4SBarry Smith }
492c733ed4SBarry Smith /* backward solve the upper triangular */
502c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) {
512c733ed4SBarry Smith v = aa + 9 * diag[i] + 9;
522c733ed4SBarry Smith vi = aj + diag[i] + 1;
532c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1;
542c733ed4SBarry Smith idt = 3 * i;
559371c9d4SSatish Balay s1 = x[idt];
569371c9d4SSatish Balay s2 = x[1 + idt];
572c733ed4SBarry Smith s3 = x[2 + idt];
582c733ed4SBarry Smith while (nz--) {
592c733ed4SBarry Smith idx = 3 * (*vi++);
609371c9d4SSatish Balay x1 = x[idx];
619371c9d4SSatish Balay x2 = x[1 + idx];
629371c9d4SSatish Balay x3 = x[2 + idx];
632c733ed4SBarry Smith s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
642c733ed4SBarry Smith s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
652c733ed4SBarry Smith s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
662c733ed4SBarry Smith v += 9;
672c733ed4SBarry Smith }
682c733ed4SBarry Smith v = aa + 9 * diag[i];
692c733ed4SBarry Smith x[idt] = v[0] * s1 + v[3] * s2 + v[6] * s3;
702c733ed4SBarry Smith x[1 + idt] = v[1] * s1 + v[4] * s2 + v[7] * s3;
712c733ed4SBarry Smith x[2 + idt] = v[2] * s1 + v[5] * s2 + v[8] * s3;
722c733ed4SBarry Smith }
732c733ed4SBarry Smith
749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b));
759566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
769566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 9 * (a->nz) - 3.0 * A->cmap->n));
77*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
782c733ed4SBarry Smith }
792c733ed4SBarry Smith
MatSolve_SeqBAIJ_3_NaturalOrdering(Mat A,Vec bb,Vec xx)80d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_3_NaturalOrdering(Mat A, Vec bb, Vec xx)
81d71ae5a4SJacob Faibussowitsch {
822c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
832c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *adiag = a->diag;
842c733ed4SBarry Smith PetscInt i, k, nz, idx, jdx, idt;
852c733ed4SBarry Smith const PetscInt bs = A->rmap->bs, bs2 = a->bs2;
862c733ed4SBarry Smith const MatScalar *aa = a->a, *v;
872c733ed4SBarry Smith PetscScalar *x;
882c733ed4SBarry Smith const PetscScalar *b;
892c733ed4SBarry Smith PetscScalar s1, s2, s3, x1, x2, x3;
902c733ed4SBarry Smith
912c733ed4SBarry Smith PetscFunctionBegin;
929566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b));
939566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
942c733ed4SBarry Smith /* forward solve the lower triangular */
952c733ed4SBarry Smith idx = 0;
969371c9d4SSatish Balay x[0] = b[idx];
979371c9d4SSatish Balay x[1] = b[1 + idx];
989371c9d4SSatish Balay x[2] = b[2 + idx];
992c733ed4SBarry Smith for (i = 1; i < n; i++) {
1002c733ed4SBarry Smith v = aa + bs2 * ai[i];
1012c733ed4SBarry Smith vi = aj + ai[i];
1022c733ed4SBarry Smith nz = ai[i + 1] - ai[i];
1032c733ed4SBarry Smith idx = bs * i;
1049371c9d4SSatish Balay s1 = b[idx];
1059371c9d4SSatish Balay s2 = b[1 + idx];
1069371c9d4SSatish Balay s3 = b[2 + idx];
1072c733ed4SBarry Smith for (k = 0; k < nz; k++) {
1082c733ed4SBarry Smith jdx = bs * vi[k];
1099371c9d4SSatish Balay x1 = x[jdx];
1109371c9d4SSatish Balay x2 = x[1 + jdx];
1119371c9d4SSatish Balay x3 = x[2 + jdx];
1122c733ed4SBarry Smith s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
1132c733ed4SBarry Smith s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
1142c733ed4SBarry Smith s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
1152c733ed4SBarry Smith
1162c733ed4SBarry Smith v += bs2;
1172c733ed4SBarry Smith }
1182c733ed4SBarry Smith
1192c733ed4SBarry Smith x[idx] = s1;
1202c733ed4SBarry Smith x[1 + idx] = s2;
1212c733ed4SBarry Smith x[2 + idx] = s3;
1222c733ed4SBarry Smith }
1232c733ed4SBarry Smith
1242c733ed4SBarry Smith /* backward solve the upper triangular */
1252c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) {
1262c733ed4SBarry Smith v = aa + bs2 * (adiag[i + 1] + 1);
1272c733ed4SBarry Smith vi = aj + adiag[i + 1] + 1;
1282c733ed4SBarry Smith nz = adiag[i] - adiag[i + 1] - 1;
1292c733ed4SBarry Smith idt = bs * i;
1309371c9d4SSatish Balay s1 = x[idt];
1319371c9d4SSatish Balay s2 = x[1 + idt];
1329371c9d4SSatish Balay s3 = x[2 + idt];
1332c733ed4SBarry Smith
1342c733ed4SBarry Smith for (k = 0; k < nz; k++) {
1352c733ed4SBarry Smith idx = bs * vi[k];
1369371c9d4SSatish Balay x1 = x[idx];
1379371c9d4SSatish Balay x2 = x[1 + idx];
1389371c9d4SSatish Balay x3 = x[2 + idx];
1392c733ed4SBarry Smith s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
1402c733ed4SBarry Smith s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
1412c733ed4SBarry Smith s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
1422c733ed4SBarry Smith
1432c733ed4SBarry Smith v += bs2;
1442c733ed4SBarry Smith }
1452c733ed4SBarry Smith /* x = inv_diagonal*x */
1462c733ed4SBarry Smith x[idt] = v[0] * s1 + v[3] * s2 + v[6] * s3;
1472c733ed4SBarry Smith x[1 + idt] = v[1] * s1 + v[4] * s2 + v[7] * s3;
1482c733ed4SBarry Smith x[2 + idt] = v[2] * s1 + v[5] * s2 + v[8] * s3;
1492c733ed4SBarry Smith }
1502c733ed4SBarry Smith
1519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b));
1529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
1539566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
154*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1552c733ed4SBarry Smith }
1562c733ed4SBarry Smith
MatForwardSolve_SeqBAIJ_3_NaturalOrdering(Mat A,Vec bb,Vec xx)157d71ae5a4SJacob Faibussowitsch PetscErrorCode MatForwardSolve_SeqBAIJ_3_NaturalOrdering(Mat A, Vec bb, Vec xx)
158d71ae5a4SJacob Faibussowitsch {
1592c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
1602c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j;
1612c733ed4SBarry Smith PetscInt i, k, nz, idx, jdx;
1622c733ed4SBarry Smith const PetscInt bs = A->rmap->bs, bs2 = a->bs2;
1632c733ed4SBarry Smith const MatScalar *aa = a->a, *v;
1642c733ed4SBarry Smith PetscScalar *x;
1652c733ed4SBarry Smith const PetscScalar *b;
1662c733ed4SBarry Smith PetscScalar s1, s2, s3, x1, x2, x3;
1672c733ed4SBarry Smith
1682c733ed4SBarry Smith PetscFunctionBegin;
1699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b));
1709566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
1712c733ed4SBarry Smith /* forward solve the lower triangular */
1722c733ed4SBarry Smith idx = 0;
1739371c9d4SSatish Balay x[0] = b[idx];
1749371c9d4SSatish Balay x[1] = b[1 + idx];
1759371c9d4SSatish Balay x[2] = b[2 + idx];
1762c733ed4SBarry Smith for (i = 1; i < n; i++) {
1772c733ed4SBarry Smith v = aa + bs2 * ai[i];
1782c733ed4SBarry Smith vi = aj + ai[i];
1792c733ed4SBarry Smith nz = ai[i + 1] - ai[i];
1802c733ed4SBarry Smith idx = bs * i;
1819371c9d4SSatish Balay s1 = b[idx];
1829371c9d4SSatish Balay s2 = b[1 + idx];
1839371c9d4SSatish Balay s3 = b[2 + idx];
1842c733ed4SBarry Smith for (k = 0; k < nz; k++) {
1852c733ed4SBarry Smith jdx = bs * vi[k];
1869371c9d4SSatish Balay x1 = x[jdx];
1879371c9d4SSatish Balay x2 = x[1 + jdx];
1889371c9d4SSatish Balay x3 = x[2 + jdx];
1892c733ed4SBarry Smith s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
1902c733ed4SBarry Smith s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
1912c733ed4SBarry Smith s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
1922c733ed4SBarry Smith
1932c733ed4SBarry Smith v += bs2;
1942c733ed4SBarry Smith }
1952c733ed4SBarry Smith
1962c733ed4SBarry Smith x[idx] = s1;
1972c733ed4SBarry Smith x[1 + idx] = s2;
1982c733ed4SBarry Smith x[2 + idx] = s3;
1992c733ed4SBarry Smith }
2002c733ed4SBarry Smith
2019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b));
2029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
2039566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.0 * bs2 * (a->nz) - bs * A->cmap->n));
204*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2052c733ed4SBarry Smith }
2062c733ed4SBarry Smith
MatBackwardSolve_SeqBAIJ_3_NaturalOrdering(Mat A,Vec bb,Vec xx)207d71ae5a4SJacob Faibussowitsch PetscErrorCode MatBackwardSolve_SeqBAIJ_3_NaturalOrdering(Mat A, Vec bb, Vec xx)
208d71ae5a4SJacob Faibussowitsch {
2092c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
2102c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *aj = a->j, *adiag = a->diag;
2112c733ed4SBarry Smith PetscInt i, k, nz, idx, idt;
2122c733ed4SBarry Smith const PetscInt bs = A->rmap->bs, bs2 = a->bs2;
2132c733ed4SBarry Smith const MatScalar *aa = a->a, *v;
2142c733ed4SBarry Smith PetscScalar *x;
2152c733ed4SBarry Smith const PetscScalar *b;
2162c733ed4SBarry Smith PetscScalar s1, s2, s3, x1, x2, x3;
2172c733ed4SBarry Smith
2182c733ed4SBarry Smith PetscFunctionBegin;
2199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b));
2209566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
2212c733ed4SBarry Smith
2222c733ed4SBarry Smith /* backward solve the upper triangular */
2232c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) {
2242c733ed4SBarry Smith v = aa + bs2 * (adiag[i + 1] + 1);
2252c733ed4SBarry Smith vi = aj + adiag[i + 1] + 1;
2262c733ed4SBarry Smith nz = adiag[i] - adiag[i + 1] - 1;
2272c733ed4SBarry Smith idt = bs * i;
2289371c9d4SSatish Balay s1 = b[idt];
2299371c9d4SSatish Balay s2 = b[1 + idt];
2309371c9d4SSatish Balay s3 = b[2 + idt];
2312c733ed4SBarry Smith
2322c733ed4SBarry Smith for (k = 0; k < nz; k++) {
2332c733ed4SBarry Smith idx = bs * vi[k];
2349371c9d4SSatish Balay x1 = x[idx];
2359371c9d4SSatish Balay x2 = x[1 + idx];
2369371c9d4SSatish Balay x3 = x[2 + idx];
2372c733ed4SBarry Smith s1 -= v[0] * x1 + v[3] * x2 + v[6] * x3;
2382c733ed4SBarry Smith s2 -= v[1] * x1 + v[4] * x2 + v[7] * x3;
2392c733ed4SBarry Smith s3 -= v[2] * x1 + v[5] * x2 + v[8] * x3;
2402c733ed4SBarry Smith
2412c733ed4SBarry Smith v += bs2;
2422c733ed4SBarry Smith }
2432c733ed4SBarry Smith /* x = inv_diagonal*x */
2442c733ed4SBarry Smith x[idt] = v[0] * s1 + v[3] * s2 + v[6] * s3;
2452c733ed4SBarry Smith x[1 + idt] = v[1] * s1 + v[4] * s2 + v[7] * s3;
2462c733ed4SBarry Smith x[2 + idt] = v[2] * s1 + v[5] * s2 + v[8] * s3;
2472c733ed4SBarry Smith }
2482c733ed4SBarry Smith
2499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b));
2509566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
2519566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
252*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2532c733ed4SBarry Smith }
254