Actual source code: bvec2.c

  1: #define PETSCVEC_DLL
  2: /*
  3:    Implements the sequential vectors.
  4: */

 6:  #include private/vecimpl.h
 7:  #include ../src/vec/vec/impls/dvecimpl.h
 8:  #include petscblaslapack.h
  9: #if defined(PETSC_HAVE_PNETCDF)
 11: #include "pnetcdf.h"
 13: #endif

 15: #include "../src/vec/vec/impls/seq/ftn-kernels/fnorm.h"
 18: PetscErrorCode VecNorm_Seq(Vec xin,NormType type,PetscReal* z)
 19: {
 20:   PetscScalar    *xx;
 22:   PetscInt       n = xin->map->n;
 23:   PetscBLASInt   one = 1, bn = PetscBLASIntCast(n);

 26:   if (type == NORM_2 || type == NORM_FROBENIUS) {
 27:     VecGetArray(xin,&xx);
 28:     /*
 29:       This is because the Fortran BLAS 1 Norm is very slow! 
 30:     */
 31: #if defined(PETSC_HAVE_SLOW_BLAS_NORM2)
 32: #if defined(PETSC_USE_FORTRAN_KERNEL_NORM)
 33:     fortrannormsqr_(xx,&n,z);
 34:     *z = sqrt(*z);
 35: #elif defined(PETSC_USE_UNROLLED_NORM)
 36:     {
 37:     PetscReal work = 0.0;
 38:     switch (n & 0x3) {
 39:       case 3: work += PetscRealPart(xx[0]*PetscConj(xx[0])); xx++;
 40:       case 2: work += PetscRealPart(xx[0]*PetscConj(xx[0])); xx++;
 41:       case 1: work += PetscRealPart(xx[0]*PetscConj(xx[0])); xx++; n -= 4;
 42:     }
 43:     while (n>0) {
 44:       work += PetscRealPart(xx[0]*PetscConj(xx[0])+xx[1]*PetscConj(xx[1])+
 45:                         xx[2]*PetscConj(xx[2])+xx[3]*PetscConj(xx[3]));
 46:       xx += 4; n -= 4;
 47:     }
 48:     *z = sqrt(work);}
 49: #else
 50:     {
 51:       PetscInt         i;
 52:       PetscScalar sum=0.0;
 53:       for (i=0; i<n; i++) {
 54:         sum += (xx[i])*(PetscConj(xx[i]));
 55:       }
 56:       *z = sqrt(PetscRealPart(sum));
 57:     }
 58: #endif
 59: #else
 60:     *z = BLASnrm2_(&bn,xx,&one);
 61: #endif
 62:     VecRestoreArray(xin,&xx);
 63:     PetscLogFlops(PetscMax(2.0*n-1,0.0));
 64:   } else if (type == NORM_INFINITY) {
 65:     PetscInt          i;
 66:     PetscReal    max = 0.0,tmp;

 68:     VecGetArray(xin,&xx);
 69:     for (i=0; i<n; i++) {
 70:       if ((tmp = PetscAbsScalar(*xx)) > max) max = tmp;
 71:       /* check special case of tmp == NaN */
 72:       if (tmp != tmp) {max = tmp; break;}
 73:       xx++;
 74:     }
 75:     VecRestoreArray(xin,&xx);
 76:     *z   = max;
 77:   } else if (type == NORM_1) {
 78:     VecGetArray(xin,&xx);
 79:     *z = BLASasum_(&bn,xx,&one);
 80:     VecRestoreArray(xin,&xx);
 81:     PetscLogFlops(PetscMax(n-1.0,0.0));
 82:   } else if (type == NORM_1_AND_2) {
 83:     VecNorm_Seq(xin,NORM_1,z);
 84:     VecNorm_Seq(xin,NORM_2,z+1);
 85:   }
 86:   return(0);
 87: }

 91: PetscErrorCode VecView_Seq_File(Vec xin,PetscViewer viewer)
 92: {
 93:   Vec_Seq           *x = (Vec_Seq *)xin->data;
 94:   PetscErrorCode    ierr;
 95:   PetscInt          i,n = xin->map->n;
 96:   const char        *name;
 97:   PetscViewerFormat format;

100:   PetscViewerGetFormat(viewer,&format);
101:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
102:     PetscObjectGetName((PetscObject)xin,&name);
103:     PetscViewerASCIIPrintf(viewer,"%s = [\n",name);
104:     for (i=0; i<n; i++) {
105: #if defined(PETSC_USE_COMPLEX)
106:       if (PetscImaginaryPart(x->array[i]) > 0.0) {
107:         PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei\n",PetscRealPart(x->array[i]),PetscImaginaryPart(x->array[i]));
108:       } else if (PetscImaginaryPart(x->array[i]) < 0.0) {
109:         PetscViewerASCIIPrintf(viewer,"%18.16e - %18.16ei\n",PetscRealPart(x->array[i]),-PetscImaginaryPart(x->array[i]));
110:       } else {
111:         PetscViewerASCIIPrintf(viewer,"%18.16e\n",PetscRealPart(x->array[i]));
112:       }
113: #else
114:       PetscViewerASCIIPrintf(viewer,"%18.16e\n",x->array[i]);
115: #endif
116:     }
117:     PetscViewerASCIIPrintf(viewer,"];\n");
118:   } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
119:     for (i=0; i<n; i++) {
120: #if defined(PETSC_USE_COMPLEX)
121:       PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e\n",PetscRealPart(x->array[i]),PetscImaginaryPart(x->array[i]));
122: #else
123:       PetscViewerASCIIPrintf(viewer,"%18.16e\n",x->array[i]);
124: #endif
125:     }
126:   } else if (format == PETSC_VIEWER_ASCII_VTK || format == PETSC_VIEWER_ASCII_VTK_CELL) {
127:     /* 
128:        state 0: No header has been output
129:        state 1: Only POINT_DATA has been output
130:        state 2: Only CELL_DATA has been output
131:        state 3: Output both, POINT_DATA last
132:        state 4: Output both, CELL_DATA last 
133:     */
134:     static PetscInt stateId = -1;
135:     int outputState = 0;
136:     PetscTruth hasState;
137:     int doOutput = 0;
138:     PetscInt bs, b;

140:     if (stateId < 0) {
141:       PetscObjectComposedDataRegister(&stateId);
142:     }
143:     PetscObjectComposedDataGetInt((PetscObject) viewer, stateId, outputState, hasState);
144:     if (!hasState) {
145:       outputState = 0;
146:     }
147:     PetscObjectGetName((PetscObject) xin, &name);
148:     VecGetBlockSize(xin, &bs);
149:     if ((bs < 1) || (bs > 3)) {
150:       SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %d", bs);
151:     }
152:     if (format == PETSC_VIEWER_ASCII_VTK) {
153:       if (outputState == 0) {
154:         outputState = 1;
155:         doOutput = 1;
156:       } else if (outputState == 1) {
157:         doOutput = 0;
158:       } else if (outputState == 2) {
159:         outputState = 3;
160:         doOutput = 1;
161:       } else if (outputState == 3) {
162:         doOutput = 0;
163:       } else if (outputState == 4) {
164:         SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Tried to output POINT_DATA again after intervening CELL_DATA");
165:       }
166:       if (doOutput) {
167:         PetscViewerASCIIPrintf(viewer, "POINT_DATA %d\n", n);
168:       }
169:     } else {
170:       if (outputState == 0) {
171:         outputState = 2;
172:         doOutput = 1;
173:       } else if (outputState == 1) {
174:         outputState = 4;
175:         doOutput = 1;
176:       } else if (outputState == 2) {
177:         doOutput = 0;
178:       } else if (outputState == 3) {
179:         SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Tried to output CELL_DATA again after intervening POINT_DATA");
180:       } else if (outputState == 4) {
181:         doOutput = 0;
182:       }
183:       if (doOutput) {
184:         PetscViewerASCIIPrintf(viewer, "CELL_DATA %d\n", n);
185:       }
186:     }
187:     PetscObjectComposedDataSetInt((PetscObject) viewer, stateId, outputState);
188:     if (name) {
189:       PetscViewerASCIIPrintf(viewer, "SCALARS %s double %d\n", name, bs);
190:     } else {
191:       PetscViewerASCIIPrintf(viewer, "SCALARS scalars double %d\n", bs);
192:     }
193:     PetscViewerASCIIPrintf(viewer, "LOOKUP_TABLE default\n");
194:     for (i=0; i<n/bs; i++) {
195:       for (b=0; b<bs; b++) {
196:         if (b > 0) {
197:           PetscViewerASCIIPrintf(viewer," ");
198:         }
199: #if !defined(PETSC_USE_COMPLEX)
200:         PetscViewerASCIIPrintf(viewer,"%G",x->array[i*bs+b]);
201: #endif
202:       }
203:       PetscViewerASCIIPrintf(viewer,"\n");
204:     }
205:   } else if (format == PETSC_VIEWER_ASCII_VTK_COORDS) {
206:     PetscInt bs, b;

208:     VecGetBlockSize(xin, &bs);
209:     if ((bs < 1) || (bs > 3)) {
210:       SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %d", bs);
211:     }
212:     for (i=0; i<n/bs; i++) {
213:       for (b=0; b<bs; b++) {
214:         if (b > 0) {
215:           PetscViewerASCIIPrintf(viewer," ");
216:         }
217: #if !defined(PETSC_USE_COMPLEX)
218:         PetscViewerASCIIPrintf(viewer,"%G",x->array[i*bs+b]);
219: #endif
220:       }
221:       for (b=bs; b<3; b++) {
222:         PetscViewerASCIIPrintf(viewer," 0.0");
223:       }
224:       PetscViewerASCIIPrintf(viewer,"\n");
225:     }
226:   } else if (format == PETSC_VIEWER_ASCII_PCICE) {
227:     PetscInt bs, b;

229:     VecGetBlockSize(xin, &bs);
230:     if ((bs < 1) || (bs > 3)) {
231:       SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "PCICE can only handle up to 3D objects, but vector dimension is %d", bs);
232:     }
233:     PetscViewerASCIIPrintf(viewer,"%D\n", xin->map->N/bs);
234:     for (i=0; i<n/bs; i++) {
235:       PetscViewerASCIIPrintf(viewer,"%7D   ", i+1);
236:       for (b=0; b<bs; b++) {
237:         if (b > 0) {
238:           PetscViewerASCIIPrintf(viewer," ");
239:         }
240: #if !defined(PETSC_USE_COMPLEX)
241:         PetscViewerASCIIPrintf(viewer,"% 12.5E",x->array[i*bs+b]);
242: #endif
243:       }
244:       PetscViewerASCIIPrintf(viewer,"\n");
245:     }
246:   } else if (format == PETSC_VIEWER_ASCII_PYLITH) {
247:     PetscInt bs, b;

249:     VecGetBlockSize(xin, &bs);
250:     if (bs != 3) {
251:       SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "PyLith can only handle 3D objects, but vector dimension is %d", bs);
252:     }
253:     PetscViewerASCIIPrintf(viewer,"#\n");
254:     PetscViewerASCIIPrintf(viewer,"#  Node      X-coord           Y-coord           Z-coord\n");
255:     PetscViewerASCIIPrintf(viewer,"#\n");
256:     for (i=0; i<n/bs; i++) {
257:       PetscViewerASCIIPrintf(viewer,"%7D ", i+1);
258:       for (b=0; b<bs; b++) {
259:         if (b > 0) {
260:           PetscViewerASCIIPrintf(viewer," ");
261:         }
262: #if !defined(PETSC_USE_COMPLEX)
263:         PetscViewerASCIIPrintf(viewer,"% 16.8E",x->array[i*bs+b]);
264: #endif
265:       }
266:       PetscViewerASCIIPrintf(viewer,"\n");
267:     }
268:   } else {
269:     for (i=0; i<n; i++) {
270:       if (format == PETSC_VIEWER_ASCII_INDEX) {
271:         PetscViewerASCIIPrintf(viewer,"%D: ",i);
272:       }
273: #if defined(PETSC_USE_COMPLEX)
274:       if (PetscImaginaryPart(x->array[i]) > 0.0) {
275:         PetscViewerASCIIPrintf(viewer,"%G + %G i\n",PetscRealPart(x->array[i]),PetscImaginaryPart(x->array[i]));
276:       } else if (PetscImaginaryPart(x->array[i]) < 0.0) {
277:         PetscViewerASCIIPrintf(viewer,"%G - %G i\n",PetscRealPart(x->array[i]),-PetscImaginaryPart(x->array[i]));
278:       } else {
279:         PetscViewerASCIIPrintf(viewer,"%G\n",PetscRealPart(x->array[i]));
280:       }
281: #else
282:       PetscViewerASCIIPrintf(viewer,"%G\n",x->array[i]);
283: #endif
284:     }
285:   }
286:   PetscViewerFlush(viewer);
287:   return(0);
288: }

292: static PetscErrorCode VecView_Seq_Draw_LG(Vec xin,PetscViewer v)
293: {
294:   Vec_Seq        *x = (Vec_Seq *)xin->data;
296:   PetscInt       i,n = xin->map->n;
297:   PetscDraw      win;
298:   PetscReal      *xx;
299:   PetscDrawLG    lg;

302:   PetscViewerDrawGetDrawLG(v,0,&lg);
303:   PetscDrawLGGetDraw(lg,&win);
304:   PetscDrawCheckResizedWindow(win);
305:   PetscDrawLGReset(lg);
306:   PetscMalloc((n+1)*sizeof(PetscReal),&xx);
307:   for (i=0; i<n; i++) {
308:     xx[i] = (PetscReal) i;
309:   }
310: #if !defined(PETSC_USE_COMPLEX)
311:   PetscDrawLGAddPoints(lg,n,&xx,&x->array);
312: #else 
313:   {
314:     PetscReal *yy;
315:     PetscMalloc((n+1)*sizeof(PetscReal),&yy);
316:     for (i=0; i<n; i++) {
317:       yy[i] = PetscRealPart(x->array[i]);
318:     }
319:     PetscDrawLGAddPoints(lg,n,&xx,&yy);
320:     PetscFree(yy);
321:   }
322: #endif
323:   PetscFree(xx);
324:   PetscDrawLGDraw(lg);
325:   PetscDrawSynchronizedFlush(win);
326:   return(0);
327: }

331: static PetscErrorCode VecView_Seq_Draw(Vec xin,PetscViewer v)
332: {
333:   PetscErrorCode    ierr;
334:   PetscDraw         draw;
335:   PetscTruth        isnull;
336:   PetscViewerFormat format;

339:   PetscViewerDrawGetDraw(v,0,&draw);
340:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
341: 
342:   PetscViewerGetFormat(v,&format);
343:   /*
344:      Currently it only supports drawing to a line graph */
345:   if (format != PETSC_VIEWER_DRAW_LG) {
346:     PetscViewerPushFormat(v,PETSC_VIEWER_DRAW_LG);
347:   }
348:   VecView_Seq_Draw_LG(xin,v);
349:   if (format != PETSC_VIEWER_DRAW_LG) {
350:     PetscViewerPopFormat(v);
351:   }

353:   return(0);
354: }

358: static PetscErrorCode VecView_Seq_Binary(Vec xin,PetscViewer viewer)
359: {
360:   Vec_Seq        *x = (Vec_Seq *)xin->data;
362:   int            fdes;
363:   PetscInt       n = xin->map->n,cookie=VEC_FILE_COOKIE;
364:   FILE           *file;
365: #if defined(PETSC_HAVE_MPIIO)
366:   PetscTruth     isMPIIO;
367: #endif


371:   /* Write vector header */
372:   PetscViewerBinaryWrite(viewer,&cookie,1,PETSC_INT,PETSC_FALSE);
373:   PetscViewerBinaryWrite(viewer,&n,1,PETSC_INT,PETSC_FALSE);

375:   /* Write vector contents */
376: #if defined(PETSC_HAVE_MPIIO)
377:   PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);
378:   if (!isMPIIO) {
379: #endif
380:     PetscViewerBinaryGetDescriptor(viewer,&fdes);
381:     PetscBinaryWrite(fdes,x->array,n,PETSC_SCALAR,PETSC_FALSE);
382: #if defined(PETSC_HAVE_MPIIO)
383:   } else {
384:     MPI_Offset   off;
385:     MPI_File     mfdes;
386:     PetscMPIInt  gsizes[1],lsizes[1],lstarts[1];
387:     MPI_Datatype view;

389:     gsizes[0]  = PetscMPIIntCast(n);
390:     lsizes[0]  = PetscMPIIntCast(n);
391:     lstarts[0] = 0;
392:     MPI_Type_create_subarray(1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);
393:     MPI_Type_commit(&view);

395:     PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);
396:     PetscViewerBinaryGetMPIIOOffset(viewer,&off);
397:     MPIU_File_write_all(mfdes,x->array,lsizes[0],MPIU_SCALAR,MPI_STATUS_IGNORE);
398:     PetscViewerBinaryAddMPIIOOffset(viewer,n*sizeof(PetscScalar));
399:     MPI_Type_free(&view);
400:   }
401: #endif

403:   PetscViewerBinaryGetInfoPointer(viewer,&file);
404:   if (file && xin->map->bs > 1) {
405:     if (((PetscObject)xin)->prefix) {
406:       PetscFPrintf(PETSC_COMM_SELF,file,"-%s_vecload_block_size %D\n",((PetscObject)xin)->prefix,xin->map->bs);
407:     } else {
408:       PetscFPrintf(PETSC_COMM_SELF,file,"-vecload_block_size %D\n",xin->map->bs);
409:     }
410:   }
411:   return(0);
412: }

414: #if defined(PETSC_HAVE_PNETCDF)
417: PetscErrorCode VecView_Seq_Netcdf(Vec xin,PetscViewer v)
418: {
420:   int            n = xin->map->n,ncid,xdim,xdim_num=1,xin_id,xstart=0;
421:   PetscScalar    *xarray;

424: #if !defined(PETSC_USE_COMPLEX)
425:   VecGetArray(xin,&xarray);
426:   PetscViewerNetcdfGetID(v,&ncid);
427:   if (ncid < 0) SETERRQ(PETSC_ERR_ORDER,"First call PetscViewerNetcdfOpen to create NetCDF dataset");
428:   /* define dimensions */
429:   ncmpi_def_dim(ncid,"PETSc_Vector_Global_Size",n,&xdim);
430:   /* define variables */
431:   ncmpi_def_var(ncid,"PETSc_Vector_Seq",NC_DOUBLE,xdim_num,&xdim,&xin_id);
432:   /* leave define mode */
433:   ncmpi_enddef(ncid);
434:   /* store the vector */
435:   VecGetOwnershipRange(xin,&xstart,PETSC_NULL);
436:   ncmpi_put_vara_double_all(ncid,xin_id,(const MPI_Offset*)&xstart,(const MPI_Offset*)&n,xarray);
437: #else 
438:     PetscPrintf(PETSC_COMM_WORLD,"NetCDF viewer not supported for complex numbers\n");
439: #endif
440:   return(0);
441: }
442: #endif

444: #if defined(PETSC_HAVE_MATLAB_ENGINE)
445: #include "mat.h"   /* Matlab include file */
449: PetscErrorCode VecView_Seq_Matlab(Vec vec,PetscViewer viewer)
450: {
452:   PetscInt       n;
453:   PetscScalar    *array;
454: 
456:   VecGetLocalSize(vec,&n);
457:   PetscObjectName((PetscObject)vec);
458:   VecGetArray(vec,&array);
459:   PetscViewerMatlabPutArray(viewer,n,1,array,((PetscObject)vec)->name);
460:   VecRestoreArray(vec,&array);
461:   return(0);
462: }
464: #endif

468: PetscErrorCode VecView_Seq(Vec xin,PetscViewer viewer)
469: {
471:   PetscTruth     isdraw,iascii,issocket,isbinary;
472: #if defined(PETSC_HAVE_MATHEMATICA)
473:   PetscTruth     ismathematica;
474: #endif
475: #if defined(PETSC_HAVE_PNETCDF)
476:   PetscTruth     isnetcdf;
477: #endif
478: #if defined(PETSC_HAVE_MATLAB_ENGINE)
479:   PetscTruth     ismatlab;
480: #endif

483:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
484:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
485:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
486:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
487: #if defined(PETSC_HAVE_MATHEMATICA)
488:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_MATHEMATICA,&ismathematica);
489: #endif
490: #if defined(PETSC_HAVE_PNETCDF)
491:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_NETCDF,&isnetcdf);
492: #endif
493: #if defined(PETSC_HAVE_MATLAB_ENGINE)
494:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_MATLAB,&ismatlab);
495: #endif

497:   if (isdraw){
498:     VecView_Seq_Draw(xin,viewer);
499:   } else if (iascii){
500:     VecView_Seq_File(xin,viewer);
501:   } else if (isbinary) {
502:     VecView_Seq_Binary(xin,viewer);
503: #if defined(PETSC_HAVE_MATHEMATICA)
504:   } else if (ismathematica) {
505:     PetscViewerMathematicaPutVector(viewer,xin);
506: #endif
507: #if defined(PETSC_HAVE_PNETCDF)
508:   } else if (isnetcdf) {
509:     VecView_Seq_Netcdf(xin,viewer);
510: #endif
511: #if defined(PETSC_HAVE_MATLAB_ENGINE)
512:   } else if (ismatlab) {
513:     VecView_Seq_Matlab(xin,viewer);
514: #endif
515:   } else {
516:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by this vector object",((PetscObject)viewer)->type_name);
517:   }
518:   return(0);
519: }

523: PetscErrorCode VecGetValues_Seq(Vec xin,PetscInt ni,const PetscInt ix[],PetscScalar y[])
524: {
525:   Vec_Seq     *x = (Vec_Seq *)xin->data;
526:   PetscScalar *xx = x->array;
527:   PetscInt    i;

530:   for (i=0; i<ni; i++) {
531:     if (xin->stash.ignorenegidx && ix[i] < 0) continue;
532: #if defined(PETSC_USE_DEBUG)
533:     if (ix[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
534:     if (ix[i] >= xin->map->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D to large maximum allowed %D",ix[i],xin->map->n);
535: #endif
536:     y[i] = xx[ix[i]];
537:   }
538:   return(0);
539: }

543: PetscErrorCode VecSetValues_Seq(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode m)
544: {
545:   Vec_Seq     *x = (Vec_Seq *)xin->data;
546:   PetscScalar *xx = x->array;
547:   PetscInt    i;

550:   if (m == INSERT_VALUES) {
551:     for (i=0; i<ni; i++) {
552:       if (xin->stash.ignorenegidx && ix[i] < 0) continue;
553: #if defined(PETSC_USE_DEBUG)
554:       if (ix[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
555:       if (ix[i] >= xin->map->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",ix[i],xin->map->n);
556: #endif
557:       xx[ix[i]] = y[i];
558:     }
559:   } else {
560:     for (i=0; i<ni; i++) {
561:       if (xin->stash.ignorenegidx && ix[i] < 0) continue;
562: #if defined(PETSC_USE_DEBUG)
563:       if (ix[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
564:       if (ix[i] >= xin->map->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",ix[i],xin->map->n);
565: #endif
566:       xx[ix[i]] += y[i];
567:     }
568:   }
569:   return(0);
570: }

574: PetscErrorCode VecSetValuesBlocked_Seq(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar yin[],InsertMode m)
575: {
576:   Vec_Seq     *x = (Vec_Seq *)xin->data;
577:   PetscScalar *xx = x->array,*y = (PetscScalar*)yin;
578:   PetscInt    i,bs = xin->map->bs,start,j;

580:   /*
581:        For optimization could treat bs = 2, 3, 4, 5 as special cases with loop unrolling
582:   */
584:   if (m == INSERT_VALUES) {
585:     for (i=0; i<ni; i++) {
586:       start = bs*ix[i];
587:       if (start < 0) continue;
588: #if defined(PETSC_USE_DEBUG)
589:       if (start >= xin->map->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",start,xin->map->n);
590: #endif
591:       for (j=0; j<bs; j++) {
592:         xx[start+j] = y[j];
593:       }
594:       y += bs;
595:     }
596:   } else {
597:     for (i=0; i<ni; i++) {
598:       start = bs*ix[i];
599:       if (start < 0) continue;
600: #if defined(PETSC_USE_DEBUG)
601:       if (start >= xin->map->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",start,xin->map->n);
602: #endif
603:       for (j=0; j<bs; j++) {
604:         xx[start+j] += y[j];
605:       }
606:       y += bs;
607:     }
608:   }
609:   return(0);
610: }


615: PetscErrorCode VecDestroy_Seq(Vec v)
616: {
617:   Vec_Seq        *vs = (Vec_Seq*)v->data;


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

625: #if defined(PETSC_USE_LOG)
626:   PetscLogObjectState((PetscObject)v,"Length=%D",v->map->n);
627: #endif
628:   PetscFree(vs->array_allocated);
629:   PetscFree(vs);

631:   return(0);
632: }

634: #if 0
637: static PetscErrorCode VecPublish_Seq(PetscObject obj)
638: {
640:   return(0);
641: }
642: #endif

644: EXTERN PetscErrorCode VecLoad_Binary(PetscViewer, const VecType, Vec*);

646: static struct _VecOps DvOps = {VecDuplicate_Seq, /* 1 */
647:             VecDuplicateVecs_Default,
648:             VecDestroyVecs_Default,
649:             VecDot_Seq,
650:             VecMDot_Seq,
651:             VecNorm_Seq,
652:             VecTDot_Seq,
653:             VecMTDot_Seq,
654:             VecScale_Seq,
655:             VecCopy_Seq, /* 10 */
656:             VecSet_Seq,
657:             VecSwap_Seq,
658:             VecAXPY_Seq,
659:             VecAXPBY_Seq,
660:             VecMAXPY_Seq,
661:             VecAYPX_Seq,
662:             VecWAXPY_Seq,
663:             VecAXPBYPCZ_Seq,
664:             VecPointwiseMult_Seq,
665:             VecPointwiseDivide_Seq,
666:             VecSetValues_Seq, /* 20 */
667:             0,0,
668:             VecGetArray_Seq,
669:             VecGetSize_Seq,
670:             VecGetSize_Seq,
671:             VecRestoreArray_Seq,
672:             VecMax_Seq,
673:             VecMin_Seq,
674:             VecSetRandom_Seq,
675:             VecSetOption_Seq, /* 30 */
676:             VecSetValuesBlocked_Seq,
677:             VecDestroy_Seq,
678:             VecView_Seq,
679:             VecPlaceArray_Seq,
680:             VecReplaceArray_Seq,
681:             VecDot_Seq,
682:             VecTDot_Seq,
683:             VecNorm_Seq,
684:             VecMDot_Seq,
685:             VecMTDot_Seq, /* 40 */
686:             VecLoadIntoVector_Default,
687:             0, /* VecLoadIntoVectorNative */
688:             VecReciprocal_Default,
689:             0, /* VecViewNative */
690:             VecConjugate_Seq,
691:             0,
692:             0,
693:             VecResetArray_Seq,
694:             0,
695:             VecMaxPointwiseDivide_Seq,
696:             VecLoad_Binary, /* 50 */
697:             VecPointwiseMax_Seq,
698:             VecPointwiseMaxAbs_Seq,
699:             VecPointwiseMin_Seq,
700:             VecGetValues_Seq
701:           };


704: /*
705:       This is called by VecCreate_Seq() (i.e. VecCreateSeq()) and VecCreateSeqWithArray()
706: */
709: static PetscErrorCode VecCreate_Seq_Private(Vec v,const PetscScalar array[])
710: {
711:   Vec_Seq        *s;

715:   PetscNewLog(v,Vec_Seq,&s);
716:   PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));
717:   v->data            = (void*)s;
718:   v->petscnative     = PETSC_TRUE;
719:   s->array           = (PetscScalar *)array;
720:   s->array_allocated = 0;

722:   if (v->map->bs == -1) v->map->bs = 1;
723:   PetscLayoutSetUp(v->map);
724:   PetscObjectChangeTypeName((PetscObject)v,VECSEQ);
725: #if defined(PETSC_HAVE_MATLAB_ENGINE)
726:   PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEnginePut_C","VecMatlabEnginePut_Default",VecMatlabEnginePut_Default);
727:   PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEngineGet_C","VecMatlabEngineGet_Default",VecMatlabEngineGet_Default);
728: #endif
729:   PetscPublishAll(v);
730:   return(0);
731: }

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

739:    Collective on MPI_Comm

741:    Input Parameter:
742: +  comm - the communicator, should be PETSC_COMM_SELF
743: .  n - the vector length 
744: -  array - memory where the vector elements are to be stored.

746:    Output Parameter:
747: .  V - the vector

749:    Notes:
750:    Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
751:    same type as an existing vector.

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

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

759:    Level: intermediate

761:    Concepts: vectors^creating with array

763: .seealso: VecCreateMPIWithArray(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), 
764:           VecCreateGhost(), VecCreateSeq(), VecPlaceArray()
765: @*/
766: PetscErrorCode  VecCreateSeqWithArray(MPI_Comm comm,PetscInt n,const PetscScalar array[],Vec *V)
767: {
769:   PetscMPIInt    size;

772:   VecCreate(comm,V);
773:   VecSetSizes(*V,n,n);
774:   MPI_Comm_size(comm,&size);
775:   if (size > 1) {
776:     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot create VECSEQ on more than one process");
777:   }
778:   VecCreate_Seq_Private(*V,array);
779:   return(0);
780: }

782: /*MC
783:    VECSEQ - VECSEQ = "seq" - The basic sequential vector

785:    Options Database Keys:
786: . -vec_type seq - sets the vector type to VECSEQ during a call to VecSetFromOptions()

788:   Level: beginner

790: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateSeqWithArray(), VECMPI, VecType, VecCreateMPI(), VecCreateSeq()
791: M*/

796: PetscErrorCode  VecCreate_Seq(Vec V)
797: {
798:   Vec_Seq        *s;
799:   PetscScalar    *array;
801:   PetscInt       n = PetscMax(V->map->n,V->map->N);
802:   PetscMPIInt    size;

805:   MPI_Comm_size(((PetscObject)V)->comm,&size);
806:   if (size > 1) {
807:     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot create VECSEQ on more than one process");
808:   }
809:   PetscMalloc(n*sizeof(PetscScalar),&array);
810:   PetscLogObjectMemory(V, n*sizeof(PetscScalar));
811:   PetscMemzero(array,n*sizeof(PetscScalar));
812:   VecCreate_Seq_Private(V,array);
813:   s    = (Vec_Seq*)V->data;
814:   s->array_allocated = array;
815:   return(0);
816: }


822: PetscErrorCode VecDuplicate_Seq(Vec win,Vec *V)
823: {

827:   VecCreateSeq(((PetscObject)win)->comm,win->map->n,V);
828:   if (win->mapping) {
829:     PetscObjectReference((PetscObject)win->mapping);
830:     (*V)->mapping = win->mapping;
831:   }
832:   if (win->bmapping) {
833:     PetscObjectReference((PetscObject)win->bmapping);
834:     (*V)->bmapping = win->bmapping;
835:   }
836:   (*V)->map->bs = win->map->bs;
837:   PetscOListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*V))->olist);
838:   PetscFListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*V))->qlist);

840:   (*V)->stash.ignorenegidx = win->stash.ignorenegidx;

842:   return(0);
843: }

847: PetscErrorCode VecSetOption_Seq(Vec v,VecOption op,PetscTruth flag)
848: {
850:   if (op == VEC_IGNORE_NEGATIVE_INDICES) {
851:     v->stash.ignorenegidx = flag;
852:   }
853:   return(0);
854: }