1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 233879625SMatthew G. Knepley 3552f7358SJed Brown /*@ 4552f7358SJed Brown DMPlexGetPointLocal - get location of point data in local Vec 5552f7358SJed Brown 6552f7358SJed Brown Not Collective 7552f7358SJed Brown 8552f7358SJed Brown Input Arguments: 9552f7358SJed Brown + dm - DM defining the topological space 10552f7358SJed Brown - point - topological point 11552f7358SJed Brown 12552f7358SJed Brown Output Arguments: 13552f7358SJed Brown + start - start of point data 14552f7358SJed Brown - end - end of point data 15552f7358SJed Brown 16a89cf0ddSMatthew G. Knepley Note: This is a half open interval [start, end) 17a89cf0ddSMatthew G. Knepley 18552f7358SJed Brown Level: intermediate 19552f7358SJed Brown 20a89cf0ddSMatthew G. Knepley .seealso: DMPlexGetPointLocalField(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef() 21552f7358SJed Brown @*/ 22552f7358SJed Brown PetscErrorCode DMPlexGetPointLocal(DM dm, PetscInt point, PetscInt *start, PetscInt *end) 23552f7358SJed Brown { 24a89cf0ddSMatthew G. Knepley PetscInt s, e; 25552f7358SJed Brown PetscErrorCode ierr; 26552f7358SJed Brown 27552f7358SJed Brown PetscFunctionBegin; 28552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 29a89cf0ddSMatthew G. Knepley if (start) PetscValidPointer(start, 3); 30a89cf0ddSMatthew G. Knepley if (end) PetscValidPointer(end, 4); 31*088fd00dSMatthew G. Knepley ierr = DMGetLocalOffset_Private(dm, point, &s, &e);CHKERRQ(ierr); 32a89cf0ddSMatthew G. Knepley if (start) *start = s; 33a89cf0ddSMatthew G. Knepley if (end) *end = e; 34552f7358SJed Brown PetscFunctionReturn(0); 35552f7358SJed Brown } 36552f7358SJed Brown 37552f7358SJed Brown /*@ 38552f7358SJed Brown DMPlexPointLocalRead - return read access to a point in local array 39552f7358SJed Brown 40552f7358SJed Brown Not Collective 41552f7358SJed Brown 42552f7358SJed Brown Input Arguments: 43552f7358SJed Brown + dm - DM defining topological space 44552f7358SJed Brown . point - topological point 45552f7358SJed Brown - array - array to index into 46552f7358SJed Brown 47552f7358SJed Brown Output Arguments: 48552f7358SJed Brown . ptr - address of read reference to point data, type generic so user can place in structure 49552f7358SJed Brown 50552f7358SJed Brown Level: intermediate 51552f7358SJed Brown 52552f7358SJed Brown Note: 53552f7358SJed Brown A common usage when data sizes are known statically: 54552f7358SJed Brown 55552f7358SJed Brown $ const struct { PetscScalar foo,bar,baz; } *ptr; 56552f7358SJed Brown $ DMPlexPointLocalRead(dm,point,array,&ptr); 57552f7358SJed Brown $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 58552f7358SJed Brown 59552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead() 60552f7358SJed Brown @*/ 61081a2d76SSatish Balay PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,void *ptr) 62552f7358SJed Brown { 63552f7358SJed Brown PetscErrorCode ierr; 64a89cf0ddSMatthew G. Knepley PetscInt start, end; 65552f7358SJed Brown 66552f7358SJed Brown PetscFunctionBegin; 67552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 68552f7358SJed Brown PetscValidScalarPointer(array,3); 69552f7358SJed Brown PetscValidPointer(ptr,4); 70*088fd00dSMatthew G. Knepley ierr = DMGetLocalOffset_Private(dm,point,&start,&end);CHKERRQ(ierr); 71a89cf0ddSMatthew G. Knepley *(const PetscScalar**)ptr = (start < end) ? array + start : NULL; 72552f7358SJed Brown PetscFunctionReturn(0); 73552f7358SJed Brown } 74552f7358SJed Brown 75552f7358SJed Brown /*@ 76552f7358SJed Brown DMPlexPointLocalRef - return read/write access to a point in local array 77552f7358SJed Brown 78552f7358SJed Brown Not Collective 79552f7358SJed Brown 80552f7358SJed Brown Input Arguments: 81552f7358SJed Brown + dm - DM defining topological space 82552f7358SJed Brown . point - topological point 83552f7358SJed Brown - array - array to index into 84552f7358SJed Brown 85552f7358SJed Brown Output Arguments: 86552f7358SJed Brown . ptr - address of reference to point data, type generic so user can place in structure 87552f7358SJed Brown 88552f7358SJed Brown Level: intermediate 89552f7358SJed Brown 90552f7358SJed Brown Note: 91552f7358SJed Brown A common usage when data sizes are known statically: 92552f7358SJed Brown 93552f7358SJed Brown $ struct { PetscScalar foo,bar,baz; } *ptr; 94552f7358SJed Brown $ DMPlexPointLocalRef(dm,point,array,&ptr); 95552f7358SJed Brown $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 96552f7358SJed Brown 97552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 98552f7358SJed Brown @*/ 99552f7358SJed Brown PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 100552f7358SJed Brown { 101552f7358SJed Brown PetscErrorCode ierr; 102a89cf0ddSMatthew G. Knepley PetscInt start, end; 103552f7358SJed Brown 104552f7358SJed Brown PetscFunctionBegin; 105552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 106552f7358SJed Brown PetscValidScalarPointer(array,3); 107552f7358SJed Brown PetscValidPointer(ptr,4); 108*088fd00dSMatthew G. Knepley ierr = DMGetLocalOffset_Private(dm,point,&start,&end);CHKERRQ(ierr); 109a89cf0ddSMatthew G. Knepley *(PetscScalar**)ptr = (start < end) ? array + start : NULL; 110a89cf0ddSMatthew G. Knepley PetscFunctionReturn(0); 111a89cf0ddSMatthew G. Knepley } 112a89cf0ddSMatthew G. Knepley 113a89cf0ddSMatthew G. Knepley /*@ 114a89cf0ddSMatthew G. Knepley DMPlexGetPointLocalField - get location of point field data in local Vec 115a89cf0ddSMatthew G. Knepley 116a89cf0ddSMatthew G. Knepley Not Collective 117a89cf0ddSMatthew G. Knepley 118a89cf0ddSMatthew G. Knepley Input Arguments: 119a89cf0ddSMatthew G. Knepley + dm - DM defining the topological space 120a89cf0ddSMatthew G. Knepley . point - topological point 121a89cf0ddSMatthew G. Knepley - field - the field number 122a89cf0ddSMatthew G. Knepley 123a89cf0ddSMatthew G. Knepley Output Arguments: 124a89cf0ddSMatthew G. Knepley + start - start of point data 125a89cf0ddSMatthew G. Knepley - end - end of point data 126a89cf0ddSMatthew G. Knepley 127a89cf0ddSMatthew G. Knepley Note: This is a half open interval [start, end) 128a89cf0ddSMatthew G. Knepley 129a89cf0ddSMatthew G. Knepley Level: intermediate 130a89cf0ddSMatthew G. Knepley 131a89cf0ddSMatthew G. Knepley .seealso: DMPlexGetPointLocal(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef() 132a89cf0ddSMatthew G. Knepley @*/ 133a89cf0ddSMatthew G. Knepley PetscErrorCode DMPlexGetPointLocalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end) 134a89cf0ddSMatthew G. Knepley { 135a89cf0ddSMatthew G. Knepley PetscInt s, e; 136a89cf0ddSMatthew G. Knepley PetscErrorCode ierr; 137a89cf0ddSMatthew G. Knepley 138a89cf0ddSMatthew G. Knepley PetscFunctionBegin; 139a89cf0ddSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 140a89cf0ddSMatthew G. Knepley if (start) PetscValidPointer(start, 4); 141a89cf0ddSMatthew G. Knepley if (end) PetscValidPointer(end, 5); 142*088fd00dSMatthew G. Knepley ierr = DMGetLocalFieldOffset_Private(dm, point, field, &s, &e);CHKERRQ(ierr); 143a89cf0ddSMatthew G. Knepley if (start) *start = s; 144a89cf0ddSMatthew G. Knepley if (end) *end = e; 145552f7358SJed Brown PetscFunctionReturn(0); 146552f7358SJed Brown } 147552f7358SJed Brown 1481ce3176fSMatthew G. Knepley /*@ 1491ce3176fSMatthew G. Knepley DMPlexPointLocalFieldRead - return read access to a field on a point in local array 1501ce3176fSMatthew G. Knepley 1511ce3176fSMatthew G. Knepley Not Collective 1521ce3176fSMatthew G. Knepley 1531ce3176fSMatthew G. Knepley Input Arguments: 1541ce3176fSMatthew G. Knepley + dm - DM defining topological space 1551ce3176fSMatthew G. Knepley . point - topological point 1561ce3176fSMatthew G. Knepley . field - field number 1571ce3176fSMatthew G. Knepley - array - array to index into 1581ce3176fSMatthew G. Knepley 1591ce3176fSMatthew G. Knepley Output Arguments: 1601ce3176fSMatthew G. Knepley . ptr - address of read reference to point data, type generic so user can place in structure 1611ce3176fSMatthew G. Knepley 1621ce3176fSMatthew G. Knepley Level: intermediate 1631ce3176fSMatthew G. Knepley 1641ce3176fSMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 1651ce3176fSMatthew G. Knepley @*/ 166081a2d76SSatish Balay PetscErrorCode DMPlexPointLocalFieldRead(DM dm, PetscInt point,PetscInt field,const PetscScalar *array,void *ptr) 1671ce3176fSMatthew G. Knepley { 1681ce3176fSMatthew G. Knepley PetscErrorCode ierr; 169a89cf0ddSMatthew G. Knepley PetscInt start, end; 1701ce3176fSMatthew G. Knepley 1711ce3176fSMatthew G. Knepley PetscFunctionBegin; 1721ce3176fSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1731ce3176fSMatthew G. Knepley PetscValidScalarPointer(array,3); 1741ce3176fSMatthew G. Knepley PetscValidPointer(ptr,4); 175*088fd00dSMatthew G. Knepley ierr = DMGetLocalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 1761ce3176fSMatthew G. Knepley *(const PetscScalar**)ptr = array + start; 1771ce3176fSMatthew G. Knepley PetscFunctionReturn(0); 1781ce3176fSMatthew G. Knepley } 1791ce3176fSMatthew G. Knepley 1804824f456SMatthew G. Knepley /*@ 1814824f456SMatthew G. Knepley DMPlexPointLocalFieldRef - return read/write access to a field on a point in local array 1824824f456SMatthew G. Knepley 1834824f456SMatthew G. Knepley Not Collective 1844824f456SMatthew G. Knepley 1854824f456SMatthew G. Knepley Input Arguments: 1864824f456SMatthew G. Knepley + dm - DM defining topological space 1874824f456SMatthew G. Knepley . point - topological point 1884824f456SMatthew G. Knepley . field - field number 1894824f456SMatthew G. Knepley - array - array to index into 1904824f456SMatthew G. Knepley 1914824f456SMatthew G. Knepley Output Arguments: 1924824f456SMatthew G. Knepley . ptr - address of reference to point data, type generic so user can place in structure 1934824f456SMatthew G. Knepley 1944824f456SMatthew G. Knepley Level: intermediate 1954824f456SMatthew G. Knepley 1964824f456SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 1974824f456SMatthew G. Knepley @*/ 1984824f456SMatthew G. Knepley PetscErrorCode DMPlexPointLocalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr) 1994824f456SMatthew G. Knepley { 2004824f456SMatthew G. Knepley PetscErrorCode ierr; 201a89cf0ddSMatthew G. Knepley PetscInt start, end; 2024824f456SMatthew G. Knepley 2034824f456SMatthew G. Knepley PetscFunctionBegin; 2044824f456SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2054824f456SMatthew G. Knepley PetscValidScalarPointer(array,3); 2064824f456SMatthew G. Knepley PetscValidPointer(ptr,4); 207*088fd00dSMatthew G. Knepley ierr = DMGetLocalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 2084824f456SMatthew G. Knepley *(PetscScalar**)ptr = array + start; 2094824f456SMatthew G. Knepley PetscFunctionReturn(0); 2104824f456SMatthew G. Knepley } 2114824f456SMatthew G. Knepley 212552f7358SJed Brown /*@ 213552f7358SJed Brown DMPlexGetPointGlobal - get location of point data in global Vec 214552f7358SJed Brown 215552f7358SJed Brown Not Collective 216552f7358SJed Brown 217552f7358SJed Brown Input Arguments: 218552f7358SJed Brown + dm - DM defining the topological space 219552f7358SJed Brown - point - topological point 220552f7358SJed Brown 221552f7358SJed Brown Output Arguments: 222a89cf0ddSMatthew G. Knepley + start - start of point data; returns -(globalStart+1) if point is not owned 223a89cf0ddSMatthew G. Knepley - end - end of point data; returns -(globalEnd+1) if point is not owned 224a89cf0ddSMatthew G. Knepley 225a89cf0ddSMatthew G. Knepley Note: This is a half open interval [start, end) 226552f7358SJed Brown 227552f7358SJed Brown Level: intermediate 228552f7358SJed Brown 229a89cf0ddSMatthew G. Knepley .seealso: DMPlexGetPointGlobalField(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef() 230552f7358SJed Brown @*/ 231552f7358SJed Brown PetscErrorCode DMPlexGetPointGlobal(DM dm, PetscInt point, PetscInt *start, PetscInt *end) 232552f7358SJed Brown { 233a89cf0ddSMatthew G. Knepley PetscInt s, e; 234552f7358SJed Brown PetscErrorCode ierr; 235552f7358SJed Brown 236552f7358SJed Brown PetscFunctionBegin; 237552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 238a89cf0ddSMatthew G. Knepley if (start) PetscValidPointer(start, 3); 239a89cf0ddSMatthew G. Knepley if (end) PetscValidPointer(end, 4); 240*088fd00dSMatthew G. Knepley ierr = DMGetGlobalOffset_Private(dm, point, &s, &e);CHKERRQ(ierr); 241a89cf0ddSMatthew G. Knepley if (start) *start = s; 242a89cf0ddSMatthew G. Knepley if (end) *end = e; 243552f7358SJed Brown PetscFunctionReturn(0); 244552f7358SJed Brown } 245552f7358SJed Brown 246552f7358SJed Brown /*@ 247552f7358SJed Brown DMPlexPointGlobalRead - return read access to a point in global array 248552f7358SJed Brown 249552f7358SJed Brown Not Collective 250552f7358SJed Brown 251552f7358SJed Brown Input Arguments: 252552f7358SJed Brown + dm - DM defining topological space 253552f7358SJed Brown . point - topological point 254552f7358SJed Brown - array - array to index into 255552f7358SJed Brown 256552f7358SJed Brown Output Arguments: 2570298fd71SBarry 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 258552f7358SJed Brown 259552f7358SJed Brown Level: intermediate 260552f7358SJed Brown 261552f7358SJed Brown Note: 262552f7358SJed Brown A common usage when data sizes are known statically: 263552f7358SJed Brown 264552f7358SJed Brown $ const struct { PetscScalar foo,bar,baz; } *ptr; 265552f7358SJed Brown $ DMPlexPointGlobalRead(dm,point,array,&ptr); 266552f7358SJed Brown $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 267552f7358SJed Brown 268552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef() 269552f7358SJed Brown @*/ 270552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 271552f7358SJed Brown { 27279532bb4SMatthew G. Knepley PetscInt start, end; 273552f7358SJed Brown PetscErrorCode ierr; 274552f7358SJed Brown 275552f7358SJed Brown PetscFunctionBegin; 276552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 277552f7358SJed Brown PetscValidScalarPointer(array, 3); 278552f7358SJed Brown PetscValidPointer(ptr, 4); 279*088fd00dSMatthew G. Knepley ierr = DMGetGlobalOffset_Private(dm, point, &start, &end);CHKERRQ(ierr); 28079532bb4SMatthew G. Knepley *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 281552f7358SJed Brown PetscFunctionReturn(0); 282552f7358SJed Brown } 283552f7358SJed Brown 284552f7358SJed Brown /*@ 285552f7358SJed Brown DMPlexPointGlobalRef - return read/write access to a point in global array 286552f7358SJed Brown 287552f7358SJed Brown Not Collective 288552f7358SJed Brown 289552f7358SJed Brown Input Arguments: 290552f7358SJed Brown + dm - DM defining topological space 291552f7358SJed Brown . point - topological point 292552f7358SJed Brown - array - array to index into 293552f7358SJed Brown 294552f7358SJed Brown Output Arguments: 2950298fd71SBarry Smith . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 296552f7358SJed Brown 297552f7358SJed Brown Level: intermediate 298552f7358SJed Brown 299552f7358SJed Brown Note: 300552f7358SJed Brown A common usage when data sizes are known statically: 301552f7358SJed Brown 302552f7358SJed Brown $ struct { PetscScalar foo,bar,baz; } *ptr; 303552f7358SJed Brown $ DMPlexPointGlobalRef(dm,point,array,&ptr); 304552f7358SJed Brown $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 305552f7358SJed Brown 306552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead() 307552f7358SJed Brown @*/ 308552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 309552f7358SJed Brown { 31079532bb4SMatthew G. Knepley PetscInt start, end; 311552f7358SJed Brown PetscErrorCode ierr; 312552f7358SJed Brown 313552f7358SJed Brown PetscFunctionBegin; 314552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 315552f7358SJed Brown PetscValidScalarPointer(array, 3); 316552f7358SJed Brown PetscValidPointer(ptr, 4); 317*088fd00dSMatthew G. Knepley ierr = DMGetGlobalOffset_Private(dm, point, &start, &end);CHKERRQ(ierr); 31879532bb4SMatthew G. Knepley *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 319552f7358SJed Brown PetscFunctionReturn(0); 320552f7358SJed Brown } 32133879625SMatthew G. Knepley 322a89cf0ddSMatthew G. Knepley /*@ 323a89cf0ddSMatthew G. Knepley DMPlexGetPointGlobalField - get location of point field data in global Vec 324a89cf0ddSMatthew G. Knepley 325a89cf0ddSMatthew G. Knepley Not Collective 326a89cf0ddSMatthew G. Knepley 327a89cf0ddSMatthew G. Knepley Input Arguments: 328a89cf0ddSMatthew G. Knepley + dm - DM defining the topological space 329a89cf0ddSMatthew G. Knepley . point - topological point 330a89cf0ddSMatthew G. Knepley - field - the field number 331a89cf0ddSMatthew G. Knepley 332a89cf0ddSMatthew G. Knepley Output Arguments: 333a89cf0ddSMatthew G. Knepley + start - start of point data; returns -(globalStart+1) if point is not owned 334a89cf0ddSMatthew G. Knepley - end - end of point data; returns -(globalEnd+1) if point is not owned 335a89cf0ddSMatthew G. Knepley 336a89cf0ddSMatthew G. Knepley Note: This is a half open interval [start, end) 337a89cf0ddSMatthew G. Knepley 338a89cf0ddSMatthew G. Knepley Level: intermediate 339a89cf0ddSMatthew G. Knepley 340a89cf0ddSMatthew G. Knepley .seealso: DMPlexGetPointGlobal(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef() 341a89cf0ddSMatthew G. Knepley @*/ 342a89cf0ddSMatthew G. Knepley PetscErrorCode DMPlexGetPointGlobalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end) 343a89cf0ddSMatthew G. Knepley { 344a89cf0ddSMatthew G. Knepley PetscInt s, e; 345a89cf0ddSMatthew G. Knepley PetscErrorCode ierr; 346a89cf0ddSMatthew G. Knepley 347a89cf0ddSMatthew G. Knepley PetscFunctionBegin; 348a89cf0ddSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 349a89cf0ddSMatthew G. Knepley if (start) PetscValidPointer(start, 4); 350a89cf0ddSMatthew G. Knepley if (end) PetscValidPointer(end, 5); 351*088fd00dSMatthew G. Knepley ierr = DMGetGlobalFieldOffset_Private(dm, point, field, &s, &e);CHKERRQ(ierr); 352a89cf0ddSMatthew G. Knepley if (start) *start = s; 353a89cf0ddSMatthew G. Knepley if (end) *end = e; 354a89cf0ddSMatthew G. Knepley PetscFunctionReturn(0); 355a89cf0ddSMatthew G. Knepley } 356a89cf0ddSMatthew G. Knepley 35733879625SMatthew G. Knepley /*@ 35833879625SMatthew G. Knepley DMPlexPointGlobalFieldRead - return read access to a field on a point in global array 35933879625SMatthew G. Knepley 36033879625SMatthew G. Knepley Not Collective 36133879625SMatthew G. Knepley 36233879625SMatthew G. Knepley Input Arguments: 36333879625SMatthew G. Knepley + dm - DM defining topological space 36433879625SMatthew G. Knepley . point - topological point 36533879625SMatthew G. Knepley . field - field number 36633879625SMatthew G. Knepley - array - array to index into 36733879625SMatthew G. Knepley 36833879625SMatthew G. Knepley Output Arguments: 36933879625SMatthew G. Knepley . ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 37033879625SMatthew G. Knepley 37133879625SMatthew G. Knepley Level: intermediate 37233879625SMatthew G. Knepley 37333879625SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef() 37433879625SMatthew G. Knepley @*/ 375081a2d76SSatish Balay PetscErrorCode DMPlexPointGlobalFieldRead(DM dm,PetscInt point,PetscInt field,const PetscScalar *array,void *ptr) 37633879625SMatthew G. Knepley { 37779532bb4SMatthew G. Knepley PetscInt start, end; 37833879625SMatthew G. Knepley PetscErrorCode ierr; 37933879625SMatthew G. Knepley 38033879625SMatthew G. Knepley PetscFunctionBegin; 38133879625SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 38233879625SMatthew G. Knepley PetscValidScalarPointer(array, 3); 38333879625SMatthew G. Knepley PetscValidPointer(ptr, 4); 384*088fd00dSMatthew G. Knepley ierr = DMGetGlobalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 38579532bb4SMatthew G. Knepley *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 38633879625SMatthew G. Knepley PetscFunctionReturn(0); 38733879625SMatthew G. Knepley } 38833879625SMatthew G. Knepley 38933879625SMatthew G. Knepley /*@ 39033879625SMatthew G. Knepley DMPlexPointGlobalFieldRef - return read/write access to a field on a point in global array 39133879625SMatthew G. Knepley 39233879625SMatthew G. Knepley Not Collective 39333879625SMatthew G. Knepley 39433879625SMatthew G. Knepley Input Arguments: 39533879625SMatthew G. Knepley + dm - DM defining topological space 39633879625SMatthew G. Knepley . point - topological point 39733879625SMatthew G. Knepley . field - field number 39833879625SMatthew G. Knepley - array - array to index into 39933879625SMatthew G. Knepley 40033879625SMatthew G. Knepley Output Arguments: 40133879625SMatthew G. Knepley . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 40233879625SMatthew G. Knepley 40333879625SMatthew G. Knepley Level: intermediate 40433879625SMatthew G. Knepley 40533879625SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead() 40633879625SMatthew G. Knepley @*/ 40733879625SMatthew G. Knepley PetscErrorCode DMPlexPointGlobalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr) 40833879625SMatthew G. Knepley { 40979532bb4SMatthew G. Knepley PetscInt start, end; 41033879625SMatthew G. Knepley PetscErrorCode ierr; 41133879625SMatthew G. Knepley 41233879625SMatthew G. Knepley PetscFunctionBegin; 41333879625SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 41433879625SMatthew G. Knepley PetscValidScalarPointer(array, 3); 41533879625SMatthew G. Knepley PetscValidPointer(ptr, 4); 416*088fd00dSMatthew G. Knepley ierr = DMGetGlobalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 41779532bb4SMatthew G. Knepley *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 41833879625SMatthew G. Knepley PetscFunctionReturn(0); 41933879625SMatthew G. Knepley } 420