DMPlexPointGlobalRead#

return read access to a point in global array

Synopsis#

#include "petscdmplex.h"   
PetscErrorCode DMPlexPointGlobalRead(DM dm, PetscInt point, const PetscScalar *array, const void *ptr)

Not Collective

Input Parameters#

  • dm - DM defining topological space

  • point - topological point

  • array - array to index into

Output Parameter#

  • ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned

Note#

A common usage when data sizes are known statically:

  const struct { PetscScalar foo,bar,baz; } *ptr;
  DMPlexPointGlobalRead(dm,point,array,&ptr);
  x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;

See Also#

DMPlex: Unstructured Grids, DM, DMPLEX, DMGetLocalSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef()

Level#

intermediate

Location#

src/dm/impls/plex/plexpoint.c

Examples#

src/ts/tutorials/ex11_sa.c
src/ts/tutorials/ex11.c


Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages