xref: /petsc/src/mat/impls/sbaij/seq/relax.h (revision 03d09022fc983d93283cb27fabb49a20e7633e00)
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__
6eeffb40dSHong Zhang #define __FUNCT__ "MatMult_SeqSBAIJ_1_Hermitian"
7018dd85eSSatish Balay #if defined(USESHORT)
8018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian_ushort(Mat A,Vec xx,Vec zz)
9018dd85eSSatish Balay #else
10018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian(Mat A,Vec xx,Vec zz)
11018dd85eSSatish Balay #endif
12eeffb40dSHong Zhang {
13eeffb40dSHong Zhang   Mat_SeqSBAIJ         *a = (Mat_SeqSBAIJ*)A->data;
14eeffb40dSHong Zhang   const PetscScalar    *x;
15eeffb40dSHong Zhang   PetscScalar          *z,x1,sum;
16eeffb40dSHong Zhang   const MatScalar      *v;
17eeffb40dSHong Zhang   MatScalar            vj;
18eeffb40dSHong Zhang   PetscErrorCode       ierr;
19eeffb40dSHong Zhang   PetscInt             mbs=a->mbs,i,j,nz;
20eeffb40dSHong Zhang   const PetscInt       *ai=a->i;
21eeffb40dSHong Zhang #if defined(USESHORT)
22eeffb40dSHong Zhang   const unsigned short *ib=a->jshort;
23eeffb40dSHong Zhang   unsigned short       ibt;
24eeffb40dSHong Zhang #else
25eeffb40dSHong Zhang   const PetscInt       *ib=a->j;
26eeffb40dSHong Zhang   PetscInt             ibt;
27eeffb40dSHong Zhang #endif
28*03d09022SShri Abhyankar   PetscInt             nonzerorow;
29eeffb40dSHong Zhang 
30eeffb40dSHong Zhang   PetscFunctionBegin;
31eeffb40dSHong Zhang   ierr = VecSet(zz,0.0);CHKERRQ(ierr);
32eeffb40dSHong Zhang   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
33eeffb40dSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
34eeffb40dSHong Zhang 
35eeffb40dSHong Zhang   v  = a->a;
36eeffb40dSHong Zhang   for (i=0; i<mbs; i++) {
37eeffb40dSHong Zhang     nz   = ai[i+1] - ai[i];  /* length of i_th row of A */
38*03d09022SShri Abhyankar     if (!nz) continue; /* Move to the next row if the current row is empty */
39*03d09022SShri Abhyankar     nonzerorow++;
40eeffb40dSHong Zhang     x1   = x[i];
41eeffb40dSHong Zhang     sum  = v[0]*x1;          /* diagonal term */
42eeffb40dSHong Zhang     for (j=1; j<nz; j++) {
43a6f0575fSBarry Smith       ibt  = ib[j];
44eeffb40dSHong Zhang       vj   = v[j];
45eeffb40dSHong Zhang       sum += vj * x[ibt];   /* (strict upper triangular part of A)*x  */
46a6f0575fSBarry Smith       z[ibt] += PetscConj(v[j]) * x1;    /* (strict lower triangular part of A)*x  */
47eeffb40dSHong Zhang     }
48eeffb40dSHong Zhang     z[i] += sum;
49eeffb40dSHong Zhang     v    += nz;
50eeffb40dSHong Zhang     ib   += nz;
51eeffb40dSHong Zhang   }
52eeffb40dSHong Zhang 
53eeffb40dSHong Zhang   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
54eeffb40dSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
55*03d09022SShri Abhyankar   ierr = PetscLogFlops(2.0*(2.0*a->nz - nonzerorow) - nonzerorow);CHKERRQ(ierr);
56eeffb40dSHong Zhang   PetscFunctionReturn(0);
57eeffb40dSHong Zhang }
58eeffb40dSHong Zhang 
59eeffb40dSHong Zhang #undef __FUNCT__
60b5b17502SBarry Smith #define __FUNCT__ "MatMult_SeqSBAIJ_1"
61018dd85eSSatish Balay #if defined(USESHORT)
62018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_ushort(Mat A,Vec xx,Vec zz)
63018dd85eSSatish Balay #else
64018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1(Mat A,Vec xx,Vec zz)
65018dd85eSSatish Balay #endif
66b5b17502SBarry Smith {
67b5b17502SBarry Smith   Mat_SeqSBAIJ         *a = (Mat_SeqSBAIJ*)A->data;
68b5b17502SBarry Smith   const PetscScalar    *x;
69b5b17502SBarry Smith   PetscScalar          *z,x1,sum;
70b5b17502SBarry Smith   const MatScalar      *v;
71b5b17502SBarry Smith   MatScalar            vj;
72b5b17502SBarry Smith   PetscErrorCode       ierr;
73b5b17502SBarry Smith   PetscInt             mbs=a->mbs,i,j,nz;
74b5b17502SBarry Smith   const PetscInt       *ai=a->i;
75b5b17502SBarry Smith #if defined(USESHORT)
76b5b17502SBarry Smith   const unsigned short *ib=a->jshort;
77b5b17502SBarry Smith   unsigned short       ibt;
78b5b17502SBarry Smith #else
79b5b17502SBarry Smith   const PetscInt       *ib=a->j;
80b5b17502SBarry Smith   PetscInt             ibt;
81b5b17502SBarry Smith #endif
82*03d09022SShri Abhyankar   PetscInt             nonzerorow=0;
83b5b17502SBarry Smith 
84b5b17502SBarry Smith   PetscFunctionBegin;
85b5b17502SBarry Smith   ierr = VecSet(zz,0.0);CHKERRQ(ierr);
86b5b17502SBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
87b5b17502SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
88b5b17502SBarry Smith 
89b5b17502SBarry Smith   v  = a->a;
90b5b17502SBarry Smith   for (i=0; i<mbs; i++) {
91b5b17502SBarry Smith     nz   = ai[i+1] - ai[i];        /* length of i_th row of A */
92*03d09022SShri Abhyankar     if (!nz) continue; /* Move to the next row if the current row is empty */
93*03d09022SShri Abhyankar     nonzerorow++;
94b5b17502SBarry Smith     x1   = x[i];
95b5b17502SBarry Smith     sum  = v[0]*x1;                /* diagonal term */
96b5b17502SBarry Smith     for (j=1; j<nz; j++) {
97b5b17502SBarry Smith       ibt = ib[j];
98b5b17502SBarry Smith       vj  = v[j];
99b5b17502SBarry Smith       z[ibt] += vj * x1;       /* (strict lower triangular part of A)*x  */
100b5b17502SBarry Smith       sum    += vj * x[ibt]; /* (strict upper triangular part of A)*x  */
101b5b17502SBarry Smith     }
102b5b17502SBarry Smith     z[i] += sum;
103b5b17502SBarry Smith     v    += nz;
104b5b17502SBarry Smith     ib   += nz;
105b5b17502SBarry Smith   }
106b5b17502SBarry Smith 
107b5b17502SBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
108b5b17502SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
109*03d09022SShri Abhyankar   ierr = PetscLogFlops(2.0*(2.0*a->nz - nonzerorow) - nonzerorow);CHKERRQ(ierr);
110b5b17502SBarry Smith   PetscFunctionReturn(0);
111b5b17502SBarry Smith }
112b5b17502SBarry Smith 
113b5b17502SBarry Smith #undef __FUNCT__
11441f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqSBAIJ"
115018dd85eSSatish Balay #if defined(USESHORT)
116018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ_ushort(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
117018dd85eSSatish Balay #else
118018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
119018dd85eSSatish Balay #endif
120b5b17502SBarry Smith {
121b5b17502SBarry Smith   Mat_SeqSBAIJ         *a = (Mat_SeqSBAIJ*)A->data;
122b5b17502SBarry Smith   const MatScalar      *aa=a->a,*v,*v1,*aidiag;
123b5b17502SBarry Smith   PetscScalar          *x,*t,sum;
124b5b17502SBarry Smith   const PetscScalar    *b;
125b5b17502SBarry Smith   MatScalar            tmp;
126b5b17502SBarry Smith   PetscErrorCode       ierr;
127b5b17502SBarry Smith   PetscInt             m=a->mbs,bs=A->rmap->bs,j;
128b5b17502SBarry Smith   const PetscInt       *ai=a->i;
129b5b17502SBarry Smith #if defined(USESHORT)
130b5b17502SBarry Smith   const unsigned short *aj=a->jshort,*vj,*vj1;
131b5b17502SBarry Smith #else
132b5b17502SBarry Smith   const PetscInt       *aj=a->j,*vj,*vj1;
133b5b17502SBarry Smith #endif
134b5b17502SBarry Smith   PetscInt             nz,nz1,i;
135b5b17502SBarry Smith 
136b5b17502SBarry Smith   PetscFunctionBegin;
137b5b17502SBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat");
138b5b17502SBarry Smith 
139b5b17502SBarry Smith   its = its*lits;
140b5b17502SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
141b5b17502SBarry Smith 
142b5b17502SBarry Smith   if (bs > 1) SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
143b5b17502SBarry Smith 
144b5b17502SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
145b5b17502SBarry Smith   if (xx != bb) {
146b5b17502SBarry Smith     ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
147b5b17502SBarry Smith   } else {
148b5b17502SBarry Smith     b = x;
149b5b17502SBarry Smith   }
150b5b17502SBarry Smith 
151b5b17502SBarry Smith   if (!a->idiagvalid) {
152b5b17502SBarry Smith     if (!a->idiag) {
153b5b17502SBarry Smith       ierr = PetscMalloc(m*sizeof(PetscScalar),&a->idiag);CHKERRQ(ierr);
154b5b17502SBarry Smith     }
155b5b17502SBarry Smith     for (i=0; i<a->mbs; i++) a->idiag[i] = 1.0/a->a[a->i[i]];
156b5b17502SBarry Smith     a->idiagvalid = PETSC_TRUE;
157b5b17502SBarry Smith   }
158b5b17502SBarry Smith 
15941f059aeSBarry Smith   if (!a->sor_work) {
16041f059aeSBarry Smith     ierr = PetscMalloc(m*sizeof(PetscScalar),&a->sor_work);CHKERRQ(ierr);
161b5b17502SBarry Smith   }
16241f059aeSBarry Smith   t = a->sor_work;
163b5b17502SBarry Smith 
164b5b17502SBarry Smith   aidiag = a->idiag;
165b5b17502SBarry Smith 
166b5b17502SBarry Smith   if (flag == SOR_APPLY_UPPER) {
167b5b17502SBarry Smith     /* apply (U + D/omega) to the vector */
168b5b17502SBarry Smith     PetscScalar d;
169b5b17502SBarry Smith     for (i=0; i<m; i++) {
170b5b17502SBarry Smith       d    = fshift + aa[ai[i]];
171b5b17502SBarry Smith       nz   = ai[i+1] - ai[i] - 1;
172b5b17502SBarry Smith       vj   = aj + ai[i] + 1;
173b5b17502SBarry Smith       v    = aa + ai[i] + 1;
174b5b17502SBarry Smith       sum  = b[i]*d/omega;
175b5b17502SBarry Smith       PetscSparseDensePlusDot(sum,b,v,vj,nz);
176b5b17502SBarry Smith       x[i] = sum;
177b5b17502SBarry Smith     }
178b5b17502SBarry Smith     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
179b5b17502SBarry Smith   }
180b5b17502SBarry Smith 
181b5b17502SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
182b5b17502SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
183b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
184b5b17502SBarry Smith 
185b5b17502SBarry Smith       v  = aa + 1;
186b5b17502SBarry Smith       vj = aj + 1;
187b5b17502SBarry Smith       for (i=0; i<m; i++){
188b5b17502SBarry Smith         nz = ai[i+1] - ai[i] - 1;
189b5b17502SBarry Smith         tmp = - (x[i] = omega*t[i]*aidiag[i]);
190b5b17502SBarry Smith         for (j=0; j<nz; j++) {
191b5b17502SBarry Smith           t[vj[j]] += tmp*v[j];
192b5b17502SBarry Smith         }
193b5b17502SBarry Smith         v  += nz + 1;
194b5b17502SBarry Smith         vj += nz + 1;
195b5b17502SBarry Smith       }
196b5b17502SBarry Smith       ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
197b5b17502SBarry Smith     }
198b5b17502SBarry Smith 
199b5b17502SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
200b5b17502SBarry Smith       int nz2;
201b5b17502SBarry Smith       if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)){
202b5b17502SBarry Smith #if defined(PETSC_USE_BACKWARD_LOOP)
203b5b17502SBarry Smith 	v  = aa + ai[m] - 1;
204b5b17502SBarry Smith 	vj = aj + ai[m] - 1;
205b5b17502SBarry Smith 	for (i=m-1; i>=0; i--){
206b5b17502SBarry Smith           sum = b[i];
207b5b17502SBarry Smith 	  nz  = ai[i+1] - ai[i] - 1;
208b5b17502SBarry Smith           {PetscInt __i;for(__i=0;__i<nz;__i++) sum -= v[-__i] * x[vj[-__i]];}
209b5b17502SBarry Smith #else
210b5b17502SBarry Smith 	v  = aa + ai[m-1] + 1;
211b5b17502SBarry Smith 	vj = aj + ai[m-1] + 1;
212b5b17502SBarry Smith 	nz = 0;
213b5b17502SBarry Smith 	for (i=m-1; i>=0; i--){
214b5b17502SBarry Smith           sum = b[i];
215b5b17502SBarry Smith 	  nz2 = ai[i] - ai[i-1] - 1;
216b5b17502SBarry Smith 	  PETSC_Prefetch(v-nz2-1,0,1);
217b5b17502SBarry Smith 	  PETSC_Prefetch(vj-nz2-1,0,1);
218b5b17502SBarry Smith 	  PetscSparseDenseMinusDot(sum,x,v,vj,nz);
219b5b17502SBarry Smith           nz   = nz2;
220b5b17502SBarry Smith #endif
221b5b17502SBarry Smith           x[i] = omega*sum*aidiag[i];
222b5b17502SBarry Smith 	  v  -= nz + 1;
223b5b17502SBarry Smith 	  vj -= nz + 1;
224b5b17502SBarry Smith 	}
225b5b17502SBarry Smith 	ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
226b5b17502SBarry Smith       } else {
227b5b17502SBarry Smith         v  = aa + ai[m-1] + 1;
228b5b17502SBarry Smith 	vj = aj + ai[m-1] + 1;
229b5b17502SBarry Smith 	nz = 0;
230b5b17502SBarry Smith 	for (i=m-1; i>=0; i--){
231b5b17502SBarry Smith           sum = t[i];
232b5b17502SBarry Smith 	  nz2 = ai[i] - ai[i-1] - 1;
233b5b17502SBarry Smith 	  PETSC_Prefetch(v-nz2-1,0,1);
234b5b17502SBarry Smith 	  PETSC_Prefetch(vj-nz2-1,0,1);
235b5b17502SBarry Smith 	  PetscSparseDenseMinusDot(sum,x,v,vj,nz);
236b5b17502SBarry Smith           x[i] = (1-omega)*x[i] + omega*sum*aidiag[i];
237b5b17502SBarry Smith 	  nz  = nz2;
238b5b17502SBarry Smith 	  v  -= nz + 1;
239b5b17502SBarry Smith 	  vj -= nz + 1;
240b5b17502SBarry Smith 	}
241b5b17502SBarry Smith 	ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
242b5b17502SBarry Smith       }
243b5b17502SBarry Smith     }
244b5b17502SBarry Smith     its--;
245b5b17502SBarry Smith   }
246b5b17502SBarry Smith 
247b5b17502SBarry Smith   while (its--) {
248b5b17502SBarry Smith     /*
249b5b17502SBarry Smith        forward sweep:
250b5b17502SBarry Smith        for i=0,...,m-1:
251b5b17502SBarry Smith          sum[i] = (b[i] - U(i,:)x )/d[i];
252b5b17502SBarry Smith          x[i]   = (1-omega)x[i] + omega*sum[i];
253b5b17502SBarry Smith          b      = b - x[i]*U^T(i,:);
254b5b17502SBarry Smith 
255b5b17502SBarry Smith     */
256b5b17502SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
257b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
258b5b17502SBarry Smith 
259b5b17502SBarry Smith       for (i=0; i<m; i++){
260b5b17502SBarry Smith         v  = aa + ai[i] + 1; v1=v;
261b5b17502SBarry Smith         vj = aj + ai[i] + 1; vj1=vj;
262b5b17502SBarry Smith         nz = ai[i+1] - ai[i] - 1; nz1=nz;
263b5b17502SBarry Smith         sum = t[i];
264b5b17502SBarry Smith         ierr = PetscLogFlops(4.0*nz-2);CHKERRQ(ierr);
265b5b17502SBarry Smith         while (nz1--) sum -= (*v1++)*x[*vj1++];
266b5b17502SBarry Smith         x[i] = (1-omega)*x[i] + omega*sum*aidiag[i];
267b5b17502SBarry Smith         while (nz--) t[*vj++] -= x[i]*(*v++);
268b5b17502SBarry Smith       }
269b5b17502SBarry Smith     }
270b5b17502SBarry Smith 
271b5b17502SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
272b5b17502SBarry Smith       /*
273b5b17502SBarry Smith        backward sweep:
274b5b17502SBarry Smith        b = b - x[i]*U^T(i,:), i=0,...,n-2
275b5b17502SBarry Smith        for i=m-1,...,0:
276b5b17502SBarry Smith          sum[i] = (b[i] - U(i,:)x )/d[i];
277b5b17502SBarry Smith          x[i]   = (1-omega)x[i] + omega*sum[i];
278b5b17502SBarry Smith       */
279b5b17502SBarry Smith       /* if there was a forward sweep done above then I thing the next two for loops are not needed */
280b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
281b5b17502SBarry Smith 
282b5b17502SBarry Smith       for (i=0; i<m-1; i++){  /* update rhs */
283b5b17502SBarry Smith         v  = aa + ai[i] + 1;
284b5b17502SBarry Smith         vj = aj + ai[i] + 1;
285b5b17502SBarry Smith         nz = ai[i+1] - ai[i] - 1;
286b5b17502SBarry Smith         ierr = PetscLogFlops(2.0*nz-1);CHKERRQ(ierr);
287b5b17502SBarry Smith         while (nz--) t[*vj++] -= x[i]*(*v++);
288b5b17502SBarry Smith       }
289b5b17502SBarry Smith       for (i=m-1; i>=0; i--){
290b5b17502SBarry Smith         v  = aa + ai[i] + 1;
291b5b17502SBarry Smith         vj = aj + ai[i] + 1;
292b5b17502SBarry Smith         nz = ai[i+1] - ai[i] - 1;
293b5b17502SBarry Smith         ierr = PetscLogFlops(2.0*nz-1);CHKERRQ(ierr);
294b5b17502SBarry Smith         sum = t[i];
295b5b17502SBarry Smith         while (nz--) sum -= x[*vj++]*(*v++);
296b5b17502SBarry Smith         x[i] =   (1-omega)*x[i] + omega*sum*aidiag[i];
297b5b17502SBarry Smith       }
298b5b17502SBarry Smith     }
299b5b17502SBarry Smith   }
300b5b17502SBarry Smith 
301b5b17502SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
302b5b17502SBarry Smith   if (bb != xx) {
303b5b17502SBarry Smith     ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);
304b5b17502SBarry Smith   }
305b5b17502SBarry Smith   PetscFunctionReturn(0);
306b5b17502SBarry Smith }
307