Actual source code: baijsolvtrannat5.c
1: #include <../src/mat/impls/baij/seq/baij.h>
3: PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx)
4: {
5: Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
6: const PetscInt *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j;
7: PetscInt i, nz, idx, idt, oidx;
8: const MatScalar *aa = a->a, *v;
9: PetscScalar s1, s2, s3, s4, s5, x1, x2, x3, x4, x5, *x;
11: PetscFunctionBegin;
12: PetscCall(VecCopy(bb, xx));
13: PetscCall(VecGetArray(xx, &x));
15: /* forward solve the U^T */
16: idx = 0;
17: for (i = 0; i < n; i++) {
18: v = aa + 25 * diag[i];
19: /* multiply by the inverse of the block diagonal */
20: x1 = x[idx];
21: x2 = x[1 + idx];
22: x3 = x[2 + idx];
23: x4 = x[3 + idx];
24: x5 = x[4 + idx];
25: s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5;
26: s2 = v[5] * x1 + v[6] * x2 + v[7] * x3 + v[8] * x4 + v[9] * x5;
27: s3 = v[10] * x1 + v[11] * x2 + v[12] * x3 + v[13] * x4 + v[14] * x5;
28: s4 = v[15] * x1 + v[16] * x2 + v[17] * x3 + v[18] * x4 + v[19] * x5;
29: s5 = v[20] * x1 + v[21] * x2 + v[22] * x3 + v[23] * x4 + v[24] * x5;
30: v += 25;
32: vi = aj + diag[i] + 1;
33: nz = ai[i + 1] - diag[i] - 1;
34: while (nz--) {
35: oidx = 5 * (*vi++);
36: x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
37: x[oidx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
38: x[oidx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
39: x[oidx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
40: x[oidx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
41: v += 25;
42: }
43: x[idx] = s1;
44: x[1 + idx] = s2;
45: x[2 + idx] = s3;
46: x[3 + idx] = s4;
47: x[4 + idx] = s5;
48: idx += 5;
49: }
50: /* backward solve the L^T */
51: for (i = n - 1; i >= 0; i--) {
52: v = aa + 25 * diag[i] - 25;
53: vi = aj + diag[i] - 1;
54: nz = diag[i] - ai[i];
55: idt = 5 * i;
56: s1 = x[idt];
57: s2 = x[1 + idt];
58: s3 = x[2 + idt];
59: s4 = x[3 + idt];
60: s5 = x[4 + idt];
61: while (nz--) {
62: idx = 5 * (*vi--);
63: x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
64: x[idx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
65: x[idx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
66: x[idx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
67: x[idx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
68: v -= 25;
69: }
70: }
71: PetscCall(VecRestoreArray(xx, &x));
72: PetscCall(PetscLogFlops(2.0 * 25 * (a->nz) - 5.0 * A->cmap->n));
73: PetscFunctionReturn(PETSC_SUCCESS);
74: }
76: PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat A, Vec bb, Vec xx)
77: {
78: Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
79: const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
80: PetscInt nz, idx, idt, j, i, oidx;
81: const PetscInt bs = A->rmap->bs, bs2 = a->bs2;
82: const MatScalar *aa = a->a, *v;
83: PetscScalar s1, s2, s3, s4, s5, x1, x2, x3, x4, x5, *x;
85: PetscFunctionBegin;
86: PetscCall(VecCopy(bb, xx));
87: PetscCall(VecGetArray(xx, &x));
89: /* forward solve the U^T */
90: idx = 0;
91: for (i = 0; i < n; i++) {
92: v = aa + bs2 * diag[i];
93: /* multiply by the inverse of the block diagonal */
94: x1 = x[idx];
95: x2 = x[1 + idx];
96: x3 = x[2 + idx];
97: x4 = x[3 + idx];
98: x5 = x[4 + idx];
99: s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5;
100: s2 = v[5] * x1 + v[6] * x2 + v[7] * x3 + v[8] * x4 + v[9] * x5;
101: s3 = v[10] * x1 + v[11] * x2 + v[12] * x3 + v[13] * x4 + v[14] * x5;
102: s4 = v[15] * x1 + v[16] * x2 + v[17] * x3 + v[18] * x4 + v[19] * x5;
103: s5 = v[20] * x1 + v[21] * x2 + v[22] * x3 + v[23] * x4 + v[24] * x5;
104: v -= bs2;
106: vi = aj + diag[i] - 1;
107: nz = diag[i] - diag[i + 1] - 1;
108: for (j = 0; j > -nz; j--) {
109: oidx = bs * vi[j];
110: x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
111: x[oidx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
112: x[oidx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
113: x[oidx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
114: x[oidx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
115: v -= bs2;
116: }
117: x[idx] = s1;
118: x[1 + idx] = s2;
119: x[2 + idx] = s3;
120: x[3 + idx] = s4;
121: x[4 + idx] = s5;
122: idx += bs;
123: }
124: /* backward solve the L^T */
125: for (i = n - 1; i >= 0; i--) {
126: v = aa + bs2 * ai[i];
127: vi = aj + ai[i];
128: nz = ai[i + 1] - ai[i];
129: idt = bs * i;
130: s1 = x[idt];
131: s2 = x[1 + idt];
132: s3 = x[2 + idt];
133: s4 = x[3 + idt];
134: s5 = x[4 + idt];
135: for (j = 0; j < nz; j++) {
136: idx = bs * vi[j];
137: x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
138: x[idx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
139: x[idx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
140: x[idx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
141: x[idx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
142: v += bs2;
143: }
144: }
145: PetscCall(VecRestoreArray(xx, &x));
146: PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
147: PetscFunctionReturn(PETSC_SUCCESS);
148: }