Actual source code: plexpoint.c
petsc-3.11.4 2019-09-28
1: #include <petsc/private/dmpleximpl.h>
3: /*@
4: DMPlexGetPointLocal - get location of point data in local Vec
6: Not Collective
8: Input Arguments:
9: + dm - DM defining the topological space
10: - point - topological point
12: Output Arguments:
13: + start - start of point data
14: - end - end of point data
16: Note: This is a half open interval [start, end)
18: Level: intermediate
20: .seealso: DMPlexGetPointLocalField(), DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef()
21: @*/
22: PetscErrorCode DMPlexGetPointLocal(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
23: {
24: PetscInt s, e;
31: DMGetLocalOffset_Private(dm, point, &s, &e);
32: if (start) *start = s;
33: if (end) *end = e;
34: return(0);
35: }
37: /*@
38: DMPlexPointLocalRead - return read access to a point in local array
40: Not Collective
42: Input Arguments:
43: + dm - DM defining topological space
44: . point - topological point
45: - array - array to index into
47: Output Arguments:
48: . ptr - address of read reference to point data, type generic so user can place in structure
50: Level: intermediate
52: Note:
53: A common usage when data sizes are known statically:
55: $ const struct { PetscScalar foo,bar,baz; } *ptr;
56: $ DMPlexPointLocalRead(dm,point,array,&ptr);
57: $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
59: .seealso: DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead()
60: @*/
61: PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,void *ptr)
62: {
64: PetscInt start, end;
70: DMGetLocalOffset_Private(dm,point,&start,&end);
71: *(const PetscScalar**)ptr = (start < end) ? array + start : NULL;
72: return(0);
73: }
75: /*@
76: DMPlexPointLocalRef - return read/write access to a point in local array
78: Not Collective
80: Input Arguments:
81: + dm - DM defining topological space
82: . point - topological point
83: - array - array to index into
85: Output Arguments:
86: . ptr - address of reference to point data, type generic so user can place in structure
88: Level: intermediate
90: Note:
91: A common usage when data sizes are known statically:
93: $ struct { PetscScalar foo,bar,baz; } *ptr;
94: $ DMPlexPointLocalRef(dm,point,array,&ptr);
95: $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
97: .seealso: DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
98: @*/
99: PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
100: {
102: PetscInt start, end;
108: DMGetLocalOffset_Private(dm,point,&start,&end);
109: *(PetscScalar**)ptr = (start < end) ? array + start : NULL;
110: return(0);
111: }
113: /*@
114: DMPlexGetPointLocalField - get location of point field data in local Vec
116: Not Collective
118: Input Arguments:
119: + dm - DM defining the topological space
120: . point - topological point
121: - field - the field number
123: Output Arguments:
124: + start - start of point data
125: - end - end of point data
127: Note: This is a half open interval [start, end)
129: Level: intermediate
131: .seealso: DMPlexGetPointLocal(), DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef()
132: @*/
133: PetscErrorCode DMPlexGetPointLocalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
134: {
135: PetscInt s, e;
142: DMGetLocalFieldOffset_Private(dm, point, field, &s, &e);
143: if (start) *start = s;
144: if (end) *end = e;
145: return(0);
146: }
148: /*@
149: DMPlexPointLocalFieldRead - return read access to a field on a point in local array
151: Not Collective
153: Input Arguments:
154: + dm - DM defining topological space
155: . point - topological point
156: . field - field number
157: - array - array to index into
159: Output Arguments:
160: . ptr - address of read reference to point data, type generic so user can place in structure
162: Level: intermediate
164: .seealso: DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
165: @*/
166: PetscErrorCode DMPlexPointLocalFieldRead(DM dm, PetscInt point,PetscInt field,const PetscScalar *array,void *ptr)
167: {
169: PetscInt start, end;
175: DMGetLocalFieldOffset_Private(dm, point, field, &start, &end);
176: *(const PetscScalar**)ptr = array + start;
177: return(0);
178: }
180: /*@
181: DMPlexPointLocalFieldRef - return read/write access to a field on a point in local array
183: Not Collective
185: Input Arguments:
186: + dm - DM defining topological space
187: . point - topological point
188: . field - field number
189: - array - array to index into
191: Output Arguments:
192: . ptr - address of reference to point data, type generic so user can place in structure
194: Level: intermediate
196: .seealso: DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
197: @*/
198: PetscErrorCode DMPlexPointLocalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr)
199: {
201: PetscInt start, end;
207: DMGetLocalFieldOffset_Private(dm, point, field, &start, &end);
208: *(PetscScalar**)ptr = array + start;
209: return(0);
210: }
212: /*@
213: DMPlexGetPointGlobal - get location of point data in global Vec
215: Not Collective
217: Input Arguments:
218: + dm - DM defining the topological space
219: - point - topological point
221: Output Arguments:
222: + start - start of point data; returns -(globalStart+1) if point is not owned
223: - end - end of point data; returns -(globalEnd+1) if point is not owned
225: Note: This is a half open interval [start, end)
227: Level: intermediate
229: .seealso: DMPlexGetPointGlobalField(), DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef()
230: @*/
231: PetscErrorCode DMPlexGetPointGlobal(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
232: {
233: PetscInt s, e;
240: DMGetGlobalOffset_Private(dm, point, &s, &e);
241: if (start) *start = s;
242: if (end) *end = e;
243: return(0);
244: }
246: /*@
247: DMPlexPointGlobalRead - return read access to a point in global array
249: Not Collective
251: Input Arguments:
252: + dm - DM defining topological space
253: . point - topological point
254: - array - array to index into
256: Output Arguments:
257: . ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
259: Level: intermediate
261: Note:
262: A common usage when data sizes are known statically:
264: $ const struct { PetscScalar foo,bar,baz; } *ptr;
265: $ DMPlexPointGlobalRead(dm,point,array,&ptr);
266: $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
268: .seealso: DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef()
269: @*/
270: PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
271: {
272: PetscInt start, end;
279: DMGetGlobalOffset_Private(dm, point, &start, &end);
280: *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL;
281: return(0);
282: }
284: /*@
285: DMPlexPointGlobalRef - return read/write access to a point in global array
287: Not Collective
289: Input Arguments:
290: + dm - DM defining topological space
291: . point - topological point
292: - array - array to index into
294: Output Arguments:
295: . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
297: Level: intermediate
299: Note:
300: A common usage when data sizes are known statically:
302: $ struct { PetscScalar foo,bar,baz; } *ptr;
303: $ DMPlexPointGlobalRef(dm,point,array,&ptr);
304: $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
306: .seealso: DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead()
307: @*/
308: PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
309: {
310: PetscInt start, end;
317: DMGetGlobalOffset_Private(dm, point, &start, &end);
318: *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL;
319: return(0);
320: }
322: /*@
323: DMPlexGetPointGlobalField - get location of point field data in global Vec
325: Not Collective
327: Input Arguments:
328: + dm - DM defining the topological space
329: . point - topological point
330: - field - the field number
332: Output Arguments:
333: + start - start of point data; returns -(globalStart+1) if point is not owned
334: - end - end of point data; returns -(globalEnd+1) if point is not owned
336: Note: This is a half open interval [start, end)
338: Level: intermediate
340: .seealso: DMPlexGetPointGlobal(), DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef()
341: @*/
342: PetscErrorCode DMPlexGetPointGlobalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
343: {
344: PetscInt s, e;
351: DMGetGlobalFieldOffset_Private(dm, point, field, &s, &e);
352: if (start) *start = s;
353: if (end) *end = e;
354: return(0);
355: }
357: /*@
358: DMPlexPointGlobalFieldRead - return read access to a field on a point in global array
360: Not Collective
362: Input Arguments:
363: + dm - DM defining topological space
364: . point - topological point
365: . field - field number
366: - array - array to index into
368: Output Arguments:
369: . ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
371: Level: intermediate
373: .seealso: DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef()
374: @*/
375: PetscErrorCode DMPlexPointGlobalFieldRead(DM dm,PetscInt point,PetscInt field,const PetscScalar *array,void *ptr)
376: {
377: PetscInt start, end;
384: DMGetGlobalFieldOffset_Private(dm, point, field, &start, &end);
385: *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL;
386: return(0);
387: }
389: /*@
390: DMPlexPointGlobalFieldRef - return read/write access to a field on a point in global array
392: Not Collective
394: Input Arguments:
395: + dm - DM defining topological space
396: . point - topological point
397: . field - field number
398: - array - array to index into
400: Output Arguments:
401: . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
403: Level: intermediate
405: .seealso: DMGetSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead()
406: @*/
407: PetscErrorCode DMPlexPointGlobalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr)
408: {
409: PetscInt start, end;
416: DMGetGlobalFieldOffset_Private(dm, point, field, &start, &end);
417: *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL;
418: return(0);
419: }