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: #if defined(PETSC_HAVE_HDF5)
 17: #endif

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

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

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

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

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

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

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

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

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

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

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

335: static PetscErrorCode VecView_Seq_Draw(Vec xin,PetscViewer v)
336: {
337:   PetscErrorCode    ierr;
338:   PetscDraw         draw;
339:   PetscTruth        isnull;
340:   PetscViewerFormat format;

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

357:   return(0);
358: }

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


375:   /* Write vector header */
376:   PetscViewerBinaryWrite(viewer,&cookie,1,PETSC_INT,PETSC_FALSE);
377:   PetscViewerBinaryWrite(viewer,&n,1,PETSC_INT,PETSC_FALSE);

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

393:     gsizes[0]  = PetscMPIIntCast(n);
394:     lsizes[0]  = PetscMPIIntCast(n);
395:     lstarts[0] = 0;
396:     MPI_Type_create_subarray(1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);
397:     MPI_Type_commit(&view);

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

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

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

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

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

472: PetscErrorCode VecView_Seq(Vec xin,PetscViewer viewer)
473: {
475:   PetscTruth     isdraw,iascii,issocket,isbinary;
476: #if defined(PETSC_HAVE_MATHEMATICA)
477:   PetscTruth     ismathematica;
478: #endif
479: #if defined(PETSC_HAVE_PNETCDF)
480:   PetscTruth     isnetcdf;
481: #endif
482: #if defined(PETSC_HAVE_MATLAB_ENGINE)
483:   PetscTruth     ismatlab;
484: #endif
485: #if defined(PETSC_HAVE_HDF5)
486:   PetscTruth     ishdf5;
487: #endif

490:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
491:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
492:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
493:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
494: #if defined(PETSC_HAVE_MATHEMATICA)
495:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_MATHEMATICA,&ismathematica);
496: #endif
497: #if defined(PETSC_HAVE_HDF5)
498:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_HDF5,&ishdf5);
499: #endif
500: #if defined(PETSC_HAVE_PNETCDF)
501:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_NETCDF,&isnetcdf);
502: #endif
503: #if defined(PETSC_HAVE_MATLAB_ENGINE)
504:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_MATLAB,&ismatlab);
505: #endif

507:   if (isdraw){
508:     VecView_Seq_Draw(xin,viewer);
509:   } else if (iascii){
510:     VecView_Seq_File(xin,viewer);
511:   } else if (isbinary) {
512:     VecView_Seq_Binary(xin,viewer);
513: #if defined(PETSC_HAVE_MATHEMATICA)
514:   } else if (ismathematica) {
515:     PetscViewerMathematicaPutVector(viewer,xin);
516: #endif
517: #if defined(PETSC_HAVE_HDF5)
518:   } else if (ishdf5) {
519:     VecView_MPI_HDF5(xin,viewer); /* Reusing VecView_MPI_HDF5 ... don't want code duplication*/
520: #endif
521: #if defined(PETSC_HAVE_PNETCDF)
522:   } else if (isnetcdf) {
523:     VecView_Seq_Netcdf(xin,viewer);
524: #endif
525: #if defined(PETSC_HAVE_MATLAB_ENGINE)
526:   } else if (ismatlab) {
527:     VecView_Seq_Matlab(xin,viewer);
528: #endif
529:   } else {
530:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by this vector object",((PetscObject)viewer)->type_name);
531:   }
532:   return(0);
533: }

537: PetscErrorCode VecGetValues_Seq(Vec xin,PetscInt ni,const PetscInt ix[],PetscScalar y[])
538: {
539:   Vec_Seq     *x = (Vec_Seq *)xin->data;
540:   PetscScalar *xx = x->array;
541:   PetscInt    i;

544:   for (i=0; i<ni; i++) {
545:     if (xin->stash.ignorenegidx && ix[i] < 0) continue;
546: #if defined(PETSC_USE_DEBUG)
547:     if (ix[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
548:     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);
549: #endif
550:     y[i] = xx[ix[i]];
551:   }
552:   return(0);
553: }

557: PetscErrorCode VecSetValues_Seq(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode m)
558: {
559:   Vec_Seq     *x = (Vec_Seq *)xin->data;
560:   PetscScalar *xx = x->array;
561:   PetscInt    i;

564:   if (m == INSERT_VALUES) {
565:     for (i=0; i<ni; i++) {
566:       if (xin->stash.ignorenegidx && ix[i] < 0) continue;
567: #if defined(PETSC_USE_DEBUG)
568:       if (ix[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
569:       if (ix[i] >= xin->map->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",ix[i],xin->map->n);
570: #endif
571:       xx[ix[i]] = y[i];
572:     }
573:   } else {
574:     for (i=0; i<ni; i++) {
575:       if (xin->stash.ignorenegidx && ix[i] < 0) continue;
576: #if defined(PETSC_USE_DEBUG)
577:       if (ix[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
578:       if (ix[i] >= xin->map->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",ix[i],xin->map->n);
579: #endif
580:       xx[ix[i]] += y[i];
581:     }
582:   }
583:   return(0);
584: }

588: PetscErrorCode VecSetValuesBlocked_Seq(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar yin[],InsertMode m)
589: {
590:   Vec_Seq     *x = (Vec_Seq *)xin->data;
591:   PetscScalar *xx = x->array,*y = (PetscScalar*)yin;
592:   PetscInt    i,bs = xin->map->bs,start,j;

594:   /*
595:        For optimization could treat bs = 2, 3, 4, 5 as special cases with loop unrolling
596:   */
598:   if (m == INSERT_VALUES) {
599:     for (i=0; i<ni; i++) {
600:       start = bs*ix[i];
601:       if (start < 0) continue;
602: #if defined(PETSC_USE_DEBUG)
603:       if (start >= xin->map->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",start,xin->map->n);
604: #endif
605:       for (j=0; j<bs; j++) {
606:         xx[start+j] = y[j];
607:       }
608:       y += bs;
609:     }
610:   } else {
611:     for (i=0; i<ni; i++) {
612:       start = bs*ix[i];
613:       if (start < 0) continue;
614: #if defined(PETSC_USE_DEBUG)
615:       if (start >= xin->map->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",start,xin->map->n);
616: #endif
617:       for (j=0; j<bs; j++) {
618:         xx[start+j] += y[j];
619:       }
620:       y += bs;
621:     }
622:   }
623:   return(0);
624: }


629: PetscErrorCode VecDestroy_Seq(Vec v)
630: {
631:   Vec_Seq        *vs = (Vec_Seq*)v->data;


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

639: #if defined(PETSC_USE_LOG)
640:   PetscLogObjectState((PetscObject)v,"Length=%D",v->map->n);
641: #endif
642:   PetscFree(vs->array_allocated);
643:   PetscFree(vs);

645:   return(0);
646: }

648: #if 0
651: static PetscErrorCode VecPublish_Seq(PetscObject obj)
652: {
654:   return(0);
655: }
656: #endif

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

660: static struct _VecOps DvOps = {VecDuplicate_Seq, /* 1 */
661:             VecDuplicateVecs_Default,
662:             VecDestroyVecs_Default,
663:             VecDot_Seq,
664:             VecMDot_Seq,
665:             VecNorm_Seq,
666:             VecTDot_Seq,
667:             VecMTDot_Seq,
668:             VecScale_Seq,
669:             VecCopy_Seq, /* 10 */
670:             VecSet_Seq,
671:             VecSwap_Seq,
672:             VecAXPY_Seq,
673:             VecAXPBY_Seq,
674:             VecMAXPY_Seq,
675:             VecAYPX_Seq,
676:             VecWAXPY_Seq,
677:             VecAXPBYPCZ_Seq,
678:             VecPointwiseMult_Seq,
679:             VecPointwiseDivide_Seq,
680:             VecSetValues_Seq, /* 20 */
681:             0,0,
682:             VecGetArray_Seq,
683:             VecGetSize_Seq,
684:             VecGetSize_Seq,
685:             VecRestoreArray_Seq,
686:             VecMax_Seq,
687:             VecMin_Seq,
688:             VecSetRandom_Seq,
689:             VecSetOption_Seq, /* 30 */
690:             VecSetValuesBlocked_Seq,
691:             VecDestroy_Seq,
692:             VecView_Seq,
693:             VecPlaceArray_Seq,
694:             VecReplaceArray_Seq,
695:             VecDot_Seq,
696:             VecTDot_Seq,
697:             VecNorm_Seq,
698:             VecMDot_Seq,
699:             VecMTDot_Seq, /* 40 */
700:             VecLoadIntoVector_Default,
701:             0, /* VecLoadIntoVectorNative */
702:             VecReciprocal_Default,
703:             0, /* VecViewNative */
704:             VecConjugate_Seq,
705:             0,
706:             0,
707:             VecResetArray_Seq,
708:             0,
709:             VecMaxPointwiseDivide_Seq,
710:             VecLoad_Binary, /* 50 */
711:             VecPointwiseMax_Seq,
712:             VecPointwiseMaxAbs_Seq,
713:             VecPointwiseMin_Seq,
714:             VecGetValues_Seq
715:           };


718: /*
719:       This is called by VecCreate_Seq() (i.e. VecCreateSeq()) and VecCreateSeqWithArray()
720: */
723: static PetscErrorCode VecCreate_Seq_Private(Vec v,const PetscScalar array[])
724: {
725:   Vec_Seq        *s;

729:   PetscNewLog(v,Vec_Seq,&s);
730:   PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));
731:   v->data            = (void*)s;
732:   v->petscnative     = PETSC_TRUE;
733:   s->array           = (PetscScalar *)array;
734:   s->array_allocated = 0;

736:   if (v->map->bs == -1) v->map->bs = 1;
737:   PetscLayoutSetUp(v->map);
738:   PetscObjectChangeTypeName((PetscObject)v,VECSEQ);
739: #if defined(PETSC_HAVE_MATLAB_ENGINE)
740:   PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEnginePut_C","VecMatlabEnginePut_Default",VecMatlabEnginePut_Default);
741:   PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEngineGet_C","VecMatlabEngineGet_Default",VecMatlabEngineGet_Default);
742: #endif
743:   PetscPublishAll(v);
744:   return(0);
745: }

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

753:    Collective on MPI_Comm

755:    Input Parameter:
756: +  comm - the communicator, should be PETSC_COMM_SELF
757: .  n - the vector length 
758: -  array - memory where the vector elements are to be stored.

760:    Output Parameter:
761: .  V - the vector

763:    Notes:
764:    Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
765:    same type as an existing vector.

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

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

773:    Level: intermediate

775:    Concepts: vectors^creating with array

777: .seealso: VecCreateMPIWithArray(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), 
778:           VecCreateGhost(), VecCreateSeq(), VecPlaceArray()
779: @*/
780: PetscErrorCode  VecCreateSeqWithArray(MPI_Comm comm,PetscInt n,const PetscScalar array[],Vec *V)
781: {
783:   PetscMPIInt    size;

786:   VecCreate(comm,V);
787:   VecSetSizes(*V,n,n);
788:   MPI_Comm_size(comm,&size);
789:   if (size > 1) {
790:     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot create VECSEQ on more than one process");
791:   }
792:   VecCreate_Seq_Private(*V,array);
793:   return(0);
794: }

796: /*MC
797:    VECSEQ - VECSEQ = "seq" - The basic sequential vector

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

802:   Level: beginner

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

810: PetscErrorCode  VecCreate_Seq(Vec V)
811: {
812:   Vec_Seq        *s;
813:   PetscScalar    *array;
815:   PetscInt       n = PetscMax(V->map->n,V->map->N);
816:   PetscMPIInt    size;

819:   MPI_Comm_size(((PetscObject)V)->comm,&size);
820:   if (size > 1) {
821:     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot create VECSEQ on more than one process");
822:   }
823:   PetscMalloc(n*sizeof(PetscScalar),&array);
824:   PetscLogObjectMemory(V, n*sizeof(PetscScalar));
825:   PetscMemzero(array,n*sizeof(PetscScalar));
826:   VecCreate_Seq_Private(V,array);
827:   s    = (Vec_Seq*)V->data;
828:   s->array_allocated = array;
829:   return(0);
830: }


836: PetscErrorCode VecDuplicate_Seq(Vec win,Vec *V)
837: {

841:   VecCreateSeq(((PetscObject)win)->comm,win->map->n,V);
842:   if (win->mapping) {
843:     PetscObjectReference((PetscObject)win->mapping);
844:     (*V)->mapping = win->mapping;
845:   }
846:   if (win->bmapping) {
847:     PetscObjectReference((PetscObject)win->bmapping);
848:     (*V)->bmapping = win->bmapping;
849:   }
850:   (*V)->map->bs = win->map->bs;
851:   PetscOListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*V))->olist);
852:   PetscFListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*V))->qlist);

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

856:   return(0);
857: }

861: PetscErrorCode VecSetOption_Seq(Vec v,VecOption op,PetscTruth flag)
862: {
864:   if (op == VEC_IGNORE_NEGATIVE_INDICES) {
865:     v->stash.ignorenegidx = flag;
866:   }
867:   return(0);
868: }