xref: /petsc/src/mat/impls/baij/seq/baijsolvtrannat6.c (revision 31d78bcd2b98084dc1368b20eb1129c8b9fb39fe)
12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith 
MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)3d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx)
4d71ae5a4SJacob Faibussowitsch {
52c733ed4SBarry Smith   Mat_SeqBAIJ     *a    = (Mat_SeqBAIJ *)A->data;
62c733ed4SBarry Smith   const PetscInt  *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j;
72c733ed4SBarry Smith   PetscInt         i, nz, idx, idt, oidx;
82c733ed4SBarry Smith   const MatScalar *aa = a->a, *v;
92c733ed4SBarry Smith   PetscScalar      s1, s2, s3, s4, s5, s6, x1, x2, x3, x4, x5, x6, *x;
102c733ed4SBarry Smith 
112c733ed4SBarry Smith   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(VecCopy(bb, xx));
139566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
142c733ed4SBarry Smith 
152c733ed4SBarry Smith   /* forward solve the U^T */
162c733ed4SBarry Smith   idx = 0;
172c733ed4SBarry Smith   for (i = 0; i < n; i++) {
182c733ed4SBarry Smith     v = aa + 36 * diag[i];
192c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
209371c9d4SSatish Balay     x1 = x[idx];
219371c9d4SSatish Balay     x2 = x[1 + idx];
229371c9d4SSatish Balay     x3 = x[2 + idx];
239371c9d4SSatish Balay     x4 = x[3 + idx];
249371c9d4SSatish Balay     x5 = x[4 + idx];
252c733ed4SBarry Smith     x6 = x[5 + idx];
262c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6;
272c733ed4SBarry Smith     s2 = v[6] * x1 + v[7] * x2 + v[8] * x3 + v[9] * x4 + v[10] * x5 + v[11] * x6;
282c733ed4SBarry Smith     s3 = v[12] * x1 + v[13] * x2 + v[14] * x3 + v[15] * x4 + v[16] * x5 + v[17] * x6;
292c733ed4SBarry Smith     s4 = v[18] * x1 + v[19] * x2 + v[20] * x3 + v[21] * x4 + v[22] * x5 + v[23] * x6;
302c733ed4SBarry Smith     s5 = v[24] * x1 + v[25] * x2 + v[26] * x3 + v[27] * x4 + v[28] * x5 + v[29] * x6;
312c733ed4SBarry Smith     s6 = v[30] * x1 + v[31] * x2 + v[32] * x3 + v[33] * x4 + v[34] * x5 + v[35] * x6;
322c733ed4SBarry Smith     v += 36;
332c733ed4SBarry Smith 
342c733ed4SBarry Smith     vi = aj + diag[i] + 1;
352c733ed4SBarry Smith     nz = ai[i + 1] - diag[i] - 1;
362c733ed4SBarry Smith     while (nz--) {
372c733ed4SBarry Smith       oidx = 6 * (*vi++);
382c733ed4SBarry Smith       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
392c733ed4SBarry Smith       x[oidx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
402c733ed4SBarry Smith       x[oidx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
412c733ed4SBarry Smith       x[oidx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
422c733ed4SBarry Smith       x[oidx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
432c733ed4SBarry Smith       x[oidx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
442c733ed4SBarry Smith       v += 36;
452c733ed4SBarry Smith     }
469371c9d4SSatish Balay     x[idx]     = s1;
479371c9d4SSatish Balay     x[1 + idx] = s2;
489371c9d4SSatish Balay     x[2 + idx] = s3;
499371c9d4SSatish Balay     x[3 + idx] = s4;
509371c9d4SSatish Balay     x[4 + idx] = s5;
512c733ed4SBarry Smith     x[5 + idx] = s6;
522c733ed4SBarry Smith     idx += 6;
532c733ed4SBarry Smith   }
542c733ed4SBarry Smith   /* backward solve the L^T */
552c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
562c733ed4SBarry Smith     v   = aa + 36 * diag[i] - 36;
572c733ed4SBarry Smith     vi  = aj + diag[i] - 1;
582c733ed4SBarry Smith     nz  = diag[i] - ai[i];
592c733ed4SBarry Smith     idt = 6 * i;
609371c9d4SSatish Balay     s1  = x[idt];
619371c9d4SSatish Balay     s2  = x[1 + idt];
629371c9d4SSatish Balay     s3  = x[2 + idt];
639371c9d4SSatish Balay     s4  = x[3 + idt];
649371c9d4SSatish Balay     s5  = x[4 + idt];
652c733ed4SBarry Smith     s6  = x[5 + idt];
662c733ed4SBarry Smith     while (nz--) {
672c733ed4SBarry Smith       idx = 6 * (*vi--);
682c733ed4SBarry Smith       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
692c733ed4SBarry Smith       x[idx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
702c733ed4SBarry Smith       x[idx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
712c733ed4SBarry Smith       x[idx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
722c733ed4SBarry Smith       x[idx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
732c733ed4SBarry Smith       x[idx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
742c733ed4SBarry Smith       v -= 36;
752c733ed4SBarry Smith     }
762c733ed4SBarry Smith   }
779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
789566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 36 * (a->nz) - 6.0 * A->cmap->n));
79*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
802c733ed4SBarry Smith }
812c733ed4SBarry Smith 
MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat A,Vec bb,Vec xx)82d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat A, Vec bb, Vec xx)
83d71ae5a4SJacob Faibussowitsch {
842c733ed4SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data;
852c733ed4SBarry Smith   const PetscInt   n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
862c733ed4SBarry Smith   PetscInt         nz, idx, idt, j, i, oidx;
872c733ed4SBarry Smith   const PetscInt   bs = A->rmap->bs, bs2 = a->bs2;
882c733ed4SBarry Smith   const MatScalar *aa = a->a, *v;
892c733ed4SBarry Smith   PetscScalar      s1, s2, s3, s4, s5, s6, x1, x2, x3, x4, x5, x6, *x;
902c733ed4SBarry Smith 
912c733ed4SBarry Smith   PetscFunctionBegin;
929566063dSJacob Faibussowitsch   PetscCall(VecCopy(bb, xx));
939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
942c733ed4SBarry Smith 
952c733ed4SBarry Smith   /* forward solve the U^T */
962c733ed4SBarry Smith   idx = 0;
972c733ed4SBarry Smith   for (i = 0; i < n; i++) {
982c733ed4SBarry Smith     v = aa + bs2 * diag[i];
992c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
1009371c9d4SSatish Balay     x1 = x[idx];
1019371c9d4SSatish Balay     x2 = x[1 + idx];
1029371c9d4SSatish Balay     x3 = x[2 + idx];
1039371c9d4SSatish Balay     x4 = x[3 + idx];
1049371c9d4SSatish Balay     x5 = x[4 + idx];
1059371c9d4SSatish Balay     x6 = x[5 + idx];
1062c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6;
1072c733ed4SBarry Smith     s2 = v[6] * x1 + v[7] * x2 + v[8] * x3 + v[9] * x4 + v[10] * x5 + v[11] * x6;
1082c733ed4SBarry Smith     s3 = v[12] * x1 + v[13] * x2 + v[14] * x3 + v[15] * x4 + v[16] * x5 + v[17] * x6;
1092c733ed4SBarry Smith     s4 = v[18] * x1 + v[19] * x2 + v[20] * x3 + v[21] * x4 + v[22] * x5 + v[23] * x6;
1102c733ed4SBarry Smith     s5 = v[24] * x1 + v[25] * x2 + v[26] * x3 + v[27] * x4 + v[28] * x5 + v[29] * x6;
1112c733ed4SBarry Smith     s6 = v[30] * x1 + v[31] * x2 + v[32] * x3 + v[33] * x4 + v[34] * x5 + v[35] * x6;
1122c733ed4SBarry Smith     v -= bs2;
1132c733ed4SBarry Smith 
1142c733ed4SBarry Smith     vi = aj + diag[i] - 1;
1152c733ed4SBarry Smith     nz = diag[i] - diag[i + 1] - 1;
1162c733ed4SBarry Smith     for (j = 0; j > -nz; j--) {
1172c733ed4SBarry Smith       oidx = bs * vi[j];
1182c733ed4SBarry Smith       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
1192c733ed4SBarry Smith       x[oidx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
1202c733ed4SBarry Smith       x[oidx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
1212c733ed4SBarry Smith       x[oidx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
1222c733ed4SBarry Smith       x[oidx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
1232c733ed4SBarry Smith       x[oidx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
1242c733ed4SBarry Smith       v -= bs2;
1252c733ed4SBarry Smith     }
1269371c9d4SSatish Balay     x[idx]     = s1;
1279371c9d4SSatish Balay     x[1 + idx] = s2;
1289371c9d4SSatish Balay     x[2 + idx] = s3;
1299371c9d4SSatish Balay     x[3 + idx] = s4;
1309371c9d4SSatish Balay     x[4 + idx] = s5;
1312c733ed4SBarry Smith     x[5 + idx] = s6;
1322c733ed4SBarry Smith     idx += bs;
1332c733ed4SBarry Smith   }
1342c733ed4SBarry Smith   /* backward solve the L^T */
1352c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1362c733ed4SBarry Smith     v   = aa + bs2 * ai[i];
1372c733ed4SBarry Smith     vi  = aj + ai[i];
1382c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1392c733ed4SBarry Smith     idt = bs * i;
1409371c9d4SSatish Balay     s1  = x[idt];
1419371c9d4SSatish Balay     s2  = x[1 + idt];
1429371c9d4SSatish Balay     s3  = x[2 + idt];
1439371c9d4SSatish Balay     s4  = x[3 + idt];
1449371c9d4SSatish Balay     s5  = x[4 + idt];
1452c733ed4SBarry Smith     s6  = x[5 + idt];
1462c733ed4SBarry Smith     for (j = 0; j < nz; j++) {
1472c733ed4SBarry Smith       idx = bs * vi[j];
1482c733ed4SBarry Smith       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6;
1492c733ed4SBarry Smith       x[idx + 1] -= v[6] * s1 + v[7] * s2 + v[8] * s3 + v[9] * s4 + v[10] * s5 + v[11] * s6;
1502c733ed4SBarry Smith       x[idx + 2] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4 + v[16] * s5 + v[17] * s6;
1512c733ed4SBarry Smith       x[idx + 3] -= v[18] * s1 + v[19] * s2 + v[20] * s3 + v[21] * s4 + v[22] * s5 + v[23] * s6;
1522c733ed4SBarry Smith       x[idx + 4] -= v[24] * s1 + v[25] * s2 + v[26] * s3 + v[27] * s4 + v[28] * s5 + v[29] * s6;
1532c733ed4SBarry Smith       x[idx + 5] -= v[30] * s1 + v[31] * s2 + v[32] * s3 + v[33] * s4 + v[34] * s5 + v[35] * s6;
1542c733ed4SBarry Smith       v += bs2;
1552c733ed4SBarry Smith     }
1562c733ed4SBarry Smith   }
1579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1589566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
159*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1602c733ed4SBarry Smith }
161