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