Actual source code: bvec2.c
1: /*$Id: bvec2.c,v 1.202 2001/09/12 03:26:24 bsmith Exp $*/
2: /*
3: Implements the sequential vectors.
4: */
6: #include src/vec/vecimpl.h
7: #include src/vec/impls/dvecimpl.h
8: #include petscblaslapack.h
9: #if defined(PETSC_HAVE_AMS)
10: EXTERN int PetscViewerAMSGetAMSComm(PetscViewer,AMS_Comm *);
11: #endif
13: int VecNorm_Seq(Vec xin,NormType type,PetscReal* z)
14: {
15: PetscScalar *xa;
16: int ierr,one = 1;
19: if (type == NORM_2) {
20: VecGetArrayFast(xin,&xa);
21: /*
22: This is because the Fortran BLAS 1 Norm is very slow!
23: */
24: #if defined(PETSC_HAVE_SLOW_NRM2)
25: {
26: int i;
27: PetscScalar sum=0.0;
28: for (i=0; i<xin->n; i++) {
29: sum += (xa[i])*(PetscConj(xa[i]));
30: }
31: *z = sqrt(PetscRealPart(sum));
32: }
33: #else
34: *z = BLnrm2_(&xin->n,xa,&one);
35: #endif
36: VecRestoreArrayFast(xin,&xa);
37: PetscLogFlops(2*xin->n-1);
38: } else if (type == NORM_INFINITY) {
39: int i,n = xin->n;
40: PetscReal max = 0.0,tmp;
42: VecGetArrayFast(xin,&xa);
43: for (i=0; i<n; i++) {
44: if ((tmp = PetscAbsScalar(*xa)) > max) max = tmp;
45: /* check special case of tmp == NaN */
46: if (tmp != tmp) {max = tmp; break;}
47: xa++;
48: }
49: VecRestoreArrayFast(xin,&xa);
50: *z = max;
51: } else if (type == NORM_1) {
52: VecGetArrayFast(xin,&xa);
53: *z = BLasum_(&xin->n,xa,&one);
54: VecRestoreArrayFast(xin,&xa);
55: PetscLogFlops(xin->n-1);
56: } else if (type == NORM_1_AND_2) {
57: VecNorm_Seq(xin,NORM_1,z);
58: VecNorm_Seq(xin,NORM_2,z+1);
59: }
60: return(0);
61: }
63: #include petscviewer.h
64: #include petscsys.h
66: int VecView_Seq_File(Vec xin,PetscViewer viewer)
67: {
68: Vec_Seq *x = (Vec_Seq *)xin->data;
69: int i,n = xin->n,ierr;
70: char *name;
71: PetscViewerFormat format;
74: PetscViewerGetFormat(viewer,&format);
75: if (format == PETSC_VIEWER_ASCII_MATLAB) {
76: PetscObjectGetName((PetscObject)xin,&name);
77: PetscViewerASCIIPrintf(viewer,"%s = [n",name);
78: for (i=0; i<n; i++) {
79: #if defined(PETSC_USE_COMPLEX)
80: if (PetscImaginaryPart(x->array[i]) > 0.0) {
81: PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ein",PetscRealPart(x->array[i]),PetscImaginaryPart(x->array[i]));
82: } else if (PetscImaginaryPart(x->array[i]) < 0.0) {
83: PetscViewerASCIIPrintf(viewer,"%18.16e - %18.16ein",PetscRealPart(x->array[i]),-PetscImaginaryPart(x->array[i]));
84: } else {
85: PetscViewerASCIIPrintf(viewer,"%18.16en",PetscRealPart(x->array[i]));
86: }
87: #else
88: PetscViewerASCIIPrintf(viewer,"%18.16en",x->array[i]);
89: #endif
90: }
91: PetscViewerASCIIPrintf(viewer,"];n");
92: } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
93: for (i=0; i<n; i++) {
94: #if defined(PETSC_USE_COMPLEX)
95: PetscViewerASCIIPrintf(viewer,"%18.16e %18.16en",PetscRealPart(x->array[i]),PetscImaginaryPart(x->array[i]));
96: #else
97: PetscViewerASCIIPrintf(viewer,"%18.16en",x->array[i]);
98: #endif
99: }
100: } else {
101: for (i=0; i<n; i++) {
102: if (format == PETSC_VIEWER_ASCII_INDEX) {
103: PetscViewerASCIIPrintf(viewer,"%d: ",i);
104: }
105: #if defined(PETSC_USE_COMPLEX)
106: if (PetscImaginaryPart(x->array[i]) > 0.0) {
107: PetscViewerASCIIPrintf(viewer,"%g + %g in",PetscRealPart(x->array[i]),PetscImaginaryPart(x->array[i]));
108: } else if (PetscImaginaryPart(x->array[i]) < 0.0) {
109: PetscViewerASCIIPrintf(viewer,"%g - %g in",PetscRealPart(x->array[i]),-PetscImaginaryPart(x->array[i]));
110: } else {
111: PetscViewerASCIIPrintf(viewer,"%gn",PetscRealPart(x->array[i]));
112: }
113: #else
114: PetscViewerASCIIPrintf(viewer,"%gn",x->array[i]);
115: #endif
116: }
117: }
118: PetscViewerFlush(viewer);
119: return(0);
120: }
122: static int VecView_Seq_Draw_LG(Vec xin,PetscViewer v)
123: {
124: Vec_Seq *x = (Vec_Seq *)xin->data;
125: int i,n = xin->n,ierr;
126: PetscDraw win;
127: PetscReal *xx;
128: PetscDrawLG lg;
131: PetscViewerDrawGetDrawLG(v,0,&lg);
132: PetscDrawLGGetDraw(lg,&win);
133: PetscDrawCheckResizedWindow(win);
134: PetscDrawLGReset(lg);
135: PetscMalloc((n+1)*sizeof(PetscReal),&xx);
136: for (i=0; i<n; i++) {
137: xx[i] = (PetscReal) i;
138: }
139: #if !defined(PETSC_USE_COMPLEX)
140: PetscDrawLGAddPoints(lg,n,&xx,&x->array);
141: #else
142: {
143: PetscReal *yy;
144: PetscMalloc((n+1)*sizeof(PetscReal),&yy);
145: for (i=0; i<n; i++) {
146: yy[i] = PetscRealPart(x->array[i]);
147: }
148: PetscDrawLGAddPoints(lg,n,&xx,&yy);
149: PetscFree(yy);
150: }
151: #endif
152: PetscFree(xx);
153: PetscDrawLGDraw(lg);
154: PetscDrawSynchronizedFlush(win);
155: return(0);
156: }
158: static int VecView_Seq_Draw(Vec xin,PetscViewer v)
159: {
160: int ierr;
161: PetscDraw draw;
162: PetscTruth isnull;
163: PetscViewerFormat format;
166: PetscViewerDrawGetDraw(v,0,&draw);
167: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
168:
169: PetscViewerGetFormat(v,&format);
170: /*
171: Currently it only supports drawing to a line graph */
172: if (format != PETSC_VIEWER_DRAW_LG) {
173: PetscViewerPushFormat(v,PETSC_VIEWER_DRAW_LG);
174: }
175: VecView_Seq_Draw_LG(xin,v);
176: if (format != PETSC_VIEWER_DRAW_LG) {
177: PetscViewerPopFormat(v);
178: }
180: return(0);
181: }
183: static int VecView_Seq_Binary(Vec xin,PetscViewer viewer)
184: {
185: Vec_Seq *x = (Vec_Seq *)xin->data;
186: int ierr,fdes,n = xin->n,cookie=VEC_FILE_COOKIE;
187: FILE *file;
190: ierr = PetscViewerBinaryGetDescriptor(viewer,&fdes);
191: /* Write vector header */
192: PetscBinaryWrite(fdes,&cookie,1,PETSC_INT,0);
193: PetscBinaryWrite(fdes,&n,1,PETSC_INT,0);
195: /* Write vector contents */
196: PetscBinaryWrite(fdes,x->array,n,PETSC_SCALAR,0);
198: PetscViewerBinaryGetInfoPointer(viewer,&file);
199: if (file && xin->bs > 1) {
200: fprintf(file,"-vecload_block_size %dn",xin->bs);
201: }
202: return(0);
203: }
206: int VecView_Seq(Vec xin,PetscViewer viewer)
207: {
208: Vec_Seq *x = (Vec_Seq *)xin->data;
209: int ierr;
210: PetscTruth isdraw,isascii,issocket,isbinary,ismathematica;
213: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
214: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
215: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
216: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
217: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_MATHEMATICA,&ismathematica);
218: if (isdraw){
219: VecView_Seq_Draw(xin,viewer);
220: } else if (isascii){
221: VecView_Seq_File(xin,viewer);
222: } else if (issocket) {
223: PetscViewerSocketPutScalar(viewer,xin->n,1,x->array);
224: } else if (isbinary) {
225: VecView_Seq_Binary(xin,viewer);
226: } else if (ismathematica) {
227: PetscViewerMathematicaPutVector(viewer,xin);
228: } else {
229: SETERRQ1(1,"Viewer type %s not supported by this vector object",((PetscObject)viewer)->type_name);
230: }
231: return(0);
232: }
234: int VecSetValues_Seq(Vec xin,int ni,const int ix[],const PetscScalar y[],InsertMode m)
235: {
236: Vec_Seq *x = (Vec_Seq *)xin->data;
237: PetscScalar *xx = x->array;
238: int i;
241: if (m == INSERT_VALUES) {
242: for (i=0; i<ni; i++) {
243: if (ix[i] < 0) continue;
244: #if defined(PETSC_USE_BOPT_g)
245: if (ix[i] >= xin->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",ix[i],xin->n);
246: #endif
247: xx[ix[i]] = y[i];
248: }
249: } else {
250: for (i=0; i<ni; i++) {
251: if (ix[i] < 0) continue;
252: #if defined(PETSC_USE_BOPT_g)
253: if (ix[i] >= xin->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",ix[i],xin->n);
254: #endif
255: xx[ix[i]] += y[i];
256: }
257: }
258: return(0);
259: }
261: int VecSetValuesBlocked_Seq(Vec xin,int ni,const int ix[],const PetscScalar yin[],InsertMode m)
262: {
263: Vec_Seq *x = (Vec_Seq *)xin->data;
264: PetscScalar *xx = x->array,*y = (PetscScalar*)yin;
265: int i,bs = xin->bs,start,j;
267: /*
268: For optimization could treat bs = 2, 3, 4, 5 as special cases with loop unrolling
269: */
271: if (m == INSERT_VALUES) {
272: for (i=0; i<ni; i++) {
273: start = bs*ix[i];
274: if (start < 0) continue;
275: #if defined(PETSC_USE_BOPT_g)
276: if (start >= xin->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",start,xin->n);
277: #endif
278: for (j=0; j<bs; j++) {
279: xx[start+j] = y[j];
280: }
281: y += bs;
282: }
283: } else {
284: for (i=0; i<ni; i++) {
285: start = bs*ix[i];
286: if (start < 0) continue;
287: #if defined(PETSC_USE_BOPT_g)
288: if (start >= xin->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",start,xin->n);
289: #endif
290: for (j=0; j<bs; j++) {
291: xx[start+j] += y[j];
292: }
293: y += bs;
294: }
295: }
296: return(0);
297: }
300: int VecDestroy_Seq(Vec v)
301: {
302: Vec_Seq *vs = (Vec_Seq*)v->data;
303: int ierr;
307: /* if memory was published with AMS then destroy it */
308: PetscObjectDepublish(v);
310: #if defined(PETSC_USE_LOG)
311: PetscLogObjectState((PetscObject)v,"Length=%d",v->n);
312: #endif
313: if (vs->array_allocated) {PetscFree(vs->array_allocated);}
314: PetscFree(vs);
316: return(0);
317: }
319: static int VecPublish_Seq(PetscObject obj)
320: {
321: #if defined(PETSC_HAVE_AMS)
322: Vec v = (Vec) obj;
323: Vec_Seq *s = (Vec_Seq*)v->data;
324: int ierr,(*f)(AMS_Memory,char *,Vec);
325: #endif
329: #if defined(PETSC_HAVE_AMS)
330: /* if it is already published then return */
331: if (v->amem >=0) return(0);
333: /* if array in vector was not allocated (for example PCSetUp_BJacobi_Singleblock()) then
334: cannot AMS publish the object*/
335: if (!s->array) return(0);
337: PetscObjectPublishBaseBegin(obj);
338: AMS_Memory_add_field((AMS_Memory)v->amem,"values",s->array,v->n,AMS_DOUBLE,AMS_READ,
339: AMS_DISTRIBUTED,AMS_REDUCT_UNDEF);
341: /* if the vector knows its "layout" let it set it*/
342: PetscObjectQueryFunction(obj,"AMSSetFieldBlock_C",(void (**)(void))&f);
343: if (f) {
344: (*f)((AMS_Memory)v->amem,"values",v);
345: }
346: PetscObjectPublishBaseEnd(obj);
347: #endif
349: return(0);
350: }
352: static struct _VecOps DvOps = {VecDuplicate_Seq,
353: VecDuplicateVecs_Default,
354: VecDestroyVecs_Default,
355: VecDot_Seq,
356: VecMDot_Seq,
357: VecNorm_Seq,
358: VecTDot_Seq,
359: VecMTDot_Seq,
360: VecScale_Seq,
361: VecCopy_Seq,
362: VecSet_Seq,
363: VecSwap_Seq,
364: VecAXPY_Seq,
365: VecAXPBY_Seq,
366: VecMAXPY_Seq,
367: VecAYPX_Seq,
368: VecWAXPY_Seq,
369: VecPointwiseMult_Seq,
370: VecPointwiseDivide_Seq,
371: VecSetValues_Seq,0,0,
372: VecGetArray_Seq,
373: VecGetSize_Seq,
374: VecGetSize_Seq,
375: VecRestoreArray_Seq,
376: VecMax_Seq,
377: VecMin_Seq,
378: VecSetRandom_Seq,0,
379: VecSetValuesBlocked_Seq,
380: VecDestroy_Seq,
381: VecView_Seq,
382: VecPlaceArray_Seq,
383: VecReplaceArray_Seq,
384: VecDot_Seq,
385: VecTDot_Seq,
386: VecNorm_Seq,
387: VecLoadIntoVector_Default,
388: VecReciprocal_Default,
389: 0, /* VecViewNative */
390: VecConjugate_Seq,
391: 0,
392: 0,
393: VecResetArray_Seq};
396: /*
397: This is called by VecCreate_Seq() (i.e. VecCreateSeq()) and VecCreateSeqWithArray()
398: */
399: static int VecCreate_Seq_Private(Vec v,const PetscScalar array[])
400: {
401: Vec_Seq *s;
402: int ierr;
405: PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));
406: PetscNew(Vec_Seq,&s);
407: PetscMemzero(s,sizeof(Vec_Seq));
408: v->data = (void*)s;
409: v->bops->publish = VecPublish_Seq;
410: v->n = PetscMax(v->n,v->N);
411: v->N = PetscMax(v->n,v->N);
412: v->bs = -1;
413: v->petscnative = PETSC_TRUE;
414: s->array = (PetscScalar *)array;
415: s->array_allocated = 0;
416: if (!v->map) {
417: PetscMapCreateMPI(v->comm,v->n,v->N,&v->map);
418: }
419: PetscObjectChangeTypeName((PetscObject)v,VECSEQ);
420: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
421: PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEnginePut_C","VecMatlabEnginePut_Default",VecMatlabEnginePut_Default);
422: PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscMatlabEngineGet_C","VecMatlabEngineGet_Default",VecMatlabEngineGet_Default);
423: #endif
424: PetscPublishAll(v);
425: return(0);
426: }
428: /*@C
429: VecCreateSeqWithArray - Creates a standard,sequential array-style vector,
430: where the user provides the array space to store the vector values.
432: Collective on MPI_Comm
434: Input Parameter:
435: + comm - the communicator, should be PETSC_COMM_SELF
436: . n - the vector length
437: - array - memory where the vector elements are to be stored.
439: Output Parameter:
440: . V - the vector
442: Notes:
443: Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
444: same type as an existing vector.
446: If the user-provided array is PETSC_NULL, then VecPlaceArray() can be used
447: at a later stage to SET the array for storing the vector values.
449: PETSc does NOT free the array when the vector is destroyed via VecDestroy().
450: The user should not free the array until the vector is destroyed.
452: Level: intermediate
454: Concepts: vectors^creating with array
456: .seealso: VecCreateMPIWithArray(), VecCreate(), VecDuplicate(), VecDuplicateVecs(),
457: VecCreateGhost(), VecCreateSeq(), VecPlaceArray()
458: @*/
459: int VecCreateSeqWithArray(MPI_Comm comm,int n,const PetscScalar array[],Vec *V)
460: {
461: int ierr;
464: VecCreate(comm,V);
465: VecSetSizes(*V,n,n);
466: VecCreate_Seq_Private(*V,array);
467: return(0);
468: }
470: EXTERN_C_BEGIN
471: int VecCreate_Seq(Vec V)
472: {
473: Vec_Seq *s;
474: PetscScalar *array;
475: int ierr,n = PetscMax(V->n,V->N);
478: PetscMalloc((n+1)*sizeof(PetscScalar),&array);
479: PetscMemzero(array,n*sizeof(PetscScalar));
480: VecCreate_Seq_Private(V,array);
481: s = (Vec_Seq*)V->data;
482: s->array_allocated = array;
483: VecSetSerializeType(V,VEC_SER_SEQ_BINARY);
484: return(0);
485: }
487: int VecSerialize_Seq(MPI_Comm comm, Vec *vec, PetscViewer viewer, PetscTruth store)
488: {
489: Vec v;
490: Vec_Seq *x;
491: PetscScalar *array;
492: int fd;
493: int vars;
494: int ierr;
497: PetscViewerBinaryGetDescriptor(viewer, &fd);
498: if (store) {
499: v = *vec;
500: x = (Vec_Seq *) v->data;
501: PetscBinaryWrite(fd, &v->n, 1, PETSC_INT, 0);
502: PetscBinaryWrite(fd, x->array, v->n, PETSC_SCALAR, 0);
503: } else {
504: PetscBinaryRead(fd, &vars, 1, PETSC_INT);
505: VecCreate(comm, &v);
506: VecSetSizes(v, vars, vars);
507: PetscMalloc(vars * sizeof(PetscScalar), &array);
508: PetscBinaryRead(fd, array, vars, PETSC_SCALAR);
509: VecCreate_Seq_Private(v, array);
510: ((Vec_Seq *) v->data)->array_allocated = array;
512: VecAssemblyBegin(v);
513: VecAssemblyEnd(v);
514: *vec = v;
515: }
517: return(0);
518: }
519: EXTERN_C_END
522: int VecDuplicate_Seq(Vec win,Vec *V)
523: {
524: int ierr;
527: VecCreateSeq(win->comm,win->n,V);
528: if (win->mapping) {
529: (*V)->mapping = win->mapping;
530: PetscObjectReference((PetscObject)win->mapping);
531: }
532: if (win->bmapping) {
533: (*V)->bmapping = win->bmapping;
534: PetscObjectReference((PetscObject)win->bmapping);
535: }
536: (*V)->bs = win->bs;
537: PetscOListDuplicate(win->olist,&(*V)->olist);
538: PetscFListDuplicate(win->qlist,&(*V)->qlist);
539: return(0);
540: }