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 18552f7358SJed Brown Level: intermediate 19552f7358SJed Brown 20552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef() 21552f7358SJed Brown @*/ 22552f7358SJed Brown PetscErrorCode DMPlexGetPointLocal(DM dm,PetscInt point,PetscInt *start,PetscInt *end) 23552f7358SJed Brown { 24552f7358SJed Brown PetscErrorCode ierr; 25552f7358SJed Brown PetscInt offset,dof; 26552f7358SJed Brown 27552f7358SJed Brown PetscFunctionBegin; 28552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 29552f7358SJed Brown ierr = PetscSectionGetOffset(dm->defaultSection,point,&offset);CHKERRQ(ierr); 30552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,point,&dof);CHKERRQ(ierr); 31552f7358SJed Brown if (start) *start = offset; 32552f7358SJed Brown if (end) *end = offset + dof; 33552f7358SJed Brown PetscFunctionReturn(0); 34552f7358SJed Brown } 35552f7358SJed Brown 36552f7358SJed Brown #undef __FUNCT__ 37552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRead" 38552f7358SJed Brown /*@ 39552f7358SJed Brown DMPlexPointLocalRead - return read access to a point in local array 40552f7358SJed Brown 41552f7358SJed Brown Not Collective 42552f7358SJed Brown 43552f7358SJed Brown Input Arguments: 44552f7358SJed Brown + dm - DM defining topological space 45552f7358SJed Brown . point - topological point 46552f7358SJed Brown - array - array to index into 47552f7358SJed Brown 48552f7358SJed Brown Output Arguments: 49552f7358SJed Brown . ptr - address of read reference to point data, type generic so user can place in structure 50552f7358SJed Brown 51552f7358SJed Brown Level: intermediate 52552f7358SJed Brown 53552f7358SJed Brown Note: 54552f7358SJed Brown A common usage when data sizes are known statically: 55552f7358SJed Brown 56552f7358SJed Brown $ const struct { PetscScalar foo,bar,baz; } *ptr; 57552f7358SJed Brown $ DMPlexPointLocalRead(dm,point,array,&ptr); 58552f7358SJed Brown $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 59552f7358SJed Brown 60552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead() 61552f7358SJed Brown @*/ 62552f7358SJed Brown PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 63552f7358SJed Brown { 64552f7358SJed Brown PetscErrorCode ierr; 65552f7358SJed Brown PetscInt start; 66552f7358SJed Brown 67552f7358SJed Brown PetscFunctionBegin; 68552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 69552f7358SJed Brown PetscValidScalarPointer(array,3); 70552f7358SJed Brown PetscValidPointer(ptr,4); 71552f7358SJed Brown ierr = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr); 72552f7358SJed Brown *(const PetscScalar**)ptr = array + start; 73552f7358SJed Brown PetscFunctionReturn(0); 74552f7358SJed Brown } 75552f7358SJed Brown 76552f7358SJed Brown #undef __FUNCT__ 77552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRef" 78552f7358SJed Brown /*@ 79552f7358SJed Brown DMPlexPointLocalRef - return read/write access to a point in local array 80552f7358SJed Brown 81552f7358SJed Brown Not Collective 82552f7358SJed Brown 83552f7358SJed Brown Input Arguments: 84552f7358SJed Brown + dm - DM defining topological space 85552f7358SJed Brown . point - topological point 86552f7358SJed Brown - array - array to index into 87552f7358SJed Brown 88552f7358SJed Brown Output Arguments: 89552f7358SJed Brown . ptr - address of reference to point data, type generic so user can place in structure 90552f7358SJed Brown 91552f7358SJed Brown Level: intermediate 92552f7358SJed Brown 93552f7358SJed Brown Note: 94552f7358SJed Brown A common usage when data sizes are known statically: 95552f7358SJed Brown 96552f7358SJed Brown $ struct { PetscScalar foo,bar,baz; } *ptr; 97552f7358SJed Brown $ DMPlexPointLocalRef(dm,point,array,&ptr); 98552f7358SJed Brown $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 99552f7358SJed Brown 100552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 101552f7358SJed Brown @*/ 102552f7358SJed Brown PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 103552f7358SJed Brown { 104552f7358SJed Brown PetscErrorCode ierr; 105552f7358SJed Brown PetscInt start; 106552f7358SJed Brown 107552f7358SJed Brown PetscFunctionBegin; 108552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 109552f7358SJed Brown PetscValidScalarPointer(array,3); 110552f7358SJed Brown PetscValidPointer(ptr,4); 111552f7358SJed Brown ierr = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr); 112552f7358SJed Brown *(PetscScalar**)ptr = array + start; 113552f7358SJed Brown PetscFunctionReturn(0); 114552f7358SJed Brown } 115552f7358SJed Brown 116552f7358SJed Brown #undef __FUNCT__ 1171ce3176fSMatthew G. Knepley #define __FUNCT__ "DMPlexPointLocalFieldRead" 1181ce3176fSMatthew G. Knepley /*@ 1191ce3176fSMatthew G. Knepley DMPlexPointLocalFieldRead - return read access to a field on a point in local array 1201ce3176fSMatthew G. Knepley 1211ce3176fSMatthew G. Knepley Not Collective 1221ce3176fSMatthew G. Knepley 1231ce3176fSMatthew G. Knepley Input Arguments: 1241ce3176fSMatthew G. Knepley + dm - DM defining topological space 1251ce3176fSMatthew G. Knepley . point - topological point 1261ce3176fSMatthew G. Knepley . field - field number 1271ce3176fSMatthew G. Knepley - array - array to index into 1281ce3176fSMatthew G. Knepley 1291ce3176fSMatthew G. Knepley Output Arguments: 1301ce3176fSMatthew G. Knepley . ptr - address of read reference to point data, type generic so user can place in structure 1311ce3176fSMatthew G. Knepley 1321ce3176fSMatthew G. Knepley Level: intermediate 1331ce3176fSMatthew G. Knepley 1341ce3176fSMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 1351ce3176fSMatthew G. Knepley @*/ 1361ce3176fSMatthew G. Knepley PetscErrorCode DMPlexPointLocalFieldRead(DM dm,PetscInt point,PetscInt field,const PetscScalar *array,const void *ptr) 1371ce3176fSMatthew G. Knepley { 1381ce3176fSMatthew G. Knepley PetscErrorCode ierr; 1391ce3176fSMatthew G. Knepley PetscInt start; 1401ce3176fSMatthew G. Knepley 1411ce3176fSMatthew G. Knepley PetscFunctionBegin; 1421ce3176fSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1431ce3176fSMatthew G. Knepley PetscValidScalarPointer(array,3); 1441ce3176fSMatthew G. Knepley PetscValidPointer(ptr,4); 1451ce3176fSMatthew G. Knepley ierr = DMPlexGetLocalFieldOffset_Private(dm,point,field,&start);CHKERRQ(ierr); 1461ce3176fSMatthew G. Knepley *(const PetscScalar**)ptr = array + start; 1471ce3176fSMatthew G. Knepley PetscFunctionReturn(0); 1481ce3176fSMatthew G. Knepley } 1491ce3176fSMatthew G. Knepley 1501ce3176fSMatthew G. Knepley #undef __FUNCT__ 1514824f456SMatthew G. Knepley #define __FUNCT__ "DMPlexPointLocalFieldRef" 1524824f456SMatthew G. Knepley /*@ 1534824f456SMatthew G. Knepley DMPlexPointLocalFieldRef - return read/write access to a field on a point in local array 1544824f456SMatthew G. Knepley 1554824f456SMatthew G. Knepley Not Collective 1564824f456SMatthew G. Knepley 1574824f456SMatthew G. Knepley Input Arguments: 1584824f456SMatthew G. Knepley + dm - DM defining topological space 1594824f456SMatthew G. Knepley . point - topological point 1604824f456SMatthew G. Knepley . field - field number 1614824f456SMatthew G. Knepley - array - array to index into 1624824f456SMatthew G. Knepley 1634824f456SMatthew G. Knepley Output Arguments: 1644824f456SMatthew G. Knepley . ptr - address of reference to point data, type generic so user can place in structure 1654824f456SMatthew G. Knepley 1664824f456SMatthew G. Knepley Level: intermediate 1674824f456SMatthew G. Knepley 1684824f456SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 1694824f456SMatthew G. Knepley @*/ 1704824f456SMatthew G. Knepley PetscErrorCode DMPlexPointLocalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr) 1714824f456SMatthew G. Knepley { 1724824f456SMatthew G. Knepley PetscErrorCode ierr; 1734824f456SMatthew G. Knepley PetscInt start; 1744824f456SMatthew G. Knepley 1754824f456SMatthew G. Knepley PetscFunctionBegin; 1764824f456SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1774824f456SMatthew G. Knepley PetscValidScalarPointer(array,3); 1784824f456SMatthew G. Knepley PetscValidPointer(ptr,4); 1794824f456SMatthew G. Knepley ierr = DMPlexGetLocalFieldOffset_Private(dm,point,field,&start);CHKERRQ(ierr); 1804824f456SMatthew G. Knepley *(PetscScalar**)ptr = array + start; 1814824f456SMatthew G. Knepley PetscFunctionReturn(0); 1824824f456SMatthew G. Knepley } 1834824f456SMatthew G. Knepley 1844824f456SMatthew G. Knepley #undef __FUNCT__ 185552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointGlobal" 186552f7358SJed Brown /*@ 187552f7358SJed Brown DMPlexGetPointGlobal - get location of point data in global Vec 188552f7358SJed Brown 189552f7358SJed Brown Not Collective 190552f7358SJed Brown 191552f7358SJed Brown Input Arguments: 192552f7358SJed Brown + dm - DM defining the topological space 193552f7358SJed Brown - point - topological point 194552f7358SJed Brown 195552f7358SJed Brown Output Arguments: 196552f7358SJed Brown + start - start of point data; returns -(global_start+1) if point is not owned 197552f7358SJed Brown - end - end of point data; returns -(global_end+1) if point is not owned 198552f7358SJed Brown 199552f7358SJed Brown Level: intermediate 200552f7358SJed Brown 201552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef() 202552f7358SJed Brown @*/ 203552f7358SJed Brown PetscErrorCode DMPlexGetPointGlobal(DM dm,PetscInt point,PetscInt *start,PetscInt *end) 204552f7358SJed Brown { 205552f7358SJed Brown PetscErrorCode ierr; 206552f7358SJed Brown PetscInt offset,dof; 207552f7358SJed Brown 208552f7358SJed Brown PetscFunctionBegin; 209552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 210552f7358SJed Brown ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,&offset);CHKERRQ(ierr); 211552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr); 212552f7358SJed Brown if (start) *start = offset; 213552f7358SJed Brown if (end) *end = offset + dof; 214552f7358SJed Brown PetscFunctionReturn(0); 215552f7358SJed Brown } 216552f7358SJed Brown 217552f7358SJed Brown #undef __FUNCT__ 218552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRead" 219552f7358SJed Brown /*@ 220552f7358SJed Brown DMPlexPointGlobalRead - return read access to a point in global array 221552f7358SJed Brown 222552f7358SJed Brown Not Collective 223552f7358SJed Brown 224552f7358SJed Brown Input Arguments: 225552f7358SJed Brown + dm - DM defining topological space 226552f7358SJed Brown . point - topological point 227552f7358SJed Brown - array - array to index into 228552f7358SJed Brown 229552f7358SJed Brown Output Arguments: 2300298fd71SBarry 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 231552f7358SJed Brown 232552f7358SJed Brown Level: intermediate 233552f7358SJed Brown 234552f7358SJed Brown Note: 235552f7358SJed Brown A common usage when data sizes are known statically: 236552f7358SJed Brown 237552f7358SJed Brown $ const struct { PetscScalar foo,bar,baz; } *ptr; 238552f7358SJed Brown $ DMPlexPointGlobalRead(dm,point,array,&ptr); 239552f7358SJed Brown $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 240552f7358SJed Brown 241552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef() 242552f7358SJed Brown @*/ 243552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 244552f7358SJed Brown { 245*79532bb4SMatthew G. Knepley PetscInt start, end; 246552f7358SJed Brown PetscErrorCode ierr; 247552f7358SJed Brown 248552f7358SJed Brown PetscFunctionBegin; 249552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 250552f7358SJed Brown PetscValidScalarPointer(array, 3); 251552f7358SJed Brown PetscValidPointer(ptr, 4); 252*79532bb4SMatthew G. Knepley ierr = DMPlexGetGlobalOffset_Private(dm, point, &start, &end);CHKERRQ(ierr); 253*79532bb4SMatthew G. Knepley *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 254552f7358SJed Brown PetscFunctionReturn(0); 255552f7358SJed Brown } 256552f7358SJed Brown 257552f7358SJed Brown #undef __FUNCT__ 258552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRef" 259552f7358SJed Brown /*@ 260552f7358SJed Brown DMPlexPointGlobalRef - return read/write access to a point in global array 261552f7358SJed Brown 262552f7358SJed Brown Not Collective 263552f7358SJed Brown 264552f7358SJed Brown Input Arguments: 265552f7358SJed Brown + dm - DM defining topological space 266552f7358SJed Brown . point - topological point 267552f7358SJed Brown - array - array to index into 268552f7358SJed Brown 269552f7358SJed Brown Output Arguments: 2700298fd71SBarry Smith . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 271552f7358SJed Brown 272552f7358SJed Brown Level: intermediate 273552f7358SJed Brown 274552f7358SJed Brown Note: 275552f7358SJed Brown A common usage when data sizes are known statically: 276552f7358SJed Brown 277552f7358SJed Brown $ struct { PetscScalar foo,bar,baz; } *ptr; 278552f7358SJed Brown $ DMPlexPointGlobalRef(dm,point,array,&ptr); 279552f7358SJed Brown $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 280552f7358SJed Brown 281552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead() 282552f7358SJed Brown @*/ 283552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 284552f7358SJed Brown { 285*79532bb4SMatthew G. Knepley PetscInt start, end; 286552f7358SJed Brown PetscErrorCode ierr; 287552f7358SJed Brown 288552f7358SJed Brown PetscFunctionBegin; 289552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 290552f7358SJed Brown PetscValidScalarPointer(array, 3); 291552f7358SJed Brown PetscValidPointer(ptr, 4); 292*79532bb4SMatthew G. Knepley ierr = DMPlexGetGlobalOffset_Private(dm, point, &start, &end);CHKERRQ(ierr); 293*79532bb4SMatthew G. Knepley *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 294552f7358SJed Brown PetscFunctionReturn(0); 295552f7358SJed Brown } 29633879625SMatthew G. Knepley 29733879625SMatthew G. Knepley #undef __FUNCT__ 29833879625SMatthew G. Knepley #define __FUNCT__ "DMPlexPointGlobalFieldRead" 29933879625SMatthew G. Knepley /*@ 30033879625SMatthew G. Knepley DMPlexPointGlobalFieldRead - return read access to a field on a point in global array 30133879625SMatthew G. Knepley 30233879625SMatthew G. Knepley Not Collective 30333879625SMatthew G. Knepley 30433879625SMatthew G. Knepley Input Arguments: 30533879625SMatthew G. Knepley + dm - DM defining topological space 30633879625SMatthew G. Knepley . point - topological point 30733879625SMatthew G. Knepley . field - field number 30833879625SMatthew G. Knepley - array - array to index into 30933879625SMatthew G. Knepley 31033879625SMatthew G. Knepley Output Arguments: 31133879625SMatthew 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 31233879625SMatthew G. Knepley 31333879625SMatthew G. Knepley Level: intermediate 31433879625SMatthew G. Knepley 31533879625SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef() 31633879625SMatthew G. Knepley @*/ 31733879625SMatthew G. Knepley PetscErrorCode DMPlexPointGlobalFieldRead(DM dm,PetscInt point,PetscInt field,const PetscScalar *array,const void *ptr) 31833879625SMatthew G. Knepley { 319*79532bb4SMatthew G. Knepley PetscInt start, end; 32033879625SMatthew G. Knepley PetscErrorCode ierr; 32133879625SMatthew G. Knepley 32233879625SMatthew G. Knepley PetscFunctionBegin; 32333879625SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 32433879625SMatthew G. Knepley PetscValidScalarPointer(array, 3); 32533879625SMatthew G. Knepley PetscValidPointer(ptr, 4); 326*79532bb4SMatthew G. Knepley ierr = DMPlexGetGlobalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 327*79532bb4SMatthew G. Knepley *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 32833879625SMatthew G. Knepley PetscFunctionReturn(0); 32933879625SMatthew G. Knepley } 33033879625SMatthew G. Knepley 33133879625SMatthew G. Knepley #undef __FUNCT__ 33233879625SMatthew G. Knepley #define __FUNCT__ "DMPlexPointGlobalFieldRef" 33333879625SMatthew G. Knepley /*@ 33433879625SMatthew G. Knepley DMPlexPointGlobalFieldRef - return read/write access to a field on a point in global array 33533879625SMatthew G. Knepley 33633879625SMatthew G. Knepley Not Collective 33733879625SMatthew G. Knepley 33833879625SMatthew G. Knepley Input Arguments: 33933879625SMatthew G. Knepley + dm - DM defining topological space 34033879625SMatthew G. Knepley . point - topological point 34133879625SMatthew G. Knepley . field - field number 34233879625SMatthew G. Knepley - array - array to index into 34333879625SMatthew G. Knepley 34433879625SMatthew G. Knepley Output Arguments: 34533879625SMatthew 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 34633879625SMatthew G. Knepley 34733879625SMatthew G. Knepley Level: intermediate 34833879625SMatthew G. Knepley 34933879625SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead() 35033879625SMatthew G. Knepley @*/ 35133879625SMatthew G. Knepley PetscErrorCode DMPlexPointGlobalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr) 35233879625SMatthew G. Knepley { 353*79532bb4SMatthew G. Knepley PetscInt start, end; 35433879625SMatthew G. Knepley PetscErrorCode ierr; 35533879625SMatthew G. Knepley 35633879625SMatthew G. Knepley PetscFunctionBegin; 35733879625SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 35833879625SMatthew G. Knepley PetscValidScalarPointer(array, 3); 35933879625SMatthew G. Knepley PetscValidPointer(ptr, 4); 360*79532bb4SMatthew G. Knepley ierr = DMPlexGetGlobalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 361*79532bb4SMatthew G. Knepley *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 36233879625SMatthew G. Knepley PetscFunctionReturn(0); 36333879625SMatthew G. Knepley } 364