xref: /petsc/src/mat/impls/baij/seq/baijsolvtran2.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_2_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,x1,x2,*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      = 2*c[i];
28*2c733ed4SBarry Smith     t[ii]   = b[ic];
29*2c733ed4SBarry Smith     t[ii+1] = b[ic+1];
30*2c733ed4SBarry Smith     ii     += 2;
31*2c733ed4SBarry Smith   }
32*2c733ed4SBarry Smith 
33*2c733ed4SBarry Smith   /* forward solve the U^T */
34*2c733ed4SBarry Smith   idx = 0;
35*2c733ed4SBarry Smith   for (i=0; i<n; i++) {
36*2c733ed4SBarry Smith 
37*2c733ed4SBarry Smith     v = aa + 4*diag[i];
38*2c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
39*2c733ed4SBarry Smith     x1 = t[idx];   x2 = t[1+idx];
40*2c733ed4SBarry Smith     s1 = v[0]*x1  +  v[1]*x2;
41*2c733ed4SBarry Smith     s2 = v[2]*x1  +  v[3]*x2;
42*2c733ed4SBarry Smith     v += 4;
43*2c733ed4SBarry Smith 
44*2c733ed4SBarry Smith     vi = aj + diag[i] + 1;
45*2c733ed4SBarry Smith     nz = ai[i+1] - diag[i] - 1;
46*2c733ed4SBarry Smith     while (nz--) {
47*2c733ed4SBarry Smith       oidx       = 2*(*vi++);
48*2c733ed4SBarry Smith       t[oidx]   -= v[0]*s1  +  v[1]*s2;
49*2c733ed4SBarry Smith       t[oidx+1] -= v[2]*s1  +  v[3]*s2;
50*2c733ed4SBarry Smith       v         += 4;
51*2c733ed4SBarry Smith     }
52*2c733ed4SBarry Smith     t[idx] = s1;t[1+idx] = s2;
53*2c733ed4SBarry Smith     idx   += 2;
54*2c733ed4SBarry Smith   }
55*2c733ed4SBarry Smith   /* backward solve the L^T */
56*2c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
57*2c733ed4SBarry Smith     v   = aa + 4*diag[i] - 4;
58*2c733ed4SBarry Smith     vi  = aj + diag[i] - 1;
59*2c733ed4SBarry Smith     nz  = diag[i] - ai[i];
60*2c733ed4SBarry Smith     idt = 2*i;
61*2c733ed4SBarry Smith     s1  = t[idt];  s2 = t[1+idt];
62*2c733ed4SBarry Smith     while (nz--) {
63*2c733ed4SBarry Smith       idx       = 2*(*vi--);
64*2c733ed4SBarry Smith       t[idx]   -=  v[0]*s1 +  v[1]*s2;
65*2c733ed4SBarry Smith       t[idx+1] -=  v[2]*s1 +  v[3]*s2;
66*2c733ed4SBarry Smith       v        -= 4;
67*2c733ed4SBarry Smith     }
68*2c733ed4SBarry Smith   }
69*2c733ed4SBarry Smith 
70*2c733ed4SBarry Smith   /* copy t into x according to permutation */
71*2c733ed4SBarry Smith   ii = 0;
72*2c733ed4SBarry Smith   for (i=0; i<n; i++) {
73*2c733ed4SBarry Smith     ir      = 2*r[i];
74*2c733ed4SBarry Smith     x[ir]   = t[ii];
75*2c733ed4SBarry Smith     x[ir+1] = t[ii+1];
76*2c733ed4SBarry Smith     ii     += 2;
77*2c733ed4SBarry Smith   }
78*2c733ed4SBarry Smith 
79*2c733ed4SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
80*2c733ed4SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
81*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
82*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
83*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*4*(a->nz) - 2.0*A->cmap->n);CHKERRQ(ierr);
84*2c733ed4SBarry Smith   PetscFunctionReturn(0);
85*2c733ed4SBarry Smith }
86*2c733ed4SBarry Smith 
87*2c733ed4SBarry Smith PetscErrorCode MatSolveTranspose_SeqBAIJ_2(Mat A,Vec bb,Vec xx)
88*2c733ed4SBarry Smith {
89*2c733ed4SBarry Smith   Mat_SeqBAIJ       *a=(Mat_SeqBAIJ*)A->data;
90*2c733ed4SBarry Smith   PetscErrorCode    ierr;
91*2c733ed4SBarry Smith   IS                iscol=a->col,isrow=a->row;
92*2c733ed4SBarry Smith   const PetscInt    n    =a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
93*2c733ed4SBarry Smith   const PetscInt    *r,*c,*rout,*cout;
94*2c733ed4SBarry Smith   PetscInt          nz,idx,idt,j,i,oidx,ii,ic,ir;
95*2c733ed4SBarry Smith   const PetscInt    bs =A->rmap->bs,bs2=a->bs2;
96*2c733ed4SBarry Smith   const MatScalar   *aa=a->a,*v;
97*2c733ed4SBarry Smith   PetscScalar       s1,s2,x1,x2,*x,*t;
98*2c733ed4SBarry Smith   const PetscScalar *b;
99*2c733ed4SBarry Smith 
100*2c733ed4SBarry Smith   PetscFunctionBegin;
101*2c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
102*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
103*2c733ed4SBarry Smith   t    = a->solve_work;
104*2c733ed4SBarry Smith 
105*2c733ed4SBarry Smith   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
106*2c733ed4SBarry Smith   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
107*2c733ed4SBarry Smith 
108*2c733ed4SBarry Smith   /* copy b into temp work space according to permutation */
109*2c733ed4SBarry Smith   for (i=0; i<n; i++) {
110*2c733ed4SBarry Smith     ii    = bs*i; ic = bs*c[i];
111*2c733ed4SBarry Smith     t[ii] = b[ic]; t[ii+1] = b[ic+1];
112*2c733ed4SBarry Smith   }
113*2c733ed4SBarry Smith 
114*2c733ed4SBarry Smith   /* forward solve the U^T */
115*2c733ed4SBarry Smith   idx = 0;
116*2c733ed4SBarry Smith   for (i=0; i<n; i++) {
117*2c733ed4SBarry Smith     v = aa + bs2*diag[i];
118*2c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
119*2c733ed4SBarry Smith     x1 = t[idx];   x2 = t[1+idx];
120*2c733ed4SBarry Smith     s1 = v[0]*x1  +  v[1]*x2;
121*2c733ed4SBarry Smith     s2 = v[2]*x1  +  v[3]*x2;
122*2c733ed4SBarry Smith     v -= bs2;
123*2c733ed4SBarry Smith 
124*2c733ed4SBarry Smith     vi = aj + diag[i] - 1;
125*2c733ed4SBarry Smith     nz = diag[i] - diag[i+1] - 1;
126*2c733ed4SBarry Smith     for (j=0; j>-nz; j--) {
127*2c733ed4SBarry Smith       oidx       = bs*vi[j];
128*2c733ed4SBarry Smith       t[oidx]   -= v[0]*s1  +  v[1]*s2;
129*2c733ed4SBarry Smith       t[oidx+1] -= v[2]*s1  +  v[3]*s2;
130*2c733ed4SBarry Smith       v         -= bs2;
131*2c733ed4SBarry Smith     }
132*2c733ed4SBarry Smith     t[idx] = s1;t[1+idx] = s2;
133*2c733ed4SBarry Smith     idx   += bs;
134*2c733ed4SBarry Smith   }
135*2c733ed4SBarry Smith   /* backward solve the L^T */
136*2c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
137*2c733ed4SBarry Smith     v   = aa + bs2*ai[i];
138*2c733ed4SBarry Smith     vi  = aj + ai[i];
139*2c733ed4SBarry Smith     nz  = ai[i+1] - ai[i];
140*2c733ed4SBarry Smith     idt = bs*i;
141*2c733ed4SBarry Smith     s1  = t[idt];  s2 = t[1+idt];
142*2c733ed4SBarry Smith     for (j=0; j<nz; j++) {
143*2c733ed4SBarry Smith       idx       = bs*vi[j];
144*2c733ed4SBarry Smith       t[idx]   -=  v[0]*s1 +  v[1]*s2;
145*2c733ed4SBarry Smith       t[idx+1] -=  v[2]*s1 +  v[3]*s2;
146*2c733ed4SBarry Smith       v        += bs2;
147*2c733ed4SBarry Smith     }
148*2c733ed4SBarry Smith   }
149*2c733ed4SBarry Smith 
150*2c733ed4SBarry Smith   /* copy t into x according to permutation */
151*2c733ed4SBarry Smith   for (i=0; i<n; i++) {
152*2c733ed4SBarry Smith     ii    = bs*i;  ir = bs*r[i];
153*2c733ed4SBarry Smith     x[ir] = t[ii];  x[ir+1] = t[ii+1];
154*2c733ed4SBarry Smith   }
155*2c733ed4SBarry Smith 
156*2c733ed4SBarry Smith   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
157*2c733ed4SBarry Smith   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
158*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
159*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
160*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr);
161*2c733ed4SBarry Smith   PetscFunctionReturn(0);
162*2c733ed4SBarry Smith }
163