Actual source code: baijsolvnat7.c
1: #include <../src/mat/impls/baij/seq/baij.h>
2: #include <petsc/private/kernels/blockinvert.h>
4: PetscErrorCode MatSolve_SeqBAIJ_7_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx)
5: {
6: Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
7: const PetscInt *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j;
8: PetscInt i, nz, idx, idt, jdx;
9: const MatScalar *aa = a->a, *v;
10: PetscScalar *x, s1, s2, s3, s4, s5, s6, s7, x1, x2, x3, x4, x5, x6, x7;
11: const PetscScalar *b;
13: VecGetArrayRead(bb, &b);
14: VecGetArray(xx, &x);
15: /* forward solve the lower triangular */
16: idx = 0;
17: x[0] = b[idx];
18: x[1] = b[1 + idx];
19: x[2] = b[2 + idx];
20: x[3] = b[3 + idx];
21: x[4] = b[4 + idx];
22: x[5] = b[5 + idx];
23: x[6] = b[6 + idx];
24: for (i = 1; i < n; i++) {
25: v = aa + 49 * ai[i];
26: vi = aj + ai[i];
27: nz = diag[i] - ai[i];
28: idx = 7 * i;
29: s1 = b[idx];
30: s2 = b[1 + idx];
31: s3 = b[2 + idx];
32: s4 = b[3 + idx];
33: s5 = b[4 + idx];
34: s6 = b[5 + idx];
35: s7 = b[6 + idx];
36: while (nz--) {
37: jdx = 7 * (*vi++);
38: x1 = x[jdx];
39: x2 = x[1 + jdx];
40: x3 = x[2 + jdx];
41: x4 = x[3 + jdx];
42: x5 = x[4 + jdx];
43: x6 = x[5 + jdx];
44: x7 = x[6 + jdx];
45: s1 -= v[0] * x1 + v[7] * x2 + v[14] * x3 + v[21] * x4 + v[28] * x5 + v[35] * x6 + v[42] * x7;
46: s2 -= v[1] * x1 + v[8] * x2 + v[15] * x3 + v[22] * x4 + v[29] * x5 + v[36] * x6 + v[43] * x7;
47: s3 -= v[2] * x1 + v[9] * x2 + v[16] * x3 + v[23] * x4 + v[30] * x5 + v[37] * x6 + v[44] * x7;
48: s4 -= v[3] * x1 + v[10] * x2 + v[17] * x3 + v[24] * x4 + v[31] * x5 + v[38] * x6 + v[45] * x7;
49: s5 -= v[4] * x1 + v[11] * x2 + v[18] * x3 + v[25] * x4 + v[32] * x5 + v[39] * x6 + v[46] * x7;
50: s6 -= v[5] * x1 + v[12] * x2 + v[19] * x3 + v[26] * x4 + v[33] * x5 + v[40] * x6 + v[47] * x7;
51: s7 -= v[6] * x1 + v[13] * x2 + v[20] * x3 + v[27] * x4 + v[34] * x5 + v[41] * x6 + v[48] * x7;
52: v += 49;
53: }
54: x[idx] = s1;
55: x[1 + idx] = s2;
56: x[2 + idx] = s3;
57: x[3 + idx] = s4;
58: x[4 + idx] = s5;
59: x[5 + idx] = s6;
60: x[6 + idx] = s7;
61: }
62: /* backward solve the upper triangular */
63: for (i = n - 1; i >= 0; i--) {
64: v = aa + 49 * diag[i] + 49;
65: vi = aj + diag[i] + 1;
66: nz = ai[i + 1] - diag[i] - 1;
67: idt = 7 * i;
68: s1 = x[idt];
69: s2 = x[1 + idt];
70: s3 = x[2 + idt];
71: s4 = x[3 + idt];
72: s5 = x[4 + idt];
73: s6 = x[5 + idt];
74: s7 = x[6 + idt];
75: while (nz--) {
76: idx = 7 * (*vi++);
77: x1 = x[idx];
78: x2 = x[1 + idx];
79: x3 = x[2 + idx];
80: x4 = x[3 + idx];
81: x5 = x[4 + idx];
82: x6 = x[5 + idx];
83: x7 = x[6 + idx];
84: s1 -= v[0] * x1 + v[7] * x2 + v[14] * x3 + v[21] * x4 + v[28] * x5 + v[35] * x6 + v[42] * x7;
85: s2 -= v[1] * x1 + v[8] * x2 + v[15] * x3 + v[22] * x4 + v[29] * x5 + v[36] * x6 + v[43] * x7;
86: s3 -= v[2] * x1 + v[9] * x2 + v[16] * x3 + v[23] * x4 + v[30] * x5 + v[37] * x6 + v[44] * x7;
87: s4 -= v[3] * x1 + v[10] * x2 + v[17] * x3 + v[24] * x4 + v[31] * x5 + v[38] * x6 + v[45] * x7;
88: s5 -= v[4] * x1 + v[11] * x2 + v[18] * x3 + v[25] * x4 + v[32] * x5 + v[39] * x6 + v[46] * x7;
89: s6 -= v[5] * x1 + v[12] * x2 + v[19] * x3 + v[26] * x4 + v[33] * x5 + v[40] * x6 + v[47] * x7;
90: s7 -= v[6] * x1 + v[13] * x2 + v[20] * x3 + v[27] * x4 + v[34] * x5 + v[41] * x6 + v[48] * x7;
91: v += 49;
92: }
93: v = aa + 49 * diag[i];
94: x[idt] = v[0] * s1 + v[7] * s2 + v[14] * s3 + v[21] * s4 + v[28] * s5 + v[35] * s6 + v[42] * s7;
95: x[1 + idt] = v[1] * s1 + v[8] * s2 + v[15] * s3 + v[22] * s4 + v[29] * s5 + v[36] * s6 + v[43] * s7;
96: x[2 + idt] = v[2] * s1 + v[9] * s2 + v[16] * s3 + v[23] * s4 + v[30] * s5 + v[37] * s6 + v[44] * s7;
97: x[3 + idt] = v[3] * s1 + v[10] * s2 + v[17] * s3 + v[24] * s4 + v[31] * s5 + v[38] * s6 + v[45] * s7;
98: x[4 + idt] = v[4] * s1 + v[11] * s2 + v[18] * s3 + v[25] * s4 + v[32] * s5 + v[39] * s6 + v[46] * s7;
99: x[5 + idt] = v[5] * s1 + v[12] * s2 + v[19] * s3 + v[26] * s4 + v[33] * s5 + v[40] * s6 + v[47] * s7;
100: x[6 + idt] = v[6] * s1 + v[13] * s2 + v[20] * s3 + v[27] * s4 + v[34] * s5 + v[41] * s6 + v[48] * s7;
101: }
103: VecRestoreArrayRead(bb, &b);
104: VecRestoreArray(xx, &x);
105: PetscLogFlops(2.0 * 36 * (a->nz) - 6.0 * A->cmap->n);
106: return 0;
107: }
109: PetscErrorCode MatSolve_SeqBAIJ_7_NaturalOrdering(Mat A, Vec bb, Vec xx)
110: {
111: Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data;
112: const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *adiag = a->diag;
113: PetscInt i, k, nz, idx, jdx, idt;
114: const PetscInt bs = A->rmap->bs, bs2 = a->bs2;
115: const MatScalar *aa = a->a, *v;
116: PetscScalar *x;
117: const PetscScalar *b;
118: PetscScalar s1, s2, s3, s4, s5, s6, s7, x1, x2, x3, x4, x5, x6, x7;
120: VecGetArrayRead(bb, &b);
121: VecGetArray(xx, &x);
122: /* forward solve the lower triangular */
123: idx = 0;
124: x[0] = b[idx];
125: x[1] = b[1 + idx];
126: x[2] = b[2 + idx];
127: x[3] = b[3 + idx];
128: x[4] = b[4 + idx];
129: x[5] = b[5 + idx];
130: x[6] = b[6 + idx];
131: for (i = 1; i < n; i++) {
132: v = aa + bs2 * ai[i];
133: vi = aj + ai[i];
134: nz = ai[i + 1] - ai[i];
135: idx = bs * i;
136: s1 = b[idx];
137: s2 = b[1 + idx];
138: s3 = b[2 + idx];
139: s4 = b[3 + idx];
140: s5 = b[4 + idx];
141: s6 = b[5 + idx];
142: s7 = b[6 + idx];
143: for (k = 0; k < nz; k++) {
144: jdx = bs * vi[k];
145: x1 = x[jdx];
146: x2 = x[1 + jdx];
147: x3 = x[2 + jdx];
148: x4 = x[3 + jdx];
149: x5 = x[4 + jdx];
150: x6 = x[5 + jdx];
151: x7 = x[6 + jdx];
152: s1 -= v[0] * x1 + v[7] * x2 + v[14] * x3 + v[21] * x4 + v[28] * x5 + v[35] * x6 + v[42] * x7;
153: s2 -= v[1] * x1 + v[8] * x2 + v[15] * x3 + v[22] * x4 + v[29] * x5 + v[36] * x6 + v[43] * x7;
154: s3 -= v[2] * x1 + v[9] * x2 + v[16] * x3 + v[23] * x4 + v[30] * x5 + v[37] * x6 + v[44] * x7;
155: s4 -= v[3] * x1 + v[10] * x2 + v[17] * x3 + v[24] * x4 + v[31] * x5 + v[38] * x6 + v[45] * x7;
156: s5 -= v[4] * x1 + v[11] * x2 + v[18] * x3 + v[25] * x4 + v[32] * x5 + v[39] * x6 + v[46] * x7;
157: s6 -= v[5] * x1 + v[12] * x2 + v[19] * x3 + v[26] * x4 + v[33] * x5 + v[40] * x6 + v[47] * x7;
158: s7 -= v[6] * x1 + v[13] * x2 + v[20] * x3 + v[27] * x4 + v[34] * x5 + v[41] * x6 + v[48] * x7;
159: v += bs2;
160: }
162: x[idx] = s1;
163: x[1 + idx] = s2;
164: x[2 + idx] = s3;
165: x[3 + idx] = s4;
166: x[4 + idx] = s5;
167: x[5 + idx] = s6;
168: x[6 + idx] = s7;
169: }
171: /* backward solve the upper triangular */
172: for (i = n - 1; i >= 0; i--) {
173: v = aa + bs2 * (adiag[i + 1] + 1);
174: vi = aj + adiag[i + 1] + 1;
175: nz = adiag[i] - adiag[i + 1] - 1;
176: idt = bs * i;
177: s1 = x[idt];
178: s2 = x[1 + idt];
179: s3 = x[2 + idt];
180: s4 = x[3 + idt];
181: s5 = x[4 + idt];
182: s6 = x[5 + idt];
183: s7 = x[6 + idt];
184: for (k = 0; k < nz; k++) {
185: idx = bs * vi[k];
186: x1 = x[idx];
187: x2 = x[1 + idx];
188: x3 = x[2 + idx];
189: x4 = x[3 + idx];
190: x5 = x[4 + idx];
191: x6 = x[5 + idx];
192: x7 = x[6 + idx];
193: s1 -= v[0] * x1 + v[7] * x2 + v[14] * x3 + v[21] * x4 + v[28] * x5 + v[35] * x6 + v[42] * x7;
194: s2 -= v[1] * x1 + v[8] * x2 + v[15] * x3 + v[22] * x4 + v[29] * x5 + v[36] * x6 + v[43] * x7;
195: s3 -= v[2] * x1 + v[9] * x2 + v[16] * x3 + v[23] * x4 + v[30] * x5 + v[37] * x6 + v[44] * x7;
196: s4 -= v[3] * x1 + v[10] * x2 + v[17] * x3 + v[24] * x4 + v[31] * x5 + v[38] * x6 + v[45] * x7;
197: s5 -= v[4] * x1 + v[11] * x2 + v[18] * x3 + v[25] * x4 + v[32] * x5 + v[39] * x6 + v[46] * x7;
198: s6 -= v[5] * x1 + v[12] * x2 + v[19] * x3 + v[26] * x4 + v[33] * x5 + v[40] * x6 + v[47] * x7;
199: s7 -= v[6] * x1 + v[13] * x2 + v[20] * x3 + v[27] * x4 + v[34] * x5 + v[41] * x6 + v[48] * x7;
200: v += bs2;
201: }
202: /* x = inv_diagonal*x */
203: x[idt] = v[0] * s1 + v[7] * s2 + v[14] * s3 + v[21] * s4 + v[28] * s5 + v[35] * s6 + v[42] * s7;
204: x[1 + idt] = v[1] * s1 + v[8] * s2 + v[15] * s3 + v[22] * s4 + v[29] * s5 + v[36] * s6 + v[43] * s7;
205: x[2 + idt] = v[2] * s1 + v[9] * s2 + v[16] * s3 + v[23] * s4 + v[30] * s5 + v[37] * s6 + v[44] * s7;
206: x[3 + idt] = v[3] * s1 + v[10] * s2 + v[17] * s3 + v[24] * s4 + v[31] * s5 + v[38] * s6 + v[45] * s7;
207: x[4 + idt] = v[4] * s1 + v[11] * s2 + v[18] * s3 + v[25] * s4 + v[32] * s5 + v[39] * s6 + v[46] * s7;
208: x[5 + idt] = v[5] * s1 + v[12] * s2 + v[19] * s3 + v[26] * s4 + v[33] * s5 + v[40] * s6 + v[47] * s7;
209: x[6 + idt] = v[6] * s1 + v[13] * s2 + v[20] * s3 + v[27] * s4 + v[34] * s5 + v[41] * s6 + v[48] * s7;
210: }
212: VecRestoreArrayRead(bb, &b);
213: VecRestoreArray(xx, &x);
214: PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n);
215: return 0;
216: }