Actual source code: veccomp.c

slepc-3.16.3 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: PETSC_STATIC_INLINE void SumNorm2(PetscReal*,PetscReal*,PetscReal*,PetscReal*);
 20: PETSC_STATIC_INLINE PetscReal GetNorm2(PetscReal,PetscReal);
 21: PETSC_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: PETSC_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: PETSC_STATIC_INLINE PetscReal GetNorm2(PetscReal ssq,PetscReal scale)
 44: {
 45:   return scale*PetscSqrtReal(ssq);
 46: }

 48: PETSC_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;

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

 86: static PetscErrorCode VecCompNormEnd(void)
 87: {

 91:   MPI_Type_free(&MPIU_NORM2);
 92:   MPI_Type_free(&MPIU_NORM1_AND_2);
 93:   MPI_Op_free(&MPIU_NORM2_SUM);
 94:   VecCompInitialized = PETSC_FALSE;
 95:   return(0);
 96: }

 98: static PetscErrorCode VecCompNormInit(void)
 99: {

103:   MPI_Type_contiguous(2,MPIU_REAL,&MPIU_NORM2);
104:   MPI_Type_commit(&MPIU_NORM2);
105:   MPI_Type_contiguous(3,MPIU_REAL,&MPIU_NORM1_AND_2);
106:   MPI_Type_commit(&MPIU_NORM1_AND_2);
107:   MPI_Op_create(SlepcSumNorm2_Local,PETSC_TRUE,&MPIU_NORM2_SUM);
108:   PetscRegisterFinalize(VecCompNormEnd);
109:   return(0);
110: }

112: PetscErrorCode VecDestroy_Comp(Vec v)
113: {
114:   Vec_Comp       *vs = (Vec_Comp*)v->data;
115:   PetscInt       i;

119: #if defined(PETSC_USE_LOG)
120:   PetscLogObjectState((PetscObject)v,"Length=%D",v->map->n);
121: #endif
122:   for (i=0;i<vs->nx;i++) {
123:     VecDestroy(&vs->x[i]);
124:   }
125:   if (--vs->n->friends <= 0) {
126:     PetscFree(vs->n);
127:   }
128:   PetscFree(vs->x);
129:   PetscFree(vs);
130:   PetscObjectComposeFunction((PetscObject)v,"VecCompSetSubVecs_C",NULL);
131:   PetscObjectComposeFunction((PetscObject)v,"VecCompGetSubVecs_C",NULL);
132:   return(0);
133: }

135: static struct _VecOps DvOps = {VecDuplicate_Comp, /* 1 */
136:             VecDuplicateVecs_Comp,
137:             VecDestroyVecs_Comp,
138:             VecDot_Comp_MPI,
139:             VecMDot_Comp_MPI,
140:             VecNorm_Comp_MPI,
141:             VecTDot_Comp_MPI,
142:             VecMTDot_Comp_MPI,
143:             VecScale_Comp,
144:             VecCopy_Comp, /* 10 */
145:             VecSet_Comp,
146:             VecSwap_Comp,
147:             VecAXPY_Comp,
148:             VecAXPBY_Comp,
149:             VecMAXPY_Comp,
150:             VecAYPX_Comp,
151:             VecWAXPY_Comp,
152:             VecAXPBYPCZ_Comp,
153:             VecPointwiseMult_Comp,
154:             VecPointwiseDivide_Comp,
155:             0, /* 20 */
156:             0,0,
157:             0 /*VecGetArray_Seq*/,
158:             VecGetSize_Comp,
159:             VecGetLocalSize_Comp,
160:             0/*VecRestoreArray_Seq*/,
161:             VecMax_Comp,
162:             VecMin_Comp,
163:             VecSetRandom_Comp,
164:             0, /* 30 */
165:             0,
166:             VecDestroy_Comp,
167:             VecView_Comp,
168:             0/*VecPlaceArray_Seq*/,
169:             0/*VecReplaceArray_Seq*/,
170:             VecDot_Comp_Seq,
171:             VecTDot_Comp_Seq,
172:             VecNorm_Comp_Seq,
173:             VecMDot_Comp_Seq,
174:             VecMTDot_Comp_Seq, /* 40 */
175:             0,
176:             VecReciprocal_Comp,
177:             VecConjugate_Comp,
178:             0,0,
179:             0/*VecResetArray_Seq*/,
180:             0,
181:             VecMaxPointwiseDivide_Comp,
182:             VecPointwiseMax_Comp,
183:             VecPointwiseMaxAbs_Comp,
184:             VecPointwiseMin_Comp,
185:             0,
186:             VecSqrtAbs_Comp,
187:             VecAbs_Comp,
188:             VecExp_Comp,
189:             VecLog_Comp,
190:             VecShift_Comp,
191:             0,
192:             0,
193:             0,
194:             VecDotNorm2_Comp_MPI
195:           };

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

205:   if (m<=0) SETERRQ1(PetscObjectComm((PetscObject)w),PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %D",m);
206:   PetscMalloc1(m,V);
207:   for (i=0;i<m;i++) { VecDuplicate(w,*V+i); }
208:   return(0);
209: }

211: PetscErrorCode VecDestroyVecs_Comp(PetscInt m,Vec v[])
212: {
214:   PetscInt       i;

218:   if (m<=0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %D",m);
219:   for (i=0;i<m;i++) { VecDestroy(&v[i]); }
220:   PetscFree(v);
221:   return(0);
222: }

224: static PetscErrorCode VecCreate_Comp_Private(Vec v,Vec *x,PetscInt nx,PetscBool x_to_me,Vec_Comp_N *n)
225: {
226:   Vec_Comp       *s;
228:   PetscInt       N=0,lN=0,i,k;

231:   if (!VecCompInitialized) {
232:     VecCompInitialized = PETSC_TRUE;
233:     VecRegister(VECCOMP,VecCreate_Comp);
234:     VecCompNormInit();
235:   }

237:   /* Allocate a new Vec_Comp */
238:   if (v->data) { PetscFree(v->data); }
239:   PetscNewLog(v,&s);
240:   PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));
241:   v->data        = (void*)s;
242:   v->petscnative = PETSC_FALSE;

244:   /* Allocate the array of Vec, if it is needed to be done */
245:   if (!x_to_me) {
246:     if (nx) { PetscMalloc1(nx,&s->x); }
247:     if (x) { PetscArraycpy(s->x,x,nx); }
248:   } else s->x = x;

250:   s->nx = nx;

252:   if (nx && x) {
253:     /* Allocate the shared structure, if it is not given */
254:     if (!n) {
255:       for (i=0;i<nx;i++) {
256:         VecGetSize(x[i],&k);
257:         N+= k;
258:         VecGetLocalSize(x[i],&k);
259:         lN+= k;
260:       }
261:       PetscNewLog(v,&n);
262:       s->n = n;
263:       n->n = nx;
264:       n->N = N;
265:       n->lN = lN;
266:       n->friends = 1;
267:     } else { /* If not, check in the vector in the shared structure */
268:       s->n = n;
269:       s->n->friends++;
270:     }

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

276:   PetscObjectChangeTypeName((PetscObject)v,VECCOMP);
277:   PetscObjectComposeFunction((PetscObject)v,"VecCompSetSubVecs_C",VecCompSetSubVecs_Comp);
278:   PetscObjectComposeFunction((PetscObject)v,"VecCompGetSubVecs_C",VecCompGetSubVecs_Comp);
279:   return(0);
280: }

282: SLEPC_EXTERN PetscErrorCode VecCreate_Comp(Vec V)
283: {

287:   VecCreate_Comp_Private(V,NULL,0,PETSC_FALSE,NULL);
288:   return(0);
289: }

291: /*@C
292:    VecCreateComp - Creates a new vector containing several subvectors,
293:    each stored separately.

295:    Collective

297:    Input Parameters:
298: +  comm - communicator for the new Vec
299: .  Nx   - array of (initial) global sizes of child vectors
300: .  n    - number of child vectors
301: .  t    - type of the child vectors
302: -  Vparent - (optional) template vector

304:    Output Parameter:
305: .  V - new vector

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

311:    Level: developer

313: .seealso: VecCreateCompWithVecs(), VecCompSetSubVecs()
314: @*/
315: PetscErrorCode VecCreateComp(MPI_Comm comm,PetscInt *Nx,PetscInt n,VecType t,Vec Vparent,Vec *V)
316: {
318:   Vec            *x;
319:   PetscInt       i;

322:   VecCreate(comm,V);
323:   PetscMalloc1(n,&x);
324:   PetscLogObjectMemory((PetscObject)*V,n*sizeof(Vec));
325:   for (i=0;i<n;i++) {
326:     VecCreate(comm,&x[i]);
327:     VecSetSizes(x[i],PETSC_DECIDE,Nx[i]);
328:     VecSetType(x[i],t);
329:   }
330:   VecCreate_Comp_Private(*V,x,n,PETSC_TRUE,Vparent?((Vec_Comp*)Vparent->data)->n:NULL);
331:   return(0);
332: }

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

338:    Collective on x

340:    Input Parameters:
341: +  x - array of Vecs
342: .  n - number of child vectors
343: -  Vparent - (optional) template vector

345:    Output Parameter:
346: .  V - new vector

348:    Level: developer

350: .seealso: VecCreateComp()
351: @*/
352: PetscErrorCode VecCreateCompWithVecs(Vec *x,PetscInt n,Vec Vparent,Vec *V)
353: {
355:   PetscInt       i;

361:   VecCreate(PetscObjectComm((PetscObject)x[0]),V);
362:   for (i=0;i<n;i++) {
363:     PetscObjectReference((PetscObject)x[i]);
364:   }
365:   VecCreate_Comp_Private(*V,x,n,PETSC_FALSE,Vparent?((Vec_Comp*)Vparent->data)->n:NULL);
366:   return(0);
367: }

369: PetscErrorCode VecDuplicate_Comp(Vec win,Vec *V)
370: {
372:   Vec            *x;
373:   PetscInt       i;
374:   Vec_Comp       *s = (Vec_Comp*)win->data;

377:   SlepcValidVecComp(win,1);
378:   VecCreate(PetscObjectComm((PetscObject)win),V);
379:   PetscMalloc1(s->nx,&x);
380:   PetscLogObjectMemory((PetscObject)*V,s->nx*sizeof(Vec));
381:   for (i=0;i<s->nx;i++) {
382:     if (s->x[i]) {
383:       VecDuplicate(s->x[i],&x[i]);
384:     } else x[i] = NULL;
385:   }
386:   VecCreate_Comp_Private(*V,x,s->nx,PETSC_TRUE,s->n);
387:   return(0);
388: }

390: static PetscErrorCode VecCompGetSubVecs_Comp(Vec win,PetscInt *n,const Vec **x)
391: {
392:   Vec_Comp *s = (Vec_Comp*)win->data;

395:   if (x) *x = s->x;
396:   if (n) *n = s->n->n;
397:   return(0);
398: }

400: /*@C
401:    VecCompGetSubVecs - Returns the entire array of vectors defining a
402:    compound vector.

404:    Collective on win

406:    Input Parameter:
407: .  win - compound vector

409:    Output Parameters:
410: +  n - number of child vectors
411: -  x - array of child vectors

413:    Level: developer

415: .seealso: VecCreateComp()
416: @*/
417: PetscErrorCode VecCompGetSubVecs(Vec win,PetscInt *n,const Vec **x)
418: {

423:   PetscUseMethod(win,"VecCompGetSubVecs_C",(Vec,PetscInt*,const Vec**),(win,n,x));
424:   return(0);
425: }

427: static PetscErrorCode VecCompSetSubVecs_Comp(Vec win,PetscInt n,Vec *x)
428: {
429:   Vec_Comp       *s = (Vec_Comp*)win->data;
430:   PetscInt       i,N,nlocal;
431:   Vec_Comp_N     *nn;

435:   if (!s) SETERRQ(PetscObjectComm((PetscObject)win),PETSC_ERR_ORDER,"Must call VecSetSizes first");
436:   if (!s->nx) {
437:     /* vector has been created via VecCreate+VecSetType+VecSetSizes, so allocate data structures */
438:     PetscMalloc1(n,&s->x);
439:     PetscLogObjectMemory((PetscObject)win,n*sizeof(Vec));
440:     VecGetSize(win,&N);
441:     if (N%n) SETERRQ2(PetscObjectComm((PetscObject)win),1,"Global dimension %D is not divisible by %D",N,n);
442:     VecGetLocalSize(win,&nlocal);
443:     if (nlocal%n) SETERRQ2(PetscObjectComm((PetscObject)win),1,"Local dimension %D is not divisible by %D",nlocal,n);
444:     s->nx = n;
445:     for (i=0;i<n;i++) {
446:       VecCreate(PetscObjectComm((PetscObject)win),&s->x[i]);
447:       VecSetSizes(s->x[i],nlocal/n,N/n);
448:       VecSetFromOptions(s->x[i]);
449:     }
450:     if (!s->n) {
451:       PetscNewLog(win,&nn);
452:       s->n = nn;
453:       nn->N = N;
454:       nn->lN = nlocal;
455:       nn->friends = 1;
456:     }
457:   } else if (n > s->nx) SETERRQ1(PetscObjectComm((PetscObject)win),PETSC_ERR_SUP,"Number of child vectors cannot be larger than %D",s->nx);
458:   if (x) {
459:     PetscArraycpy(s->x,x,n);
460:   }
461:   s->n->n = n;
462:   return(0);
463: }

465: /*@C
466:    VecCompSetSubVecs - Resets the number of subvectors defining a compound vector,
467:    or replaces the subvectors.

469:    Collective on win

471:    Input Parameters:
472: +  win - compound vector
473: .  n - number of child vectors
474: -  x - array of child vectors

476:    Note:
477:    It is not possible to increase the number of subvectors with respect to the
478:    number set at its creation.

480:    Level: developer

482: .seealso: VecCreateComp(), VecCompGetSubVecs()
483: @*/
484: PetscErrorCode VecCompSetSubVecs(Vec win,PetscInt n,Vec *x)
485: {

491:   PetscTryMethod(win,"VecCompSetSubVecs_C",(Vec,PetscInt,Vec*),(win,n,x));
492:   return(0);
493: }

495: PetscErrorCode VecAXPY_Comp(Vec v,PetscScalar alpha,Vec w)
496: {
498:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
499:   PetscInt       i;

502:   SlepcValidVecComp(v,1);
503:   SlepcValidVecComp(w,3);
504:   for (i=0;i<vs->n->n;i++) {
505:     VecAXPY(vs->x[i],alpha,ws->x[i]);
506:   }
507:   return(0);
508: }

510: PetscErrorCode VecAYPX_Comp(Vec v,PetscScalar alpha,Vec w)
511: {
513:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
514:   PetscInt       i;

517:   SlepcValidVecComp(v,1);
518:   SlepcValidVecComp(w,3);
519:   for (i=0;i<vs->n->n;i++) {
520:     VecAYPX(vs->x[i],alpha,ws->x[i]);
521:   }
522:   return(0);
523: }

525: PetscErrorCode VecAXPBY_Comp(Vec v,PetscScalar alpha,PetscScalar beta,Vec w)
526: {
528:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
529:   PetscInt       i;

532:   SlepcValidVecComp(v,1);
533:   SlepcValidVecComp(w,4);
534:   for (i=0;i<vs->n->n;i++) {
535:     VecAXPBY(vs->x[i],alpha,beta,ws->x[i]);
536:   }
537:   return(0);
538: }

540: PetscErrorCode VecMAXPY_Comp(Vec v,PetscInt n,const PetscScalar *alpha,Vec *w)
541: {
543:   Vec_Comp       *vs = (Vec_Comp*)v->data;
544:   Vec            *wx;
545:   PetscInt       i,j;

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

551:   PetscMalloc1(n,&wx);

553:   for (j=0;j<vs->n->n;j++) {
554:     for (i=0;i<n;i++) wx[i] = ((Vec_Comp*)w[i]->data)->x[j];
555:     VecMAXPY(vs->x[j],n,alpha,wx);
556:   }

558:   PetscFree(wx);
559:   return(0);
560: }

562: PetscErrorCode VecWAXPY_Comp(Vec v,PetscScalar alpha,Vec w,Vec z)
563: {
565:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data,*zs = (Vec_Comp*)z->data;
566:   PetscInt       i;

569:   SlepcValidVecComp(v,1);
570:   SlepcValidVecComp(w,3);
571:   SlepcValidVecComp(z,4);
572:   for (i=0;i<vs->n->n;i++) {
573:     VecWAXPY(vs->x[i],alpha,ws->x[i],zs->x[i]);
574:   }
575:   return(0);
576: }

578: PetscErrorCode VecAXPBYPCZ_Comp(Vec v,PetscScalar alpha,PetscScalar beta,PetscScalar gamma,Vec w,Vec z)
579: {
580:   PetscErrorCode  ierr;
581:   Vec_Comp        *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data,*zs = (Vec_Comp*)z->data;
582:   PetscInt        i;

585:   SlepcValidVecComp(v,1);
586:   SlepcValidVecComp(w,5);
587:   SlepcValidVecComp(z,6);
588:   for (i=0;i<vs->n->n;i++) {
589:     VecAXPBYPCZ(vs->x[i],alpha,beta,gamma,ws->x[i],zs->x[i]);
590:   }
591:   return(0);
592: }

594: PetscErrorCode VecGetSize_Comp(Vec v,PetscInt *size)
595: {
596:   Vec_Comp *vs = (Vec_Comp*)v->data;

600:   if (vs->n) {
601:     SlepcValidVecComp(v,1);
602:     *size = vs->n->N;
603:   } else *size = v->map->N;
604:   return(0);
605: }

607: PetscErrorCode VecGetLocalSize_Comp(Vec v,PetscInt *size)
608: {
609:   Vec_Comp *vs = (Vec_Comp*)v->data;

613:   if (vs->n) {
614:     SlepcValidVecComp(v,1);
615:     *size = vs->n->lN;
616:   } else *size = v->map->n;
617:   return(0);
618: }

620: PetscErrorCode VecMax_Comp(Vec v,PetscInt *idx,PetscReal *z)
621: {
623:   Vec_Comp       *vs = (Vec_Comp*)v->data;
624:   PetscInt       idxp,s=0,s0;
625:   PetscReal      zp,z0;
626:   PetscInt       i;

629:   SlepcValidVecComp(v,1);
630:   if (!idx && !z) return(0);

632:   if (vs->n->n > 0) {
633:     VecMax(vs->x[0],idx?&idxp:NULL,&zp);
634:   } else {
635:     zp = PETSC_MIN_REAL;
636:     if (idx) idxp = -1;
637:   }
638:   for (i=1;i<vs->n->n;i++) {
639:     VecGetSize(vs->x[i-1],&s0);
640:     s += s0;
641:     VecMax(vs->x[i],idx?&idxp:NULL,&z0);
642:     if (zp < z0) {
643:       if (idx) *idx = s+idxp;
644:       zp = z0;
645:     }
646:   }
647:   if (z) *z = zp;
648:   return(0);
649: }

651: PetscErrorCode VecMin_Comp(Vec v,PetscInt *idx,PetscReal *z)
652: {
654:   Vec_Comp       *vs = (Vec_Comp*)v->data;
655:   PetscInt       idxp,s=0,s0;
656:   PetscReal      zp,z0;
657:   PetscInt       i;

660:   SlepcValidVecComp(v,1);
661:   if (!idx && !z) return(0);

663:   if (vs->n->n > 0) {
664:     VecMin(vs->x[0],idx?&idxp:NULL,&zp);
665:   } else {
666:     zp = PETSC_MAX_REAL;
667:     if (idx) idxp = -1;
668:   }
669:   for (i=1;i<vs->n->n;i++) {
670:     VecGetSize(vs->x[i-1],&s0);
671:     s += s0;
672:     VecMin(vs->x[i],idx?&idxp:NULL,&z0);
673:     if (zp > z0) {
674:       if (idx) *idx = s+idxp;
675:       zp = z0;
676:     }
677:   }
678:   if (z) *z = zp;
679:   return(0);
680: }

682: PetscErrorCode VecMaxPointwiseDivide_Comp(Vec v,Vec w,PetscReal *m)
683: {
685:   Vec_Comp       *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
686:   PetscReal      work;
687:   PetscInt       i;

690:   SlepcValidVecComp(v,1);
691:   SlepcValidVecComp(w,2);
692:   if (!m || vs->n->n == 0) return(0);
693:   VecMaxPointwiseDivide(vs->x[0],ws->x[0],m);
694:   for (i=1;i<vs->n->n;i++) {
695:     VecMaxPointwiseDivide(vs->x[i],ws->x[i],&work);
696:     *m = PetscMax(*m,work);
697:   }
698:   return(0);
699: }


706: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v) \
707: { \
708:   PetscErrorCode  ierr; \
709:   Vec_Comp        *vs = (Vec_Comp*)v->data; \
710:   PetscInt        i; \
711: \
713:   SlepcValidVecComp(v,1); \
714:   for (i=0;i<vs->n->n;i++) { \
715:     __COMPOSE2__(Vec,NAME)(vs->x[i]); \
716:   } \
717:   return(0);\
718: }

720: __FUNC_TEMPLATE1__(Conjugate)
721: __FUNC_TEMPLATE1__(Reciprocal)
722: __FUNC_TEMPLATE1__(SqrtAbs)
723: __FUNC_TEMPLATE1__(Abs)
724: __FUNC_TEMPLATE1__(Exp)
725: __FUNC_TEMPLATE1__(Log)

728: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v,T0 __a) \
729: { \
730:   PetscErrorCode  ierr; \
731:   Vec_Comp        *vs = (Vec_Comp*)v->data; \
732:   PetscInt        i; \
733: \
735:   SlepcValidVecComp(v,1); \
736:   for (i=0;i<vs->n->n;i++) { \
737:     __COMPOSE2__(Vec,NAME)(vs->x[i],__a); \
738:   } \
739:   return(0);\
740: }

742: __FUNC_TEMPLATE2__(Set,PetscScalar)
743: __FUNC_TEMPLATE2__(View,PetscViewer)
744: __FUNC_TEMPLATE2__(Scale,PetscScalar)
745: __FUNC_TEMPLATE2__(SetRandom,PetscRandom)
746: __FUNC_TEMPLATE2__(Shift,PetscScalar)

749: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v,Vec w) \
750: { \
751:   PetscErrorCode  ierr; \
752:   Vec_Comp        *vs = (Vec_Comp*)v->data,\
753:                   *ws = (Vec_Comp*)w->data; \
754:   PetscInt        i; \
755: \
757:   SlepcValidVecComp(v,1); \
758:   SlepcValidVecComp(w,2); \
759:   for (i=0;i<vs->n->n;i++) { \
760:     __COMPOSE2__(Vec,NAME)(vs->x[i],ws->x[i]); \
761:   } \
762:   return(0);\
763: }

765: __FUNC_TEMPLATE3__(Copy)
766: __FUNC_TEMPLATE3__(Swap)

769: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v,Vec w,Vec z) \
770: { \
771:   PetscErrorCode  ierr; \
772:   Vec_Comp        *vs = (Vec_Comp*)v->data, \
773:                   *ws = (Vec_Comp*)w->data, \
774:                   *zs = (Vec_Comp*)z->data; \
775:   PetscInt        i; \
776: \
778:   SlepcValidVecComp(v,1); \
779:   SlepcValidVecComp(w,2); \
780:   SlepcValidVecComp(z,3); \
781:   for (i=0;i<vs->n->n;i++) { \
782:     __COMPOSE2__(Vec,NAME)(vs->x[i],ws->x[i],zs->x[i]); \
783:   } \
784:   return(0);\
785: }

787: __FUNC_TEMPLATE4__(PointwiseMax)
788: __FUNC_TEMPLATE4__(PointwiseMaxAbs)
789: __FUNC_TEMPLATE4__(PointwiseMin)
790: __FUNC_TEMPLATE4__(PointwiseMult)
791: __FUNC_TEMPLATE4__(PointwiseDivide)