xref: /petsc/src/mat/impls/sbaij/seq/relax.h (revision 54e8760d7c85eb9a87c81c504c093cab8d2a6da8)
1b5b17502SBarry Smith 
2b5b17502SBarry Smith /*
3b5b17502SBarry Smith     This is included by sbaij.c to generate unsigned short and regular versions of these two functions
4b5b17502SBarry Smith */
5*54e8760dSRichard Tran Mills 
6*54e8760dSRichard Tran Mills /* We cut-and-past below from aij.h to make a "no_function" version of PetscSparseDensePlusDot().
7*54e8760dSRichard Tran Mills  * This is necessary because the USESHORT case cannot use the inlined functions that may be employed. */
8*54e8760dSRichard Tran Mills 
9*54e8760dSRichard Tran Mills #if defined(PETSC_KERNEL_USE_UNROLL_4)
10*54e8760dSRichard Tran Mills #define PetscSparseDensePlusDot_no_function(sum,r,xv,xi,nnz) { \
11*54e8760dSRichard Tran Mills     if (nnz > 0) { \
12*54e8760dSRichard Tran Mills       switch (nnz & 0x3) { \
13*54e8760dSRichard Tran Mills       case 3: sum += *xv++ *r[*xi++]; \
14*54e8760dSRichard Tran Mills       case 2: sum += *xv++ *r[*xi++]; \
15*54e8760dSRichard Tran Mills       case 1: sum += *xv++ *r[*xi++]; \
16*54e8760dSRichard Tran Mills         nnz       -= 4;} \
17*54e8760dSRichard Tran Mills       while (nnz > 0) { \
18*54e8760dSRichard Tran Mills         sum +=  xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \
19*54e8760dSRichard Tran Mills                xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \
20*54e8760dSRichard Tran Mills         xv += 4; xi += 4; nnz -= 4; }}}
21*54e8760dSRichard Tran Mills 
22*54e8760dSRichard Tran Mills #elif defined(PETSC_KERNEL_USE_UNROLL_2)
23*54e8760dSRichard Tran Mills #define PetscSparseDensePlusDot_no_function(sum,r,xv,xi,nnz) { \
24*54e8760dSRichard Tran Mills     PetscInt __i,__i1,__i2; \
25*54e8760dSRichard Tran Mills     for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \
26*54e8760dSRichard Tran Mills                                     sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \
27*54e8760dSRichard Tran Mills     if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];}
28*54e8760dSRichard Tran Mills 
29*54e8760dSRichard Tran Mills #else
30*54e8760dSRichard Tran Mills #define PetscSparseDensePlusDot_no_function(sum,r,xv,xi,nnz) { \
31*54e8760dSRichard Tran Mills     PetscInt __i; \
32*54e8760dSRichard Tran Mills     for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];}
33*54e8760dSRichard Tran Mills #endif
34*54e8760dSRichard Tran Mills 
35*54e8760dSRichard Tran Mills 
36018dd85eSSatish Balay #if defined(USESHORT)
37018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian_ushort(Mat A,Vec xx,Vec zz)
38018dd85eSSatish Balay #else
39018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian(Mat A,Vec xx,Vec zz)
40018dd85eSSatish Balay #endif
41eeffb40dSHong Zhang {
42eeffb40dSHong Zhang   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
43eeffb40dSHong Zhang   const PetscScalar *x;
44eeffb40dSHong Zhang   PetscScalar       *z,x1,sum;
45eeffb40dSHong Zhang   const MatScalar   *v;
46eeffb40dSHong Zhang   MatScalar         vj;
47eeffb40dSHong Zhang   PetscErrorCode    ierr;
48eeffb40dSHong Zhang   PetscInt          mbs=a->mbs,i,j,nz;
49eeffb40dSHong Zhang   const PetscInt    *ai=a->i;
50eeffb40dSHong Zhang #if defined(USESHORT)
51eeffb40dSHong Zhang   const unsigned short *ib=a->jshort;
52eeffb40dSHong Zhang   unsigned short       ibt;
53eeffb40dSHong Zhang #else
54eeffb40dSHong Zhang   const PetscInt *ib=a->j;
55eeffb40dSHong Zhang   PetscInt       ibt;
56eeffb40dSHong Zhang #endif
578dca527aSShri Abhyankar   PetscInt nonzerorow = 0,jmin;
58eeffb40dSHong Zhang 
59eeffb40dSHong Zhang   PetscFunctionBegin;
60eeffb40dSHong Zhang   ierr = VecSet(zz,0.0);CHKERRQ(ierr);
613649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
62eeffb40dSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
63eeffb40dSHong Zhang 
64eeffb40dSHong Zhang   v = a->a;
65eeffb40dSHong Zhang   for (i=0; i<mbs; i++) {
66eeffb40dSHong Zhang     nz = ai[i+1] - ai[i];    /* length of i_th row of A */
6703d09022SShri Abhyankar     if (!nz) continue; /* Move to the next row if the current row is empty */
6803d09022SShri Abhyankar     nonzerorow++;
69eeffb40dSHong Zhang     x1   = x[i];
708dca527aSShri Abhyankar     sum  = 0.0;
718dca527aSShri Abhyankar     jmin = 0;
728dca527aSShri Abhyankar     if (ib[0] == i) {
73eeffb40dSHong Zhang       sum = v[0]*x1;           /* diagonal term */
748dca527aSShri Abhyankar       jmin++;
758dca527aSShri Abhyankar     }
768dca527aSShri Abhyankar     for (j=jmin; j<nz; j++) {
77a6f0575fSBarry Smith       ibt     = ib[j];
78eeffb40dSHong Zhang       vj      = v[j];
79eeffb40dSHong Zhang       sum    += vj * x[ibt]; /* (strict upper triangular part of A)*x  */
80a6f0575fSBarry Smith       z[ibt] += PetscConj(v[j]) * x1;    /* (strict lower triangular part of A)*x  */
81eeffb40dSHong Zhang     }
82eeffb40dSHong Zhang     z[i] += sum;
83eeffb40dSHong Zhang     v    +=    nz;
84eeffb40dSHong Zhang     ib   += nz;
85eeffb40dSHong Zhang   }
86eeffb40dSHong Zhang 
873649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
88eeffb40dSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
8903d09022SShri Abhyankar   ierr = PetscLogFlops(2.0*(2.0*a->nz - nonzerorow) - nonzerorow);CHKERRQ(ierr);
90eeffb40dSHong Zhang   PetscFunctionReturn(0);
91eeffb40dSHong Zhang }
92eeffb40dSHong Zhang 
93018dd85eSSatish Balay #if defined(USESHORT)
94018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_ushort(Mat A,Vec xx,Vec zz)
95018dd85eSSatish Balay #else
96018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1(Mat A,Vec xx,Vec zz)
97018dd85eSSatish Balay #endif
98b5b17502SBarry Smith {
99b5b17502SBarry Smith   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
100b5b17502SBarry Smith   const PetscScalar *x;
101b5b17502SBarry Smith   PetscScalar       *z,x1,sum;
102b5b17502SBarry Smith   const MatScalar   *v;
103b5b17502SBarry Smith   MatScalar         vj;
104b5b17502SBarry Smith   PetscErrorCode    ierr;
105b5b17502SBarry Smith   PetscInt          mbs=a->mbs,i,j,nz;
106b5b17502SBarry Smith   const PetscInt    *ai=a->i;
107b5b17502SBarry Smith #if defined(USESHORT)
108b5b17502SBarry Smith   const unsigned short *ib=a->jshort;
109b5b17502SBarry Smith   unsigned short       ibt;
110b5b17502SBarry Smith #else
111b5b17502SBarry Smith   const PetscInt *ib=a->j;
112b5b17502SBarry Smith   PetscInt       ibt;
113b5b17502SBarry Smith #endif
1148dca527aSShri Abhyankar   PetscInt nonzerorow=0,jmin;
115b5b17502SBarry Smith 
116b5b17502SBarry Smith   PetscFunctionBegin;
117b5b17502SBarry Smith   ierr = VecSet(zz,0.0);CHKERRQ(ierr);
1183649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
119b5b17502SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
120b5b17502SBarry Smith 
121b5b17502SBarry Smith   v = a->a;
122b5b17502SBarry Smith   for (i=0; i<mbs; i++) {
123b5b17502SBarry Smith     nz = ai[i+1] - ai[i];          /* length of i_th row of A */
12403d09022SShri Abhyankar     if (!nz) continue; /* Move to the next row if the current row is empty */
12503d09022SShri Abhyankar     nonzerorow++;
1268dca527aSShri Abhyankar     sum  = 0.0;
1278dca527aSShri Abhyankar     jmin = 0;
128b5b17502SBarry Smith     x1   = x[i];
1298dca527aSShri Abhyankar     if (ib[0] == i) {
130b5b17502SBarry Smith       sum = v[0]*x1;                 /* diagonal term */
1318dca527aSShri Abhyankar       jmin++;
1328dca527aSShri Abhyankar     }
133444d8c10SJed Brown     PetscPrefetchBlock(ib+nz,nz,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
134444d8c10SJed Brown     PetscPrefetchBlock(v+nz,nz,0,PETSC_PREFETCH_HINT_NTA);  /* Entries for the next row */
1358dca527aSShri Abhyankar     for (j=jmin; j<nz; j++) {
136b5b17502SBarry Smith       ibt     = ib[j];
137b5b17502SBarry Smith       vj      = v[j];
138b5b17502SBarry Smith       z[ibt] += vj * x1;       /* (strict lower triangular part of A)*x  */
139b5b17502SBarry Smith       sum    += vj * x[ibt];   /* (strict upper triangular part of A)*x  */
140b5b17502SBarry Smith     }
141b5b17502SBarry Smith     z[i] += sum;
142b5b17502SBarry Smith     v    += nz;
143b5b17502SBarry Smith     ib   += nz;
144b5b17502SBarry Smith   }
145b5b17502SBarry Smith 
1463649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
147b5b17502SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
14803d09022SShri Abhyankar   ierr = PetscLogFlops(2.0*(2.0*a->nz - nonzerorow) - nonzerorow);CHKERRQ(ierr);
149b5b17502SBarry Smith   PetscFunctionReturn(0);
150b5b17502SBarry Smith }
151b5b17502SBarry Smith 
152018dd85eSSatish Balay #if defined(USESHORT)
153018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ_ushort(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
154018dd85eSSatish Balay #else
155018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
156018dd85eSSatish Balay #endif
157b5b17502SBarry Smith {
158b5b17502SBarry Smith   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
159b5b17502SBarry Smith   const MatScalar   *aa=a->a,*v,*v1,*aidiag;
160b5b17502SBarry Smith   PetscScalar       *x,*t,sum;
161b5b17502SBarry Smith   const PetscScalar *b;
162b5b17502SBarry Smith   MatScalar         tmp;
163b5b17502SBarry Smith   PetscErrorCode    ierr;
164b5b17502SBarry Smith   PetscInt          m  =a->mbs,bs=A->rmap->bs,j;
165b5b17502SBarry Smith   const PetscInt    *ai=a->i;
166b5b17502SBarry Smith #if defined(USESHORT)
167b5b17502SBarry Smith   const unsigned short *aj=a->jshort,*vj,*vj1;
168b5b17502SBarry Smith #else
169b5b17502SBarry Smith   const PetscInt *aj=a->j,*vj,*vj1;
170b5b17502SBarry Smith #endif
171b5b17502SBarry Smith   PetscInt nz,nz1,i;
172b5b17502SBarry Smith 
173b5b17502SBarry Smith   PetscFunctionBegin;
174422a814eSBarry Smith   if (fshift == -1.0) fshift = 0.0; /* negative fshift indicates do not error on zero diagonal; this code never errors on zero diagonal */
175e32f2f54SBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
176b5b17502SBarry Smith 
177b5b17502SBarry Smith   its = its*lits;
178e32f2f54SBarry Smith   if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
179b5b17502SBarry Smith 
180e32f2f54SBarry Smith   if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
181b5b17502SBarry Smith 
182b5b17502SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1833649974fSBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
184b5b17502SBarry Smith 
185b5b17502SBarry Smith   if (!a->idiagvalid) {
186b5b17502SBarry Smith     if (!a->idiag) {
187785e854fSJed Brown       ierr = PetscMalloc1(m,&a->idiag);CHKERRQ(ierr);
188b5b17502SBarry Smith     }
189b5b17502SBarry Smith     for (i=0; i<a->mbs; i++) a->idiag[i] = 1.0/a->a[a->i[i]];
190b5b17502SBarry Smith     a->idiagvalid = PETSC_TRUE;
191b5b17502SBarry Smith   }
192b5b17502SBarry Smith 
19341f059aeSBarry Smith   if (!a->sor_work) {
194785e854fSJed Brown     ierr = PetscMalloc1(m,&a->sor_work);CHKERRQ(ierr);
195b5b17502SBarry Smith   }
19641f059aeSBarry Smith   t = a->sor_work;
197b5b17502SBarry Smith 
198b5b17502SBarry Smith   aidiag = a->idiag;
199b5b17502SBarry Smith 
200b5b17502SBarry Smith   if (flag == SOR_APPLY_UPPER) {
201b5b17502SBarry Smith     /* apply (U + D/omega) to the vector */
202b5b17502SBarry Smith     PetscScalar d;
203b5b17502SBarry Smith     for (i=0; i<m; i++) {
204b5b17502SBarry Smith       d   = fshift + aa[ai[i]];
205b5b17502SBarry Smith       nz  = ai[i+1] - ai[i] - 1;
206b5b17502SBarry Smith       vj  = aj + ai[i] + 1;
207b5b17502SBarry Smith       v   = aa + ai[i] + 1;
208b5b17502SBarry Smith       sum = b[i]*d/omega;
209*54e8760dSRichard Tran Mills #ifdef USESHORT
210*54e8760dSRichard Tran Mills       PetscSparseDensePlusDot_no_function(sum,b,v,vj,nz);
211*54e8760dSRichard Tran Mills #else
212b5b17502SBarry Smith       PetscSparseDensePlusDot(sum,b,v,vj,nz);
213*54e8760dSRichard Tran Mills #endif
214b5b17502SBarry Smith       x[i] = sum;
215b5b17502SBarry Smith     }
216b5b17502SBarry Smith     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
217b5b17502SBarry Smith   }
218b5b17502SBarry Smith 
219b5b17502SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
220b5b17502SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
221b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
222b5b17502SBarry Smith 
223b5b17502SBarry Smith       v  = aa + 1;
224b5b17502SBarry Smith       vj = aj + 1;
225b5b17502SBarry Smith       for (i=0; i<m; i++) {
226b5b17502SBarry Smith         nz  = ai[i+1] - ai[i] - 1;
227b5b17502SBarry Smith         tmp = -(x[i] = omega*t[i]*aidiag[i]);
22826fbe8dcSKarl Rupp         for (j=0; j<nz; j++) t[vj[j]] += tmp*v[j];
229b5b17502SBarry Smith         v  += nz + 1;
230b5b17502SBarry Smith         vj += nz + 1;
231b5b17502SBarry Smith       }
232b5b17502SBarry Smith       ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
233b5b17502SBarry Smith     }
234b5b17502SBarry Smith 
235b5b17502SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
236b5b17502SBarry Smith       int nz2;
237b5b17502SBarry Smith       if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) {
238b5b17502SBarry Smith #if defined(PETSC_USE_BACKWARD_LOOP)
239b5b17502SBarry Smith         v  = aa + ai[m] - 1;
240b5b17502SBarry Smith         vj = aj + ai[m] - 1;
241b5b17502SBarry Smith         for (i=m-1; i>=0; i--) {
242b5b17502SBarry Smith           sum = b[i];
243b5b17502SBarry Smith           nz  = ai[i+1] - ai[i] - 1;
244b5b17502SBarry Smith           {PetscInt __i;for (__i=0; __i<nz; __i++) sum -= v[-__i] * x[vj[-__i]];}
245b5b17502SBarry Smith #else
246b5b17502SBarry Smith         v  = aa + ai[m-1] + 1;
247b5b17502SBarry Smith         vj = aj + ai[m-1] + 1;
248b5b17502SBarry Smith         nz = 0;
249b5b17502SBarry Smith         for (i=m-1; i>=0; i--) {
250b5b17502SBarry Smith           sum = b[i];
2516c11f900SJed Brown           nz2 = ai[i] - ai[PetscMax(i-1,0)] - 1; /* avoid referencing ai[-1], nonsense nz2 is okay on last iteration */
25250d8bf02SJed Brown           PETSC_Prefetch(v-nz2-1,0,PETSC_PREFETCH_HINT_NTA);
25350d8bf02SJed Brown           PETSC_Prefetch(vj-nz2-1,0,PETSC_PREFETCH_HINT_NTA);
254b5b17502SBarry Smith           PetscSparseDenseMinusDot(sum,x,v,vj,nz);
255b5b17502SBarry Smith           nz = nz2;
256b5b17502SBarry Smith #endif
257b5b17502SBarry Smith           x[i] = omega*sum*aidiag[i];
258b5b17502SBarry Smith           v   -= nz + 1;
259b5b17502SBarry Smith           vj  -= nz + 1;
260b5b17502SBarry Smith         }
261b5b17502SBarry Smith         ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
262b5b17502SBarry Smith       } else {
263b5b17502SBarry Smith         v  = aa + ai[m-1] + 1;
264b5b17502SBarry Smith         vj = aj + ai[m-1] + 1;
265b5b17502SBarry Smith         nz = 0;
266b5b17502SBarry Smith         for (i=m-1; i>=0; i--) {
267b5b17502SBarry Smith           sum = t[i];
2686c11f900SJed Brown           nz2 = ai[i] - ai[PetscMax(i-1,0)] - 1; /* avoid referencing ai[-1], nonsense nz2 is okay on last iteration */
26950d8bf02SJed Brown           PETSC_Prefetch(v-nz2-1,0,PETSC_PREFETCH_HINT_NTA);
27050d8bf02SJed Brown           PETSC_Prefetch(vj-nz2-1,0,PETSC_PREFETCH_HINT_NTA);
271b5b17502SBarry Smith           PetscSparseDenseMinusDot(sum,x,v,vj,nz);
272b5b17502SBarry Smith           x[i] = (1-omega)*x[i] + omega*sum*aidiag[i];
273b5b17502SBarry Smith           nz   = nz2;
274b5b17502SBarry Smith           v   -= nz + 1;
275b5b17502SBarry Smith           vj  -= nz + 1;
276b5b17502SBarry Smith         }
277b5b17502SBarry Smith         ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
278b5b17502SBarry Smith       }
279b5b17502SBarry Smith     }
280b5b17502SBarry Smith     its--;
281b5b17502SBarry Smith   }
282b5b17502SBarry Smith 
283b5b17502SBarry Smith   while (its--) {
284b5b17502SBarry Smith     /*
285b5b17502SBarry Smith        forward sweep:
286b5b17502SBarry Smith        for i=0,...,m-1:
287b5b17502SBarry Smith          sum[i] = (b[i] - U(i,:)x)/d[i];
288b5b17502SBarry Smith          x[i]   = (1-omega)x[i] + omega*sum[i];
289b5b17502SBarry Smith          b      = b - x[i]*U^T(i,:);
290b5b17502SBarry Smith 
291b5b17502SBarry Smith     */
292b5b17502SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
293b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
294b5b17502SBarry Smith 
295b5b17502SBarry Smith       for (i=0; i<m; i++) {
296b5b17502SBarry Smith         v    = aa + ai[i] + 1; v1=v;
297b5b17502SBarry Smith         vj   = aj + ai[i] + 1; vj1=vj;
298b5b17502SBarry Smith         nz   = ai[i+1] - ai[i] - 1; nz1=nz;
299b5b17502SBarry Smith         sum  = t[i];
300b5b17502SBarry Smith         while (nz1--) sum -= (*v1++)*x[*vj1++];
301b5b17502SBarry Smith         x[i] = (1-omega)*x[i] + omega*sum*aidiag[i];
302b5b17502SBarry Smith         while (nz--) t[*vj++] -= x[i]*(*v++);
303b5b17502SBarry Smith       }
304dcca0238SPatrick Lacasse       ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr);
305b5b17502SBarry Smith     }
306b5b17502SBarry Smith 
307b5b17502SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
308b5b17502SBarry Smith       /*
309b5b17502SBarry Smith        backward sweep:
310b5b17502SBarry Smith        b = b - x[i]*U^T(i,:), i=0,...,n-2
311b5b17502SBarry Smith        for i=m-1,...,0:
312b5b17502SBarry Smith          sum[i] = (b[i] - U(i,:)x)/d[i];
313b5b17502SBarry Smith          x[i]   = (1-omega)x[i] + omega*sum[i];
314b5b17502SBarry Smith       */
315b5b17502SBarry Smith       /* if there was a forward sweep done above then I thing the next two for loops are not needed */
316b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
317b5b17502SBarry Smith 
318b5b17502SBarry Smith       for (i=0; i<m-1; i++) {  /* update rhs */
319b5b17502SBarry Smith         v    = aa + ai[i] + 1;
320b5b17502SBarry Smith         vj   = aj + ai[i] + 1;
321b5b17502SBarry Smith         nz   = ai[i+1] - ai[i] - 1;
322b5b17502SBarry Smith         while (nz--) t[*vj++] -= x[i]*(*v++);
323b5b17502SBarry Smith       }
324dcca0238SPatrick Lacasse       ierr = PetscLogFlops(2.0*(a->nz - m));CHKERRQ(ierr);
325b5b17502SBarry Smith       for (i=m-1; i>=0; i--) {
326b5b17502SBarry Smith         v    = aa + ai[i] + 1;
327b5b17502SBarry Smith         vj   = aj + ai[i] + 1;
328b5b17502SBarry Smith         nz   = ai[i+1] - ai[i] - 1;
329b5b17502SBarry Smith         sum  = t[i];
330b5b17502SBarry Smith         while (nz--) sum -= x[*vj++]*(*v++);
331b5b17502SBarry Smith         x[i] =   (1-omega)*x[i] + omega*sum*aidiag[i];
332b5b17502SBarry Smith       }
333dcca0238SPatrick Lacasse       ierr = PetscLogFlops(2.0*(a->nz + m));CHKERRQ(ierr);
334b5b17502SBarry Smith     }
335b5b17502SBarry Smith   }
336b5b17502SBarry Smith 
337b5b17502SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3383649974fSBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
339b5b17502SBarry Smith   PetscFunctionReturn(0);
340b5b17502SBarry Smith }
341