xref: /petsc/src/dm/impls/plex/plexpoint.c (revision 1ce3176fd9e5f4f81eecfcb7b09f4042f482252a)
134541f0dSBarry Smith #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2ad7c26b0SJed Brown #include <petsc-private/isimpl.h>     /* for inline access to atlasOff */
3552f7358SJed Brown 
4552f7358SJed Brown #undef __FUNCT__
5552f7358SJed Brown #define __FUNCT__ "DMPlexGetLocalOffset_Private"
6552f7358SJed Brown PETSC_STATIC_INLINE PetscErrorCode DMPlexGetLocalOffset_Private(DM dm,PetscInt point,PetscInt *offset)
7552f7358SJed Brown {
8552f7358SJed Brown   PetscFunctionBegin;
9552f7358SJed Brown #if defined(PETSC_USE_DEBUG)
10552f7358SJed Brown   {
11552f7358SJed Brown     PetscErrorCode ierr;
12552f7358SJed Brown     ierr = PetscSectionGetOffset(dm->defaultSection,point,offset);CHKERRQ(ierr);
13552f7358SJed Brown   }
14552f7358SJed Brown #else
15552f7358SJed Brown   {
16552f7358SJed Brown     PetscSection s = dm->defaultSection;
174c0d72c2SMatthew G. Knepley     *offset = s->atlasOff[point - s->pStart];
18552f7358SJed Brown   }
19552f7358SJed Brown #endif
20552f7358SJed Brown   PetscFunctionReturn(0);
21552f7358SJed Brown }
22552f7358SJed Brown 
23552f7358SJed Brown #undef __FUNCT__
244824f456SMatthew G. Knepley #define __FUNCT__ "DMPlexGetLocalFieldOffset_Private"
254824f456SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMPlexGetLocalFieldOffset_Private(DM dm,PetscInt point,PetscInt field,PetscInt *offset)
264824f456SMatthew G. Knepley {
274824f456SMatthew G. Knepley   PetscFunctionBegin;
284824f456SMatthew G. Knepley #if defined(PETSC_USE_DEBUG)
294824f456SMatthew G. Knepley   {
304824f456SMatthew G. Knepley     PetscErrorCode ierr;
314824f456SMatthew G. Knepley     ierr = PetscSectionGetFieldOffset(dm->defaultSection,point,field,offset);CHKERRQ(ierr);
324824f456SMatthew G. Knepley   }
334824f456SMatthew G. Knepley #else
344824f456SMatthew G. Knepley   {
354824f456SMatthew G. Knepley     PetscSection s = dm->defaultSection->field[field];
364824f456SMatthew G. Knepley     *offset = s->atlasOff[point - s->pStart];
374824f456SMatthew G. Knepley   }
384824f456SMatthew G. Knepley #endif
394824f456SMatthew G. Knepley   PetscFunctionReturn(0);
404824f456SMatthew G. Knepley }
414824f456SMatthew G. Knepley 
424824f456SMatthew G. Knepley #undef __FUNCT__
43552f7358SJed Brown #define __FUNCT__ "DMPlexGetGlobalOffset_Private"
44552f7358SJed Brown PETSC_STATIC_INLINE PetscErrorCode DMPlexGetGlobalOffset_Private(DM dm,PetscInt point,PetscInt *offset)
45552f7358SJed Brown {
46552f7358SJed Brown   PetscFunctionBegin;
47552f7358SJed Brown #if defined(PETSC_USE_DEBUG)
48552f7358SJed Brown   {
49552f7358SJed Brown     PetscErrorCode ierr;
50552f7358SJed Brown     PetscInt       dof,cdof;
51552f7358SJed Brown     ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,offset);CHKERRQ(ierr);
52552f7358SJed Brown     ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr);
53552f7358SJed Brown     ierr = PetscSectionGetConstraintDof(dm->defaultGlobalSection,point,&cdof);CHKERRQ(ierr);
54552f7358SJed Brown     if (dof-cdof <= 0) *offset = -1; /* Indicates no data */
55552f7358SJed Brown   }
56552f7358SJed Brown #else
57552f7358SJed Brown   {
58552f7358SJed Brown     PetscSection s = dm->defaultGlobalSection;
59552f7358SJed Brown     PetscInt     dof,cdof;
604c0d72c2SMatthew G. Knepley     *offset = s->atlasOff[point - s->pStart];
614c0d72c2SMatthew G. Knepley     dof     = s->atlasDof[point - s->pStart];
624c0d72c2SMatthew G. Knepley     cdof    = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
63552f7358SJed Brown     if (dof-cdof <= 0) *offset = -1;
64552f7358SJed Brown   }
65552f7358SJed Brown #endif
66552f7358SJed Brown   PetscFunctionReturn(0);
67552f7358SJed Brown }
68552f7358SJed Brown 
69552f7358SJed Brown #undef __FUNCT__
70552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointLocal"
71552f7358SJed Brown /*@
72552f7358SJed Brown    DMPlexGetPointLocal - get location of point data in local Vec
73552f7358SJed Brown 
74552f7358SJed Brown    Not Collective
75552f7358SJed Brown 
76552f7358SJed Brown    Input Arguments:
77552f7358SJed Brown +  dm - DM defining the topological space
78552f7358SJed Brown -  point - topological point
79552f7358SJed Brown 
80552f7358SJed Brown    Output Arguments:
81552f7358SJed Brown +  start - start of point data
82552f7358SJed Brown -  end - end of point data
83552f7358SJed Brown 
84552f7358SJed Brown    Level: intermediate
85552f7358SJed Brown 
86552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef()
87552f7358SJed Brown @*/
88552f7358SJed Brown PetscErrorCode DMPlexGetPointLocal(DM dm,PetscInt point,PetscInt *start,PetscInt *end)
89552f7358SJed Brown {
90552f7358SJed Brown   PetscErrorCode ierr;
91552f7358SJed Brown   PetscInt       offset,dof;
92552f7358SJed Brown 
93552f7358SJed Brown   PetscFunctionBegin;
94552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
95552f7358SJed Brown   ierr = PetscSectionGetOffset(dm->defaultSection,point,&offset);CHKERRQ(ierr);
96552f7358SJed Brown   ierr = PetscSectionGetDof(dm->defaultSection,point,&dof);CHKERRQ(ierr);
97552f7358SJed Brown   if (start) *start = offset;
98552f7358SJed Brown   if (end) *end = offset + dof;
99552f7358SJed Brown   PetscFunctionReturn(0);
100552f7358SJed Brown }
101552f7358SJed Brown 
102552f7358SJed Brown #undef __FUNCT__
103552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRead"
104552f7358SJed Brown /*@
105552f7358SJed Brown    DMPlexPointLocalRead - return read access to a point in local array
106552f7358SJed Brown 
107552f7358SJed Brown    Not Collective
108552f7358SJed Brown 
109552f7358SJed Brown    Input Arguments:
110552f7358SJed Brown +  dm - DM defining topological space
111552f7358SJed Brown .  point - topological point
112552f7358SJed Brown -  array - array to index into
113552f7358SJed Brown 
114552f7358SJed Brown    Output Arguments:
115552f7358SJed Brown .  ptr - address of read reference to point data, type generic so user can place in structure
116552f7358SJed Brown 
117552f7358SJed Brown    Level: intermediate
118552f7358SJed Brown 
119552f7358SJed Brown    Note:
120552f7358SJed Brown    A common usage when data sizes are known statically:
121552f7358SJed Brown 
122552f7358SJed Brown $  const struct { PetscScalar foo,bar,baz; } *ptr;
123552f7358SJed Brown $  DMPlexPointLocalRead(dm,point,array,&ptr);
124552f7358SJed Brown $  x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
125552f7358SJed Brown 
126552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead()
127552f7358SJed Brown @*/
128552f7358SJed Brown PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
129552f7358SJed Brown {
130552f7358SJed Brown   PetscErrorCode ierr;
131552f7358SJed Brown   PetscInt       start;
132552f7358SJed Brown 
133552f7358SJed Brown   PetscFunctionBegin;
134552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
135552f7358SJed Brown   PetscValidScalarPointer(array,3);
136552f7358SJed Brown   PetscValidPointer(ptr,4);
137552f7358SJed Brown   ierr                      = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr);
138552f7358SJed Brown   *(const PetscScalar**)ptr = array + start;
139552f7358SJed Brown   PetscFunctionReturn(0);
140552f7358SJed Brown }
141552f7358SJed Brown 
142552f7358SJed Brown #undef __FUNCT__
143552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRef"
144552f7358SJed Brown /*@
145552f7358SJed Brown    DMPlexPointLocalRef - return read/write access to a point in local array
146552f7358SJed Brown 
147552f7358SJed Brown    Not Collective
148552f7358SJed Brown 
149552f7358SJed Brown    Input Arguments:
150552f7358SJed Brown +  dm - DM defining topological space
151552f7358SJed Brown .  point - topological point
152552f7358SJed Brown -  array - array to index into
153552f7358SJed Brown 
154552f7358SJed Brown    Output Arguments:
155552f7358SJed Brown .  ptr - address of reference to point data, type generic so user can place in structure
156552f7358SJed Brown 
157552f7358SJed Brown    Level: intermediate
158552f7358SJed Brown 
159552f7358SJed Brown    Note:
160552f7358SJed Brown    A common usage when data sizes are known statically:
161552f7358SJed Brown 
162552f7358SJed Brown $  struct { PetscScalar foo,bar,baz; } *ptr;
163552f7358SJed Brown $  DMPlexPointLocalRef(dm,point,array,&ptr);
164552f7358SJed Brown $  ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
165552f7358SJed Brown 
166552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
167552f7358SJed Brown @*/
168552f7358SJed Brown PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
169552f7358SJed Brown {
170552f7358SJed Brown   PetscErrorCode ierr;
171552f7358SJed Brown   PetscInt       start;
172552f7358SJed Brown 
173552f7358SJed Brown   PetscFunctionBegin;
174552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
175552f7358SJed Brown   PetscValidScalarPointer(array,3);
176552f7358SJed Brown   PetscValidPointer(ptr,4);
177552f7358SJed Brown   ierr                = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr);
178552f7358SJed Brown   *(PetscScalar**)ptr = array + start;
179552f7358SJed Brown   PetscFunctionReturn(0);
180552f7358SJed Brown }
181552f7358SJed Brown 
182552f7358SJed Brown #undef __FUNCT__
183*1ce3176fSMatthew G. Knepley #define __FUNCT__ "DMPlexPointLocalFieldRead"
184*1ce3176fSMatthew G. Knepley /*@
185*1ce3176fSMatthew G. Knepley    DMPlexPointLocalFieldRead - return read access to a field on a point in local array
186*1ce3176fSMatthew G. Knepley 
187*1ce3176fSMatthew G. Knepley    Not Collective
188*1ce3176fSMatthew G. Knepley 
189*1ce3176fSMatthew G. Knepley    Input Arguments:
190*1ce3176fSMatthew G. Knepley +  dm - DM defining topological space
191*1ce3176fSMatthew G. Knepley .  point - topological point
192*1ce3176fSMatthew G. Knepley .  field - field number
193*1ce3176fSMatthew G. Knepley -  array - array to index into
194*1ce3176fSMatthew G. Knepley 
195*1ce3176fSMatthew G. Knepley    Output Arguments:
196*1ce3176fSMatthew G. Knepley .  ptr - address of read reference to point data, type generic so user can place in structure
197*1ce3176fSMatthew G. Knepley 
198*1ce3176fSMatthew G. Knepley    Level: intermediate
199*1ce3176fSMatthew G. Knepley 
200*1ce3176fSMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
201*1ce3176fSMatthew G. Knepley @*/
202*1ce3176fSMatthew G. Knepley PetscErrorCode DMPlexPointLocalFieldRead(DM dm,PetscInt point,PetscInt field,const PetscScalar *array,const void *ptr)
203*1ce3176fSMatthew G. Knepley {
204*1ce3176fSMatthew G. Knepley   PetscErrorCode ierr;
205*1ce3176fSMatthew G. Knepley   PetscInt       start;
206*1ce3176fSMatthew G. Knepley 
207*1ce3176fSMatthew G. Knepley   PetscFunctionBegin;
208*1ce3176fSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
209*1ce3176fSMatthew G. Knepley   PetscValidScalarPointer(array,3);
210*1ce3176fSMatthew G. Knepley   PetscValidPointer(ptr,4);
211*1ce3176fSMatthew G. Knepley   ierr                      = DMPlexGetLocalFieldOffset_Private(dm,point,field,&start);CHKERRQ(ierr);
212*1ce3176fSMatthew G. Knepley   *(const PetscScalar**)ptr = array + start;
213*1ce3176fSMatthew G. Knepley   PetscFunctionReturn(0);
214*1ce3176fSMatthew G. Knepley }
215*1ce3176fSMatthew G. Knepley 
216*1ce3176fSMatthew G. Knepley #undef __FUNCT__
2174824f456SMatthew G. Knepley #define __FUNCT__ "DMPlexPointLocalFieldRef"
2184824f456SMatthew G. Knepley /*@
2194824f456SMatthew G. Knepley    DMPlexPointLocalFieldRef - return read/write access to a field on a point in local array
2204824f456SMatthew G. Knepley 
2214824f456SMatthew G. Knepley    Not Collective
2224824f456SMatthew G. Knepley 
2234824f456SMatthew G. Knepley    Input Arguments:
2244824f456SMatthew G. Knepley +  dm - DM defining topological space
2254824f456SMatthew G. Knepley .  point - topological point
2264824f456SMatthew G. Knepley .  field - field number
2274824f456SMatthew G. Knepley -  array - array to index into
2284824f456SMatthew G. Knepley 
2294824f456SMatthew G. Knepley    Output Arguments:
2304824f456SMatthew G. Knepley .  ptr - address of reference to point data, type generic so user can place in structure
2314824f456SMatthew G. Knepley 
2324824f456SMatthew G. Knepley    Level: intermediate
2334824f456SMatthew G. Knepley 
2344824f456SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
2354824f456SMatthew G. Knepley @*/
2364824f456SMatthew G. Knepley PetscErrorCode DMPlexPointLocalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr)
2374824f456SMatthew G. Knepley {
2384824f456SMatthew G. Knepley   PetscErrorCode ierr;
2394824f456SMatthew G. Knepley   PetscInt       start;
2404824f456SMatthew G. Knepley 
2414824f456SMatthew G. Knepley   PetscFunctionBegin;
2424824f456SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2434824f456SMatthew G. Knepley   PetscValidScalarPointer(array,3);
2444824f456SMatthew G. Knepley   PetscValidPointer(ptr,4);
2454824f456SMatthew G. Knepley   ierr                = DMPlexGetLocalFieldOffset_Private(dm,point,field,&start);CHKERRQ(ierr);
2464824f456SMatthew G. Knepley   *(PetscScalar**)ptr = array + start;
2474824f456SMatthew G. Knepley   PetscFunctionReturn(0);
2484824f456SMatthew G. Knepley }
2494824f456SMatthew G. Knepley 
2504824f456SMatthew G. Knepley #undef __FUNCT__
251552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointGlobal"
252552f7358SJed Brown /*@
253552f7358SJed Brown    DMPlexGetPointGlobal - get location of point data in global Vec
254552f7358SJed Brown 
255552f7358SJed Brown    Not Collective
256552f7358SJed Brown 
257552f7358SJed Brown    Input Arguments:
258552f7358SJed Brown +  dm - DM defining the topological space
259552f7358SJed Brown -  point - topological point
260552f7358SJed Brown 
261552f7358SJed Brown    Output Arguments:
262552f7358SJed Brown +  start - start of point data; returns -(global_start+1) if point is not owned
263552f7358SJed Brown -  end - end of point data; returns -(global_end+1) if point is not owned
264552f7358SJed Brown 
265552f7358SJed Brown    Level: intermediate
266552f7358SJed Brown 
267552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef()
268552f7358SJed Brown @*/
269552f7358SJed Brown PetscErrorCode DMPlexGetPointGlobal(DM dm,PetscInt point,PetscInt *start,PetscInt *end)
270552f7358SJed Brown {
271552f7358SJed Brown   PetscErrorCode ierr;
272552f7358SJed Brown   PetscInt       offset,dof;
273552f7358SJed Brown 
274552f7358SJed Brown   PetscFunctionBegin;
275552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
276552f7358SJed Brown   ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,&offset);CHKERRQ(ierr);
277552f7358SJed Brown   ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr);
278552f7358SJed Brown   if (start) *start = offset;
279552f7358SJed Brown   if (end) *end = offset + dof;
280552f7358SJed Brown   PetscFunctionReturn(0);
281552f7358SJed Brown }
282552f7358SJed Brown 
283552f7358SJed Brown #undef __FUNCT__
284552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRead"
285552f7358SJed Brown /*@
286552f7358SJed Brown    DMPlexPointGlobalRead - return read access to a point in global array
287552f7358SJed Brown 
288552f7358SJed Brown    Not Collective
289552f7358SJed Brown 
290552f7358SJed Brown    Input Arguments:
291552f7358SJed Brown +  dm - DM defining topological space
292552f7358SJed Brown .  point - topological point
293552f7358SJed Brown -  array - array to index into
294552f7358SJed Brown 
295552f7358SJed Brown    Output Arguments:
2960298fd71SBarry Smith .  ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
297552f7358SJed Brown 
298552f7358SJed Brown    Level: intermediate
299552f7358SJed Brown 
300552f7358SJed Brown    Note:
301552f7358SJed Brown    A common usage when data sizes are known statically:
302552f7358SJed Brown 
303552f7358SJed Brown $  const struct { PetscScalar foo,bar,baz; } *ptr;
304552f7358SJed Brown $  DMPlexPointGlobalRead(dm,point,array,&ptr);
305552f7358SJed Brown $  x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
306552f7358SJed Brown 
307552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef()
308552f7358SJed Brown @*/
309552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
310552f7358SJed Brown {
311552f7358SJed Brown   PetscErrorCode ierr;
312552f7358SJed Brown   PetscInt       start;
313552f7358SJed Brown 
314552f7358SJed Brown   PetscFunctionBegin;
315552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
316552f7358SJed Brown   PetscValidScalarPointer(array,3);
317552f7358SJed Brown   PetscValidPointer(ptr,4);
318552f7358SJed Brown   ierr                      = DMPlexGetGlobalOffset_Private(dm,point,&start);CHKERRQ(ierr);
3190298fd71SBarry Smith   *(const PetscScalar**)ptr = (start >= 0) ? array + start - dm->map->rstart : NULL;
320552f7358SJed Brown   PetscFunctionReturn(0);
321552f7358SJed Brown }
322552f7358SJed Brown 
323552f7358SJed Brown #undef __FUNCT__
324552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRef"
325552f7358SJed Brown /*@
326552f7358SJed Brown    DMPlexPointGlobalRef - return read/write access to a point in global array
327552f7358SJed Brown 
328552f7358SJed Brown    Not Collective
329552f7358SJed Brown 
330552f7358SJed Brown    Input Arguments:
331552f7358SJed Brown +  dm - DM defining topological space
332552f7358SJed Brown .  point - topological point
333552f7358SJed Brown -  array - array to index into
334552f7358SJed Brown 
335552f7358SJed Brown    Output Arguments:
3360298fd71SBarry Smith .  ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
337552f7358SJed Brown 
338552f7358SJed Brown    Level: intermediate
339552f7358SJed Brown 
340552f7358SJed Brown    Note:
341552f7358SJed Brown    A common usage when data sizes are known statically:
342552f7358SJed Brown 
343552f7358SJed Brown $  struct { PetscScalar foo,bar,baz; } *ptr;
344552f7358SJed Brown $  DMPlexPointGlobalRef(dm,point,array,&ptr);
345552f7358SJed Brown $  ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
346552f7358SJed Brown 
347552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead()
348552f7358SJed Brown @*/
349552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
350552f7358SJed Brown {
351552f7358SJed Brown   PetscErrorCode ierr;
352552f7358SJed Brown   PetscInt       start;
353552f7358SJed Brown 
354552f7358SJed Brown   PetscFunctionBegin;
355552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
356552f7358SJed Brown   PetscValidScalarPointer(array,3);
357552f7358SJed Brown   PetscValidPointer(ptr,4);
358552f7358SJed Brown   ierr                = DMPlexGetGlobalOffset_Private(dm,point,&start);CHKERRQ(ierr);
3590298fd71SBarry Smith   *(PetscScalar**)ptr = (start >= 0) ? array + start - dm->map->rstart : NULL;
360552f7358SJed Brown   PetscFunctionReturn(0);
361552f7358SJed Brown }
362