xref: /petsc/src/mat/impls/baij/seq/baijsolvnat1.c (revision 31d78bcd2b98084dc1368b20eb1129c8b9fb39fe)
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_1_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)8d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_1_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;
142c733ed4SBarry Smith   const PetscScalar *b;
152c733ed4SBarry Smith   PetscScalar        s1, x1;
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;
242c733ed4SBarry Smith   x[0] = b[0];
252c733ed4SBarry Smith   for (i = 1; i < n; i++) {
262c733ed4SBarry Smith     v  = aa + ai[i];
272c733ed4SBarry Smith     vi = aj + ai[i];
282c733ed4SBarry Smith     nz = diag[i] - ai[i];
292c733ed4SBarry Smith     idx += 1;
302c733ed4SBarry Smith     s1 = b[idx];
312c733ed4SBarry Smith     while (nz--) {
322c733ed4SBarry Smith       jdx = *vi++;
332c733ed4SBarry Smith       x1  = x[jdx];
342c733ed4SBarry Smith       s1 -= v[0] * x1;
352c733ed4SBarry Smith       v += 1;
362c733ed4SBarry Smith     }
372c733ed4SBarry Smith     x[idx] = s1;
382c733ed4SBarry Smith   }
392c733ed4SBarry Smith   /* backward solve the upper triangular */
402c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
412c733ed4SBarry Smith     v   = aa + diag[i] + 1;
422c733ed4SBarry Smith     vi  = aj + diag[i] + 1;
432c733ed4SBarry Smith     nz  = ai[i + 1] - diag[i] - 1;
442c733ed4SBarry Smith     idt = i;
452c733ed4SBarry Smith     s1  = x[idt];
462c733ed4SBarry Smith     while (nz--) {
472c733ed4SBarry Smith       idx = *vi++;
482c733ed4SBarry Smith       x1  = x[idx];
492c733ed4SBarry Smith       s1 -= v[0] * x1;
502c733ed4SBarry Smith       v += 1;
512c733ed4SBarry Smith     }
522c733ed4SBarry Smith     v      = aa + diag[i];
532c733ed4SBarry Smith     x[idt] = v[0] * s1;
542c733ed4SBarry Smith   }
559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
579566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * (a->nz) - A->cmap->n));
58*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
592c733ed4SBarry Smith }
602c733ed4SBarry Smith 
MatForwardSolve_SeqBAIJ_1_NaturalOrdering(Mat A,Vec bb,Vec xx)61d71ae5a4SJacob Faibussowitsch PetscErrorCode MatForwardSolve_SeqBAIJ_1_NaturalOrdering(Mat A, Vec bb, Vec xx)
62d71ae5a4SJacob Faibussowitsch {
632c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
642c733ed4SBarry Smith   const PetscInt     n = a->mbs, *ai = a->i, *aj = a->j, *vi;
652c733ed4SBarry Smith   PetscScalar       *x, sum;
662c733ed4SBarry Smith   const PetscScalar *b;
672c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
682c733ed4SBarry Smith   PetscInt           i, nz;
692c733ed4SBarry Smith 
702c733ed4SBarry Smith   PetscFunctionBegin;
71*3ba16761SJacob Faibussowitsch   if (!n) PetscFunctionReturn(PETSC_SUCCESS);
722c733ed4SBarry Smith 
739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
752c733ed4SBarry Smith 
762c733ed4SBarry Smith   /* forward solve the lower triangular */
772c733ed4SBarry Smith   x[0] = b[0];
782c733ed4SBarry Smith   v    = aa;
792c733ed4SBarry Smith   vi   = aj;
802c733ed4SBarry Smith   for (i = 1; i < n; i++) {
812c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
822c733ed4SBarry Smith     sum = b[i];
832c733ed4SBarry Smith     PetscSparseDenseMinusDot(sum, x, v, vi, nz);
842c733ed4SBarry Smith     v += nz;
852c733ed4SBarry Smith     vi += nz;
862c733ed4SBarry Smith     x[i] = sum;
872c733ed4SBarry Smith   }
889566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(a->nz - A->cmap->n));
899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
91*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
922c733ed4SBarry Smith }
932c733ed4SBarry Smith 
MatBackwardSolve_SeqBAIJ_1_NaturalOrdering(Mat A,Vec bb,Vec xx)94d71ae5a4SJacob Faibussowitsch PetscErrorCode MatBackwardSolve_SeqBAIJ_1_NaturalOrdering(Mat A, Vec bb, Vec xx)
95d71ae5a4SJacob Faibussowitsch {
962c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
972c733ed4SBarry Smith   const PetscInt     n = a->mbs, *aj = a->j, *adiag = a->diag, *vi;
982c733ed4SBarry Smith   PetscScalar       *x, sum;
992c733ed4SBarry Smith   const PetscScalar *b;
1002c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
1012c733ed4SBarry Smith   PetscInt           i, nz;
1022c733ed4SBarry Smith 
1032c733ed4SBarry Smith   PetscFunctionBegin;
104*3ba16761SJacob Faibussowitsch   if (!n) PetscFunctionReturn(PETSC_SUCCESS);
1052c733ed4SBarry Smith 
1069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1082c733ed4SBarry Smith 
1092c733ed4SBarry Smith   /* backward solve the upper triangular */
1102c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1112c733ed4SBarry Smith     v   = aa + adiag[i + 1] + 1;
1122c733ed4SBarry Smith     vi  = aj + adiag[i + 1] + 1;
1132c733ed4SBarry Smith     nz  = adiag[i] - adiag[i + 1] - 1;
1142c733ed4SBarry Smith     sum = b[i];
1152c733ed4SBarry Smith     PetscSparseDenseMinusDot(sum, x, v, vi, nz);
1162c733ed4SBarry Smith     x[i] = sum * v[nz]; /* x[i]=aa[adiag[i]]*sum; v++; */
1172c733ed4SBarry Smith   }
1182c733ed4SBarry Smith 
1199566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * a->nz - A->cmap->n));
1209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
122*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1232c733ed4SBarry Smith }
1242c733ed4SBarry Smith 
MatSolve_SeqBAIJ_1_NaturalOrdering(Mat A,Vec bb,Vec xx)125d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_1_NaturalOrdering(Mat A, Vec bb, Vec xx)
126d71ae5a4SJacob Faibussowitsch {
1272c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
1282c733ed4SBarry Smith   const PetscInt     n = a->mbs, *ai = a->i, *aj = a->j, *adiag = a->diag, *vi;
1292c733ed4SBarry Smith   PetscScalar       *x, sum;
1302c733ed4SBarry Smith   const PetscScalar *b;
1312c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
1322c733ed4SBarry Smith   PetscInt           i, nz;
1332c733ed4SBarry Smith 
1342c733ed4SBarry Smith   PetscFunctionBegin;
135*3ba16761SJacob Faibussowitsch   if (!n) PetscFunctionReturn(PETSC_SUCCESS);
1362c733ed4SBarry Smith 
1379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1392c733ed4SBarry Smith 
1402c733ed4SBarry Smith   /* forward solve the lower triangular */
1412c733ed4SBarry Smith   x[0] = b[0];
1422c733ed4SBarry Smith   v    = aa;
1432c733ed4SBarry Smith   vi   = aj;
1442c733ed4SBarry Smith   for (i = 1; i < n; i++) {
1452c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1462c733ed4SBarry Smith     sum = b[i];
1472c733ed4SBarry Smith     PetscSparseDenseMinusDot(sum, x, v, vi, nz);
1482c733ed4SBarry Smith     v += nz;
1492c733ed4SBarry Smith     vi += nz;
1502c733ed4SBarry Smith     x[i] = sum;
1512c733ed4SBarry Smith   }
1522c733ed4SBarry Smith 
1532c733ed4SBarry Smith   /* backward solve the upper triangular */
1542c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1552c733ed4SBarry Smith     v   = aa + adiag[i + 1] + 1;
1562c733ed4SBarry Smith     vi  = aj + adiag[i + 1] + 1;
1572c733ed4SBarry Smith     nz  = adiag[i] - adiag[i + 1] - 1;
1582c733ed4SBarry Smith     sum = x[i];
1592c733ed4SBarry Smith     PetscSparseDenseMinusDot(sum, x, v, vi, nz);
1602c733ed4SBarry Smith     x[i] = sum * v[nz]; /* x[i]=aa[adiag[i]]*sum; v++; */
1612c733ed4SBarry Smith   }
1622c733ed4SBarry Smith 
1639566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * a->nz - A->cmap->n));
1649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
166*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1672c733ed4SBarry Smith }
168