xref: /petsc/src/dm/impls/plex/plexpoint.c (revision 661006243e0904bc232e219bc5c02c8d300d1d6d)
1 #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2 
3 /*@
4    DMPlexGetPointLocal - get location of point data in local Vec
5 
6    Not Collective
7 
8    Input Parameters:
9 +  dm - DM defining the topological space
10 -  point - topological point
11 
12    Output Parameters:
13 +  start - start of point data
14 -  end - end of point data
15 
16    Note: This is a half open interval [start, end)
17 
18    Level: intermediate
19 
20 .seealso: `DMPlexGetPointLocalField()`, `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexPointLocalRead()`, `DMPlexPointLocalRead()`, `DMPlexPointLocalRef()`
21 @*/
22 PetscErrorCode DMPlexGetPointLocal(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
23 {
24   PetscInt s, e;
25 
26   PetscFunctionBegin;
27   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28   if (start) PetscValidIntPointer(start, 3);
29   if (end) PetscValidIntPointer(end, 4);
30   PetscCall(DMGetLocalOffset_Private(dm, point, &s, &e));
31   if (start) *start = s;
32   if (end) *end = e;
33   PetscFunctionReturn(0);
34 }
35 
36 /*@
37    DMPlexPointLocalRead - return read access to a point in local array
38 
39    Not Collective
40 
41    Input Parameters:
42 +  dm - DM defining topological space
43 .  point - topological point
44 -  array - array to index into
45 
46    Output Parameter:
47 .  ptr - address of read reference to point data, type generic so user can place in structure
48 
49    Level: intermediate
50 
51    Note:
52    A common usage when data sizes are known statically:
53 
54 $  const struct { PetscScalar foo,bar,baz; } *ptr;
55 $  DMPlexPointLocalRead(dm,point,array,&ptr);
56 $  x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
57 
58 .seealso: `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexGetPointLocal()`, `DMPlexPointGlobalRead()`
59 @*/
60 PetscErrorCode DMPlexPointLocalRead(DM dm, PetscInt point, const PetscScalar *array, void *ptr)
61 {
62   PetscInt start, end;
63 
64   PetscFunctionBegin;
65   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
66   PetscValidScalarPointer(array, 3);
67   PetscValidPointer(ptr, 4);
68   PetscCall(DMGetLocalOffset_Private(dm, point, &start, &end));
69   *(const PetscScalar **)ptr = (start < end) ? array + start : NULL;
70   PetscFunctionReturn(0);
71 }
72 
73 /*@
74    DMPlexPointLocalRef - return read/write access to a point in local array
75 
76    Not Collective
77 
78    Input Parameters:
79 +  dm - DM defining topological space
80 .  point - topological point
81 -  array - array to index into
82 
83    Output Parameter:
84 .  ptr - address of reference to point data, type generic so user can place in structure
85 
86    Level: intermediate
87 
88    Note:
89    A common usage when data sizes are known statically:
90 
91 $  struct { PetscScalar foo,bar,baz; } *ptr;
92 $  DMPlexPointLocalRef(dm,point,array,&ptr);
93 $  ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
94 
95 .seealso: `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexGetPointLocal()`, `DMPlexPointGlobalRef()`
96 @*/
97 PetscErrorCode DMPlexPointLocalRef(DM dm, PetscInt point, PetscScalar *array, void *ptr)
98 {
99   PetscInt start, end;
100 
101   PetscFunctionBegin;
102   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
103   PetscValidScalarPointer(array, 3);
104   PetscValidPointer(ptr, 4);
105   PetscCall(DMGetLocalOffset_Private(dm, point, &start, &end));
106   *(PetscScalar **)ptr = (start < end) ? array + start : NULL;
107   PetscFunctionReturn(0);
108 }
109 
110 /*@
111   DMPlexGetPointLocalField - get location of point field data in local Vec
112 
113   Not Collective
114 
115   Input Parameters:
116 + dm - DM defining the topological space
117 . point - topological point
118 - field - the field number
119 
120   Output Parameters:
121 + start - start of point data
122 - end - end of point data
123 
124   Note: This is a half open interval [start, end)
125 
126   Level: intermediate
127 
128 .seealso: `DMPlexGetPointLocal()`, `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexPointLocalRead()`, `DMPlexPointLocalRead()`, `DMPlexPointLocalRef()`
129 @*/
130 PetscErrorCode DMPlexGetPointLocalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
131 {
132   PetscInt s, e;
133 
134   PetscFunctionBegin;
135   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
136   if (start) PetscValidIntPointer(start, 4);
137   if (end) PetscValidIntPointer(end, 5);
138   PetscCall(DMGetLocalFieldOffset_Private(dm, point, field, &s, &e));
139   if (start) *start = s;
140   if (end) *end = e;
141   PetscFunctionReturn(0);
142 }
143 
144 /*@
145    DMPlexPointLocalFieldRead - return read access to a field on a point in local array
146 
147    Not Collective
148 
149    Input Parameters:
150 +  dm - DM defining topological space
151 .  point - topological point
152 .  field - field number
153 -  array - array to index into
154 
155    Output Parameter:
156 .  ptr - address of read reference to point data, type generic so user can place in structure
157 
158    Level: intermediate
159 
160 .seealso: `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexGetPointLocal()`, `DMPlexPointGlobalRef()`
161 @*/
162 PetscErrorCode DMPlexPointLocalFieldRead(DM dm, PetscInt point, PetscInt field, const PetscScalar *array, void *ptr)
163 {
164   PetscInt start, end;
165 
166   PetscFunctionBegin;
167   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
168   PetscValidScalarPointer(array, 4);
169   PetscValidPointer(ptr, 5);
170   PetscCall(DMGetLocalFieldOffset_Private(dm, point, field, &start, &end));
171   *(const PetscScalar **)ptr = array + start;
172   PetscFunctionReturn(0);
173 }
174 
175 /*@
176    DMPlexPointLocalFieldRef - return read/write access to a field on a point in local array
177 
178    Not Collective
179 
180    Input Parameters:
181 +  dm - DM defining topological space
182 .  point - topological point
183 .  field - field number
184 -  array - array to index into
185 
186    Output Parameter:
187 .  ptr - address of reference to point data, type generic so user can place in structure
188 
189    Level: intermediate
190 
191 .seealso: `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexGetPointLocal()`, `DMPlexPointGlobalRef()`
192 @*/
193 PetscErrorCode DMPlexPointLocalFieldRef(DM dm, PetscInt point, PetscInt field, PetscScalar *array, void *ptr)
194 {
195   PetscInt start, end;
196 
197   PetscFunctionBegin;
198   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
199   PetscValidScalarPointer(array, 4);
200   PetscValidPointer(ptr, 5);
201   PetscCall(DMGetLocalFieldOffset_Private(dm, point, field, &start, &end));
202   *(PetscScalar **)ptr = array + start;
203   PetscFunctionReturn(0);
204 }
205 
206 /*@
207   DMPlexGetPointGlobal - get location of point data in global Vec
208 
209   Not Collective
210 
211   Input Parameters:
212 + dm - DM defining the topological space
213 - point - topological point
214 
215   Output Parameters:
216 + start - start of point data; returns -(globalStart+1) if point is not owned
217 - end - end of point data; returns -(globalEnd+1) if point is not owned
218 
219   Note: This is a half open interval [start, end)
220 
221   Level: intermediate
222 
223 .seealso: `DMPlexGetPointGlobalField()`, `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexPointGlobalRead()`, `DMPlexGetPointLocal()`, `DMPlexPointGlobalRead()`, `DMPlexPointGlobalRef()`
224 @*/
225 PetscErrorCode DMPlexGetPointGlobal(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
226 {
227   PetscInt s, e;
228 
229   PetscFunctionBegin;
230   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
231   if (start) PetscValidIntPointer(start, 3);
232   if (end) PetscValidIntPointer(end, 4);
233   PetscCall(DMGetGlobalOffset_Private(dm, point, &s, &e));
234   if (start) *start = s;
235   if (end) *end = e;
236   PetscFunctionReturn(0);
237 }
238 
239 /*@
240    DMPlexPointGlobalRead - return read access to a point in global array
241 
242    Not Collective
243 
244    Input Parameters:
245 +  dm - DM defining topological space
246 .  point - topological point
247 -  array - array to index into
248 
249    Output Parameter:
250 .  ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
251 
252    Level: intermediate
253 
254    Note:
255    A common usage when data sizes are known statically:
256 
257 $  const struct { PetscScalar foo,bar,baz; } *ptr;
258 $  DMPlexPointGlobalRead(dm,point,array,&ptr);
259 $  x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
260 
261 .seealso: `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexGetPointGlobal()`, `DMPlexPointLocalRead()`, `DMPlexPointGlobalRef()`
262 @*/
263 PetscErrorCode DMPlexPointGlobalRead(DM dm, PetscInt point, const PetscScalar *array, const void *ptr)
264 {
265   PetscInt start, end;
266 
267   PetscFunctionBegin;
268   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
269   PetscValidScalarPointer(array, 3);
270   PetscValidPointer(ptr, 4);
271   PetscCall(DMGetGlobalOffset_Private(dm, point, &start, &end));
272   *(const PetscScalar **)ptr = (start < end) ? array + start - dm->map->rstart : NULL;
273   PetscFunctionReturn(0);
274 }
275 
276 /*@
277    DMPlexPointGlobalRef - return read/write access to a point in global array
278 
279    Not Collective
280 
281    Input Parameters:
282 +  dm - DM defining topological space
283 .  point - topological point
284 -  array - array to index into
285 
286    Output Parameter:
287 .  ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
288 
289    Level: intermediate
290 
291    Note:
292    A common usage when data sizes are known statically:
293 
294 $  struct { PetscScalar foo,bar,baz; } *ptr;
295 $  DMPlexPointGlobalRef(dm,point,array,&ptr);
296 $  ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
297 
298 .seealso: `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexGetPointGlobal()`, `DMPlexPointLocalRef()`, `DMPlexPointGlobalRead()`
299 @*/
300 PetscErrorCode DMPlexPointGlobalRef(DM dm, PetscInt point, PetscScalar *array, void *ptr)
301 {
302   PetscInt start, end;
303 
304   PetscFunctionBegin;
305   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
306   PetscValidScalarPointer(array, 3);
307   PetscValidPointer(ptr, 4);
308   PetscCall(DMGetGlobalOffset_Private(dm, point, &start, &end));
309   *(PetscScalar **)ptr = (start < end) ? array + start - dm->map->rstart : NULL;
310   PetscFunctionReturn(0);
311 }
312 
313 /*@
314   DMPlexGetPointGlobalField - get location of point field data in global Vec
315 
316   Not Collective
317 
318   Input Parameters:
319 + dm - DM defining the topological space
320 . point - topological point
321 - field - the field number
322 
323   Output Parameters:
324 + start - start of point data; returns -(globalStart+1) if point is not owned
325 - end - end of point data; returns -(globalEnd+1) if point is not owned
326 
327   Note: This is a half open interval [start, end)
328 
329   Level: intermediate
330 
331 .seealso: `DMPlexGetPointGlobal()`, `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexPointGlobalRead()`, `DMPlexGetPointLocal()`, `DMPlexPointGlobalRead()`, `DMPlexPointGlobalRef()`
332 @*/
333 PetscErrorCode DMPlexGetPointGlobalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
334 {
335   PetscInt s, e;
336 
337   PetscFunctionBegin;
338   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
339   if (start) PetscValidIntPointer(start, 4);
340   if (end) PetscValidIntPointer(end, 5);
341   PetscCall(DMGetGlobalFieldOffset_Private(dm, point, field, &s, &e));
342   if (start) *start = s;
343   if (end) *end = e;
344   PetscFunctionReturn(0);
345 }
346 
347 /*@
348    DMPlexPointGlobalFieldRead - return read access to a field on a point in global array
349 
350    Not Collective
351 
352    Input Parameters:
353 +  dm - DM defining topological space
354 .  point - topological point
355 .  field - field number
356 -  array - array to index into
357 
358    Output Parameter:
359 .  ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
360 
361    Level: intermediate
362 
363 .seealso: `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexGetPointGlobal()`, `DMPlexPointLocalRead()`, `DMPlexPointGlobalRef()`
364 @*/
365 PetscErrorCode DMPlexPointGlobalFieldRead(DM dm, PetscInt point, PetscInt field, const PetscScalar *array, void *ptr)
366 {
367   PetscInt start, end;
368 
369   PetscFunctionBegin;
370   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
371   PetscValidScalarPointer(array, 4);
372   PetscValidPointer(ptr, 5);
373   PetscCall(DMGetGlobalFieldOffset_Private(dm, point, field, &start, &end));
374   *(const PetscScalar **)ptr = (start < end) ? array + start - dm->map->rstart : NULL;
375   PetscFunctionReturn(0);
376 }
377 
378 /*@
379    DMPlexPointGlobalFieldRef - return read/write access to a field on a point in global array
380 
381    Not Collective
382 
383    Input Parameters:
384 +  dm - DM defining topological space
385 .  point - topological point
386 .  field - field number
387 -  array - array to index into
388 
389    Output Parameter:
390 .  ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
391 
392    Level: intermediate
393 
394 .seealso: `DMGetLocalSection()`, `PetscSectionGetOffset()`, `PetscSectionGetDof()`, `DMPlexGetPointGlobal()`, `DMPlexPointLocalRef()`, `DMPlexPointGlobalRead()`
395 @*/
396 PetscErrorCode DMPlexPointGlobalFieldRef(DM dm, PetscInt point, PetscInt field, PetscScalar *array, void *ptr)
397 {
398   PetscInt start, end;
399 
400   PetscFunctionBegin;
401   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
402   PetscValidScalarPointer(array, 4);
403   PetscValidPointer(ptr, 5);
404   PetscCall(DMGetGlobalFieldOffset_Private(dm, point, field, &start, &end));
405   *(PetscScalar **)ptr = (start < end) ? array + start - dm->map->rstart : NULL;
406   PetscFunctionReturn(0);
407 }
408