xref: /petsc/src/mat/impls/baij/seq/baijsolvtrannat4.c (revision 31d78bcd2b98084dc1368b20eb1129c8b9fb39fe)
12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith 
MatSolveTranspose_SeqBAIJ_4_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)3d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_4_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, x1, x2, x3, x4, *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 + 16 * 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];
242c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4;
252c733ed4SBarry Smith     s2 = v[4] * x1 + v[5] * x2 + v[6] * x3 + v[7] * x4;
262c733ed4SBarry Smith     s3 = v[8] * x1 + v[9] * x2 + v[10] * x3 + v[11] * x4;
272c733ed4SBarry Smith     s4 = v[12] * x1 + v[13] * x2 + v[14] * x3 + v[15] * x4;
282c733ed4SBarry Smith     v += 16;
292c733ed4SBarry Smith 
302c733ed4SBarry Smith     vi = aj + diag[i] + 1;
312c733ed4SBarry Smith     nz = ai[i + 1] - diag[i] - 1;
322c733ed4SBarry Smith     while (nz--) {
332c733ed4SBarry Smith       oidx = 4 * (*vi++);
342c733ed4SBarry Smith       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4;
352c733ed4SBarry Smith       x[oidx + 1] -= v[4] * s1 + v[5] * s2 + v[6] * s3 + v[7] * s4;
362c733ed4SBarry Smith       x[oidx + 2] -= v[8] * s1 + v[9] * s2 + v[10] * s3 + v[11] * s4;
372c733ed4SBarry Smith       x[oidx + 3] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4;
382c733ed4SBarry Smith       v += 16;
392c733ed4SBarry Smith     }
409371c9d4SSatish Balay     x[idx]     = s1;
419371c9d4SSatish Balay     x[1 + idx] = s2;
429371c9d4SSatish Balay     x[2 + idx] = s3;
439371c9d4SSatish Balay     x[3 + idx] = s4;
442c733ed4SBarry Smith     idx += 4;
452c733ed4SBarry Smith   }
462c733ed4SBarry Smith   /* backward solve the L^T */
472c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
482c733ed4SBarry Smith     v   = aa + 16 * diag[i] - 16;
492c733ed4SBarry Smith     vi  = aj + diag[i] - 1;
502c733ed4SBarry Smith     nz  = diag[i] - ai[i];
512c733ed4SBarry Smith     idt = 4 * i;
529371c9d4SSatish Balay     s1  = x[idt];
539371c9d4SSatish Balay     s2  = x[1 + idt];
549371c9d4SSatish Balay     s3  = x[2 + idt];
559371c9d4SSatish Balay     s4  = x[3 + idt];
562c733ed4SBarry Smith     while (nz--) {
572c733ed4SBarry Smith       idx = 4 * (*vi--);
582c733ed4SBarry Smith       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4;
592c733ed4SBarry Smith       x[idx + 1] -= v[4] * s1 + v[5] * s2 + v[6] * s3 + v[7] * s4;
602c733ed4SBarry Smith       x[idx + 2] -= v[8] * s1 + v[9] * s2 + v[10] * s3 + v[11] * s4;
612c733ed4SBarry Smith       x[idx + 3] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4;
622c733ed4SBarry Smith       v -= 16;
632c733ed4SBarry Smith     }
642c733ed4SBarry Smith   }
659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
669566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * 16 * (a->nz) - 4.0 * A->cmap->n));
67*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
682c733ed4SBarry Smith }
692c733ed4SBarry Smith 
MatSolveTranspose_SeqBAIJ_4_NaturalOrdering(Mat A,Vec bb,Vec xx)70d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolveTranspose_SeqBAIJ_4_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, *diag = a->diag;
742c733ed4SBarry Smith   PetscInt         nz, idx, idt, j, i, oidx;
752c733ed4SBarry Smith   const PetscInt   bs = A->rmap->bs, bs2 = a->bs2;
762c733ed4SBarry Smith   const MatScalar *aa = a->a, *v;
772c733ed4SBarry Smith   PetscScalar      s1, s2, s3, s4, x1, x2, x3, x4, *x;
782c733ed4SBarry Smith 
792c733ed4SBarry Smith   PetscFunctionBegin;
809566063dSJacob Faibussowitsch   PetscCall(VecCopy(bb, xx));
819566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
822c733ed4SBarry Smith 
832c733ed4SBarry Smith   /* forward solve the U^T */
842c733ed4SBarry Smith   idx = 0;
852c733ed4SBarry Smith   for (i = 0; i < n; i++) {
862c733ed4SBarry Smith     v = aa + bs2 * diag[i];
872c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
889371c9d4SSatish Balay     x1 = x[idx];
899371c9d4SSatish Balay     x2 = x[1 + idx];
909371c9d4SSatish Balay     x3 = x[2 + idx];
919371c9d4SSatish Balay     x4 = x[3 + idx];
922c733ed4SBarry Smith     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4;
932c733ed4SBarry Smith     s2 = v[4] * x1 + v[5] * x2 + v[6] * x3 + v[7] * x4;
942c733ed4SBarry Smith     s3 = v[8] * x1 + v[9] * x2 + v[10] * x3 + v[11] * x4;
952c733ed4SBarry Smith     s4 = v[12] * x1 + v[13] * x2 + v[14] * x3 + v[15] * x4;
962c733ed4SBarry Smith     v -= bs2;
972c733ed4SBarry Smith 
982c733ed4SBarry Smith     vi = aj + diag[i] - 1;
992c733ed4SBarry Smith     nz = diag[i] - diag[i + 1] - 1;
1002c733ed4SBarry Smith     for (j = 0; j > -nz; j--) {
1012c733ed4SBarry Smith       oidx = bs * vi[j];
1022c733ed4SBarry Smith       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4;
1032c733ed4SBarry Smith       x[oidx + 1] -= v[4] * s1 + v[5] * s2 + v[6] * s3 + v[7] * s4;
1042c733ed4SBarry Smith       x[oidx + 2] -= v[8] * s1 + v[9] * s2 + v[10] * s3 + v[11] * s4;
1052c733ed4SBarry Smith       x[oidx + 3] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4;
1062c733ed4SBarry Smith       v -= bs2;
1072c733ed4SBarry Smith     }
1089371c9d4SSatish Balay     x[idx]     = s1;
1099371c9d4SSatish Balay     x[1 + idx] = s2;
1109371c9d4SSatish Balay     x[2 + idx] = s3;
1119371c9d4SSatish Balay     x[3 + idx] = s4;
1122c733ed4SBarry Smith     idx += bs;
1132c733ed4SBarry Smith   }
1142c733ed4SBarry Smith   /* backward solve the L^T */
1152c733ed4SBarry Smith   for (i = n - 1; i >= 0; i--) {
1162c733ed4SBarry Smith     v   = aa + bs2 * ai[i];
1172c733ed4SBarry Smith     vi  = aj + ai[i];
1182c733ed4SBarry Smith     nz  = ai[i + 1] - ai[i];
1192c733ed4SBarry Smith     idt = bs * i;
1209371c9d4SSatish Balay     s1  = x[idt];
1219371c9d4SSatish Balay     s2  = x[1 + idt];
1229371c9d4SSatish Balay     s3  = x[2 + idt];
1239371c9d4SSatish Balay     s4  = x[3 + idt];
1242c733ed4SBarry Smith     for (j = 0; j < nz; j++) {
1252c733ed4SBarry Smith       idx = bs * vi[j];
1262c733ed4SBarry Smith       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4;
1272c733ed4SBarry Smith       x[idx + 1] -= v[4] * s1 + v[5] * s2 + v[6] * s3 + v[7] * s4;
1282c733ed4SBarry Smith       x[idx + 2] -= v[8] * s1 + v[9] * s2 + v[10] * s3 + v[11] * s4;
1292c733ed4SBarry Smith       x[idx + 3] -= v[12] * s1 + v[13] * s2 + v[14] * s3 + v[15] * s4;
1302c733ed4SBarry Smith       v += bs2;
1312c733ed4SBarry Smith     }
1322c733ed4SBarry Smith   }
1339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
1349566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
135*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1362c733ed4SBarry Smith }
137