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