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