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: #undef __FUNCT__  
  9: int VecDestroy_MPI(Vec v)
 10: {
 11:   Vec_MPI *x = (Vec_MPI*)v->data;
 12:   int     ierr;


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

 19: #if defined(PETSC_USE_LOG)
 20:   PetscLogObjectState((PetscObject)v,"Length=%d",v->N);
 21: #endif  
 22:   if (x->array_allocated) {PetscFree(x->array_allocated);}

 24:   /* Destroy local representation of vector if it exists */
 25:   if (x->localrep) {
 26:     VecDestroy(x->localrep);
 27:     if (x->localupdate) {VecScatterDestroy(x->localupdate);}
 28:   }
 29:   /* Destroy the stashes: note the order - so that the tags are freed properly */
 30:   VecStashDestroy_Private(&v->bstash);
 31:   VecStashDestroy_Private(&v->stash);
 32:   PetscFree(x);
 33:   return(0);
 34: }

 36: #undef __FUNCT__  
 38: int VecView_MPI_ASCII(Vec xin,PetscViewer viewer)
 39: {
 40:   int               i,rank,len,work = xin->n,n,j,size,ierr,cnt,tag = ((PetscObject)viewer)->tag;
 41:   MPI_Status        status;
 42:   PetscScalar       *values,*xarray;
 43:   char              *name;
 44:   PetscViewerFormat format;

 47:   VecGetArrayFast(xin,&xarray);
 48:   /* determine maximum message to arrive */
 49:   MPI_Comm_rank(xin->comm,&rank);
 50:   MPI_Reduce(&work,&len,1,MPI_INT,MPI_MAX,0,xin->comm);
 51:   MPI_Comm_size(xin->comm,&size);

 53:   if (!rank) {
 54:     PetscMalloc((len+1)*sizeof(PetscScalar),&values);
 55:     PetscViewerGetFormat(viewer,&format);
 56:     /*
 57:         Matlab format and ASCII format are very similar except 
 58:         Matlab uses %18.16e format while ASCII uses %g
 59:     */
 60:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
 61:       PetscObjectGetName((PetscObject)xin,&name);
 62:       PetscViewerASCIIPrintf(viewer,"%s = [n",name);
 63:       for (i=0; i<xin->n; i++) {
 64: #if defined(PETSC_USE_COMPLEX)
 65:         if (PetscImaginaryPart(xarray[i]) > 0.0) {
 66:           PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ein",PetscRealPart(xarray[i]),PetscImaginaryPart(xarray[i]));
 67:         } else if (PetscImaginaryPart(xarray[i]) < 0.0) {
 68:           PetscViewerASCIIPrintf(viewer,"%18.16e - %18.16ein",PetscRealPart(xarray[i]),-PetscImaginaryPart(xarray[i]));
 69:         } else {
 70:           PetscViewerASCIIPrintf(viewer,"%18.16en",PetscRealPart(xarray[i]));
 71:         }
 72: #else
 73:         PetscViewerASCIIPrintf(viewer,"%18.16en",xarray[i]);
 74: #endif
 75:       }
 76:       /* receive and print messages */
 77:       for (j=1; j<size; j++) {
 78:         MPI_Recv(values,len,MPIU_SCALAR,j,tag,xin->comm,&status);
 79:         MPI_Get_count(&status,MPIU_SCALAR,&n);
 80:         for (i=0; i<n; i++) {
 81: #if defined(PETSC_USE_COMPLEX)
 82:           if (PetscImaginaryPart(values[i]) > 0.0) {
 83:             PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e in",PetscRealPart(values[i]),PetscImaginaryPart(values[i]));
 84:           } else if (PetscImaginaryPart(values[i]) < 0.0) {
 85:             PetscViewerASCIIPrintf(viewer,"%18.16e - %18.16e in",PetscRealPart(values[i]),-PetscImaginaryPart(values[i]));
 86:           } else {
 87:             PetscViewerASCIIPrintf(viewer,"%18.16en",PetscRealPart(values[i]));
 88:           }
 89: #else
 90:           PetscViewerASCIIPrintf(viewer,"%18.16en",values[i]);
 91: #endif
 92:         }
 93:       }
 94:       PetscViewerASCIIPrintf(viewer,"];n");

 96:     } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
 97:       for (i=0; i<xin->n; i++) {
 98: #if defined(PETSC_USE_COMPLEX)
 99:         PetscViewerASCIIPrintf(viewer,"%18.16e %18.16en",PetscRealPart(xarray[i]),PetscImaginaryPart(xarray[i]));
100: #else
101:         PetscViewerASCIIPrintf(viewer,"%18.16en",xarray[i]);
102: #endif
103:       }
104:       /* receive and print messages */
105:       for (j=1; j<size; j++) {
106:         MPI_Recv(values,len,MPIU_SCALAR,j,tag,xin->comm,&status);
107:         MPI_Get_count(&status,MPIU_SCALAR,&n);
108:         for (i=0; i<n; i++) {
109: #if defined(PETSC_USE_COMPLEX)
110:           PetscViewerASCIIPrintf(viewer,"%18.16e %18.16en",PetscRealPart(values[i]),PetscImaginaryPart(values[i]));
111: #else
112:           PetscViewerASCIIPrintf(viewer,"%18.16en",values[i]);
113: #endif
114:         }
115:       }

117:     } else {
118:       if (format != PETSC_VIEWER_ASCII_COMMON) {PetscViewerASCIIPrintf(viewer,"Processor [%d]n",rank);}
119:       cnt = 0;
120:       for (i=0; i<xin->n; i++) {
121:         if (format == PETSC_VIEWER_ASCII_INDEX) {
122:           PetscViewerASCIIPrintf(viewer,"%d: ",cnt++);
123:         }
124: #if defined(PETSC_USE_COMPLEX)
125:         if (PetscImaginaryPart(xarray[i]) > 0.0) {
126:           PetscViewerASCIIPrintf(viewer,"%g + %g in",PetscRealPart(xarray[i]),PetscImaginaryPart(xarray[i]));
127:         } else if (PetscImaginaryPart(xarray[i]) < 0.0) {
128:           PetscViewerASCIIPrintf(viewer,"%g - %g in",PetscRealPart(xarray[i]),-PetscImaginaryPart(xarray[i]));
129:         } else {
130:           PetscViewerASCIIPrintf(viewer,"%gn",PetscRealPart(xarray[i]));
131:         }
132: #else
133:         PetscViewerASCIIPrintf(viewer,"%gn",xarray[i]);
134: #endif
135:       }
136:       /* receive and print messages */
137:       for (j=1; j<size; j++) {
138:         MPI_Recv(values,len,MPIU_SCALAR,j,tag,xin->comm,&status);
139:         MPI_Get_count(&status,MPIU_SCALAR,&n);
140:         if (format != PETSC_VIEWER_ASCII_COMMON) {
141:           PetscViewerASCIIPrintf(viewer,"Processor [%d]n",j);
142:         }
143:         for (i=0; i<n; i++) {
144:           if (format == PETSC_VIEWER_ASCII_INDEX) {
145:             PetscViewerASCIIPrintf(viewer,"%d: ",cnt++);
146:           }
147: #if defined(PETSC_USE_COMPLEX)
148:           if (PetscImaginaryPart(values[i]) > 0.0) {
149:             PetscViewerASCIIPrintf(viewer,"%g + %g in",PetscRealPart(values[i]),PetscImaginaryPart(values[i]));
150:           } else if (PetscImaginaryPart(values[i]) < 0.0) {
151:             PetscViewerASCIIPrintf(viewer,"%g - %g in",PetscRealPart(values[i]),-PetscImaginaryPart(values[i]));
152:           } else {
153:             PetscViewerASCIIPrintf(viewer,"%gn",PetscRealPart(values[i]));
154:           }
155: #else
156:           PetscViewerASCIIPrintf(viewer,"%gn",values[i]);
157: #endif
158:         }
159:       }
160:     }
161:     PetscFree(values);
162:   } else {
163:     /* send values */
164:     MPI_Send(xarray,xin->n,MPIU_SCALAR,0,tag,xin->comm);
165:   }
166:   PetscViewerFlush(viewer);
167:   VecRestoreArrayFast(xin,&xarray);
168:   return(0);
169: }

171: #undef __FUNCT__  
173: int VecView_MPI_Binary(Vec xin,PetscViewer viewer)
174: {
175:   int         rank,ierr,len,work = xin->n,n,j,size,fdes,tag = ((PetscObject)viewer)->tag;
176:   MPI_Status  status;
177:   PetscScalar *values,*xarray;
178:   FILE        *file;

181:   VecGetArrayFast(xin,&xarray);
182:   PetscViewerBinaryGetDescriptor(viewer,&fdes);

184:   /* determine maximum message to arrive */
185:   MPI_Comm_rank(xin->comm,&rank);
186:   MPI_Reduce(&work,&len,1,MPI_INT,MPI_MAX,0,xin->comm);
187:   MPI_Comm_size(xin->comm,&size);

189:   if (!rank) {
190:     int cookie = VEC_FILE_COOKIE;
191:     PetscBinaryWrite(fdes,&cookie,1,PETSC_INT,0);
192:     PetscBinaryWrite(fdes,&xin->N,1,PETSC_INT,0);
193:     PetscBinaryWrite(fdes,xarray,xin->n,PETSC_SCALAR,0);

195:     PetscMalloc((len+1)*sizeof(PetscScalar),&values);
196:     /* receive and print messages */
197:     for (j=1; j<size; j++) {
198:       MPI_Recv(values,len,MPIU_SCALAR,j,tag,xin->comm,&status);
199:       MPI_Get_count(&status,MPIU_SCALAR,&n);
200:       PetscBinaryWrite(fdes,values,n,PETSC_SCALAR,0);
201:     }
202:     PetscFree(values);
203:     PetscViewerBinaryGetInfoPointer(viewer,&file);
204:     if (file && xin->bs > 1) {
205:       fprintf(file,"-vecload_block_size %dn",xin->bs);
206:     }
207:   } else {
208:     /* send values */
209:     MPI_Send(xarray,xin->n,MPIU_SCALAR,0,tag,xin->comm);
210:   }
211:   VecRestoreArrayFast(xin,&xarray);
212:   return(0);
213: }

215: #undef __FUNCT__  
217: int VecView_MPI_Draw_LG(Vec xin,PetscViewer viewer)
218: {
219:   int         i,rank,size,N = xin->N,*lens,ierr;
220:   PetscDraw   draw;
221:   PetscReal   *xx,*yy;
222:   PetscDrawLG lg;
223:   PetscTruth  isnull;
224:   PetscScalar *xarray;

227:   PetscViewerDrawGetDraw(viewer,0,&draw);
228:   PetscDrawIsNull(draw,&isnull);
229:   if (isnull) return(0);

231:   VecGetArrayFast(xin,&xarray);
232:   PetscViewerDrawGetDrawLG(viewer,0,&lg);
233:   PetscDrawCheckResizedWindow(draw);
234:   MPI_Comm_rank(xin->comm,&rank);
235:   MPI_Comm_size(xin->comm,&size);
236:   if (!rank) {
237:     PetscDrawLGReset(lg);
238:     PetscMalloc(2*(N+1)*sizeof(PetscReal),&xx);
239:     for (i=0; i<N; i++) {xx[i] = (PetscReal) i;}
240:     yy   = xx + N;
241:     PetscMalloc(size*sizeof(int),&lens);
242:     for (i=0; i<size; i++) {
243:       lens[i] = xin->map->range[i+1] - xin->map->range[i];
244:     }
245: #if !defined(PETSC_USE_COMPLEX)
246:     MPI_Gatherv(xarray,xin->n,MPIU_REAL,yy,lens,xin->map->range,MPIU_REAL,0,xin->comm);
247: #else
248:     {
249:       PetscReal *xr;
250:       PetscMalloc((xin->n+1)*sizeof(PetscReal),&xr);
251:       for (i=0; i<xin->n; i++) {
252:         xr[i] = PetscRealPart(xarray[i]);
253:       }
254:       MPI_Gatherv(xr,xin->n,MPIU_REAL,yy,lens,xin->map->range,MPIU_REAL,0,xin->comm);
255:       PetscFree(xr);
256:     }
257: #endif
258:     PetscFree(lens);
259:     PetscDrawLGAddPoints(lg,N,&xx,&yy);
260:     PetscFree(xx);
261:   } else {
262: #if !defined(PETSC_USE_COMPLEX)
263:     MPI_Gatherv(xarray,xin->n,MPIU_REAL,0,0,0,MPIU_REAL,0,xin->comm);
264: #else
265:     {
266:       PetscReal *xr;
267:       PetscMalloc((xin->n+1)*sizeof(PetscReal),&xr);
268:       for (i=0; i<xin->n; i++) {
269:         xr[i] = PetscRealPart(xarray[i]);
270:       }
271:       MPI_Gatherv(xr,xin->n,MPIU_REAL,0,0,0,MPIU_REAL,0,xin->comm);
272:       PetscFree(xr);
273:     }
274: #endif
275:   }
276:   PetscDrawLGDraw(lg);
277:   PetscDrawSynchronizedFlush(draw);
278:   VecRestoreArrayFast(xin,&xarray);
279:   return(0);
280: }

282: EXTERN_C_BEGIN
283: #undef __FUNCT__  
285: int VecView_MPI_Draw(Vec xin,PetscViewer viewer)
286: {
287:   int           i,rank,size,ierr,start,end,tag = ((PetscObject)viewer)->tag;
288:   MPI_Status    status;
289:   PetscReal     coors[4],ymin,ymax,xmin,xmax,tmp;
290:   PetscDraw     draw;
291:   PetscTruth    isnull;
292:   PetscDrawAxis axis;
293:   PetscScalar   *xarray;

296:   PetscViewerDrawGetDraw(viewer,0,&draw);
297:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

299:   VecGetArrayFast(xin,&xarray);
300:   PetscDrawCheckResizedWindow(draw);
301:   xmin = 1.e20; xmax = -1.e20;
302:   for (i=0; i<xin->n; i++) {
303: #if defined(PETSC_USE_COMPLEX)
304:     if (PetscRealPart(xarray[i]) < xmin) xmin = PetscRealPart(xarray[i]);
305:     if (PetscRealPart(xarray[i]) > xmax) xmax = PetscRealPart(xarray[i]);
306: #else
307:     if (xarray[i] < xmin) xmin = xarray[i];
308:     if (xarray[i] > xmax) xmax = xarray[i];
309: #endif
310:   }
311:   if (xmin + 1.e-10 > xmax) {
312:     xmin -= 1.e-5;
313:     xmax += 1.e-5;
314:   }
315:   MPI_Reduce(&xmin,&ymin,1,MPIU_REAL,MPI_MIN,0,xin->comm);
316:   MPI_Reduce(&xmax,&ymax,1,MPIU_REAL,MPI_MAX,0,xin->comm);
317:   MPI_Comm_size(xin->comm,&size);
318:   MPI_Comm_rank(xin->comm,&rank);
319:   PetscDrawAxisCreate(draw,&axis);
320:   PetscLogObjectParent(draw,axis);
321:   if (!rank) {
322:     PetscDrawClear(draw);
323:     PetscDrawFlush(draw);
324:     PetscDrawAxisSetLimits(axis,0.0,(double)xin->N,ymin,ymax);
325:     PetscDrawAxisDraw(axis);
326:     PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);
327:   }
328:   PetscDrawAxisDestroy(axis);
329:   MPI_Bcast(coors,4,MPIU_REAL,0,xin->comm);
330:   if (rank) {PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);}
331:   /* draw local part of vector */
332:   VecGetOwnershipRange(xin,&start,&end);
333:   if (rank < size-1) { /*send value to right */
334:     MPI_Send(&xarray[xin->n-1],1,MPIU_REAL,rank+1,tag,xin->comm);
335:   }
336:   for (i=1; i<xin->n; i++) {
337: #if !defined(PETSC_USE_COMPLEX)
338:     PetscDrawLine(draw,(PetscReal)(i-1+start),xarray[i-1],(PetscReal)(i+start),
339:                    xarray[i],PETSC_DRAW_RED);
340: #else
341:     PetscDrawLine(draw,(PetscReal)(i-1+start),PetscRealPart(xarray[i-1]),(PetscReal)(i+start),
342:                    PetscRealPart(xarray[i]),PETSC_DRAW_RED);
343: #endif
344:   }
345:   if (rank) { /* receive value from right */
346:     MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag,xin->comm,&status);
347: #if !defined(PETSC_USE_COMPLEX)
348:     PetscDrawLine(draw,(PetscReal)start-1,tmp,(PetscReal)start,xarray[0],PETSC_DRAW_RED);
349: #else
350:     PetscDrawLine(draw,(PetscReal)start-1,tmp,(PetscReal)start,PetscRealPart(xarray[0]),PETSC_DRAW_RED);
351: #endif
352:   }
353:   PetscDrawSynchronizedFlush(draw);
354:   PetscDrawPause(draw);
355:   VecRestoreArrayFast(xin,&xarray);
356:   return(0);
357: }
358: EXTERN_C_END

360: #undef __FUNCT__  
362: int VecView_MPI_Socket(Vec xin,PetscViewer viewer)
363: {
364:   int         i,rank,size,N = xin->N,*lens,ierr;
365:   PetscScalar *xx,*xarray;

368:   VecGetArrayFast(xin,&xarray);
369:   MPI_Comm_rank(xin->comm,&rank);
370:   MPI_Comm_size(xin->comm,&size);
371:   if (!rank) {
372:     PetscMalloc((N+1)*sizeof(PetscScalar),&xx);
373:     PetscMalloc(size*sizeof(int),&lens);
374:     for (i=0; i<size; i++) {
375:       lens[i] = xin->map->range[i+1] - xin->map->range[i];
376:     }
377:     MPI_Gatherv(xarray,xin->n,MPIU_SCALAR,xx,lens,xin->map->range,MPIU_SCALAR,0,xin->comm);
378:     PetscFree(lens);
379:     PetscViewerSocketPutScalar(viewer,N,1,xx);
380:     PetscFree(xx);
381:   } else {
382:     MPI_Gatherv(xarray,xin->n,MPIU_REAL,0,0,0,MPIU_REAL,0,xin->comm);
383:   }
384:   VecRestoreArrayFast(xin,&xarray);
385:   return(0);
386: }

388: #undef __FUNCT__  
390: int VecView_MPI(Vec xin,PetscViewer viewer)
391: {
392:   int        ierr;
393:   PetscTruth isascii,issocket,isbinary,isdraw,ismathematica;

396:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
397:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
398:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
399:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
400:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_MATHEMATICA,&ismathematica);
401:   if (isascii){
402:     VecView_MPI_ASCII(xin,viewer);
403:   } else if (issocket) {
404:     VecView_MPI_Socket(xin,viewer);
405:   } else if (isbinary) {
406:     VecView_MPI_Binary(xin,viewer);
407:   } else if (isdraw) {
408:     PetscViewerFormat format;

410:     PetscViewerGetFormat(viewer,&format);
411:     if (format == PETSC_VIEWER_DRAW_LG) {
412:       VecView_MPI_Draw_LG(xin,viewer);
413:     } else {
414:       VecView_MPI_Draw(xin,viewer);
415:     }
416:   } else if (ismathematica) {
417:     PetscViewerMathematicaPutVector(viewer,xin);
418:   } else {
419:     SETERRQ1(1,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
420:   }
421:   return(0);
422: }

424: #undef __FUNCT__  
426: int VecGetSize_MPI(Vec xin,int *N)
427: {
429:   *N = xin->N;
430:   return(0);
431: }

433: #undef __FUNCT__  
435: int VecSetValues_MPI(Vec xin,int ni,const int ix[],const PetscScalar y[],InsertMode addv)
436: {
437:   int           rank = xin->stash.rank, *owners = xin->map->range,start = owners[rank];
438:   int           end = owners[rank+1],i,row,ierr;
439:   PetscScalar   *xx;

442: #if defined(PETSC_USE_BOPT_g)
443:   if (xin->stash.insertmode == INSERT_VALUES && addv == ADD_VALUES) {
444:    SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You have already inserted values; you cannot now add");
445:   } else if (xin->stash.insertmode == ADD_VALUES && addv == INSERT_VALUES) {
446:    SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You have already added values; you cannot now insert");
447:   }
448: #endif
449:   VecGetArrayFast(xin,&xx);
450:   xin->stash.insertmode = addv;

452:   if (addv == INSERT_VALUES) {
453:     for (i=0; i<ni; i++) {
454:       if ((row = ix[i]) >= start && row < end) {
455:         xx[row-start] = y[i];
456:       } else if (!xin->stash.donotstash) {
457:         if (ix[i] < 0) continue;
458: #if defined(PETSC_USE_BOPT_g)
459:         if (ix[i] >= xin->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",ix[i],xin->N);
460: #endif
461:         VecStashValue_Private(&xin->stash,row,y[i]);
462:       }
463:     }
464:   } else {
465:     for (i=0; i<ni; i++) {
466:       if ((row = ix[i]) >= start && row < end) {
467:         xx[row-start] += y[i];
468:       } else if (!xin->stash.donotstash) {
469:         if (ix[i] < 0) continue;
470: #if defined(PETSC_USE_BOPT_g)
471:         if (ix[i] > xin->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d maximum %d",ix[i],xin->N);
472: #endif        
473:         VecStashValue_Private(&xin->stash,row,y[i]);
474:       }
475:     }
476:   }
477:   VecRestoreArrayFast(xin,&xx);
478:   return(0);
479: }

481: #undef __FUNCT__  
483: int VecSetValuesBlocked_MPI(Vec xin,int ni,const int ix[],const PetscScalar yin[],InsertMode addv)
484: {
485:   int           rank = xin->stash.rank,*owners = xin->map->range,start = owners[rank];
486:   int           end = owners[rank+1],i,row,bs = xin->bs,j,ierr;
487:   PetscScalar   *xx,*y = (PetscScalar*)yin;

490:   VecGetArrayFast(xin,&xx);
491: #if defined(PETSC_USE_BOPT_g)
492:   if (xin->stash.insertmode == INSERT_VALUES && addv == ADD_VALUES) {
493:    SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You have already inserted values; you cannot now add");
494:   }
495:   else if (xin->stash.insertmode == ADD_VALUES && addv == INSERT_VALUES) {
496:    SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You have already added values; you cannot now insert");
497:   }
498: #endif
499:   xin->stash.insertmode = addv;

501:   if (addv == INSERT_VALUES) {
502:     for (i=0; i<ni; i++) {
503:       if ((row = bs*ix[i]) >= start && row < end) {
504:         for (j=0; j<bs; j++) {
505:           xx[row-start+j] = y[j];
506:         }
507:       } else if (!xin->stash.donotstash) {
508:         if (ix[i] < 0) continue;
509: #if defined(PETSC_USE_BOPT_g)
510:         if (ix[i] >= xin->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d max %d",ix[i],xin->N);
511: #endif
512:         VecStashValuesBlocked_Private(&xin->bstash,ix[i],y);
513:       }
514:       y += bs;
515:     }
516:   } else {
517:     for (i=0; i<ni; i++) {
518:       if ((row = bs*ix[i]) >= start && row < end) {
519:         for (j=0; j<bs; j++) {
520:           xx[row-start+j] += y[j];
521:         }
522:       } else if (!xin->stash.donotstash) {
523:         if (ix[i] < 0) continue;
524: #if defined(PETSC_USE_BOPT_g)
525:         if (ix[i] > xin->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %d max %d",ix[i],xin->N);
526: #endif
527:         VecStashValuesBlocked_Private(&xin->bstash,ix[i],y);
528:       }
529:       y += bs;
530:     }
531:   }
532:   VecRestoreArrayFast(xin,&xx);
533:   return(0);
534: }

536: /*
537:    Since nsends or nreceives may be zero we add 1 in certain mallocs
538: to make sure we never malloc an empty one.      
539: */
540: #undef __FUNCT__  
542: int VecAssemblyBegin_MPI(Vec xin)
543: {
544:   int         *owners = xin->map->range,*bowners,ierr,size,i,bs,nstash,reallocs;
545:   InsertMode  addv;
546:   MPI_Comm    comm = xin->comm;

549:   if (xin->stash.donotstash) {
550:     return(0);
551:   }

553:   MPI_Allreduce(&xin->stash.insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
554:   if (addv == (ADD_VALUES|INSERT_VALUES)) {
555:     SETERRQ(PETSC_ERR_ARG_NOTSAMETYPE,"Some processors inserted values while others added");
556:   }
557:   xin->stash.insertmode = addv; /* in case this processor had no cache */
558: 
559:   bs = xin->bs;
560:   MPI_Comm_size(xin->comm,&size);
561:   if (!xin->bstash.bowners && xin->bs != -1) {
562:     PetscMalloc((size+1)*sizeof(int),&bowners);
563:     for (i=0; i<size+1; i++){ bowners[i] = owners[i]/bs;}
564:     xin->bstash.bowners = bowners;
565:   } else {
566:     bowners = xin->bstash.bowners;
567:   }
568:   VecStashScatterBegin_Private(&xin->stash,owners);
569:   VecStashScatterBegin_Private(&xin->bstash,bowners);
570:   VecStashGetInfo_Private(&xin->stash,&nstash,&reallocs);
571:   PetscLogInfo(0,"VecAssemblyBegin_MPI:Stash has %d entries, uses %d mallocs.n",nstash,reallocs);
572:   VecStashGetInfo_Private(&xin->bstash,&nstash,&reallocs);
573:   PetscLogInfo(0,"VecAssemblyBegin_MPI:Block-Stash has %d entries, uses %d mallocs.n",nstash,reallocs);

575:   return(0);
576: }

578: #undef __FUNCT__  
580: int VecAssemblyEnd_MPI(Vec vec)
581: {
582:   int         ierr,base,i,j,n,*row,flg,bs;
583:   PetscScalar *val,*vv,*array,*xarray;

586:   if (!vec->stash.donotstash) {
587:     VecGetArrayFast(vec,&xarray);
588:     base = vec->map->range[vec->stash.rank];
589:     bs   = vec->bs;

591:     /* Process the stash */
592:     while (1) {
593:       VecStashScatterGetMesg_Private(&vec->stash,&n,&row,&val,&flg);
594:       if (!flg) break;
595:       if (vec->stash.insertmode == ADD_VALUES) {
596:         for (i=0; i<n; i++) { xarray[row[i] - base] += val[i]; }
597:       } else if (vec->stash.insertmode == INSERT_VALUES) {
598:         for (i=0; i<n; i++) { xarray[row[i] - base] = val[i]; }
599:       } else {
600:         SETERRQ(PETSC_ERR_ARG_CORRUPT,"Insert mode is not set correctly; corrupted vector");
601:       }
602:     }
603:     VecStashScatterEnd_Private(&vec->stash);

605:     /* now process the block-stash */
606:     while (1) {
607:       VecStashScatterGetMesg_Private(&vec->bstash,&n,&row,&val,&flg);
608:       if (!flg) break;
609:       for (i=0; i<n; i++) {
610:         array = xarray+row[i]*bs-base;
611:         vv    = val+i*bs;
612:         if (vec->stash.insertmode == ADD_VALUES) {
613:           for (j=0; j<bs; j++) { array[j] += vv[j];}
614:         } else if (vec->stash.insertmode == INSERT_VALUES) {
615:           for (j=0; j<bs; j++) { array[j] = vv[j]; }
616:         } else {
617:           SETERRQ(PETSC_ERR_ARG_CORRUPT,"Insert mode is not set correctly; corrupted vector");
618:         }
619:       }
620:     }
621:     VecStashScatterEnd_Private(&vec->bstash);
622:     VecRestoreArrayFast(vec,&xarray);
623:   }
624:   vec->stash.insertmode = NOT_SET_VALUES;
625:   return(0);
626: }