xref: /petsc/src/dm/impls/plex/plexpoint.c (revision 0298fd7132830bec7daee99a80be0eddb2b310a5)
1552f7358SJed Brown #include <petsc-private/pleximpl.h>   /*I      "petscdmplex.h"   I*/
2552f7358SJed Brown 
3552f7358SJed Brown #undef __FUNCT__
4552f7358SJed Brown #define __FUNCT__ "DMPlexGetLocalOffset_Private"
5552f7358SJed Brown PETSC_STATIC_INLINE PetscErrorCode DMPlexGetLocalOffset_Private(DM dm,PetscInt point,PetscInt *offset)
6552f7358SJed Brown {
7552f7358SJed Brown   PetscFunctionBegin;
8552f7358SJed Brown #if defined(PETSC_USE_DEBUG)
9552f7358SJed Brown   {
10552f7358SJed Brown     PetscErrorCode ierr;
11552f7358SJed Brown     ierr = PetscSectionGetOffset(dm->defaultSection,point,offset);CHKERRQ(ierr);
12552f7358SJed Brown   }
13552f7358SJed Brown #else
14552f7358SJed Brown   {
15552f7358SJed Brown     PetscSection s = dm->defaultSection;
16552f7358SJed Brown     *offset = s->atlasOff[point - s->atlasLayout.pStart];
17552f7358SJed Brown   }
18552f7358SJed Brown #endif
19552f7358SJed Brown   PetscFunctionReturn(0);
20552f7358SJed Brown }
21552f7358SJed Brown 
22552f7358SJed Brown #undef __FUNCT__
23552f7358SJed Brown #define __FUNCT__ "DMPlexGetGlobalOffset_Private"
24552f7358SJed Brown PETSC_STATIC_INLINE PetscErrorCode DMPlexGetGlobalOffset_Private(DM dm,PetscInt point,PetscInt *offset)
25552f7358SJed Brown {
26552f7358SJed Brown   PetscFunctionBegin;
27552f7358SJed Brown #if defined(PETSC_USE_DEBUG)
28552f7358SJed Brown   {
29552f7358SJed Brown     PetscErrorCode ierr;
30552f7358SJed Brown     PetscInt       dof,cdof;
31552f7358SJed Brown     ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,offset);CHKERRQ(ierr);
32552f7358SJed Brown     ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr);
33552f7358SJed Brown     ierr = PetscSectionGetConstraintDof(dm->defaultGlobalSection,point,&cdof);CHKERRQ(ierr);
34552f7358SJed Brown     if (dof-cdof <= 0) *offset = -1; /* Indicates no data */
35552f7358SJed Brown   }
36552f7358SJed Brown #else
37552f7358SJed Brown   {
38552f7358SJed Brown     PetscSection s = dm->defaultGlobalSection;
39552f7358SJed Brown     PetscInt     dof,cdof;
40552f7358SJed Brown     *offset = s->atlasOff[point - s->atlasLayout.pStart];
41552f7358SJed Brown     dof     = s->atlasDof[point - s->atlasLayout.pStart];
42552f7358SJed Brown     cdof    = s->bc ? s->bc->atlasDof[point - s->bc->atlasLayout.pStart] : 0;
43552f7358SJed Brown     if (dof-cdof <= 0) *offset = -1;
44552f7358SJed Brown   }
45552f7358SJed Brown #endif
46552f7358SJed Brown   PetscFunctionReturn(0);
47552f7358SJed Brown }
48552f7358SJed Brown 
49552f7358SJed Brown #undef __FUNCT__
50552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointLocal"
51552f7358SJed Brown /*@
52552f7358SJed Brown    DMPlexGetPointLocal - get location of point data in local Vec
53552f7358SJed Brown 
54552f7358SJed Brown    Not Collective
55552f7358SJed Brown 
56552f7358SJed Brown    Input Arguments:
57552f7358SJed Brown +  dm - DM defining the topological space
58552f7358SJed Brown -  point - topological point
59552f7358SJed Brown 
60552f7358SJed Brown    Output Arguments:
61552f7358SJed Brown +  start - start of point data
62552f7358SJed Brown -  end - end of point data
63552f7358SJed Brown 
64552f7358SJed Brown    Level: intermediate
65552f7358SJed Brown 
66552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef()
67552f7358SJed Brown @*/
68552f7358SJed Brown PetscErrorCode DMPlexGetPointLocal(DM dm,PetscInt point,PetscInt *start,PetscInt *end)
69552f7358SJed Brown {
70552f7358SJed Brown   PetscErrorCode ierr;
71552f7358SJed Brown   PetscInt       offset,dof;
72552f7358SJed Brown 
73552f7358SJed Brown   PetscFunctionBegin;
74552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
75552f7358SJed Brown   ierr = PetscSectionGetOffset(dm->defaultSection,point,&offset);CHKERRQ(ierr);
76552f7358SJed Brown   ierr = PetscSectionGetDof(dm->defaultSection,point,&dof);CHKERRQ(ierr);
77552f7358SJed Brown   if (start) *start = offset;
78552f7358SJed Brown   if (end) *end = offset + dof;
79552f7358SJed Brown   PetscFunctionReturn(0);
80552f7358SJed Brown }
81552f7358SJed Brown 
82552f7358SJed Brown #undef __FUNCT__
83552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRead"
84552f7358SJed Brown /*@
85552f7358SJed Brown    DMPlexPointLocalRead - return read access to a point in local array
86552f7358SJed Brown 
87552f7358SJed Brown    Not Collective
88552f7358SJed Brown 
89552f7358SJed Brown    Input Arguments:
90552f7358SJed Brown +  dm - DM defining topological space
91552f7358SJed Brown .  point - topological point
92552f7358SJed Brown -  array - array to index into
93552f7358SJed Brown 
94552f7358SJed Brown    Output Arguments:
95552f7358SJed Brown .  ptr - address of read reference to point data, type generic so user can place in structure
96552f7358SJed Brown 
97552f7358SJed Brown    Level: intermediate
98552f7358SJed Brown 
99552f7358SJed Brown    Note:
100552f7358SJed Brown    A common usage when data sizes are known statically:
101552f7358SJed Brown 
102552f7358SJed Brown $  const struct { PetscScalar foo,bar,baz; } *ptr;
103552f7358SJed Brown $  DMPlexPointLocalRead(dm,point,array,&ptr);
104552f7358SJed Brown $  x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
105552f7358SJed Brown 
106552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead()
107552f7358SJed Brown @*/
108552f7358SJed Brown PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
109552f7358SJed Brown {
110552f7358SJed Brown   PetscErrorCode ierr;
111552f7358SJed Brown   PetscInt       start;
112552f7358SJed Brown 
113552f7358SJed Brown   PetscFunctionBegin;
114552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
115552f7358SJed Brown   PetscValidScalarPointer(array,3);
116552f7358SJed Brown   PetscValidPointer(ptr,4);
117552f7358SJed Brown   ierr                      = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr);
118552f7358SJed Brown   *(const PetscScalar**)ptr = array + start;
119552f7358SJed Brown   PetscFunctionReturn(0);
120552f7358SJed Brown }
121552f7358SJed Brown 
122552f7358SJed Brown #undef __FUNCT__
123552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRef"
124552f7358SJed Brown /*@
125552f7358SJed Brown    DMPlexPointLocalRef - return read/write access to a point in local array
126552f7358SJed Brown 
127552f7358SJed Brown    Not Collective
128552f7358SJed Brown 
129552f7358SJed Brown    Input Arguments:
130552f7358SJed Brown +  dm - DM defining topological space
131552f7358SJed Brown .  point - topological point
132552f7358SJed Brown -  array - array to index into
133552f7358SJed Brown 
134552f7358SJed Brown    Output Arguments:
135552f7358SJed Brown .  ptr - address of reference to point data, type generic so user can place in structure
136552f7358SJed Brown 
137552f7358SJed Brown    Level: intermediate
138552f7358SJed Brown 
139552f7358SJed Brown    Note:
140552f7358SJed Brown    A common usage when data sizes are known statically:
141552f7358SJed Brown 
142552f7358SJed Brown $  struct { PetscScalar foo,bar,baz; } *ptr;
143552f7358SJed Brown $  DMPlexPointLocalRef(dm,point,array,&ptr);
144552f7358SJed Brown $  ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
145552f7358SJed Brown 
146552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
147552f7358SJed Brown @*/
148552f7358SJed Brown PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
149552f7358SJed Brown {
150552f7358SJed Brown   PetscErrorCode ierr;
151552f7358SJed Brown   PetscInt       start;
152552f7358SJed Brown 
153552f7358SJed Brown   PetscFunctionBegin;
154552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
155552f7358SJed Brown   PetscValidScalarPointer(array,3);
156552f7358SJed Brown   PetscValidPointer(ptr,4);
157552f7358SJed Brown   ierr                = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr);
158552f7358SJed Brown   *(PetscScalar**)ptr = array + start;
159552f7358SJed Brown   PetscFunctionReturn(0);
160552f7358SJed Brown }
161552f7358SJed Brown 
162552f7358SJed Brown #undef __FUNCT__
163552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointGlobal"
164552f7358SJed Brown /*@
165552f7358SJed Brown    DMPlexGetPointGlobal - get location of point data in global Vec
166552f7358SJed Brown 
167552f7358SJed Brown    Not Collective
168552f7358SJed Brown 
169552f7358SJed Brown    Input Arguments:
170552f7358SJed Brown +  dm - DM defining the topological space
171552f7358SJed Brown -  point - topological point
172552f7358SJed Brown 
173552f7358SJed Brown    Output Arguments:
174552f7358SJed Brown +  start - start of point data; returns -(global_start+1) if point is not owned
175552f7358SJed Brown -  end - end of point data; returns -(global_end+1) if point is not owned
176552f7358SJed Brown 
177552f7358SJed Brown    Level: intermediate
178552f7358SJed Brown 
179552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef()
180552f7358SJed Brown @*/
181552f7358SJed Brown PetscErrorCode DMPlexGetPointGlobal(DM dm,PetscInt point,PetscInt *start,PetscInt *end)
182552f7358SJed Brown {
183552f7358SJed Brown   PetscErrorCode ierr;
184552f7358SJed Brown   PetscInt       offset,dof;
185552f7358SJed Brown 
186552f7358SJed Brown   PetscFunctionBegin;
187552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
188552f7358SJed Brown   ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,&offset);CHKERRQ(ierr);
189552f7358SJed Brown   ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr);
190552f7358SJed Brown   if (start) *start = offset;
191552f7358SJed Brown   if (end) *end = offset + dof;
192552f7358SJed Brown   PetscFunctionReturn(0);
193552f7358SJed Brown }
194552f7358SJed Brown 
195552f7358SJed Brown #undef __FUNCT__
196552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRead"
197552f7358SJed Brown /*@
198552f7358SJed Brown    DMPlexPointGlobalRead - return read access to a point in global array
199552f7358SJed Brown 
200552f7358SJed Brown    Not Collective
201552f7358SJed Brown 
202552f7358SJed Brown    Input Arguments:
203552f7358SJed Brown +  dm - DM defining topological space
204552f7358SJed Brown .  point - topological point
205552f7358SJed Brown -  array - array to index into
206552f7358SJed Brown 
207552f7358SJed Brown    Output Arguments:
208*0298fd71SBarry Smith .  ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
209552f7358SJed Brown 
210552f7358SJed Brown    Level: intermediate
211552f7358SJed Brown 
212552f7358SJed Brown    Note:
213552f7358SJed Brown    A common usage when data sizes are known statically:
214552f7358SJed Brown 
215552f7358SJed Brown $  const struct { PetscScalar foo,bar,baz; } *ptr;
216552f7358SJed Brown $  DMPlexPointGlobalRead(dm,point,array,&ptr);
217552f7358SJed Brown $  x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
218552f7358SJed Brown 
219552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef()
220552f7358SJed Brown @*/
221552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
222552f7358SJed Brown {
223552f7358SJed Brown   PetscErrorCode ierr;
224552f7358SJed Brown   PetscInt       start;
225552f7358SJed Brown 
226552f7358SJed Brown   PetscFunctionBegin;
227552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
228552f7358SJed Brown   PetscValidScalarPointer(array,3);
229552f7358SJed Brown   PetscValidPointer(ptr,4);
230552f7358SJed Brown   ierr                      = DMPlexGetGlobalOffset_Private(dm,point,&start);CHKERRQ(ierr);
231*0298fd71SBarry Smith   *(const PetscScalar**)ptr = (start >= 0) ? array + start - dm->map->rstart : NULL;
232552f7358SJed Brown   PetscFunctionReturn(0);
233552f7358SJed Brown }
234552f7358SJed Brown 
235552f7358SJed Brown #undef __FUNCT__
236552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRef"
237552f7358SJed Brown /*@
238552f7358SJed Brown    DMPlexPointGlobalRef - return read/write access to a point in global array
239552f7358SJed Brown 
240552f7358SJed Brown    Not Collective
241552f7358SJed Brown 
242552f7358SJed Brown    Input Arguments:
243552f7358SJed Brown +  dm - DM defining topological space
244552f7358SJed Brown .  point - topological point
245552f7358SJed Brown -  array - array to index into
246552f7358SJed Brown 
247552f7358SJed Brown    Output Arguments:
248*0298fd71SBarry Smith .  ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
249552f7358SJed Brown 
250552f7358SJed Brown    Level: intermediate
251552f7358SJed Brown 
252552f7358SJed Brown    Note:
253552f7358SJed Brown    A common usage when data sizes are known statically:
254552f7358SJed Brown 
255552f7358SJed Brown $  struct { PetscScalar foo,bar,baz; } *ptr;
256552f7358SJed Brown $  DMPlexPointGlobalRef(dm,point,array,&ptr);
257552f7358SJed Brown $  ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
258552f7358SJed Brown 
259552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead()
260552f7358SJed Brown @*/
261552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
262552f7358SJed Brown {
263552f7358SJed Brown   PetscErrorCode ierr;
264552f7358SJed Brown   PetscInt       start;
265552f7358SJed Brown 
266552f7358SJed Brown   PetscFunctionBegin;
267552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
268552f7358SJed Brown   PetscValidScalarPointer(array,3);
269552f7358SJed Brown   PetscValidPointer(ptr,4);
270552f7358SJed Brown   ierr                = DMPlexGetGlobalOffset_Private(dm,point,&start);CHKERRQ(ierr);
271*0298fd71SBarry Smith   *(PetscScalar**)ptr = (start >= 0) ? array + start - dm->map->rstart : NULL;
272552f7358SJed Brown   PetscFunctionReturn(0);
273552f7358SJed Brown }
274