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