xref: /petsc/src/mat/impls/baij/seq/baijsolvnat3.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 /*
5*2c733ed4SBarry Smith       Special case where the matrix was ILU(0) factored in the natural
6*2c733ed4SBarry Smith    ordering. This eliminates the need for the column and row permutation.
7*2c733ed4SBarry Smith */
8*2c733ed4SBarry Smith PetscErrorCode MatSolve_SeqBAIJ_3_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
9*2c733ed4SBarry Smith {
10*2c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
11*2c733ed4SBarry Smith   const PetscInt    n  =a->mbs,*ai=a->i,*aj=a->j;
12*2c733ed4SBarry Smith   PetscErrorCode    ierr;
13*2c733ed4SBarry Smith   const PetscInt    *diag = a->diag,*vi;
14*2c733ed4SBarry Smith   const MatScalar   *aa   =a->a,*v;
15*2c733ed4SBarry Smith   PetscScalar       *x,s1,s2,s3,x1,x2,x3;
16*2c733ed4SBarry Smith   const PetscScalar *b;
17*2c733ed4SBarry Smith   PetscInt          jdx,idt,idx,nz,i;
18*2c733ed4SBarry Smith 
19*2c733ed4SBarry Smith   PetscFunctionBegin;
20*2c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
21*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
22*2c733ed4SBarry Smith 
23*2c733ed4SBarry Smith   /* forward solve the lower triangular */
24*2c733ed4SBarry Smith   idx  = 0;
25*2c733ed4SBarry Smith   x[0] = b[0]; x[1] = b[1]; x[2] = b[2];
26*2c733ed4SBarry Smith   for (i=1; i<n; i++) {
27*2c733ed4SBarry Smith     v    =  aa      + 9*ai[i];
28*2c733ed4SBarry Smith     vi   =  aj      + ai[i];
29*2c733ed4SBarry Smith     nz   =  diag[i] - ai[i];
30*2c733ed4SBarry Smith     idx +=  3;
31*2c733ed4SBarry Smith     s1   =  b[idx];s2 = b[1+idx];s3 = b[2+idx];
32*2c733ed4SBarry Smith     while (nz--) {
33*2c733ed4SBarry Smith       jdx = 3*(*vi++);
34*2c733ed4SBarry Smith       x1  = x[jdx];x2 = x[1+jdx];x3 = x[2+jdx];
35*2c733ed4SBarry Smith       s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3;
36*2c733ed4SBarry Smith       s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3;
37*2c733ed4SBarry Smith       s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3;
38*2c733ed4SBarry Smith       v  += 9;
39*2c733ed4SBarry Smith     }
40*2c733ed4SBarry Smith     x[idx]   = s1;
41*2c733ed4SBarry Smith     x[1+idx] = s2;
42*2c733ed4SBarry Smith     x[2+idx] = s3;
43*2c733ed4SBarry Smith   }
44*2c733ed4SBarry Smith   /* backward solve the upper triangular */
45*2c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
46*2c733ed4SBarry Smith     v   = aa + 9*diag[i] + 9;
47*2c733ed4SBarry Smith     vi  = aj + diag[i] + 1;
48*2c733ed4SBarry Smith     nz  = ai[i+1] - diag[i] - 1;
49*2c733ed4SBarry Smith     idt = 3*i;
50*2c733ed4SBarry Smith     s1  = x[idt];  s2 = x[1+idt];
51*2c733ed4SBarry Smith     s3  = x[2+idt];
52*2c733ed4SBarry Smith     while (nz--) {
53*2c733ed4SBarry Smith       idx = 3*(*vi++);
54*2c733ed4SBarry Smith       x1  = x[idx];   x2 = x[1+idx];x3    = x[2+idx];
55*2c733ed4SBarry Smith       s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3;
56*2c733ed4SBarry Smith       s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3;
57*2c733ed4SBarry Smith       s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3;
58*2c733ed4SBarry Smith       v  += 9;
59*2c733ed4SBarry Smith     }
60*2c733ed4SBarry Smith     v        = aa +  9*diag[i];
61*2c733ed4SBarry Smith     x[idt]   = v[0]*s1 + v[3]*s2 + v[6]*s3;
62*2c733ed4SBarry Smith     x[1+idt] = v[1]*s1 + v[4]*s2 + v[7]*s3;
63*2c733ed4SBarry Smith     x[2+idt] = v[2]*s1 + v[5]*s2 + v[8]*s3;
64*2c733ed4SBarry Smith   }
65*2c733ed4SBarry Smith 
66*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
67*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
68*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*9*(a->nz) - 3.0*A->cmap->n);CHKERRQ(ierr);
69*2c733ed4SBarry Smith   PetscFunctionReturn(0);
70*2c733ed4SBarry Smith }
71*2c733ed4SBarry Smith 
72*2c733ed4SBarry Smith PetscErrorCode MatSolve_SeqBAIJ_3_NaturalOrdering(Mat A,Vec bb,Vec xx)
73*2c733ed4SBarry Smith {
74*2c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
75*2c733ed4SBarry Smith   const PetscInt    n  =a->mbs,*vi,*ai=a->i,*aj=a->j,*adiag=a->diag;
76*2c733ed4SBarry Smith   PetscErrorCode    ierr;
77*2c733ed4SBarry Smith   PetscInt          i,k,nz,idx,jdx,idt;
78*2c733ed4SBarry Smith   const PetscInt    bs = A->rmap->bs,bs2 = a->bs2;
79*2c733ed4SBarry Smith   const MatScalar   *aa=a->a,*v;
80*2c733ed4SBarry Smith   PetscScalar       *x;
81*2c733ed4SBarry Smith   const PetscScalar *b;
82*2c733ed4SBarry Smith   PetscScalar       s1,s2,s3,x1,x2,x3;
83*2c733ed4SBarry Smith 
84*2c733ed4SBarry Smith   PetscFunctionBegin;
85*2c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
86*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
87*2c733ed4SBarry Smith   /* forward solve the lower triangular */
88*2c733ed4SBarry Smith   idx  = 0;
89*2c733ed4SBarry Smith   x[0] = b[idx]; x[1] = b[1+idx];x[2] = b[2+idx];
90*2c733ed4SBarry Smith   for (i=1; i<n; i++) {
91*2c733ed4SBarry Smith     v   = aa + bs2*ai[i];
92*2c733ed4SBarry Smith     vi  = aj + ai[i];
93*2c733ed4SBarry Smith     nz  = ai[i+1] - ai[i];
94*2c733ed4SBarry Smith     idx = bs*i;
95*2c733ed4SBarry Smith     s1  = b[idx];s2 = b[1+idx];s3 = b[2+idx];
96*2c733ed4SBarry Smith     for (k=0; k<nz; k++) {
97*2c733ed4SBarry Smith       jdx = bs*vi[k];
98*2c733ed4SBarry Smith       x1  = x[jdx];x2 = x[1+jdx]; x3 =x[2+jdx];
99*2c733ed4SBarry Smith       s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3;
100*2c733ed4SBarry Smith       s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3;
101*2c733ed4SBarry Smith       s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3;
102*2c733ed4SBarry Smith 
103*2c733ed4SBarry Smith       v +=  bs2;
104*2c733ed4SBarry Smith     }
105*2c733ed4SBarry Smith 
106*2c733ed4SBarry Smith     x[idx]   = s1;
107*2c733ed4SBarry Smith     x[1+idx] = s2;
108*2c733ed4SBarry Smith     x[2+idx] = s3;
109*2c733ed4SBarry Smith   }
110*2c733ed4SBarry Smith 
111*2c733ed4SBarry Smith   /* backward solve the upper triangular */
112*2c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
113*2c733ed4SBarry Smith     v   = aa + bs2*(adiag[i+1]+1);
114*2c733ed4SBarry Smith     vi  = aj + adiag[i+1]+1;
115*2c733ed4SBarry Smith     nz  = adiag[i] - adiag[i+1]-1;
116*2c733ed4SBarry Smith     idt = bs*i;
117*2c733ed4SBarry Smith     s1  = x[idt];  s2 = x[1+idt];s3 = x[2+idt];
118*2c733ed4SBarry Smith 
119*2c733ed4SBarry Smith     for (k=0; k<nz; k++) {
120*2c733ed4SBarry Smith       idx = bs*vi[k];
121*2c733ed4SBarry Smith       x1  = x[idx];   x2 = x[1+idx]; x3 = x[2+idx];
122*2c733ed4SBarry Smith       s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3;
123*2c733ed4SBarry Smith       s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3;
124*2c733ed4SBarry Smith       s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3;
125*2c733ed4SBarry Smith 
126*2c733ed4SBarry Smith       v +=  bs2;
127*2c733ed4SBarry Smith     }
128*2c733ed4SBarry Smith     /* x = inv_diagonal*x */
129*2c733ed4SBarry Smith     x[idt]   = v[0]*s1 + v[3]*s2 + v[6]*s3;
130*2c733ed4SBarry Smith     x[1+idt] = v[1]*s1 + v[4]*s2 + v[7]*s3;
131*2c733ed4SBarry Smith     x[2+idt] = v[2]*s1 + v[5]*s2 + v[8]*s3;
132*2c733ed4SBarry Smith 
133*2c733ed4SBarry Smith   }
134*2c733ed4SBarry Smith 
135*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
136*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
137*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr);
138*2c733ed4SBarry Smith   PetscFunctionReturn(0);
139*2c733ed4SBarry Smith }
140*2c733ed4SBarry Smith 
141*2c733ed4SBarry Smith PetscErrorCode MatForwardSolve_SeqBAIJ_3_NaturalOrdering(Mat A,Vec bb,Vec xx)
142*2c733ed4SBarry Smith {
143*2c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
144*2c733ed4SBarry Smith   const PetscInt    n  =a->mbs,*vi,*ai=a->i,*aj=a->j;
145*2c733ed4SBarry Smith   PetscErrorCode    ierr;
146*2c733ed4SBarry Smith   PetscInt          i,k,nz,idx,jdx;
147*2c733ed4SBarry Smith   const PetscInt    bs = A->rmap->bs,bs2 = a->bs2;
148*2c733ed4SBarry Smith   const MatScalar   *aa=a->a,*v;
149*2c733ed4SBarry Smith   PetscScalar       *x;
150*2c733ed4SBarry Smith   const PetscScalar *b;
151*2c733ed4SBarry Smith   PetscScalar       s1,s2,s3,x1,x2,x3;
152*2c733ed4SBarry Smith 
153*2c733ed4SBarry Smith   PetscFunctionBegin;
154*2c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
155*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
156*2c733ed4SBarry Smith   /* forward solve the lower triangular */
157*2c733ed4SBarry Smith   idx  = 0;
158*2c733ed4SBarry Smith   x[0] = b[idx]; x[1] = b[1+idx];x[2] = b[2+idx];
159*2c733ed4SBarry Smith   for (i=1; i<n; i++) {
160*2c733ed4SBarry Smith     v   = aa + bs2*ai[i];
161*2c733ed4SBarry Smith     vi  = aj + ai[i];
162*2c733ed4SBarry Smith     nz  = ai[i+1] - ai[i];
163*2c733ed4SBarry Smith     idx = bs*i;
164*2c733ed4SBarry Smith     s1  = b[idx];s2 = b[1+idx];s3 = b[2+idx];
165*2c733ed4SBarry Smith     for (k=0; k<nz; k++) {
166*2c733ed4SBarry Smith       jdx = bs*vi[k];
167*2c733ed4SBarry Smith       x1  = x[jdx];x2 = x[1+jdx]; x3 =x[2+jdx];
168*2c733ed4SBarry Smith       s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3;
169*2c733ed4SBarry Smith       s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3;
170*2c733ed4SBarry Smith       s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3;
171*2c733ed4SBarry Smith 
172*2c733ed4SBarry Smith       v +=  bs2;
173*2c733ed4SBarry Smith     }
174*2c733ed4SBarry Smith 
175*2c733ed4SBarry Smith     x[idx]   = s1;
176*2c733ed4SBarry Smith     x[1+idx] = s2;
177*2c733ed4SBarry Smith     x[2+idx] = s3;
178*2c733ed4SBarry Smith   }
179*2c733ed4SBarry Smith 
180*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
181*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
182*2c733ed4SBarry Smith   ierr = PetscLogFlops(1.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr);
183*2c733ed4SBarry Smith   PetscFunctionReturn(0);
184*2c733ed4SBarry Smith }
185*2c733ed4SBarry Smith 
186*2c733ed4SBarry Smith 
187*2c733ed4SBarry Smith PetscErrorCode MatBackwardSolve_SeqBAIJ_3_NaturalOrdering(Mat A,Vec bb,Vec xx)
188*2c733ed4SBarry Smith {
189*2c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
190*2c733ed4SBarry Smith   const PetscInt    n  =a->mbs,*vi,*aj=a->j,*adiag=a->diag;
191*2c733ed4SBarry Smith   PetscErrorCode    ierr;
192*2c733ed4SBarry Smith   PetscInt          i,k,nz,idx,idt;
193*2c733ed4SBarry Smith   const PetscInt    bs = A->rmap->bs,bs2 = a->bs2;
194*2c733ed4SBarry Smith   const MatScalar   *aa=a->a,*v;
195*2c733ed4SBarry Smith   PetscScalar       *x;
196*2c733ed4SBarry Smith   const PetscScalar *b;
197*2c733ed4SBarry Smith   PetscScalar       s1,s2,s3,x1,x2,x3;
198*2c733ed4SBarry Smith 
199*2c733ed4SBarry Smith   PetscFunctionBegin;
200*2c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
201*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
202*2c733ed4SBarry Smith 
203*2c733ed4SBarry Smith   /* backward solve the upper triangular */
204*2c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
205*2c733ed4SBarry Smith     v   = aa + bs2*(adiag[i+1]+1);
206*2c733ed4SBarry Smith     vi  = aj + adiag[i+1]+1;
207*2c733ed4SBarry Smith     nz  = adiag[i] - adiag[i+1]-1;
208*2c733ed4SBarry Smith     idt = bs*i;
209*2c733ed4SBarry Smith     s1  = b[idt];  s2 = b[1+idt];s3 = b[2+idt];
210*2c733ed4SBarry Smith 
211*2c733ed4SBarry Smith     for (k=0; k<nz; k++) {
212*2c733ed4SBarry Smith       idx = bs*vi[k];
213*2c733ed4SBarry Smith       x1  = x[idx];   x2 = x[1+idx]; x3 = x[2+idx];
214*2c733ed4SBarry Smith       s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3;
215*2c733ed4SBarry Smith       s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3;
216*2c733ed4SBarry Smith       s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3;
217*2c733ed4SBarry Smith 
218*2c733ed4SBarry Smith       v +=  bs2;
219*2c733ed4SBarry Smith     }
220*2c733ed4SBarry Smith     /* x = inv_diagonal*x */
221*2c733ed4SBarry Smith     x[idt]   = v[0]*s1 + v[3]*s2 + v[6]*s3;
222*2c733ed4SBarry Smith     x[1+idt] = v[1]*s1 + v[4]*s2 + v[7]*s3;
223*2c733ed4SBarry Smith     x[2+idt] = v[2]*s1 + v[5]*s2 + v[8]*s3;
224*2c733ed4SBarry Smith 
225*2c733ed4SBarry Smith   }
226*2c733ed4SBarry Smith 
227*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
228*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
229*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr);
230*2c733ed4SBarry Smith   PetscFunctionReturn(0);
231*2c733ed4SBarry Smith }
232