Actual source code: pdvec.c

petsc-dev 2014-02-02
Report Typos and Errors
  2: /*
  3:      Code for some of the parallel vector primatives.
  4: */
  5: #include <../src/vec/vec/impls/mpi/pvecimpl.h>   /*I  "petscvec.h"   I*/
  6: #include <petscviewerhdf5.h>

 10: PetscErrorCode VecDestroy_MPI(Vec v)
 11: {
 12:   Vec_MPI        *x = (Vec_MPI*)v->data;

 16: #if defined(PETSC_USE_LOG)
 17:   PetscLogObjectState((PetscObject)v,"Length=%D",v->map->N);
 18: #endif
 19:   if (!x) return(0);
 20:   PetscFree(x->array_allocated);

 22:   /* Destroy local representation of vector if it exists */
 23:   if (x->localrep) {
 24:     VecDestroy(&x->localrep);
 25:     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(v->data);
 31:   return(0);
 32: }

 36: PetscErrorCode VecView_MPI_ASCII(Vec xin,PetscViewer viewer)
 37: {
 38:   PetscErrorCode    ierr;
 39:   PetscInt          i,work = xin->map->n,cnt,len,nLen;
 40:   PetscMPIInt       j,n = 0,size,rank,tag = ((PetscObject)viewer)->tag;
 41:   MPI_Status        status;
 42:   PetscScalar       *values;
 43:   const PetscScalar *xarray;
 44:   const char        *name;
 45:   PetscViewerFormat format;

 48:   VecGetArrayRead(xin,&xarray);
 49:   /* determine maximum message to arrive */
 50:   MPI_Comm_rank(PetscObjectComm((PetscObject)xin),&rank);
 51:   MPI_Reduce(&work,&len,1,MPIU_INT,MPI_MAX,0,PetscObjectComm((PetscObject)xin));
 52:   MPI_Comm_size(PetscObjectComm((PetscObject)xin),&size);

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

 97:     } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
 98:       for (i=0; i<xin->map->n; i++) {
 99: #if defined(PETSC_USE_COMPLEX)
100:         PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e\n",(double)PetscRealPart(xarray[i]),(double)PetscImaginaryPart(xarray[i]));
101: #else
102:         PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)xarray[i]);
103: #endif
104:       }
105:       /* receive and print messages */
106:       for (j=1; j<size; j++) {
107:         MPI_Recv(values,(PetscMPIInt)len,MPIU_SCALAR,j,tag,PetscObjectComm((PetscObject)xin),&status);
108:         MPI_Get_count(&status,MPIU_SCALAR,&n);
109:         for (i=0; i<n; i++) {
110: #if defined(PETSC_USE_COMPLEX)
111:           PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e\n",(double)PetscRealPart(values[i]),(double)PetscImaginaryPart(values[i]));
112: #else
113:           PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)values[i]);
114: #endif
115:         }
116:       }
117:     } else if (format == PETSC_VIEWER_ASCII_VTK || format == PETSC_VIEWER_ASCII_VTK_CELL) {
118:       /*
119:         state 0: No header has been output
120:         state 1: Only POINT_DATA has been output
121:         state 2: Only CELL_DATA has been output
122:         state 3: Output both, POINT_DATA last
123:         state 4: Output both, CELL_DATA last
124:       */
125:       static PetscInt stateId     = -1;
126:       int             outputState = 0;
127:       int             doOutput    = 0;
128:       PetscBool       hasState;
129:       PetscInt        bs, b;

131:       if (stateId < 0) {
132:         PetscObjectComposedDataRegister(&stateId);
133:       }
134:       PetscObjectComposedDataGetInt((PetscObject) viewer, stateId, outputState, hasState);
135:       if (!hasState) outputState = 0;

137:       PetscObjectGetName((PetscObject) xin, &name);
138:       VecGetLocalSize(xin, &nLen);
139:       PetscMPIIntCast(nLen,&n);
140:       VecGetBlockSize(xin, &bs);
141:       if (format == PETSC_VIEWER_ASCII_VTK) {
142:         if (outputState == 0) {
143:           outputState = 1;
144:           doOutput    = 1;
145:         } else if (outputState == 1) doOutput = 0;
146:         else if (outputState == 2) {
147:           outputState = 3;
148:           doOutput    = 1;
149:         } else if (outputState == 3) doOutput = 0;
150:         else if (outputState == 4) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Tried to output POINT_DATA again after intervening CELL_DATA");

152:         if (doOutput) {
153:           PetscViewerASCIIPrintf(viewer, "POINT_DATA %d\n", xin->map->N/bs);
154:         }
155:       } else {
156:         if (outputState == 0) {
157:           outputState = 2;
158:           doOutput    = 1;
159:         } else if (outputState == 1) {
160:           outputState = 4;
161:           doOutput    = 1;
162:         } else if (outputState == 2) doOutput = 0;
163:         else if (outputState == 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Tried to output CELL_DATA again after intervening POINT_DATA");
164:         else if (outputState == 4) doOutput = 0;

166:         if (doOutput) {
167:           PetscViewerASCIIPrintf(viewer, "CELL_DATA %d\n", xin->map->N/bs);
168:         }
169:       }
170:       PetscObjectComposedDataSetInt((PetscObject) viewer, stateId, outputState);
171:       if (name) {
172:         if (bs == 3) {
173:           PetscViewerASCIIPrintf(viewer, "VECTORS %s double\n", name);
174:         } else {
175:           PetscViewerASCIIPrintf(viewer, "SCALARS %s double %d\n", name, bs);
176:         }
177:       } else {
178:         PetscViewerASCIIPrintf(viewer, "SCALARS scalars double %d\n", bs);
179:       }
180:       if (bs != 3) {
181:         PetscViewerASCIIPrintf(viewer, "LOOKUP_TABLE default\n");
182:       }
183:       for (i=0; i<n/bs; i++) {
184:         for (b=0; b<bs; b++) {
185:           if (b > 0) {
186:             PetscViewerASCIIPrintf(viewer," ");
187:           }
188:           PetscViewerASCIIPrintf(viewer,"%g",(double)PetscRealPart(xarray[i*bs+b]));
189:         }
190:         PetscViewerASCIIPrintf(viewer,"\n");
191:       }
192:       for (j=1; j<size; j++) {
193:         MPI_Recv(values,(PetscMPIInt)len,MPIU_SCALAR,j,tag,PetscObjectComm((PetscObject)xin),&status);
194:         MPI_Get_count(&status,MPIU_SCALAR,&n);
195:         for (i=0; i<n/bs; i++) {
196:           for (b=0; b<bs; b++) {
197:             if (b > 0) {
198:               PetscViewerASCIIPrintf(viewer," ");
199:             }
200:             PetscViewerASCIIPrintf(viewer,"%g",(double)PetscRealPart(values[i*bs+b]));
201:           }
202:           PetscViewerASCIIPrintf(viewer,"\n");
203:         }
204:       }
205:     } else if (format == PETSC_VIEWER_ASCII_VTK_COORDS) {
206:       PetscInt bs, b;

208:       VecGetLocalSize(xin, &nLen);
209:       PetscMPIIntCast(nLen,&n);
210:       VecGetBlockSize(xin, &bs);
211:       if ((bs < 1) || (bs > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %d", bs);

213:       for (i=0; i<n/bs; i++) {
214:         for (b=0; b<bs; b++) {
215:           if (b > 0) {
216:             PetscViewerASCIIPrintf(viewer," ");
217:           }
218:           PetscViewerASCIIPrintf(viewer,"%g",(double)PetscRealPart(xarray[i*bs+b]));
219:         }
220:         for (b=bs; b<3; b++) {
221:           PetscViewerASCIIPrintf(viewer," 0.0");
222:         }
223:         PetscViewerASCIIPrintf(viewer,"\n");
224:       }
225:       for (j=1; j<size; j++) {
226:         MPI_Recv(values,(PetscMPIInt)len,MPIU_SCALAR,j,tag,PetscObjectComm((PetscObject)xin),&status);
227:         MPI_Get_count(&status,MPIU_SCALAR,&n);
228:         for (i=0; i<n/bs; i++) {
229:           for (b=0; b<bs; b++) {
230:             if (b > 0) {
231:               PetscViewerASCIIPrintf(viewer," ");
232:             }
233:             PetscViewerASCIIPrintf(viewer,"%g",(double)PetscRealPart(values[i*bs+b]));
234:           }
235:           for (b=bs; b<3; b++) {
236:             PetscViewerASCIIPrintf(viewer," 0.0");
237:           }
238:           PetscViewerASCIIPrintf(viewer,"\n");
239:         }
240:       }
241:     } else if (format == PETSC_VIEWER_ASCII_PCICE) {
242:       PetscInt bs, b, vertexCount = 1;

244:       VecGetLocalSize(xin, &nLen);
245:       PetscMPIIntCast(nLen,&n);
246:       VecGetBlockSize(xin, &bs);
247:       if ((bs < 1) || (bs > 3)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "PCICE can only handle up to 3D objects, but vector dimension is %d", bs);

249:       PetscViewerASCIIPrintf(viewer,"%D\n", xin->map->N/bs);
250:       for (i=0; i<n/bs; i++) {
251:         PetscViewerASCIIPrintf(viewer,"%7D   ", vertexCount++);
252:         for (b=0; b<bs; b++) {
253:           if (b > 0) {
254:             PetscViewerASCIIPrintf(viewer," ");
255:           }
256: #if !defined(PETSC_USE_COMPLEX)
257:           PetscViewerASCIIPrintf(viewer,"% 12.5E",(double)xarray[i*bs+b]);
258: #endif
259:         }
260:         PetscViewerASCIIPrintf(viewer,"\n");
261:       }
262:       for (j=1; j<size; j++) {
263:         MPI_Recv(values,(PetscMPIInt)len,MPIU_SCALAR,j,tag,PetscObjectComm((PetscObject)xin),&status);
264:         MPI_Get_count(&status,MPIU_SCALAR,&n);
265:         for (i=0; i<n/bs; i++) {
266:           PetscViewerASCIIPrintf(viewer,"%7D   ", vertexCount++);
267:           for (b=0; b<bs; b++) {
268:             if (b > 0) {
269:               PetscViewerASCIIPrintf(viewer," ");
270:             }
271: #if !defined(PETSC_USE_COMPLEX)
272:             PetscViewerASCIIPrintf(viewer,"% 12.5E",(double)values[i*bs+b]);
273: #endif
274:           }
275:           PetscViewerASCIIPrintf(viewer,"\n");
276:         }
277:       }
278:     } else {
279:       PetscObjectPrintClassNamePrefixType((PetscObject)xin,viewer);
280:       if (format != PETSC_VIEWER_ASCII_COMMON) {PetscViewerASCIIPrintf(viewer,"Process [%d]\n",rank);}
281:       cnt = 0;
282:       for (i=0; i<xin->map->n; i++) {
283:         if (format == PETSC_VIEWER_ASCII_INDEX) {
284:           PetscViewerASCIIPrintf(viewer,"%D: ",cnt++);
285:         }
286: #if defined(PETSC_USE_COMPLEX)
287:         if (PetscImaginaryPart(xarray[i]) > 0.0) {
288:           PetscViewerASCIIPrintf(viewer,"%g + %g i\n",(double)PetscRealPart(xarray[i]),(double)PetscImaginaryPart(xarray[i]));
289:         } else if (PetscImaginaryPart(xarray[i]) < 0.0) {
290:           PetscViewerASCIIPrintf(viewer,"%g - %g i\n",(double)PetscRealPart(xarray[i]),-(double)PetscImaginaryPart(xarray[i]));
291:         } else {
292:           PetscViewerASCIIPrintf(viewer,"%g\n",(double)PetscRealPart(xarray[i]));
293:         }
294: #else
295:         PetscViewerASCIIPrintf(viewer,"%g\n",(double)xarray[i]);
296: #endif
297:       }
298:       /* receive and print messages */
299:       for (j=1; j<size; j++) {
300:         MPI_Recv(values,(PetscMPIInt)len,MPIU_SCALAR,j,tag,PetscObjectComm((PetscObject)xin),&status);
301:         MPI_Get_count(&status,MPIU_SCALAR,&n);
302:         if (format != PETSC_VIEWER_ASCII_COMMON) {
303:           PetscViewerASCIIPrintf(viewer,"Process [%d]\n",j);
304:         }
305:         for (i=0; i<n; i++) {
306:           if (format == PETSC_VIEWER_ASCII_INDEX) {
307:             PetscViewerASCIIPrintf(viewer,"%D: ",cnt++);
308:           }
309: #if defined(PETSC_USE_COMPLEX)
310:           if (PetscImaginaryPart(values[i]) > 0.0) {
311:             PetscViewerASCIIPrintf(viewer,"%g + %g i\n",(double)PetscRealPart(values[i]),(double)PetscImaginaryPart(values[i]));
312:           } else if (PetscImaginaryPart(values[i]) < 0.0) {
313:             PetscViewerASCIIPrintf(viewer,"%g - %g i\n",(double)PetscRealPart(values[i]),-(double)PetscImaginaryPart(values[i]));
314:           } else {
315:             PetscViewerASCIIPrintf(viewer,"%g\n",(double)PetscRealPart(values[i]));
316:           }
317: #else
318:           PetscViewerASCIIPrintf(viewer,"%g\n",(double)values[i]);
319: #endif
320:         }
321:       }
322:     }
323:     PetscFree(values);
324:   } else {
325:     PetscViewerGetFormat(viewer,&format);
326:     if (format == PETSC_VIEWER_ASCII_MATLAB) {
327:       /* this may be a collective operation so make sure everyone calls it */
328:       PetscObjectGetName((PetscObject)xin,&name);
329:     }
330:     /* send values */
331:     MPI_Send((void*)xarray,xin->map->n,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)xin));
332:   }
333:   PetscViewerFlush(viewer);
334:   VecRestoreArrayRead(xin,&xarray);
335:   return(0);
336: }

340: PetscErrorCode VecView_MPI_Binary(Vec xin,PetscViewer viewer)
341: {
342:   PetscErrorCode    ierr;
343:   PetscMPIInt       rank,size,mesgsize,tag = ((PetscObject)viewer)->tag, mesglen;
344:   PetscInt          len,n = xin->map->n,j,tr[2];
345:   int               fdes;
346:   MPI_Status        status;
347:   PetscScalar       *values;
348:   const PetscScalar *xarray;
349:   FILE              *file;
350: #if defined(PETSC_HAVE_MPIIO)
351:   PetscBool         isMPIIO;
352: #endif
353:   PetscBool         skipHeader;
354:   PetscInt          message_count,flowcontrolcount;
355:   PetscViewerFormat format;

358:   VecGetArrayRead(xin,&xarray);
359:   PetscViewerBinaryGetDescriptor(viewer,&fdes);
360:   PetscViewerBinaryGetSkipHeader(viewer,&skipHeader);

362:   /* determine maximum message to arrive */
363:   MPI_Comm_rank(PetscObjectComm((PetscObject)xin),&rank);
364:   MPI_Comm_size(PetscObjectComm((PetscObject)xin),&size);

366:   if (!skipHeader) {
367:     tr[0] = VEC_FILE_CLASSID;
368:     tr[1] = xin->map->N;
369:     PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_FALSE);
370:   }

372: #if defined(PETSC_HAVE_MPIIO)
373:   PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);
374:   if (!isMPIIO) {
375: #endif
376:     PetscViewerFlowControlStart(viewer,&message_count,&flowcontrolcount);
377:     if (!rank) {
378:       PetscBinaryWrite(fdes,(void*)xarray,xin->map->n,PETSC_SCALAR,PETSC_FALSE);

380:       len = 0;
381:       for (j=1; j<size; j++) len = PetscMax(len,xin->map->range[j+1]-xin->map->range[j]);
382:       PetscMalloc1(len,&values);
383:       PetscMPIIntCast(len,&mesgsize);
384:       /* receive and save messages */
385:       for (j=1; j<size; j++) {
386:         PetscViewerFlowControlStepMaster(viewer,j,&message_count,flowcontrolcount);
387:         MPI_Recv(values,mesgsize,MPIU_SCALAR,j,tag,PetscObjectComm((PetscObject)xin),&status);
388:         MPI_Get_count(&status,MPIU_SCALAR,&mesglen);
389:         n    = (PetscInt)mesglen;
390:         PetscBinaryWrite(fdes,values,n,PETSC_SCALAR,PETSC_TRUE);
391:       }
392:       PetscViewerFlowControlEndMaster(viewer,&message_count);
393:       PetscFree(values);
394:     } else {
395:       PetscViewerFlowControlStepWorker(viewer,rank,&message_count);
396:       PetscMPIIntCast(xin->map->n,&mesgsize);
397:       MPI_Send((void*)xarray,mesgsize,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)xin));
398:       PetscViewerFlowControlEndWorker(viewer,&message_count);
399:     }
400:     PetscViewerGetFormat(viewer,&format);
401:     if (format == PETSC_VIEWER_BINARY_MATLAB) {
402:       MPI_Comm   comm;
403:       FILE       *info;
404:       const char *name;

406:       PetscObjectGetName((PetscObject)xin,&name);
407:       PetscObjectGetComm((PetscObject)viewer,&comm);
408:       PetscViewerBinaryGetInfoPointer(viewer,&info);
409:       PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");
410:       PetscFPrintf(comm,info,"#$$ Set.%s = PetscBinaryRead(fd);\n",name);
411:       PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");
412:     }
413: #if defined(PETSC_HAVE_MPIIO)
414:   } else {
415:     MPI_Offset   off;
416:     MPI_File     mfdes;
417:     PetscMPIInt  lsize;

419:     PetscMPIIntCast(n,&lsize);
420:     PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);
421:     PetscViewerBinaryGetMPIIOOffset(viewer,&off);
422:     off += xin->map->rstart*sizeof(PetscScalar); /* off is MPI_Offset, not PetscMPIInt */
423:     MPI_File_set_view(mfdes,off,MPIU_SCALAR,MPIU_SCALAR,(char*)"native",MPI_INFO_NULL);
424:     MPIU_File_write_all(mfdes,(void*)xarray,lsize,MPIU_SCALAR,MPI_STATUS_IGNORE);
425:     PetscViewerBinaryAddMPIIOOffset(viewer,xin->map->N*sizeof(PetscScalar));
426:   }
427: #endif

429:   VecRestoreArrayRead(xin,&xarray);
430:   if (!rank) {
431:     PetscViewerBinaryGetInfoPointer(viewer,&file);
432:     if (file) {
433:       if (((PetscObject)xin)->prefix) {
434:         PetscFPrintf(PETSC_COMM_SELF,file,"-%svecload_block_size %D\n",((PetscObject)xin)->prefix,xin->map->bs);
435:       } else {
436:         PetscFPrintf(PETSC_COMM_SELF,file,"-vecload_block_size %D\n",xin->map->bs);
437:       }
438:     }
439:   }
440:   return(0);
441: }

443: #include <petscdraw.h>
446: PetscErrorCode VecView_MPI_Draw_LG(Vec xin,PetscViewer viewer)
447: {
448:   PetscDraw      draw;
449:   PetscBool      isnull;

452: #if defined(PETSC_USE_64BIT_INDICES)
454:   PetscViewerDrawGetDraw(viewer,0,&draw);
455:   PetscDrawIsNull(draw,&isnull);
456:   if (isnull) return(0);
457:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported with 64 bit integers");
458: #else
459:   PetscMPIInt       size,rank;
460:   int               i,N = xin->map->N,*lens;
461:   PetscReal         *xx,*yy;
462:   PetscDrawLG       lg;
463:   const PetscScalar *xarray;

466:   PetscViewerDrawGetDraw(viewer,0,&draw);
467:   PetscDrawIsNull(draw,&isnull);
468:   if (isnull) return(0);

470:   VecGetArrayRead(xin,&xarray);
471:   PetscViewerDrawGetDrawLG(viewer,0,&lg);
472:   PetscDrawCheckResizedWindow(draw);
473:   MPI_Comm_rank(PetscObjectComm((PetscObject)xin),&rank);
474:   MPI_Comm_size(PetscObjectComm((PetscObject)xin),&size);
475:   if (!rank) {
476:     PetscDrawLGReset(lg);
477:     PetscMalloc2(N,&xx,N,&yy);
478:     for (i=0; i<N; i++) xx[i] = (PetscReal) i;
479:     PetscMalloc1(size,&lens);
480:     for (i=0; i<size; i++) lens[i] = xin->map->range[i+1] - xin->map->range[i];

482: #if !defined(PETSC_USE_COMPLEX)
483:     MPI_Gatherv((void*)xarray,xin->map->n,MPIU_REAL,yy,lens,xin->map->range,MPIU_REAL,0,PetscObjectComm((PetscObject)xin));
484: #else
485:     {
486:       PetscReal *xr;
487:       PetscMalloc1((xin->map->n+1),&xr);
488:       for (i=0; i<xin->map->n; i++) xr[i] = PetscRealPart(xarray[i]);

490:       MPI_Gatherv(xr,xin->map->n,MPIU_REAL,yy,lens,xin->map->range,MPIU_REAL,0,PetscObjectComm((PetscObject)xin));
491:       PetscFree(xr);
492:     }
493: #endif
494:     PetscFree(lens);
495:     PetscDrawLGAddPoints(lg,N,&xx,&yy);
496:     PetscFree2(xx,yy);
497:   } else {
498: #if !defined(PETSC_USE_COMPLEX)
499:     MPI_Gatherv((void*)xarray,xin->map->n,MPIU_REAL,0,0,0,MPIU_REAL,0,PetscObjectComm((PetscObject)xin));
500: #else
501:     {
502:       PetscReal *xr;
503:       PetscMalloc1((xin->map->n+1),&xr);
504:       for (i=0; i<xin->map->n; i++) xr[i] = PetscRealPart(xarray[i]);

506:       MPI_Gatherv(xr,xin->map->n,MPIU_REAL,0,0,0,MPIU_REAL,0,PetscObjectComm((PetscObject)xin));
507:       PetscFree(xr);
508:     }
509: #endif
510:   }
511:   PetscDrawLGDraw(lg);
512:   PetscDrawSynchronizedFlush(draw);
513:   VecRestoreArrayRead(xin,&xarray);
514: #endif
515:   return(0);
516: }

520: PetscErrorCode  VecView_MPI_Draw(Vec xin,PetscViewer viewer)
521: {
522:   PetscErrorCode    ierr;
523:   PetscMPIInt       rank,size,tag = ((PetscObject)viewer)->tag;
524:   PetscInt          i,start,end;
525:   MPI_Status        status;
526:   PetscReal         coors[4],ymin,ymax,xmin,xmax,tmp = 0.0;
527:   PetscDraw         draw;
528:   PetscBool         isnull;
529:   PetscDrawAxis     axis;
530:   const PetscScalar *xarray;

533:   PetscViewerDrawGetDraw(viewer,0,&draw);
534:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

536:   VecGetArrayRead(xin,&xarray);
537:   PetscDrawCheckResizedWindow(draw);
538:   xmin = 1.e20; xmax = -1.e20;
539:   for (i=0; i<xin->map->n; i++) {
540: #if defined(PETSC_USE_COMPLEX)
541:     if (PetscRealPart(xarray[i]) < xmin) xmin = PetscRealPart(xarray[i]);
542:     if (PetscRealPart(xarray[i]) > xmax) xmax = PetscRealPart(xarray[i]);
543: #else
544:     if (xarray[i] < xmin) xmin = xarray[i];
545:     if (xarray[i] > xmax) xmax = xarray[i];
546: #endif
547:   }
548:   if (xmin + 1.e-10 > xmax) {
549:     xmin -= 1.e-5;
550:     xmax += 1.e-5;
551:   }
552:   MPI_Reduce(&xmin,&ymin,1,MPIU_REAL,MPIU_MIN,0,PetscObjectComm((PetscObject)xin));
553:   MPI_Reduce(&xmax,&ymax,1,MPIU_REAL,MPIU_MAX,0,PetscObjectComm((PetscObject)xin));
554:   MPI_Comm_size(PetscObjectComm((PetscObject)xin),&size);
555:   MPI_Comm_rank(PetscObjectComm((PetscObject)xin),&rank);
556:   PetscDrawAxisCreate(draw,&axis);
557:   PetscLogObjectParent((PetscObject)draw,(PetscObject)axis);
558:   if (!rank) {
559:     PetscDrawClear(draw);
560:     PetscDrawFlush(draw);
561:     PetscDrawAxisSetLimits(axis,0.0,(double)xin->map->N,ymin,ymax);
562:     PetscDrawAxisDraw(axis);
563:     PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);
564:   }
565:   PetscDrawAxisDestroy(&axis);
566:   MPI_Bcast(coors,4,MPIU_REAL,0,PetscObjectComm((PetscObject)xin));
567:   if (rank) {PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);}
568:   /* draw local part of vector */
569:   VecGetOwnershipRange(xin,&start,&end);
570:   if (rank < size-1) { /*send value to right */
571:     MPI_Send((void*)&xarray[xin->map->n-1],1,MPIU_REAL,rank+1,tag,PetscObjectComm((PetscObject)xin));
572:   }
573:   for (i=1; i<xin->map->n; i++) {
574: #if !defined(PETSC_USE_COMPLEX)
575:     PetscDrawLine(draw,(PetscReal)(i-1+start),xarray[i-1],(PetscReal)(i+start),xarray[i],PETSC_DRAW_RED);
576: #else
577:     PetscDrawLine(draw,(PetscReal)(i-1+start),PetscRealPart(xarray[i-1]),(PetscReal)(i+start),PetscRealPart(xarray[i]),PETSC_DRAW_RED);
578: #endif
579:   }
580:   if (rank) { /* receive value from right */
581:     MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag,PetscObjectComm((PetscObject)xin),&status);
582: #if !defined(PETSC_USE_COMPLEX)
583:     PetscDrawLine(draw,(PetscReal)start-1,tmp,(PetscReal)start,xarray[0],PETSC_DRAW_RED);
584: #else
585:     PetscDrawLine(draw,(PetscReal)start-1,tmp,(PetscReal)start,PetscRealPart(xarray[0]),PETSC_DRAW_RED);
586: #endif
587:   }
588:   PetscDrawSynchronizedFlush(draw);
589:   PetscDrawPause(draw);
590:   VecRestoreArrayRead(xin,&xarray);
591:   return(0);
592: }

594: #if defined(PETSC_HAVE_MATLAB_ENGINE)
597: PetscErrorCode VecView_MPI_Matlab(Vec xin,PetscViewer viewer)
598: {
599:   PetscErrorCode    ierr;
600:   PetscMPIInt       rank,size,*lens;
601:   PetscInt          i,N = xin->map->N;
602:   const PetscScalar *xarray;
603:   PetscScalar       *xx;

606:   VecGetArrayRead(xin,&xarray);
607:   MPI_Comm_rank(PetscObjectComm((PetscObject)xin),&rank);
608:   MPI_Comm_size(PetscObjectComm((PetscObject)xin),&size);
609:   if (!rank) {
610:     PetscMalloc1(N,&xx);
611:     PetscMalloc1(size,&lens);
612:     for (i=0; i<size; i++) lens[i] = xin->map->range[i+1] - xin->map->range[i];

614:     MPI_Gatherv((void*)xarray,xin->map->n,MPIU_SCALAR,xx,lens,xin->map->range,MPIU_SCALAR,0,PetscObjectComm((PetscObject)xin));
615:     PetscFree(lens);

617:     PetscObjectName((PetscObject)xin);
618:     PetscViewerMatlabPutArray(viewer,N,1,xx,((PetscObject)xin)->name);

620:     PetscFree(xx);
621:   } else {
622:     MPI_Gatherv((void*)xarray,xin->map->n,MPIU_SCALAR,0,0,0,MPIU_SCALAR,0,PetscObjectComm((PetscObject)xin));
623:   }
624:   VecRestoreArrayRead(xin,&xarray);
625:   return(0);
626: }
627: #endif

629: #if defined(PETSC_HAVE_HDF5)
632: PetscErrorCode VecView_MPI_HDF5(Vec xin, PetscViewer viewer)
633: {
634:   /* TODO: It looks like we can remove the H5Sclose(filespace) and H5Dget_space(dset_id). Why do we do this? */
635:   hid_t             filespace; /* file dataspace identifier */
636:   hid_t             chunkspace; /* chunk dataset property identifier */
637:   hid_t             plist_id;  /* property list identifier */
638:   hid_t             dset_id;   /* dataset identifier */
639:   hid_t             memspace;  /* memory dataspace identifier */
640:   hid_t             file_id;
641:   hid_t             group;
642:   hid_t             scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
643:   herr_t            status;
644:   PetscInt          bs = xin->map->bs > 0 ? xin->map->bs : 1;
645:   hsize_t           dim;
646:   hsize_t           maxDims[4], dims[4], chunkDims[4], count[4],offset[4];
647:   PetscInt          timestep;
648:   PetscInt          low;
649:   const PetscScalar *x;
650:   const char        *vecname;
651:   PetscErrorCode    ierr;

654:   PetscViewerHDF5OpenGroup(viewer, &file_id, &group);
655:   PetscViewerHDF5GetTimestep(viewer, &timestep);

657:   /* Create the dataspace for the dataset.
658:    *
659:    * dims - holds the current dimensions of the dataset
660:    *
661:    * maxDims - holds the maximum dimensions of the dataset (unlimited
662:    * for the number of time steps with the current dimensions for the
663:    * other dimensions; so only additional time steps can be added).
664:    *
665:    * chunkDims - holds the size of a single time step (required to
666:    * permit extending dataset).
667:    */
668:   dim = 0;
669:   if (timestep >= 0) {
670:     dims[dim]      = timestep+1;
671:     maxDims[dim]   = H5S_UNLIMITED;
672:     chunkDims[dim] = 1;
673:     ++dim;
674:   }
675:   PetscHDF5IntCast(xin->map->N/bs,dims + dim);

677:   maxDims[dim]   = dims[dim];
678:   chunkDims[dim] = dims[dim];
679:   ++dim;
680:   if (bs >= 1) {
681:     dims[dim]      = bs;
682:     maxDims[dim]   = dims[dim];
683:     chunkDims[dim] = dims[dim];
684:     ++dim;
685:   }
686: #if defined(PETSC_USE_COMPLEX)
687:   dims[dim]      = 2;
688:   maxDims[dim]   = dims[dim];
689:   chunkDims[dim] = dims[dim];
690:   ++dim;
691: #endif
692:   filespace = H5Screate_simple(dim, dims, maxDims);
693:   if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()");

695: #if defined(PETSC_USE_REAL_SINGLE)
696:   scalartype = H5T_NATIVE_FLOAT;
697: #elif defined(PETSC_USE_REAL___FLOAT128)
698: #error "HDF5 output with 128 bit floats not supported."
699: #else
700:   scalartype = H5T_NATIVE_DOUBLE;
701: #endif

703:   /* Create the dataset with default properties and close filespace */
704:   PetscObjectGetName((PetscObject) xin, &vecname);
705:   if (!H5Lexists(group, vecname, H5P_DEFAULT)) {
706:     /* Create chunk */
707:     chunkspace = H5Pcreate(H5P_DATASET_CREATE);
708:     if (chunkspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Pcreate()");
709:     status = H5Pset_chunk(chunkspace, dim, chunkDims);CHKERRQ(status);

711: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
712:     dset_id = H5Dcreate2(group, vecname, scalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT);
713: #else
714:     dset_id = H5Dcreate(group, vecname, scalartype, filespace, H5P_DEFAULT);
715: #endif
716:     if (dset_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dcreate2()");
717:     status = H5Pclose(chunkspace);CHKERRQ(status);
718:   } else {
719:     dset_id = H5Dopen2(group, vecname, H5P_DEFAULT);
720:     status  = H5Dset_extent(dset_id, dims);CHKERRQ(status);
721:   }
722:   status = H5Sclose(filespace);CHKERRQ(status);

724:   /* Each process defines a dataset and writes it to the hyperslab in the file */
725:   dim = 0;
726:   if (timestep >= 0) {
727:     count[dim] = 1;
728:     ++dim;
729:   }
730:   PetscHDF5IntCast(xin->map->n/bs,count + dim);
731:   ++dim;
732:   if (bs >= 1) {
733:     count[dim] = bs;
734:     ++dim;
735:   }
736: #if defined(PETSC_USE_COMPLEX)
737:   count[dim] = 2;
738:   ++dim;
739: #endif
740:   if (xin->map->n > 0) {
741:     memspace = H5Screate_simple(dim, count, NULL);
742:     if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()");
743:   } else {
744:     /* Can't create dataspace with zero for any dimension, so create null dataspace. */
745:     memspace = H5Screate(H5S_NULL);
746:     if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate()");
747:   }

749:   /* Select hyperslab in the file */
750:   VecGetOwnershipRange(xin, &low, NULL);
751:   dim  = 0;
752:   if (timestep >= 0) {
753:     offset[dim] = timestep;
754:     ++dim;
755:   }
756:   PetscHDF5IntCast(low/bs,offset + dim);
757:   ++dim;
758:   if (bs >= 1) {
759:     offset[dim] = 0;
760:     ++dim;
761:   }
762: #if defined(PETSC_USE_COMPLEX)
763:   offset[dim] = 0;
764:   ++dim;
765: #endif
766:   if (xin->map->n > 0) {
767:     filespace = H5Dget_space(dset_id);
768:     if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dget_space()");
769:     status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status);
770:   } else {
771:     /* Create null filespace to match null memspace. */
772:     filespace = H5Screate(H5S_NULL);
773:     if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate(H5S_NULL)");
774:   }

776:   /* Create property list for collective dataset write */
777:   plist_id = H5Pcreate(H5P_DATASET_XFER);
778:   if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Pcreate()");
779: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
780:   status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
781: #endif
782:   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */

784:   VecGetArrayRead(xin, &x);
785:   status = H5Dwrite(dset_id, scalartype, memspace, filespace, plist_id, x);CHKERRQ(status);
786:   status = H5Fflush(file_id, H5F_SCOPE_GLOBAL);CHKERRQ(status);
787:   VecRestoreArrayRead(xin, &x);

789:   /* Close/release resources */
790:   if (group != file_id) {
791:     status = H5Gclose(group);CHKERRQ(status);
792:   }
793:   status = H5Pclose(plist_id);CHKERRQ(status);
794:   status = H5Sclose(filespace);CHKERRQ(status);
795:   status = H5Sclose(memspace);CHKERRQ(status);
796:   status = H5Dclose(dset_id);CHKERRQ(status);
797:   PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);
798:   return(0);
799: }
800: #endif

804: PetscErrorCode VecView_MPI(Vec xin,PetscViewer viewer)
805: {
807:   PetscBool      iascii,isbinary,isdraw;
808: #if defined(PETSC_HAVE_MATHEMATICA)
809:   PetscBool      ismathematica;
810: #endif
811: #if defined(PETSC_HAVE_HDF5)
812:   PetscBool      ishdf5;
813: #endif
814: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
815:   PetscBool      ismatlab;
816: #endif

819:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
820:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
821:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
822: #if defined(PETSC_HAVE_MATHEMATICA)
823:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATHEMATICA,&ismathematica);
824: #endif
825: #if defined(PETSC_HAVE_HDF5)
826:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
827: #endif
828: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
829:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERMATLAB,&ismatlab);
830: #endif
831:   if (iascii) {
832:     VecView_MPI_ASCII(xin,viewer);
833:   } else if (isbinary) {
834:     VecView_MPI_Binary(xin,viewer);
835:   } else if (isdraw) {
836:     PetscViewerFormat format;

838:     PetscViewerGetFormat(viewer,&format);
839:     if (format == PETSC_VIEWER_DRAW_LG) {
840:       VecView_MPI_Draw_LG(xin,viewer);
841:     } else {
842:       VecView_MPI_Draw(xin,viewer);
843:     }
844: #if defined(PETSC_HAVE_MATHEMATICA)
845:   } else if (ismathematica) {
846:     PetscViewerMathematicaPutVector(viewer,xin);
847: #endif
848: #if defined(PETSC_HAVE_HDF5)
849:   } else if (ishdf5) {
850:     VecView_MPI_HDF5(xin,viewer);
851: #endif
852: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
853:   } else if (ismatlab) {
854:     VecView_MPI_Matlab(xin,viewer);
855: #endif
856:   }
857:   return(0);
858: }

862: PetscErrorCode VecGetSize_MPI(Vec xin,PetscInt *N)
863: {
865:   *N = xin->map->N;
866:   return(0);
867: }

871: PetscErrorCode VecGetValues_MPI(Vec xin,PetscInt ni,const PetscInt ix[],PetscScalar y[])
872: {
873:   const PetscScalar *xx;
874:   PetscInt          i,tmp,start = xin->map->range[xin->stash.rank];
875:   PetscErrorCode    ierr;

878:   VecGetArrayRead(xin,&xx);
879:   for (i=0; i<ni; i++) {
880:     if (xin->stash.ignorenegidx && ix[i] < 0) continue;
881:     tmp = ix[i] - start;
882: #if defined(PETSC_USE_DEBUG)
883:     if (tmp < 0 || tmp >= xin->map->n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Can only get local values, trying %D",ix[i]);
884: #endif
885:     y[i] = xx[tmp];
886:   }
887:   VecRestoreArrayRead(xin,&xx);
888:   return(0);
889: }

893: PetscErrorCode VecSetValues_MPI(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode addv)
894: {
896:   PetscMPIInt    rank    = xin->stash.rank;
897:   PetscInt       *owners = xin->map->range,start = owners[rank];
898:   PetscInt       end     = owners[rank+1],i,row;
899:   PetscScalar    *xx;

902: #if defined(PETSC_USE_DEBUG)
903:   if (xin->stash.insertmode == INSERT_VALUES && addv == ADD_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"You have already inserted values; you cannot now add");
904:   else if (xin->stash.insertmode == ADD_VALUES && addv == INSERT_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"You have already added values; you cannot now insert");
905: #endif
906:   VecGetArray(xin,&xx);
907:   xin->stash.insertmode = addv;

909:   if (addv == INSERT_VALUES) {
910:     for (i=0; i<ni; i++) {
911:       if (xin->stash.ignorenegidx && ix[i] < 0) continue;
912: #if defined(PETSC_USE_DEBUG)
913:       if (ix[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
914: #endif
915:       if ((row = ix[i]) >= start && row < end) {
916:         xx[row-start] = y[i];
917:       } else if (!xin->stash.donotstash) {
918: #if defined(PETSC_USE_DEBUG)
919:         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);
920: #endif
921:         VecStashValue_Private(&xin->stash,row,y[i]);
922:       }
923:     }
924:   } else {
925:     for (i=0; i<ni; i++) {
926:       if (xin->stash.ignorenegidx && ix[i] < 0) continue;
927: #if defined(PETSC_USE_DEBUG)
928:       if (ix[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D cannot be negative",ix[i]);
929: #endif
930:       if ((row = ix[i]) >= start && row < end) {
931:         xx[row-start] += y[i];
932:       } else if (!xin->stash.donotstash) {
933: #if defined(PETSC_USE_DEBUG)
934:         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);
935: #endif
936:         VecStashValue_Private(&xin->stash,row,y[i]);
937:       }
938:     }
939:   }
940:   VecRestoreArray(xin,&xx);
941:   return(0);
942: }

946: PetscErrorCode VecSetValuesBlocked_MPI(Vec xin,PetscInt ni,const PetscInt ix[],const PetscScalar yin[],InsertMode addv)
947: {
948:   PetscMPIInt    rank    = xin->stash.rank;
949:   PetscInt       *owners = xin->map->range,start = owners[rank];
951:   PetscInt       end = owners[rank+1],i,row,bs = xin->map->bs,j;
952:   PetscScalar    *xx,*y = (PetscScalar*)yin;

955:   VecGetArray(xin,&xx);
956: #if defined(PETSC_USE_DEBUG)
957:   if (xin->stash.insertmode == INSERT_VALUES && addv == ADD_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"You have already inserted values; you cannot now add");
958:   else if (xin->stash.insertmode == ADD_VALUES && addv == INSERT_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"You have already added values; you cannot now insert");
959: #endif
960:   xin->stash.insertmode = addv;

962:   if (addv == INSERT_VALUES) {
963:     for (i=0; i<ni; i++) {
964:       if ((row = bs*ix[i]) >= start && row < end) {
965:         for (j=0; j<bs; j++) xx[row-start+j] = y[j];
966:       } else if (!xin->stash.donotstash) {
967:         if (ix[i] < 0) continue;
968: #if defined(PETSC_USE_DEBUG)
969:         if (ix[i] >= xin->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D max %D",ix[i],xin->map->N);
970: #endif
971:         VecStashValuesBlocked_Private(&xin->bstash,ix[i],y);
972:       }
973:       y += bs;
974:     }
975:   } else {
976:     for (i=0; i<ni; i++) {
977:       if ((row = bs*ix[i]) >= start && row < end) {
978:         for (j=0; j<bs; j++) xx[row-start+j] += y[j];
979:       } else if (!xin->stash.donotstash) {
980:         if (ix[i] < 0) continue;
981: #if defined(PETSC_USE_DEBUG)
982:         if (ix[i] > xin->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D max %D",ix[i],xin->map->N);
983: #endif
984:         VecStashValuesBlocked_Private(&xin->bstash,ix[i],y);
985:       }
986:       y += bs;
987:     }
988:   }
989:   VecRestoreArray(xin,&xx);
990:   return(0);
991: }

993: /*
994:    Since nsends or nreceives may be zero we add 1 in certain mallocs
995: to make sure we never malloc an empty one.
996: */
999: PetscErrorCode VecAssemblyBegin_MPI(Vec xin)
1000: {
1002:   PetscInt       *owners = xin->map->range,*bowners,i,bs,nstash,reallocs;
1003:   PetscMPIInt    size;
1004:   InsertMode     addv;
1005:   MPI_Comm       comm;

1008:   PetscObjectGetComm((PetscObject)xin,&comm);
1009:   if (xin->stash.donotstash) return(0);

1011:   MPI_Allreduce((PetscEnum*)&xin->stash.insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);
1012:   if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(comm,PETSC_ERR_ARG_NOTSAMETYPE,"Some processors inserted values while others added");
1013:   xin->stash.insertmode = addv; /* in case this processor had no cache */

1015:   bs   = xin->map->bs;
1016:   MPI_Comm_size(PetscObjectComm((PetscObject)xin),&size);
1017:   if (!xin->bstash.bowners && xin->map->bs != -1) {
1018:     PetscMalloc1((size+1),&bowners);
1019:     for (i=0; i<size+1; i++) bowners[i] = owners[i]/bs;
1020:     xin->bstash.bowners = bowners;
1021:   } else bowners = xin->bstash.bowners;

1023:   VecStashScatterBegin_Private(&xin->stash,owners);
1024:   VecStashScatterBegin_Private(&xin->bstash,bowners);
1025:   VecStashGetInfo_Private(&xin->stash,&nstash,&reallocs);
1026:   PetscInfo2(xin,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
1027:   VecStashGetInfo_Private(&xin->bstash,&nstash,&reallocs);
1028:   PetscInfo2(xin,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
1029:   return(0);
1030: }

1034: PetscErrorCode VecAssemblyEnd_MPI(Vec vec)
1035: {
1037:   PetscInt       base,i,j,*row,flg,bs;
1038:   PetscMPIInt    n;
1039:   PetscScalar    *val,*vv,*array,*xarray;

1042:   if (!vec->stash.donotstash) {
1043:     VecGetArray(vec,&xarray);
1044:     base = vec->map->range[vec->stash.rank];
1045:     bs   = vec->map->bs;

1047:     /* Process the stash */
1048:     while (1) {
1049:       VecStashScatterGetMesg_Private(&vec->stash,&n,&row,&val,&flg);
1050:       if (!flg) break;
1051:       if (vec->stash.insertmode == ADD_VALUES) {
1052:         for (i=0; i<n; i++) xarray[row[i] - base] += val[i];
1053:       } else if (vec->stash.insertmode == INSERT_VALUES) {
1054:         for (i=0; i<n; i++) xarray[row[i] - base] = val[i];
1055:       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Insert mode is not set correctly; corrupted vector");
1056:     }
1057:     VecStashScatterEnd_Private(&vec->stash);

1059:     /* now process the block-stash */
1060:     while (1) {
1061:       VecStashScatterGetMesg_Private(&vec->bstash,&n,&row,&val,&flg);
1062:       if (!flg) break;
1063:       for (i=0; i<n; i++) {
1064:         array = xarray+row[i]*bs-base;
1065:         vv    = val+i*bs;
1066:         if (vec->stash.insertmode == ADD_VALUES) {
1067:           for (j=0; j<bs; j++) array[j] += vv[j];
1068:         } else if (vec->stash.insertmode == INSERT_VALUES) {
1069:           for (j=0; j<bs; j++) array[j] = vv[j];
1070:         } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Insert mode is not set correctly; corrupted vector");
1071:       }
1072:     }
1073:     VecStashScatterEnd_Private(&vec->bstash);
1074:     VecRestoreArray(vec,&xarray);
1075:   }
1076:   vec->stash.insertmode = NOT_SET_VALUES;
1077:   return(0);
1078: }