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