xref: /petsc/src/mat/impls/baij/seq/baijsolvnat2.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 /*
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 */
8d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_2_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, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
122c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
132c733ed4SBarry Smith   PetscScalar       *x, s1, s2, x1, x2;
142c733ed4SBarry Smith   const PetscScalar *b;
152c733ed4SBarry Smith   PetscInt           jdx, idt, idx, nz, i;
162c733ed4SBarry Smith 
172c733ed4SBarry Smith   PetscFunctionBegin;
189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
199566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
202c733ed4SBarry Smith 
212c733ed4SBarry Smith   /* forward solve the lower triangular */
222c733ed4SBarry Smith   idx  = 0;
239371c9d4SSatish Balay   x[0] = b[0];
249371c9d4SSatish Balay   x[1] = b[1];
252c733ed4SBarry Smith   for (i = 1; i < n; i++) {
262c733ed4SBarry Smith     v  = aa + 4 * ai[i];
272c733ed4SBarry Smith     vi = aj + ai[i];
282c733ed4SBarry Smith     nz = diag[i] - ai[i];
292c733ed4SBarry Smith     idx += 2;
309371c9d4SSatish Balay     s1 = b[idx];
319371c9d4SSatish Balay     s2 = b[1 + idx];
322c733ed4SBarry Smith     while (nz--) {
332c733ed4SBarry Smith       jdx = 2 * (*vi++);
349371c9d4SSatish Balay       x1  = x[jdx];
359371c9d4SSatish Balay       x2  = x[1 + jdx];
362c733ed4SBarry Smith       s1 -= v[0] * x1 + v[2] * x2;
372c733ed4SBarry Smith       s2 -= v[1] * x1 + v[3] * x2;
382c733ed4SBarry Smith       v += 4;
392c733ed4SBarry Smith     }
402c733ed4SBarry Smith     x[idx]     = s1;
412c733ed4SBarry Smith     x[1 + idx] = s2;
422c733ed4SBarry Smith   }
432c733ed4SBarry Smith   /* backward solve the upper triangular */
442c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
452c733ed4SBarry Smith     v   = aa + 4 * diag[i] + 4;
462c733ed4SBarry Smith     vi  = aj + diag[i] + 1;
472c733ed4SBarry Smith     nz  = ai[i + 1] - diag[i] - 1;
482c733ed4SBarry Smith     idt = 2 * i;
499371c9d4SSatish Balay     s1  = x[idt];
509371c9d4SSatish Balay     s2  = x[1 + idt];
512c733ed4SBarry Smith     while (nz--) {
522c733ed4SBarry Smith       idx = 2 * (*vi++);
539371c9d4SSatish Balay       x1  = x[idx];
549371c9d4SSatish Balay       x2  = x[1 + idx];
552c733ed4SBarry Smith       s1 -= v[0] * x1 + v[2] * x2;
562c733ed4SBarry Smith       s2 -= v[1] * x1 + v[3] * x2;
572c733ed4SBarry Smith       v += 4;
582c733ed4SBarry Smith     }
592c733ed4SBarry Smith     v          = aa + 4 * diag[i];
602c733ed4SBarry Smith     x[idt]     = v[0] * s1 + v[2] * s2;
612c733ed4SBarry Smith     x[1 + idt] = v[1] * s1 + v[3] * s2;
622c733ed4SBarry Smith   }
632c733ed4SBarry Smith 
649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
669566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 4 * (a->nz) - 2.0 * A->cmap->n));
67*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
682c733ed4SBarry Smith }
692c733ed4SBarry Smith 
70d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_2_NaturalOrdering(Mat A, Vec bb, Vec xx)
71d71ae5a4SJacob Faibussowitsch {
722c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
732c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *ai = a->i, *aj = a->j, *adiag = a->diag;
742c733ed4SBarry Smith   PetscInt           i, k, nz, idx, idt, jdx;
752c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
762c733ed4SBarry Smith   PetscScalar       *x, s1, s2, x1, x2;
772c733ed4SBarry Smith   const PetscScalar *b;
782c733ed4SBarry Smith 
792c733ed4SBarry Smith   PetscFunctionBegin;
809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
819566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
822c733ed4SBarry Smith   /* forward solve the lower triangular */
832c733ed4SBarry Smith   idx  = 0;
849371c9d4SSatish Balay   x[0] = b[idx];
859371c9d4SSatish Balay   x[1] = b[1 + idx];
862c733ed4SBarry Smith   for (i = 1; i < n; i++) {
872c733ed4SBarry Smith     v   = aa + 4 * ai[i];
882c733ed4SBarry Smith     vi  = aj + ai[i];
892c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
902c733ed4SBarry Smith     idx = 2 * i;
919371c9d4SSatish Balay     s1  = b[idx];
929371c9d4SSatish Balay     s2  = b[1 + idx];
932c733ed4SBarry Smith     PetscPrefetchBlock(vi + nz, nz, 0, PETSC_PREFETCH_HINT_NTA);
942c733ed4SBarry Smith     PetscPrefetchBlock(v + 4 * nz, 4 * nz, 0, PETSC_PREFETCH_HINT_NTA);
952c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
962c733ed4SBarry Smith       jdx = 2 * vi[k];
979371c9d4SSatish Balay       x1  = x[jdx];
989371c9d4SSatish Balay       x2  = x[1 + jdx];
992c733ed4SBarry Smith       s1 -= v[0] * x1 + v[2] * x2;
1002c733ed4SBarry Smith       s2 -= v[1] * x1 + v[3] * x2;
1012c733ed4SBarry Smith       v += 4;
1022c733ed4SBarry Smith     }
1032c733ed4SBarry Smith     x[idx]     = s1;
1042c733ed4SBarry Smith     x[1 + idx] = s2;
1052c733ed4SBarry Smith   }
1062c733ed4SBarry Smith 
1072c733ed4SBarry Smith   /* backward solve the upper triangular */
1082c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1092c733ed4SBarry Smith     v   = aa + 4 * (adiag[i + 1] + 1);
1102c733ed4SBarry Smith     vi  = aj + adiag[i + 1] + 1;
1112c733ed4SBarry Smith     nz  = adiag[i] - adiag[i + 1] - 1;
1122c733ed4SBarry Smith     idt = 2 * i;
1139371c9d4SSatish Balay     s1  = x[idt];
1149371c9d4SSatish Balay     s2  = x[1 + idt];
1152c733ed4SBarry Smith     PetscPrefetchBlock(vi + nz, nz, 0, PETSC_PREFETCH_HINT_NTA);
1162c733ed4SBarry Smith     PetscPrefetchBlock(v + 4 * nz, 4 * nz, 0, PETSC_PREFETCH_HINT_NTA);
1172c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
1182c733ed4SBarry Smith       idx = 2 * vi[k];
1199371c9d4SSatish Balay       x1  = x[idx];
1209371c9d4SSatish Balay       x2  = x[1 + idx];
1212c733ed4SBarry Smith       s1 -= v[0] * x1 + v[2] * x2;
1222c733ed4SBarry Smith       s2 -= v[1] * x1 + v[3] * x2;
1232c733ed4SBarry Smith       v += 4;
1242c733ed4SBarry Smith     }
1252c733ed4SBarry Smith     /* x = inv_diagonal*x */
1262c733ed4SBarry Smith     x[idt]     = v[0] * s1 + v[2] * s2;
1272c733ed4SBarry Smith     x[1 + idt] = v[1] * s1 + v[3] * s2;
1282c733ed4SBarry Smith   }
1292c733ed4SBarry Smith 
1309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1329566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 4 * (a->nz) - 2.0 * A->cmap->n));
133*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1342c733ed4SBarry Smith }
1352c733ed4SBarry Smith 
136d71ae5a4SJacob Faibussowitsch PetscErrorCode MatForwardSolve_SeqBAIJ_2_NaturalOrdering(Mat A, Vec bb, Vec xx)
137d71ae5a4SJacob Faibussowitsch {
1382c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
1392c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *ai = a->i, *aj = a->j;
1402c733ed4SBarry Smith   PetscInt           i, k, nz, idx, jdx;
1412c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
1422c733ed4SBarry Smith   PetscScalar       *x, s1, s2, x1, x2;
1432c733ed4SBarry Smith   const PetscScalar *b;
1442c733ed4SBarry Smith 
1452c733ed4SBarry Smith   PetscFunctionBegin;
1469566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1482c733ed4SBarry Smith   /* forward solve the lower triangular */
1492c733ed4SBarry Smith   idx  = 0;
1509371c9d4SSatish Balay   x[0] = b[idx];
1519371c9d4SSatish Balay   x[1] = b[1 + idx];
1522c733ed4SBarry Smith   for (i = 1; i < n; i++) {
1532c733ed4SBarry Smith     v   = aa + 4 * ai[i];
1542c733ed4SBarry Smith     vi  = aj + ai[i];
1552c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1562c733ed4SBarry Smith     idx = 2 * i;
1579371c9d4SSatish Balay     s1  = b[idx];
1589371c9d4SSatish Balay     s2  = b[1 + idx];
1592c733ed4SBarry Smith     PetscPrefetchBlock(vi + nz, nz, 0, PETSC_PREFETCH_HINT_NTA);
1602c733ed4SBarry Smith     PetscPrefetchBlock(v + 4 * nz, 4 * nz, 0, PETSC_PREFETCH_HINT_NTA);
1612c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
1622c733ed4SBarry Smith       jdx = 2 * vi[k];
1639371c9d4SSatish Balay       x1  = x[jdx];
1649371c9d4SSatish Balay       x2  = x[1 + jdx];
1652c733ed4SBarry Smith       s1 -= v[0] * x1 + v[2] * x2;
1662c733ed4SBarry Smith       s2 -= v[1] * x1 + v[3] * x2;
1672c733ed4SBarry Smith       v += 4;
1682c733ed4SBarry Smith     }
1692c733ed4SBarry Smith     x[idx]     = s1;
1702c733ed4SBarry Smith     x[1 + idx] = s2;
1712c733ed4SBarry Smith   }
1722c733ed4SBarry Smith 
1739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1759566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(4.0 * (a->nz) - A->cmap->n));
176*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1772c733ed4SBarry Smith }
1782c733ed4SBarry Smith 
179d71ae5a4SJacob Faibussowitsch PetscErrorCode MatBackwardSolve_SeqBAIJ_2_NaturalOrdering(Mat A, Vec bb, Vec xx)
180d71ae5a4SJacob Faibussowitsch {
1812c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
1822c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *aj = a->j, *adiag = a->diag;
1832c733ed4SBarry Smith   PetscInt           i, k, nz, idx, idt;
1842c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
1852c733ed4SBarry Smith   PetscScalar       *x, s1, s2, x1, x2;
1862c733ed4SBarry Smith   const PetscScalar *b;
1872c733ed4SBarry Smith 
1882c733ed4SBarry Smith   PetscFunctionBegin;
1899566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1912c733ed4SBarry Smith 
1922c733ed4SBarry Smith   /* backward solve the upper triangular */
1932c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1942c733ed4SBarry Smith     v   = aa + 4 * (adiag[i + 1] + 1);
1952c733ed4SBarry Smith     vi  = aj + adiag[i + 1] + 1;
1962c733ed4SBarry Smith     nz  = adiag[i] - adiag[i + 1] - 1;
1972c733ed4SBarry Smith     idt = 2 * i;
1989371c9d4SSatish Balay     s1  = b[idt];
1999371c9d4SSatish Balay     s2  = b[1 + idt];
2002c733ed4SBarry Smith     PetscPrefetchBlock(vi + nz, nz, 0, PETSC_PREFETCH_HINT_NTA);
2012c733ed4SBarry Smith     PetscPrefetchBlock(v + 4 * nz, 4 * nz, 0, PETSC_PREFETCH_HINT_NTA);
2022c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
2032c733ed4SBarry Smith       idx = 2 * vi[k];
2049371c9d4SSatish Balay       x1  = x[idx];
2059371c9d4SSatish Balay       x2  = x[1 + idx];
2062c733ed4SBarry Smith       s1 -= v[0] * x1 + v[2] * x2;
2072c733ed4SBarry Smith       s2 -= v[1] * x1 + v[3] * x2;
2082c733ed4SBarry Smith       v += 4;
2092c733ed4SBarry Smith     }
2102c733ed4SBarry Smith     /* x = inv_diagonal*x */
2112c733ed4SBarry Smith     x[idt]     = v[0] * s1 + v[2] * s2;
2122c733ed4SBarry Smith     x[1 + idt] = v[1] * s1 + v[3] * s2;
2132c733ed4SBarry Smith   }
2142c733ed4SBarry Smith 
2159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
2169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
2179566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(4.0 * a->nz - A->cmap->n));
218*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2192c733ed4SBarry Smith }
220