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