<html><head><meta http-equiv="Content-Type" content="text/html; charset=us-ascii"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class=""><br class=""></div>   You need to first convert the Vec to the "natural" ordering and then bring that vector down to one rank. Something like<div class=""><br class=""></div><div class="">DMDACreateNaturalVector(dm,&n);</div><div class="">DMDAGlobalToNaturalBegin/End(dm,   n);</div><div class=""><blockquote type="cite" class=""><div id="divtagdefaultwrapper" dir="ltr" class="" style="font-size: 12pt; font-family: Calibri, Helvetica, sans-serif, Helvetica, EmojiFont, "Apple Color Emoji", "Segoe UI Emoji", NotoColorEmoji, "Segoe UI Symbol", "Android Emoji", EmojiSymbols;"><div class="" style="margin-top: 0px; margin-bottom: 0px;">VecScatterCreateToZero(n, &scat, &Xseq);</div><div class="" style="margin-top: 0px; margin-bottom: 0px;">  VecScatterBegin(scat, n Xseq, INSERT_VALUES, SCATTER_FORWARD);</div><div class="" style="margin-top: 0px; margin-bottom: 0px;">  VecScatterEnd(scat, n,  Xseq, INSERT_VALUES, SCATTER_FORWARD);</div></div></blockquote><div class=""><br class=""></div><div class=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On Oct 11, 2022, at 7:45 PM, Jorti, Zakariae via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><meta charset="UTF-8" class=""><div id="divtagdefaultwrapper" dir="ltr" style="font-style: normal; font-variant-caps: normal; font-weight: 400; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none; font-size: 12pt; font-family: Calibri, Helvetica, sans-serif, Helvetica, EmojiFont, "Apple Color Emoji", "Segoe UI Emoji", NotoColorEmoji, "Segoe UI Symbol", "Android Emoji", EmojiSymbols;" class=""><div style="margin-top: 0px; margin-bottom: 0px;" class="">Hello, </div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">I have a code that handles a PETSc Vec on many procs and I would like to use VecScatterCreateToZero to have all elements of this vector on a single proc, call VecGetArrayRead on this proc to get the corresponding array to carry out some diagnostics.</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">Unfortunately, what I noticed is that the ordering of the initial Vec is not preserved after the VecScatterCreateToZero call. Is there a way to have the same ordering for both Vecs? </div><div style="margin-top: 0px; margin-bottom: 0px;" class="">You will find below a minimal example that demonstrates the issue.</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">Many thanks,</div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">Zakariae   </div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">-------------------------------------------------------------</div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">static char help[] = "Demonstrates ordering change after VecScatterCreateToZero call.\n\n";</div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">#include <petscpf.h></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">#include <petscdm.h></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">#include <petscdmda.h></div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">int main(int argc,char **argv)</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">{</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  Vec            xy;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  DM             da;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  PetscErrorCode ierr;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  PetscInt       m = 11, n = 11, dof = 2;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  PetscMPIInt rank;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  DMDACoor2d **coors;</div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,m,n,PETSC_DECIDE,PETSC_DECIDE,dof,1,0,0,&da);CHKERRQ(ierr);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  ierr = DMSetFromOptions(da);CHKERRQ(ierr);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  ierr = DMSetUp(da);CHKERRQ(ierr);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  ierr = DMGetCoordinates(da,&xy);CHKERRQ(ierr);</div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  PetscInt      i, j, ixs, ixm, iys, iym;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);</div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  DMDAGetCorners(da, &ixs, &iys, 0, &ixm, &iym, 0);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  DMDAVecGetArray(da, xy, &coors);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  for (j = iys; j < iys + iym; j++) {</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">    for (i = ixs; i < ixs + ixm; i++) {</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">        PetscPrintf(PETSC_COMM_SELF, "rank=%d, %d, %d, (%g, %g)\n",rank, </div><div style="margin-top: 0px; margin-bottom: 0px;" class="">                    i, j,coors[j][i].x,coors[j][i].y);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">    }</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  }</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  DMDAVecRestoreArray(da, xy, &coors);</div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  VecScatter  scat;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  Vec         Xseq;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  const PetscScalar *array; </div><p style="margin-top: 0px; margin-bottom: 0px;" class="">   </p><div style="margin-top: 0px; margin-bottom: 0px;" class="">  /* create scater to zero */</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  VecScatterCreateToZero(xy, &scat, &Xseq);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  VecScatterBegin(scat, xy, Xseq, INSERT_VALUES, SCATTER_FORWARD);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  VecScatterEnd(scat, xy, Xseq, INSERT_VALUES, SCATTER_FORWARD);</div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  if (rank == 0) {</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">    PetscInt sizeX;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">    VecGetSize(Xseq, &sizeX); </div><div style="margin-top: 0px; margin-bottom: 0px;" class="">    PetscPrintf(PETSC_COMM_SELF,"The size of Xseq is %d, and the grid size is %d\n",sizeX,11*11);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">    VecGetArrayRead(Xseq, &array);</div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">    for (j = 0; j < 11; j++) {</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">      for (i = 0; i < 11; i++) {</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">        PetscPrintf(PETSC_COMM_SELF, "%d,%d, (%g,%g)\n", i, j, (double)array[2*(j*11+i)], (double)array[1+2*(j*11+i)]);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">      }</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">    }</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">    VecRestoreArrayRead(Xseq, &array);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  }</div><div style="margin-top: 0px; margin-bottom: 0px;" class=""><br class=""></div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  /*</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">     Free work space.  All PETSc objects should be destroyed when they</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">     are no longer needed.</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  */</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  ierr = DMDestroy(&da);CHKERRQ(ierr);</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  ierr = PetscFinalize();</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">  return ierr;</div><div style="margin-top: 0px; margin-bottom: 0px;" class="">}</div></div></div></blockquote></div><br class=""></div></div></body></html>