Actual source code: plexpoint.c
petsc-3.6.1 2015-08-06
1: #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
5: /*@
6: DMPlexGetPointLocal - get location of point data in local Vec
8: Not Collective
10: Input Arguments:
11: + dm - DM defining the topological space
12: - point - topological point
14: Output Arguments:
15: + start - start of point data
16: - end - end of point data
18: Note: This is a half open interval [start, end)
20: Level: intermediate
22: .seealso: DMPlexGetPointLocalField(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef()
23: @*/
24: PetscErrorCode DMPlexGetPointLocal(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
25: {
26: PetscInt s, e;
33: DMPlexGetLocalOffset_Private(dm, point, &s, &e);
34: if (start) *start = s;
35: if (end) *end = e;
36: return(0);
37: }
41: /*@
42: DMPlexPointLocalRead - return read access to a point in local array
44: Not Collective
46: Input Arguments:
47: + dm - DM defining topological space
48: . point - topological point
49: - array - array to index into
51: Output Arguments:
52: . ptr - address of read reference to point data, type generic so user can place in structure
54: Level: intermediate
56: Note:
57: A common usage when data sizes are known statically:
59: $ const struct { PetscScalar foo,bar,baz; } *ptr;
60: $ DMPlexPointLocalRead(dm,point,array,&ptr);
61: $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
63: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead()
64: @*/
65: PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
66: {
68: PetscInt start, end;
74: DMPlexGetLocalOffset_Private(dm,point,&start,&end);
75: *(const PetscScalar**)ptr = (start < end) ? array + start : NULL;
76: return(0);
77: }
81: /*@
82: DMPlexPointLocalRef - return read/write access to a point in local array
84: Not Collective
86: Input Arguments:
87: + dm - DM defining topological space
88: . point - topological point
89: - array - array to index into
91: Output Arguments:
92: . ptr - address of reference to point data, type generic so user can place in structure
94: Level: intermediate
96: Note:
97: A common usage when data sizes are known statically:
99: $ struct { PetscScalar foo,bar,baz; } *ptr;
100: $ DMPlexPointLocalRef(dm,point,array,&ptr);
101: $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
103: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
104: @*/
105: PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
106: {
108: PetscInt start, end;
114: DMPlexGetLocalOffset_Private(dm,point,&start,&end);
115: *(PetscScalar**)ptr = (start < end) ? array + start : NULL;
116: return(0);
117: }
121: /*@
122: DMPlexGetPointLocalField - get location of point field data in local Vec
124: Not Collective
126: Input Arguments:
127: + dm - DM defining the topological space
128: . point - topological point
129: - field - the field number
131: Output Arguments:
132: + start - start of point data
133: - end - end of point data
135: Note: This is a half open interval [start, end)
137: Level: intermediate
139: .seealso: DMPlexGetPointLocal(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef()
140: @*/
141: PetscErrorCode DMPlexGetPointLocalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
142: {
143: PetscInt s, e;
150: DMPlexGetLocalFieldOffset_Private(dm, point, field, &s, &e);
151: if (start) *start = s;
152: if (end) *end = e;
153: return(0);
154: }
158: /*@
159: DMPlexPointLocalFieldRead - return read access to a field on a point in local array
161: Not Collective
163: Input Arguments:
164: + dm - DM defining topological space
165: . point - topological point
166: . field - field number
167: - array - array to index into
169: Output Arguments:
170: . ptr - address of read reference to point data, type generic so user can place in structure
172: Level: intermediate
174: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
175: @*/
176: PetscErrorCode DMPlexPointLocalFieldRead(DM dm, PetscInt point,PetscInt field,const PetscScalar *array,const void *ptr)
177: {
179: PetscInt start, end;
185: DMPlexGetLocalFieldOffset_Private(dm, point, field, &start, &end);
186: *(const PetscScalar**)ptr = array + start;
187: return(0);
188: }
192: /*@
193: DMPlexPointLocalFieldRef - return read/write access to a field on a point in local array
195: Not Collective
197: Input Arguments:
198: + dm - DM defining topological space
199: . point - topological point
200: . field - field number
201: - array - array to index into
203: Output Arguments:
204: . ptr - address of reference to point data, type generic so user can place in structure
206: Level: intermediate
208: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
209: @*/
210: PetscErrorCode DMPlexPointLocalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr)
211: {
213: PetscInt start, end;
219: DMPlexGetLocalFieldOffset_Private(dm, point, field, &start, &end);
220: *(PetscScalar**)ptr = array + start;
221: return(0);
222: }
226: /*@
227: DMPlexGetPointGlobal - get location of point data in global Vec
229: Not Collective
231: Input Arguments:
232: + dm - DM defining the topological space
233: - point - topological point
235: Output Arguments:
236: + start - start of point data; returns -(globalStart+1) if point is not owned
237: - end - end of point data; returns -(globalEnd+1) if point is not owned
239: Note: This is a half open interval [start, end)
241: Level: intermediate
243: .seealso: DMPlexGetPointGlobalField(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef()
244: @*/
245: PetscErrorCode DMPlexGetPointGlobal(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
246: {
247: PetscInt s, e;
254: DMPlexGetGlobalOffset_Private(dm, point, &s, &e);
255: if (start) *start = s;
256: if (end) *end = e;
257: return(0);
258: }
262: /*@
263: DMPlexPointGlobalRead - return read access to a point in global array
265: Not Collective
267: Input Arguments:
268: + dm - DM defining topological space
269: . point - topological point
270: - array - array to index into
272: Output Arguments:
273: . ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
275: Level: intermediate
277: Note:
278: A common usage when data sizes are known statically:
280: $ const struct { PetscScalar foo,bar,baz; } *ptr;
281: $ DMPlexPointGlobalRead(dm,point,array,&ptr);
282: $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
284: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef()
285: @*/
286: PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
287: {
288: PetscInt start, end;
295: DMPlexGetGlobalOffset_Private(dm, point, &start, &end);
296: *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL;
297: return(0);
298: }
302: /*@
303: DMPlexPointGlobalRef - return read/write access to a point in global array
305: Not Collective
307: Input Arguments:
308: + dm - DM defining topological space
309: . point - topological point
310: - array - array to index into
312: Output Arguments:
313: . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
315: Level: intermediate
317: Note:
318: A common usage when data sizes are known statically:
320: $ struct { PetscScalar foo,bar,baz; } *ptr;
321: $ DMPlexPointGlobalRef(dm,point,array,&ptr);
322: $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
324: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead()
325: @*/
326: PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
327: {
328: PetscInt start, end;
335: DMPlexGetGlobalOffset_Private(dm, point, &start, &end);
336: *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL;
337: return(0);
338: }
342: /*@
343: DMPlexGetPointGlobalField - get location of point field data in global Vec
345: Not Collective
347: Input Arguments:
348: + dm - DM defining the topological space
349: . point - topological point
350: - field - the field number
352: Output Arguments:
353: + start - start of point data; returns -(globalStart+1) if point is not owned
354: - end - end of point data; returns -(globalEnd+1) if point is not owned
356: Note: This is a half open interval [start, end)
358: Level: intermediate
360: .seealso: DMPlexGetPointGlobal(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef()
361: @*/
362: PetscErrorCode DMPlexGetPointGlobalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
363: {
364: PetscInt s, e;
371: DMPlexGetGlobalFieldOffset_Private(dm, point, field, &s, &e);
372: if (start) *start = s;
373: if (end) *end = e;
374: return(0);
375: }
379: /*@
380: DMPlexPointGlobalFieldRead - return read access to a field on a point in global array
382: Not Collective
384: Input Arguments:
385: + dm - DM defining topological space
386: . point - topological point
387: . field - field number
388: - array - array to index into
390: Output Arguments:
391: . ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
393: Level: intermediate
395: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef()
396: @*/
397: PetscErrorCode DMPlexPointGlobalFieldRead(DM dm,PetscInt point,PetscInt field,const PetscScalar *array,const void *ptr)
398: {
399: PetscInt start, end;
406: DMPlexGetGlobalFieldOffset_Private(dm, point, field, &start, &end);
407: *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL;
408: return(0);
409: }
413: /*@
414: DMPlexPointGlobalFieldRef - return read/write access to a field on a point in global array
416: Not Collective
418: Input Arguments:
419: + dm - DM defining topological space
420: . point - topological point
421: . field - field number
422: - array - array to index into
424: Output Arguments:
425: . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
427: Level: intermediate
429: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead()
430: @*/
431: PetscErrorCode DMPlexPointGlobalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr)
432: {
433: PetscInt start, end;
440: DMPlexGetGlobalFieldOffset_Private(dm, point, field, &start, &end);
441: *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL;
442: return(0);
443: }