Actual source code: veccomp.c

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

  6:    This file is part of SLEPc.
  7:       
  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 "petscmat.h"
 23: #include "private/vecimpl.h"          /*I  "petscvec.h"   I*/
 24: #include "veccomp_private.h"
 25:  #include veccomp.h

 27: typedef struct {
 28:   PetscInt      n,        /* number of active subvectors */
 29:                 N,        /* virtual global size */
 30:                 lN,       /* virtual local size */
 31:                 friends;  /* number of vectors sharing this structure */
 32: } Vec_Comp_N;

 34: typedef struct {
 35:   Vec           *x;       /* the vectors */
 36:   PetscInt      nx;       /* number of available subvectors */
 37:   Vec_Comp_N    *n;       /* structure shared by friend vectors */
 38: } Vec_Comp;

 40: #if defined(PETSC_USE_DEBUG)
 42:   if (((Vec_Comp*)(x)->data)->nx < ((Vec_Comp*)(x)->data)->n->n) { \
 43:     return PetscError(__LINE__,__FUNCT__,__FILE__,__SDIR__,1,0,"Invalid number of subvectors required!");}
 44: #else
 46: #endif

 48: static PetscErrorCode VecCreate_Comp_Private(Vec v, Vec *x, PetscInt nx,
 49:                                              PetscTruth x_to_me, Vec_Comp_N* n);

 51: #include "veccomp0.h"

 54: #include "veccomp0.h"

 58: PetscErrorCode PETSCVEC_DLLEXPORT VecRegister_Comp(const char path[])
 59: {


 64:   VecRegisterDynamic(VECCOMP, path, "VecCreate_Comp", VecCreate_Comp);
 65: 

 67:   return(0);
 68: }


 73: PetscErrorCode VecDestroy_Comp(Vec v)
 74: {
 75:   Vec_Comp       *vs = (Vec_Comp*)v->data;
 76:   PetscInt       i;


 81:   /* if memory was published with AMS then destroy it */
 82:   PetscObjectDepublish(v);

 84: #if defined(PETSC_USE_LOG)
 85:   PetscLogObjectState((PetscObject)v,"Length=%D",v->map->n);
 86: #endif
 87:   for(i=0; i<vs->nx; i++) {
 88:     VecDestroy(vs->x[i]);
 89:   }
 90:   if(--vs->n->friends <= 0) {
 91:     PetscFree(vs->n);
 92:   }
 93:   PetscFree(vs->x);
 94:   PetscFree(vs);

 96:   return(0);
 97: }


100: static struct _VecOps DvOps = {VecDuplicate_Comp, /* 1 */
101:             VecDuplicateVecs_Default,
102:             VecDestroyVecs_Default,
103:             VecDot_Comp_MPI,
104:             VecMDot_Comp_MPI,
105:             VecNorm_Comp_MPI,
106:             0,
107:             0,
108:             VecScale_Comp,
109:             VecCopy_Comp, /* 10 */
110:             VecSet_Comp,
111:             VecSwap_Comp,
112:             VecAXPY_Comp,
113:             VecAXPBY_Comp,
114:             VecMAXPY_Comp,
115:             VecAXPY_Comp/*VecAYPX_Comp*/,
116:             VecWAXPY_Comp,
117:             VecAXPBYPCZ_Comp,
118:             VecPointwiseMult_Comp,
119:             VecPointwiseDivide_Comp,
120:             0, /* 20 */
121:             0,0,
122:             0 /*VecGetArray_Seq*/,
123:             VecGetSize_Comp,
124:             VecGetLocalSize_Comp,
125:             0/*VecRestoreArray_Seq*/,
126:             VecMax_Comp,
127:             VecMin_Comp,
128:             VecSetRandom_Comp,
129:             0, /* 30 */
130:             0,
131:             VecDestroy_Comp,
132:             VecView_Comp,
133:             0/*VecPlaceArray_Seq*/,
134:             0/*VecReplaceArray_Seq*/,
135:             VecDot_Comp_Seq,
136:             0,
137:             VecNorm_Comp_Seq,
138:             VecMDot_Comp_Seq,
139:             0, /* 40 */
140:             0,
141:             0, /* VecLoadIntoVectorNative */
142:             VecReciprocal_Comp,
143:             0, /* VecViewNative */
144:             VecConjugate_Comp,
145:             0,0,
146:             0/*VecResetArray_Seq*/,
147:             0,
148:             VecMaxPointwiseDivide_Comp,
149:             0/*VecLoad_Binary*/, /* 50 */
150:             VecPointwiseMax_Comp,
151:             VecPointwiseMaxAbs_Comp,
152:             VecPointwiseMin_Comp,
153:             0,
154:             VecSqrt_Comp,
155:             VecAbs_Comp,
156:             VecExp_Comp,
157:             VecLog_Comp,
158:             0/*VecShift_Comp*/,
159:             0,
160:             VecDotNorm2_Comp_MPI
161:           };

165: static PetscErrorCode VecCreate_Comp_Private(Vec v, Vec *x, PetscInt nx,
166:                                              PetscTruth x_to_me, Vec_Comp_N *n)
167: {
168:   Vec_Comp        *s;
169:   PetscErrorCode  ierr;
170:   PetscInt        N=0, lN=0, i, k;


174:   /* Allocate a new Vec_Comp */
175:   if (v->data) { PetscFree(v->data);  }
176:   PetscNewLog(v, Vec_Comp, &s);
177:   PetscMemcpy(v->ops, &DvOps, sizeof(DvOps));
178:   v->data  = (void*)s;
179:   v->petscnative     = PETSC_FALSE;

181:   /* Allocate the array of Vec, if it is needed to be done */
182:   if (x_to_me != PETSC_TRUE) {
183:     PetscMalloc(sizeof(Vec)*nx, &s->x);
184:     PetscMemcpy(s->x, x, sizeof(Vec)*nx);
185:   } else
186:     s->x = x;

188:   s->nx = nx;
189:   for(i=0; i<nx; i++) {
190:     VecGetSize(x[i], &k);  N+= k;
191:     VecGetLocalSize(x[i], &k);  lN+= k;
192:   }
193: 
194:   /* Allocate the shared structure, if it is not given */
195:   if (!n) {
196:     PetscNewLog(v, Vec_Comp_N, &n);
197:     s->n = n;
198:     n->n = nx;
199:     n->N = N;
200:     n->lN = lN;
201:     n->friends = 1;
202:   }

204:   /* If not, check in the vector in the shared structure */
205:   else {
206:     s->n = n;
207:     s->n->friends++;
208:     s->n->n = nx;
209:   }

211:   /* Set the virtual sizes as the real sizes of the vector */
212:   VecSetSizes(v, s->n->lN, s->n->N);

214:   PetscObjectChangeTypeName((PetscObject)v,VECCOMP);
215:   PetscPublishAll(v);

217:   return(0);
218: }


224: PetscErrorCode PETSCVEC_DLLEXPORT VecCreate_Comp(Vec V)
225: {
226:   PetscErrorCode  ierr;


230:   VecCreate_Comp_Private(V, PETSC_NULL, 0, PETSC_FALSE, PETSC_NULL);
231: 

233:   return(0);
234: }


240: PetscErrorCode PETSCVEC_DLLEXPORT VecCreateComp(MPI_Comm comm, PetscInt *Nx,
241:                                                 PetscInt n, const VecType t,
242:                                                 Vec Vparent, Vec *V)
243: {
244:   PetscErrorCode  ierr;
245:   Vec             *x;
246:   PetscInt        i;


250:   VecCreate(comm, V);
251:   PetscMalloc(sizeof(Vec)*n, &x);
252:   for(i=0; i<n; i++) {
253:     VecCreate(comm, &x[i]);
254:     VecSetSizes(x[i], PETSC_DECIDE, Nx[i]);
255:     VecSetType(x[i], t);
256:   }
257:   VecCreate_Comp_Private(*V, x, n, PETSC_TRUE,
258:                            Vparent?((Vec_Comp*)Vparent->data)->n:PETSC_NULL);
259: 

261:   return(0);
262: }

266: PetscErrorCode PETSCVEC_DLLEXPORT VecCreateCompWithVecs(Vec *x, PetscInt n,
267:                                                         Vec Vparent, Vec *V)
268: {
269:   PetscErrorCode  ierr;
270:   PetscInt        i;


274:   VecCreate(((PetscObject)x[0])->comm, V);
275:   for(i=0; i<n; i++) {
276:     PetscObjectReference((PetscObject)x[i]);
277:   }
278:   VecCreate_Comp_Private(*V, x, n, PETSC_FALSE,
279:                            Vparent?((Vec_Comp*)Vparent->data)->n:PETSC_NULL);
280: 

282:   return(0);
283: }


288: PetscErrorCode VecDuplicate_Comp(Vec win, Vec *V)
289: {
290:   PetscErrorCode  ierr;
291:   Vec             *x;
292:   PetscInt        i;
293:   Vec_Comp        *s = (Vec_Comp*)win->data;



299:   VecCreate(((PetscObject)win)->comm, V);
300:   PetscMalloc(sizeof(Vec)*s->nx, &x);
301:   for(i=0; i<s->nx; i++) {
302:     VecDuplicate(s->x[i], &x[i]);
303:   }
304:   VecCreate_Comp_Private(*V, x, s->nx, PETSC_TRUE, s->n);

306:   return(0);
307: }


312: PetscErrorCode VecCompGetVecs(Vec win, const Vec **x, PetscInt *n)
313: {
314:   Vec_Comp        *s = (Vec_Comp*)win->data;



320:   if(x) *x = s->x;
321:   if(n) *n = s->nx;

323:   return(0);
324: }


329: PetscErrorCode VecCompSetVecs(Vec win, Vec *x, PetscInt n)
330: {
331:   Vec_Comp        *s = (Vec_Comp*)win->data;
332:   PetscErrorCode  ierr;



338:   if(x) {
339:     if (n > s->nx) {
340:       PetscFree(s->x);
341:       PetscMalloc(sizeof(Vec)*n, &s->x);
342:     }
343:     PetscMemcpy(s->x, x, sizeof(Vec)*n);
344:     s->nx = n;
345:   }
346:   s->n->n = n;

348:   return(0);
349: }


354: PetscErrorCode VecAXPY_Comp(Vec v, PetscScalar alpha, Vec w)
355: {
356:   PetscErrorCode  ierr;
357:   Vec_Comp        *vs = (Vec_Comp*)v->data,
358:                   *ws = (Vec_Comp*)w->data;
359:   PetscInt        i;



366:   for(i=0; i<vs->n->n; i++) {
367:     VecAXPY(vs->x[i], alpha, ws->x[i]);
368:   }

370:   return(0);
371: }


376: PetscErrorCode VecAXPBY_Comp(Vec v, PetscScalar alpha, PetscScalar beta, Vec w)
377: {
378:   PetscErrorCode  ierr;
379:   Vec_Comp        *vs = (Vec_Comp*)v->data,
380:                   *ws = (Vec_Comp*)w->data;
381:   PetscInt        i;



388:   for(i=0; i<vs->n->n; i++) {
389:     VecAXPBY(vs->x[i], alpha, beta, ws->x[i]);
390:   }

392:   return(0);
393: }


398: PetscErrorCode VecMAXPY_Comp(Vec v, PetscInt n, const PetscScalar *alpha,
399:                              Vec *w)
400: {
401:   PetscErrorCode  ierr;
402:   Vec_Comp        *vs = (Vec_Comp*)v->data;
403:   Vec             *wx;
404:   PetscInt        i, j;



411:   PetscMalloc(sizeof(Vec)*n, &wx);

413:   for(j=0; j<vs->n->n; j++) {
414:     for(i=0; i<n; i++) wx[i] = ((Vec_Comp*)w[i]->data)->x[j];
415:     VecMAXPY(vs->x[j], n, alpha, wx);
416:   }

418:   PetscFree(wx);

420:   return(0);
421: }


426: PetscErrorCode VecWAXPY_Comp(Vec v, PetscScalar alpha, Vec w, Vec z)
427: {
428:   PetscErrorCode  ierr;
429:   Vec_Comp        *vs = (Vec_Comp*)v->data,
430:                   *ws = (Vec_Comp*)w->data,
431:                   *zs = (Vec_Comp*)z->data;
432:   PetscInt        i;



440:   for(i=0; i<vs->n->n; i++) {
441:     VecWAXPY(vs->x[i], alpha, ws->x[i], zs->x[i]);
442:   }

444:   return(0);
445: }


450: PetscErrorCode VecAXPBYPCZ_Comp(Vec v, PetscScalar alpha, PetscScalar beta,
451:                                 PetscScalar gamma, Vec w, Vec z)
452: {
453:   PetscErrorCode  ierr;
454:   Vec_Comp        *vs = (Vec_Comp*)v->data,
455:                   *ws = (Vec_Comp*)w->data,
456:                   *zs = (Vec_Comp*)z->data;
457:   PetscInt        i;



465:   for(i=0; i<vs->n->n; i++) {
466:     VecAXPBYPCZ(vs->x[i], alpha, beta, gamma, ws->x[i], zs->x[i]);
467: 
468:   }

470:   return(0);
471: }


476: PetscErrorCode VecGetSize_Comp(Vec v, PetscInt *size)
477: {
478:   Vec_Comp        *vs = (Vec_Comp*)v->data;



484:   if (size) *size = vs->n->N;

486:   return(0);
487: }


492: PetscErrorCode VecGetLocalSize_Comp(Vec v, PetscInt *size)
493: {
494:   Vec_Comp        *vs = (Vec_Comp*)v->data;



500:   if (size) *size = vs->n->lN;

502:   return(0);
503: }



509: PetscErrorCode VecMax_Comp(Vec v, PetscInt *idx, PetscReal *z)
510: {
511:   PetscErrorCode  ierr;
512:   Vec_Comp        *vs = (Vec_Comp*)v->data;
513:   PetscInt        idxp, s=0, s0;
514:   PetscReal       zp, z0;
515:   PetscInt        i;



521:   if(!idx && !z)
522:     return(0);

524:   if (vs->n->n > 0) {
525:     VecMax(vs->x[0], idx?&idxp:PETSC_NULL, &zp);
526:   }
527:   for (i=1; i<vs->n->n; i++) {
528:     VecGetSize(vs->x[i-1], &s0);  s+= s0;
529:     VecMax(vs->x[i], idx?&idxp:PETSC_NULL, &z0);
530:     if(zp < z0) {
531:       if(idx) *idx = s+idxp;
532:       zp = z0;
533:     }
534:   }
535:   if (z) *z = zp;

537:   return(0);
538: }


543: PetscErrorCode VecMin_Comp(Vec v, PetscInt *idx, PetscReal *z)
544: {
545:   PetscErrorCode  ierr;
546:   Vec_Comp        *vs = (Vec_Comp*)v->data;
547:   PetscInt        idxp, s=0, s0;
548:   PetscReal       zp, z0;
549:   PetscInt        i;



555:   if(!idx && !z)
556:     return(0);

558:   if (vs->n->n > 0) {
559:     VecMin(vs->x[0], idx?&idxp:PETSC_NULL, &zp);
560:   }
561:   for (i=1; i<vs->n->n; i++) {
562:     VecGetSize(vs->x[i-1], &s0);  s+= s0;
563:     VecMin(vs->x[i], idx?&idxp:PETSC_NULL, &z0);
564:     if(zp > z0) {
565:       if(idx) *idx = s+idxp;
566:       zp = z0;
567:     }
568:   }
569:   if (z) *z = zp;

571:   return(0);
572: }


577: PetscErrorCode VecMaxPointwiseDivide_Comp(Vec v, Vec w, PetscReal *m)
578: {
579:   PetscErrorCode  ierr;
580:   Vec_Comp        *vs = (Vec_Comp*)v->data,
581:                   *ws = (Vec_Comp*)w->data;
582:   PetscReal       work;
583:   PetscInt        i;



590:   if(!m || vs->n->n == 0)
591:     return(0);
592:   VecMaxPointwiseDivide(vs->x[0], ws->x[0], m);
593:   for (i=1; i<vs->n->n; i++) {
594:     VecMaxPointwiseDivide(vs->x[i], ws->x[i], &work);
595:     *m = PetscMax(*m, work);
596:   }

598:   return(0);
599: }



607: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v) \
608: { \
609:   PetscErrorCode  ierr; \
610:   Vec_Comp        *vs = (Vec_Comp*)v->data; \
611:   PetscInt        i; \
612: \
614: \
616: \
617:   for(i=0; i<vs->n->n; i++) { \
618:     __COMPOSE2__(Vec,NAME)(vs->x[i]);  \
619:   } \
620: \
621:   return(0);\
622: } \

626: __FUNC_TEMPLATE1__(Conjugate)

630: __FUNC_TEMPLATE1__(Reciprocal)

634: __FUNC_TEMPLATE1__(Sqrt)

638: __FUNC_TEMPLATE1__(Abs)

642: __FUNC_TEMPLATE1__(Exp)

646: __FUNC_TEMPLATE1__(Log)


650: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v, T0 __a) \
651: { \
652:   PetscErrorCode  ierr; \
653:   Vec_Comp        *vs = (Vec_Comp*)v->data; \
654:   PetscInt        i; \
655: \
657: \
659: \
660:   for(i=0; i<vs->n->n; i++) { \
661:     __COMPOSE2__(Vec,NAME)(vs->x[i], __a);  \
662:   } \
663: \
664:   return(0);\
665: } \

669: __FUNC_TEMPLATE2__(Set, PetscScalar)

673: __FUNC_TEMPLATE2__(View, PetscViewer)

677: __FUNC_TEMPLATE2__(Scale, PetscScalar)

681: __FUNC_TEMPLATE2__(SetRandom, PetscRandom)

685: __FUNC_TEMPLATE2__(Shift, PetscScalar)


689: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v, Vec w) \
690: { \
691:   PetscErrorCode  ierr; \
692:   Vec_Comp        *vs = (Vec_Comp*)v->data, \
693:                   *ws = (Vec_Comp*)w->data; \
694:   PetscInt        i; \
695: \
697: \
700: \
701:   for(i=0; i<vs->n->n; i++) { \
702:     __COMPOSE2__(Vec,NAME)(vs->x[i], ws->x[i]);  \
703:   } \
704: \
705:   return(0);\
706: } \

710: __FUNC_TEMPLATE3__(Copy)

714: __FUNC_TEMPLATE3__(Swap)


718: PetscErrorCode __COMPOSE3__(Vec,NAME,_Comp)(Vec v, Vec w, Vec z) \
719: { \
720:   PetscErrorCode  ierr; \
721:   Vec_Comp        *vs = (Vec_Comp*)v->data, \
722:                   *ws = (Vec_Comp*)w->data, \
723:                   *zs = (Vec_Comp*)z->data; \
724:   PetscInt        i; \
725: \
727: \
731: \
732:   for(i=0; i<vs->n->n; i++) { \
733:     __COMPOSE2__(Vec,NAME)(vs->x[i], ws->x[i], zs->x[i]); \
734:      \
735:   } \
736: \
737:   return(0);\
738: } \

742: __FUNC_TEMPLATE4__(PointwiseMax)

746: __FUNC_TEMPLATE4__(PointwiseMaxAbs)

750: __FUNC_TEMPLATE4__(PointwiseMin)

754: __FUNC_TEMPLATE4__(PointwiseMult)

758: __FUNC_TEMPLATE4__(PointwiseDivide)