xref: /petsc/src/mat/impls/baij/seq/baijsolvnat4.c (revision a02648fdf9ec0d41d7b5ca02cb70ddcfa0e65728)
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_4_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)8d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx)
9d71ae5a4SJacob Faibussowitsch {
102c733ed4SBarry Smith   Mat_SeqBAIJ       *a  = (Mat_SeqBAIJ *)A->data;
112c733ed4SBarry Smith   PetscInt           n  = a->mbs;
122c733ed4SBarry Smith   const PetscInt    *ai = a->i, *aj = a->j;
132c733ed4SBarry Smith   const PetscInt    *diag = a->diag;
142c733ed4SBarry Smith   const MatScalar   *aa   = a->a;
152c733ed4SBarry Smith   PetscScalar       *x;
162c733ed4SBarry Smith   const PetscScalar *b;
172c733ed4SBarry Smith 
182c733ed4SBarry Smith   PetscFunctionBegin;
199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
209566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
212c733ed4SBarry Smith 
222c733ed4SBarry Smith #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJ)
232c733ed4SBarry Smith   {
242c733ed4SBarry Smith     static PetscScalar w[2000]; /* very BAD need to fix */
252c733ed4SBarry Smith     fortransolvebaij4_(&n, x, ai, aj, diag, aa, b, w);
262c733ed4SBarry Smith   }
272c733ed4SBarry Smith #elif defined(PETSC_USE_FORTRAN_KERNEL_SOLVEBAIJUNROLL)
282c733ed4SBarry Smith   fortransolvebaij4unroll_(&n, x, ai, aj, diag, aa, b);
292c733ed4SBarry Smith #else
302c733ed4SBarry Smith   {
312c733ed4SBarry Smith     PetscScalar      s1, s2, s3, s4, x1, x2, x3, x4;
322c733ed4SBarry Smith     const MatScalar *v;
332c733ed4SBarry Smith     PetscInt         jdx, idt, idx, nz, i, ai16;
342c733ed4SBarry Smith     const PetscInt  *vi;
352c733ed4SBarry Smith 
362c733ed4SBarry Smith     /* forward solve the lower triangular */
372c733ed4SBarry Smith     idx  = 0;
389371c9d4SSatish Balay     x[0] = b[0];
399371c9d4SSatish Balay     x[1] = b[1];
409371c9d4SSatish Balay     x[2] = b[2];
419371c9d4SSatish Balay     x[3] = b[3];
422c733ed4SBarry Smith     for (i = 1; i < n; i++) {
432c733ed4SBarry Smith       v  = aa + 16 * ai[i];
442c733ed4SBarry Smith       vi = aj + ai[i];
452c733ed4SBarry Smith       nz = diag[i] - ai[i];
462c733ed4SBarry Smith       idx += 4;
479371c9d4SSatish Balay       s1 = b[idx];
489371c9d4SSatish Balay       s2 = b[1 + idx];
499371c9d4SSatish Balay       s3 = b[2 + idx];
509371c9d4SSatish Balay       s4 = b[3 + idx];
512c733ed4SBarry Smith       while (nz--) {
522c733ed4SBarry Smith         jdx = 4 * (*vi++);
539371c9d4SSatish Balay         x1  = x[jdx];
549371c9d4SSatish Balay         x2  = x[1 + jdx];
559371c9d4SSatish Balay         x3  = x[2 + jdx];
569371c9d4SSatish Balay         x4  = x[3 + jdx];
572c733ed4SBarry Smith         s1 -= v[0] * x1 + v[4] * x2 + v[8] * x3 + v[12] * x4;
582c733ed4SBarry Smith         s2 -= v[1] * x1 + v[5] * x2 + v[9] * x3 + v[13] * x4;
592c733ed4SBarry Smith         s3 -= v[2] * x1 + v[6] * x2 + v[10] * x3 + v[14] * x4;
602c733ed4SBarry Smith         s4 -= v[3] * x1 + v[7] * x2 + v[11] * x3 + v[15] * x4;
612c733ed4SBarry Smith         v += 16;
622c733ed4SBarry Smith       }
632c733ed4SBarry Smith       x[idx]     = s1;
642c733ed4SBarry Smith       x[1 + idx] = s2;
652c733ed4SBarry Smith       x[2 + idx] = s3;
662c733ed4SBarry Smith       x[3 + idx] = s4;
672c733ed4SBarry Smith     }
682c733ed4SBarry Smith     /* backward solve the upper triangular */
692c733ed4SBarry Smith     idt = 4 * (n - 1);
702c733ed4SBarry Smith     for (i = n - 1; i >= 0; i--) {
712c733ed4SBarry Smith       ai16 = 16 * diag[i];
722c733ed4SBarry Smith       v    = aa + ai16 + 16;
732c733ed4SBarry Smith       vi   = aj + diag[i] + 1;
742c733ed4SBarry Smith       nz   = ai[i + 1] - diag[i] - 1;
759371c9d4SSatish Balay       s1   = x[idt];
769371c9d4SSatish Balay       s2   = x[1 + idt];
779371c9d4SSatish Balay       s3   = x[2 + idt];
789371c9d4SSatish Balay       s4   = x[3 + idt];
792c733ed4SBarry Smith       while (nz--) {
802c733ed4SBarry Smith         idx = 4 * (*vi++);
819371c9d4SSatish Balay         x1  = x[idx];
829371c9d4SSatish Balay         x2  = x[1 + idx];
839371c9d4SSatish Balay         x3  = x[2 + idx];
849371c9d4SSatish Balay         x4  = x[3 + idx];
852c733ed4SBarry Smith         s1 -= v[0] * x1 + v[4] * x2 + v[8] * x3 + v[12] * x4;
862c733ed4SBarry Smith         s2 -= v[1] * x1 + v[5] * x2 + v[9] * x3 + v[13] * x4;
872c733ed4SBarry Smith         s3 -= v[2] * x1 + v[6] * x2 + v[10] * x3 + v[14] * x4;
882c733ed4SBarry Smith         s4 -= v[3] * x1 + v[7] * x2 + v[11] * x3 + v[15] * x4;
892c733ed4SBarry Smith         v += 16;
902c733ed4SBarry Smith       }
912c733ed4SBarry Smith       v          = aa + ai16;
922c733ed4SBarry Smith       x[idt]     = v[0] * s1 + v[4] * s2 + v[8] * s3 + v[12] * s4;
932c733ed4SBarry Smith       x[1 + idt] = v[1] * s1 + v[5] * s2 + v[9] * s3 + v[13] * s4;
942c733ed4SBarry Smith       x[2 + idt] = v[2] * s1 + v[6] * s2 + v[10] * s3 + v[14] * s4;
952c733ed4SBarry Smith       x[3 + idt] = v[3] * s1 + v[7] * s2 + v[11] * s3 + v[15] * s4;
962c733ed4SBarry Smith       idt -= 4;
972c733ed4SBarry Smith     }
982c733ed4SBarry Smith   }
992c733ed4SBarry Smith #endif
1002c733ed4SBarry Smith 
1019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1039566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 16 * (a->nz) - 4.0 * A->cmap->n));
104*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1052c733ed4SBarry Smith }
1062c733ed4SBarry Smith 
MatSolve_SeqBAIJ_4_NaturalOrdering(Mat A,Vec bb,Vec xx)107d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_4_NaturalOrdering(Mat A, Vec bb, Vec xx)
108d71ae5a4SJacob Faibussowitsch {
1092c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data;
1102c733ed4SBarry Smith   const PetscInt     n = a->mbs, *vi, *ai = a->i, *aj = a->j, *adiag = a->diag;
1112c733ed4SBarry Smith   PetscInt           i, k, nz, idx, jdx, idt;
1122c733ed4SBarry Smith   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
1132c733ed4SBarry Smith   const MatScalar   *aa = a->a, *v;
1142c733ed4SBarry Smith   PetscScalar       *x;
1152c733ed4SBarry Smith   const PetscScalar *b;
1162c733ed4SBarry Smith   PetscScalar        s1, s2, s3, s4, x1, x2, x3, x4;
1172c733ed4SBarry Smith 
1182c733ed4SBarry Smith   PetscFunctionBegin;
1199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1209566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
1212c733ed4SBarry Smith   /* forward solve the lower triangular */
1222c733ed4SBarry Smith   idx  = 0;
1239371c9d4SSatish Balay   x[0] = b[idx];
1249371c9d4SSatish Balay   x[1] = b[1 + idx];
1259371c9d4SSatish Balay   x[2] = b[2 + idx];
1269371c9d4SSatish Balay   x[3] = b[3 + idx];
1272c733ed4SBarry Smith   for (i = 1; i < n; i++) {
1282c733ed4SBarry Smith     v   = aa + bs2 * ai[i];
1292c733ed4SBarry Smith     vi  = aj + ai[i];
1302c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1312c733ed4SBarry Smith     idx = bs * i;
1329371c9d4SSatish Balay     s1  = b[idx];
1339371c9d4SSatish Balay     s2  = b[1 + idx];
1349371c9d4SSatish Balay     s3  = b[2 + idx];
1359371c9d4SSatish Balay     s4  = b[3 + idx];
1362c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
1372c733ed4SBarry Smith       jdx = bs * vi[k];
1389371c9d4SSatish Balay       x1  = x[jdx];
1399371c9d4SSatish Balay       x2  = x[1 + jdx];
1409371c9d4SSatish Balay       x3  = x[2 + jdx];
1419371c9d4SSatish Balay       x4  = x[3 + jdx];
1422c733ed4SBarry Smith       s1 -= v[0] * x1 + v[4] * x2 + v[8] * x3 + v[12] * x4;
1432c733ed4SBarry Smith       s2 -= v[1] * x1 + v[5] * x2 + v[9] * x3 + v[13] * x4;
1442c733ed4SBarry Smith       s3 -= v[2] * x1 + v[6] * x2 + v[10] * x3 + v[14] * x4;
1452c733ed4SBarry Smith       s4 -= v[3] * x1 + v[7] * x2 + v[11] * x3 + v[15] * x4;
1462c733ed4SBarry Smith 
1472c733ed4SBarry Smith       v += bs2;
1482c733ed4SBarry Smith     }
1492c733ed4SBarry Smith 
1502c733ed4SBarry Smith     x[idx]     = s1;
1512c733ed4SBarry Smith     x[1 + idx] = s2;
1522c733ed4SBarry Smith     x[2 + idx] = s3;
1532c733ed4SBarry Smith     x[3 + idx] = s4;
1542c733ed4SBarry Smith   }
1552c733ed4SBarry Smith 
1562c733ed4SBarry Smith   /* backward solve the upper triangular */
1572c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1582c733ed4SBarry Smith     v   = aa + bs2 * (adiag[i + 1] + 1);
1592c733ed4SBarry Smith     vi  = aj + adiag[i + 1] + 1;
1602c733ed4SBarry Smith     nz  = adiag[i] - adiag[i + 1] - 1;
1612c733ed4SBarry Smith     idt = bs * i;
1629371c9d4SSatish Balay     s1  = x[idt];
1639371c9d4SSatish Balay     s2  = x[1 + idt];
1649371c9d4SSatish Balay     s3  = x[2 + idt];
1659371c9d4SSatish Balay     s4  = x[3 + idt];
1662c733ed4SBarry Smith 
1672c733ed4SBarry Smith     for (k = 0; k < nz; k++) {
1682c733ed4SBarry Smith       idx = bs * vi[k];
1699371c9d4SSatish Balay       x1  = x[idx];
1709371c9d4SSatish Balay       x2  = x[1 + idx];
1719371c9d4SSatish Balay       x3  = x[2 + idx];
1729371c9d4SSatish Balay       x4  = x[3 + idx];
1732c733ed4SBarry Smith       s1 -= v[0] * x1 + v[4] * x2 + v[8] * x3 + v[12] * x4;
1742c733ed4SBarry Smith       s2 -= v[1] * x1 + v[5] * x2 + v[9] * x3 + v[13] * x4;
1752c733ed4SBarry Smith       s3 -= v[2] * x1 + v[6] * x2 + v[10] * x3 + v[14] * x4;
1762c733ed4SBarry Smith       s4 -= v[3] * x1 + v[7] * x2 + v[11] * x3 + v[15] * x4;
1772c733ed4SBarry Smith 
1782c733ed4SBarry Smith       v += bs2;
1792c733ed4SBarry Smith     }
1802c733ed4SBarry Smith     /* x = inv_diagonal*x */
1812c733ed4SBarry Smith     x[idt]     = v[0] * s1 + v[4] * s2 + v[8] * s3 + v[12] * s4;
182feb237baSPierre Jolivet     x[1 + idt] = v[1] * s1 + v[5] * s2 + v[9] * s3 + v[13] * s4;
1832c733ed4SBarry Smith     x[2 + idt] = v[2] * s1 + v[6] * s2 + v[10] * s3 + v[14] * s4;
1842c733ed4SBarry Smith     x[3 + idt] = v[3] * s1 + v[7] * s2 + v[11] * s3 + v[15] * s4;
1852c733ed4SBarry Smith   }
1862c733ed4SBarry Smith 
1879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
1889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1899566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
190*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1912c733ed4SBarry Smith }
192