[petsc-dev] VecScatter problem

Matthew Knepley knepley at gmail.com
Thu Aug 25 14:56:31 CDT 2016


On Thu, Aug 25, 2016 at 2:38 PM, Robert Hager <rhager at pppl.gov> wrote:

> Isn’t xxvec(0) basically a pointer to n1 after the call to
> VecCreateSeqWithArray? If so, the VecScatter should update the values in n1
> and setting a_n1=n1 in the end makes sense.
>

xxvec(0) holds a pointer (n1) to the data. That pointer does not change
when the scatter writes in to the Vec, so
at the end the pointer n1 is still the same, namely a_n1.

   Matt


> Robert
>
>
>
> On Aug 25, 2016, at 3:16 PM, Matthew Knepley <knepley at gmail.com> wrote:
>
> On Thu, Aug 25, 2016 at 1:59 PM, Mark Adams <mfadams at lbl.gov> wrote:
>
>> 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.
>>
>> 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?
>>
>> Thanks,
>> Mark
>>
>> subroutine scatter_to_xgc(a_ts,a_XX,a_n1,a_apar,a_phi,ierr)
>>   use petscts
>>   use sml_module,only:sml_mype
>>   use xgc_ts_module
>>   implicit none
>>   type(xgc_ts),intent(in)::a_ts
>>   Vec,intent(in)::a_XX
>>   real (kind=8),dimension(a_ts%nnode)::a_n1,a_apar,a_phi
>>   PetscErrorCode,intent(out)::ierr
>>   ! locals
>>   PetscInt,parameter::ione=1
>>   PetscScalar,dimension(a_ts%nnode)::n1,apar,phi
>>   Vec::xxvec(0:2)
>>
>>   ! scatter solution back - n1
>>   n1 = a_n1
>>
>
> Let me be clearer. Here n1 is equal to a_n1, and used to back xxvec(0). It
> does not change
> in the course of the function, so why would you have
>
>   a_n1 = n
>
> later on? What am I missing?
>
>    Matt
>
>
>>   call VecCreateSeqWithArray(PETSC_COMM_SELF,ione,a_ts%nnode,n1,xxv
>> ec(0),ierr);CHKERRQ(ierr)
>>   call VecScatterBegin(a_ts%from_petsc(0),a_XX,xxvec(0),INSERT_VALU
>> ES,SCATTER_FORWARD,ierr)
>>
>>   ! scatter solution back - apar
>>   apar = a_apar
>>   call VecCreateSeqWithArray(PETSC_COMM_SELF,ione,a_ts%nnode,apar,
>> xxvec(1),ierr);CHKERRQ(ierr)
>>   call VecScatterBegin(a_ts%from_petsc(1),a_XX,xxvec(1),INSERT_VALU
>> ES,SCATTER_FORWARD,ierr)
>>
>>   ! scatter solution back - phi
>>   phi = a_phi
>>   call VecCreateSeqWithArray(PETSC_COMM_SELF,ione,a_ts%nnode,phi,xx
>> vec(2),ierr);CHKERRQ(ierr)
>>   call VecScatterBegin(a_ts%from_petsc(2),a_XX,xxvec(2),INSERT_VALU
>> ES,SCATTER_FORWARD,ierr)
>>
>>   ! end
>>   call VecScatterEnd(  a_ts%from_petsc(0),a_XX,xxvec
>> (0),INSERT_VALUES,SCATTER_FORWARD,ierr)
>>   a_n1 = n1
>>   call VecDestroy(xxvec(0),ierr)
>>
>>   call VecScatterEnd(  a_ts%from_petsc(1),a_XX,xxvec
>> (1),INSERT_VALUES,SCATTER_FORWARD,ierr)
>>   a_apar = apar
>>   call VecDestroy(xxvec(1),ierr)
>>
>>   call VecScatterEnd(  a_ts%from_petsc(2),a_XX,xxvec
>> (2),INSERT_VALUES,SCATTER_FORWARD,ierr)
>>   a_phi = phi
>>   call VecDestroy(xxvec(2),ierr)
>>
>>   return
>> end subroutine scatter_to_xgc
>>
>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
>
>


-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20160825/e4d583d4/attachment.html>


More information about the petsc-dev mailing list