Actual source code: pdvec.c
1: /* $Id: pdvec.c,v 1.154 2001/09/11 16:32:01 bsmith Exp $*/
2: /*
3: Code for some of the parallel vector primatives.
4: */
5: #include src/vec/impls/mpi/pvecimpl.h
7: int VecDestroy_MPI(Vec v)
8: {
9: Vec_MPI *x = (Vec_MPI*)v->data;
10: int ierr;
14: /* if memory was published with AMS then destroy it */
15: PetscObjectDepublish(v);
17: #if defined(PETSC_USE_LOG)
18: PetscLogObjectState((PetscObject)v,"Length=%d",v->N);
19: #endif
20: if (x->array_allocated) {PetscFree(x->array_allocated);}
22: /* Destroy local representation of vector if it exists */
23: if (x->localrep) {
24: VecDestroy(x->localrep);
25: if (x->localupdate) {VecScatterDestroy(x->localupdate);}
26: }
27: /* Destroy the stashes: note the order - so that the tags are freed properly */
28: VecStashDestroy_Private(&v->bstash);
29: VecStashDestroy_Private(&v->stash);
30: PetscFree(x);
31: return(0);
32: }
34: int VecView_MPI_ASCII(Vec xin,PetscViewer viewer)
35: {
36: int i,rank,len,work = xin->n,n,j,size,ierr,cnt,tag = ((PetscObject)viewer)->tag;
37: MPI_Status status;
38: PetscScalar *values,*xarray;
39: char *name;
40: PetscViewerFormat format;
43: VecGetArrayFast(xin,&xarray);
44: /* determine maximum message to arrive */
45: MPI_Comm_rank(xin->comm,&rank);
46: MPI_Reduce(&work,&len,1,MPI_INT,MPI_MAX,0,xin->comm);
47: MPI_Comm_size(xin->comm,&size);
49: if (!rank) {
50: PetscMalloc((len+1)*sizeof(PetscScalar),&values);
51: PetscViewerGetFormat(viewer,&format);
52: /*
53: Matlab format and ASCII format are very similar except
54: Matlab uses %18.16e format while ASCII uses %g
55: */
56: if (format == PETSC_VIEWER_ASCII_MATLAB) {
57: PetscObjectGetName((PetscObject)xin,&name);
58: PetscViewerASCIIPrintf(viewer,"%s = [n",name);
59: for (i=0; i<xin->n; i++) {
60: #if defined(PETSC_USE_COMPLEX)
61: if (PetscImaginaryPart(xarray[i]) > 0.0) {
62: PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ein",PetscRealPart(xarray[i]),PetscImaginaryPart(xarray[i]));
63: } else if (PetscImaginaryPart(xarray[i]) < 0.0) {
64: PetscViewerASCIIPrintf(viewer,"%18.16e - %18.16ein",PetscRealPart(xarray[i]),-PetscImaginaryPart(xarray[i]));
65: } else {
66: PetscViewerASCIIPrintf(viewer,"%18.16en",PetscRealPart(xarray[i]));
67: }
68: #else
69: PetscViewerASCIIPrintf(viewer,"%18.16en",xarray[i]);
70: #endif
71: }
72: /* receive and print messages */
73: for (j=1; j<size; j++) {
74: MPI_Recv(values,len,MPIU_SCALAR,j,tag,xin->comm,&status);
75: MPI_Get_count(&status,MPIU_SCALAR,&n);
76: for (i=0; i<n; i++) {
77: #if defined(PETSC_USE_COMPLEX)
78: if (PetscImaginaryPart(values[i]) > 0.0) {
79: PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e in",PetscRealPart(values[i]),PetscImaginaryPart(values[i]));
80: } else if (PetscImaginaryPart(values[i]) < 0.0) {
81: PetscViewerASCIIPrintf(viewer,"%18.16e - %18.16e in",PetscRealPart(values[i]),-PetscImaginaryPart(values[i]));
82: } else {
83: PetscViewerASCIIPrintf(viewer,"%18.16en",PetscRealPart(values[i]));
84: }
85: #else
86: PetscViewerASCIIPrintf(viewer,"%18.16en",values[i]);
87: #endif
88: }
89: }
90: PetscViewerASCIIPrintf(viewer,"];n");
92: } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
93: for (i=0; i<xin->n; i++) {
94: #if defined(PETSC_USE_COMPLEX)
95: PetscViewerASCIIPrintf(viewer,"%18.16e %18.16en",PetscRealPart(xarray[i]),PetscImaginaryPart(xarray[i]));
96: #else
97: PetscViewerASCIIPrintf(viewer,"%18.16en",xarray[i]);
98: #endif
99: }
100: /* receive and print messages */
101: for (j=1; j<size; j++) {
102: MPI_Recv(values,len,MPIU_SCALAR,j,tag,xin->comm,&status);
103: MPI_Get_count(&status,MPIU_SCALAR,&n);
104: for (i=0; i<n; i++) {
105: #if defined(PETSC_USE_COMPLEX)
106: PetscViewerASCIIPrintf(viewer,"%18.16e %18.16en",PetscRealPart(values[i]),PetscImaginaryPart(values[i]));
107: #else
108: PetscViewerASCIIPrintf(viewer,"%18.16en",values[i]);
109: #endif
110: }
111: }
113: } else {
114: if (format != PETSC_VIEWER_ASCII_COMMON) {PetscViewerASCIIPrintf(viewer,"Processor [%d]n",rank);}
115: cnt = 0;
116: for (i=0; i<xin->n; i++) {
117: if (format == PETSC_VIEWER_ASCII_INDEX) {
118: PetscViewerASCIIPrintf(viewer,"%d: ",cnt++);
119: }
120: #if defined(PETSC_USE_COMPLEX)
121: if (PetscImaginaryPart(xarray[i]) > 0.0) {
122: PetscViewerASCIIPrintf(viewer,"%g + %g in",PetscRealPart(xarray[i]),PetscImaginaryPart(xarray[i]));
123: } else if (PetscImaginaryPart(xarray[i]) < 0.0) {
124: PetscViewerASCIIPrintf(viewer,"%g - %g in",PetscRealPart(xarray[i]),-PetscImaginaryPart(xarray[i]));
125: } else {
126: PetscViewerASCIIPrintf(viewer,"%gn",PetscRealPart(xarray[i]));
127: }
128: #else
129: PetscViewerASCIIPrintf(viewer,"%gn",xarray[i]);
130: #endif
131: }
132: /* receive and print messages */
133: for (j=1; j<size; j++) {
134: MPI_Recv(values,len,MPIU_SCALAR,j,tag,xin->comm,&status);
135: MPI_Get_count(&status,MPIU_SCALAR,&n);
136: if (format != PETSC_VIEWER_ASCII_COMMON) {
137: PetscViewerASCIIPrintf(viewer,"Processor [%d]n",j);
138: }
139: for (i=0; i<n; i++) {
140: if (format == PETSC_VIEWER_ASCII_INDEX) {
141: PetscViewerASCIIPrintf(viewer,"%d: ",cnt++);
142: }
143: #if defined(PETSC_USE_COMPLEX)
144: if (PetscImaginaryPart(values[i]) > 0.0) {
145: PetscViewerASCIIPrintf(viewer,"%g + %g in",PetscRealPart(values[i]),PetscImaginaryPart(values[i]));
146: } else if (PetscImaginaryPart(values[i]) < 0.0) {
147: PetscViewerASCIIPrintf(viewer,"%g - %g in",PetscRealPart(values[i]),-PetscImaginaryPart(values[i]));
148: } else {
149: PetscViewerASCIIPrintf(viewer,"%gn",PetscRealPart(values[i]));
150: }
151: #else
152: PetscViewerASCIIPrintf(viewer,"%gn",values[i]);
153: #endif
154: }
155: }
156: }
157: PetscFree(values);
158: } else {
159: /* send values */
160: MPI_Send(xarray,xin->n,MPIU_SCALAR,0,tag,xin->comm);
161: }
162: PetscViewerFlush(viewer);
163: VecRestoreArrayFast(xin,&xarray);
164: return(0);
165: }
167: int VecView_MPI_Binary(Vec xin,PetscViewer viewer)
168: {
169: int rank,ierr,len,work = xin->n,n,j,size,fdes,tag = ((PetscObject)viewer)->tag;
170: MPI_Status status;
171: PetscScalar *values,*xarray;
172: FILE *file;
175: VecGetArrayFast(xin,&xarray);
176: PetscViewerBinaryGetDescriptor(viewer,&fdes);
178: /* determine maximum message to arrive */
179: MPI_Comm_rank(xin->comm,&rank);
180: MPI_Reduce(&work,&len,1,MPI_INT,MPI_MAX,0,xin->comm);
181: MPI_Comm_size(xin->comm,&size);
183: if (!rank) {
184: int cookie = VEC_FILE_COOKIE;
185: PetscBinaryWrite(fdes,&cookie,1,PETSC_INT,0);
186: PetscBinaryWrite(fdes,&xin->N,1,PETSC_INT,0);
187: PetscBinaryWrite(fdes,xarray,xin->n,PETSC_SCALAR,0);
189: PetscMalloc((len+1)*sizeof(PetscScalar),&values);
190: /* receive and print messages */
191: for (j=1; j<size; j++) {
192: MPI_Recv(values,len,MPIU_SCALAR,j,tag,xin->comm,&status);
193: MPI_Get_count(&status,MPIU_SCALAR,&n);
194: PetscBinaryWrite(fdes,values,n,PETSC_SCALAR,0);
195: }
196: PetscFree(values);
197: PetscViewerBinaryGetInfoPointer(viewer,&file);
198: if (file && xin->bs > 1) {
199: fprintf(file,"-vecload_block_size %dn",xin->bs);
200: }
201: } else {
202: /* send values */
203: MPI_Send(xarray,xin->n,MPIU_SCALAR,0,tag,xin->comm);
204: }
205: VecRestoreArrayFast(xin,&xarray);
206: return(0);
207: }
209: int VecView_MPI_Draw_LG(Vec xin,PetscViewer viewer)
210: {
211: int i,rank,size,N = xin->N,*lens,ierr;
212: PetscDraw draw;
213: PetscReal *xx,*yy;
214: PetscDrawLG lg;
215: PetscTruth isnull;
216: PetscScalar *xarray;
219: PetscViewerDrawGetDraw(viewer,0,&draw);
220: PetscDrawIsNull(draw,&isnull);
221: if (isnull) return(0);
223: VecGetArrayFast(xin,&xarray);
224: PetscViewerDrawGetDrawLG(viewer,0,&lg);
225: PetscDrawCheckResizedWindow(draw);
226: MPI_Comm_rank(xin->comm,&rank);
227: MPI_Comm_size(xin->comm,&size);
228: if (!rank) {
229: PetscDrawLGReset(lg);
230: PetscMalloc(2*(N+1)*sizeof(PetscReal),&xx);
231: for (i=0; i<N; i++) {xx[i] = (PetscReal) i;}
232: yy = xx + N;
233: PetscMalloc(size*sizeof(int),&lens);
234: for (i=0; i<size; i++) {
235: lens[i] = xin->map->range[i+1] - xin->map->range[i];
236: }
237: #if !defined(PETSC_USE_COMPLEX)
238: MPI_Gatherv(xarray,xin->n,MPIU_REAL,yy,lens,xin->map->range,MPIU_REAL,0,xin->comm);
239: #else
240: {
241: PetscReal *xr;
242: PetscMalloc((xin->n+1)*sizeof(PetscReal),&xr);
243: for (i=0; i<xin->n; i++) {
244: xr[i] = PetscRealPart(xarray[i]);
245: }
246: MPI_Gatherv(xr,xin->n,MPIU_REAL,yy,lens,xin->map->range,MPIU_REAL,0,xin->comm);
247: PetscFree(xr);
248: }
249: #endif
250: PetscFree(lens);
251: PetscDrawLGAddPoints(lg,N,&xx,&yy);
252: PetscFree(xx);
253: } else {
254: #if !defined(PETSC_USE_COMPLEX)
255: MPI_Gatherv(xarray,xin->n,MPIU_REAL,0,0,0,MPIU_REAL,0,xin->comm);
256: #else
257: {
258: PetscReal *xr;
259: PetscMalloc((xin->n+1)*sizeof(PetscReal),&xr);
260: for (i=0; i<xin->n; i++) {
261: xr[i] = PetscRealPart(xarray[i]);
262: }
263: MPI_Gatherv(xr,xin->n,MPIU_REAL,0,0,0,MPIU_REAL,0,xin->comm);
264: PetscFree(xr);
265: }
266: #endif
267: }
268: PetscDrawLGDraw(lg);
269: PetscDrawSynchronizedFlush(draw);
270: VecRestoreArrayFast(xin,&xarray);
271: return(0);
272: }
274: EXTERN_C_BEGIN
275: int VecView_MPI_Draw(Vec xin,PetscViewer viewer)
276: {
277: int i,rank,size,ierr,start,end,tag = ((PetscObject)viewer)->tag;
278: MPI_Status status;
279: PetscReal coors[4],ymin,ymax,xmin,xmax,tmp;
280: PetscDraw draw;
281: PetscTruth isnull;
282: PetscDrawAxis axis;
283: PetscScalar *xarray;
286: PetscViewerDrawGetDraw(viewer,0,&draw);
287: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
289: VecGetArrayFast(xin,&xarray);
290: PetscDrawCheckResizedWindow(draw);
291: xmin = 1.e20; xmax = -1.e20;
292: for (i=0; i<xin->n; i++) {
293: #if defined(PETSC_USE_COMPLEX)
294: if (PetscRealPart(xarray[i]) < xmin) xmin = PetscRealPart(xarray[i]);
295: if (PetscRealPart(xarray[i]) > xmax) xmax = PetscRealPart(xarray[i]);
296: #else
297: if (xarray[i] < xmin) xmin = xarray[i];
298: if (xarray[i] > xmax) xmax = xarray[i];
299: #endif
300: }
301: if (xmin + 1.e-10 > xmax) {
302: xmin -= 1.e-5;
303: xmax += 1.e-5;
304: }
305: MPI_Reduce(&xmin,&ymin,1,MPIU_REAL,MPI_MIN,0,xin->comm);
306: MPI_Reduce(&xmax,&ymax,1,MPIU_REAL,MPI_MAX,0,xin->comm);
307: MPI_Comm_size(xin->comm,&size);
308: MPI_Comm_rank(xin->comm,&rank);
309: PetscDrawAxisCreate(draw,&axis);
310: PetscLogObjectParent(draw,axis);
311: if (!rank) {
312: PetscDrawClear(draw);
313: PetscDrawFlush(draw);
314: PetscDrawAxisSetLimits(axis,0.0,(double)xin->N,ymin,ymax);
315: PetscDrawAxisDraw(axis);
316: PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);
317: }
318: PetscDrawAxisDestroy(axis);
319: MPI_Bcast(coors,4,MPIU_REAL,0,xin->comm);
320: if (rank) {PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);}
321: /* draw local part of vector */
322: VecGetOwnershipRange(xin,&start,&end);
323: if (rank < size-1) { /*send value to right */
324: MPI_Send(&xarray[xin->n-1],1,MPIU_REAL,rank+1,tag,xin->comm);
325: }
326: for (i=1; i<xin->n; i++) {
327: #if !defined(PETSC_USE_COMPLEX)
328: PetscDrawLine(draw,(PetscReal)(i-1+start),xarray[i-1],(PetscReal)(i+start),
329: xarray[i],PETSC_DRAW_RED);
330: #else
331: PetscDrawLine(draw,(PetscReal)(i-1+start),PetscRealPart(xarray[i-1]),(PetscReal)(i+start),
332: PetscRealPart(xarray[i]),PETSC_DRAW_RED);
333: #endif
334: }
335: if (rank) { /* receive value from right */
336: MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag,xin->comm,&status);
337: #if !defined(PETSC_USE_COMPLEX)
338: PetscDrawLine(draw,(PetscReal)start-1,tmp,(PetscReal)start,xarray[0],PETSC_DRAW_RED);
339: #else
340: PetscDrawLine(draw,(PetscReal)start-1,tmp,(PetscReal)start,PetscRealPart(xarray[0]),PETSC_DRAW_RED);
341: #endif
342: }
343: PetscDrawSynchronizedFlush(draw);
344: PetscDrawPause(draw);
345: VecRestoreArrayFast(xin,&xarray);
346: return(0);
347: }
348: EXTERN_C_END
350: int VecView_MPI_Socket(Vec xin,PetscViewer viewer)
351: {
352: int i,rank,size,N = xin->N,*lens,ierr;
353: PetscScalar *xx,*xarray;
356: VecGetArrayFast(xin,&xarray);
357: MPI_Comm_rank(xin->comm,&rank);
358: MPI_Comm_size(xin->comm,&size);
359: if (!rank) {
360: PetscMalloc((N+1)*sizeof(PetscScalar),&xx);
361: PetscMalloc(size*sizeof(int),&lens);
362: for (i=0; i<size; i++) {
363: lens[i] = xin->map->range[i+1] - xin->map->range[i];
364: }
365: MPI_Gatherv(xarray,xin->n,MPIU_SCALAR,xx,lens,xin->map->range,MPIU_SCALAR,0,xin->comm);
366: PetscFree(lens);
367: PetscViewerSocketPutScalar(viewer,N,1,xx);
368: PetscFree(xx);
369: } else {
370: MPI_Gatherv(xarray,xin->n,MPIU_REAL,0,0,0,MPIU_REAL,0,xin->comm);
371: }
372: VecRestoreArrayFast(xin,&xarray);
373: return(0);
374: }
376: int VecView_MPI(Vec xin,PetscViewer viewer)
377: {
378: int ierr;
379: PetscTruth isascii,issocket,isbinary,isdraw,ismathematica;
382: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
383: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
384: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
385: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
386: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_MATHEMATICA,&ismathematica);
387: if (isascii){
388: VecView_MPI_ASCII(xin,viewer);
389: } else if (issocket) {
390: VecView_MPI_Socket(xin,viewer);
391: } else if (isbinary) {
392: VecView_MPI_Binary(xin,viewer);
393: } else if (isdraw) {
394: PetscViewerFormat format;
396: PetscViewerGetFormat(viewer,&format);
397: if (format == PETSC_VIEWER_DRAW_LG) {
398: VecView_MPI_Draw_LG(xin,viewer);
399: } else {
400: VecView_MPI_Draw(xin,viewer);
401: }
402: } else if (ismathematica) {
403: PetscViewerMathematicaPutVector(viewer,xin);
404: } else {
405: SETERRQ1(1,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
406: }
407: return(0);
408: }
410: int VecGetSize_MPI(Vec xin,int *N)
411: {
413: *N = xin->N;
414: return(0);
415: }
417: int VecSetValues_MPI(Vec xin,int ni,const int ix[],const PetscScalar y[],InsertMode addv)
418: {
419: int rank = xin->stash.rank, *owners = xin->map->range,start = owners[rank];
420: int end = owners[rank+1],i,row,ierr;
421: PetscScalar *xx;
424: #if defined(PETSC_USE_BOPT_g)
425: if (xin->stash.insertmode == INSERT_VALUES && addv == ADD_VALUES) {
426: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You have already inserted values; you cannot now add");
427: } else if (xin->stash.insertmode == ADD_VALUES && addv == INSERT_VALUES) {
428: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You have already added values; you cannot now insert");
429: }
430: #endif
431: VecGetArrayFast(xin,&xx);
432: xin->stash.insertmode = addv;
434: if (addv == INSERT_VALUES) {
435: for (i=0; i<ni; i++) {
436: if ((row = ix[i]) >= start && row < end) {
437: xx[row-start] = y[i];
438: } else if (!xin->stash.donotstash) {
439: if (ix[i] < 0) continue;
440: #if defined(PETSC_USE_BOPT_g)
441: if (ix[i] >= xin->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",ix[i],xin->N);
442: #endif
443: VecStashValue_Private(&xin->stash,row,y[i]);
444: }
445: }
446: } else {
447: for (i=0; i<ni; i++) {
448: if ((row = ix[i]) >= start && row < end) {
449: xx[row-start] += y[i];
450: } else if (!xin->stash.donotstash) {
451: if (ix[i] < 0) continue;
452: #if defined(PETSC_USE_BOPT_g)
453: if (ix[i] > xin->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",ix[i],xin->N);
454: #endif
455: VecStashValue_Private(&xin->stash,row,y[i]);
456: }
457: }
458: }
459: VecRestoreArrayFast(xin,&xx);
460: return(0);
461: }
463: int VecSetValuesBlocked_MPI(Vec xin,int ni,const int ix[],const PetscScalar yin[],InsertMode addv)
464: {
465: int rank = xin->stash.rank,*owners = xin->map->range,start = owners[rank];
466: int end = owners[rank+1],i,row,bs = xin->bs,j,ierr;
467: PetscScalar *xx,*y = (PetscScalar*)yin;
470: VecGetArrayFast(xin,&xx);
471: #if defined(PETSC_USE_BOPT_g)
472: if (xin->stash.insertmode == INSERT_VALUES && addv == ADD_VALUES) {
473: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You have already inserted values; you cannot now add");
474: }
475: else if (xin->stash.insertmode == ADD_VALUES && addv == INSERT_VALUES) {
476: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You have already added values; you cannot now insert");
477: }
478: #endif
479: xin->stash.insertmode = addv;
481: if (addv == INSERT_VALUES) {
482: for (i=0; i<ni; i++) {
483: if ((row = bs*ix[i]) >= start && row < end) {
484: for (j=0; j<bs; j++) {
485: xx[row-start+j] = y[j];
486: }
487: } else if (!xin->stash.donotstash) {
488: if (ix[i] < 0) continue;
489: #if defined(PETSC_USE_BOPT_g)
490: if (ix[i] >= xin->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d max %d",ix[i],xin->N);
491: #endif
492: VecStashValuesBlocked_Private(&xin->bstash,ix[i],y);
493: }
494: y += bs;
495: }
496: } else {
497: for (i=0; i<ni; i++) {
498: if ((row = bs*ix[i]) >= start && row < end) {
499: for (j=0; j<bs; j++) {
500: xx[row-start+j] += y[j];
501: }
502: } else if (!xin->stash.donotstash) {
503: if (ix[i] < 0) continue;
504: #if defined(PETSC_USE_BOPT_g)
505: if (ix[i] > xin->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d max %d",ix[i],xin->N);
506: #endif
507: VecStashValuesBlocked_Private(&xin->bstash,ix[i],y);
508: }
509: y += bs;
510: }
511: }
512: VecRestoreArrayFast(xin,&xx);
513: return(0);
514: }
516: /*
517: Since nsends or nreceives may be zero we add 1 in certain mallocs
518: to make sure we never malloc an empty one.
519: */
520: int VecAssemblyBegin_MPI(Vec xin)
521: {
522: int *owners = xin->map->range,*bowners,ierr,size,i,bs,nstash,reallocs;
523: InsertMode addv;
524: MPI_Comm comm = xin->comm;
527: if (xin->stash.donotstash) {
528: return(0);
529: }
531: MPI_Allreduce(&xin->stash.insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
532: if (addv == (ADD_VALUES|INSERT_VALUES)) {
533: SETERRQ(PETSC_ERR_ARG_NOTSAMETYPE,"Some processors inserted values while others added");
534: }
535: xin->stash.insertmode = addv; /* in case this processor had no cache */
536:
537: bs = xin->bs;
538: MPI_Comm_size(xin->comm,&size);
539: if (!xin->bstash.bowners && xin->bs != -1) {
540: PetscMalloc((size+1)*sizeof(int),&bowners);
541: for (i=0; i<size+1; i++){ bowners[i] = owners[i]/bs;}
542: xin->bstash.bowners = bowners;
543: } else {
544: bowners = xin->bstash.bowners;
545: }
546: VecStashScatterBegin_Private(&xin->stash,owners);
547: VecStashScatterBegin_Private(&xin->bstash,bowners);
548: VecStashGetInfo_Private(&xin->stash,&nstash,&reallocs);
549: PetscLogInfo(0,"VecAssemblyBegin_MPI:Stash has %d entries, uses %d mallocs.n",nstash,reallocs);
550: VecStashGetInfo_Private(&xin->bstash,&nstash,&reallocs);
551: PetscLogInfo(0,"VecAssemblyBegin_MPI:Block-Stash has %d entries, uses %d mallocs.n",nstash,reallocs);
553: return(0);
554: }
556: int VecAssemblyEnd_MPI(Vec vec)
557: {
558: int ierr,base,i,j,n,*row,flg,bs;
559: PetscScalar *val,*vv,*array,*xarray;
562: if (!vec->stash.donotstash) {
563: VecGetArrayFast(vec,&xarray);
564: base = vec->map->range[vec->stash.rank];
565: bs = vec->bs;
567: /* Process the stash */
568: while (1) {
569: VecStashScatterGetMesg_Private(&vec->stash,&n,&row,&val,&flg);
570: if (!flg) break;
571: if (vec->stash.insertmode == ADD_VALUES) {
572: for (i=0; i<n; i++) { xarray[row[i] - base] += val[i]; }
573: } else if (vec->stash.insertmode == INSERT_VALUES) {
574: for (i=0; i<n; i++) { xarray[row[i] - base] = val[i]; }
575: } else {
576: SETERRQ(PETSC_ERR_ARG_CORRUPT,"Insert mode is not set correctly; corrupted vector");
577: }
578: }
579: VecStashScatterEnd_Private(&vec->stash);
581: /* now process the block-stash */
582: while (1) {
583: VecStashScatterGetMesg_Private(&vec->bstash,&n,&row,&val,&flg);
584: if (!flg) break;
585: for (i=0; i<n; i++) {
586: array = xarray+row[i]*bs-base;
587: vv = val+i*bs;
588: if (vec->stash.insertmode == ADD_VALUES) {
589: for (j=0; j<bs; j++) { array[j] += vv[j];}
590: } else if (vec->stash.insertmode == INSERT_VALUES) {
591: for (j=0; j<bs; j++) { array[j] = vv[j]; }
592: } else {
593: SETERRQ(PETSC_ERR_ARG_CORRUPT,"Insert mode is not set correctly; corrupted vector");
594: }
595: }
596: }
597: VecStashScatterEnd_Private(&vec->bstash);
598: VecRestoreArrayFast(vec,&xarray);
599: }
600: vec->stash.insertmode = NOT_SET_VALUES;
601: return(0);
602: }