xref: /petsc/src/dm/impls/plex/plexpoint.c (revision 552f735842aa6127d09b62f740a903cdd0631614)
1*552f7358SJed Brown #include <petsc-private/pleximpl.h>   /*I      "petscdmplex.h"   I*/
2*552f7358SJed Brown 
3*552f7358SJed Brown #undef __FUNCT__
4*552f7358SJed Brown #define __FUNCT__ "DMPlexGetLocalOffset_Private"
5*552f7358SJed Brown PETSC_STATIC_INLINE PetscErrorCode DMPlexGetLocalOffset_Private(DM dm,PetscInt point,PetscInt *offset)
6*552f7358SJed Brown {
7*552f7358SJed Brown   PetscFunctionBegin;
8*552f7358SJed Brown #if defined(PETSC_USE_DEBUG)
9*552f7358SJed Brown   {
10*552f7358SJed Brown     PetscErrorCode ierr;
11*552f7358SJed Brown     ierr = PetscSectionGetOffset(dm->defaultSection,point,offset);CHKERRQ(ierr);
12*552f7358SJed Brown   }
13*552f7358SJed Brown #else
14*552f7358SJed Brown   {
15*552f7358SJed Brown     PetscSection s = dm->defaultSection;
16*552f7358SJed Brown     *offset = s->atlasOff[point - s->atlasLayout.pStart];
17*552f7358SJed Brown   }
18*552f7358SJed Brown #endif
19*552f7358SJed Brown   PetscFunctionReturn(0);
20*552f7358SJed Brown }
21*552f7358SJed Brown 
22*552f7358SJed Brown #undef __FUNCT__
23*552f7358SJed Brown #define __FUNCT__ "DMPlexGetGlobalOffset_Private"
24*552f7358SJed Brown PETSC_STATIC_INLINE PetscErrorCode DMPlexGetGlobalOffset_Private(DM dm,PetscInt point,PetscInt *offset)
25*552f7358SJed Brown {
26*552f7358SJed Brown   PetscFunctionBegin;
27*552f7358SJed Brown #if defined(PETSC_USE_DEBUG)
28*552f7358SJed Brown   {
29*552f7358SJed Brown     PetscErrorCode ierr;
30*552f7358SJed Brown     PetscInt dof,cdof;
31*552f7358SJed Brown     ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,offset);CHKERRQ(ierr);
32*552f7358SJed Brown     ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr);
33*552f7358SJed Brown     ierr = PetscSectionGetConstraintDof(dm->defaultGlobalSection,point,&cdof);CHKERRQ(ierr);
34*552f7358SJed Brown     if (dof-cdof <= 0) *offset = -1; /* Indicates no data */
35*552f7358SJed Brown   }
36*552f7358SJed Brown #else
37*552f7358SJed Brown   {
38*552f7358SJed Brown     PetscSection s = dm->defaultGlobalSection;
39*552f7358SJed Brown     PetscInt dof,cdof;
40*552f7358SJed Brown     *offset = s->atlasOff[point - s->atlasLayout.pStart];
41*552f7358SJed Brown     dof = s->atlasDof[point - s->atlasLayout.pStart];
42*552f7358SJed Brown     cdof = s->bc ? s->bc->atlasDof[point - s->bc->atlasLayout.pStart] : 0;
43*552f7358SJed Brown     if (dof-cdof <= 0) *offset = -1;
44*552f7358SJed Brown   }
45*552f7358SJed Brown #endif
46*552f7358SJed Brown   PetscFunctionReturn(0);
47*552f7358SJed Brown }
48*552f7358SJed Brown 
49*552f7358SJed Brown #undef __FUNCT__
50*552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointLocal"
51*552f7358SJed Brown /*@
52*552f7358SJed Brown    DMPlexGetPointLocal - get location of point data in local Vec
53*552f7358SJed Brown 
54*552f7358SJed Brown    Not Collective
55*552f7358SJed Brown 
56*552f7358SJed Brown    Input Arguments:
57*552f7358SJed Brown +  dm - DM defining the topological space
58*552f7358SJed Brown -  point - topological point
59*552f7358SJed Brown 
60*552f7358SJed Brown    Output Arguments:
61*552f7358SJed Brown +  start - start of point data
62*552f7358SJed Brown -  end - end of point data
63*552f7358SJed Brown 
64*552f7358SJed Brown    Level: intermediate
65*552f7358SJed Brown 
66*552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef()
67*552f7358SJed Brown @*/
68*552f7358SJed Brown PetscErrorCode DMPlexGetPointLocal(DM dm,PetscInt point,PetscInt *start,PetscInt *end)
69*552f7358SJed Brown {
70*552f7358SJed Brown   PetscErrorCode ierr;
71*552f7358SJed Brown   PetscInt offset,dof;
72*552f7358SJed Brown 
73*552f7358SJed Brown   PetscFunctionBegin;
74*552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
75*552f7358SJed Brown   ierr = PetscSectionGetOffset(dm->defaultSection,point,&offset);CHKERRQ(ierr);
76*552f7358SJed Brown   ierr = PetscSectionGetDof(dm->defaultSection,point,&dof);CHKERRQ(ierr);
77*552f7358SJed Brown   if (start) *start = offset;
78*552f7358SJed Brown   if (end) *end = offset + dof;
79*552f7358SJed Brown   PetscFunctionReturn(0);
80*552f7358SJed Brown }
81*552f7358SJed Brown 
82*552f7358SJed Brown #undef __FUNCT__
83*552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRead"
84*552f7358SJed Brown /*@
85*552f7358SJed Brown    DMPlexPointLocalRead - return read access to a point in local array
86*552f7358SJed Brown 
87*552f7358SJed Brown    Not Collective
88*552f7358SJed Brown 
89*552f7358SJed Brown    Input Arguments:
90*552f7358SJed Brown +  dm - DM defining topological space
91*552f7358SJed Brown .  point - topological point
92*552f7358SJed Brown -  array - array to index into
93*552f7358SJed Brown 
94*552f7358SJed Brown    Output Arguments:
95*552f7358SJed Brown .  ptr - address of read reference to point data, type generic so user can place in structure
96*552f7358SJed Brown 
97*552f7358SJed Brown    Level: intermediate
98*552f7358SJed Brown 
99*552f7358SJed Brown    Note:
100*552f7358SJed Brown    A common usage when data sizes are known statically:
101*552f7358SJed Brown 
102*552f7358SJed Brown $  const struct { PetscScalar foo,bar,baz; } *ptr;
103*552f7358SJed Brown $  DMPlexPointLocalRead(dm,point,array,&ptr);
104*552f7358SJed Brown $  x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
105*552f7358SJed Brown 
106*552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead()
107*552f7358SJed Brown @*/
108*552f7358SJed Brown PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
109*552f7358SJed Brown {
110*552f7358SJed Brown   PetscErrorCode ierr;
111*552f7358SJed Brown   PetscInt start;
112*552f7358SJed Brown 
113*552f7358SJed Brown   PetscFunctionBegin;
114*552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
115*552f7358SJed Brown   PetscValidScalarPointer(array,3);
116*552f7358SJed Brown   PetscValidPointer(ptr,4);
117*552f7358SJed Brown   ierr = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr);
118*552f7358SJed Brown   *(const PetscScalar **)ptr = array + start;
119*552f7358SJed Brown   PetscFunctionReturn(0);
120*552f7358SJed Brown }
121*552f7358SJed Brown 
122*552f7358SJed Brown #undef __FUNCT__
123*552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRef"
124*552f7358SJed Brown /*@
125*552f7358SJed Brown    DMPlexPointLocalRef - return read/write access to a point in local array
126*552f7358SJed Brown 
127*552f7358SJed Brown    Not Collective
128*552f7358SJed Brown 
129*552f7358SJed Brown    Input Arguments:
130*552f7358SJed Brown +  dm - DM defining topological space
131*552f7358SJed Brown .  point - topological point
132*552f7358SJed Brown -  array - array to index into
133*552f7358SJed Brown 
134*552f7358SJed Brown    Output Arguments:
135*552f7358SJed Brown .  ptr - address of reference to point data, type generic so user can place in structure
136*552f7358SJed Brown 
137*552f7358SJed Brown    Level: intermediate
138*552f7358SJed Brown 
139*552f7358SJed Brown    Note:
140*552f7358SJed Brown    A common usage when data sizes are known statically:
141*552f7358SJed Brown 
142*552f7358SJed Brown $  struct { PetscScalar foo,bar,baz; } *ptr;
143*552f7358SJed Brown $  DMPlexPointLocalRef(dm,point,array,&ptr);
144*552f7358SJed Brown $  ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
145*552f7358SJed Brown 
146*552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
147*552f7358SJed Brown @*/
148*552f7358SJed Brown PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
149*552f7358SJed Brown {
150*552f7358SJed Brown   PetscErrorCode ierr;
151*552f7358SJed Brown   PetscInt start;
152*552f7358SJed Brown 
153*552f7358SJed Brown   PetscFunctionBegin;
154*552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
155*552f7358SJed Brown   PetscValidScalarPointer(array,3);
156*552f7358SJed Brown   PetscValidPointer(ptr,4);
157*552f7358SJed Brown   ierr = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr);
158*552f7358SJed Brown   *(PetscScalar **)ptr = array + start;
159*552f7358SJed Brown   PetscFunctionReturn(0);
160*552f7358SJed Brown }
161*552f7358SJed Brown 
162*552f7358SJed Brown #undef __FUNCT__
163*552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointGlobal"
164*552f7358SJed Brown /*@
165*552f7358SJed Brown    DMPlexGetPointGlobal - get location of point data in global Vec
166*552f7358SJed Brown 
167*552f7358SJed Brown    Not Collective
168*552f7358SJed Brown 
169*552f7358SJed Brown    Input Arguments:
170*552f7358SJed Brown +  dm - DM defining the topological space
171*552f7358SJed Brown -  point - topological point
172*552f7358SJed Brown 
173*552f7358SJed Brown    Output Arguments:
174*552f7358SJed Brown +  start - start of point data; returns -(global_start+1) if point is not owned
175*552f7358SJed Brown -  end - end of point data; returns -(global_end+1) if point is not owned
176*552f7358SJed Brown 
177*552f7358SJed Brown    Level: intermediate
178*552f7358SJed Brown 
179*552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef()
180*552f7358SJed Brown @*/
181*552f7358SJed Brown PetscErrorCode DMPlexGetPointGlobal(DM dm,PetscInt point,PetscInt *start,PetscInt *end)
182*552f7358SJed Brown {
183*552f7358SJed Brown   PetscErrorCode ierr;
184*552f7358SJed Brown   PetscInt offset,dof;
185*552f7358SJed Brown 
186*552f7358SJed Brown   PetscFunctionBegin;
187*552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
188*552f7358SJed Brown   ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,&offset);CHKERRQ(ierr);
189*552f7358SJed Brown   ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr);
190*552f7358SJed Brown   if (start) *start = offset;
191*552f7358SJed Brown   if (end) *end = offset + dof;
192*552f7358SJed Brown   PetscFunctionReturn(0);
193*552f7358SJed Brown }
194*552f7358SJed Brown 
195*552f7358SJed Brown #undef __FUNCT__
196*552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRead"
197*552f7358SJed Brown /*@
198*552f7358SJed Brown    DMPlexPointGlobalRead - return read access to a point in global array
199*552f7358SJed Brown 
200*552f7358SJed Brown    Not Collective
201*552f7358SJed Brown 
202*552f7358SJed Brown    Input Arguments:
203*552f7358SJed Brown +  dm - DM defining topological space
204*552f7358SJed Brown .  point - topological point
205*552f7358SJed Brown -  array - array to index into
206*552f7358SJed Brown 
207*552f7358SJed Brown    Output Arguments:
208*552f7358SJed Brown .  ptr - address of read reference to point data, type generic so user can place in structure; returns PETSC_NULL if global point is not owned
209*552f7358SJed Brown 
210*552f7358SJed Brown    Level: intermediate
211*552f7358SJed Brown 
212*552f7358SJed Brown    Note:
213*552f7358SJed Brown    A common usage when data sizes are known statically:
214*552f7358SJed Brown 
215*552f7358SJed Brown $  const struct { PetscScalar foo,bar,baz; } *ptr;
216*552f7358SJed Brown $  DMPlexPointGlobalRead(dm,point,array,&ptr);
217*552f7358SJed Brown $  x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
218*552f7358SJed Brown 
219*552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef()
220*552f7358SJed Brown @*/
221*552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
222*552f7358SJed Brown {
223*552f7358SJed Brown   PetscErrorCode ierr;
224*552f7358SJed Brown   PetscInt start;
225*552f7358SJed Brown 
226*552f7358SJed Brown   PetscFunctionBegin;
227*552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
228*552f7358SJed Brown   PetscValidScalarPointer(array,3);
229*552f7358SJed Brown   PetscValidPointer(ptr,4);
230*552f7358SJed Brown   ierr = DMPlexGetGlobalOffset_Private(dm,point,&start);CHKERRQ(ierr);
231*552f7358SJed Brown   *(const PetscScalar **)ptr = (start >= 0) ? array + start - dm->map->rstart : PETSC_NULL;
232*552f7358SJed Brown   PetscFunctionReturn(0);
233*552f7358SJed Brown }
234*552f7358SJed Brown 
235*552f7358SJed Brown #undef __FUNCT__
236*552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRef"
237*552f7358SJed Brown /*@
238*552f7358SJed Brown    DMPlexPointGlobalRef - return read/write access to a point in global array
239*552f7358SJed Brown 
240*552f7358SJed Brown    Not Collective
241*552f7358SJed Brown 
242*552f7358SJed Brown    Input Arguments:
243*552f7358SJed Brown +  dm - DM defining topological space
244*552f7358SJed Brown .  point - topological point
245*552f7358SJed Brown -  array - array to index into
246*552f7358SJed Brown 
247*552f7358SJed Brown    Output Arguments:
248*552f7358SJed Brown .  ptr - address of reference to point data, type generic so user can place in structure; returns PETSC_NULL if global point is not owned
249*552f7358SJed Brown 
250*552f7358SJed Brown    Level: intermediate
251*552f7358SJed Brown 
252*552f7358SJed Brown    Note:
253*552f7358SJed Brown    A common usage when data sizes are known statically:
254*552f7358SJed Brown 
255*552f7358SJed Brown $  struct { PetscScalar foo,bar,baz; } *ptr;
256*552f7358SJed Brown $  DMPlexPointGlobalRef(dm,point,array,&ptr);
257*552f7358SJed Brown $  ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
258*552f7358SJed Brown 
259*552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead()
260*552f7358SJed Brown @*/
261*552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
262*552f7358SJed Brown {
263*552f7358SJed Brown   PetscErrorCode ierr;
264*552f7358SJed Brown   PetscInt start;
265*552f7358SJed Brown 
266*552f7358SJed Brown   PetscFunctionBegin;
267*552f7358SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
268*552f7358SJed Brown   PetscValidScalarPointer(array,3);
269*552f7358SJed Brown   PetscValidPointer(ptr,4);
270*552f7358SJed Brown   ierr = DMPlexGetGlobalOffset_Private(dm,point,&start);CHKERRQ(ierr);
271*552f7358SJed Brown   *(PetscScalar **)ptr = (start >= 0) ? array + start - dm->map->rstart : PETSC_NULL;
272*552f7358SJed Brown   PetscFunctionReturn(0);
273*552f7358SJed Brown }
274