xref: /petsc/src/mat/impls/baij/seq/baijsolvnat1.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 /*
6*2c733ed4SBarry Smith       Special case where the matrix was ILU(0) factored in the natural
7*2c733ed4SBarry Smith    ordering. This eliminates the need for the column and row permutation.
8*2c733ed4SBarry Smith */
9*2c733ed4SBarry Smith PetscErrorCode MatSolve_SeqBAIJ_1_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
10*2c733ed4SBarry Smith {
11*2c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
12*2c733ed4SBarry Smith   const PetscInt    n  = a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
13*2c733ed4SBarry Smith   PetscErrorCode    ierr;
14*2c733ed4SBarry Smith   const MatScalar   *aa=a->a,*v;
15*2c733ed4SBarry Smith   PetscScalar       *x;
16*2c733ed4SBarry Smith   const PetscScalar *b;
17*2c733ed4SBarry Smith   PetscScalar       s1,x1;
18*2c733ed4SBarry Smith   PetscInt          jdx,idt,idx,nz,i;
19*2c733ed4SBarry Smith 
20*2c733ed4SBarry Smith   PetscFunctionBegin;
21*2c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
22*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
23*2c733ed4SBarry Smith 
24*2c733ed4SBarry Smith   /* forward solve the lower triangular */
25*2c733ed4SBarry Smith   idx  = 0;
26*2c733ed4SBarry Smith   x[0] = b[0];
27*2c733ed4SBarry Smith   for (i=1; i<n; i++) {
28*2c733ed4SBarry Smith     v    =  aa      + ai[i];
29*2c733ed4SBarry Smith     vi   =  aj      + ai[i];
30*2c733ed4SBarry Smith     nz   =  diag[i] - ai[i];
31*2c733ed4SBarry Smith     idx +=  1;
32*2c733ed4SBarry Smith     s1   =  b[idx];
33*2c733ed4SBarry Smith     while (nz--) {
34*2c733ed4SBarry Smith       jdx = *vi++;
35*2c733ed4SBarry Smith       x1  = x[jdx];
36*2c733ed4SBarry Smith       s1 -= v[0]*x1;
37*2c733ed4SBarry Smith       v  += 1;
38*2c733ed4SBarry Smith     }
39*2c733ed4SBarry Smith     x[idx] = s1;
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 + diag[i] + 1;
44*2c733ed4SBarry Smith     vi  = aj + diag[i] + 1;
45*2c733ed4SBarry Smith     nz  = ai[i+1] - diag[i] - 1;
46*2c733ed4SBarry Smith     idt = i;
47*2c733ed4SBarry Smith     s1  = x[idt];
48*2c733ed4SBarry Smith     while (nz--) {
49*2c733ed4SBarry Smith       idx = *vi++;
50*2c733ed4SBarry Smith       x1  = x[idx];
51*2c733ed4SBarry Smith       s1 -= v[0]*x1;
52*2c733ed4SBarry Smith       v  += 1;
53*2c733ed4SBarry Smith     }
54*2c733ed4SBarry Smith     v      = aa +  diag[i];
55*2c733ed4SBarry Smith     x[idt] = v[0]*s1;
56*2c733ed4SBarry Smith   }
57*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
58*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
59*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*(a->nz) - A->cmap->n);CHKERRQ(ierr);
60*2c733ed4SBarry Smith   PetscFunctionReturn(0);
61*2c733ed4SBarry Smith }
62*2c733ed4SBarry Smith 
63*2c733ed4SBarry Smith 
64*2c733ed4SBarry Smith PetscErrorCode MatForwardSolve_SeqBAIJ_1_NaturalOrdering(Mat A,Vec bb,Vec xx)
65*2c733ed4SBarry Smith {
66*2c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
67*2c733ed4SBarry Smith   PetscErrorCode    ierr;
68*2c733ed4SBarry Smith   const PetscInt    n = a->mbs,*ai = a->i,*aj = a->j,*vi;
69*2c733ed4SBarry Smith   PetscScalar       *x,sum;
70*2c733ed4SBarry Smith   const PetscScalar *b;
71*2c733ed4SBarry Smith   const MatScalar   *aa = a->a,*v;
72*2c733ed4SBarry Smith   PetscInt          i,nz;
73*2c733ed4SBarry Smith 
74*2c733ed4SBarry Smith   PetscFunctionBegin;
75*2c733ed4SBarry Smith   if (!n) PetscFunctionReturn(0);
76*2c733ed4SBarry Smith 
77*2c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
78*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
79*2c733ed4SBarry Smith 
80*2c733ed4SBarry Smith   /* forward solve the lower triangular */
81*2c733ed4SBarry Smith   x[0] = b[0];
82*2c733ed4SBarry Smith   v    = aa;
83*2c733ed4SBarry Smith   vi   = aj;
84*2c733ed4SBarry Smith   for (i=1; i<n; i++) {
85*2c733ed4SBarry Smith     nz  = ai[i+1] - ai[i];
86*2c733ed4SBarry Smith     sum = b[i];
87*2c733ed4SBarry Smith     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
88*2c733ed4SBarry Smith     v   += nz;
89*2c733ed4SBarry Smith     vi  += nz;
90*2c733ed4SBarry Smith     x[i] = sum;
91*2c733ed4SBarry Smith   }
92*2c733ed4SBarry Smith   ierr = PetscLogFlops(a->nz - A->cmap->n);CHKERRQ(ierr);
93*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
94*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
95*2c733ed4SBarry Smith   PetscFunctionReturn(0);
96*2c733ed4SBarry Smith }
97*2c733ed4SBarry Smith 
98*2c733ed4SBarry Smith PetscErrorCode MatBackwardSolve_SeqBAIJ_1_NaturalOrdering(Mat A,Vec bb,Vec xx)
99*2c733ed4SBarry Smith {
100*2c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
101*2c733ed4SBarry Smith   PetscErrorCode    ierr;
102*2c733ed4SBarry Smith   const PetscInt    n = a->mbs,*aj = a->j,*adiag = a->diag,*vi;
103*2c733ed4SBarry Smith   PetscScalar       *x,sum;
104*2c733ed4SBarry Smith   const PetscScalar *b;
105*2c733ed4SBarry Smith   const MatScalar   *aa = a->a,*v;
106*2c733ed4SBarry Smith   PetscInt          i,nz;
107*2c733ed4SBarry Smith 
108*2c733ed4SBarry Smith   PetscFunctionBegin;
109*2c733ed4SBarry Smith   if (!n) PetscFunctionReturn(0);
110*2c733ed4SBarry Smith 
111*2c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
112*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
113*2c733ed4SBarry Smith 
114*2c733ed4SBarry Smith   /* backward solve the upper triangular */
115*2c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
116*2c733ed4SBarry Smith     v   = aa + adiag[i+1] + 1;
117*2c733ed4SBarry Smith     vi  = aj + adiag[i+1] + 1;
118*2c733ed4SBarry Smith     nz  = adiag[i] - adiag[i+1]-1;
119*2c733ed4SBarry Smith     sum = b[i];
120*2c733ed4SBarry Smith     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
121*2c733ed4SBarry Smith     x[i] = sum*v[nz]; /* x[i]=aa[adiag[i]]*sum; v++; */
122*2c733ed4SBarry Smith   }
123*2c733ed4SBarry Smith 
124*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
125*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
126*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
127*2c733ed4SBarry Smith   PetscFunctionReturn(0);
128*2c733ed4SBarry Smith }
129*2c733ed4SBarry Smith 
130*2c733ed4SBarry Smith PetscErrorCode MatSolve_SeqBAIJ_1_NaturalOrdering(Mat A,Vec bb,Vec xx)
131*2c733ed4SBarry Smith {
132*2c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
133*2c733ed4SBarry Smith   PetscErrorCode    ierr;
134*2c733ed4SBarry Smith   const PetscInt    n = a->mbs,*ai = a->i,*aj = a->j,*adiag = a->diag,*vi;
135*2c733ed4SBarry Smith   PetscScalar       *x,sum;
136*2c733ed4SBarry Smith   const PetscScalar *b;
137*2c733ed4SBarry Smith   const MatScalar   *aa = a->a,*v;
138*2c733ed4SBarry Smith   PetscInt          i,nz;
139*2c733ed4SBarry Smith 
140*2c733ed4SBarry Smith   PetscFunctionBegin;
141*2c733ed4SBarry Smith   if (!n) PetscFunctionReturn(0);
142*2c733ed4SBarry Smith 
143*2c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
144*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
145*2c733ed4SBarry Smith 
146*2c733ed4SBarry Smith   /* forward solve the lower triangular */
147*2c733ed4SBarry Smith   x[0] = b[0];
148*2c733ed4SBarry Smith   v    = aa;
149*2c733ed4SBarry Smith   vi   = aj;
150*2c733ed4SBarry Smith   for (i=1; i<n; i++) {
151*2c733ed4SBarry Smith     nz  = ai[i+1] - ai[i];
152*2c733ed4SBarry Smith     sum = b[i];
153*2c733ed4SBarry Smith     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
154*2c733ed4SBarry Smith     v   += nz;
155*2c733ed4SBarry Smith     vi  += nz;
156*2c733ed4SBarry Smith     x[i] = sum;
157*2c733ed4SBarry Smith   }
158*2c733ed4SBarry Smith 
159*2c733ed4SBarry Smith   /* backward solve the upper triangular */
160*2c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
161*2c733ed4SBarry Smith     v   = aa + adiag[i+1] + 1;
162*2c733ed4SBarry Smith     vi  = aj + adiag[i+1] + 1;
163*2c733ed4SBarry Smith     nz  = adiag[i] - adiag[i+1]-1;
164*2c733ed4SBarry Smith     sum = x[i];
165*2c733ed4SBarry Smith     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
166*2c733ed4SBarry Smith     x[i] = sum*v[nz]; /* x[i]=aa[adiag[i]]*sum; v++; */
167*2c733ed4SBarry Smith   }
168*2c733ed4SBarry Smith 
169*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr);
170*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
171*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
172*2c733ed4SBarry Smith   PetscFunctionReturn(0);
173*2c733ed4SBarry Smith }
174*2c733ed4SBarry Smith 
175