[petsc-users] Accessing Vector's ghost values

Mohammad Mirzadeh mirzadeh at gmail.com
Fri Feb 24 04:43:43 CST 2012


Thanks Bojan. I suppose that would work. Problem arises, however, when the
actual application numbering matters and you have to preserve it somehow.
One example, is where you have a library that is currently sequential and a
lot of functionality is directly based on the data structure used for the
grid. For example when computing derivatives or performing interpolation
and so on.

Naturally you would want to reuse as much of the code as possible. So, you
would think you can have a function, or a class, that, for example, takes
the ghosted vector and the original global numbering and returns the value
so that you could still use the old functions, say to compute gradients.
Now, when you are constructing the ghosted vector, you know order you are
putting the ghost nodes, lets say the first one has global index 42, next
has 21, next 56 and so on ... However, when you call the function with the
global index 21, you need to know what local ghost index it corresponds to
in the local form of the ghosted vector. I can think of the following
options to do this:

1- Search the ghost nodes indecies that you created the ghosted vector
with: here, search int ghost [] = {42, 21, 56} for 21 and return the index
which is 2. This, you can do at O(log n) e.g with binary search
2- Make a vector consisting of _all_ nodes on each processor that returns
the ghost index, i.e. have something like ghostIndex[21] = 2,
ghostIndex[42] = 1, ghostIndex[56] = 3. Using this you get the index in
O(1) but at the cost of extra memory
3- Make an ao that maps {1,2,3} to {42, 21, 56}. This seems to be the best
of both worlds

So, in short, I'm interested in something that wont alter (at
least permanently) the application global numbering so that I can use my
old pieces of code and at the same time something that can compute the
ghost location fast and efficient from a given global numbering.

Mohammad

On Fri, Feb 24, 2012 at 1:54 AM, Niceno Bojan <bojan.niceno at psi.ch> wrote:

> I do re-numbering during the domain decomposition stage.  I re-number
> nodes (unknows) in the way I believe PETSc does internaly anyhow.  It was a
> guess-work, I admit, but my communication patterns indeed work.  Here is
> what I do:
>
>
>  /*-------------------------------------+
>  |    Assign global numbers to nodes    |
>  +-------------------------------------*/
>  Int global_number = 0;
>  for(Int p=0; p<number_of_partitions; p++) {
>    for(Int n=0; n<size(); n++) {
>      if( nodes[n].partition == p) {
>        nodes[n].global_number = global_number;
>        global_number++;
>      }
>    }
>  }
>  assert( global_number == nodes.size());
>
>
> -----Original Message-----
> From: petsc-users-bounces at mcs.anl.gov on behalf of Mohammad Mirzadeh
> Sent: Fri 2/24/2012 10:42 AM
> To: PETSc users list
> Subject: Re: [petsc-users] Accessing Vector's ghost values
>
> Actually now that we are at it, what would you guys recommend as the best
> way to find local indecies from global ones when accessing the _ghost_
> values in the local form? This becomes even more complicated when your
> global petsc ordering is different from global application ordering.
>
> Its quite easy to access the local indecies for the local nodes -- you use
> ao to convert application global to petsc global and from petsc global, it
> is easy to get petsc local. However, when you have ghost nodes, there is no
> such simple relationship between petsc global and the _ghosted_ petsc local
> (Even though there is a simple one from petsc local to petsc global). The
> best thing I have came up with so far is to use a sequential second ao to
> map local indecies of ghost points to petsc global indecies. Is there a
> better way of doing this?
>
> Mohammad
>
> On Thu, Feb 23, 2012 at 7:44 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
>
> >
> > I have added the following to the VecGhostGetLocalForm() manual page to
> > make this clearer for everyone:
> >
> >    To update the ghost values from the locations on the other processes
> > one must call
> >    VecGhostUpdateBegin() and VecGhostUpdateEnd() before accessing the
> > ghost values. Thus normal
> >    usage is
> > $     VecGhostUpdateBegin(x,INSERT_VALUES,SCATTER_FORWARD);
> > $     VecGhostUpdateEnd(x,INSERT_VALUES,SCATTER_FORWARD);
> > $     VecGhostGetLocalForm(x,&xlocal);
> > $     VecGetArray(xlocal,&xvalues);
> > $        /* access the non-ghost values in locations xvalues[0:n-1] and
> > ghost values in locations xvalues[n:n+nghost]; */
> > $     VecRestoreArray(xlocal,&xvalues);
> > $     VecGhostRestoreLocalForm(x,&xlocal);
> >
> >
> > On Feb 23, 2012, at 5:38 PM, Bojan Niceno wrote:
> >
> > > Yeeee-ha!
> > >
> > > VecGhostUpdateBegin() / VecGhostUpdateEnd() was indeed what I was
> > missing.
> > >
> > > Thank you Mohammad, thank you all who helped me today!
> > >
> > >
> > >     Cheers
> > >
> > >
> > >     Bojan
> > >
> > >
> > >
> > > On 2/24/2012 12:27 AM, Mohammad Mirzadeh wrote:
> > >> You also need calls to VecGhostUpdateBegin()/VecGhostUpdateEnd()
> > functions to update the ghost values if you change them in the global
> > representation. See Petsc Manual 3.2 pp 55-56
> > >>
> > >> Mohammad
> > >>
> > >> On Thu, Feb 23, 2012 at 3:18 PM, Bojan Niceno <bojan.niceno at psi.ch>
> > wrote:
> > >> On 2/23/2012 11:49 PM, Mohammad Mirzadeh wrote:
> > >>> based on,
> > >>>
> > >>> VecCreateGhost(PETSC_COMM_WORLD, n, PETSC_DECIDE, nghost, &ghosts[0],
> > &x);
> > >>>
> > >>>
> > >>> it seems to me that x is your actual ghosted vector. If this is true,
> > then you need to get its "local" form via VecGhostGetLocalForm(). Once
> you
> > have done, you should be able to access the ghosted nodes. Are you
> calling
> > this function anywhere?
> > >>
> > >> I tried that before.  I did:
> > >>
> > >> Vec lx;
> > >> VecGhostGetLocalForm(x, &lx)
> > >>
> > >> then I copied "lx" to my variable, like
> > >>
> > >> for(Int i=1; i<N; i++)    /* N includes buffer cells */
> > >>   unk[i] = lx[i]
> > >>
> > >> but ghost values were also zero.   I am thinking that PETSc somehow
> > clears the ghost values after a call to KSP.  Is it the case?
> > >>
> > >>
> > >>     Kind regards,
> > >>
> > >>
> > >>     Bojan
> > >>
> > >>
> > >>>
> > >>> On Thu, Feb 23, 2012 at 2:32 PM, Bojan Niceno <bojan.niceno at psi.ch>
> > wrote:
> > >>> Dear Mohammad,
> > >>>
> > >>>
> > >>> it doesn't help me, or I did not understand your explanation.
> > >>>
> > >>> If I do this:
> > >>>
> > >>>   /* copy internal values (THIS WORKS, BUT COPIES NO BUFFER VALUES)
> */
> > >>>
> > >>>   for(Int i=0; i<n; i++) {
> > >>>     Int gi = mesh.nodes[i].global_number;
> > >>>     VecGetValues(x, 1, &gi, &unk[i]);
> > >>>   }
> > >>>
> > >>>   /* copy ghost values (CREATES MANY WARNINGS */
> > >>>
> > >>>   for(Int i=n; i<N; i++) {
> > >>>     VecGetValues(x, 1, &i, &unk[i]);
> > >>>   }
> > >>>
> > >>> I get arnings are like this.
> > >>>
> > >>> [0]PETSC ERROR: --------------------- Error Message
> > ------------------------------------
> > >>> [0]PETSC ERROR: Argument out of range!
> > >>> [0]PETSC ERROR: Can only get local values, trying 3518!
> > >>> [3]PETSC ERROR: --------------------- Error Message
> > ------------------------------------
> > >>> [3]PETSC ERROR: Argument out of range!
> > >>> [3]PETSC ERROR: Can only get local values, trying 3511!
> > >>> [3]PETSC ERROR:
> > ------------------------------------------------------------------------
> > >>>
> > >>> What am I doing wrong here?
> > >>>
> > >>>
> > >>>     Cheers,
> > >>>
> > >>>
> > >>>     Bojan
> > >>>
> > >>>
> > >>>
> > >>> On 2/23/2012 11:23 PM, Mohammad Mirzadeh wrote:
> > >>>> just index x with the local numberings. if you have 'n' local nodes
> > and 'g' ghost nodes in x, ghost nodes indecies run from 'n' to 'n+g-1'
> > >>>>
> > >>>> On Feb 23, 2012 1:16 PM, "Bojan Niceno" <bojan.niceno at psi.ch>
> wrote:
> > >>>> Dear Jed,
> > >>>>
> > >>>> thanks.
> > >>>>
> > >>>> Now I have the following:
> > >>>>
> > >>>> - Array unk, which should hold values inside the partition and in
> > ghost cells.  It is big enough to hold both
> > >>>> - Vec x, created by command VecCreateGhost, with proper padding for
> > ghost cells
> > >>>> - Successful call to linear solver in parallel.
> > >>>>
> > >>>> But I need to copy ghost values from x to my array unk.
> > >>>>
> > >>>> How can I do it?
> > >>>>
> > >>>>
> > >>>>     Kind regards,
> > >>>>
> > >>>>
> > >>>>     Bojan
> > >>>>
> > >>>> On 2/23/2012 10:10 PM, Jed Brown wrote:
> > >>>>> On Thu, Feb 23, 2012 at 15:05, Bojan Niceno <bojan.niceno at psi.ch>
> > wrote:
> > >>>>> No, I use global.
> > >>>>>
> > >>>>> The local form is just a local vector. It doesn't even know that a
> > global problem exists. You can't index into it using global indices. (In
> > general, there is no efficient way to look up information in the local
> > vector (includes ghost points) using global indices.)
> > >>>>>
> > >>>>>
> > >>>>>   for(Int i=0; i<n; i++) {
> > >>>>>     Int gi = mesh.nodes[i].global_number;
> > >>>>>     VecGetValues(x, 1, &gi, &unk[i]);
> > >>>>>   }
> > >>>>>
> > >>>>> "n" is defined as the number of cells inside, i.e. without buffers.
> >  "unk" is my external array.  If I try to access buffer values, I use:
> > >>>>>
> > >>>>>   for(Int i=0; i<N; i++) {
> > >>>>>     Int gi = mesh.nodes[i].global_number;
> > >>>>>     VecGetValues(x, 1, &gi, &unk[i]);
> > >>>>>   }
> > >>>>>
> > >>>>> But then I end up with tons of warnings, presumably because I am
> > going beyond "n".  Vector x was created with VecCreateGhost.
> > >>>>>
> > >>>>
> > >>>>
> > >>>> --
> > >>>> <Mail Attachment.png>
> > >>>
> > >>>
> > >>> --
> > >>> <Mail Attachment.png>
> > >>>
> > >>
> > >>
> > >> --
> > >> <Mail Attachment.png>
> > >>
> > >
> > >
> > > --
> > > <Signature.png>
> >
> >
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20120224/e7d86e12/attachment-0001.htm>


More information about the petsc-users mailing list