12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith
MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)3d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_5_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, x1, x2, x3, x4, x5, *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 + 25 * 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 s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5;
262c733ed4SBarry Smith s2 = v[5] * x1 + v[6] * x2 + v[7] * x3 + v[8] * x4 + v[9] * x5;
272c733ed4SBarry Smith s3 = v[10] * x1 + v[11] * x2 + v[12] * x3 + v[13] * x4 + v[14] * x5;
282c733ed4SBarry Smith s4 = v[15] * x1 + v[16] * x2 + v[17] * x3 + v[18] * x4 + v[19] * x5;
292c733ed4SBarry Smith s5 = v[20] * x1 + v[21] * x2 + v[22] * x3 + v[23] * x4 + v[24] * x5;
302c733ed4SBarry Smith v += 25;
312c733ed4SBarry Smith
322c733ed4SBarry Smith vi = aj + diag[i] + 1;
332c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1;
342c733ed4SBarry Smith while (nz--) {
352c733ed4SBarry Smith oidx = 5 * (*vi++);
362c733ed4SBarry Smith x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
372c733ed4SBarry Smith x[oidx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
382c733ed4SBarry Smith x[oidx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
392c733ed4SBarry Smith x[oidx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
402c733ed4SBarry Smith x[oidx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
412c733ed4SBarry Smith v += 25;
422c733ed4SBarry Smith }
439371c9d4SSatish Balay x[idx] = s1;
449371c9d4SSatish Balay x[1 + idx] = s2;
459371c9d4SSatish Balay x[2 + idx] = s3;
469371c9d4SSatish Balay x[3 + idx] = s4;
479371c9d4SSatish Balay x[4 + idx] = s5;
482c733ed4SBarry Smith idx += 5;
492c733ed4SBarry Smith }
502c733ed4SBarry Smith /* backward solve the L^T */
512c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) {
522c733ed4SBarry Smith v = aa + 25 * diag[i] - 25;
532c733ed4SBarry Smith vi = aj + diag[i] - 1;
542c733ed4SBarry Smith nz = diag[i] - ai[i];
552c733ed4SBarry Smith idt = 5 * i;
569371c9d4SSatish Balay s1 = x[idt];
579371c9d4SSatish Balay s2 = x[1 + idt];
589371c9d4SSatish Balay s3 = x[2 + idt];
599371c9d4SSatish Balay s4 = x[3 + idt];
609371c9d4SSatish Balay s5 = x[4 + idt];
612c733ed4SBarry Smith while (nz--) {
622c733ed4SBarry Smith idx = 5 * (*vi--);
632c733ed4SBarry Smith x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
642c733ed4SBarry Smith x[idx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
652c733ed4SBarry Smith x[idx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
662c733ed4SBarry Smith x[idx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
672c733ed4SBarry Smith x[idx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
682c733ed4SBarry Smith v -= 25;
692c733ed4SBarry Smith }
702c733ed4SBarry Smith }
719566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
729566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 25 * (a->nz) - 5.0 * A->cmap->n));
73*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
742c733ed4SBarry Smith }
752c733ed4SBarry Smith
MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat A,Vec bb,Vec xx)76d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat A, Vec bb, Vec xx)
77d71ae5a4SJacob Faibussowitsch {
782c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
792c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
802c733ed4SBarry Smith PetscInt nz, idx, idt, j, i, oidx;
812c733ed4SBarry Smith const PetscInt bs = A->rmap->bs, bs2 = a->bs2;
822c733ed4SBarry Smith const MatScalar *aa = a->a, *v;
832c733ed4SBarry Smith PetscScalar s1, s2, s3, s4, s5, x1, x2, x3, x4, x5, *x;
842c733ed4SBarry Smith
852c733ed4SBarry Smith PetscFunctionBegin;
869566063dSJacob Faibussowitsch PetscCall(VecCopy(bb, xx));
879566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x));
882c733ed4SBarry Smith
892c733ed4SBarry Smith /* forward solve the U^T */
902c733ed4SBarry Smith idx = 0;
912c733ed4SBarry Smith for (i = 0; i < n; i++) {
922c733ed4SBarry Smith v = aa + bs2 * diag[i];
932c733ed4SBarry Smith /* multiply by the inverse of the block diagonal */
949371c9d4SSatish Balay x1 = x[idx];
959371c9d4SSatish Balay x2 = x[1 + idx];
969371c9d4SSatish Balay x3 = x[2 + idx];
979371c9d4SSatish Balay x4 = x[3 + idx];
982c733ed4SBarry Smith x5 = x[4 + idx];
992c733ed4SBarry Smith s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5;
1002c733ed4SBarry Smith s2 = v[5] * x1 + v[6] * x2 + v[7] * x3 + v[8] * x4 + v[9] * x5;
1012c733ed4SBarry Smith s3 = v[10] * x1 + v[11] * x2 + v[12] * x3 + v[13] * x4 + v[14] * x5;
1022c733ed4SBarry Smith s4 = v[15] * x1 + v[16] * x2 + v[17] * x3 + v[18] * x4 + v[19] * x5;
1032c733ed4SBarry Smith s5 = v[20] * x1 + v[21] * x2 + v[22] * x3 + v[23] * x4 + v[24] * x5;
1042c733ed4SBarry Smith v -= bs2;
1052c733ed4SBarry Smith
1062c733ed4SBarry Smith vi = aj + diag[i] - 1;
1072c733ed4SBarry Smith nz = diag[i] - diag[i + 1] - 1;
1082c733ed4SBarry Smith for (j = 0; j > -nz; j--) {
1092c733ed4SBarry Smith oidx = bs * vi[j];
1102c733ed4SBarry Smith x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
1112c733ed4SBarry Smith x[oidx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
1122c733ed4SBarry Smith x[oidx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
1132c733ed4SBarry Smith x[oidx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
1142c733ed4SBarry Smith x[oidx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
1152c733ed4SBarry Smith v -= bs2;
1162c733ed4SBarry Smith }
1179371c9d4SSatish Balay x[idx] = s1;
1189371c9d4SSatish Balay x[1 + idx] = s2;
1199371c9d4SSatish Balay x[2 + idx] = s3;
1209371c9d4SSatish Balay x[3 + idx] = s4;
1219371c9d4SSatish Balay x[4 + idx] = s5;
1222c733ed4SBarry Smith idx += bs;
1232c733ed4SBarry Smith }
1242c733ed4SBarry Smith /* backward solve the L^T */
1252c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) {
1262c733ed4SBarry Smith v = aa + bs2 * ai[i];
1272c733ed4SBarry Smith vi = aj + ai[i];
1282c733ed4SBarry Smith nz = ai[i + 1] - ai[i];
1292c733ed4SBarry Smith idt = bs * i;
1309371c9d4SSatish Balay s1 = x[idt];
1319371c9d4SSatish Balay s2 = x[1 + idt];
1329371c9d4SSatish Balay s3 = x[2 + idt];
1339371c9d4SSatish Balay s4 = x[3 + idt];
1349371c9d4SSatish Balay s5 = x[4 + idt];
1352c733ed4SBarry Smith for (j = 0; j < nz; j++) {
1362c733ed4SBarry Smith idx = bs * vi[j];
1372c733ed4SBarry Smith x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
1382c733ed4SBarry Smith x[idx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
1392c733ed4SBarry Smith x[idx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
1402c733ed4SBarry Smith x[idx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
1412c733ed4SBarry Smith x[idx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
1422c733ed4SBarry Smith v += bs2;
1432c733ed4SBarry Smith }
1442c733ed4SBarry Smith }
1459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x));
1469566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
147*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1482c733ed4SBarry Smith }
149