Actual source code: bvec2.c

  1: /*$Id: bvec2.c,v 1.202 2001/09/12 03:26:24 bsmith Exp $*/
  2: /*
  3:    Implements the sequential vectors.
  4: */

 6:  #include src/vec/vecimpl.h
 7:  #include src/vec/impls/dvecimpl.h
 8:  #include petscblaslapack.h
  9: #if defined(PETSC_HAVE_AMS)
 10: EXTERN int PetscViewerAMSGetAMSComm(PetscViewer,AMS_Comm *);
 11: #endif

 13: int VecNorm_Seq(Vec xin,NormType type,PetscReal* z)
 14: {
 15:   PetscScalar *xa;
 16:   int         ierr,one = 1;

 19:   if (type == NORM_2) {
 20:     VecGetArrayFast(xin,&xa);
 21:     /*
 22:       This is because the Fortran BLAS 1 Norm is very slow! 
 23:     */
 24: #if defined(PETSC_HAVE_SLOW_NRM2)
 25:     {
 26:       int         i;
 27:       PetscScalar sum=0.0;
 28:       for (i=0; i<xin->n; i++) {
 29:         sum += (xa[i])*(PetscConj(xa[i]));
 30:       }
 31:       *z = sqrt(PetscRealPart(sum));
 32:     }
 33: #else
 34:     *z = BLnrm2_(&xin->n,xa,&one);
 35: #endif
 36:     VecRestoreArrayFast(xin,&xa);
 37:     PetscLogFlops(2*xin->n-1);
 38:   } else if (type == NORM_INFINITY) {
 39:     int          i,n = xin->n;
 40:     PetscReal    max = 0.0,tmp;

 42:     VecGetArrayFast(xin,&xa);
 43:     for (i=0; i<n; i++) {
 44:       if ((tmp = PetscAbsScalar(*xa)) > max) max = tmp;
 45:       /* check special case of tmp == NaN */
 46:       if (tmp != tmp) {max = tmp; break;}
 47:       xa++;
 48:     }
 49:     VecRestoreArrayFast(xin,&xa);
 50:     *z   = max;
 51:   } else if (type == NORM_1) {
 52:     VecGetArrayFast(xin,&xa);
 53:     *z = BLasum_(&xin->n,xa,&one);
 54:     VecRestoreArrayFast(xin,&xa);
 55:     PetscLogFlops(xin->n-1);
 56:   } else if (type == NORM_1_AND_2) {
 57:     VecNorm_Seq(xin,NORM_1,z);
 58:     VecNorm_Seq(xin,NORM_2,z+1);
 59:   }
 60:   return(0);
 61: }

 63:  #include petscviewer.h
 64:  #include petscsys.h

 66: int VecView_Seq_File(Vec xin,PetscViewer viewer)
 67: {
 68:   Vec_Seq           *x = (Vec_Seq *)xin->data;
 69:   int               i,n = xin->n,ierr;
 70:   char              *name;
 71:   PetscViewerFormat format;

 74:   PetscViewerGetFormat(viewer,&format);
 75:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
 76:     PetscObjectGetName((PetscObject)xin,&name);
 77:     PetscViewerASCIIPrintf(viewer,"%s = [n",name);
 78:     for (i=0; i<n; i++) {
 79: #if defined(PETSC_USE_COMPLEX)
 80:       if (PetscImaginaryPart(x->array[i]) > 0.0) {
 81:         PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ein",PetscRealPart(x->array[i]),PetscImaginaryPart(x->array[i]));
 82:       } else if (PetscImaginaryPart(x->array[i]) < 0.0) {
 83:         PetscViewerASCIIPrintf(viewer,"%18.16e - %18.16ein",PetscRealPart(x->array[i]),-PetscImaginaryPart(x->array[i]));
 84:       } else {
 85:         PetscViewerASCIIPrintf(viewer,"%18.16en",PetscRealPart(x->array[i]));
 86:       }
 87: #else
 88:       PetscViewerASCIIPrintf(viewer,"%18.16en",x->array[i]);
 89: #endif
 90:     }
 91:     PetscViewerASCIIPrintf(viewer,"];n");
 92:   } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
 93:     for (i=0; i<n; i++) {
 94: #if defined(PETSC_USE_COMPLEX)
 95:       PetscViewerASCIIPrintf(viewer,"%18.16e %18.16en",PetscRealPart(x->array[i]),PetscImaginaryPart(x->array[i]));
 96: #else
 97:       PetscViewerASCIIPrintf(viewer,"%18.16en",x->array[i]);
 98: #endif
 99:     }
100:   } else {
101:     for (i=0; i<n; i++) {
102:       if (format == PETSC_VIEWER_ASCII_INDEX) {
103:         PetscViewerASCIIPrintf(viewer,"%d: ",i);
104:       }
105: #if defined(PETSC_USE_COMPLEX)
106:       if (PetscImaginaryPart(x->array[i]) > 0.0) {
107:         PetscViewerASCIIPrintf(viewer,"%g + %g in",PetscRealPart(x->array[i]),PetscImaginaryPart(x->array[i]));
108:       } else if (PetscImaginaryPart(x->array[i]) < 0.0) {
109:         PetscViewerASCIIPrintf(viewer,"%g - %g in",PetscRealPart(x->array[i]),-PetscImaginaryPart(x->array[i]));
110:       } else {
111:         PetscViewerASCIIPrintf(viewer,"%gn",PetscRealPart(x->array[i]));
112:       }
113: #else
114:       PetscViewerASCIIPrintf(viewer,"%gn",x->array[i]);
115: #endif
116:     }
117:   }
118:   PetscViewerFlush(viewer);
119:   return(0);
120: }

122: static int VecView_Seq_Draw_LG(Vec xin,PetscViewer v)
123: {
124:   Vec_Seq     *x = (Vec_Seq *)xin->data;
125:   int         i,n = xin->n,ierr;
126:   PetscDraw   win;
127:   PetscReal   *xx;
128:   PetscDrawLG lg;

131:   PetscViewerDrawGetDrawLG(v,0,&lg);
132:   PetscDrawLGGetDraw(lg,&win);
133:   PetscDrawCheckResizedWindow(win);
134:   PetscDrawLGReset(lg);
135:   PetscMalloc((n+1)*sizeof(PetscReal),&xx);
136:   for (i=0; i<n; i++) {
137:     xx[i] = (PetscReal) i;
138:   }
139: #if !defined(PETSC_USE_COMPLEX)
140:   PetscDrawLGAddPoints(lg,n,&xx,&x->array);
141: #else 
142:   {
143:     PetscReal *yy;
144:     PetscMalloc((n+1)*sizeof(PetscReal),&yy);
145:     for (i=0; i<n; i++) {
146:       yy[i] = PetscRealPart(x->array[i]);
147:     }
148:     PetscDrawLGAddPoints(lg,n,&xx,&yy);
149:     PetscFree(yy);
150:   }
151: #endif
152:   PetscFree(xx);
153:   PetscDrawLGDraw(lg);
154:   PetscDrawSynchronizedFlush(win);
155:   return(0);
156: }

158: static int VecView_Seq_Draw(Vec xin,PetscViewer v)
159: {
160:   int               ierr;
161:   PetscDraw         draw;
162:   PetscTruth        isnull;
163:   PetscViewerFormat format;

166:   PetscViewerDrawGetDraw(v,0,&draw);
167:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
168: 
169:   PetscViewerGetFormat(v,&format);
170:   /*
171:      Currently it only supports drawing to a line graph */
172:   if (format != PETSC_VIEWER_DRAW_LG) {
173:     PetscViewerPushFormat(v,PETSC_VIEWER_DRAW_LG);
174:   }
175:   VecView_Seq_Draw_LG(xin,v);
176:   if (format != PETSC_VIEWER_DRAW_LG) {
177:     PetscViewerPopFormat(v);
178:   }

180:   return(0);
181: }

183: static int VecView_Seq_Binary(Vec xin,PetscViewer viewer)
184: {
185:   Vec_Seq  *x = (Vec_Seq *)xin->data;
186:   int      ierr,fdes,n = xin->n,cookie=VEC_FILE_COOKIE;
187:   FILE     *file;

190:   ierr  = PetscViewerBinaryGetDescriptor(viewer,&fdes);
191:   /* Write vector header */
192:   PetscBinaryWrite(fdes,&cookie,1,PETSC_INT,0);
193:   PetscBinaryWrite(fdes,&n,1,PETSC_INT,0);

195:   /* Write vector contents */
196:   PetscBinaryWrite(fdes,x->array,n,PETSC_SCALAR,0);

198:   PetscViewerBinaryGetInfoPointer(viewer,&file);
199:   if (file && xin->bs > 1) {
200:     fprintf(file,"-vecload_block_size %dn",xin->bs);
201:   }
202:   return(0);
203: }


206: int VecView_Seq(Vec xin,PetscViewer viewer)
207: {
208:   Vec_Seq     *x = (Vec_Seq *)xin->data;
209:   int         ierr;
210:   PetscTruth  isdraw,isascii,issocket,isbinary;

213:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
214:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
215:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
216:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
217:   if (isdraw){
218:     VecView_Seq_Draw(xin,viewer);
219:   } else if (isascii){
220:     VecView_Seq_File(xin,viewer);
221:   } else if (issocket) {
222:     PetscViewerSocketPutScalar(viewer,xin->n,1,x->array);
223:   } else if (isbinary) {
224:     VecView_Seq_Binary(xin,viewer);
225:   } else {
226:     SETERRQ1(1,"Viewer type %s not supported by this vector object",((PetscObject)viewer)->type_name);
227:   }
228:   return(0);
229: }

231: int VecSetValues_Seq(Vec xin,int ni,const int ix[],const PetscScalar y[],InsertMode m)
232: {
233:   Vec_Seq     *x = (Vec_Seq *)xin->data;
234:   PetscScalar *xx = x->array;
235:   int         i;

238:   if (m == INSERT_VALUES) {
239:     for (i=0; i<ni; i++) {
240:       if (ix[i] < 0) continue;
241: #if defined(PETSC_USE_BOPT_g)
242:       if (ix[i] >= xin->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",ix[i],xin->n);
243: #endif
244:       xx[ix[i]] = y[i];
245:     }
246:   } else {
247:     for (i=0; i<ni; i++) {
248:       if (ix[i] < 0) continue;
249: #if defined(PETSC_USE_BOPT_g)
250:       if (ix[i] >= xin->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",ix[i],xin->n);
251: #endif
252:       xx[ix[i]] += y[i];
253:     }
254:   }
255:   return(0);
256: }

258: int VecSetValuesBlocked_Seq(Vec xin,int ni,const int ix[],const PetscScalar yin[],InsertMode m)
259: {
260:   Vec_Seq     *x = (Vec_Seq *)xin->data;
261:   PetscScalar *xx = x->array,*y = (PetscScalar*)yin;
262:   int         i,bs = xin->bs,start,j;

264:   /*
265:        For optimization could treat bs = 2, 3, 4, 5 as special cases with loop unrolling
266:   */
268:   if (m == INSERT_VALUES) {
269:     for (i=0; i<ni; i++) {
270:       start = bs*ix[i];
271:       if (start < 0) continue;
272: #if defined(PETSC_USE_BOPT_g)
273:       if (start >= xin->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",start,xin->n);
274: #endif
275:       for (j=0; j<bs; j++) {
276:         xx[start+j] = y[j];
277:       }
278:       y += bs;
279:     }
280:   } else {
281:     for (i=0; i<ni; i++) {
282:       start = bs*ix[i];
283:       if (start < 0) continue;
284: #if defined(PETSC_USE_BOPT_g)
285:       if (start >= xin->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",start,xin->n);
286: #endif
287:       for (j=0; j<bs; j++) {
288:         xx[start+j] += y[j];
289:       }
290:       y += bs;
291:     }
292:   }
293:   return(0);
294: }


297: int VecDestroy_Seq(Vec v)
298: {
299:   Vec_Seq *vs = (Vec_Seq*)v->data;
300:   int     ierr;


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

307: #if defined(PETSC_USE_LOG)
308:   PetscLogObjectState((PetscObject)v,"Length=%d",v->n);
309: #endif
310:   if (vs->array_allocated) {PetscFree(vs->array_allocated);}
311:   PetscFree(vs);

313:   return(0);
314: }

316: static int VecPublish_Seq(PetscObject obj)
317: {
318: #if defined(PETSC_HAVE_AMS)
319:   Vec          v = (Vec) obj;
320:   Vec_Seq      *s = (Vec_Seq*)v->data;
321:   int          ierr,(*f)(AMS_Memory,char *,Vec);
322: #endif


326: #if defined(PETSC_HAVE_AMS)
327:   /* if it is already published then return */
328:   if (v->amem >=0) return(0);

330:   /* if array in vector was not allocated (for example PCSetUp_BJacobi_Singleblock()) then
331:      cannot AMS publish the object*/
332:   if (!s->array) return(0);

334:   PetscObjectPublishBaseBegin(obj);
335:   AMS_Memory_add_field((AMS_Memory)v->amem,"values",s->array,v->n,AMS_DOUBLE,AMS_READ,
336:                                 AMS_DISTRIBUTED,AMS_REDUCT_UNDEF);

338:   /* if the vector knows its "layout" let it set it*/
339:   PetscObjectQueryFunction(obj,"AMSSetFieldBlock_C",(void (**)(void))&f);
340:   if (f) {
341:     (*f)((AMS_Memory)v->amem,"values",v);
342:   }
343:   PetscObjectPublishBaseEnd(obj);
344: #endif

346:   return(0);
347: }

349: static struct _VecOps DvOps = {VecDuplicate_Seq,
350:             VecDuplicateVecs_Default,
351:             VecDestroyVecs_Default,
352:             VecDot_Seq,
353:             VecMDot_Seq,
354:             VecNorm_Seq,
355:             VecTDot_Seq,
356:             VecMTDot_Seq,
357:             VecScale_Seq,
358:             VecCopy_Seq,
359:             VecSet_Seq,
360:             VecSwap_Seq,
361:             VecAXPY_Seq,
362:             VecAXPBY_Seq,
363:             VecMAXPY_Seq,
364:             VecAYPX_Seq,
365:             VecWAXPY_Seq,
366:             VecPointwiseMult_Seq,
367:             VecPointwiseDivide_Seq,
368:             VecSetValues_Seq,0,0,
369:             VecGetArray_Seq,
370:             VecGetSize_Seq,
371:             VecGetSize_Seq,
372:             VecRestoreArray_Seq,
373:             VecMax_Seq,
374:             VecMin_Seq,
375:             VecSetRandom_Seq,0,
376:             VecSetValuesBlocked_Seq,
377:             VecDestroy_Seq,
378:             VecView_Seq,
379:             VecPlaceArray_Seq,
380:             VecReplaceArray_Seq,
381:             VecDot_Seq,
382:             VecTDot_Seq,
383:             VecNorm_Seq,
384:             VecLoadIntoVector_Default,
385:             VecReciprocal_Default,
386:             0, /* VecViewNative */
387:             VecConjugate_Seq,
388:             0,
389:             0,
390:             VecResetArray_Seq};


393: /*
394:       This is called by VecCreate_Seq() (i.e. VecCreateSeq()) and VecCreateSeqWithArray()
395: */
396: static int VecCreate_Seq_Private(Vec v,const PetscScalar array[])
397: {
398:   Vec_Seq *s;
399:   int     ierr;

402:   PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));
403:   PetscNew(Vec_Seq,&s);
404:   PetscMemzero(s,sizeof(Vec_Seq));
405:   v->data            = (void*)s;
406:   v->bops->publish   = VecPublish_Seq;
407:   v->n               = PetscMax(v->n,v->N);
408:   v->N               = PetscMax(v->n,v->N);
409:   v->bs              = -1;
410:   v->petscnative     = PETSC_TRUE;
411:   s->array           = (PetscScalar *)array;
412:   s->array_allocated = 0;
413:   if (!v->map) {
414:     PetscMapCreateMPI(v->comm,v->n,v->N,&v->map);
415:   }
416:   PetscObjectChangeTypeName((PetscObject)v,VECSEQ);
417: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
418:   PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEnginePut_C","VecMatlabEnginePut_Default",VecMatlabEnginePut_Default);
419:   PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEngineGet_C","VecMatlabEngineGet_Default",VecMatlabEngineGet_Default);
420: #endif
421:   PetscPublishAll(v);
422:   return(0);
423: }

425: /*@C
426:    VecCreateSeqWithArray - Creates a standard,sequential array-style vector,
427:    where the user provides the array space to store the vector values.

429:    Collective on MPI_Comm

431:    Input Parameter:
432: +  comm - the communicator, should be PETSC_COMM_SELF
433: .  n - the vector length 
434: -  array - memory where the vector elements are to be stored.

436:    Output Parameter:
437: .  V - the vector

439:    Notes:
440:    Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
441:    same type as an existing vector.

443:    If the user-provided array is PETSC_NULL, then VecPlaceArray() can be used
444:    at a later stage to SET the array for storing the vector values.

446:    PETSc does NOT free the array when the vector is destroyed via VecDestroy().
447:    The user should not free the array until the vector is destroyed.

449:    Level: intermediate

451:    Concepts: vectors^creating with array

453: .seealso: VecCreateMPIWithArray(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), 
454:           VecCreateGhost(), VecCreateSeq(), VecPlaceArray()
455: @*/
456: int VecCreateSeqWithArray(MPI_Comm comm,int n,const PetscScalar array[],Vec *V)
457: {
458:   int  ierr;

461:   VecCreate(comm,V);
462:   VecSetSizes(*V,n,n);
463:   VecCreate_Seq_Private(*V,array);
464:   return(0);
465: }

467: EXTERN_C_BEGIN
468: int VecCreate_Seq(Vec V)
469: {
470:   Vec_Seq      *s;
471:   PetscScalar  *array;
472:   int          ierr,n = PetscMax(V->n,V->N);

475:   PetscMalloc((n+1)*sizeof(PetscScalar),&array);
476:   PetscMemzero(array,n*sizeof(PetscScalar));
477:   VecCreate_Seq_Private(V,array);
478:   s    = (Vec_Seq*)V->data;
479:   s->array_allocated = array;
480:   VecSetSerializeType(V,VEC_SER_SEQ_BINARY);
481:   return(0);
482: }

484: int VecSerialize_Seq(MPI_Comm comm, Vec *vec, PetscViewer viewer, PetscTruth store)
485: {
486:   Vec          v;
487:   Vec_Seq     *x;
488:   PetscScalar *array;
489:   int          fd;
490:   int          vars;
491:   int          ierr;

494:   PetscViewerBinaryGetDescriptor(viewer, &fd);
495:   if (store) {
496:     v    = *vec;
497:     x    = (Vec_Seq *) v->data;
498:     PetscBinaryWrite(fd, &v->n,     1,    PETSC_INT,    0);
499:     PetscBinaryWrite(fd,  x->array, v->n, PETSC_SCALAR, 0);
500:   } else {
501:     PetscBinaryRead(fd, &vars,  1,    PETSC_INT);
502:     VecCreate(comm, &v);
503:     VecSetSizes(v, vars, vars);
504:     PetscMalloc(vars * sizeof(PetscScalar), &array);
505:     PetscBinaryRead(fd,  array, vars, PETSC_SCALAR);
506:     VecCreate_Seq_Private(v, array);
507:     ((Vec_Seq *) v->data)->array_allocated = array;

509:     VecAssemblyBegin(v);
510:     VecAssemblyEnd(v);
511:     *vec = v;
512:   }

514:   return(0);
515: }
516: EXTERN_C_END


519: int VecDuplicate_Seq(Vec win,Vec *V)
520: {
521:   int     ierr;

524:   VecCreateSeq(win->comm,win->n,V);
525:   if (win->mapping) {
526:     (*V)->mapping = win->mapping;
527:     PetscObjectReference((PetscObject)win->mapping);
528:   }
529:   if (win->bmapping) {
530:     (*V)->bmapping = win->bmapping;
531:     PetscObjectReference((PetscObject)win->bmapping);
532:   }
533:   (*V)->bs = win->bs;
534:   PetscOListDuplicate(win->olist,&(*V)->olist);
535:   PetscFListDuplicate(win->qlist,&(*V)->qlist);
536:   return(0);
537: }