1 #include <../src/mat/impls/baij/seq/baij.h>
2
MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)3 PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
4 {
5 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
6 PetscErrorCode ierr;
7 const PetscInt *diag=a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j;
8 PetscInt i,nz,idx,idt,oidx;
9 const MatScalar *aa=a->a,*v;
10 PetscScalar s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6,*x;
11
12 PetscFunctionBegin;
13 ierr = VecCopy(bb,xx);CHKERRQ(ierr);
14 ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
15
16 /* forward solve the U^T */
17 idx = 0;
18 for (i=0; i<n; i++) {
19
20 v = aa + 36*diag[i];
21 /* multiply by the inverse of the block diagonal */
22 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; x4 = x[3+idx]; x5 = x[4+idx];
23 x6 = x[5+idx];
24 s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5 + v[5]*x6;
25 s2 = v[6]*x1 + v[7]*x2 + v[8]*x3 + v[9]*x4 + v[10]*x5 + v[11]*x6;
26 s3 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4 + v[16]*x5 + v[17]*x6;
27 s4 = v[18]*x1 + v[19]*x2 + v[20]*x3 + v[21]*x4 + v[22]*x5 + v[23]*x6;
28 s5 = v[24]*x1 + v[25]*x2 + v[26]*x3 + v[27]*x4 + v[28]*x5 + v[29]*x6;
29 s6 = v[30]*x1 + v[31]*x2 + v[32]*x3 + v[33]*x4 + v[34]*x5 + v[35]*x6;
30 v += 36;
31
32 vi = aj + diag[i] + 1;
33 nz = ai[i+1] - diag[i] - 1;
34 while (nz--) {
35 oidx = 6*(*vi++);
36 x[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6;
37 x[oidx+1] -= v[6]*s1 + v[7]*s2 + v[8]*s3 + v[9]*s4 + v[10]*s5 + v[11]*s6;
38 x[oidx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
39 x[oidx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
40 x[oidx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
41 x[oidx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
42 v += 36;
43 }
44 x[idx] = s1;x[1+idx] = s2; x[2+idx] = s3;x[3+idx] = s4; x[4+idx] = s5;
45 x[5+idx] = s6;
46 idx += 6;
47 }
48 /* backward solve the L^T */
49 for (i=n-1; i>=0; i--) {
50 v = aa + 36*diag[i] - 36;
51 vi = aj + diag[i] - 1;
52 nz = diag[i] - ai[i];
53 idt = 6*i;
54 s1 = x[idt]; s2 = x[1+idt]; s3 = x[2+idt];s4 = x[3+idt]; s5 = x[4+idt];
55 s6 = x[5+idt];
56 while (nz--) {
57 idx = 6*(*vi--);
58 x[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6;
59 x[idx+1] -= v[6]*s1 + v[7]*s2 + v[8]*s3 + v[9]*s4 + v[10]*s5 + v[11]*s6;
60 x[idx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
61 x[idx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
62 x[idx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
63 x[idx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
64 v -= 36;
65 }
66 }
67 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
68 ierr = PetscLogFlops(2.0*36*(a->nz) - 6.0*A->cmap->n);CHKERRQ(ierr);
69 PetscFunctionReturn(0);
70 }
71
MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat A,Vec bb,Vec xx)72 PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat A,Vec bb,Vec xx)
73 {
74 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
75 PetscErrorCode ierr;
76 const PetscInt n=a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
77 PetscInt nz,idx,idt,j,i,oidx;
78 const PetscInt bs =A->rmap->bs,bs2=a->bs2;
79 const MatScalar *aa=a->a,*v;
80 PetscScalar s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6,*x;
81
82 PetscFunctionBegin;
83 ierr = VecCopy(bb,xx);CHKERRQ(ierr);
84 ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
85
86 /* forward solve the U^T */
87 idx = 0;
88 for (i=0; i<n; i++) {
89 v = aa + bs2*diag[i];
90 /* multiply by the inverse of the block diagonal */
91 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; x4 = x[3+idx];
92 x5 = x[4+idx]; x6 = x[5+idx];
93 s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5 + v[5]*x6;
94 s2 = v[6]*x1 + v[7]*x2 + v[8]*x3 + v[9]*x4 + v[10]*x5 + v[11]*x6;
95 s3 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4 + v[16]*x5 + v[17]*x6;
96 s4 = v[18]*x1 + v[19]*x2 + v[20]*x3 + v[21]*x4 + v[22]*x5 + v[23]*x6;
97 s5 = v[24]*x1 + v[25]*x2 + v[26]*x3 + v[27]*x4 + v[28]*x5 + v[29]*x6;
98 s6 = v[30]*x1 + v[31]*x2 + v[32]*x3 + v[33]*x4 + v[34]*x5 + v[35]*x6;
99 v -= bs2;
100
101 vi = aj + diag[i] - 1;
102 nz = diag[i] - diag[i+1] - 1;
103 for (j=0; j>-nz; j--) {
104 oidx = bs*vi[j];
105 x[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6;
106 x[oidx+1] -= v[6]*s1 + v[7]*s2 + v[8]*s3 + v[9]*s4 + v[10]*s5 + v[11]*s6;
107 x[oidx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
108 x[oidx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
109 x[oidx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
110 x[oidx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
111 v -= bs2;
112 }
113 x[idx] = s1;x[1+idx] = s2; x[2+idx] = s3; x[3+idx] = s4; x[4+idx] = s5;
114 x[5+idx] = s6;
115 idx += bs;
116 }
117 /* backward solve the L^T */
118 for (i=n-1; i>=0; i--) {
119 v = aa + bs2*ai[i];
120 vi = aj + ai[i];
121 nz = ai[i+1] - ai[i];
122 idt = bs*i;
123 s1 = x[idt]; s2 = x[1+idt]; s3 = x[2+idt]; s4 = x[3+idt]; s5 = x[4+idt];
124 s6 = x[5+idt];
125 for (j=0; j<nz; j++) {
126 idx = bs*vi[j];
127 x[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6;
128 x[idx+1] -= v[6]*s1 + v[7]*s2 + v[8]*s3 + v[9]*s4 + v[10]*s5 + v[11]*s6;
129 x[idx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
130 x[idx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
131 x[idx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
132 x[idx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
133 v += bs2;
134 }
135 }
136 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
137 ierr = PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr);
138 PetscFunctionReturn(0);
139 }
140