Actual source code: baijsolvtran7.c
1: #include <../src/mat/impls/baij/seq/baij.h>
2: #include <petsc/private/kernels/blockinvert.h>
4: PetscErrorCode MatSolveTranspose_SeqBAIJ_7_inplace(Mat A,Vec bb,Vec xx)
5: {
6: Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data;
7: IS iscol=a->col,isrow=a->row;
8: const PetscInt *r,*c,*rout,*cout;
9: const PetscInt *diag=a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j;
10: PetscInt i,nz,idx,idt,ii,ic,ir,oidx;
11: const MatScalar *aa=a->a,*v;
12: PetscScalar s1,s2,s3,s4,s5,s6,s7,x1,x2,x3,x4,x5,x6,x7,*x,*t;
13: const PetscScalar *b;
15: VecGetArrayRead(bb,&b);
16: VecGetArray(xx,&x);
17: t = a->solve_work;
19: ISGetIndices(isrow,&rout); r = rout;
20: ISGetIndices(iscol,&cout); c = cout;
22: /* copy the b into temp work space according to permutation */
23: ii = 0;
24: for (i=0; i<n; i++) {
25: ic = 7*c[i];
26: t[ii] = b[ic];
27: t[ii+1] = b[ic+1];
28: t[ii+2] = b[ic+2];
29: t[ii+3] = b[ic+3];
30: t[ii+4] = b[ic+4];
31: t[ii+5] = b[ic+5];
32: t[ii+6] = b[ic+6];
33: ii += 7;
34: }
36: /* forward solve the U^T */
37: idx = 0;
38: for (i=0; i<n; i++) {
40: v = aa + 49*diag[i];
41: /* multiply by the inverse of the block diagonal */
42: x1 = t[idx]; x2 = t[1+idx]; x3 = t[2+idx]; x4 = t[3+idx]; x5 = t[4+idx];
43: x6 = t[5+idx]; x7 = t[6+idx];
44: s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5 + v[5]*x6 + v[6]*x7;
45: s2 = v[7]*x1 + v[8]*x2 + v[9]*x3 + v[10]*x4 + v[11]*x5 + v[12]*x6 + v[13]*x7;
46: s3 = v[14]*x1 + v[15]*x2 + v[16]*x3 + v[17]*x4 + v[18]*x5 + v[19]*x6 + v[20]*x7;
47: s4 = v[21]*x1 + v[22]*x2 + v[23]*x3 + v[24]*x4 + v[25]*x5 + v[26]*x6 + v[27]*x7;
48: s5 = v[28]*x1 + v[29]*x2 + v[30]*x3 + v[31]*x4 + v[32]*x5 + v[33]*x6 + v[34]*x7;
49: s6 = v[35]*x1 + v[36]*x2 + v[37]*x3 + v[38]*x4 + v[39]*x5 + v[40]*x6 + v[41]*x7;
50: s7 = v[42]*x1 + v[43]*x2 + v[44]*x3 + v[45]*x4 + v[46]*x5 + v[47]*x6 + v[48]*x7;
51: v += 49;
53: vi = aj + diag[i] + 1;
54: nz = ai[i+1] - diag[i] - 1;
55: while (nz--) {
56: oidx = 7*(*vi++);
57: t[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6 + v[6]*s7;
58: t[oidx+1] -= v[7]*s1 + v[8]*s2 + v[9]*s3 + v[10]*s4 + v[11]*s5 + v[12]*s6 + v[13]*s7;
59: t[oidx+2] -= v[14]*s1 + v[15]*s2 + v[16]*s3 + v[17]*s4 + v[18]*s5 + v[19]*s6 + v[20]*s7;
60: t[oidx+3] -= v[21]*s1 + v[22]*s2 + v[23]*s3 + v[24]*s4 + v[25]*s5 + v[26]*s6 + v[27]*s7;
61: t[oidx+4] -= v[28]*s1 + v[29]*s2 + v[30]*s3 + v[31]*s4 + v[32]*s5 + v[33]*s6 + v[34]*s7;
62: t[oidx+5] -= v[35]*s1 + v[36]*s2 + v[37]*s3 + v[38]*s4 + v[39]*s5 + v[40]*s6 + v[41]*s7;
63: t[oidx+6] -= v[42]*s1 + v[43]*s2 + v[44]*s3 + v[45]*s4 + v[46]*s5 + v[47]*s6 + v[48]*s7;
64: v += 49;
65: }
66: t[idx] = s1;t[1+idx] = s2; t[2+idx] = s3;t[3+idx] = s4; t[4+idx] = s5;
67: t[5+idx] = s6;t[6+idx] = s7;
68: idx += 7;
69: }
70: /* backward solve the L^T */
71: for (i=n-1; i>=0; i--) {
72: v = aa + 49*diag[i] - 49;
73: vi = aj + diag[i] - 1;
74: nz = diag[i] - ai[i];
75: idt = 7*i;
76: s1 = t[idt]; s2 = t[1+idt]; s3 = t[2+idt];s4 = t[3+idt]; s5 = t[4+idt];
77: s6 = t[5+idt];s7 = t[6+idt];
78: while (nz--) {
79: idx = 7*(*vi--);
80: t[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6 + v[6]*s7;
81: t[idx+1] -= v[7]*s1 + v[8]*s2 + v[9]*s3 + v[10]*s4 + v[11]*s5 + v[12]*s6 + v[13]*s7;
82: t[idx+2] -= v[14]*s1 + v[15]*s2 + v[16]*s3 + v[17]*s4 + v[18]*s5 + v[19]*s6 + v[20]*s7;
83: t[idx+3] -= v[21]*s1 + v[22]*s2 + v[23]*s3 + v[24]*s4 + v[25]*s5 + v[26]*s6 + v[27]*s7;
84: t[idx+4] -= v[28]*s1 + v[29]*s2 + v[30]*s3 + v[31]*s4 + v[32]*s5 + v[33]*s6 + v[34]*s7;
85: t[idx+5] -= v[35]*s1 + v[36]*s2 + v[37]*s3 + v[38]*s4 + v[39]*s5 + v[40]*s6 + v[41]*s7;
86: t[idx+6] -= v[42]*s1 + v[43]*s2 + v[44]*s3 + v[45]*s4 + v[46]*s5 + v[47]*s6 + v[48]*s7;
87: v -= 49;
88: }
89: }
91: /* copy t into x according to permutation */
92: ii = 0;
93: for (i=0; i<n; i++) {
94: ir = 7*r[i];
95: x[ir] = t[ii];
96: x[ir+1] = t[ii+1];
97: x[ir+2] = t[ii+2];
98: x[ir+3] = t[ii+3];
99: x[ir+4] = t[ii+4];
100: x[ir+5] = t[ii+5];
101: x[ir+6] = t[ii+6];
102: ii += 7;
103: }
105: ISRestoreIndices(isrow,&rout);
106: ISRestoreIndices(iscol,&cout);
107: VecRestoreArrayRead(bb,&b);
108: VecRestoreArray(xx,&x);
109: PetscLogFlops(2.0*49*(a->nz) - 7.0*A->cmap->n);
110: return 0;
111: }
112: PetscErrorCode MatSolveTranspose_SeqBAIJ_7(Mat A,Vec bb,Vec xx)
113: {
114: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
115: IS iscol=a->col,isrow=a->row;
116: const PetscInt n =a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
117: const PetscInt *r,*c,*rout,*cout;
118: PetscInt nz,idx,idt,j,i,oidx,ii,ic,ir;
119: const PetscInt bs =A->rmap->bs,bs2=a->bs2;
120: const MatScalar *aa=a->a,*v;
121: PetscScalar s1,s2,s3,s4,s5,s6,s7,x1,x2,x3,x4,x5,x6,x7,*x,*t;
122: const PetscScalar *b;
124: VecGetArrayRead(bb,&b);
125: VecGetArray(xx,&x);
126: t = a->solve_work;
128: ISGetIndices(isrow,&rout); r = rout;
129: ISGetIndices(iscol,&cout); c = cout;
131: /* copy b into temp work space according to permutation */
132: for (i=0; i<n; i++) {
133: ii = bs*i; ic = bs*c[i];
134: t[ii] = b[ic]; t[ii+1] = b[ic+1]; t[ii+2] = b[ic+2]; t[ii+3] = b[ic+3];
135: t[ii+4] = b[ic+4]; t[ii+5] = b[ic+5]; t[ii+6] = b[ic+6];
136: }
138: /* forward solve the U^T */
139: idx = 0;
140: for (i=0; i<n; i++) {
141: v = aa + bs2*diag[i];
142: /* multiply by the inverse of the block diagonal */
143: x1 = t[idx]; x2 = t[1+idx]; x3 = t[2+idx]; x4 = t[3+idx]; x5 = t[4+idx];
144: x6 = t[5+idx]; x7 = t[6+idx];
145: s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5 + v[5]*x6 + v[6]*x7;
146: s2 = v[7]*x1 + v[8]*x2 + v[9]*x3 + v[10]*x4 + v[11]*x5 + v[12]*x6 + v[13]*x7;
147: s3 = v[14]*x1 + v[15]*x2 + v[16]*x3 + v[17]*x4 + v[18]*x5 + v[19]*x6 + v[20]*x7;
148: s4 = v[21]*x1 + v[22]*x2 + v[23]*x3 + v[24]*x4 + v[25]*x5 + v[26]*x6 + v[27]*x7;
149: s5 = v[28]*x1 + v[29]*x2 + v[30]*x3 + v[31]*x4 + v[32]*x5 + v[33]*x6 + v[34]*x7;
150: s6 = v[35]*x1 + v[36]*x2 + v[37]*x3 + v[38]*x4 + v[39]*x5 + v[40]*x6 + v[41]*x7;
151: s7 = v[42]*x1 + v[43]*x2 + v[44]*x3 + v[45]*x4 + v[46]*x5 + v[47]*x6 + v[48]*x7;
152: v -= bs2;
154: vi = aj + diag[i] - 1;
155: nz = diag[i] - diag[i+1] - 1;
156: for (j=0; j>-nz; j--) {
157: oidx = bs*vi[j];
158: t[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6 + v[6]*s7;
159: t[oidx+1] -= v[7]*s1 + v[8]*s2 + v[9]*s3 + v[10]*s4 + v[11]*s5 + v[12]*s6 + v[13]*s7;
160: t[oidx+2] -= v[14]*s1 + v[15]*s2 + v[16]*s3 + v[17]*s4 + v[18]*s5 + v[19]*s6 + v[20]*s7;
161: t[oidx+3] -= v[21]*s1 + v[22]*s2 + v[23]*s3 + v[24]*s4 + v[25]*s5 + v[26]*s6 + v[27]*s7;
162: t[oidx+4] -= v[28]*s1 + v[29]*s2 + v[30]*s3 + v[31]*s4 + v[32]*s5 + v[33]*s6 + v[34]*s7;
163: t[oidx+5] -= v[35]*s1 + v[36]*s2 + v[37]*s3 + v[38]*s4 + v[39]*s5 + v[40]*s6 + v[41]*s7;
164: t[oidx+6] -= v[42]*s1 + v[43]*s2 + v[44]*s3 + v[45]*s4 + v[46]*s5 + v[47]*s6 + v[48]*s7;
165: v -= bs2;
166: }
167: t[idx] = s1;t[1+idx] = s2; t[2+idx] = s3; t[3+idx] = s4; t[4+idx] =s5;
168: t[5+idx] = s6; t[6+idx] = s7;
169: idx += bs;
170: }
171: /* backward solve the L^T */
172: for (i=n-1; i>=0; i--) {
173: v = aa + bs2*ai[i];
174: vi = aj + ai[i];
175: nz = ai[i+1] - ai[i];
176: idt = bs*i;
177: s1 = t[idt]; s2 = t[1+idt]; s3 = t[2+idt]; s4 = t[3+idt]; s5 = t[4+idt];
178: s6 = t[5+idt]; s7 = t[6+idt];
179: for (j=0; j<nz; j++) {
180: idx = bs*vi[j];
181: t[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6 + v[6]*s7;
182: t[idx+1] -= v[7]*s1 + v[8]*s2 + v[9]*s3 + v[10]*s4 + v[11]*s5 + v[12]*s6 + v[13]*s7;
183: t[idx+2] -= v[14]*s1 + v[15]*s2 + v[16]*s3 + v[17]*s4 + v[18]*s5 + v[19]*s6 + v[20]*s7;
184: t[idx+3] -= v[21]*s1 + v[22]*s2 + v[23]*s3 + v[24]*s4 + v[25]*s5 + v[26]*s6 + v[27]*s7;
185: t[idx+4] -= v[28]*s1 + v[29]*s2 + v[30]*s3 + v[31]*s4 + v[32]*s5 + v[33]*s6 + v[34]*s7;
186: t[idx+5] -= v[35]*s1 + v[36]*s2 + v[37]*s3 + v[38]*s4 + v[39]*s5 + v[40]*s6 + v[41]*s7;
187: t[idx+6] -= v[42]*s1 + v[43]*s2 + v[44]*s3 + v[45]*s4 + v[46]*s5 + v[47]*s6 + v[48]*s7;
188: v += bs2;
189: }
190: }
192: /* copy t into x according to permutation */
193: for (i=0; i<n; i++) {
194: ii = bs*i; ir = bs*r[i];
195: x[ir] = t[ii]; x[ir+1] = t[ii+1]; x[ir+2] = t[ii+2]; x[ir+3] = t[ii+3];
196: x[ir+4] = t[ii+4]; x[ir+5] = t[ii+5]; x[ir+6] = t[ii+6];
197: }
199: ISRestoreIndices(isrow,&rout);
200: ISRestoreIndices(iscol,&cout);
201: VecRestoreArrayRead(bb,&b);
202: VecRestoreArray(xx,&x);
203: PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);
204: return 0;
205: }