Actual source code: baijsolvtrannat3.c
1: #include <../src/mat/impls/baij/seq/baij.h>
3: PetscErrorCode MatSolveTranspose_SeqBAIJ_3_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
4: {
5: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
6: const PetscInt n=a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
7: PetscInt i,nz,idx,idt,oidx;
8: const MatScalar *aa=a->a,*v;
9: PetscScalar s1,s2,s3,x1,x2,x3,*x;
11: VecCopy(bb,xx);
12: VecGetArray(xx,&x);
14: /* forward solve the U^T */
15: idx = 0;
16: for (i=0; i<n; i++) {
18: v = aa + 9*diag[i];
19: /* multiply by the inverse of the block diagonal */
20: x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx];
21: s1 = v[0]*x1 + v[1]*x2 + v[2]*x3;
22: s2 = v[3]*x1 + v[4]*x2 + v[5]*x3;
23: s3 = v[6]*x1 + v[7]*x2 + v[8]*x3;
24: v += 9;
26: vi = aj + diag[i] + 1;
27: nz = ai[i+1] - diag[i] - 1;
28: while (nz--) {
29: oidx = 3*(*vi++);
30: x[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3;
31: x[oidx+1] -= v[3]*s1 + v[4]*s2 + v[5]*s3;
32: x[oidx+2] -= v[6]*s1 + v[7]*s2 + v[8]*s3;
33: v += 9;
34: }
35: x[idx] = s1;x[1+idx] = s2; x[2+idx] = s3;
36: idx += 3;
37: }
38: /* backward solve the L^T */
39: for (i=n-1; i>=0; i--) {
40: v = aa + 9*diag[i] - 9;
41: vi = aj + diag[i] - 1;
42: nz = diag[i] - ai[i];
43: idt = 3*i;
44: s1 = x[idt]; s2 = x[1+idt]; s3 = x[2+idt];
45: while (nz--) {
46: idx = 3*(*vi--);
47: x[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3;
48: x[idx+1] -= v[3]*s1 + v[4]*s2 + v[5]*s3;
49: x[idx+2] -= v[6]*s1 + v[7]*s2 + v[8]*s3;
50: v -= 9;
51: }
52: }
53: VecRestoreArray(xx,&x);
54: PetscLogFlops(2.0*9.0*(a->nz) - 3.0*A->cmap->n);
55: return 0;
56: }
58: PetscErrorCode MatSolveTranspose_SeqBAIJ_3_NaturalOrdering(Mat A,Vec bb,Vec xx)
59: {
60: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
61: const PetscInt n=a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
62: PetscInt nz,idx,idt,j,i,oidx;
63: const PetscInt bs =A->rmap->bs,bs2=a->bs2;
64: const MatScalar *aa=a->a,*v;
65: PetscScalar s1,s2,s3,x1,x2,x3,*x;
67: VecCopy(bb,xx);
68: VecGetArray(xx,&x);
70: /* forward solve the U^T */
71: idx = 0;
72: for (i=0; i<n; i++) {
73: v = aa + bs2*diag[i];
74: /* multiply by the inverse of the block diagonal */
75: x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx];
76: s1 = v[0]*x1 + v[1]*x2 + v[2]*x3;
77: s2 = v[3]*x1 + v[4]*x2 + v[5]*x3;
78: s3 = v[6]*x1 + v[7]*x2 + v[8]*x3;
79: v -= bs2;
81: vi = aj + diag[i] - 1;
82: nz = diag[i] - diag[i+1] - 1;
83: for (j=0; j>-nz; j--) {
84: oidx = bs*vi[j];
85: x[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3;
86: x[oidx+1] -= v[3]*s1 + v[4]*s2 + v[5]*s3;
87: x[oidx+2] -= v[6]*s1 + v[7]*s2 + v[8]*s3;
88: v -= bs2;
89: }
90: x[idx] = s1;x[1+idx] = s2; x[2+idx] = s3;
91: idx += bs;
92: }
93: /* backward solve the L^T */
94: for (i=n-1; i>=0; i--) {
95: v = aa + bs2*ai[i];
96: vi = aj + ai[i];
97: nz = ai[i+1] - ai[i];
98: idt = bs*i;
99: s1 = x[idt]; s2 = x[1+idt]; s3 = x[2+idt];
100: for (j=0; j<nz; j++) {
101: idx = bs*vi[j];
102: x[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3;
103: x[idx+1] -= v[3]*s1 + v[4]*s2 + v[5]*s3;
104: x[idx+2] -= v[6]*s1 + v[7]*s2 + v[8]*s3;
105: v += bs2;
106: }
107: }
108: VecRestoreArray(xx,&x);
109: PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);
110: return 0;
111: }