Actual source code: veccomp.c
slepc-3.6.1 2015-09-03
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, 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 <slepc/private/vecimplslepc.h> /*I "slepcvec.h" I*/
24: /* Private MPI datatypes and operators */
25: static MPI_Datatype MPIU_NORM2=0, MPIU_NORM1_AND_2=0;
26: static MPI_Op MPIU_NORM2_SUM=0;
27: static PetscBool VecCompInitialized = PETSC_FALSE;
29: /* Private inline functions */
30: PETSC_STATIC_INLINE void SumNorm2(PetscReal *,PetscReal *,PetscReal *,PetscReal *);
31: PETSC_STATIC_INLINE PetscReal GetNorm2(PetscReal,PetscReal);
32: PETSC_STATIC_INLINE void AddNorm2(PetscReal *,PetscReal *,PetscReal);
34: #include "veccomp0.h"
37: #include "veccomp0.h"
39: PETSC_STATIC_INLINE void SumNorm2(PetscReal *ssq0,PetscReal *scale0,PetscReal *ssq1,PetscReal *scale1)
40: {
41: PetscReal q;
42: if (*scale0 > *scale1) {
43: q = *scale1/(*scale0);
44: *ssq1 = *ssq0 + q*q*(*ssq1);
45: *scale1 = *scale0;
46: } else {
47: q = *scale0/(*scale1);
48: *ssq1 += q*q*(*ssq0);
49: }
50: }
52: PETSC_STATIC_INLINE PetscReal GetNorm2(PetscReal ssq,PetscReal scale)
53: {
54: return scale*PetscSqrtReal(ssq);
55: }
57: PETSC_STATIC_INLINE void AddNorm2(PetscReal *ssq,PetscReal *scale,PetscReal x)
58: {
59: PetscReal absx,q;
60: if (x != 0.0) {
61: absx = PetscAbs(x);
62: if (*scale < absx) {
63: q = *scale/absx;
64: *ssq = 1.0 + *ssq*q*q;
65: *scale = absx;
66: } else {
67: q = absx/(*scale);
68: *ssq += q*q;
69: }
70: }
71: }
75: static void SlepcSumNorm2_Local(void *in,void *out,PetscMPIInt *cnt,MPI_Datatype *datatype)
76: {
77: PetscInt i,count = *cnt;
80: if (*datatype == MPIU_NORM2) {
81: PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
82: for (i=0; i<count; i++) {
83: SumNorm2(&xin[i*2],&xin[i*2+1],&xout[i*2],&xout[i*2+1]);
84: }
85: } else if (*datatype == MPIU_NORM1_AND_2) {
86: PetscReal *xin = (PetscReal*)in,*xout = (PetscReal*)out;
87: for (i=0; i<count; i++) {
88: xout[i*3]+= xin[i*3];
89: SumNorm2(&xin[i*3+1],&xin[i*3+2],&xout[i*3+1],&xout[i*3+2]);
90: }
91: } else {
92: (*PetscErrorPrintf)("Can only handle MPIU_NORM* data types");
93: MPI_Abort(MPI_COMM_WORLD,1);
94: }
95: PetscFunctionReturnVoid();
96: }
100: static PetscErrorCode VecNormCompEnd(void)
101: {
105: MPI_Type_free(&MPIU_NORM2);
106: MPI_Type_free(&MPIU_NORM1_AND_2);
107: MPI_Op_free(&MPIU_NORM2_SUM);
108: VecCompInitialized = PETSC_FALSE;
109: return(0);
110: }
114: static PetscErrorCode VecNormCompInit()
115: {
119: MPI_Type_contiguous(sizeof(PetscReal)*2,MPI_BYTE,&MPIU_NORM2);
120: MPI_Type_commit(&MPIU_NORM2);
121: MPI_Type_contiguous(sizeof(PetscReal)*3,MPI_BYTE,&MPIU_NORM1_AND_2);
122: MPI_Type_commit(&MPIU_NORM1_AND_2);
123: MPI_Op_create(SlepcSumNorm2_Local,1,&MPIU_NORM2_SUM);
124: PetscRegisterFinalize(VecNormCompEnd);
125: return(0);
126: }
130: PetscErrorCode VecDestroy_Comp(Vec v)
131: {
132: Vec_Comp *vs = (Vec_Comp*)v->data;
133: PetscInt i;
137: #if defined(PETSC_USE_LOG)
138: PetscLogObjectState((PetscObject)v,"Length=%D",v->map->n);
139: #endif
140: for (i=0;i<vs->nx;i++) {
141: VecDestroy(&vs->x[i]);
142: }
143: if (--vs->n->friends <= 0) {
144: PetscFree(vs->n);
145: }
146: PetscFree(vs->x);
147: PetscFree(vs);
148: return(0);
149: }
151: static struct _VecOps DvOps = {VecDuplicate_Comp, /* 1 */
152: VecDuplicateVecs_Comp,
153: VecDestroyVecs_Comp,
154: VecDot_Comp_MPI,
155: VecMDot_Comp_MPI,
156: VecNorm_Comp_MPI,
157: VecTDot_Comp_MPI,
158: VecMTDot_Comp_MPI,
159: VecScale_Comp,
160: VecCopy_Comp, /* 10 */
161: VecSet_Comp,
162: VecSwap_Comp,
163: VecAXPY_Comp,
164: VecAXPBY_Comp,
165: VecMAXPY_Comp,
166: VecAYPX_Comp,
167: VecWAXPY_Comp,
168: VecAXPBYPCZ_Comp,
169: VecPointwiseMult_Comp,
170: VecPointwiseDivide_Comp,
171: 0, /* 20 */
172: 0,0,
173: 0 /*VecGetArray_Seq*/,
174: VecGetSize_Comp,
175: VecGetLocalSize_Comp,
176: 0/*VecRestoreArray_Seq*/,
177: VecMax_Comp,
178: VecMin_Comp,
179: VecSetRandom_Comp,
180: 0, /* 30 */
181: 0,
182: VecDestroy_Comp,
183: VecView_Comp,
184: 0/*VecPlaceArray_Seq*/,
185: 0/*VecReplaceArray_Seq*/,
186: VecDot_Comp_Seq,
187: 0,
188: VecNorm_Comp_Seq,
189: VecMDot_Comp_Seq,
190: 0, /* 40 */
191: 0,
192: VecReciprocal_Comp,
193: VecConjugate_Comp,
194: 0,0,
195: 0/*VecResetArray_Seq*/,
196: 0,
197: VecMaxPointwiseDivide_Comp,
198: VecPointwiseMax_Comp,
199: VecPointwiseMaxAbs_Comp,
200: VecPointwiseMin_Comp,
201: 0,
202: VecSqrtAbs_Comp,
203: VecAbs_Comp,
204: VecExp_Comp,
205: VecLog_Comp,
206: 0/*VecShift_Comp*/,
207: 0,
208: 0,
209: 0,
210: VecDotNorm2_Comp_MPI
211: };
215: PetscErrorCode VecDuplicateVecs_Comp(Vec w,PetscInt m,Vec *V[])
216: {
218: PetscInt i;
223: if (m<=0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %D",m);
224: PetscMalloc1(m,V);
225: for (i=0;i<m;i++) { VecDuplicate(w,*V+i); }
226: return(0);
227: }
231: PetscErrorCode VecDestroyVecs_Comp(PetscInt m,Vec v[])
232: {
234: PetscInt i;
238: if (m<=0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %D",m);
239: for (i=0;i<m;i++) if (v[i]) { VecDestroy(&v[i]); }
240: PetscFree(v);
241: return(0);
242: }
246: static PetscErrorCode VecCreate_Comp_Private(Vec v,Vec *x,PetscInt nx,PetscBool x_to_me,Vec_Comp_N *n)
247: {
248: Vec_Comp *s;
250: PetscInt N=0,lN=0,i,k;
253: if (!VecCompInitialized) {
254: VecCompInitialized = PETSC_TRUE;
255: VecRegister(VECCOMP,VecCreate_Comp);
256: VecNormCompInit();
257: }
259: /* Allocate a new Vec_Comp */
260: if (v->data) { PetscFree(v->data); }
261: PetscNewLog(v,&s);
262: PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));
263: v->data = (void*)s;
264: v->petscnative = PETSC_FALSE;
266: /* Allocate the array of Vec, if it is needed to be done */
267: if (x_to_me != PETSC_TRUE) {
268: PetscMalloc(sizeof(Vec)*nx,&s->x);
269: PetscMemcpy(s->x,x,sizeof(Vec)*nx);
270: } else s->x = x;
272: s->nx = nx;
274: /* Allocate the shared structure, if it is not given */
275: if (!n) {
276: for (i=0;i<nx;i++) {
277: VecGetSize(x[i],&k);
278: N+= k;
279: VecGetLocalSize(x[i],&k);
280: lN+= k;
281: }
282: PetscNewLog(v,&n);
283: s->n = n;
284: n->n = nx;
285: n->N = N;
286: n->lN = lN;
287: n->friends = 1;
288: } else { /* If not, check in the vector in the shared structure */
289: N = n->N;
290: lN = n->lN;
291: s->n = n;
292: s->n->friends++;
293: }
295: /* Set the virtual sizes as the real sizes of the vector */
296: VecSetSizes(v,s->n->lN,s->n->N);
298: PetscObjectChangeTypeName((PetscObject)v,VECCOMP);
299: return(0);
300: }
304: PETSC_EXTERN PetscErrorCode VecCreate_Comp(Vec V)
305: {
309: VecCreate_Comp_Private(V,NULL,0,PETSC_FALSE,NULL);
310: return(0);
311: }
315: /*@C
316: VecCreateComp - Creates a new vector containing several subvectors,
317: each stored separately.
319: Collective on Vec
321: Input Parameters:
322: + comm - communicator for the new Vec
323: . Nx - array of (initial) global sizes of child vectors
324: . n - number of child vectors
325: . t - type of the child vectors
326: - Vparent - (optional) template vector
328: Output Parameter:
329: . V - new vector
331: Notes:
332: This is similar to PETSc's VecNest but customized for SLEPc's needs. In particular,
333: the number of child vectors can be modified dynamically, with VecCompSetSubVecs().
335: Level: developer
337: .seealso: VecCreateCompWithVecs(), VecCompSetSubVecs()
338: @*/
339: PetscErrorCode VecCreateComp(MPI_Comm comm,PetscInt *Nx,PetscInt n,VecType t,Vec Vparent,Vec *V)
340: {
342: Vec *x;
343: PetscInt i;
346: VecCreate(comm,V);
347: PetscMalloc1(n,&x);
348: PetscLogObjectMemory((PetscObject)*V,n*sizeof(Vec));
349: for (i=0;i<n;i++) {
350: VecCreate(comm,&x[i]);
351: VecSetSizes(x[i],PETSC_DECIDE,Nx[i]);
352: VecSetType(x[i],t);
353: }
354: VecCreate_Comp_Private(*V,x,n,PETSC_TRUE,Vparent?((Vec_Comp*)Vparent->data)->n:NULL);
355: return(0);
356: }
360: /*@C
361: VecCreateCompWithVecs - Creates a new vector containing several subvectors,
362: each stored separately, from an array of Vecs.
364: Collective on Vec
366: Input Parameters:
367: + x - array of Vecs
368: . n - number of child vectors
369: - Vparent - (optional) template vector
371: Output Parameter:
372: . V - new vector
374: Level: developer
376: .seealso: VecCreateComp()
377: @*/
378: PetscErrorCode VecCreateCompWithVecs(Vec *x,PetscInt n,Vec Vparent,Vec *V)
379: {
381: PetscInt i;
384: VecCreate(PetscObjectComm((PetscObject)x[0]),V);
385: for (i=0;i<n;i++) {
386: PetscObjectReference((PetscObject)x[i]);
387: }
388: VecCreate_Comp_Private(*V,x,n,PETSC_FALSE,Vparent?((Vec_Comp*)Vparent->data)->n:NULL);
389: return(0);
390: }
394: PetscErrorCode VecDuplicate_Comp(Vec win,Vec *V)
395: {
397: Vec *x;
398: PetscInt i;
399: Vec_Comp *s = (Vec_Comp*)win->data;
402: SlepcValidVecComp(win);
403: VecCreate(PetscObjectComm((PetscObject)win),V);
404: PetscMalloc1(s->nx,&x);
405: PetscLogObjectMemory((PetscObject)*V,s->nx*sizeof(Vec));
406: for (i=0;i<s->nx;i++) {
407: if (s->x[i]) {
408: VecDuplicate(s->x[i],&x[i]);
409: } else x[i] = NULL;
410: }
411: VecCreate_Comp_Private(*V,x,s->nx,PETSC_TRUE,s->n);
412: return(0);
413: }
417: /*@C
418: VecCompGetSubVecs - Returns the entire array of vectors defining a
419: compound vector.
421: Collective on Vec
423: Input Parameter:
424: . win - compound vector
426: Output Parameters:
427: + n - number of child vectors
428: - x - array of child vectors
430: Level: developer
432: .seealso: VecCreateComp()
433: @*/
434: PetscErrorCode VecCompGetSubVecs(Vec win,PetscInt *n,const Vec **x)
435: {
436: Vec_Comp *s = (Vec_Comp*)win->data;
439: SlepcValidVecComp(win);
440: if (x) *x = s->x;
441: if (n) *n = s->nx;
442: return(0);
443: }
447: /*@C
448: VecCompSetSubVecs - Resets the number of subvectors defining a compound vector,
449: of replaces the subvectors.
451: Collective on Vec
453: Input Parameters:
454: + win - compound vector
455: . n - number of child vectors
456: - x - array of child vectors
458: Level: developer
460: .seealso: VecCreateComp()
461: @*/
462: PetscErrorCode VecCompSetSubVecs(Vec win,PetscInt n,Vec *x)
463: {
464: Vec_Comp *s = (Vec_Comp*)win->data;
468: if (x) {
469: if (n > s->nx) {
470: PetscFree(s->x);
471: PetscMalloc(sizeof(Vec)*n,&s->x);
472: }
473: PetscMemcpy(s->x,x,sizeof(Vec)*n);
474: s->nx = n;
475: }
476: s->n->n = n;
477: return(0);
478: }
482: PetscErrorCode VecAXPY_Comp(Vec v,PetscScalar alpha,Vec w)
483: {
485: Vec_Comp *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
486: PetscInt i;
489: SlepcValidVecComp(v);
490: SlepcValidVecComp(w);
491: for (i=0;i<vs->n->n;i++) {
492: VecAXPY(vs->x[i],alpha,ws->x[i]);
493: }
494: return(0);
495: }
499: PetscErrorCode VecAYPX_Comp(Vec v,PetscScalar alpha,Vec w)
500: {
502: Vec_Comp *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
503: PetscInt i;
506: SlepcValidVecComp(v);
507: SlepcValidVecComp(w);
508: for (i=0;i<vs->n->n;i++) {
509: VecAYPX(vs->x[i],alpha,ws->x[i]);
510: }
511: return(0);
512: }
516: PetscErrorCode VecAXPBY_Comp(Vec v,PetscScalar alpha,PetscScalar beta,Vec w)
517: {
519: Vec_Comp *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
520: PetscInt i;
523: SlepcValidVecComp(v);
524: SlepcValidVecComp(w);
525: for (i=0;i<vs->n->n;i++) {
526: VecAXPBY(vs->x[i],alpha,beta,ws->x[i]);
527: }
528: return(0);
529: }
533: PetscErrorCode VecMAXPY_Comp(Vec v,PetscInt n,const PetscScalar *alpha,Vec *w)
534: {
536: Vec_Comp *vs = (Vec_Comp*)v->data;
537: Vec *wx;
538: PetscInt i,j;
541: SlepcValidVecComp(v);
542: for (i=0;i<n;i++) SlepcValidVecComp(w[i]);
544: PetscMalloc(sizeof(Vec)*n,&wx);
546: for (j=0;j<vs->n->n;j++) {
547: for (i=0;i<n;i++) wx[i] = ((Vec_Comp*)w[i]->data)->x[j];
548: VecMAXPY(vs->x[j],n,alpha,wx);
549: }
551: PetscFree(wx);
552: return(0);
553: }
557: PetscErrorCode VecWAXPY_Comp(Vec v,PetscScalar alpha,Vec w,Vec z)
558: {
560: Vec_Comp *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data,*zs = (Vec_Comp*)z->data;
561: PetscInt i;
564: SlepcValidVecComp(v);
565: SlepcValidVecComp(w);
566: SlepcValidVecComp(z);
567: for (i=0;i<vs->n->n;i++) {
568: VecWAXPY(vs->x[i],alpha,ws->x[i],zs->x[i]);
569: }
570: return(0);
571: }
575: PetscErrorCode VecAXPBYPCZ_Comp(Vec v,PetscScalar alpha,PetscScalar beta,PetscScalar gamma,Vec w,Vec z)
576: {
577: PetscErrorCode ierr;
578: Vec_Comp *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data,*zs = (Vec_Comp*)z->data;
579: PetscInt i;
582: SlepcValidVecComp(v);
583: SlepcValidVecComp(w);
584: SlepcValidVecComp(z);
585: for (i=0;i<vs->n->n;i++) {
586: VecAXPBYPCZ(vs->x[i],alpha,beta,gamma,ws->x[i],zs->x[i]);
587: }
588: return(0);
589: }
593: PetscErrorCode VecGetSize_Comp(Vec v,PetscInt *size)
594: {
595: Vec_Comp *vs = (Vec_Comp*)v->data;
598: SlepcValidVecComp(v);
600: *size = vs->n->N;
601: return(0);
602: }
606: PetscErrorCode VecGetLocalSize_Comp(Vec v,PetscInt *size)
607: {
608: Vec_Comp *vs = (Vec_Comp*)v->data;
611: SlepcValidVecComp(v);
613: *size = vs->n->lN;
614: return(0);
615: }
619: PetscErrorCode VecMax_Comp(Vec v,PetscInt *idx,PetscReal *z)
620: {
622: Vec_Comp *vs = (Vec_Comp*)v->data;
623: PetscInt idxp,s=0,s0;
624: PetscReal zp,z0;
625: PetscInt i;
628: SlepcValidVecComp(v);
629: if (!idx && !z) return(0);
631: if (vs->n->n > 0) {
632: VecMax(vs->x[0],idx?&idxp:NULL,&zp);
633: }
634: for (i=1;i<vs->n->n;i++) {
635: VecGetSize(vs->x[i-1],&s0);
636: s+= s0;
637: VecMax(vs->x[i],idx?&idxp:NULL,&z0);
638: if (zp < z0) {
639: if (idx) *idx = s+idxp;
640: zp = z0;
641: }
642: }
643: if (z) *z = zp;
644: return(0);
645: }
649: PetscErrorCode VecMin_Comp(Vec v,PetscInt *idx,PetscReal *z)
650: {
652: Vec_Comp *vs = (Vec_Comp*)v->data;
653: PetscInt idxp,s=0,s0;
654: PetscReal zp,z0;
655: PetscInt i;
658: SlepcValidVecComp(v);
659: if (!idx && !z) return(0);
661: if (vs->n->n > 0) {
662: VecMin(vs->x[0],idx?&idxp:NULL,&zp);
663: }
664: for (i=1;i<vs->n->n;i++) {
665: VecGetSize(vs->x[i-1],&s0);
666: s+= s0;
667: VecMin(vs->x[i],idx?&idxp:NULL,&z0);
668: if (zp > z0) {
669: if (idx) *idx = s+idxp;
670: zp = z0;
671: }
672: }
673: if (z) *z = zp;
674: return(0);
675: }
679: PetscErrorCode VecMaxPointwiseDivide_Comp(Vec v,Vec w,PetscReal *m)
680: {
682: Vec_Comp *vs = (Vec_Comp*)v->data,*ws = (Vec_Comp*)w->data;
683: PetscReal work;
684: PetscInt i;
687: SlepcValidVecComp(v);
688: SlepcValidVecComp(w);
689: if (!m || vs->n->n == 0) return(0);
690: VecMaxPointwiseDivide(vs->x[0],ws->x[0],m);
691: for (i=1;i<vs->n->n;i++) {
692: VecMaxPointwiseDivide(vs->x[i],ws->x[i],&work);
693: *m = PetscMax(*m,work);
694: }
695: return(0);
696: }
704: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v) \
705: { \
706: PetscErrorCode ierr; \
707: Vec_Comp *vs = (Vec_Comp*)v->data; \
708: PetscInt i; \
709: \
711: SlepcValidVecComp(v); \
712: for (i=0;i<vs->n->n;i++) { \
713: __COMPOSE2__(Vec,NAME)(vs->x[i]); \
714: } \
715: return(0);\
716: }
720: __FUNC_TEMPLATE1__(Conjugate)
724: __FUNC_TEMPLATE1__(Reciprocal)
728: __FUNC_TEMPLATE1__(SqrtAbs)
732: __FUNC_TEMPLATE1__(Abs)
736: __FUNC_TEMPLATE1__(Exp)
740: __FUNC_TEMPLATE1__(Log)
744: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v,T0 __a) \
745: { \
746: PetscErrorCode ierr; \
747: Vec_Comp *vs = (Vec_Comp*)v->data; \
748: PetscInt i; \
749: \
751: SlepcValidVecComp(v); \
752: for (i=0;i<vs->n->n;i++) { \
753: __COMPOSE2__(Vec,NAME)(vs->x[i],__a); \
754: } \
755: return(0);\
756: }
760: __FUNC_TEMPLATE2__(Set,PetscScalar)
764: __FUNC_TEMPLATE2__(View,PetscViewer)
768: __FUNC_TEMPLATE2__(Scale,PetscScalar)
772: __FUNC_TEMPLATE2__(SetRandom,PetscRandom)
776: __FUNC_TEMPLATE2__(Shift,PetscScalar)
780: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v,Vec w) \
781: { \
782: PetscErrorCode ierr; \
783: Vec_Comp *vs = (Vec_Comp*)v->data,\
784: *ws = (Vec_Comp*)w->data; \
785: PetscInt i; \
786: \
788: SlepcValidVecComp(v); \
789: SlepcValidVecComp(w); \
790: for (i=0;i<vs->n->n;i++) { \
791: __COMPOSE2__(Vec,NAME)(vs->x[i],ws->x[i]); \
792: } \
793: return(0);\
794: }
798: __FUNC_TEMPLATE3__(Copy)
802: __FUNC_TEMPLATE3__(Swap)
806: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v,Vec w,Vec z) \
807: { \
808: PetscErrorCode ierr; \
809: Vec_Comp *vs = (Vec_Comp*)v->data, \
810: *ws = (Vec_Comp*)w->data, \
811: *zs = (Vec_Comp*)z->data; \
812: PetscInt i; \
813: \
815: SlepcValidVecComp(v); \
816: SlepcValidVecComp(w); \
817: SlepcValidVecComp(z); \
818: for (i=0;i<vs->n->n;i++) { \
819: __COMPOSE2__(Vec,NAME)(vs->x[i],ws->x[i],zs->x[i]); \
820: } \
821: return(0);\
822: }
826: __FUNC_TEMPLATE4__(PointwiseMax)
830: __FUNC_TEMPLATE4__(PointwiseMaxAbs)
834: __FUNC_TEMPLATE4__(PointwiseMin)
838: __FUNC_TEMPLATE4__(PointwiseMult)
842: __FUNC_TEMPLATE4__(PointwiseDivide)