Actual source code: veccomp0.h

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: #include <petsc-private/vecimpl.h>     /*I  "petsvec.h"  I*/

 24: #if defined(__WITH_MPI__)
 26: #else
 28: #endif


 36: PetscErrorCode __SUF__(VecDot_Comp)(Vec a,Vec b,PetscScalar *z)
 37: {
 38:   PetscScalar    sum = 0.0,work;
 39:   PetscInt       i;
 41:   Vec_Comp       *as = (Vec_Comp*)a->data,*bs = (Vec_Comp*)b->data;

 44:   SlepcValidVecComp(a);
 45:   SlepcValidVecComp(b);
 46:   if (as->x[0]->ops->dot_local) {
 47:     for (i=0,sum=0.0;i<as->n->n;i++) {
 48:       as->x[i]->ops->dot_local(as->x[i],bs->x[i],&work);
 49:       sum += work;
 50:     }
 51: #if defined(__WITH_MPI__)
 52:     work = sum;
 53:     MPI_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)a));
 54: #endif
 55:   } else {
 56:     for (i=0,sum=0.0;i<as->n->n;i++) {
 57:       VecDot(as->x[i],bs->x[i],&work);
 58:       sum += work;
 59:     }
 60:   }
 61:   *z = sum;
 62:   return(0);
 63: }

 67: PetscErrorCode __SUF__(VecMDot_Comp)(Vec a,PetscInt n,const Vec b[],PetscScalar *z)
 68: {
 69:   PetscScalar    *work,*work0,*r;
 71:   Vec_Comp       *as = (Vec_Comp*)a->data;
 72:   Vec            *bx;
 73:   PetscInt       i,j;

 76:   SlepcValidVecComp(a);
 77:   for (i=0;i<n;i++) SlepcValidVecComp(b[i]);

 79:   if (as->n->n == 0) {
 80:     *z = 0;
 81:     return(0);
 82:   }

 84:   PetscMalloc(sizeof(PetscScalar)*n,&work0);
 85:   PetscMalloc(sizeof(Vec)*n,&bx);

 87: #if defined(__WITH_MPI__)
 88:   if (as->x[0]->ops->mdot_local) {
 89:     r = work0; work = z;
 90:   } else
 91: #endif
 92:   {
 93:     r = z; work = work0;
 94:   }

 96:   /* z[i] <- a.x' * b[i].x */
 97:   for (i=0;i<n;i++) bx[i] = ((Vec_Comp*)b[i]->data)->x[0];
 98:   if (as->x[0]->ops->mdot_local) {
 99:     as->x[0]->ops->mdot_local(as->x[0],n,bx,r);
100:   } else {
101:     VecMDot(as->x[0],n,bx,r);
102:   }
103:   for (j=0;j<as->n->n;j++) {
104:     for (i=0;i<n;i++) bx[i] = ((Vec_Comp*)b[i]->data)->x[j];
105:     if (as->x[0]->ops->mdot_local) {
106:       as->x[j]->ops->mdot_local(as->x[j],n,bx,work);
107:     } else {
108:       VecMDot(as->x[j],n,bx,work);
109:     }
110:     for (i=0;i<n;i++) r[i] += work[i];
111:   }

113:   /* If def(__WITH_MPI__) and exists mdot_local */
114: #if defined(__WITH_MPI__)
115:   if (as->x[0]->ops->mdot_local) {
116:     /* z[i] <- Allreduce(work[i]) */
117:     MPI_Allreduce(r,z,n,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)a));
118:   }
119: #endif

121:   PetscFree(work0);
122:   PetscFree(bx);
123:   return(0);
124: }

128: PetscErrorCode __SUF__(VecTDot_Comp)(Vec a,Vec b,PetscScalar *z)
129: {
130:   PetscScalar    sum = 0.0,work;
131:   PetscInt       i;
133:   Vec_Comp       *as = (Vec_Comp*)a->data,*bs = (Vec_Comp*)b->data;

136:   SlepcValidVecComp(a);
137:   SlepcValidVecComp(b);
138:   if (as->x[0]->ops->tdot_local) {
139:     for (i=0,sum=0.0;i<as->n->n;i++) {
140:       as->x[i]->ops->tdot_local(as->x[i],bs->x[i],&work);
141:       sum += work;
142:     }
143: #if defined(__WITH_MPI__)
144:     work = sum;
145:     MPI_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)a));
146: #endif
147:   } else {
148:     for (i=0,sum=0.0;i<as->n->n;i++) {
149:       VecTDot(as->x[i],bs->x[i],&work);
150:       sum += work;
151:     }
152:   }
153:   *z = sum;
154:   return(0);
155: }

159: PetscErrorCode __SUF__(VecMTDot_Comp)(Vec a,PetscInt n,const Vec b[],PetscScalar *z)
160: {
161:   PetscScalar    *work,*work0,*r;
163:   Vec_Comp       *as = (Vec_Comp*)a->data;
164:   Vec            *bx;
165:   PetscInt       i,j;

168:   SlepcValidVecComp(a);
169:   for (i=0;i<n;i++) SlepcValidVecComp(b[i]);

171:   if (as->n->n == 0) {
172:     *z = 0;
173:     return(0);
174:   }

176:   PetscMalloc(sizeof(PetscScalar)*n,&work0);
177:   PetscMalloc(sizeof(Vec)*n,&bx);

179: #if defined(__WITH_MPI__)
180:   if (as->x[0]->ops->mtdot_local) {
181:     r = work0; work = z;
182:   } else
183: #endif
184:   {
185:     r = z; work = work0;
186:   }

188:   /* z[i] <- a.x' * b[i].x */
189:   for (i=0;i<n;i++) bx[i] = ((Vec_Comp*)b[i]->data)->x[0];
190:   if (as->x[0]->ops->mtdot_local) {
191:     as->x[0]->ops->mtdot_local(as->x[0],n,bx,r);
192:   } else {
193:     VecMTDot(as->x[0],n,bx,r);
194:   }
195:   for (j=0;j<as->n->n;j++) {
196:     for (i=0;i<n;i++) bx[i] = ((Vec_Comp*)b[i]->data)->x[j];
197:     if (as->x[0]->ops->mtdot_local) {
198:       as->x[j]->ops->mtdot_local(as->x[j],n,bx,work);
199:     } else {
200:       VecMTDot(as->x[j],n,bx,work);
201:     }
202:     for (i=0;i<n;i++) r[i] += work[i];
203:   }

205:   /* If def(__WITH_MPI__) and exists mtdot_local */
206: #if defined(__WITH_MPI__)
207:   if (as->x[0]->ops->mtdot_local) {
208:     /* z[i] <- Allreduce(work[i]) */
209:     MPI_Allreduce(r,z,n,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)a));
210:   }
211: #endif

213:   PetscFree(work0);
214:   PetscFree(bx);
215:   return(0);
216: }

220: PetscErrorCode __SUF__(VecNorm_Comp)(Vec a,NormType t,PetscReal *norm)
221: {
222:   PetscReal      work[3],s=0.0;
224:   Vec_Comp       *as = (Vec_Comp*)a->data;
225:   PetscInt       i;

228:   SlepcValidVecComp(a);
229:   /* Initialize norm */
230:   switch (t) {
231:     case NORM_1: case NORM_INFINITY: *norm = 0.0; break;
232:     case NORM_2: case NORM_FROBENIUS: *norm = 1.0; s = 0.0; break;
233:     case NORM_1_AND_2: norm[0] = 0.0; norm[1] = 1.0; s = 0.0; break;
234:   }
235:   for (i=0;i<as->n->n;i++) {
236:     if (as->x[0]->ops->norm_local) {
237:       as->x[0]->ops->norm_local(as->x[i],t,work);
238:     } else {
239:       VecNorm(as->x[i],t,work);
240:     }
241:     /* norm+= work */
242:     switch (t) {
243:       case NORM_1: *norm+= *work; break;
244:       case NORM_2: case NORM_FROBENIUS: AddNorm2(norm,&s,*work); break;
245:       case NORM_1_AND_2: norm[0]+= work[0]; AddNorm2(&norm[1],&s,work[1]); break;
246:       case NORM_INFINITY: *norm = PetscMax(*norm,*work); break;
247:     }
248:   }

250:   /* If def(__WITH_MPI__) and exists norm_local */
251: #if defined(__WITH_MPI__)
252:   if (as->x[0]->ops->norm_local) {
253:     PetscReal work0[3];
254:     /* norm <- Allreduce(work) */
255:     switch (t) {
256:     case NORM_1:
257:       work[0] = *norm;
258:       MPI_Allreduce(work,norm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)a));
259:       break;
260:     case NORM_2: case NORM_FROBENIUS:
261:       work[0] = *norm; work[1] = s;
262:       MPI_Allreduce(work,work0,1,MPIU_NORM2,MPIU_NORM2_SUM,PetscObjectComm((PetscObject)a));
263:       *norm = GetNorm2(work0[0],work0[1]);
264:       break;
265:     case NORM_1_AND_2:
266:       work[0] = norm[0]; work[1] = norm[1]; work[2] = s;
267:       MPI_Allreduce(work,work0,1,MPIU_NORM1_AND_2,MPIU_NORM2_SUM,PetscObjectComm((PetscObject)a));
268:       norm[0] = work0[0];
269:       norm[1] = GetNorm2(work0[1],work0[2]);
270:       break;
271:     case NORM_INFINITY:
272:       work[0] = *norm;
273:       MPI_Allreduce(work,norm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)a));
274:       break;
275:     }
276:   }
277: #else
278:   /* Norm correction */
279:   switch (t) {
280:     case NORM_2: case NORM_FROBENIUS: *norm = GetNorm2(*norm,s); break;
281:     case NORM_1_AND_2: norm[1] = GetNorm2(norm[1],s); break;
282:     default: ;
283:   }
284: #endif
285:   return(0);
286: }

290: PetscErrorCode __SUF__(VecDotNorm2_Comp)(Vec v,Vec w,PetscScalar *dp,PetscScalar *nm)
291: {
292:   PetscScalar    *vx,*wx,dp0,nm0,dp1,nm1;
294:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
295:   PetscInt       i,n;
296:   PetscBool      t0,t1;
297: #if defined(__WITH_MPI__)
298:   PetscScalar    work[4];
299: #endif

302:   /* Compute recursively the local part */
303:   dp0 = nm0 = 0.0;
304:   PetscObjectTypeCompare((PetscObject)v,VECCOMP,&t0);
305:   PetscObjectTypeCompare((PetscObject)w,VECCOMP,&t1);
306:   if (t0 && t1) {
307:     SlepcValidVecComp(v);
308:     SlepcValidVecComp(w);
309:     for (i=0;i<vs->n->n;i++) {
310:       VecDotNorm2_Comp_Seq(vs->x[i],ws->x[i],&dp1,&nm1);
311:       dp0 += dp1;
312:       nm0 += nm1;
313:     }
314:   } else if (!t0 && !t1) {
315:     VecGetLocalSize(v,&n);
316:     VecGetArray(v,&vx);
317:     VecGetArray(w,&wx);
318:     for (i=0;i<n;i++) {
319:       dp0 += vx[i]*PetscConj(wx[i]);
320:       nm0 += wx[i]*PetscConj(wx[i]);
321:     }
322:     VecRestoreArray(v,&vx);
323:     VecRestoreArray(w,&wx);
324:   } else SETERRQ(PetscObjectComm((PetscObject)v),PETSC_ERR_ARG_INCOMP,"Incompatible vector types");

326: #if defined(__WITH_MPI__)
327:     /* [dp, nm] <- Allreduce([dp0, nm0]) */
328:     work[0] = dp0; work[1] = nm0;
329:     MPI_Allreduce(work,&work[2],2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)v));
330:     *dp = work[2]; *nm = work[3];
331: #else
332:     *dp = dp0; *nm = nm0;
333: #endif
334:   return(0);
335: }