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