xref: /petsc/src/mat/impls/sbaij/seq/relax.h (revision 3649974ff2fe2c90f327f29aac9f68ec978ade9e)
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
28eeffb40dSHong Zhang 
29eeffb40dSHong Zhang   PetscFunctionBegin;
30eeffb40dSHong Zhang   ierr = VecSet(zz,0.0);CHKERRQ(ierr);
31*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
32eeffb40dSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
33eeffb40dSHong Zhang 
34eeffb40dSHong Zhang   v  = a->a;
35eeffb40dSHong Zhang   for (i=0; i<mbs; i++) {
36eeffb40dSHong Zhang     nz   = ai[i+1] - ai[i];  /* length of i_th row of A */
37eeffb40dSHong Zhang     x1   = x[i];
38eeffb40dSHong Zhang     sum  = v[0]*x1;          /* diagonal term */
39eeffb40dSHong Zhang     for (j=1; j<nz; j++) {
40a6f0575fSBarry Smith       ibt  = ib[j];
41eeffb40dSHong Zhang       vj   = v[j];
42eeffb40dSHong Zhang       sum += vj * x[ibt];   /* (strict upper triangular part of A)*x  */
43a6f0575fSBarry Smith       z[ibt] += PetscConj(v[j]) * x1;    /* (strict lower triangular part of A)*x  */
44eeffb40dSHong Zhang     }
45eeffb40dSHong Zhang     z[i] += sum;
46eeffb40dSHong Zhang     v    += nz;
47eeffb40dSHong Zhang     ib   += nz;
48eeffb40dSHong Zhang   }
49eeffb40dSHong Zhang 
50*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
51eeffb40dSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
52eeffb40dSHong Zhang   ierr = PetscLogFlops(2.0*(2.0*a->nz - mbs) - mbs);CHKERRQ(ierr);
53eeffb40dSHong Zhang   PetscFunctionReturn(0);
54eeffb40dSHong Zhang }
55eeffb40dSHong Zhang 
56eeffb40dSHong Zhang #undef __FUNCT__
57b5b17502SBarry Smith #define __FUNCT__ "MatMult_SeqSBAIJ_1"
58018dd85eSSatish Balay #if defined(USESHORT)
59018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_ushort(Mat A,Vec xx,Vec zz)
60018dd85eSSatish Balay #else
61018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1(Mat A,Vec xx,Vec zz)
62018dd85eSSatish Balay #endif
63b5b17502SBarry Smith {
64b5b17502SBarry Smith   Mat_SeqSBAIJ         *a = (Mat_SeqSBAIJ*)A->data;
65b5b17502SBarry Smith   const PetscScalar    *x;
66b5b17502SBarry Smith   PetscScalar          *z,x1,sum;
67b5b17502SBarry Smith   const MatScalar      *v;
68b5b17502SBarry Smith   MatScalar            vj;
69b5b17502SBarry Smith   PetscErrorCode       ierr;
70b5b17502SBarry Smith   PetscInt             mbs=a->mbs,i,j,nz;
71b5b17502SBarry Smith   const PetscInt       *ai=a->i;
72b5b17502SBarry Smith #if defined(USESHORT)
73b5b17502SBarry Smith   const unsigned short *ib=a->jshort;
74b5b17502SBarry Smith   unsigned short       ibt;
75b5b17502SBarry Smith #else
76b5b17502SBarry Smith   const PetscInt       *ib=a->j;
77b5b17502SBarry Smith   PetscInt             ibt;
78b5b17502SBarry Smith #endif
79b5b17502SBarry Smith 
80b5b17502SBarry Smith   PetscFunctionBegin;
81b5b17502SBarry Smith   ierr = VecSet(zz,0.0);CHKERRQ(ierr);
82*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
83b5b17502SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
84b5b17502SBarry Smith 
85b5b17502SBarry Smith   v  = a->a;
86b5b17502SBarry Smith   for (i=0; i<mbs; i++) {
87b5b17502SBarry Smith     nz   = ai[i+1] - ai[i];        /* length of i_th row of A */
88b5b17502SBarry Smith     x1   = x[i];
89b5b17502SBarry Smith     sum  = v[0]*x1;                /* diagonal term */
90b5b17502SBarry Smith     for (j=1; j<nz; j++) {
91b5b17502SBarry Smith       ibt = ib[j];
92b5b17502SBarry Smith       vj  = v[j];
93b5b17502SBarry Smith       z[ibt] += vj * x1;       /* (strict lower triangular part of A)*x  */
94b5b17502SBarry Smith       sum    += vj * x[ibt]; /* (strict upper triangular part of A)*x  */
95b5b17502SBarry Smith     }
96b5b17502SBarry Smith     z[i] += sum;
97b5b17502SBarry Smith     v    += nz;
98b5b17502SBarry Smith     ib   += nz;
99b5b17502SBarry Smith   }
100b5b17502SBarry Smith 
101*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
102b5b17502SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
103b5b17502SBarry Smith   ierr = PetscLogFlops(2.0*(2.0*a->nz - mbs) - mbs);CHKERRQ(ierr);
104b5b17502SBarry Smith   PetscFunctionReturn(0);
105b5b17502SBarry Smith }
106b5b17502SBarry Smith 
107b5b17502SBarry Smith #undef __FUNCT__
10841f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqSBAIJ"
109018dd85eSSatish Balay #if defined(USESHORT)
110018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ_ushort(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
111018dd85eSSatish Balay #else
112018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
113018dd85eSSatish Balay #endif
114b5b17502SBarry Smith {
115b5b17502SBarry Smith   Mat_SeqSBAIJ         *a = (Mat_SeqSBAIJ*)A->data;
116b5b17502SBarry Smith   const MatScalar      *aa=a->a,*v,*v1,*aidiag;
117b5b17502SBarry Smith   PetscScalar          *x,*t,sum;
118b5b17502SBarry Smith   const PetscScalar    *b;
119b5b17502SBarry Smith   MatScalar            tmp;
120b5b17502SBarry Smith   PetscErrorCode       ierr;
121b5b17502SBarry Smith   PetscInt             m=a->mbs,bs=A->rmap->bs,j;
122b5b17502SBarry Smith   const PetscInt       *ai=a->i;
123b5b17502SBarry Smith #if defined(USESHORT)
124b5b17502SBarry Smith   const unsigned short *aj=a->jshort,*vj,*vj1;
125b5b17502SBarry Smith #else
126b5b17502SBarry Smith   const PetscInt       *aj=a->j,*vj,*vj1;
127b5b17502SBarry Smith #endif
128b5b17502SBarry Smith   PetscInt             nz,nz1,i;
129b5b17502SBarry Smith 
130b5b17502SBarry Smith   PetscFunctionBegin;
131e32f2f54SBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
132b5b17502SBarry Smith 
133b5b17502SBarry Smith   its = its*lits;
134e32f2f54SBarry 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);
135b5b17502SBarry Smith 
136e32f2f54SBarry Smith   if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
137b5b17502SBarry Smith 
138b5b17502SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
139*3649974fSBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
140b5b17502SBarry Smith 
141b5b17502SBarry Smith   if (!a->idiagvalid) {
142b5b17502SBarry Smith     if (!a->idiag) {
143b5b17502SBarry Smith       ierr = PetscMalloc(m*sizeof(PetscScalar),&a->idiag);CHKERRQ(ierr);
144b5b17502SBarry Smith     }
145b5b17502SBarry Smith     for (i=0; i<a->mbs; i++) a->idiag[i] = 1.0/a->a[a->i[i]];
146b5b17502SBarry Smith     a->idiagvalid = PETSC_TRUE;
147b5b17502SBarry Smith   }
148b5b17502SBarry Smith 
14941f059aeSBarry Smith   if (!a->sor_work) {
15041f059aeSBarry Smith     ierr = PetscMalloc(m*sizeof(PetscScalar),&a->sor_work);CHKERRQ(ierr);
151b5b17502SBarry Smith   }
15241f059aeSBarry Smith   t = a->sor_work;
153b5b17502SBarry Smith 
154b5b17502SBarry Smith   aidiag = a->idiag;
155b5b17502SBarry Smith 
156b5b17502SBarry Smith   if (flag == SOR_APPLY_UPPER) {
157b5b17502SBarry Smith     /* apply (U + D/omega) to the vector */
158b5b17502SBarry Smith     PetscScalar d;
159b5b17502SBarry Smith     for (i=0; i<m; i++) {
160b5b17502SBarry Smith       d    = fshift + aa[ai[i]];
161b5b17502SBarry Smith       nz   = ai[i+1] - ai[i] - 1;
162b5b17502SBarry Smith       vj   = aj + ai[i] + 1;
163b5b17502SBarry Smith       v    = aa + ai[i] + 1;
164b5b17502SBarry Smith       sum  = b[i]*d/omega;
165b5b17502SBarry Smith       PetscSparseDensePlusDot(sum,b,v,vj,nz);
166b5b17502SBarry Smith       x[i] = sum;
167b5b17502SBarry Smith     }
168b5b17502SBarry Smith     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
169b5b17502SBarry Smith   }
170b5b17502SBarry Smith 
171b5b17502SBarry Smith   if (flag & SOR_ZERO_INITIAL_GUESS) {
172b5b17502SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
173b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
174b5b17502SBarry Smith 
175b5b17502SBarry Smith       v  = aa + 1;
176b5b17502SBarry Smith       vj = aj + 1;
177b5b17502SBarry Smith       for (i=0; i<m; i++){
178b5b17502SBarry Smith         nz = ai[i+1] - ai[i] - 1;
179b5b17502SBarry Smith         tmp = - (x[i] = omega*t[i]*aidiag[i]);
180b5b17502SBarry Smith         for (j=0; j<nz; j++) {
181b5b17502SBarry Smith           t[vj[j]] += tmp*v[j];
182b5b17502SBarry Smith         }
183b5b17502SBarry Smith         v  += nz + 1;
184b5b17502SBarry Smith         vj += nz + 1;
185b5b17502SBarry Smith       }
186b5b17502SBarry Smith       ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
187b5b17502SBarry Smith     }
188b5b17502SBarry Smith 
189b5b17502SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
190b5b17502SBarry Smith       int nz2;
191b5b17502SBarry Smith       if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)){
192b5b17502SBarry Smith #if defined(PETSC_USE_BACKWARD_LOOP)
193b5b17502SBarry Smith 	v  = aa + ai[m] - 1;
194b5b17502SBarry Smith 	vj = aj + ai[m] - 1;
195b5b17502SBarry Smith 	for (i=m-1; i>=0; i--){
196b5b17502SBarry Smith           sum = b[i];
197b5b17502SBarry Smith 	  nz  = ai[i+1] - ai[i] - 1;
198b5b17502SBarry Smith           {PetscInt __i;for(__i=0;__i<nz;__i++) sum -= v[-__i] * x[vj[-__i]];}
199b5b17502SBarry Smith #else
200b5b17502SBarry Smith 	v  = aa + ai[m-1] + 1;
201b5b17502SBarry Smith 	vj = aj + ai[m-1] + 1;
202b5b17502SBarry Smith 	nz = 0;
203b5b17502SBarry Smith 	for (i=m-1; i>=0; i--){
204b5b17502SBarry Smith           sum = b[i];
205b5b17502SBarry Smith 	  nz2 = ai[i] - ai[i-1] - 1;
206b5b17502SBarry Smith 	  PETSC_Prefetch(v-nz2-1,0,1);
207b5b17502SBarry Smith 	  PETSC_Prefetch(vj-nz2-1,0,1);
208b5b17502SBarry Smith 	  PetscSparseDenseMinusDot(sum,x,v,vj,nz);
209b5b17502SBarry Smith           nz   = nz2;
210b5b17502SBarry Smith #endif
211b5b17502SBarry Smith           x[i] = omega*sum*aidiag[i];
212b5b17502SBarry Smith 	  v  -= nz + 1;
213b5b17502SBarry Smith 	  vj -= nz + 1;
214b5b17502SBarry Smith 	}
215b5b17502SBarry Smith 	ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
216b5b17502SBarry Smith       } else {
217b5b17502SBarry Smith         v  = aa + ai[m-1] + 1;
218b5b17502SBarry Smith 	vj = aj + ai[m-1] + 1;
219b5b17502SBarry Smith 	nz = 0;
220b5b17502SBarry Smith 	for (i=m-1; i>=0; i--){
221b5b17502SBarry Smith           sum = t[i];
222b5b17502SBarry Smith 	  nz2 = ai[i] - ai[i-1] - 1;
223b5b17502SBarry Smith 	  PETSC_Prefetch(v-nz2-1,0,1);
224b5b17502SBarry Smith 	  PETSC_Prefetch(vj-nz2-1,0,1);
225b5b17502SBarry Smith 	  PetscSparseDenseMinusDot(sum,x,v,vj,nz);
226b5b17502SBarry Smith           x[i] = (1-omega)*x[i] + omega*sum*aidiag[i];
227b5b17502SBarry Smith 	  nz  = nz2;
228b5b17502SBarry Smith 	  v  -= nz + 1;
229b5b17502SBarry Smith 	  vj -= nz + 1;
230b5b17502SBarry Smith 	}
231b5b17502SBarry Smith 	ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
232b5b17502SBarry Smith       }
233b5b17502SBarry Smith     }
234b5b17502SBarry Smith     its--;
235b5b17502SBarry Smith   }
236b5b17502SBarry Smith 
237b5b17502SBarry Smith   while (its--) {
238b5b17502SBarry Smith     /*
239b5b17502SBarry Smith        forward sweep:
240b5b17502SBarry Smith        for i=0,...,m-1:
241b5b17502SBarry Smith          sum[i] = (b[i] - U(i,:)x )/d[i];
242b5b17502SBarry Smith          x[i]   = (1-omega)x[i] + omega*sum[i];
243b5b17502SBarry Smith          b      = b - x[i]*U^T(i,:);
244b5b17502SBarry Smith 
245b5b17502SBarry Smith     */
246b5b17502SBarry Smith     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
247b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
248b5b17502SBarry Smith 
249b5b17502SBarry Smith       for (i=0; i<m; i++){
250b5b17502SBarry Smith         v  = aa + ai[i] + 1; v1=v;
251b5b17502SBarry Smith         vj = aj + ai[i] + 1; vj1=vj;
252b5b17502SBarry Smith         nz = ai[i+1] - ai[i] - 1; nz1=nz;
253b5b17502SBarry Smith         sum = t[i];
254b5b17502SBarry Smith         ierr = PetscLogFlops(4.0*nz-2);CHKERRQ(ierr);
255b5b17502SBarry Smith         while (nz1--) sum -= (*v1++)*x[*vj1++];
256b5b17502SBarry Smith         x[i] = (1-omega)*x[i] + omega*sum*aidiag[i];
257b5b17502SBarry Smith         while (nz--) t[*vj++] -= x[i]*(*v++);
258b5b17502SBarry Smith       }
259b5b17502SBarry Smith     }
260b5b17502SBarry Smith 
261b5b17502SBarry Smith     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
262b5b17502SBarry Smith       /*
263b5b17502SBarry Smith        backward sweep:
264b5b17502SBarry Smith        b = b - x[i]*U^T(i,:), i=0,...,n-2
265b5b17502SBarry Smith        for i=m-1,...,0:
266b5b17502SBarry Smith          sum[i] = (b[i] - U(i,:)x )/d[i];
267b5b17502SBarry Smith          x[i]   = (1-omega)x[i] + omega*sum[i];
268b5b17502SBarry Smith       */
269b5b17502SBarry Smith       /* if there was a forward sweep done above then I thing the next two for loops are not needed */
270b5b17502SBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
271b5b17502SBarry Smith 
272b5b17502SBarry Smith       for (i=0; i<m-1; i++){  /* update rhs */
273b5b17502SBarry Smith         v  = aa + ai[i] + 1;
274b5b17502SBarry Smith         vj = aj + ai[i] + 1;
275b5b17502SBarry Smith         nz = ai[i+1] - ai[i] - 1;
276b5b17502SBarry Smith         ierr = PetscLogFlops(2.0*nz-1);CHKERRQ(ierr);
277b5b17502SBarry Smith         while (nz--) t[*vj++] -= x[i]*(*v++);
278b5b17502SBarry Smith       }
279b5b17502SBarry Smith       for (i=m-1; i>=0; i--){
280b5b17502SBarry Smith         v  = aa + ai[i] + 1;
281b5b17502SBarry Smith         vj = aj + ai[i] + 1;
282b5b17502SBarry Smith         nz = ai[i+1] - ai[i] - 1;
283b5b17502SBarry Smith         ierr = PetscLogFlops(2.0*nz-1);CHKERRQ(ierr);
284b5b17502SBarry Smith         sum = t[i];
285b5b17502SBarry Smith         while (nz--) sum -= x[*vj++]*(*v++);
286b5b17502SBarry Smith         x[i] =   (1-omega)*x[i] + omega*sum*aidiag[i];
287b5b17502SBarry Smith       }
288b5b17502SBarry Smith     }
289b5b17502SBarry Smith   }
290b5b17502SBarry Smith 
291b5b17502SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
292*3649974fSBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
293b5b17502SBarry Smith   PetscFunctionReturn(0);
294b5b17502SBarry Smith }
295