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