xref: /petsc/src/mat/impls/baij/seq/baijsolvtran6.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
1 #include <../src/mat/impls/baij/seq/baij.h>
2 #include <petsc/private/kernels/blockinvert.h>
3 
4 PetscErrorCode MatSolveTranspose_SeqBAIJ_6_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,s5,s6,x1,x2,x3,x4,x5,x6,*x,*t;
13   const PetscScalar *b;
14 
15   PetscFunctionBegin;
16   PetscCall(VecGetArrayRead(bb,&b));
17   PetscCall(VecGetArray(xx,&x));
18   t    = a->solve_work;
19 
20   PetscCall(ISGetIndices(isrow,&rout)); r = rout;
21   PetscCall(ISGetIndices(iscol,&cout)); c = cout;
22 
23   /* copy the b into temp work space according to permutation */
24   ii = 0;
25   for (i=0; i<n; i++) {
26     ic      = 6*c[i];
27     t[ii]   = b[ic];
28     t[ii+1] = b[ic+1];
29     t[ii+2] = b[ic+2];
30     t[ii+3] = b[ic+3];
31     t[ii+4] = b[ic+4];
32     t[ii+5] = b[ic+5];
33     ii     += 6;
34   }
35 
36   /* forward solve the U^T */
37   idx = 0;
38   for (i=0; i<n; i++) {
39 
40     v = aa + 36*diag[i];
41     /* multiply by the inverse of the block diagonal */
42     x1 = t[idx];   x2 = t[1+idx]; x3    = t[2+idx]; x4 = t[3+idx]; x5 = t[4+idx];
43     x6 = t[5+idx];
44     s1 = v[0]*x1  +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5 +  v[5]*x6;
45     s2 = v[6]*x1  +  v[7]*x2 +  v[8]*x3 +  v[9]*x4 + v[10]*x5 + v[11]*x6;
46     s3 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4 + v[16]*x5 + v[17]*x6;
47     s4 = v[18]*x1 + v[19]*x2 + v[20]*x3 + v[21]*x4 + v[22]*x5 + v[23]*x6;
48     s5 = v[24]*x1 + v[25]*x2 + v[26]*x3 + v[27]*x4 + v[28]*x5 + v[29]*x6;
49     s6 = v[30]*x1 + v[31]*x2 + v[32]*x3 + v[33]*x4 + v[34]*x5 + v[35]*x6;
50     v += 36;
51 
52     vi = aj + diag[i] + 1;
53     nz = ai[i+1] - diag[i] - 1;
54     while (nz--) {
55       oidx       = 6*(*vi++);
56       t[oidx]   -= v[0]*s1  +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5 +  v[5]*s6;
57       t[oidx+1] -= v[6]*s1  +  v[7]*s2 +  v[8]*s3 +  v[9]*s4 + v[10]*s5 + v[11]*s6;
58       t[oidx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
59       t[oidx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
60       t[oidx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
61       t[oidx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
62       v         += 36;
63     }
64     t[idx]   = s1;t[1+idx] = s2; t[2+idx] = s3;t[3+idx] = s4; t[4+idx] = s5;
65     t[5+idx] = s6;
66     idx     += 6;
67   }
68   /* backward solve the L^T */
69   for (i=n-1; i>=0; i--) {
70     v   = aa + 36*diag[i] - 36;
71     vi  = aj + diag[i] - 1;
72     nz  = diag[i] - ai[i];
73     idt = 6*i;
74     s1  = t[idt];  s2 = t[1+idt]; s3 = t[2+idt];s4 = t[3+idt]; s5 = t[4+idt];
75     s6  = t[5+idt];
76     while (nz--) {
77       idx       = 6*(*vi--);
78       t[idx]   -=  v[0]*s1 +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5 +  v[5]*s6;
79       t[idx+1] -=  v[6]*s1 +  v[7]*s2 +  v[8]*s3 +  v[9]*s4 + v[10]*s5 + v[11]*s6;
80       t[idx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
81       t[idx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
82       t[idx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
83       t[idx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
84       v        -= 36;
85     }
86   }
87 
88   /* copy t into x according to permutation */
89   ii = 0;
90   for (i=0; i<n; i++) {
91     ir      = 6*r[i];
92     x[ir]   = t[ii];
93     x[ir+1] = t[ii+1];
94     x[ir+2] = t[ii+2];
95     x[ir+3] = t[ii+3];
96     x[ir+4] = t[ii+4];
97     x[ir+5] = t[ii+5];
98     ii     += 6;
99   }
100 
101   PetscCall(ISRestoreIndices(isrow,&rout));
102   PetscCall(ISRestoreIndices(iscol,&cout));
103   PetscCall(VecRestoreArrayRead(bb,&b));
104   PetscCall(VecRestoreArray(xx,&x));
105   PetscCall(PetscLogFlops(2.0*36*(a->nz) - 6.0*A->cmap->n));
106   PetscFunctionReturn(0);
107 }
108 
109 PetscErrorCode MatSolveTranspose_SeqBAIJ_6(Mat A,Vec bb,Vec xx)
110 {
111   Mat_SeqBAIJ       *a=(Mat_SeqBAIJ*)A->data;
112   IS                iscol=a->col,isrow=a->row;
113   const PetscInt    n    =a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
114   const PetscInt    *r,*c,*rout,*cout;
115   PetscInt          nz,idx,idt,j,i,oidx,ii,ic,ir;
116   const PetscInt    bs =A->rmap->bs,bs2=a->bs2;
117   const MatScalar   *aa=a->a,*v;
118   PetscScalar       s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6,*x,*t;
119   const PetscScalar *b;
120 
121   PetscFunctionBegin;
122   PetscCall(VecGetArrayRead(bb,&b));
123   PetscCall(VecGetArray(xx,&x));
124   t    = a->solve_work;
125 
126   PetscCall(ISGetIndices(isrow,&rout)); r = rout;
127   PetscCall(ISGetIndices(iscol,&cout)); c = cout;
128 
129   /* copy b into temp work space according to permutation */
130   for (i=0; i<n; i++) {
131     ii      = bs*i; ic = bs*c[i];
132     t[ii]   = b[ic]; t[ii+1] = b[ic+1]; t[ii+2] = b[ic+2]; t[ii+3] = b[ic+3];
133     t[ii+4] = b[ic+4];  t[ii+5] = b[ic+5];
134   }
135 
136   /* forward solve the U^T */
137   idx = 0;
138   for (i=0; i<n; i++) {
139     v = aa + bs2*diag[i];
140     /* multiply by the inverse of the block diagonal */
141     x1 = t[idx];   x2 = t[1+idx]; x3    = t[2+idx]; x4 = t[3+idx]; x5 = t[4+idx];
142     x6 = t[5+idx];
143     s1 = v[0]*x1  +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5 +  v[5]*x6;
144     s2 = v[6]*x1  +  v[7]*x2 +  v[8]*x3 +  v[9]*x4 + v[10]*x5 + v[11]*x6;
145     s3 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4 + v[16]*x5 + v[17]*x6;
146     s4 = v[18]*x1 + v[19]*x2 + v[20]*x3 + v[21]*x4 + v[22]*x5 + v[23]*x6;
147     s5 = v[24]*x1 + v[25]*x2 + v[26]*x3 + v[27]*x4 + v[28]*x5 + v[29]*x6;
148     s6 = v[30]*x1 + v[31]*x2 + v[32]*x3 + v[33]*x4 + v[34]*x5 + v[35]*x6;
149     v -= bs2;
150 
151     vi = aj + diag[i] - 1;
152     nz = diag[i] - diag[i+1] - 1;
153     for (j=0; j>-nz; j--) {
154       oidx       = bs*vi[j];
155       t[oidx]   -= v[0]*s1  +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5 +  v[5]*s6;
156       t[oidx+1] -= v[6]*s1  +  v[7]*s2 +  v[8]*s3 +  v[9]*s4 + v[10]*s5 + v[11]*s6;
157       t[oidx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
158       t[oidx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
159       t[oidx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
160       t[oidx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
161       v         -= bs2;
162     }
163     t[idx]   = s1;t[1+idx] = s2;  t[2+idx] = s3;  t[3+idx] = s4; t[4+idx] =s5;
164     t[5+idx] = s6;
165     idx     += bs;
166   }
167   /* backward solve the L^T */
168   for (i=n-1; i>=0; i--) {
169     v   = aa + bs2*ai[i];
170     vi  = aj + ai[i];
171     nz  = ai[i+1] - ai[i];
172     idt = bs*i;
173     s1  = t[idt];  s2 = t[1+idt];  s3 = t[2+idt];  s4 = t[3+idt]; s5 = t[4+idt];
174     s6  = t[5+idt];
175     for (j=0; j<nz; j++) {
176       idx       = bs*vi[j];
177       t[idx]   -= v[0]*s1  +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5 +  v[5]*s6;
178       t[idx+1] -= v[6]*s1  +  v[7]*s2 +  v[8]*s3 +  v[9]*s4 + v[10]*s5 + v[11]*s6;
179       t[idx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
180       t[idx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
181       t[idx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
182       t[idx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
183       v        += bs2;
184     }
185   }
186 
187   /* copy t into x according to permutation */
188   for (i=0; i<n; i++) {
189     ii      = bs*i;  ir = bs*r[i];
190     x[ir]   = t[ii];  x[ir+1] = t[ii+1]; x[ir+2] = t[ii+2];  x[ir+3] = t[ii+3];
191     x[ir+4] = t[ii+4];  x[ir+5] = t[ii+5];
192   }
193 
194   PetscCall(ISRestoreIndices(isrow,&rout));
195   PetscCall(ISRestoreIndices(iscol,&cout));
196   PetscCall(VecRestoreArrayRead(bb,&b));
197   PetscCall(VecRestoreArray(xx,&x));
198   PetscCall(PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n));
199   PetscFunctionReturn(0);
200 }
201