<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Thu, Aug 25, 2016 at 1:59 PM, Mark Adams <span dir="ltr"><<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div>We have the subroutine below that scatters three vectors.  We have used this code on many machines and it works fine but on one machine data does not get scattered correctly. The first scatter looks OK, but it looks like the other two are missing data.</div><div><br></div><div>Am I using this correctly?  Should we use VecGetArray in here instead of just using the pointer used for construction? Is there a race condition here that I'm missing?</div><div><br></div><div>Thanks,</div><div>Mark</div><div><br></div><div>subroutine scatter_to_xgc(a_ts,a_XX,a_n1,<wbr>a_apar,a_phi,ierr)</div><div>  use petscts</div><div>  use sml_module,only:sml_mype</div><div>  use xgc_ts_module</div><div>  implicit none</div><div>  type(xgc_ts),intent(in)::a_ts</div><div>  Vec,intent(in)::a_XX</div><div>  real (kind=8),dimension(a_ts%nnode)<wbr>::a_n1,a_apar,a_phi</div><div>  PetscErrorCode,intent(out)::<wbr>ierr</div><div>  ! locals</div><div>  PetscInt,parameter::ione=1</div><div>  PetscScalar,dimension(a_ts%<wbr>nnode)::n1,apar,phi</div><div>  Vec::xxvec(0:2)</div><div><br></div><div>  ! scatter solution back - n1</div><div>  n1 = a_n1</div></div></blockquote><div><br></div><div>Let me be clearer. Here n1 is equal to a_n1, and used to back xxvec(0). It does not change</div><div>in the course of the function, so why would you have</div><div><br></div><div>  a_n1 = n</div><div><br></div><div>later on? What am I missing?</div><div><br></div><div>   Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div>  call VecCreateSeqWithArray(PETSC_<wbr>COMM_SELF,ione,a_ts%nnode,n1,<wbr>xxvec(0),ierr);CHKERRQ(ierr)</div><div>  call VecScatterBegin(a_ts%from_<wbr>petsc(0),a_XX,xxvec(0),INSERT_<wbr>VALUES,SCATTER_FORWARD,ierr)</div><div><br></div><div>  ! scatter solution back - apar</div><div>  apar = a_apar</div><div>  call VecCreateSeqWithArray(PETSC_<wbr>COMM_SELF,ione,a_ts%nnode,<wbr>apar,xxvec(1),ierr);CHKERRQ(<wbr>ierr)</div><div>  call VecScatterBegin(a_ts%from_<wbr>petsc(1),a_XX,xxvec(1),INSERT_<wbr>VALUES,SCATTER_FORWARD,ierr)</div><div><br></div><div>  ! scatter solution back - phi</div><div>  phi = a_phi</div><div>  call VecCreateSeqWithArray(PETSC_<wbr>COMM_SELF,ione,a_ts%nnode,phi,<wbr>xxvec(2),ierr);CHKERRQ(ierr)</div><div>  call VecScatterBegin(a_ts%from_<wbr>petsc(2),a_XX,xxvec(2),INSERT_<wbr>VALUES,SCATTER_FORWARD,ierr)</div><div><br></div><div>  ! end</div><div>  call VecScatterEnd(  a_ts%from_petsc(0),a_XX,<wbr>xxvec(0),INSERT_VALUES,<wbr>SCATTER_FORWARD,ierr)</div><div>  a_n1 = n1</div><div>  call VecDestroy(xxvec(0),ierr)</div><div><br></div><div>  call VecScatterEnd(  a_ts%from_petsc(1),a_XX,<wbr>xxvec(1),INSERT_VALUES,<wbr>SCATTER_FORWARD,ierr)</div><div>  a_apar = apar</div><div>  call VecDestroy(xxvec(1),ierr)</div><div><br></div><div>  call VecScatterEnd(  a_ts%from_petsc(2),a_XX,<wbr>xxvec(2),INSERT_VALUES,<wbr>SCATTER_FORWARD,ierr)</div><div>  a_phi = phi</div><div>  call VecDestroy(xxvec(2),ierr)</div><div><br></div><div>  return</div><div>end subroutine scatter_to_xgc</div></div>
</blockquote></div><br><br clear="all"><div><br></div>-- <br><div class="gmail_signature" data-smartmail="gmail_signature">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div>
</div></div>