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 src/inline/dot.h
9: #include petscblaslapack.h
10: #if defined(PETSC_HAVE_PNETCDF)
12: #include "pnetcdf.h"
14: #endif
18: PetscErrorCode VecNorm_Seq(Vec xin,NormType type,PetscReal* z)
19: {
20: PetscScalar *xx;
22: PetscInt n = xin->map.n;
23: PetscBLASInt bn = (PetscBLASInt)n,one = 1;
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(2*n-1);
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(n-1);
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: PetscInt bs, b;
129: VecGetBlockSize(xin, &bs);
130: if ((bs < 1) || (bs > 3)) {
131: SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %d", bs);
132: }
133: if (format == PETSC_VIEWER_ASCII_VTK) {
134: PetscViewerASCIIPrintf(viewer, "POINT_DATA %d\n", n);
135: } else {
136: PetscViewerASCIIPrintf(viewer, "CELL_DATA %d\n", n);
137: }
138: PetscViewerASCIIPrintf(viewer, "SCALARS scalars double %d\n", bs);
139: PetscViewerASCIIPrintf(viewer, "LOOKUP_TABLE default\n");
140: for (i=0; i<n/bs; i++) {
141: for (b=0; b<bs; b++) {
142: if (b > 0) {
143: PetscViewerASCIIPrintf(viewer," ");
144: }
145: #if !defined(PETSC_USE_COMPLEX)
146: PetscViewerASCIIPrintf(viewer,"%G",x->array[i*bs+b]);
147: #endif
148: }
149: PetscViewerASCIIPrintf(viewer,"\n");
150: }
151: } else if (format == PETSC_VIEWER_ASCII_VTK_COORDS) {
152: PetscInt bs, b;
154: VecGetBlockSize(xin, &bs);
155: if ((bs < 1) || (bs > 3)) {
156: SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %d", bs);
157: }
158: for (i=0; i<n/bs; i++) {
159: for (b=0; b<bs; b++) {
160: if (b > 0) {
161: PetscViewerASCIIPrintf(viewer," ");
162: }
163: #if !defined(PETSC_USE_COMPLEX)
164: PetscViewerASCIIPrintf(viewer,"%G",x->array[i*bs+b]);
165: #endif
166: }
167: for (b=bs; b<3; b++) {
168: PetscViewerASCIIPrintf(viewer," 0.0");
169: }
170: PetscViewerASCIIPrintf(viewer,"\n");
171: }
172: } else if (format == PETSC_VIEWER_ASCII_PCICE) {
173: PetscInt bs, b;
175: VecGetBlockSize(xin, &bs);
176: if ((bs < 1) || (bs > 3)) {
177: SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "PCICE can only handle up to 3D objects, but vector dimension is %d", bs);
178: }
179: PetscViewerASCIIPrintf(viewer,"%D\n", xin->map.N/bs);
180: for (i=0; i<n/bs; i++) {
181: PetscViewerASCIIPrintf(viewer,"%7D ", i+1);
182: for (b=0; b<bs; b++) {
183: if (b > 0) {
184: PetscViewerASCIIPrintf(viewer," ");
185: }
186: #if !defined(PETSC_USE_COMPLEX)
187: PetscViewerASCIIPrintf(viewer,"% 12.5E",x->array[i*bs+b]);
188: #endif
189: }
190: PetscViewerASCIIPrintf(viewer,"\n");
191: }
192: } else if (format == PETSC_VIEWER_ASCII_PYLITH) {
193: PetscInt bs, b;
195: VecGetBlockSize(xin, &bs);
196: if (bs != 3) {
197: SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "PyLith can only handle 3D objects, but vector dimension is %d", bs);
198: }
199: PetscViewerASCIIPrintf(viewer,"#\n");
200: PetscViewerASCIIPrintf(viewer,"# Node X-coord Y-coord Z-coord\n");
201: PetscViewerASCIIPrintf(viewer,"#\n");
202: for (i=0; i<n/bs; i++) {
203: PetscViewerASCIIPrintf(viewer,"%7D ", i+1);
204: for (b=0; b<bs; b++) {
205: if (b > 0) {
206: PetscViewerASCIIPrintf(viewer," ");
207: }
208: #if !defined(PETSC_USE_COMPLEX)
209: PetscViewerASCIIPrintf(viewer,"% 16.8E",x->array[i*bs+b]);
210: #endif
211: }
212: PetscViewerASCIIPrintf(viewer,"\n");
213: }
214: } else {
215: for (i=0; i<n; i++) {
216: if (format == PETSC_VIEWER_ASCII_INDEX) {
217: PetscViewerASCIIPrintf(viewer,"%D: ",i);
218: }
219: #if defined(PETSC_USE_COMPLEX)
220: if (PetscImaginaryPart(x->array[i]) > 0.0) {
221: PetscViewerASCIIPrintf(viewer,"%G + %G i\n",PetscRealPart(x->array[i]),PetscImaginaryPart(x->array[i]));
222: } else if (PetscImaginaryPart(x->array[i]) < 0.0) {
223: PetscViewerASCIIPrintf(viewer,"%G - %G i\n",PetscRealPart(x->array[i]),-PetscImaginaryPart(x->array[i]));
224: } else {
225: PetscViewerASCIIPrintf(viewer,"%G\n",PetscRealPart(x->array[i]));
226: }
227: #else
228: PetscViewerASCIIPrintf(viewer,"%G\n",x->array[i]);
229: #endif
230: }
231: }
232: PetscViewerFlush(viewer);
233: return(0);
234: }
238: static PetscErrorCode VecView_Seq_Draw_LG(Vec xin,PetscViewer v)
239: {
240: Vec_Seq *x = (Vec_Seq *)xin->data;
242: PetscInt i,n = xin->map.n;
243: PetscDraw win;
244: PetscReal *xx;
245: PetscDrawLG lg;
248: PetscViewerDrawGetDrawLG(v,0,&lg);
249: PetscDrawLGGetDraw(lg,&win);
250: PetscDrawCheckResizedWindow(win);
251: PetscDrawLGReset(lg);
252: PetscMalloc((n+1)*sizeof(PetscReal),&xx);
253: for (i=0; i<n; i++) {
254: xx[i] = (PetscReal) i;
255: }
256: #if !defined(PETSC_USE_COMPLEX)
257: PetscDrawLGAddPoints(lg,n,&xx,&x->array);
258: #else
259: {
260: PetscReal *yy;
261: PetscMalloc((n+1)*sizeof(PetscReal),&yy);
262: for (i=0; i<n; i++) {
263: yy[i] = PetscRealPart(x->array[i]);
264: }
265: PetscDrawLGAddPoints(lg,n,&xx,&yy);
266: PetscFree(yy);
267: }
268: #endif
269: PetscFree(xx);
270: PetscDrawLGDraw(lg);
271: PetscDrawSynchronizedFlush(win);
272: return(0);
273: }
277: static PetscErrorCode VecView_Seq_Draw(Vec xin,PetscViewer v)
278: {
279: PetscErrorCode ierr;
280: PetscDraw draw;
281: PetscTruth isnull;
282: PetscViewerFormat format;
285: PetscViewerDrawGetDraw(v,0,&draw);
286: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
287:
288: PetscViewerGetFormat(v,&format);
289: /*
290: Currently it only supports drawing to a line graph */
291: if (format != PETSC_VIEWER_DRAW_LG) {
292: PetscViewerPushFormat(v,PETSC_VIEWER_DRAW_LG);
293: }
294: VecView_Seq_Draw_LG(xin,v);
295: if (format != PETSC_VIEWER_DRAW_LG) {
296: PetscViewerPopFormat(v);
297: }
299: return(0);
300: }
304: static PetscErrorCode VecView_Seq_Binary(Vec xin,PetscViewer viewer)
305: {
306: Vec_Seq *x = (Vec_Seq *)xin->data;
308: int fdes;
309: PetscInt n = xin->map.n,cookie=VEC_FILE_COOKIE;
310: FILE *file;
313: PetscViewerBinaryGetDescriptor(viewer,&fdes);
314: /* Write vector header */
315: PetscBinaryWrite(fdes,&cookie,1,PETSC_INT,PETSC_FALSE);
316: PetscBinaryWrite(fdes,&n,1,PETSC_INT,PETSC_FALSE);
318: /* Write vector contents */
319: PetscBinaryWrite(fdes,x->array,n,PETSC_SCALAR,PETSC_FALSE);
321: PetscViewerBinaryGetInfoPointer(viewer,&file);
322: if (file && xin->map.bs > 1) {
323: if (xin->prefix) {
324: PetscFPrintf(PETSC_COMM_SELF,file,"-%s_vecload_block_size %D\n",xin->prefix,xin->map.bs);
325: } else {
326: PetscFPrintf(PETSC_COMM_SELF,file,"-vecload_block_size %D\n",xin->map.bs);
327: }
328: }
329: return(0);
330: }
332: #if defined(PETSC_HAVE_PNETCDF)
335: PetscErrorCode VecView_Seq_Netcdf(Vec xin,PetscViewer v)
336: {
338: int n = xin->map.n,ncid,xdim,xdim_num=1,xin_id,xstart=0;
339: PetscScalar *xarray;
342: #if !defined(PETSC_USE_COMPLEX)
343: VecGetArray(xin,&xarray);
344: PetscViewerNetcdfGetID(v,&ncid);
345: if (ncid < 0) SETERRQ(PETSC_ERR_ORDER,"First call PetscViewerNetcdfOpen to create NetCDF dataset");
346: /* define dimensions */
347: ncmpi_def_dim(ncid,"PETSc_Vector_Global_Size",n,&xdim);
348: /* define variables */
349: ncmpi_def_var(ncid,"PETSc_Vector_Seq",NC_DOUBLE,xdim_num,&xdim,&xin_id);
350: /* leave define mode */
351: ncmpi_enddef(ncid);
352: /* store the vector */
353: VecGetOwnershipRange(xin,&xstart,PETSC_NULL);
354: ncmpi_put_vara_double_all(ncid,xin_id,(const MPI_Offset*)&xstart,(const MPI_Offset*)&n,xarray);
355: #else
356: PetscPrintf(PETSC_COMM_WORLD,"NetCDF viewer not supported for complex numbers\n");
357: #endif
358: return(0);
359: }
360: #endif
362: #if defined(PETSC_HAVE_MATLAB)
363: #include "mat.h" /* Matlab include file */
367: PetscErrorCode VecView_Seq_Matlab(Vec vec,PetscViewer viewer)
368: {
370: PetscInt n;
371: PetscScalar *array;
372:
374: VecGetLocalSize(vec,&n);
375: PetscObjectName((PetscObject)vec);
376: VecGetArray(vec,&array);
377: PetscViewerMatlabPutArray(viewer,n,1,array,vec->name);
378: VecRestoreArray(vec,&array);
379: return(0);
380: }
382: #endif
386: PetscErrorCode VecView_Seq(Vec xin,PetscViewer viewer)
387: {
389: PetscTruth isdraw,iascii,issocket,isbinary;
390: #if defined(PETSC_HAVE_MATHEMATICA)
391: PetscTruth ismathematica;
392: #endif
393: #if defined(PETSC_HAVE_PNETCDF)
394: PetscTruth isnetcdf;
395: #endif
396: #if defined(PETSC_HAVE_MATLAB)
397: PetscTruth ismatlab;
398: #endif
401: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
402: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
403: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
404: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
405: #if defined(PETSC_HAVE_MATHEMATICA)
406: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_MATHEMATICA,&ismathematica);
407: #endif
408: #if defined(PETSC_HAVE_PNETCDF)
409: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_NETCDF,&isnetcdf);
410: #endif
411: #if defined(PETSC_HAVE_MATLAB)
412: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_MATLAB,&ismatlab);
413: #endif
415: if (isdraw){
416: VecView_Seq_Draw(xin,viewer);
417: } else if (iascii){
418: VecView_Seq_File(xin,viewer);
419: #if defined(PETSC_USE_SOCKET_VIEWER)
420: } else if (issocket) {
421: Vec_Seq *x = (Vec_Seq *)xin->data;
422: PetscViewerSocketPutScalar(viewer,xin->map.n,1,x->array);
423: #endif
424: } else if (isbinary) {
425: VecView_Seq_Binary(xin,viewer);
426: #if defined(PETSC_HAVE_MATHEMATICA)
427: } else if (ismathematica) {
428: PetscViewerMathematicaPutVector(viewer,xin);
429: #endif
430: #if defined(PETSC_HAVE_PNETCDF)
431: } else if (isnetcdf) {
432: VecView_Seq_Netcdf(xin,viewer);
433: #endif
434: #if defined(PETSC_HAVE_MATLAB)
435: } else if (ismatlab) {
436: VecView_Seq_Matlab(xin,viewer);
437: #endif
438: } else {
439: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by this vector object",((PetscObject)viewer)->type_name);
440: }
441: return(0);
442: }
446: PetscErrorCode VecGetValues_Seq(Vec xin,PetscInt ni,const PetscInt ix[],PetscScalar y[])
447: {
448: Vec_Seq *x = (Vec_Seq *)xin->data;
449: PetscScalar *xx = x->array;
450: PetscInt i;
453: for (i=0; i<ni; i++) {
454: #if defined(PETSC_USE_DEBUG)
455: if (ix[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
456: 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);
457: #endif
458: y[i] = xx[ix[i]];
459: }
460: return(0);
461: }
465: PetscErrorCode VecSetValues_Seq(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode m)
466: {
467: Vec_Seq *x = (Vec_Seq *)xin->data;
468: PetscScalar *xx = x->array;
469: PetscInt i;
472: if (m == INSERT_VALUES) {
473: for (i=0; i<ni; i++) {
474: if (ix[i] < 0) continue;
475: #if defined(PETSC_USE_DEBUG)
476: if (ix[i] >= xin->map.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",ix[i],xin->map.n);
477: #endif
478: xx[ix[i]] = y[i];
479: }
480: } else {
481: for (i=0; i<ni; i++) {
482: if (ix[i] < 0) continue;
483: #if defined(PETSC_USE_DEBUG)
484: if (ix[i] >= xin->map.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",ix[i],xin->map.n);
485: #endif
486: xx[ix[i]] += y[i];
487: }
488: }
489: return(0);
490: }
494: PetscErrorCode VecSetValuesBlocked_Seq(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar yin[],InsertMode m)
495: {
496: Vec_Seq *x = (Vec_Seq *)xin->data;
497: PetscScalar *xx = x->array,*y = (PetscScalar*)yin;
498: PetscInt i,bs = xin->map.bs,start,j;
500: /*
501: For optimization could treat bs = 2, 3, 4, 5 as special cases with loop unrolling
502: */
504: if (m == INSERT_VALUES) {
505: for (i=0; i<ni; i++) {
506: start = bs*ix[i];
507: if (start < 0) continue;
508: #if defined(PETSC_USE_DEBUG)
509: if (start >= xin->map.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",start,xin->map.n);
510: #endif
511: for (j=0; j<bs; j++) {
512: xx[start+j] = y[j];
513: }
514: y += bs;
515: }
516: } else {
517: for (i=0; i<ni; i++) {
518: start = bs*ix[i];
519: if (start < 0) continue;
520: #if defined(PETSC_USE_DEBUG)
521: if (start >= xin->map.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",start,xin->map.n);
522: #endif
523: for (j=0; j<bs; j++) {
524: xx[start+j] += y[j];
525: }
526: y += bs;
527: }
528: }
529: return(0);
530: }
535: PetscErrorCode VecDestroy_Seq(Vec v)
536: {
537: Vec_Seq *vs = (Vec_Seq*)v->data;
542: /* if memory was published with AMS then destroy it */
543: PetscObjectDepublish(v);
545: #if defined(PETSC_USE_LOG)
546: PetscLogObjectState((PetscObject)v,"Length=%D",v->map.n);
547: #endif
548: PetscFree(vs->array_allocated);
549: PetscFree(vs);
551: return(0);
552: }
556: static PetscErrorCode VecPublish_Seq(PetscObject obj)
557: {
559: return(0);
560: }
562: EXTERN PetscErrorCode VecLoad_Binary(PetscViewer, VecType, Vec*);
564: static struct _VecOps DvOps = {VecDuplicate_Seq,
565: VecDuplicateVecs_Default,
566: VecDestroyVecs_Default,
567: VecDot_Seq,
568: VecMDot_Seq,
569: VecNorm_Seq,
570: VecTDot_Seq,
571: VecMTDot_Seq,
572: VecScale_Seq,
573: VecCopy_Seq,
574: VecSet_Seq,
575: VecSwap_Seq,
576: VecAXPY_Seq,
577: VecAXPBY_Seq,
578: VecMAXPY_Seq,
579: VecAYPX_Seq,
580: VecWAXPY_Seq,
581: VecPointwiseMult_Seq,
582: VecPointwiseDivide_Seq,
583: VecSetValues_Seq,0,0,
584: VecGetArray_Seq,
585: VecGetSize_Seq,
586: VecGetSize_Seq,
587: VecRestoreArray_Seq,
588: VecMax_Seq,
589: VecMin_Seq,
590: VecSetRandom_Seq,0,
591: VecSetValuesBlocked_Seq,
592: VecDestroy_Seq,
593: VecView_Seq,
594: VecPlaceArray_Seq,
595: VecReplaceArray_Seq,
596: VecDot_Seq,
597: VecTDot_Seq,
598: VecNorm_Seq,
599: VecMDot_Seq,
600: VecMTDot_Seq,
601: VecLoadIntoVector_Default,
602: VecReciprocal_Default,
603: 0, /* VecViewNative */
604: VecConjugate_Seq,
605: 0,
606: 0,
607: VecResetArray_Seq,
608: 0,
609: VecMaxPointwiseDivide_Seq,
610: VecLoad_Binary,
611: VecPointwiseMax_Seq,
612: VecPointwiseMaxAbs_Seq,
613: VecPointwiseMin_Seq,
614: VecGetValues_Seq
615: };
618: /*
619: This is called by VecCreate_Seq() (i.e. VecCreateSeq()) and VecCreateSeqWithArray()
620: */
623: static PetscErrorCode VecCreate_Seq_Private(Vec v,const PetscScalar array[])
624: {
625: Vec_Seq *s;
629: PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));
630: PetscNew(Vec_Seq,&s);
631: v->precision = PETSC_SCALAR;
632: v->data = (void*)s;
633: v->bops->publish = VecPublish_Seq;
634: v->petscnative = PETSC_TRUE;
635: s->array = (PetscScalar *)array;
636: s->array_allocated = 0;
638: if (v->map.bs == -1) v->map.bs = 1;
639: PetscMapInitialize(v->comm,&v->map);
640: PetscObjectChangeTypeName((PetscObject)v,VECSEQ);
641: #if defined(PETSC_HAVE_MATLAB)
642: PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEnginePut_C","VecMatlabEnginePut_Default",VecMatlabEnginePut_Default);
643: PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEngineGet_C","VecMatlabEngineGet_Default",VecMatlabEngineGet_Default);
644: #endif
645: PetscPublishAll(v);
646: return(0);
647: }
651: /*@C
652: VecCreateSeqWithArray - Creates a standard,sequential array-style vector,
653: where the user provides the array space to store the vector values.
655: Collective on MPI_Comm
657: Input Parameter:
658: + comm - the communicator, should be PETSC_COMM_SELF
659: . n - the vector length
660: - array - memory where the vector elements are to be stored.
662: Output Parameter:
663: . V - the vector
665: Notes:
666: Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
667: same type as an existing vector.
669: If the user-provided array is PETSC_NULL, then VecPlaceArray() can be used
670: at a later stage to SET the array for storing the vector values.
672: PETSc does NOT free the array when the vector is destroyed via VecDestroy().
673: The user should not free the array until the vector is destroyed.
675: Level: intermediate
677: Concepts: vectors^creating with array
679: .seealso: VecCreateMPIWithArray(), VecCreate(), VecDuplicate(), VecDuplicateVecs(),
680: VecCreateGhost(), VecCreateSeq(), VecPlaceArray()
681: @*/
682: PetscErrorCode PETSCVEC_DLLEXPORT VecCreateSeqWithArray(MPI_Comm comm,PetscInt n,const PetscScalar array[],Vec *V)
683: {
685: PetscMPIInt size;
688: VecCreate(comm,V);
689: VecSetSizes(*V,n,n);
690: MPI_Comm_size(comm,&size);
691: if (size > 1) {
692: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot create VECSEQ on more than one process");
693: }
694: VecCreate_Seq_Private(*V,array);
695: return(0);
696: }
698: /*MC
699: VECSEQ - VECSEQ = "seq" - The basic sequential vector
701: Options Database Keys:
702: . -vec_type seq - sets the vector type to VECSEQ during a call to VecSetFromOptions()
704: Level: beginner
706: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateSeqWithArray(), VECMPI, VecType, VecCreateMPI(), VecCreateSeq()
707: M*/
712: PetscErrorCode PETSCVEC_DLLEXPORT VecCreate_Seq(Vec V)
713: {
714: Vec_Seq *s;
715: PetscScalar *array;
717: PetscInt n = PetscMax(V->map.n,V->map.N);
718: PetscMPIInt size;
721: MPI_Comm_size(V->comm,&size);
722: if (size > 1) {
723: SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot create VECSEQ on more than one process");
724: }
725: PetscMalloc( n*sizeof(PetscScalar),&array);
726: PetscMemzero(array,n*sizeof(PetscScalar));
727: VecCreate_Seq_Private(V,array);
728: s = (Vec_Seq*)V->data;
729: s->array_allocated = array;
730: return(0);
731: }
737: PetscErrorCode VecDuplicate_Seq(Vec win,Vec *V)
738: {
742: VecCreateSeq(win->comm,win->map.n,V);
743: if (win->mapping) {
744: (*V)->mapping = win->mapping;
745: PetscObjectReference((PetscObject)win->mapping);
746: }
747: if (win->bmapping) {
748: (*V)->bmapping = win->bmapping;
749: PetscObjectReference((PetscObject)win->bmapping);
750: }
751: (*V)->map.bs = win->map.bs;
752: PetscOListDuplicate(win->olist,&(*V)->olist);
753: PetscFListDuplicate(win->qlist,&(*V)->qlist);
754: return(0);
755: }