Actual source code: plexpoint.c
petsc-3.4.5 2014-06-29
1: #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2: #include <petsc-private/isimpl.h> /* for inline access to atlasOff */
6: PETSC_STATIC_INLINE PetscErrorCode DMPlexGetLocalOffset_Private(DM dm,PetscInt point,PetscInt *offset)
7: {
9: #if defined(PETSC_USE_DEBUG)
10: {
12: PetscSectionGetOffset(dm->defaultSection,point,offset);
13: }
14: #else
15: {
16: PetscSection s = dm->defaultSection;
17: *offset = s->atlasOff[point - s->atlasLayout.pStart];
18: }
19: #endif
20: return(0);
21: }
25: PETSC_STATIC_INLINE PetscErrorCode DMPlexGetGlobalOffset_Private(DM dm,PetscInt point,PetscInt *offset)
26: {
28: #if defined(PETSC_USE_DEBUG)
29: {
31: PetscInt dof,cdof;
32: PetscSectionGetOffset(dm->defaultGlobalSection,point,offset);
33: PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);
34: PetscSectionGetConstraintDof(dm->defaultGlobalSection,point,&cdof);
35: if (dof-cdof <= 0) *offset = -1; /* Indicates no data */
36: }
37: #else
38: {
39: PetscSection s = dm->defaultGlobalSection;
40: PetscInt dof,cdof;
41: *offset = s->atlasOff[point - s->atlasLayout.pStart];
42: dof = s->atlasDof[point - s->atlasLayout.pStart];
43: cdof = s->bc ? s->bc->atlasDof[point - s->bc->atlasLayout.pStart] : 0;
44: if (dof-cdof <= 0) *offset = -1;
45: }
46: #endif
47: return(0);
48: }
52: /*@
53: DMPlexGetPointLocal - get location of point data in local Vec
55: Not Collective
57: Input Arguments:
58: + dm - DM defining the topological space
59: - point - topological point
61: Output Arguments:
62: + start - start of point data
63: - end - end of point data
65: Level: intermediate
67: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef()
68: @*/
69: PetscErrorCode DMPlexGetPointLocal(DM dm,PetscInt point,PetscInt *start,PetscInt *end)
70: {
72: PetscInt offset,dof;
76: PetscSectionGetOffset(dm->defaultSection,point,&offset);
77: PetscSectionGetDof(dm->defaultSection,point,&dof);
78: if (start) *start = offset;
79: if (end) *end = offset + dof;
80: return(0);
81: }
85: /*@
86: DMPlexPointLocalRead - return read access to a point in local array
88: Not Collective
90: Input Arguments:
91: + dm - DM defining topological space
92: . point - topological point
93: - array - array to index into
95: Output Arguments:
96: . ptr - address of read reference to point data, type generic so user can place in structure
98: Level: intermediate
100: Note:
101: A common usage when data sizes are known statically:
103: $ const struct { PetscScalar foo,bar,baz; } *ptr;
104: $ DMPlexPointLocalRead(dm,point,array,&ptr);
105: $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
107: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead()
108: @*/
109: PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
110: {
112: PetscInt start;
118: DMPlexGetLocalOffset_Private(dm,point,&start);
119: *(const PetscScalar**)ptr = array + start;
120: return(0);
121: }
125: /*@
126: DMPlexPointLocalRef - return read/write access to a point in local array
128: Not Collective
130: Input Arguments:
131: + dm - DM defining topological space
132: . point - topological point
133: - array - array to index into
135: Output Arguments:
136: . ptr - address of reference to point data, type generic so user can place in structure
138: Level: intermediate
140: Note:
141: A common usage when data sizes are known statically:
143: $ struct { PetscScalar foo,bar,baz; } *ptr;
144: $ DMPlexPointLocalRef(dm,point,array,&ptr);
145: $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
147: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef()
148: @*/
149: PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
150: {
152: PetscInt start;
158: DMPlexGetLocalOffset_Private(dm,point,&start);
159: *(PetscScalar**)ptr = array + start;
160: return(0);
161: }
165: /*@
166: DMPlexGetPointGlobal - get location of point data in global Vec
168: Not Collective
170: Input Arguments:
171: + dm - DM defining the topological space
172: - point - topological point
174: Output Arguments:
175: + start - start of point data; returns -(global_start+1) if point is not owned
176: - end - end of point data; returns -(global_end+1) if point is not owned
178: Level: intermediate
180: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef()
181: @*/
182: PetscErrorCode DMPlexGetPointGlobal(DM dm,PetscInt point,PetscInt *start,PetscInt *end)
183: {
185: PetscInt offset,dof;
189: PetscSectionGetOffset(dm->defaultGlobalSection,point,&offset);
190: PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);
191: if (start) *start = offset;
192: if (end) *end = offset + dof;
193: return(0);
194: }
198: /*@
199: DMPlexPointGlobalRead - return read access to a point in global array
201: Not Collective
203: Input Arguments:
204: + dm - DM defining topological space
205: . point - topological point
206: - array - array to index into
208: Output Arguments:
209: . ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
211: Level: intermediate
213: Note:
214: A common usage when data sizes are known statically:
216: $ const struct { PetscScalar foo,bar,baz; } *ptr;
217: $ DMPlexPointGlobalRead(dm,point,array,&ptr);
218: $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz;
220: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef()
221: @*/
222: PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr)
223: {
225: PetscInt start;
231: DMPlexGetGlobalOffset_Private(dm,point,&start);
232: *(const PetscScalar**)ptr = (start >= 0) ? array + start - dm->map->rstart : NULL;
233: return(0);
234: }
238: /*@
239: DMPlexPointGlobalRef - return read/write access to a point in global array
241: Not Collective
243: Input Arguments:
244: + dm - DM defining topological space
245: . point - topological point
246: - array - array to index into
248: Output Arguments:
249: . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned
251: Level: intermediate
253: Note:
254: A common usage when data sizes are known statically:
256: $ struct { PetscScalar foo,bar,baz; } *ptr;
257: $ DMPlexPointGlobalRef(dm,point,array,&ptr);
258: $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5;
260: .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead()
261: @*/
262: PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr)
263: {
265: PetscInt start;
271: DMPlexGetGlobalOffset_Private(dm,point,&start);
272: *(PetscScalar**)ptr = (start >= 0) ? array + start - dm->map->rstart : NULL;
273: return(0);
274: }