Actual source code: veccomp.c

slepc-3.18.0 2022-10-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: #include <slepc/private/vecimplslepc.h>

 13: /* Private MPI datatypes and operators */
 14: static MPI_Datatype MPIU_NORM2=0, MPIU_NORM1_AND_2=0;
 15: static PetscBool VecCompInitialized = PETSC_FALSE;
 16: MPI_Op MPIU_NORM2_SUM=0;

 18: /* Private functions */
 19: static inline void SumNorm2(PetscReal*,PetscReal*,PetscReal*,PetscReal*);
 20: static inline PetscReal GetNorm2(PetscReal,PetscReal);
 21: static inline void AddNorm2(PetscReal*,PetscReal*,PetscReal);
 22: static PetscErrorCode VecCompSetSubVecs_Comp(Vec,PetscInt,Vec*);
 23: static PetscErrorCode VecCompGetSubVecs_Comp(Vec,PetscInt*,const Vec**);

 25: #include "veccomp0.h"

 28: #include "veccomp0.h"

 30: static inline void SumNorm2(PetscReal *ssq0,PetscReal *scale0,PetscReal *ssq1,PetscReal *scale1)
 31: {
 32:   PetscReal q;
 33:   if (*scale0 > *scale1) {
 34:     q = *scale1/(*scale0);
 35:     *ssq1 = *ssq0 + q*q*(*ssq1);
 36:     *scale1 = *scale0;
 37:   } else {
 38:     q = *scale0/(*scale1);
 39:     *ssq1 += q*q*(*ssq0);
 40:   }
 41: }

 43: static inline PetscReal GetNorm2(PetscReal ssq,PetscReal scale)
 44: {
 45:   return scale*PetscSqrtReal(ssq);
 46: }

 48: static inline void AddNorm2(PetscReal *ssq,PetscReal *scale,PetscReal x)
 49: {
 50:   PetscReal absx,q;
 51:   if (x != 0.0) {
 52:     absx = PetscAbs(x);
 53:     if (*scale < absx) {
 54:       q = *scale/absx;
 55:       *ssq = 1.0 + *ssq*q*q;
 56:       *scale = absx;
 57:     } else {
 58:       q = absx/(*scale);
 59:       *ssq += q*q;
 60:     }
 61:   }
 62: }

 64: SLEPC_EXTERN void MPIAPI SlepcSumNorm2_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
 65: {
 66:   PetscInt  i,count = *cnt;
 67:   PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;

 69:   if (*datatype == MPIU_NORM2) {
 70:     for (i=0;i<count;i++) {
 71:       SumNorm2(&xin[i*2],&xin[i*2+1],&xout[i*2],&xout[i*2+1]);
 72:     }
 73:   } else if (*datatype == MPIU_NORM1_AND_2) {
 74:     for (i=0;i<count;i++) {
 75:       xout[i*3] += xin[i*3];
 76:       SumNorm2(&xin[i*3+1],&xin[i*3+2],&xout[i*3+1],&xout[i*3+2]);
 77:     }
 78:   } else {
 79:     (*PetscErrorPrintf)("Can only handle MPIU_NORM* data types");
 80:     MPI_Abort(PETSC_COMM_WORLD,1);
 81:   }
 82:   return;
 83: }

 85: static PetscErrorCode VecCompNormEnd(void)
 86: {
 87:   MPI_Type_free(&MPIU_NORM2);
 88:   MPI_Type_free(&MPIU_NORM1_AND_2);
 89:   MPI_Op_free(&MPIU_NORM2_SUM);
 90:   VecCompInitialized = PETSC_FALSE;
 91:   return 0;
 92: }

 94: static PetscErrorCode VecCompNormInit(void)
 95: {
 96:   MPI_Type_contiguous(2,MPIU_REAL,&MPIU_NORM2);
 97:   MPI_Type_commit(&MPIU_NORM2);
 98:   MPI_Type_contiguous(3,MPIU_REAL,&MPIU_NORM1_AND_2);
 99:   MPI_Type_commit(&MPIU_NORM1_AND_2);
100:   MPI_Op_create(SlepcSumNorm2_Local,PETSC_TRUE,&MPIU_NORM2_SUM);
101:   PetscRegisterFinalize(VecCompNormEnd);
102:   return 0;
103: }

105: PetscErrorCode VecDestroy_Comp(Vec v)
106: {
107:   Vec_Comp       *vs = (Vec_Comp*)v->data;
108:   PetscInt       i;

110: #if defined(PETSC_USE_LOG)
111:   PetscLogObjectState((PetscObject)v,"Length=%" PetscInt_FMT,v->map->n);
112: #endif
113:   for (i=0;i<vs->nx;i++) VecDestroy(&vs->x[i]);
114:   if (--vs->n->friends <= 0) PetscFree(vs->n);
115:   PetscFree(vs->x);
116:   PetscFree(vs);
117:   PetscObjectComposeFunction((PetscObject)v,"VecCompSetSubVecs_C",NULL);
118:   PetscObjectComposeFunction((PetscObject)v,"VecCompGetSubVecs_C",NULL);
119:   return 0;
120: }

122: static struct _VecOps DvOps = {
123:   PetscDesignatedInitializer(duplicate,VecDuplicate_Comp),
124:   PetscDesignatedInitializer(duplicatevecs,VecDuplicateVecs_Comp),
125:   PetscDesignatedInitializer(destroyvecs,VecDestroyVecs_Comp),
126:   PetscDesignatedInitializer(dot,VecDot_Comp_MPI),
127:   PetscDesignatedInitializer(mdot,VecMDot_Comp_MPI),
128:   PetscDesignatedInitializer(norm,VecNorm_Comp_MPI),
129:   PetscDesignatedInitializer(tdot,VecTDot_Comp_MPI),
130:   PetscDesignatedInitializer(mtdot,VecMTDot_Comp_MPI),
131:   PetscDesignatedInitializer(scale,VecScale_Comp),
132:   PetscDesignatedInitializer(copy,VecCopy_Comp),
133:   PetscDesignatedInitializer(set,VecSet_Comp),
134:   PetscDesignatedInitializer(swap,VecSwap_Comp),
135:   PetscDesignatedInitializer(axpy,VecAXPY_Comp),
136:   PetscDesignatedInitializer(axpby,VecAXPBY_Comp),
137:   PetscDesignatedInitializer(maxpy,VecMAXPY_Comp),
138:   PetscDesignatedInitializer(aypx,VecAYPX_Comp),
139:   PetscDesignatedInitializer(waxpy,VecWAXPY_Comp),
140:   PetscDesignatedInitializer(axpbypcz,VecAXPBYPCZ_Comp),
141:   PetscDesignatedInitializer(pointwisemult,VecPointwiseMult_Comp),
142:   PetscDesignatedInitializer(pointwisedivide,VecPointwiseDivide_Comp),
143:   PetscDesignatedInitializer(setvalues,NULL),
144:   PetscDesignatedInitializer(assemblybegin,NULL),
145:   PetscDesignatedInitializer(assemblyend,NULL),
146:   PetscDesignatedInitializer(getarray,NULL),
147:   PetscDesignatedInitializer(getsize,VecGetSize_Comp),
148:   PetscDesignatedInitializer(getlocalsize,VecGetLocalSize_Comp),
149:   PetscDesignatedInitializer(restorearray,NULL),
150:   PetscDesignatedInitializer(max,VecMax_Comp),
151:   PetscDesignatedInitializer(min,VecMin_Comp),
152:   PetscDesignatedInitializer(setrandom,VecSetRandom_Comp),
153:   PetscDesignatedInitializer(setoption,NULL),
154:   PetscDesignatedInitializer(setvaluesblocked,NULL),
155:   PetscDesignatedInitializer(destroy,VecDestroy_Comp),
156:   PetscDesignatedInitializer(view,VecView_Comp),
157:   PetscDesignatedInitializer(placearray,NULL),
158:   PetscDesignatedInitializer(replacearray,NULL),
159:   PetscDesignatedInitializer(dot_local,VecDot_Comp_Seq),
160:   PetscDesignatedInitializer(tdot_local,VecTDot_Comp_Seq),
161:   PetscDesignatedInitializer(norm_local,VecNorm_Comp_Seq),
162:   PetscDesignatedInitializer(mdot_local,VecMDot_Comp_Seq),
163:   PetscDesignatedInitializer(mtdot_local,VecMTDot_Comp_Seq),
164:   PetscDesignatedInitializer(load,NULL),
165:   PetscDesignatedInitializer(reciprocal,VecReciprocal_Comp),
166:   PetscDesignatedInitializer(conjugate,VecConjugate_Comp),
167:   PetscDesignatedInitializer(setlocaltoglobalmapping,NULL),
168:   PetscDesignatedInitializer(setvalueslocal,NULL),
169:   PetscDesignatedInitializer(resetarray,NULL),
170:   PetscDesignatedInitializer(setfromoptions,NULL),
171:   PetscDesignatedInitializer(maxpointwisedivide,VecMaxPointwiseDivide_Comp),
172:   PetscDesignatedInitializer(pointwisemax,VecPointwiseMax_Comp),
173:   PetscDesignatedInitializer(pointwisemaxabs,VecPointwiseMaxAbs_Comp),
174:   PetscDesignatedInitializer(pointwisemin,VecPointwiseMin_Comp),
175:   PetscDesignatedInitializer(getvalues,NULL),
176:   PetscDesignatedInitializer(sqrt,VecSqrtAbs_Comp),
177:   PetscDesignatedInitializer(abs,VecAbs_Comp),
178:   PetscDesignatedInitializer(exp,VecExp_Comp),
179:   PetscDesignatedInitializer(log,VecLog_Comp),
180:   PetscDesignatedInitializer(shift,VecShift_Comp),
181:   PetscDesignatedInitializer(create,NULL),
182:   PetscDesignatedInitializer(stridegather,NULL),
183:   PetscDesignatedInitializer(stridescatter,NULL),
184:   PetscDesignatedInitializer(dotnorm2,VecDotNorm2_Comp_MPI),
185:   PetscDesignatedInitializer(getsubvector,NULL),
186:   PetscDesignatedInitializer(restoresubvector,NULL),
187:   PetscDesignatedInitializer(getarrayread,NULL),
188:   PetscDesignatedInitializer(restorearrayread,NULL),
189:   PetscDesignatedInitializer(stridesubsetgather,NULL),
190:   PetscDesignatedInitializer(stridesubsetscatter,NULL),
191:   PetscDesignatedInitializer(viewnative,NULL),
192:   PetscDesignatedInitializer(loadnative,NULL),
193:   PetscDesignatedInitializer(getlocalvector,NULL)
194: };

196: PetscErrorCode VecDuplicateVecs_Comp(Vec w,PetscInt m,Vec *V[])
197: {
198:   PetscInt       i;

203:   PetscMalloc1(m,V);
204:   for (i=0;i<m;i++) VecDuplicate(w,*V+i);
205:   return 0;
206: }

208: PetscErrorCode VecDestroyVecs_Comp(PetscInt m,Vec v[])
209: {
210:   PetscInt       i;

214:   for (i=0;i<m;i++) VecDestroy(&v[i]);
215:   PetscFree(v);
216:   return 0;
217: }

219: static PetscErrorCode VecCreate_Comp_Private(Vec v,Vec *x,PetscInt nx,PetscBool x_to_me,Vec_Comp_N *n)
220: {
221:   Vec_Comp       *s;
222:   PetscInt       N=0,lN=0,i,k;

224:   if (!VecCompInitialized) {
225:     VecCompInitialized = PETSC_TRUE;
226:     VecRegister(VECCOMP,VecCreate_Comp);
227:     VecCompNormInit();
228:   }

230:   /* Allocate a new Vec_Comp */
231:   if (v->data) PetscFree(v->data);
232:   PetscNew(&s);
233:   PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));
234:   v->data        = (void*)s;
235:   v->petscnative = PETSC_FALSE;

237:   /* Allocate the array of Vec, if it is needed to be done */
238:   if (!x_to_me) {
239:     if (nx) PetscMalloc1(nx,&s->x);
240:     if (x) PetscArraycpy(s->x,x,nx);
241:   } else s->x = x;

243:   s->nx = nx;

245:   if (nx && x) {
246:     /* Allocate the shared structure, if it is not given */
247:     if (!n) {
248:       for (i=0;i<nx;i++) {
249:         VecGetSize(x[i],&k);
250:         N+= k;
251:         VecGetLocalSize(x[i],&k);
252:         lN+= k;
253:       }
254:       PetscNew(&n);
255:       s->n = n;
256:       n->n = nx;
257:       n->N = N;
258:       n->lN = lN;
259:       n->friends = 1;
260:     } else { /* If not, check in the vector in the shared structure */
261:       s->n = n;
262:       s->n->friends++;
263:     }

265:     /* Set the virtual sizes as the real sizes of the vector */
266:     VecSetSizes(v,s->n->lN,s->n->N);
267:   }

269:   PetscObjectChangeTypeName((PetscObject)v,VECCOMP);
270:   PetscObjectComposeFunction((PetscObject)v,"VecCompSetSubVecs_C",VecCompSetSubVecs_Comp);
271:   PetscObjectComposeFunction((PetscObject)v,"VecCompGetSubVecs_C",VecCompGetSubVecs_Comp);
272:   return 0;
273: }

275: SLEPC_EXTERN PetscErrorCode VecCreate_Comp(Vec V)
276: {
277:   VecCreate_Comp_Private(V,NULL,0,PETSC_FALSE,NULL);
278:   return 0;
279: }

281: /*@C
282:    VecCreateComp - Creates a new vector containing several subvectors,
283:    each stored separately.

285:    Collective

287:    Input Parameters:
288: +  comm - communicator for the new Vec
289: .  Nx   - array of (initial) global sizes of child vectors
290: .  n    - number of child vectors
291: .  t    - type of the child vectors
292: -  Vparent - (optional) template vector

294:    Output Parameter:
295: .  V - new vector

297:    Notes:
298:    This is similar to PETSc's VecNest but customized for SLEPc's needs. In particular,
299:    the number of child vectors can be modified dynamically, with VecCompSetSubVecs().

301:    Level: developer

303: .seealso: VecCreateCompWithVecs(), VecCompSetSubVecs()
304: @*/
305: PetscErrorCode VecCreateComp(MPI_Comm comm,PetscInt *Nx,PetscInt n,VecType t,Vec Vparent,Vec *V)
306: {
307:   Vec            *x;
308:   PetscInt       i;

310:   VecCreate(comm,V);
311:   PetscMalloc1(n,&x);
312:   for (i=0;i<n;i++) {
313:     VecCreate(comm,&x[i]);
314:     VecSetSizes(x[i],PETSC_DECIDE,Nx[i]);
315:     VecSetType(x[i],t);
316:   }
317:   VecCreate_Comp_Private(*V,x,n,PETSC_TRUE,Vparent?((Vec_Comp*)Vparent->data)->n:NULL);
318:   return 0;
319: }

321: /*@C
322:    VecCreateCompWithVecs - Creates a new vector containing several subvectors,
323:    each stored separately, from an array of Vecs.

325:    Collective on x

327:    Input Parameters:
328: +  x - array of Vecs
329: .  n - number of child vectors
330: -  Vparent - (optional) template vector

332:    Output Parameter:
333: .  V - new vector

335:    Level: developer

337: .seealso: VecCreateComp()
338: @*/
339: PetscErrorCode VecCreateCompWithVecs(Vec *x,PetscInt n,Vec Vparent,Vec *V)
340: {
341:   PetscInt       i;

346:   VecCreate(PetscObjectComm((PetscObject)x[0]),V);
347:   for (i=0;i<n;i++) PetscObjectReference((PetscObject)x[i]);
348:   VecCreate_Comp_Private(*V,x,n,PETSC_FALSE,Vparent?((Vec_Comp*)Vparent->data)->n:NULL);
349:   return 0;
350: }

352: PetscErrorCode VecDuplicate_Comp(Vec win,Vec *V)
353: {
354:   Vec            *x;
355:   PetscInt       i;
356:   Vec_Comp       *s = (Vec_Comp*)win->data;

358:   SlepcValidVecComp(win,1);
359:   VecCreate(PetscObjectComm((PetscObject)win),V);
360:   PetscMalloc1(s->nx,&x);
361:   for (i=0;i<s->nx;i++) {
362:     if (s->x[i]) VecDuplicate(s->x[i],&x[i]);
363:     else x[i] = NULL;
364:   }
365:   VecCreate_Comp_Private(*V,x,s->nx,PETSC_TRUE,s->n);
366:   return 0;
367: }

369: static PetscErrorCode VecCompGetSubVecs_Comp(Vec win,PetscInt *n,const Vec **x)
370: {
371:   Vec_Comp *s = (Vec_Comp*)win->data;

373:   if (x) *x = s->x;
374:   if (n) *n = s->n->n;
375:   return 0;
376: }

378: /*@C
379:    VecCompGetSubVecs - Returns the entire array of vectors defining a
380:    compound vector.

382:    Collective on win

384:    Input Parameter:
385: .  win - compound vector

387:    Output Parameters:
388: +  n - number of child vectors
389: -  x - array of child vectors

391:    Level: developer

393: .seealso: VecCreateComp()
394: @*/
395: PetscErrorCode VecCompGetSubVecs(Vec win,PetscInt *n,const Vec **x)
396: {
398:   PetscUseMethod(win,"VecCompGetSubVecs_C",(Vec,PetscInt*,const Vec**),(win,n,x));
399:   return 0;
400: }

402: static PetscErrorCode VecCompSetSubVecs_Comp(Vec win,PetscInt n,Vec *x)
403: {
404:   Vec_Comp       *s = (Vec_Comp*)win->data;
405:   PetscInt       i,N,nlocal;
406:   Vec_Comp_N     *nn;

409:   if (!s->nx) {
410:     /* vector has been created via VecCreate+VecSetType+VecSetSizes, so allocate data structures */
411:     PetscMalloc1(n,&s->x);
412:     VecGetSize(win,&N);
414:     VecGetLocalSize(win,&nlocal);
416:     s->nx = n;
417:     for (i=0;i<n;i++) {
418:       VecCreate(PetscObjectComm((PetscObject)win),&s->x[i]);
419:       VecSetSizes(s->x[i],nlocal/n,N/n);
420:       VecSetFromOptions(s->x[i]);
421:     }
422:     if (!s->n) {
423:       PetscNew(&nn);
424:       s->n = nn;
425:       nn->N = N;
426:       nn->lN = nlocal;
427:       nn->friends = 1;
428:     }
430:   if (x) PetscArraycpy(s->x,x,n);
431:   s->n->n = n;
432:   return 0;
433: }

435: /*@C
436:    VecCompSetSubVecs - Resets the number of subvectors defining a compound vector,
437:    or replaces the subvectors.

439:    Collective on win

441:    Input Parameters:
442: +  win - compound vector
443: .  n - number of child vectors
444: -  x - array of child vectors

446:    Note:
447:    It is not possible to increase the number of subvectors with respect to the
448:    number set at its creation.

450:    Level: developer

452: .seealso: VecCreateComp(), VecCompGetSubVecs()
453: @*/
454: PetscErrorCode VecCompSetSubVecs(Vec win,PetscInt n,Vec *x)
455: {
458:   PetscTryMethod(win,"VecCompSetSubVecs_C",(Vec,PetscInt,Vec*),(win,n,x));
459:   return 0;
460: }

462: PetscErrorCode VecAXPY_Comp(Vec v,PetscScalar alpha,Vec w)
463: {
464:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
465:   PetscInt       i;

467:   SlepcValidVecComp(v,1);
468:   SlepcValidVecComp(w,3);
469:   for (i=0;i<vs->n->n;i++) VecAXPY(vs->x[i],alpha,ws->x[i]);
470:   return 0;
471: }

473: PetscErrorCode VecAYPX_Comp(Vec v,PetscScalar alpha,Vec w)
474: {
475:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
476:   PetscInt       i;

478:   SlepcValidVecComp(v,1);
479:   SlepcValidVecComp(w,3);
480:   for (i=0;i<vs->n->n;i++) VecAYPX(vs->x[i],alpha,ws->x[i]);
481:   return 0;
482: }

484: PetscErrorCode VecAXPBY_Comp(Vec v,PetscScalar alpha,PetscScalar beta,Vec w)
485: {
486:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
487:   PetscInt       i;

489:   SlepcValidVecComp(v,1);
490:   SlepcValidVecComp(w,4);
491:   for (i=0;i<vs->n->n;i++) VecAXPBY(vs->x[i],alpha,beta,ws->x[i]);
492:   return 0;
493: }

495: PetscErrorCode VecMAXPY_Comp(Vec v,PetscInt n,const PetscScalar *alpha,Vec *w)
496: {
497:   Vec_Comp       *vs = (Vec_Comp*)v->data;
498:   Vec            *wx;
499:   PetscInt       i,j;

501:   SlepcValidVecComp(v,1);
502:   for (i=0;i<n;i++) SlepcValidVecComp(w[i],4);

504:   PetscMalloc1(n,&wx);

506:   for (j=0;j<vs->n->n;j++) {
507:     for (i=0;i<n;i++) wx[i] = ((Vec_Comp*)w[i]->data)->x[j];
508:     VecMAXPY(vs->x[j],n,alpha,wx);
509:   }

511:   PetscFree(wx);
512:   return 0;
513: }

515: PetscErrorCode VecWAXPY_Comp(Vec v,PetscScalar alpha,Vec w,Vec z)
516: {
517:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data,*zs = (Vec_Comp*)z->data;
518:   PetscInt       i;

520:   SlepcValidVecComp(v,1);
521:   SlepcValidVecComp(w,3);
522:   SlepcValidVecComp(z,4);
523:   for (i=0;i<vs->n->n;i++) VecWAXPY(vs->x[i],alpha,ws->x[i],zs->x[i]);
524:   return 0;
525: }

527: PetscErrorCode VecAXPBYPCZ_Comp(Vec v,PetscScalar alpha,PetscScalar beta,PetscScalar gamma,Vec w,Vec z)
528: {
529:   Vec_Comp        *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data,*zs = (Vec_Comp*)z->data;
530:   PetscInt        i;

532:   SlepcValidVecComp(v,1);
533:   SlepcValidVecComp(w,5);
534:   SlepcValidVecComp(z,6);
535:   for (i=0;i<vs->n->n;i++) VecAXPBYPCZ(vs->x[i],alpha,beta,gamma,ws->x[i],zs->x[i]);
536:   return 0;
537: }

539: PetscErrorCode VecGetSize_Comp(Vec v,PetscInt *size)
540: {
541:   Vec_Comp *vs = (Vec_Comp*)v->data;

544:   if (vs->n) {
545:     SlepcValidVecComp(v,1);
546:     *size = vs->n->N;
547:   } else *size = v->map->N;
548:   return 0;
549: }

551: PetscErrorCode VecGetLocalSize_Comp(Vec v,PetscInt *size)
552: {
553:   Vec_Comp *vs = (Vec_Comp*)v->data;

556:   if (vs->n) {
557:     SlepcValidVecComp(v,1);
558:     *size = vs->n->lN;
559:   } else *size = v->map->n;
560:   return 0;
561: }

563: PetscErrorCode VecMax_Comp(Vec v,PetscInt *idx,PetscReal *z)
564: {
565:   Vec_Comp       *vs = (Vec_Comp*)v->data;
566:   PetscInt       idxp,s=0,s0;
567:   PetscReal      zp,z0;
568:   PetscInt       i;

570:   SlepcValidVecComp(v,1);
571:   if (!idx && !z) return 0;

573:   if (vs->n->n > 0) VecMax(vs->x[0],idx?&idxp:NULL,&zp);
574:   else {
575:     zp = PETSC_MIN_REAL;
576:     if (idx) idxp = -1;
577:   }
578:   for (i=1;i<vs->n->n;i++) {
579:     VecGetSize(vs->x[i-1],&s0);
580:     s += s0;
581:     VecMax(vs->x[i],idx?&idxp:NULL,&z0);
582:     if (zp < z0) {
583:       if (idx) *idx = s+idxp;
584:       zp = z0;
585:     }
586:   }
587:   if (z) *z = zp;
588:   return 0;
589: }

591: PetscErrorCode VecMin_Comp(Vec v,PetscInt *idx,PetscReal *z)
592: {
593:   Vec_Comp       *vs = (Vec_Comp*)v->data;
594:   PetscInt       idxp,s=0,s0;
595:   PetscReal      zp,z0;
596:   PetscInt       i;

598:   SlepcValidVecComp(v,1);
599:   if (!idx && !z) return 0;

601:   if (vs->n->n > 0) VecMin(vs->x[0],idx?&idxp:NULL,&zp);
602:   else {
603:     zp = PETSC_MAX_REAL;
604:     if (idx) idxp = -1;
605:   }
606:   for (i=1;i<vs->n->n;i++) {
607:     VecGetSize(vs->x[i-1],&s0);
608:     s += s0;
609:     VecMin(vs->x[i],idx?&idxp:NULL,&z0);
610:     if (zp > z0) {
611:       if (idx) *idx = s+idxp;
612:       zp = z0;
613:     }
614:   }
615:   if (z) *z = zp;
616:   return 0;
617: }

619: PetscErrorCode VecMaxPointwiseDivide_Comp(Vec v,Vec w,PetscReal *m)
620: {
621:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
622:   PetscReal      work;
623:   PetscInt       i;

625:   SlepcValidVecComp(v,1);
626:   SlepcValidVecComp(w,2);
627:   if (!m || vs->n->n == 0) return 0;
628:   VecMaxPointwiseDivide(vs->x[0],ws->x[0],m);
629:   for (i=1;i<vs->n->n;i++) {
630:     VecMaxPointwiseDivide(vs->x[i],ws->x[i],&work);
631:     *m = PetscMax(*m,work);
632:   }
633:   return 0;
634: }


641: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v) \
642: { \
643:   Vec_Comp        *vs = (Vec_Comp*)v->data; \
644:   PetscInt        i; \
645: \
646:   SlepcValidVecComp(v,1); \
647:   for (i=0;i<vs->n->n;i++) { \
648:     __COMPOSE2__(Vec,NAME)(vs->x[i]); \
649:   } \
650:   return 0;\
651: }

653: __FUNC_TEMPLATE1__(Conjugate)
654: __FUNC_TEMPLATE1__(Reciprocal)
655: __FUNC_TEMPLATE1__(SqrtAbs)
656: __FUNC_TEMPLATE1__(Abs)
657: __FUNC_TEMPLATE1__(Exp)
658: __FUNC_TEMPLATE1__(Log)

661: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v,T0 __a) \
662: { \
663:   Vec_Comp        *vs = (Vec_Comp*)v->data; \
664:   PetscInt        i; \
665: \
666:   SlepcValidVecComp(v,1); \
667:   for (i=0;i<vs->n->n;i++) { \
668:     __COMPOSE2__(Vec,NAME)(vs->x[i],__a); \
669:   } \
670:   return 0;\
671: }

673: __FUNC_TEMPLATE2__(Set,PetscScalar)
674: __FUNC_TEMPLATE2__(View,PetscViewer)
675: __FUNC_TEMPLATE2__(Scale,PetscScalar)
676: __FUNC_TEMPLATE2__(SetRandom,PetscRandom)
677: __FUNC_TEMPLATE2__(Shift,PetscScalar)

680: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v,Vec w) \
681: { \
682:   Vec_Comp        *vs = (Vec_Comp*)v->data,\
683:                   *ws = (Vec_Comp*)w->data; \
684:   PetscInt        i; \
685: \
686:   SlepcValidVecComp(v,1); \
687:   SlepcValidVecComp(w,2); \
688:   for (i=0;i<vs->n->n;i++) { \
689:     __COMPOSE2__(Vec,NAME)(vs->x[i],ws->x[i]); \
690:   } \
691:   return 0;\
692: }

694: __FUNC_TEMPLATE3__(Copy)
695: __FUNC_TEMPLATE3__(Swap)

698: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v,Vec w,Vec z) \
699: { \
700:   Vec_Comp        *vs = (Vec_Comp*)v->data, \
701:                   *ws = (Vec_Comp*)w->data, \
702:                   *zs = (Vec_Comp*)z->data; \
703:   PetscInt        i; \
704: \
705:   SlepcValidVecComp(v,1); \
706:   SlepcValidVecComp(w,2); \
707:   SlepcValidVecComp(z,3); \
708:   for (i=0;i<vs->n->n;i++) { \
709:     __COMPOSE2__(Vec,NAME)(vs->x[i],ws->x[i],zs->x[i]); \
710:   } \
711:   return 0;\
712: }

714: __FUNC_TEMPLATE4__(PointwiseMax)
715: __FUNC_TEMPLATE4__(PointwiseMaxAbs)
716: __FUNC_TEMPLATE4__(PointwiseMin)
717: __FUNC_TEMPLATE4__(PointwiseMult)
718: __FUNC_TEMPLATE4__(PointwiseDivide)