xref: /petsc/src/mat/impls/sbaij/seq/relax.h (revision 26fbe8dc1a3c99fb8dddfa572c8c6b3b4ce3ca53)
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 */
5b5b17502SBarry Smith #undef __FUNCT__
6018dd85eSSatish Balay #if defined(USESHORT)
7b583852cSJed Brown #define __FUNCT__ "MatMult_SeqSBAIJ_1_Hermitian_ushort"
8018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian_ushort(Mat A,Vec xx,Vec zz)
9018dd85eSSatish Balay #else
10b583852cSJed Brown #define __FUNCT__ "MatMult_SeqSBAIJ_1_Hermitian"
11018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian(Mat A,Vec xx,Vec zz)
12018dd85eSSatish Balay #endif
13eeffb40dSHong Zhang {
14eeffb40dSHong Zhang   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
15eeffb40dSHong Zhang   const PetscScalar *x;
16eeffb40dSHong Zhang   PetscScalar       *z,x1,sum;
17eeffb40dSHong Zhang   const MatScalar   *v;
18eeffb40dSHong Zhang   MatScalar         vj;
19eeffb40dSHong Zhang   PetscErrorCode    ierr;
20eeffb40dSHong Zhang   PetscInt          mbs=a->mbs,i,j,nz;
21eeffb40dSHong Zhang   const PetscInt    *ai=a->i;
22eeffb40dSHong Zhang #if defined(USESHORT)
23eeffb40dSHong Zhang   const unsigned short *ib=a->jshort;
24eeffb40dSHong Zhang   unsigned short       ibt;
25eeffb40dSHong Zhang #else
26eeffb40dSHong Zhang   const PetscInt *ib=a->j;
27eeffb40dSHong Zhang   PetscInt       ibt;
28eeffb40dSHong Zhang #endif
298dca527aSShri Abhyankar   PetscInt nonzerorow = 0,jmin;
30eeffb40dSHong Zhang 
31eeffb40dSHong Zhang   PetscFunctionBegin;
32eeffb40dSHong Zhang   ierr = VecSet(zz,0.0);CHKERRQ(ierr);
333649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
34eeffb40dSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
35eeffb40dSHong Zhang 
36eeffb40dSHong Zhang   v = a->a;
37eeffb40dSHong Zhang   for (i=0; i<mbs; i++) {
38eeffb40dSHong Zhang     nz = ai[i+1] - ai[i];    /* length of i_th row of A */
3903d09022SShri Abhyankar     if (!nz) continue; /* Move to the next row if the current row is empty */
4003d09022SShri Abhyankar     nonzerorow++;
41eeffb40dSHong Zhang     x1   = x[i];
428dca527aSShri Abhyankar     sum  = 0.0;
438dca527aSShri Abhyankar     jmin = 0;
448dca527aSShri Abhyankar     if (ib[0] == i) {
45eeffb40dSHong Zhang       sum = v[0]*x1;           /* diagonal term */
468dca527aSShri Abhyankar       jmin++;
478dca527aSShri Abhyankar     }
488dca527aSShri Abhyankar     for (j=jmin; j<nz; j++) {
49a6f0575fSBarry Smith       ibt     = ib[j];
50eeffb40dSHong Zhang       vj      = v[j];
51eeffb40dSHong Zhang       sum    += vj * x[ibt]; /* (strict upper triangular part of A)*x  */
52a6f0575fSBarry Smith       z[ibt] += PetscConj(v[j]) * x1;    /* (strict lower triangular part of A)*x  */
53eeffb40dSHong Zhang     }
54eeffb40dSHong Zhang     z[i] += sum;
55eeffb40dSHong Zhang     v    +=    nz;
56eeffb40dSHong Zhang     ib   += nz;
57eeffb40dSHong Zhang   }
58eeffb40dSHong Zhang 
593649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
60eeffb40dSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
6103d09022SShri Abhyankar   ierr = PetscLogFlops(2.0*(2.0*a->nz - nonzerorow) - nonzerorow);CHKERRQ(ierr);
62eeffb40dSHong Zhang   PetscFunctionReturn(0);
63eeffb40dSHong Zhang }
64eeffb40dSHong Zhang 
65eeffb40dSHong Zhang #undef __FUNCT__
66018dd85eSSatish Balay #if defined(USESHORT)
67b583852cSJed Brown #define __FUNCT__ "MatMult_SeqSBAIJ_1_ushort"
68018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_ushort(Mat A,Vec xx,Vec zz)
69018dd85eSSatish Balay #else
70b583852cSJed Brown #define __FUNCT__ "MatMult_SeqSBAIJ_1"
71018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1(Mat A,Vec xx,Vec zz)
72018dd85eSSatish Balay #endif
73b5b17502SBarry Smith {
74b5b17502SBarry Smith   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
75b5b17502SBarry Smith   const PetscScalar *x;
76b5b17502SBarry Smith   PetscScalar       *z,x1,sum;
77b5b17502SBarry Smith   const MatScalar   *v;
78b5b17502SBarry Smith   MatScalar         vj;
79b5b17502SBarry Smith   PetscErrorCode    ierr;
80b5b17502SBarry Smith   PetscInt          mbs=a->mbs,i,j,nz;
81b5b17502SBarry Smith   const PetscInt    *ai=a->i;
82b5b17502SBarry Smith #if defined(USESHORT)
83b5b17502SBarry Smith   const unsigned short *ib=a->jshort;
84b5b17502SBarry Smith   unsigned short       ibt;
85b5b17502SBarry Smith #else
86b5b17502SBarry Smith   const PetscInt *ib=a->j;
87b5b17502SBarry Smith   PetscInt       ibt;
88b5b17502SBarry Smith #endif
898dca527aSShri Abhyankar   PetscInt nonzerorow=0,jmin;
90b5b17502SBarry Smith 
91b5b17502SBarry Smith   PetscFunctionBegin;
92b5b17502SBarry Smith   ierr = VecSet(zz,0.0);CHKERRQ(ierr);
933649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
94b5b17502SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
95b5b17502SBarry Smith 
96b5b17502SBarry Smith   v = a->a;
97b5b17502SBarry Smith   for (i=0; i<mbs; i++) {
98b5b17502SBarry Smith     nz = ai[i+1] - ai[i];          /* length of i_th row of A */
9903d09022SShri Abhyankar     if (!nz) continue; /* Move to the next row if the current row is empty */
10003d09022SShri Abhyankar     nonzerorow++;
1018dca527aSShri Abhyankar     sum  = 0.0;
1028dca527aSShri Abhyankar     jmin = 0;
103b5b17502SBarry Smith     x1   = x[i];
1048dca527aSShri Abhyankar     if (ib[0] == i) {
105b5b17502SBarry Smith       sum = v[0]*x1;                 /* diagonal term */
1068dca527aSShri Abhyankar       jmin++;
1078dca527aSShri Abhyankar     }
108444d8c10SJed Brown     PetscPrefetchBlock(ib+nz,nz,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
109444d8c10SJed Brown     PetscPrefetchBlock(v+nz,nz,0,PETSC_PREFETCH_HINT_NTA);  /* Entries for the next row */
1108dca527aSShri Abhyankar     for (j=jmin; j<nz; j++) {
111b5b17502SBarry Smith       ibt     = ib[j];
112b5b17502SBarry Smith       vj      = v[j];
113b5b17502SBarry Smith       z[ibt] += vj * x1;       /* (strict lower triangular part of A)*x  */
114b5b17502SBarry Smith       sum    += vj * x[ibt];   /* (strict upper triangular part of A)*x  */
115b5b17502SBarry Smith     }
116b5b17502SBarry Smith     z[i] += sum;
117b5b17502SBarry Smith     v    += nz;
118b5b17502SBarry Smith     ib   += nz;
119b5b17502SBarry Smith   }
120b5b17502SBarry Smith 
1213649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
122b5b17502SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
12303d09022SShri Abhyankar   ierr = PetscLogFlops(2.0*(2.0*a->nz - nonzerorow) - nonzerorow);CHKERRQ(ierr);
124b5b17502SBarry Smith   PetscFunctionReturn(0);
125b5b17502SBarry Smith }
126b5b17502SBarry Smith 
127b5b17502SBarry Smith #undef __FUNCT__
128018dd85eSSatish Balay #if defined(USESHORT)
129b583852cSJed Brown #define __FUNCT__ "MatSOR_SeqSBAIJ_ushort"
130018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ_ushort(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
131018dd85eSSatish Balay #else
132b583852cSJed Brown #define __FUNCT__ "MatSOR_SeqSBAIJ"
133018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
134018dd85eSSatish Balay #endif
135b5b17502SBarry Smith {
136b5b17502SBarry Smith   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
137b5b17502SBarry Smith   const MatScalar   *aa=a->a,*v,*v1,*aidiag;
138b5b17502SBarry Smith   PetscScalar       *x,*t,sum;
139b5b17502SBarry Smith   const PetscScalar *b;
140b5b17502SBarry Smith   MatScalar         tmp;
141b5b17502SBarry Smith   PetscErrorCode    ierr;
142b5b17502SBarry Smith   PetscInt          m  =a->mbs,bs=A->rmap->bs,j;
143b5b17502SBarry Smith   const PetscInt    *ai=a->i;
144b5b17502SBarry Smith #if defined(USESHORT)
145b5b17502SBarry Smith   const unsigned short *aj=a->jshort,*vj,*vj1;
146b5b17502SBarry Smith #else
147b5b17502SBarry Smith   const PetscInt *aj=a->j,*vj,*vj1;
148b5b17502SBarry Smith #endif
149b5b17502SBarry Smith   PetscInt nz,nz1,i;
150b5b17502SBarry Smith 
151b5b17502SBarry Smith   PetscFunctionBegin;
152e32f2f54SBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
153b5b17502SBarry Smith 
154b5b17502SBarry Smith   its = its*lits;
155e32f2f54SBarry 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);
156b5b17502SBarry Smith 
157e32f2f54SBarry Smith   if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
158b5b17502SBarry Smith 
159b5b17502SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1603649974fSBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
161b5b17502SBarry Smith 
162b5b17502SBarry Smith   if (!a->idiagvalid) {
163b5b17502SBarry Smith     if (!a->idiag) {
164b5b17502SBarry Smith       ierr = PetscMalloc(m*sizeof(PetscScalar),&a->idiag);CHKERRQ(ierr);
165b5b17502SBarry Smith     }
166b5b17502SBarry Smith     for (i=0; i<a->mbs; i++) a->idiag[i] = 1.0/a->a[a->i[i]];
167b5b17502SBarry Smith     a->idiagvalid = PETSC_TRUE;
168b5b17502SBarry Smith   }
169b5b17502SBarry Smith 
17041f059aeSBarry Smith   if (!a->sor_work) {
17141f059aeSBarry Smith     ierr = PetscMalloc(m*sizeof(PetscScalar),&a->sor_work);CHKERRQ(ierr);
172b5b17502SBarry Smith   }
17341f059aeSBarry Smith   t = a->sor_work;
174b5b17502SBarry Smith 
175b5b17502SBarry Smith   aidiag = a->idiag;
176b5b17502SBarry Smith 
177b5b17502SBarry Smith   if (flag == SOR_APPLY_UPPER) {
178b5b17502SBarry Smith     /* apply (U + D/omega) to the vector */
179b5b17502SBarry Smith     PetscScalar d;
180b5b17502SBarry Smith     for (i=0; i<m; i++) {
181b5b17502SBarry Smith       d   = fshift + aa[ai[i]];
182b5b17502SBarry Smith       nz  = ai[i+1] - ai[i] - 1;
183b5b17502SBarry Smith       vj  = aj + ai[i] + 1;
184b5b17502SBarry Smith       v   = aa + ai[i] + 1;
185b5b17502SBarry Smith       sum = b[i]*d/omega;
186b5b17502SBarry Smith       PetscSparseDensePlusDot(sum,b,v,vj,nz);
187b5b17502SBarry Smith       x[i] = sum;
188b5b17502SBarry Smith     }
189b5b17502SBarry Smith     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
190b5b17502SBarry Smith   }
191b5b17502SBarry Smith 
192b5b17502SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
193b5b17502SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
194b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
195b5b17502SBarry Smith 
196b5b17502SBarry Smith       v  = aa + 1;
197b5b17502SBarry Smith       vj = aj + 1;
198b5b17502SBarry Smith       for (i=0; i<m; i++) {
199b5b17502SBarry Smith         nz  = ai[i+1] - ai[i] - 1;
200b5b17502SBarry Smith         tmp = -(x[i] = omega*t[i]*aidiag[i]);
201*26fbe8dcSKarl Rupp         for (j=0; j<nz; j++) t[vj[j]] += tmp*v[j];
202b5b17502SBarry Smith         v  += nz + 1;
203b5b17502SBarry Smith         vj += nz + 1;
204b5b17502SBarry Smith       }
205b5b17502SBarry Smith       ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
206b5b17502SBarry Smith     }
207b5b17502SBarry Smith 
208b5b17502SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
209b5b17502SBarry Smith       int nz2;
210b5b17502SBarry Smith       if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) {
211b5b17502SBarry Smith #if defined(PETSC_USE_BACKWARD_LOOP)
212b5b17502SBarry Smith         v  = aa + ai[m] - 1;
213b5b17502SBarry Smith         vj = aj + ai[m] - 1;
214b5b17502SBarry Smith         for (i=m-1; i>=0; i--) {
215b5b17502SBarry Smith           sum = b[i];
216b5b17502SBarry Smith           nz  = ai[i+1] - ai[i] - 1;
217b5b17502SBarry Smith           {PetscInt __i;for (__i=0; __i<nz; __i++) sum -= v[-__i] * x[vj[-__i]];}
218b5b17502SBarry Smith #else
219b5b17502SBarry Smith         v  = aa + ai[m-1] + 1;
220b5b17502SBarry Smith         vj = aj + ai[m-1] + 1;
221b5b17502SBarry Smith         nz = 0;
222b5b17502SBarry Smith         for (i=m-1; i>=0; i--) {
223b5b17502SBarry Smith           sum = b[i];
224b5b17502SBarry Smith           nz2 = ai[i] - ai[i-1] - 1;
22550d8bf02SJed Brown           PETSC_Prefetch(v-nz2-1,0,PETSC_PREFETCH_HINT_NTA);
22650d8bf02SJed Brown           PETSC_Prefetch(vj-nz2-1,0,PETSC_PREFETCH_HINT_NTA);
227b5b17502SBarry Smith           PetscSparseDenseMinusDot(sum,x,v,vj,nz);
228b5b17502SBarry Smith           nz = nz2;
229b5b17502SBarry Smith #endif
230b5b17502SBarry Smith           x[i] = omega*sum*aidiag[i];
231b5b17502SBarry Smith           v   -= nz + 1;
232b5b17502SBarry Smith           vj  -= nz + 1;
233b5b17502SBarry Smith         }
234b5b17502SBarry Smith         ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
235b5b17502SBarry Smith       } else {
236b5b17502SBarry Smith         v  = aa + ai[m-1] + 1;
237b5b17502SBarry Smith         vj = aj + ai[m-1] + 1;
238b5b17502SBarry Smith         nz = 0;
239b5b17502SBarry Smith         for (i=m-1; i>=0; i--) {
240b5b17502SBarry Smith           sum = t[i];
241b5b17502SBarry Smith           nz2 = ai[i] - ai[i-1] - 1;
24250d8bf02SJed Brown           PETSC_Prefetch(v-nz2-1,0,PETSC_PREFETCH_HINT_NTA);
24350d8bf02SJed Brown           PETSC_Prefetch(vj-nz2-1,0,PETSC_PREFETCH_HINT_NTA);
244b5b17502SBarry Smith           PetscSparseDenseMinusDot(sum,x,v,vj,nz);
245b5b17502SBarry Smith           x[i] = (1-omega)*x[i] + omega*sum*aidiag[i];
246b5b17502SBarry Smith           nz   = nz2;
247b5b17502SBarry Smith           v   -= nz + 1;
248b5b17502SBarry Smith           vj  -= nz + 1;
249b5b17502SBarry Smith         }
250b5b17502SBarry Smith         ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
251b5b17502SBarry Smith       }
252b5b17502SBarry Smith     }
253b5b17502SBarry Smith     its--;
254b5b17502SBarry Smith   }
255b5b17502SBarry Smith 
256b5b17502SBarry Smith   while (its--) {
257b5b17502SBarry Smith     /*
258b5b17502SBarry Smith        forward sweep:
259b5b17502SBarry Smith        for i=0,...,m-1:
260b5b17502SBarry Smith          sum[i] = (b[i] - U(i,:)x)/d[i];
261b5b17502SBarry Smith          x[i]   = (1-omega)x[i] + omega*sum[i];
262b5b17502SBarry Smith          b      = b - x[i]*U^T(i,:);
263b5b17502SBarry Smith 
264b5b17502SBarry Smith     */
265b5b17502SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
266b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
267b5b17502SBarry Smith 
268b5b17502SBarry Smith       for (i=0; i<m; i++) {
269b5b17502SBarry Smith         v    = aa + ai[i] + 1; v1=v;
270b5b17502SBarry Smith         vj   = aj + ai[i] + 1; vj1=vj;
271b5b17502SBarry Smith         nz   = ai[i+1] - ai[i] - 1; nz1=nz;
272b5b17502SBarry Smith         sum  = t[i];
273b5b17502SBarry Smith         ierr = PetscLogFlops(4.0*nz-2);CHKERRQ(ierr);
274b5b17502SBarry Smith         while (nz1--) sum -= (*v1++)*x[*vj1++];
275b5b17502SBarry Smith         x[i] = (1-omega)*x[i] + omega*sum*aidiag[i];
276b5b17502SBarry Smith         while (nz--) t[*vj++] -= x[i]*(*v++);
277b5b17502SBarry Smith       }
278b5b17502SBarry Smith     }
279b5b17502SBarry Smith 
280b5b17502SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
281b5b17502SBarry Smith       /*
282b5b17502SBarry Smith        backward sweep:
283b5b17502SBarry Smith        b = b - x[i]*U^T(i,:), i=0,...,n-2
284b5b17502SBarry Smith        for i=m-1,...,0:
285b5b17502SBarry Smith          sum[i] = (b[i] - U(i,:)x)/d[i];
286b5b17502SBarry Smith          x[i]   = (1-omega)x[i] + omega*sum[i];
287b5b17502SBarry Smith       */
288b5b17502SBarry Smith       /* if there was a forward sweep done above then I thing the next two for loops are not needed */
289b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
290b5b17502SBarry Smith 
291b5b17502SBarry Smith       for (i=0; i<m-1; i++) {  /* update rhs */
292b5b17502SBarry Smith         v    = aa + ai[i] + 1;
293b5b17502SBarry Smith         vj   = aj + ai[i] + 1;
294b5b17502SBarry Smith         nz   = ai[i+1] - ai[i] - 1;
295b5b17502SBarry Smith         ierr = PetscLogFlops(2.0*nz-1);CHKERRQ(ierr);
296b5b17502SBarry Smith         while (nz--) t[*vj++] -= x[i]*(*v++);
297b5b17502SBarry Smith       }
298b5b17502SBarry Smith       for (i=m-1; i>=0; i--) {
299b5b17502SBarry Smith         v    = aa + ai[i] + 1;
300b5b17502SBarry Smith         vj   = aj + ai[i] + 1;
301b5b17502SBarry Smith         nz   = ai[i+1] - ai[i] - 1;
302b5b17502SBarry Smith         ierr = PetscLogFlops(2.0*nz-1);CHKERRQ(ierr);
303b5b17502SBarry Smith         sum  = t[i];
304b5b17502SBarry Smith         while (nz--) sum -= x[*vj++]*(*v++);
305b5b17502SBarry Smith         x[i] =   (1-omega)*x[i] + omega*sum*aidiag[i];
306b5b17502SBarry Smith       }
307b5b17502SBarry Smith     }
308b5b17502SBarry Smith   }
309b5b17502SBarry Smith 
310b5b17502SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3113649974fSBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
312b5b17502SBarry Smith   PetscFunctionReturn(0);
313b5b17502SBarry Smith }
314