[petsc-users] Accessing Global Vector Data with
Daniel Mckinnell
d_mckinnell at aol.co.uk
Mon Aug 19 09:24:12 CDT 2019
Sorry I'm still a bit confused about this. I changed DMPlexPointLocalRef for DMPlexPointGlobalRef as shown in my attached code but got the following error:
[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation, probably memory access out of range
[0]PETSC ERROR: Try option -start_in_debugger or -on_error_attach_debugger
[0]PETSC ERROR: or see http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind
[0]PETSC ERROR: [1]PETSC ERROR: ------------------------------------------------------------------------
[1]PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation, probably memory access out of range
[1]PETSC ERROR: Try option -start_in_debugger or -on_error_attach_debugger
[1]PETSC ERROR: or see http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind
[1]PETSC ERROR: or try http://valgrind.org on GNU/linux and Apple Mac OS X to find memory corruption errors
or try http://valgrind.org on GNU/linux and Apple Mac OS X to find memory corruption errors
[0]PETSC ERROR: likely location of problem given in stack below
[0]PETSC ERROR: --------------------- Stack Frames ------------------------------------
[0]PETSC ERROR: Note: The EXACT line numbers in the stack are not available,
[0]PETSC ERROR: INSTEAD the line number of the start of the function
[0]PETSC ERROR: is given.
[0]PETSC ERROR: [0] DMPlexPointGlobalRef line 313 /home/daniel/petsc/src/dm/impls/plex/plexpoint.c
[1]PETSC ERROR: likely location of problem given in stack below
[1]PETSC ERROR: --------------------- Stack Frames ------------------------------------
[1]PETSC ERROR: Note: The EXACT line numbers in the stack are not available,
[1]PETSC ERROR: INSTEAD the line number of the start of the function
[1]PETSC ERROR: is given.
[1]PETSC ERROR: [1] DMPlexPointGlobalRef line 313 /home/daniel/petsc/src/dm/impls/plex/plexpoint.c
[1]PETSC ERROR: --------------------------------------------------------------------------
I am guessing that this is because the array is not a global one but I'm not sure how you can produce a global array from a global Vec?Thanks,Daniel
-----Original Message-----
From: Matthew Knepley <knepley at gmail.com>
To: Daniel Mckinnell <d_mckinnell at aol.co.uk>
CC: PETSc <petsc-maint at mcs.anl.gov>
Sent: Mon, 19 Aug 2019 12:13
Subject: Re: [petsc-users] Accessing Global Vector Data with
On Mon, Aug 19, 2019 at 4:46 AM Daniel Mckinnell <d_mckinnell at aol.co.uk> wrote:
Thank you for the help, the "if (voff >= 0)" line was to prevent errors from looping over ghost vertices but it is better to just loop from the start of the vertices to the end of the physical vertices.
I am really just trying to create an example file of a method to access the associated field data of a sample vertex in the global Vec, alter the field data and then set it back into the Vec object.
The easiest thing to do is use
https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexPointGlobalRef.html
Thanks,
Matt
Best Regards,Daniel
-----Original Message-----
From: Matthew Knepley <knepley at gmail.com>
To: Daniel Mckinnell <d_mckinnell at aol.co.uk>
CC: PETSc <petsc-maint at mcs.anl.gov>
Sent: Thu, 15 Aug 2019 17:19
Subject: Re: [petsc-users] Accessing Global Vector Data with
On Thu, Aug 15, 2019 at 12:05 PM Daniel Mckinnell <d_mckinnell at aol.co.uk> wrote:
Thank you, it seemed to work nicely for the fields associated with cells but I'm having problems with the field on vertices. The code for altering the field is as follows:
for (int i = vStart; i < vEnd; i++)
{
PetscScalar *values;
PetscInt vsize, voff;
PetscSectionGetDof(globalsection, i, &vsize);
PetscSectionGetOffset(globalsection, i, &voff);
You do not need the Section calls I think.
VecGetArray(u, &values);
Hoist GetArray out of v loop
if (voff >= 0)
This check is not needed I think. What you are doing here is a mix of local and global.What do you actually want to do?
{
PetscScalar *p;
DMPlexPointLocalRef(*dm, i, values, &p);
p[0] = i + 1;
}
VecRestoreArray(u, &values);
Hoist RestoreArray out of v loop.
Thanks,
Matt
The Global PetscSection is:
PetscSection Object: 2 MPI processes
type not yet set
Process 0:
( 0) dim 2 offset 0
( 1) dim 2 offset 2
( 2) dim -3 offset -13
( 3) dim -3 offset -15
( 4) dim 1 offset 4
( 5) dim 1 offset 5
( 6) dim 1 offset 6
( 7) dim -2 offset -17
( 8) dim -2 offset -18
( 9) dim -2 offset -19
( 10) dim -2 offset -20
( 11) dim -2 offset -21
( 12) dim -2 offset -22
( 13) dim 1 offset 7
( 14) dim 1 offset 8
( 15) dim 1 offset 9
( 16) dim 1 offset 10
( 17) dim 1 offset 11
( 18) dim -2 offset -23
( 19) dim -2 offset -24
( 20) dim -2 offset -25
( 21) dim -2 offset -26
( 22) dim -2 offset -27
( 23) dim -2 offset -28
( 24) dim -2 offset -29
Process 1:
( 0) dim 2 offset 12
( 1) dim 2 offset 14
( 2) dim -3 offset -1
( 3) dim -3 offset -3
( 4) dim 1 offset 16
( 5) dim 1 offset 17
( 6) dim 1 offset 18
( 7) dim 1 offset 19
( 8) dim 1 offset 20
( 9) dim 1 offset 21
( 10) dim -2 offset -5
( 11) dim -2 offset -6
( 12) dim -2 offset -7
( 13) dim 1 offset 22
( 14) dim 1 offset 23
( 15) dim 1 offset 24
( 16) dim 1 offset 25
( 17) dim 1 offset 26
( 18) dim 1 offset 27
( 19) dim 1 offset 28
( 20) dim -2 offset -8
( 21) dim -2 offset -9
( 22) dim -2 offset -10
( 23) dim -2 offset -11
( 24) dim -2 offset -12
where points 4-6 on proc 0 and points 4-9 on proc 1 are vertices.
The altered vec is:
Vec Object: 2 MPI processes
type: mpi
Process [0]
0.
0.
0.
0.
0.
0.
0.
0.
5.
6.
7.
0.
Process [1]
0.
0.
0.
0.
0.
0.
0.
0.
5.
6.
7.
8.
9.
10.
0.
0.
0.
However the altered values do not correspond with the offsets in the PetscSection and I'm not sure what is going wrong.Thanks for all of the help.Daniel
-----Original Message-----
From: Matthew Knepley via petsc-users <petsc-users at mcs.anl.gov>
To: Daniel Mckinnell <d_mckinnell at aol.co.uk>
CC: PETSc <petsc-users at mcs.anl.gov>
Sent: Thu, 15 Aug 2019 16:11
Subject: Re: [petsc-users] Accessing Global Vector Data with
On Thu, Aug 15, 2019 at 10:59 AM Daniel Mckinnell via petsc-users <petsc-users at mcs.anl.gov> wrote:
Hi,Attached is a way I came up with to access data in a global vector, is this the best way to do this or are there other ways? It would seem intuitive to use the global PetscSection and VecGetValuesSection but this doesn't seem to work on global vectors.
Instead I have used VecGetValues and VecSetValues, however I also have a problem with these when extracting more than one value, I initialise the output of VecGetValues as PetscScalar *values; and then call VecGetValues(Other stuff... , values). This seems to work some times and not others and I can't find any rhyme or reason to it?
I guess I should write something about this. I like to think of it as a sort of decision tree.
1) Get just local values, meaning those owned by this process
These can be obtained from either a local or global vector.
2) Get ghosted values, meaning those values lying on unowned points that exist in the default PetscSF
These can only get obtained from a local vector
3) Get arbitrary values
You must use a VecScatter or custom PetscSF to get these
For 1) and 2), I think the best way to get values normally is to use
https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexPointLocalRef.html https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DMPLEX/DMPlexPointGlobalRef.html
These also have Read versions, and Field version to split off a particular field in the Section.
Does this help?
Thanks,
Matt
Finally I was wondering if there is a good reference code base on Github including DMPlex that would be helpful in viewing applications of the DMPlex functions?Thanks,Daniel Mckinnell
--
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
https://www.cse.buffalo.edu/~knepley/
--
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
https://www.cse.buffalo.edu/~knepley/
--
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
https://www.cse.buffalo.edu/~knepley/
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20190819/ec63db9b/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: for_petsc_forum.cpp
Type: text/x-c++src
Size: 2543 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20190819/ec63db9b/attachment-0001.bin>
More information about the petsc-users
mailing list