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