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