134541f0dSBarry Smith #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 233879625SMatthew G. Knepley 333879625SMatthew G. Knepley #undef __FUNCT__ 4552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointLocal" 5552f7358SJed Brown /*@ 6552f7358SJed Brown DMPlexGetPointLocal - get location of point data in local Vec 7552f7358SJed Brown 8552f7358SJed Brown Not Collective 9552f7358SJed Brown 10552f7358SJed Brown Input Arguments: 11552f7358SJed Brown + dm - DM defining the topological space 12552f7358SJed Brown - point - topological point 13552f7358SJed Brown 14552f7358SJed Brown Output Arguments: 15552f7358SJed Brown + start - start of point data 16552f7358SJed Brown - end - end of point data 17552f7358SJed Brown 18*a89cf0ddSMatthew G. Knepley Note: This is a half open interval [start, end) 19*a89cf0ddSMatthew G. Knepley 20552f7358SJed Brown Level: intermediate 21552f7358SJed Brown 22*a89cf0ddSMatthew G. Knepley .seealso: DMPlexGetPointLocalField(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef() 23552f7358SJed Brown @*/ 24552f7358SJed Brown PetscErrorCode DMPlexGetPointLocal(DM dm, PetscInt point, PetscInt *start, PetscInt *end) 25552f7358SJed Brown { 26*a89cf0ddSMatthew G. Knepley PetscInt s, e; 27552f7358SJed Brown PetscErrorCode ierr; 28552f7358SJed Brown 29552f7358SJed Brown PetscFunctionBegin; 30552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 31*a89cf0ddSMatthew G. Knepley if (start) PetscValidPointer(start, 3); 32*a89cf0ddSMatthew G. Knepley if (end) PetscValidPointer(end, 4); 33*a89cf0ddSMatthew G. Knepley ierr = DMPlexGetLocalOffset_Private(dm, point, &s, &e);CHKERRQ(ierr); 34*a89cf0ddSMatthew G. Knepley if (start) *start = s; 35*a89cf0ddSMatthew G. Knepley if (end) *end = e; 36552f7358SJed Brown PetscFunctionReturn(0); 37552f7358SJed Brown } 38552f7358SJed Brown 39552f7358SJed Brown #undef __FUNCT__ 40552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRead" 41552f7358SJed Brown /*@ 42552f7358SJed Brown DMPlexPointLocalRead - return read access to a point in local array 43552f7358SJed Brown 44552f7358SJed Brown Not Collective 45552f7358SJed Brown 46552f7358SJed Brown Input Arguments: 47552f7358SJed Brown + dm - DM defining topological space 48552f7358SJed Brown . point - topological point 49552f7358SJed Brown - array - array to index into 50552f7358SJed Brown 51552f7358SJed Brown Output Arguments: 52552f7358SJed Brown . ptr - address of read reference to point data, type generic so user can place in structure 53552f7358SJed Brown 54552f7358SJed Brown Level: intermediate 55552f7358SJed Brown 56552f7358SJed Brown Note: 57552f7358SJed Brown A common usage when data sizes are known statically: 58552f7358SJed Brown 59552f7358SJed Brown $ const struct { PetscScalar foo,bar,baz; } *ptr; 60552f7358SJed Brown $ DMPlexPointLocalRead(dm,point,array,&ptr); 61552f7358SJed Brown $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 62552f7358SJed Brown 63552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead() 64552f7358SJed Brown @*/ 65552f7358SJed Brown PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 66552f7358SJed Brown { 67552f7358SJed Brown PetscErrorCode ierr; 68*a89cf0ddSMatthew G. Knepley PetscInt start, end; 69552f7358SJed Brown 70552f7358SJed Brown PetscFunctionBegin; 71552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 72552f7358SJed Brown PetscValidScalarPointer(array,3); 73552f7358SJed Brown PetscValidPointer(ptr,4); 74*a89cf0ddSMatthew G. Knepley ierr = DMPlexGetLocalOffset_Private(dm,point,&start,&end);CHKERRQ(ierr); 75*a89cf0ddSMatthew G. Knepley *(const PetscScalar**)ptr = (start < end) ? array + start : NULL; 76552f7358SJed Brown PetscFunctionReturn(0); 77552f7358SJed Brown } 78552f7358SJed Brown 79552f7358SJed Brown #undef __FUNCT__ 80552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRef" 81552f7358SJed Brown /*@ 82552f7358SJed Brown DMPlexPointLocalRef - return read/write access to a point in local array 83552f7358SJed Brown 84552f7358SJed Brown Not Collective 85552f7358SJed Brown 86552f7358SJed Brown Input Arguments: 87552f7358SJed Brown + dm - DM defining topological space 88552f7358SJed Brown . point - topological point 89552f7358SJed Brown - array - array to index into 90552f7358SJed Brown 91552f7358SJed Brown Output Arguments: 92552f7358SJed Brown . ptr - address of reference to point data, type generic so user can place in structure 93552f7358SJed Brown 94552f7358SJed Brown Level: intermediate 95552f7358SJed Brown 96552f7358SJed Brown Note: 97552f7358SJed Brown A common usage when data sizes are known statically: 98552f7358SJed Brown 99552f7358SJed Brown $ struct { PetscScalar foo,bar,baz; } *ptr; 100552f7358SJed Brown $ DMPlexPointLocalRef(dm,point,array,&ptr); 101552f7358SJed Brown $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 102552f7358SJed Brown 103552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 104552f7358SJed Brown @*/ 105552f7358SJed Brown PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 106552f7358SJed Brown { 107552f7358SJed Brown PetscErrorCode ierr; 108*a89cf0ddSMatthew G. Knepley PetscInt start, end; 109552f7358SJed Brown 110552f7358SJed Brown PetscFunctionBegin; 111552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 112552f7358SJed Brown PetscValidScalarPointer(array,3); 113552f7358SJed Brown PetscValidPointer(ptr,4); 114*a89cf0ddSMatthew G. Knepley ierr = DMPlexGetLocalOffset_Private(dm,point,&start,&end);CHKERRQ(ierr); 115*a89cf0ddSMatthew G. Knepley *(PetscScalar**)ptr = (start < end) ? array + start : NULL; 116*a89cf0ddSMatthew G. Knepley PetscFunctionReturn(0); 117*a89cf0ddSMatthew G. Knepley } 118*a89cf0ddSMatthew G. Knepley 119*a89cf0ddSMatthew G. Knepley #undef __FUNCT__ 120*a89cf0ddSMatthew G. Knepley #define __FUNCT__ "DMPlexGetPointLocalField" 121*a89cf0ddSMatthew G. Knepley /*@ 122*a89cf0ddSMatthew G. Knepley DMPlexGetPointLocalField - get location of point field data in local Vec 123*a89cf0ddSMatthew G. Knepley 124*a89cf0ddSMatthew G. Knepley Not Collective 125*a89cf0ddSMatthew G. Knepley 126*a89cf0ddSMatthew G. Knepley Input Arguments: 127*a89cf0ddSMatthew G. Knepley + dm - DM defining the topological space 128*a89cf0ddSMatthew G. Knepley . point - topological point 129*a89cf0ddSMatthew G. Knepley - field - the field number 130*a89cf0ddSMatthew G. Knepley 131*a89cf0ddSMatthew G. Knepley Output Arguments: 132*a89cf0ddSMatthew G. Knepley + start - start of point data 133*a89cf0ddSMatthew G. Knepley - end - end of point data 134*a89cf0ddSMatthew G. Knepley 135*a89cf0ddSMatthew G. Knepley Note: This is a half open interval [start, end) 136*a89cf0ddSMatthew G. Knepley 137*a89cf0ddSMatthew G. Knepley Level: intermediate 138*a89cf0ddSMatthew G. Knepley 139*a89cf0ddSMatthew G. Knepley .seealso: DMPlexGetPointLocal(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef() 140*a89cf0ddSMatthew G. Knepley @*/ 141*a89cf0ddSMatthew G. Knepley PetscErrorCode DMPlexGetPointLocalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end) 142*a89cf0ddSMatthew G. Knepley { 143*a89cf0ddSMatthew G. Knepley PetscInt s, e; 144*a89cf0ddSMatthew G. Knepley PetscErrorCode ierr; 145*a89cf0ddSMatthew G. Knepley 146*a89cf0ddSMatthew G. Knepley PetscFunctionBegin; 147*a89cf0ddSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 148*a89cf0ddSMatthew G. Knepley if (start) PetscValidPointer(start, 4); 149*a89cf0ddSMatthew G. Knepley if (end) PetscValidPointer(end, 5); 150*a89cf0ddSMatthew G. Knepley ierr = DMPlexGetLocalFieldOffset_Private(dm, point, field, &s, &e);CHKERRQ(ierr); 151*a89cf0ddSMatthew G. Knepley if (start) *start = s; 152*a89cf0ddSMatthew G. Knepley if (end) *end = e; 153552f7358SJed Brown PetscFunctionReturn(0); 154552f7358SJed Brown } 155552f7358SJed Brown 156552f7358SJed Brown #undef __FUNCT__ 1571ce3176fSMatthew G. Knepley #define __FUNCT__ "DMPlexPointLocalFieldRead" 1581ce3176fSMatthew G. Knepley /*@ 1591ce3176fSMatthew G. Knepley DMPlexPointLocalFieldRead - return read access to a field on a point in local array 1601ce3176fSMatthew G. Knepley 1611ce3176fSMatthew G. Knepley Not Collective 1621ce3176fSMatthew G. Knepley 1631ce3176fSMatthew G. Knepley Input Arguments: 1641ce3176fSMatthew G. Knepley + dm - DM defining topological space 1651ce3176fSMatthew G. Knepley . point - topological point 1661ce3176fSMatthew G. Knepley . field - field number 1671ce3176fSMatthew G. Knepley - array - array to index into 1681ce3176fSMatthew G. Knepley 1691ce3176fSMatthew G. Knepley Output Arguments: 1701ce3176fSMatthew G. Knepley . ptr - address of read reference to point data, type generic so user can place in structure 1711ce3176fSMatthew G. Knepley 1721ce3176fSMatthew G. Knepley Level: intermediate 1731ce3176fSMatthew G. Knepley 1741ce3176fSMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 1751ce3176fSMatthew G. Knepley @*/ 1761ce3176fSMatthew G. Knepley PetscErrorCode DMPlexPointLocalFieldRead(DM dm, PetscInt point,PetscInt field,const PetscScalar *array,const void *ptr) 1771ce3176fSMatthew G. Knepley { 1781ce3176fSMatthew G. Knepley PetscErrorCode ierr; 179*a89cf0ddSMatthew G. Knepley PetscInt start, end; 1801ce3176fSMatthew G. Knepley 1811ce3176fSMatthew G. Knepley PetscFunctionBegin; 1821ce3176fSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1831ce3176fSMatthew G. Knepley PetscValidScalarPointer(array,3); 1841ce3176fSMatthew G. Knepley PetscValidPointer(ptr,4); 185*a89cf0ddSMatthew G. Knepley ierr = DMPlexGetLocalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 1861ce3176fSMatthew G. Knepley *(const PetscScalar**)ptr = array + start; 1871ce3176fSMatthew G. Knepley PetscFunctionReturn(0); 1881ce3176fSMatthew G. Knepley } 1891ce3176fSMatthew G. Knepley 1901ce3176fSMatthew G. Knepley #undef __FUNCT__ 1914824f456SMatthew G. Knepley #define __FUNCT__ "DMPlexPointLocalFieldRef" 1924824f456SMatthew G. Knepley /*@ 1934824f456SMatthew G. Knepley DMPlexPointLocalFieldRef - return read/write access to a field on a point in local array 1944824f456SMatthew G. Knepley 1954824f456SMatthew G. Knepley Not Collective 1964824f456SMatthew G. Knepley 1974824f456SMatthew G. Knepley Input Arguments: 1984824f456SMatthew G. Knepley + dm - DM defining topological space 1994824f456SMatthew G. Knepley . point - topological point 2004824f456SMatthew G. Knepley . field - field number 2014824f456SMatthew G. Knepley - array - array to index into 2024824f456SMatthew G. Knepley 2034824f456SMatthew G. Knepley Output Arguments: 2044824f456SMatthew G. Knepley . ptr - address of reference to point data, type generic so user can place in structure 2054824f456SMatthew G. Knepley 2064824f456SMatthew G. Knepley Level: intermediate 2074824f456SMatthew G. Knepley 2084824f456SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 2094824f456SMatthew G. Knepley @*/ 2104824f456SMatthew G. Knepley PetscErrorCode DMPlexPointLocalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr) 2114824f456SMatthew G. Knepley { 2124824f456SMatthew G. Knepley PetscErrorCode ierr; 213*a89cf0ddSMatthew G. Knepley PetscInt start, end; 2144824f456SMatthew G. Knepley 2154824f456SMatthew G. Knepley PetscFunctionBegin; 2164824f456SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2174824f456SMatthew G. Knepley PetscValidScalarPointer(array,3); 2184824f456SMatthew G. Knepley PetscValidPointer(ptr,4); 219*a89cf0ddSMatthew G. Knepley ierr = DMPlexGetLocalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 2204824f456SMatthew G. Knepley *(PetscScalar**)ptr = array + start; 2214824f456SMatthew G. Knepley PetscFunctionReturn(0); 2224824f456SMatthew G. Knepley } 2234824f456SMatthew G. Knepley 2244824f456SMatthew G. Knepley #undef __FUNCT__ 225552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointGlobal" 226552f7358SJed Brown /*@ 227552f7358SJed Brown DMPlexGetPointGlobal - get location of point data in global Vec 228552f7358SJed Brown 229552f7358SJed Brown Not Collective 230552f7358SJed Brown 231552f7358SJed Brown Input Arguments: 232552f7358SJed Brown + dm - DM defining the topological space 233552f7358SJed Brown - point - topological point 234552f7358SJed Brown 235552f7358SJed Brown Output Arguments: 236*a89cf0ddSMatthew G. Knepley + start - start of point data; returns -(globalStart+1) if point is not owned 237*a89cf0ddSMatthew G. Knepley - end - end of point data; returns -(globalEnd+1) if point is not owned 238*a89cf0ddSMatthew G. Knepley 239*a89cf0ddSMatthew G. Knepley Note: This is a half open interval [start, end) 240552f7358SJed Brown 241552f7358SJed Brown Level: intermediate 242552f7358SJed Brown 243*a89cf0ddSMatthew G. Knepley .seealso: DMPlexGetPointGlobalField(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef() 244552f7358SJed Brown @*/ 245552f7358SJed Brown PetscErrorCode DMPlexGetPointGlobal(DM dm, PetscInt point, PetscInt *start, PetscInt *end) 246552f7358SJed Brown { 247*a89cf0ddSMatthew G. Knepley PetscInt s, e; 248552f7358SJed Brown PetscErrorCode ierr; 249552f7358SJed Brown 250552f7358SJed Brown PetscFunctionBegin; 251552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 252*a89cf0ddSMatthew G. Knepley if (start) PetscValidPointer(start, 3); 253*a89cf0ddSMatthew G. Knepley if (end) PetscValidPointer(end, 4); 254*a89cf0ddSMatthew G. Knepley ierr = DMPlexGetGlobalOffset_Private(dm, point, &s, &e);CHKERRQ(ierr); 255*a89cf0ddSMatthew G. Knepley if (start) *start = s; 256*a89cf0ddSMatthew G. Knepley if (end) *end = e; 257552f7358SJed Brown PetscFunctionReturn(0); 258552f7358SJed Brown } 259552f7358SJed Brown 260552f7358SJed Brown #undef __FUNCT__ 261552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRead" 262552f7358SJed Brown /*@ 263552f7358SJed Brown DMPlexPointGlobalRead - return read access to a point in global array 264552f7358SJed Brown 265552f7358SJed Brown Not Collective 266552f7358SJed Brown 267552f7358SJed Brown Input Arguments: 268552f7358SJed Brown + dm - DM defining topological space 269552f7358SJed Brown . point - topological point 270552f7358SJed Brown - array - array to index into 271552f7358SJed Brown 272552f7358SJed Brown Output Arguments: 2730298fd71SBarry 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 274552f7358SJed Brown 275552f7358SJed Brown Level: intermediate 276552f7358SJed Brown 277552f7358SJed Brown Note: 278552f7358SJed Brown A common usage when data sizes are known statically: 279552f7358SJed Brown 280552f7358SJed Brown $ const struct { PetscScalar foo,bar,baz; } *ptr; 281552f7358SJed Brown $ DMPlexPointGlobalRead(dm,point,array,&ptr); 282552f7358SJed Brown $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 283552f7358SJed Brown 284552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef() 285552f7358SJed Brown @*/ 286552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 287552f7358SJed Brown { 28879532bb4SMatthew G. Knepley PetscInt start, end; 289552f7358SJed Brown PetscErrorCode ierr; 290552f7358SJed Brown 291552f7358SJed Brown PetscFunctionBegin; 292552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 293552f7358SJed Brown PetscValidScalarPointer(array, 3); 294552f7358SJed Brown PetscValidPointer(ptr, 4); 29579532bb4SMatthew G. Knepley ierr = DMPlexGetGlobalOffset_Private(dm, point, &start, &end);CHKERRQ(ierr); 29679532bb4SMatthew G. Knepley *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 297552f7358SJed Brown PetscFunctionReturn(0); 298552f7358SJed Brown } 299552f7358SJed Brown 300552f7358SJed Brown #undef __FUNCT__ 301552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRef" 302552f7358SJed Brown /*@ 303552f7358SJed Brown DMPlexPointGlobalRef - return read/write access to a point in global array 304552f7358SJed Brown 305552f7358SJed Brown Not Collective 306552f7358SJed Brown 307552f7358SJed Brown Input Arguments: 308552f7358SJed Brown + dm - DM defining topological space 309552f7358SJed Brown . point - topological point 310552f7358SJed Brown - array - array to index into 311552f7358SJed Brown 312552f7358SJed Brown Output Arguments: 3130298fd71SBarry Smith . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 314552f7358SJed Brown 315552f7358SJed Brown Level: intermediate 316552f7358SJed Brown 317552f7358SJed Brown Note: 318552f7358SJed Brown A common usage when data sizes are known statically: 319552f7358SJed Brown 320552f7358SJed Brown $ struct { PetscScalar foo,bar,baz; } *ptr; 321552f7358SJed Brown $ DMPlexPointGlobalRef(dm,point,array,&ptr); 322552f7358SJed Brown $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 323552f7358SJed Brown 324552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead() 325552f7358SJed Brown @*/ 326552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 327552f7358SJed Brown { 32879532bb4SMatthew G. Knepley PetscInt start, end; 329552f7358SJed Brown PetscErrorCode ierr; 330552f7358SJed Brown 331552f7358SJed Brown PetscFunctionBegin; 332552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 333552f7358SJed Brown PetscValidScalarPointer(array, 3); 334552f7358SJed Brown PetscValidPointer(ptr, 4); 33579532bb4SMatthew G. Knepley ierr = DMPlexGetGlobalOffset_Private(dm, point, &start, &end);CHKERRQ(ierr); 33679532bb4SMatthew G. Knepley *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 337552f7358SJed Brown PetscFunctionReturn(0); 338552f7358SJed Brown } 33933879625SMatthew G. Knepley 34033879625SMatthew G. Knepley #undef __FUNCT__ 341*a89cf0ddSMatthew G. Knepley #define __FUNCT__ "DMPlexGetPointGlobalField" 342*a89cf0ddSMatthew G. Knepley /*@ 343*a89cf0ddSMatthew G. Knepley DMPlexGetPointGlobalField - get location of point field data in global Vec 344*a89cf0ddSMatthew G. Knepley 345*a89cf0ddSMatthew G. Knepley Not Collective 346*a89cf0ddSMatthew G. Knepley 347*a89cf0ddSMatthew G. Knepley Input Arguments: 348*a89cf0ddSMatthew G. Knepley + dm - DM defining the topological space 349*a89cf0ddSMatthew G. Knepley . point - topological point 350*a89cf0ddSMatthew G. Knepley - field - the field number 351*a89cf0ddSMatthew G. Knepley 352*a89cf0ddSMatthew G. Knepley Output Arguments: 353*a89cf0ddSMatthew G. Knepley + start - start of point data; returns -(globalStart+1) if point is not owned 354*a89cf0ddSMatthew G. Knepley - end - end of point data; returns -(globalEnd+1) if point is not owned 355*a89cf0ddSMatthew G. Knepley 356*a89cf0ddSMatthew G. Knepley Note: This is a half open interval [start, end) 357*a89cf0ddSMatthew G. Knepley 358*a89cf0ddSMatthew G. Knepley Level: intermediate 359*a89cf0ddSMatthew G. Knepley 360*a89cf0ddSMatthew G. Knepley .seealso: DMPlexGetPointGlobal(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef() 361*a89cf0ddSMatthew G. Knepley @*/ 362*a89cf0ddSMatthew G. Knepley PetscErrorCode DMPlexGetPointGlobalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end) 363*a89cf0ddSMatthew G. Knepley { 364*a89cf0ddSMatthew G. Knepley PetscInt s, e; 365*a89cf0ddSMatthew G. Knepley PetscErrorCode ierr; 366*a89cf0ddSMatthew G. Knepley 367*a89cf0ddSMatthew G. Knepley PetscFunctionBegin; 368*a89cf0ddSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 369*a89cf0ddSMatthew G. Knepley if (start) PetscValidPointer(start, 4); 370*a89cf0ddSMatthew G. Knepley if (end) PetscValidPointer(end, 5); 371*a89cf0ddSMatthew G. Knepley ierr = DMPlexGetGlobalFieldOffset_Private(dm, point, field, &s, &e);CHKERRQ(ierr); 372*a89cf0ddSMatthew G. Knepley if (start) *start = s; 373*a89cf0ddSMatthew G. Knepley if (end) *end = e; 374*a89cf0ddSMatthew G. Knepley PetscFunctionReturn(0); 375*a89cf0ddSMatthew G. Knepley } 376*a89cf0ddSMatthew G. Knepley 377*a89cf0ddSMatthew G. Knepley #undef __FUNCT__ 37833879625SMatthew G. Knepley #define __FUNCT__ "DMPlexPointGlobalFieldRead" 37933879625SMatthew G. Knepley /*@ 38033879625SMatthew G. Knepley DMPlexPointGlobalFieldRead - return read access to a field on a point in global array 38133879625SMatthew G. Knepley 38233879625SMatthew G. Knepley Not Collective 38333879625SMatthew G. Knepley 38433879625SMatthew G. Knepley Input Arguments: 38533879625SMatthew G. Knepley + dm - DM defining topological space 38633879625SMatthew G. Knepley . point - topological point 38733879625SMatthew G. Knepley . field - field number 38833879625SMatthew G. Knepley - array - array to index into 38933879625SMatthew G. Knepley 39033879625SMatthew G. Knepley Output Arguments: 39133879625SMatthew 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 39233879625SMatthew G. Knepley 39333879625SMatthew G. Knepley Level: intermediate 39433879625SMatthew G. Knepley 39533879625SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef() 39633879625SMatthew G. Knepley @*/ 39733879625SMatthew G. Knepley PetscErrorCode DMPlexPointGlobalFieldRead(DM dm,PetscInt point,PetscInt field,const PetscScalar *array,const void *ptr) 39833879625SMatthew G. Knepley { 39979532bb4SMatthew G. Knepley PetscInt start, end; 40033879625SMatthew G. Knepley PetscErrorCode ierr; 40133879625SMatthew G. Knepley 40233879625SMatthew G. Knepley PetscFunctionBegin; 40333879625SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 40433879625SMatthew G. Knepley PetscValidScalarPointer(array, 3); 40533879625SMatthew G. Knepley PetscValidPointer(ptr, 4); 40679532bb4SMatthew G. Knepley ierr = DMPlexGetGlobalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 40779532bb4SMatthew G. Knepley *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 40833879625SMatthew G. Knepley PetscFunctionReturn(0); 40933879625SMatthew G. Knepley } 41033879625SMatthew G. Knepley 41133879625SMatthew G. Knepley #undef __FUNCT__ 41233879625SMatthew G. Knepley #define __FUNCT__ "DMPlexPointGlobalFieldRef" 41333879625SMatthew G. Knepley /*@ 41433879625SMatthew G. Knepley DMPlexPointGlobalFieldRef - return read/write access to a field on a point in global array 41533879625SMatthew G. Knepley 41633879625SMatthew G. Knepley Not Collective 41733879625SMatthew G. Knepley 41833879625SMatthew G. Knepley Input Arguments: 41933879625SMatthew G. Knepley + dm - DM defining topological space 42033879625SMatthew G. Knepley . point - topological point 42133879625SMatthew G. Knepley . field - field number 42233879625SMatthew G. Knepley - array - array to index into 42333879625SMatthew G. Knepley 42433879625SMatthew G. Knepley Output Arguments: 42533879625SMatthew 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 42633879625SMatthew G. Knepley 42733879625SMatthew G. Knepley Level: intermediate 42833879625SMatthew G. Knepley 42933879625SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead() 43033879625SMatthew G. Knepley @*/ 43133879625SMatthew G. Knepley PetscErrorCode DMPlexPointGlobalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr) 43233879625SMatthew G. Knepley { 43379532bb4SMatthew G. Knepley PetscInt start, end; 43433879625SMatthew G. Knepley PetscErrorCode ierr; 43533879625SMatthew G. Knepley 43633879625SMatthew G. Knepley PetscFunctionBegin; 43733879625SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 43833879625SMatthew G. Knepley PetscValidScalarPointer(array, 3); 43933879625SMatthew G. Knepley PetscValidPointer(ptr, 4); 44079532bb4SMatthew G. Knepley ierr = DMPlexGetGlobalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 44179532bb4SMatthew G. Knepley *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 44233879625SMatthew G. Knepley PetscFunctionReturn(0); 44333879625SMatthew G. Knepley } 444