Actual source code: bvec2.c
2: /*
3: Implements the sequential vectors.
4: */
6: #include <private/vecimpl.h> /*I "petscvec.h" I*/
7: #include <../src/vec/vec/impls/dvecimpl.h>
8: #include <../src/vec/vec/impls/mpi/pvecimpl.h> /* For VecView_MPI_HDF5 */
9: #include <petscblaslapack.h>
11: #if defined(PETSC_HAVE_HDF5)
13: #endif
17: static PetscErrorCode VecPointwiseMax_Seq(Vec win,Vec xin,Vec yin)
18: {
20: PetscInt n = win->map->n,i;
21: PetscScalar *ww,*xx,*yy; /* cannot make xx or yy const since might be ww */
24: VecGetArrayRead(xin,(const PetscScalar**)&xx);
25: VecGetArrayRead(yin,(const PetscScalar**)&yy);
26: VecGetArray(win,&ww);
27: for (i=0; i<n; i++) {
28: ww[i] = PetscMax(PetscRealPart(xx[i]),PetscRealPart(yy[i]));
29: }
30: VecRestoreArrayRead(xin,(const PetscScalar**)&xx);
31: VecRestoreArrayRead(yin,(const PetscScalar**)&yy);
32: VecRestoreArray(win,&ww);
33: PetscLogFlops(n);
34: return(0);
35: }
39: static PetscErrorCode VecPointwiseMin_Seq(Vec win,Vec xin,Vec yin)
40: {
42: PetscInt n = win->map->n,i;
43: PetscScalar *ww,*xx,*yy; /* cannot make xx or yy const since might be ww */
46: VecGetArrayRead(xin,(const PetscScalar**)&xx);
47: VecGetArrayRead(yin,(const PetscScalar**)&yy);
48: VecGetArray(win,&ww);
49: for (i=0; i<n; i++) {
50: ww[i] = PetscMin(PetscRealPart(xx[i]),PetscRealPart(yy[i]));
51: }
52: VecRestoreArrayRead(xin,(const PetscScalar**)&xx);
53: VecRestoreArrayRead(yin,(const PetscScalar**)&yy);
54: VecRestoreArray(win,&ww);
55: PetscLogFlops(n);
56: return(0);
57: }
61: static PetscErrorCode VecPointwiseMaxAbs_Seq(Vec win,Vec xin,Vec yin)
62: {
64: PetscInt n = win->map->n,i;
65: PetscScalar *ww,*xx,*yy; /* cannot make xx or yy const since might be ww */
68: VecGetArrayRead(xin,(const PetscScalar**)&xx);
69: VecGetArrayRead(yin,(const PetscScalar**)&yy);
70: VecGetArray(win,&ww);
71: for (i=0; i<n; i++) {
72: ww[i] = PetscMax(PetscAbsScalar(xx[i]),PetscAbsScalar(yy[i]));
73: }
74: PetscLogFlops(n);
75: VecRestoreArrayRead(xin,(const PetscScalar**)&xx);
76: VecRestoreArrayRead(yin,(const PetscScalar**)&yy);
77: VecRestoreArray(win,&ww);
78: return(0);
79: }
81: #include <../src/vec/vec/impls/seq/ftn-kernels/fxtimesy.h>
84: static PetscErrorCode VecPointwiseMult_Seq(Vec win,Vec xin,Vec yin)
85: {
87: PetscInt n = win->map->n,i;
88: PetscScalar *ww,*xx,*yy; /* cannot make xx or yy const since might be ww */
91: VecGetArrayRead(xin,(const PetscScalar**)&xx);
92: VecGetArrayRead(yin,(const PetscScalar**)&yy);
93: VecGetArray(win,&ww);
94: if (ww == xx) {
95: for (i=0; i<n; i++) ww[i] *= yy[i];
96: } else if (ww == yy) {
97: for (i=0; i<n; i++) ww[i] *= xx[i];
98: } else {
99: #if defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
100: fortranxtimesy_(xx,yy,ww,&n);
101: #else
102: for (i=0; i<n; i++) ww[i] = xx[i] * yy[i];
103: #endif
104: }
105: VecRestoreArrayRead(xin,(const PetscScalar**)&xx);
106: VecRestoreArrayRead(yin,(const PetscScalar**)&yy);
107: VecRestoreArray(win,&ww);
108: PetscLogFlops(n);
109: return(0);
110: }
114: static PetscErrorCode VecPointwiseDivide_Seq(Vec win,Vec xin,Vec yin)
115: {
117: PetscInt n = win->map->n,i;
118: PetscScalar *ww,*xx,*yy; /* cannot make xx or yy const since might be ww */
121: VecGetArrayRead(xin,(const PetscScalar**)&xx);
122: VecGetArrayRead(yin,(const PetscScalar**)&yy);
123: VecGetArray(win,&ww);
124: for (i=0; i<n; i++) {
125: ww[i] = xx[i] / yy[i];
126: }
127: PetscLogFlops(n);
128: VecRestoreArrayRead(xin,(const PetscScalar**)&xx);
129: VecRestoreArrayRead(yin,(const PetscScalar**)&yy);
130: VecRestoreArray(win,&ww);
131: return(0);
132: }
136: static PetscErrorCode VecSetRandom_Seq(Vec xin,PetscRandom r)
137: {
139: PetscInt n = xin->map->n,i;
140: PetscScalar *xx;
143: VecGetArray(xin,&xx);
144: for (i=0; i<n; i++) {PetscRandomGetValue(r,&xx[i]);}
145: VecRestoreArray(xin,&xx);
146: return(0);
147: }
151: static PetscErrorCode VecGetSize_Seq(Vec vin,PetscInt *size)
152: {
154: *size = vin->map->n;
155: return(0);
156: }
160: static PetscErrorCode VecConjugate_Seq(Vec xin)
161: {
162: PetscScalar *x;
163: PetscInt n = xin->map->n;
167: VecGetArray(xin,&x);
168: while (n-->0) {
169: *x = PetscConj(*x);
170: x++;
171: }
172: VecRestoreArray(xin,&x);
173: return(0);
174: }
178: static PetscErrorCode VecResetArray_Seq(Vec vin)
179: {
180: Vec_Seq *v = (Vec_Seq *)vin->data;
183: v->array = v->unplacedarray;
184: v->unplacedarray = 0;
185: return(0);
186: }
190: static PetscErrorCode VecCopy_Seq(Vec xin,Vec yin)
191: {
192: PetscScalar *ya;
193: const PetscScalar *xa;
194: PetscErrorCode ierr;
197: if (xin != yin) {
198: VecGetArrayRead(xin,&xa);
199: VecGetArray(yin,&ya);
200: PetscMemcpy(ya,xa,xin->map->n*sizeof(PetscScalar));
201: VecRestoreArrayRead(xin,&xa);
202: VecRestoreArray(yin,&ya);
203: }
204: return(0);
205: }
209: static PetscErrorCode VecSwap_Seq(Vec xin,Vec yin)
210: {
211: PetscScalar *ya, *xa;
213: PetscBLASInt one = 1,bn = PetscBLASIntCast(xin->map->n);
216: if (xin != yin) {
217: VecGetArray(xin,&xa);
218: VecGetArray(yin,&ya);
219: BLASswap_(&bn,xa,&one,ya,&one);
220: VecRestoreArray(xin,&xa);
221: VecRestoreArray(yin,&ya);
222: }
223: return(0);
224: }
226: #include <../src/vec/vec/impls/seq/ftn-kernels/fnorm.h>
229: PetscErrorCode VecNorm_Seq(Vec xin,NormType type,PetscReal* z)
230: {
231: const PetscScalar *xx;
232: PetscErrorCode ierr;
233: PetscInt n = xin->map->n;
234: PetscBLASInt one = 1, bn = PetscBLASIntCast(n);
237: if (type == NORM_2 || type == NORM_FROBENIUS) {
238: VecGetArrayRead(xin,&xx);
239: *z = BLASnrm2_(&bn,xx,&one);
240: VecRestoreArrayRead(xin,&xx);
241: PetscLogFlops(PetscMax(2.0*n-1,0.0));
242: } else if (type == NORM_INFINITY) {
243: PetscInt i;
244: PetscReal max = 0.0,tmp;
246: VecGetArrayRead(xin,&xx);
247: for (i=0; i<n; i++) {
248: if ((tmp = PetscAbsScalar(*xx)) > max) max = tmp;
249: /* check special case of tmp == NaN */
250: if (tmp != tmp) {max = tmp; break;}
251: xx++;
252: }
253: VecRestoreArrayRead(xin,&xx);
254: *z = max;
255: } else if (type == NORM_1) {
256: VecGetArrayRead(xin,&xx);
257: *z = BLASasum_(&bn,xx,&one);
258: VecRestoreArrayRead(xin,&xx);
259: PetscLogFlops(PetscMax(n-1.0,0.0));
260: } else if (type == NORM_1_AND_2) {
261: VecNorm_Seq(xin,NORM_1,z);
262: VecNorm_Seq(xin,NORM_2,z+1);
263: }
264: return(0);
265: }
269: static PetscErrorCode VecView_Seq_ASCII(Vec xin,PetscViewer viewer)
270: {
271: PetscErrorCode ierr;
272: PetscInt i,n = xin->map->n;
273: const char *name;
274: PetscViewerFormat format;
275: const PetscScalar *xv;
278: VecGetArrayRead(xin,&xv);
279: PetscViewerGetFormat(viewer,&format);
280: if (format == PETSC_VIEWER_ASCII_MATLAB) {
281: PetscObjectGetName((PetscObject)xin,&name);
282: PetscViewerASCIIPrintf(viewer,"%s = [\n",name);
283: for (i=0; i<n; i++) {
284: #if defined(PETSC_USE_COMPLEX)
285: if (PetscImaginaryPart(xv[i]) > 0.0) {
286: PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei\n",PetscRealPart(xv[i]),PetscImaginaryPart(xv[i]));
287: } else if (PetscImaginaryPart(xv[i]) < 0.0) {
288: PetscViewerASCIIPrintf(viewer,"%18.16e - %18.16ei\n",PetscRealPart(xv[i]),-PetscImaginaryPart(xv[i]));
289: } else {
290: PetscViewerASCIIPrintf(viewer,"%18.16e\n",PetscRealPart(xv[i]));
291: }
292: #else
293: PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double) xv[i]);
294: #endif
295: }
296: PetscViewerASCIIPrintf(viewer,"];\n");
297: } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
298: for (i=0; i<n; i++) {
299: #if defined(PETSC_USE_COMPLEX)
300: PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e\n",PetscRealPart(xv[i]),PetscImaginaryPart(xv[i]));
301: #else
302: PetscViewerASCIIPrintf(viewer,"%18.16e\n",xv[i]);
303: #endif
304: }
305: } else if (format == PETSC_VIEWER_ASCII_VTK || format == PETSC_VIEWER_ASCII_VTK_CELL) {
306: /*
307: state 0: No header has been output
308: state 1: Only POINT_DATA has been output
309: state 2: Only CELL_DATA has been output
310: state 3: Output both, POINT_DATA last
311: state 4: Output both, CELL_DATA last
312: */
313: static PetscInt stateId = -1;
314: int outputState = 0;
315: PetscBool hasState;
316: int doOutput = 0;
317: PetscInt bs, b;
319: if (stateId < 0) {
320: PetscObjectComposedDataRegister(&stateId);
321: }
322: PetscObjectComposedDataGetInt((PetscObject) viewer, stateId, outputState, hasState);
323: if (!hasState) {
324: outputState = 0;
325: }
326: PetscObjectGetName((PetscObject) xin, &name);
327: VecGetBlockSize(xin, &bs);
328: if ((bs < 1) || (bs > 3)) {
329: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %d", bs);
330: }
331: if (format == PETSC_VIEWER_ASCII_VTK) {
332: if (outputState == 0) {
333: outputState = 1;
334: doOutput = 1;
335: } else if (outputState == 1) {
336: doOutput = 0;
337: } else if (outputState == 2) {
338: outputState = 3;
339: doOutput = 1;
340: } else if (outputState == 3) {
341: doOutput = 0;
342: } else if (outputState == 4) {
343: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Tried to output POINT_DATA again after intervening CELL_DATA");
344: }
345: if (doOutput) {
346: PetscViewerASCIIPrintf(viewer, "POINT_DATA %d\n", n/bs);
347: }
348: } else {
349: if (outputState == 0) {
350: outputState = 2;
351: doOutput = 1;
352: } else if (outputState == 1) {
353: outputState = 4;
354: doOutput = 1;
355: } else if (outputState == 2) {
356: doOutput = 0;
357: } else if (outputState == 3) {
358: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Tried to output CELL_DATA again after intervening POINT_DATA");
359: } else if (outputState == 4) {
360: doOutput = 0;
361: }
362: if (doOutput) {
363: PetscViewerASCIIPrintf(viewer, "CELL_DATA %d\n", n);
364: }
365: }
366: PetscObjectComposedDataSetInt((PetscObject) viewer, stateId, outputState);
367: if (name) {
368: if (bs == 3) {
369: PetscViewerASCIIPrintf(viewer, "VECTORS %s double\n", name);
370: } else {
371: PetscViewerASCIIPrintf(viewer, "SCALARS %s double %d\n", name, bs);
372: }
373: } else {
374: PetscViewerASCIIPrintf(viewer, "SCALARS scalars double %d\n", bs);
375: }
376: if (bs != 3) {
377: PetscViewerASCIIPrintf(viewer, "LOOKUP_TABLE default\n");
378: }
379: for (i=0; i<n/bs; i++) {
380: for (b=0; b<bs; b++) {
381: if (b > 0) {
382: PetscViewerASCIIPrintf(viewer," ");
383: }
384: #if !defined(PETSC_USE_COMPLEX)
385: PetscViewerASCIIPrintf(viewer,"%G",xv[i*bs+b]);
386: #endif
387: }
388: PetscViewerASCIIPrintf(viewer,"\n");
389: }
390: } else if (format == PETSC_VIEWER_ASCII_VTK_COORDS) {
391: PetscInt bs, b;
393: VecGetBlockSize(xin, &bs);
394: if ((bs < 1) || (bs > 3)) {
395: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %d", bs);
396: }
397: for (i=0; i<n/bs; i++) {
398: for (b=0; b<bs; b++) {
399: if (b > 0) {
400: PetscViewerASCIIPrintf(viewer," ");
401: }
402: #if !defined(PETSC_USE_COMPLEX)
403: PetscViewerASCIIPrintf(viewer,"%G",xv[i*bs+b]);
404: #endif
405: }
406: for (b=bs; b<3; b++) {
407: PetscViewerASCIIPrintf(viewer," 0.0");
408: }
409: PetscViewerASCIIPrintf(viewer,"\n");
410: }
411: } else if (format == PETSC_VIEWER_ASCII_PCICE) {
412: PetscInt bs, b;
414: VecGetBlockSize(xin, &bs);
415: if ((bs < 1) || (bs > 3)) {
416: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "PCICE can only handle up to 3D objects, but vector dimension is %d", bs);
417: }
418: PetscViewerASCIIPrintf(viewer,"%D\n", xin->map->N/bs);
419: for (i=0; i<n/bs; i++) {
420: PetscViewerASCIIPrintf(viewer,"%7D ", i+1);
421: for (b=0; b<bs; b++) {
422: if (b > 0) {
423: PetscViewerASCIIPrintf(viewer," ");
424: }
425: #if !defined(PETSC_USE_COMPLEX)
426: PetscViewerASCIIPrintf(viewer,"% 12.5E",xv[i*bs+b]);
427: #endif
428: }
429: PetscViewerASCIIPrintf(viewer,"\n");
430: }
431: } else {
432: PetscObjectPrintClassNamePrefixType((PetscObject)xin,viewer,"Vector Object");
433: for (i=0; i<n; i++) {
434: if (format == PETSC_VIEWER_ASCII_INDEX) {
435: PetscViewerASCIIPrintf(viewer,"%D: ",i);
436: }
437: #if defined(PETSC_USE_COMPLEX)
438: if (PetscImaginaryPart(xv[i]) > 0.0) {
439: PetscViewerASCIIPrintf(viewer,"%G + %G i\n",PetscRealPart(xv[i]),PetscImaginaryPart(xv[i]));
440: } else if (PetscImaginaryPart(xv[i]) < 0.0) {
441: PetscViewerASCIIPrintf(viewer,"%G - %G i\n",PetscRealPart(xv[i]),-PetscImaginaryPart(xv[i]));
442: } else {
443: PetscViewerASCIIPrintf(viewer,"%G\n",PetscRealPart(xv[i]));
444: }
445: #else
446: PetscViewerASCIIPrintf(viewer,"%G\n",(double) xv[i]);
447: #endif
448: }
449: }
450: PetscViewerFlush(viewer);
451: VecRestoreArrayRead(xin,&xv);
452: return(0);
453: }
457: static PetscErrorCode VecView_Seq_Draw_LG(Vec xin,PetscViewer v)
458: {
459: PetscErrorCode ierr;
460: PetscInt i,n = xin->map->n;
461: PetscDraw win;
462: PetscReal *xx;
463: PetscDrawLG lg;
464: const PetscScalar *xv;
467: PetscViewerDrawGetDrawLG(v,0,&lg);
468: PetscDrawLGGetDraw(lg,&win);
469: PetscDrawCheckResizedWindow(win);
470: PetscDrawLGReset(lg);
471: PetscMalloc((n+1)*sizeof(PetscReal),&xx);
472: for (i=0; i<n; i++) {
473: xx[i] = (PetscReal) i;
474: }
475: VecGetArrayRead(xin,&xv);
476: #if !defined(PETSC_USE_COMPLEX)
477: PetscDrawLGAddPoints(lg,n,&xx,(PetscReal**)&xv);
478: #else
479: {
480: PetscReal *yy;
481: PetscMalloc((n+1)*sizeof(PetscReal),&yy);
482: for (i=0; i<n; i++) {
483: yy[i] = PetscRealPart(xv[i]);
484: }
485: PetscDrawLGAddPoints(lg,n,&xx,&yy);
486: PetscFree(yy);
487: }
488: #endif
489: VecRestoreArrayRead(xin,&xv);
490: PetscFree(xx);
491: PetscDrawLGDraw(lg);
492: PetscDrawSynchronizedFlush(win);
493: return(0);
494: }
498: static PetscErrorCode VecView_Seq_Draw(Vec xin,PetscViewer v)
499: {
500: PetscErrorCode ierr;
501: PetscDraw draw;
502: PetscBool isnull;
503: PetscViewerFormat format;
506: PetscViewerDrawGetDraw(v,0,&draw);
507: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
508:
509: PetscViewerGetFormat(v,&format);
510: /*
511: Currently it only supports drawing to a line graph */
512: if (format != PETSC_VIEWER_DRAW_LG) {
513: PetscViewerPushFormat(v,PETSC_VIEWER_DRAW_LG);
514: }
515: VecView_Seq_Draw_LG(xin,v);
516: if (format != PETSC_VIEWER_DRAW_LG) {
517: PetscViewerPopFormat(v);
518: }
520: return(0);
521: }
525: static PetscErrorCode VecView_Seq_Binary(Vec xin,PetscViewer viewer)
526: {
527: PetscErrorCode ierr;
528: int fdes;
529: PetscInt n = xin->map->n,classid=VEC_FILE_CLASSID;
530: FILE *file;
531: const PetscScalar *xv;
532: #if defined(PETSC_HAVE_MPIIO)
533: PetscBool isMPIIO;
534: #endif
535: PetscBool skipHeader;
539: /* Write vector header */
540: PetscViewerBinaryGetSkipHeader(viewer,&skipHeader);
541: if (!skipHeader) {
542: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
543: PetscViewerBinaryWrite(viewer,&n,1,PETSC_INT,PETSC_FALSE);
544: }
546: /* Write vector contents */
547: #if defined(PETSC_HAVE_MPIIO)
548: PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);
549: if (!isMPIIO) {
550: #endif
551: PetscViewerBinaryGetDescriptor(viewer,&fdes);
552: VecGetArrayRead(xin,&xv);
553: PetscBinaryWrite(fdes,(void*)xv,n,PETSC_SCALAR,PETSC_FALSE);
554: VecRestoreArrayRead(xin,&xv);
555: #if defined(PETSC_HAVE_MPIIO)
556: } else {
557: MPI_Offset off;
558: MPI_File mfdes;
559: PetscMPIInt gsizes[1],lsizes[1],lstarts[1];
560: MPI_Datatype view;
562: gsizes[0] = PetscMPIIntCast(n);
563: lsizes[0] = PetscMPIIntCast(n);
564: lstarts[0] = 0;
565: MPI_Type_create_subarray(1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);
566: MPI_Type_commit(&view);
568: PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);
569: PetscViewerBinaryGetMPIIOOffset(viewer,&off);
570: VecGetArrayRead(xin,&xv);
571: MPIU_File_write_all(mfdes,(void*)xv,lsizes[0],MPIU_SCALAR,MPI_STATUS_IGNORE);
572: VecRestoreArrayRead(xin,&xv);
573: PetscViewerBinaryAddMPIIOOffset(viewer,n*sizeof(PetscScalar));
574: MPI_Type_free(&view);
575: }
576: #endif
578: PetscViewerBinaryGetInfoPointer(viewer,&file);
579: if (file) {
580: if (((PetscObject)xin)->prefix) {
581: PetscFPrintf(PETSC_COMM_SELF,file,"-%s_vecload_block_size %D\n",((PetscObject)xin)->prefix,xin->map->bs);
582: } else {
583: PetscFPrintf(PETSC_COMM_SELF,file,"-vecload_block_size %D\n",xin->map->bs);
584: }
585: }
586: return(0);
587: }
589: #if defined(PETSC_HAVE_MATLAB_ENGINE)
590: #include <mat.h> /* MATLAB include file */
594: static PetscErrorCode VecView_Seq_Matlab(Vec vec,PetscViewer viewer)
595: {
596: PetscErrorCode ierr;
597: PetscInt n;
598: const PetscScalar *array;
599:
601: VecGetLocalSize(vec,&n);
602: PetscObjectName((PetscObject)vec);
603: VecGetArrayRead(vec,&array);
604: PetscViewerMatlabPutArray(viewer,n,1,array,((PetscObject)vec)->name);
605: VecRestoreArrayRead(vec,&array);
606: return(0);
607: }
609: #endif
613: static PetscErrorCode VecView_Seq(Vec xin,PetscViewer viewer)
614: {
616: PetscBool isdraw,iascii,issocket,isbinary;
617: #if defined(PETSC_HAVE_MATHEMATICA)
618: PetscBool ismathematica;
619: #endif
620: #if defined(PETSC_HAVE_MATLAB_ENGINE)
621: PetscBool ismatlab;
622: #endif
623: #if defined(PETSC_HAVE_HDF5)
624: PetscBool ishdf5;
625: #endif
628: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
629: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
630: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
631: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
632: #if defined(PETSC_HAVE_MATHEMATICA)
633: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERMATHEMATICA,&ismathematica);
634: #endif
635: #if defined(PETSC_HAVE_HDF5)
636: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
637: #endif
638: #if defined(PETSC_HAVE_MATLAB_ENGINE)
639: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);
640: #endif
642: if (isdraw){
643: VecView_Seq_Draw(xin,viewer);
644: } else if (iascii){
645: VecView_Seq_ASCII(xin,viewer);
646: } else if (isbinary) {
647: VecView_Seq_Binary(xin,viewer);
648: #if defined(PETSC_HAVE_MATHEMATICA)
649: } else if (ismathematica) {
650: PetscViewerMathematicaPutVector(viewer,xin);
651: #endif
652: #if defined(PETSC_HAVE_HDF5)
653: } else if (ishdf5) {
654: VecView_MPI_HDF5(xin,viewer); /* Reusing VecView_MPI_HDF5 ... don't want code duplication*/
655: #endif
656: #if defined(PETSC_HAVE_MATLAB_ENGINE)
657: } else if (ismatlab) {
658: VecView_Seq_Matlab(xin,viewer);
659: #endif
660: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by this vector object",((PetscObject)viewer)->type_name);
661: return(0);
662: }
666: PetscErrorCode VecGetValues_Seq(Vec xin,PetscInt ni,const PetscInt ix[],PetscScalar y[])
667: {
668: const PetscScalar *xx;
669: PetscInt i;
670: PetscErrorCode ierr;
673: VecGetArrayRead(xin,&xx);
674: for (i=0; i<ni; i++) {
675: if (xin->stash.ignorenegidx && ix[i] < 0) continue;
676: #if defined(PETSC_USE_DEBUG)
677: if (ix[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
678: if (ix[i] >= xin->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D to large maximum allowed %D",ix[i],xin->map->n);
679: #endif
680: y[i] = xx[ix[i]];
681: }
682: VecRestoreArrayRead(xin,&xx);
683: return(0);
684: }
688: PetscErrorCode VecSetValues_Seq(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode m)
689: {
690: PetscScalar *xx;
691: PetscInt i;
695: VecGetArray(xin,&xx);
696: if (m == INSERT_VALUES) {
697: for (i=0; i<ni; i++) {
698: if (xin->stash.ignorenegidx && ix[i] < 0) continue;
699: #if defined(PETSC_USE_DEBUG)
700: if (ix[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
701: if (ix[i] >= xin->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",ix[i],xin->map->n);
702: #endif
703: xx[ix[i]] = y[i];
704: }
705: } else {
706: for (i=0; i<ni; i++) {
707: if (xin->stash.ignorenegidx && ix[i] < 0) continue;
708: #if defined(PETSC_USE_DEBUG)
709: if (ix[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
710: if (ix[i] >= xin->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",ix[i],xin->map->n);
711: #endif
712: xx[ix[i]] += y[i];
713: }
714: }
715: VecRestoreArray(xin,&xx);
716: return(0);
717: }
721: PetscErrorCode VecSetValuesBlocked_Seq(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar yin[],InsertMode m)
722: {
723: PetscScalar *xx,*y = (PetscScalar*)yin;
724: PetscInt i,bs = xin->map->bs,start,j;
727: /*
728: For optimization could treat bs = 2, 3, 4, 5 as special cases with loop unrolling
729: */
731: VecGetArray(xin,&xx);
732: if (m == INSERT_VALUES) {
733: for (i=0; i<ni; i++) {
734: start = bs*ix[i];
735: if (start < 0) continue;
736: #if defined(PETSC_USE_DEBUG)
737: if (start >= xin->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",start,xin->map->n);
738: #endif
739: for (j=0; j<bs; j++) {
740: xx[start+j] = y[j];
741: }
742: y += bs;
743: }
744: } else {
745: for (i=0; i<ni; i++) {
746: start = bs*ix[i];
747: if (start < 0) continue;
748: #if defined(PETSC_USE_DEBUG)
749: if (start >= xin->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",start,xin->map->n);
750: #endif
751: for (j=0; j<bs; j++) {
752: xx[start+j] += y[j];
753: }
754: y += bs;
755: }
756: }
757: VecRestoreArray(xin,&xx);
758: return(0);
759: }
763: PetscErrorCode VecDestroy_Seq(Vec v)
764: {
765: Vec_Seq *vs = (Vec_Seq*)v->data;
769: PetscObjectDepublish(v);
771: #if defined(PETSC_USE_LOG)
772: PetscLogObjectState((PetscObject)v,"Length=%D",v->map->n);
773: #endif
774: PetscFree(vs->array_allocated);
775: PetscFree(vs);
776: return(0);
777: }
781: static PetscErrorCode VecSetOption_Seq(Vec v,VecOption op,PetscBool flag)
782: {
784: if (op == VEC_IGNORE_NEGATIVE_INDICES) {
785: v->stash.ignorenegidx = flag;
786: }
787: return(0);
788: }
792: static PetscErrorCode VecDuplicate_Seq(Vec win,Vec *V)
793: {
797: VecCreate(((PetscObject)win)->comm,V);
798: PetscObjectSetPrecision((PetscObject)*V,((PetscObject)win)->precision);
799: VecSetSizes(*V,win->map->n,win->map->n);
800: VecSetType(*V,((PetscObject)win)->type_name);
801: PetscLayoutReference(win->map,&(*V)->map);
802: PetscOListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*V))->olist);
803: PetscFListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*V))->qlist);
805: (*V)->stash.ignorenegidx = win->stash.ignorenegidx;
807: return(0);
808: }
810: static struct _VecOps DvOps = {VecDuplicate_Seq, /* 1 */
811: VecDuplicateVecs_Default,
812: VecDestroyVecs_Default,
813: VecDot_Seq,
814: VecMDot_Seq,
815: VecNorm_Seq,
816: VecTDot_Seq,
817: VecMTDot_Seq,
818: VecScale_Seq,
819: VecCopy_Seq, /* 10 */
820: VecSet_Seq,
821: VecSwap_Seq,
822: VecAXPY_Seq,
823: VecAXPBY_Seq,
824: VecMAXPY_Seq,
825: VecAYPX_Seq,
826: VecWAXPY_Seq,
827: VecAXPBYPCZ_Seq,
828: VecPointwiseMult_Seq,
829: VecPointwiseDivide_Seq,
830: VecSetValues_Seq, /* 20 */
831: 0,0,
832: 0,
833: VecGetSize_Seq,
834: VecGetSize_Seq,
835: 0,
836: VecMax_Seq,
837: VecMin_Seq,
838: VecSetRandom_Seq,
839: VecSetOption_Seq, /* 30 */
840: VecSetValuesBlocked_Seq,
841: VecDestroy_Seq,
842: VecView_Seq,
843: VecPlaceArray_Seq,
844: VecReplaceArray_Seq,
845: VecDot_Seq,
846: VecTDot_Seq,
847: VecNorm_Seq,
848: VecMDot_Seq,
849: VecMTDot_Seq, /* 40 */
850: VecLoad_Default,
851: VecReciprocal_Default,
852: VecConjugate_Seq,
853: 0,
854: 0,
855: VecResetArray_Seq,
856: 0,
857: VecMaxPointwiseDivide_Seq,
858: VecPointwiseMax_Seq,
859: VecPointwiseMaxAbs_Seq,
860: VecPointwiseMin_Seq,
861: VecGetValues_Seq,
862: 0,
863: 0,
864: 0,
865: 0,
866: 0,
867: 0,
868: VecStrideGather_Default,
869: VecStrideScatter_Default
870: };
873: /*
874: This is called by VecCreate_Seq() (i.e. VecCreateSeq()) and VecCreateSeqWithArray()
875: */
878: PetscErrorCode VecCreate_Seq_Private(Vec v,const PetscScalar array[])
879: {
880: Vec_Seq *s;
884: PetscNewLog(v,Vec_Seq,&s);
885: PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));
886: v->data = (void*)s;
887: v->petscnative = PETSC_TRUE;
888: s->array = (PetscScalar *)array;
889: s->array_allocated = 0;
891: if (v->map->bs == -1) v->map->bs = 1;
892: PetscLayoutSetUp(v->map);
893: PetscObjectChangeTypeName((PetscObject)v,VECSEQ);
894: #if defined(PETSC_HAVE_MATLAB_ENGINE)
895: PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEnginePut_C","VecMatlabEnginePut_Default",VecMatlabEnginePut_Default);
896: PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEngineGet_C","VecMatlabEngineGet_Default",VecMatlabEngineGet_Default);
897: #endif
898: #if defined(PETSC_USE_MIXED_PRECISION)
899: ((PetscObject)v)->precision = (PetscPrecision)sizeof(PetscScalar);
900: #endif
901: return(0);
902: }
906: /*@C
907: VecCreateSeqWithArray - Creates a standard,sequential array-style vector,
908: where the user provides the array space to store the vector values.
910: Collective on MPI_Comm
912: Input Parameter:
913: + comm - the communicator, should be PETSC_COMM_SELF
914: . n - the vector length
915: - array - memory where the vector elements are to be stored.
917: Output Parameter:
918: . V - the vector
920: Notes:
921: Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
922: same type as an existing vector.
924: If the user-provided array is PETSC_NULL, then VecPlaceArray() can be used
925: at a later stage to SET the array for storing the vector values.
927: PETSc does NOT free the array when the vector is destroyed via VecDestroy().
928: The user should not free the array until the vector is destroyed.
930: Level: intermediate
932: Concepts: vectors^creating with array
934: .seealso: VecCreateMPIWithArray(), VecCreate(), VecDuplicate(), VecDuplicateVecs(),
935: VecCreateGhost(), VecCreateSeq(), VecPlaceArray()
936: @*/
937: PetscErrorCode VecCreateSeqWithArray(MPI_Comm comm,PetscInt n,const PetscScalar array[],Vec *V)
938: {
940: PetscMPIInt size;
943: VecCreate(comm,V);
944: VecSetSizes(*V,n,n);
945: MPI_Comm_size(comm,&size);
946: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot create VECSEQ on more than one process");
947: VecCreate_Seq_Private(*V,array);
948: return(0);
949: }