Actual source code: baijsolvtran4.c
1: #include <../src/mat/impls/baij/seq/baij.h>
2: #include <petsc/private/kernels/blockinvert.h>
4: PetscErrorCode MatSolveTranspose_SeqBAIJ_4_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,x1,x2,x3,x4,*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 = 4*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: ii += 4;
31: }
33: /* forward solve the U^T */
34: idx = 0;
35: for (i=0; i<n; i++) {
37: v = aa + 16*diag[i];
38: /* multiply by the inverse of the block diagonal */
39: x1 = t[idx]; x2 = t[1+idx]; x3 = t[2+idx]; x4 = t[3+idx];
40: s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
41: s2 = v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
42: s3 = v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
43: s4 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
44: v += 16;
46: vi = aj + diag[i] + 1;
47: nz = ai[i+1] - diag[i] - 1;
48: while (nz--) {
49: oidx = 4*(*vi++);
50: t[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4;
51: t[oidx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4;
52: t[oidx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4;
53: t[oidx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4;
54: v += 16;
55: }
56: t[idx] = s1;t[1+idx] = s2; t[2+idx] = s3;t[3+idx] = s4;
57: idx += 4;
58: }
59: /* backward solve the L^T */
60: for (i=n-1; i>=0; i--) {
61: v = aa + 16*diag[i] - 16;
62: vi = aj + diag[i] - 1;
63: nz = diag[i] - ai[i];
64: idt = 4*i;
65: s1 = t[idt]; s2 = t[1+idt]; s3 = t[2+idt];s4 = t[3+idt];
66: while (nz--) {
67: idx = 4*(*vi--);
68: t[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4;
69: t[idx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4;
70: t[idx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4;
71: t[idx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4;
72: v -= 16;
73: }
74: }
76: /* copy t into x according to permutation */
77: ii = 0;
78: for (i=0; i<n; i++) {
79: ir = 4*r[i];
80: x[ir] = t[ii];
81: x[ir+1] = t[ii+1];
82: x[ir+2] = t[ii+2];
83: x[ir+3] = t[ii+3];
84: ii += 4;
85: }
87: ISRestoreIndices(isrow,&rout);
88: ISRestoreIndices(iscol,&cout);
89: VecRestoreArrayRead(bb,&b);
90: VecRestoreArray(xx,&x);
91: PetscLogFlops(2.0*16*(a->nz) - 4.0*A->cmap->n);
92: return 0;
93: }
95: PetscErrorCode MatSolveTranspose_SeqBAIJ_4(Mat A,Vec bb,Vec xx)
96: {
97: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
98: IS iscol=a->col,isrow=a->row;
99: const PetscInt n =a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
100: const PetscInt *r,*c,*rout,*cout;
101: PetscInt nz,idx,idt,j,i,oidx,ii,ic,ir;
102: const PetscInt bs =A->rmap->bs,bs2=a->bs2;
103: const MatScalar *aa=a->a,*v;
104: PetscScalar s1,s2,s3,s4,x1,x2,x3,x4,*x,*t;
105: const PetscScalar *b;
107: VecGetArrayRead(bb,&b);
108: VecGetArray(xx,&x);
109: t = a->solve_work;
111: ISGetIndices(isrow,&rout); r = rout;
112: ISGetIndices(iscol,&cout); c = cout;
114: /* copy b into temp work space according to permutation */
115: for (i=0; i<n; i++) {
116: ii = bs*i; ic = bs*c[i];
117: t[ii] = b[ic]; t[ii+1] = b[ic+1]; t[ii+2] = b[ic+2]; t[ii+3] = b[ic+3];
118: }
120: /* forward solve the U^T */
121: idx = 0;
122: for (i=0; i<n; i++) {
123: v = aa + bs2*diag[i];
124: /* multiply by the inverse of the block diagonal */
125: x1 = t[idx]; x2 = t[1+idx]; x3 = t[2+idx]; x4 = t[3+idx];
126: s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
127: s2 = v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
128: s3 = v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
129: s4 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
130: v -= bs2;
132: vi = aj + diag[i] - 1;
133: nz = diag[i] - diag[i+1] - 1;
134: for (j=0; j>-nz; j--) {
135: oidx = bs*vi[j];
136: t[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4;
137: t[oidx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4;
138: t[oidx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4;
139: t[oidx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4;
140: v -= bs2;
141: }
142: t[idx] = s1;t[1+idx] = s2; t[2+idx] = s3; t[3+idx] = s4;
143: idx += bs;
144: }
145: /* backward solve the L^T */
146: for (i=n-1; i>=0; i--) {
147: v = aa + bs2*ai[i];
148: vi = aj + ai[i];
149: nz = ai[i+1] - ai[i];
150: idt = bs*i;
151: s1 = t[idt]; s2 = t[1+idt]; s3 = t[2+idt]; s4 = t[3+idt];
152: for (j=0; j<nz; j++) {
153: idx = bs*vi[j];
154: t[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4;
155: t[idx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4;
156: t[idx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4;
157: t[idx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4;
158: v += bs2;
159: }
160: }
162: /* copy t into x according to permutation */
163: for (i=0; i<n; i++) {
164: ii = bs*i; ir = bs*r[i];
165: x[ir] = t[ii]; x[ir+1] = t[ii+1]; x[ir+2] = t[ii+2]; x[ir+3] = t[ii+3];
166: }
168: ISRestoreIndices(isrow,&rout);
169: ISRestoreIndices(iscol,&cout);
170: VecRestoreArrayRead(bb,&b);
171: VecRestoreArray(xx,&x);
172: PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);
173: return 0;
174: }