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