1*34541f0dSBarry Smith #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2ad7c26b0SJed Brown #include <petsc-private/isimpl.h> /* for inline access to atlasOff */ 3552f7358SJed Brown 4552f7358SJed Brown #undef __FUNCT__ 5552f7358SJed Brown #define __FUNCT__ "DMPlexGetLocalOffset_Private" 6552f7358SJed Brown PETSC_STATIC_INLINE PetscErrorCode DMPlexGetLocalOffset_Private(DM dm,PetscInt point,PetscInt *offset) 7552f7358SJed Brown { 8552f7358SJed Brown PetscFunctionBegin; 9552f7358SJed Brown #if defined(PETSC_USE_DEBUG) 10552f7358SJed Brown { 11552f7358SJed Brown PetscErrorCode ierr; 12552f7358SJed Brown ierr = PetscSectionGetOffset(dm->defaultSection,point,offset);CHKERRQ(ierr); 13552f7358SJed Brown } 14552f7358SJed Brown #else 15552f7358SJed Brown { 16552f7358SJed Brown PetscSection s = dm->defaultSection; 17552f7358SJed Brown *offset = s->atlasOff[point - s->atlasLayout.pStart]; 18552f7358SJed Brown } 19552f7358SJed Brown #endif 20552f7358SJed Brown PetscFunctionReturn(0); 21552f7358SJed Brown } 22552f7358SJed Brown 23552f7358SJed Brown #undef __FUNCT__ 24552f7358SJed Brown #define __FUNCT__ "DMPlexGetGlobalOffset_Private" 25552f7358SJed Brown PETSC_STATIC_INLINE PetscErrorCode DMPlexGetGlobalOffset_Private(DM dm,PetscInt point,PetscInt *offset) 26552f7358SJed Brown { 27552f7358SJed Brown PetscFunctionBegin; 28552f7358SJed Brown #if defined(PETSC_USE_DEBUG) 29552f7358SJed Brown { 30552f7358SJed Brown PetscErrorCode ierr; 31552f7358SJed Brown PetscInt dof,cdof; 32552f7358SJed Brown ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,offset);CHKERRQ(ierr); 33552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr); 34552f7358SJed Brown ierr = PetscSectionGetConstraintDof(dm->defaultGlobalSection,point,&cdof);CHKERRQ(ierr); 35552f7358SJed Brown if (dof-cdof <= 0) *offset = -1; /* Indicates no data */ 36552f7358SJed Brown } 37552f7358SJed Brown #else 38552f7358SJed Brown { 39552f7358SJed Brown PetscSection s = dm->defaultGlobalSection; 40552f7358SJed Brown PetscInt dof,cdof; 41552f7358SJed Brown *offset = s->atlasOff[point - s->atlasLayout.pStart]; 42552f7358SJed Brown dof = s->atlasDof[point - s->atlasLayout.pStart]; 43552f7358SJed Brown cdof = s->bc ? s->bc->atlasDof[point - s->bc->atlasLayout.pStart] : 0; 44552f7358SJed Brown if (dof-cdof <= 0) *offset = -1; 45552f7358SJed Brown } 46552f7358SJed Brown #endif 47552f7358SJed Brown PetscFunctionReturn(0); 48552f7358SJed Brown } 49552f7358SJed Brown 50552f7358SJed Brown #undef __FUNCT__ 51552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointLocal" 52552f7358SJed Brown /*@ 53552f7358SJed Brown DMPlexGetPointLocal - get location of point data in local Vec 54552f7358SJed Brown 55552f7358SJed Brown Not Collective 56552f7358SJed Brown 57552f7358SJed Brown Input Arguments: 58552f7358SJed Brown + dm - DM defining the topological space 59552f7358SJed Brown - point - topological point 60552f7358SJed Brown 61552f7358SJed Brown Output Arguments: 62552f7358SJed Brown + start - start of point data 63552f7358SJed Brown - end - end of point data 64552f7358SJed Brown 65552f7358SJed Brown Level: intermediate 66552f7358SJed Brown 67552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef() 68552f7358SJed Brown @*/ 69552f7358SJed Brown PetscErrorCode DMPlexGetPointLocal(DM dm,PetscInt point,PetscInt *start,PetscInt *end) 70552f7358SJed Brown { 71552f7358SJed Brown PetscErrorCode ierr; 72552f7358SJed Brown PetscInt offset,dof; 73552f7358SJed Brown 74552f7358SJed Brown PetscFunctionBegin; 75552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 76552f7358SJed Brown ierr = PetscSectionGetOffset(dm->defaultSection,point,&offset);CHKERRQ(ierr); 77552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,point,&dof);CHKERRQ(ierr); 78552f7358SJed Brown if (start) *start = offset; 79552f7358SJed Brown if (end) *end = offset + dof; 80552f7358SJed Brown PetscFunctionReturn(0); 81552f7358SJed Brown } 82552f7358SJed Brown 83552f7358SJed Brown #undef __FUNCT__ 84552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRead" 85552f7358SJed Brown /*@ 86552f7358SJed Brown DMPlexPointLocalRead - return read access to a point in local array 87552f7358SJed Brown 88552f7358SJed Brown Not Collective 89552f7358SJed Brown 90552f7358SJed Brown Input Arguments: 91552f7358SJed Brown + dm - DM defining topological space 92552f7358SJed Brown . point - topological point 93552f7358SJed Brown - array - array to index into 94552f7358SJed Brown 95552f7358SJed Brown Output Arguments: 96552f7358SJed Brown . ptr - address of read reference to point data, type generic so user can place in structure 97552f7358SJed Brown 98552f7358SJed Brown Level: intermediate 99552f7358SJed Brown 100552f7358SJed Brown Note: 101552f7358SJed Brown A common usage when data sizes are known statically: 102552f7358SJed Brown 103552f7358SJed Brown $ const struct { PetscScalar foo,bar,baz; } *ptr; 104552f7358SJed Brown $ DMPlexPointLocalRead(dm,point,array,&ptr); 105552f7358SJed Brown $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 106552f7358SJed Brown 107552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead() 108552f7358SJed Brown @*/ 109552f7358SJed Brown PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 110552f7358SJed Brown { 111552f7358SJed Brown PetscErrorCode ierr; 112552f7358SJed Brown PetscInt start; 113552f7358SJed Brown 114552f7358SJed Brown PetscFunctionBegin; 115552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 116552f7358SJed Brown PetscValidScalarPointer(array,3); 117552f7358SJed Brown PetscValidPointer(ptr,4); 118552f7358SJed Brown ierr = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr); 119552f7358SJed Brown *(const PetscScalar**)ptr = array + start; 120552f7358SJed Brown PetscFunctionReturn(0); 121552f7358SJed Brown } 122552f7358SJed Brown 123552f7358SJed Brown #undef __FUNCT__ 124552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRef" 125552f7358SJed Brown /*@ 126552f7358SJed Brown DMPlexPointLocalRef - return read/write access to a point in local array 127552f7358SJed Brown 128552f7358SJed Brown Not Collective 129552f7358SJed Brown 130552f7358SJed Brown Input Arguments: 131552f7358SJed Brown + dm - DM defining topological space 132552f7358SJed Brown . point - topological point 133552f7358SJed Brown - array - array to index into 134552f7358SJed Brown 135552f7358SJed Brown Output Arguments: 136552f7358SJed Brown . ptr - address of reference to point data, type generic so user can place in structure 137552f7358SJed Brown 138552f7358SJed Brown Level: intermediate 139552f7358SJed Brown 140552f7358SJed Brown Note: 141552f7358SJed Brown A common usage when data sizes are known statically: 142552f7358SJed Brown 143552f7358SJed Brown $ struct { PetscScalar foo,bar,baz; } *ptr; 144552f7358SJed Brown $ DMPlexPointLocalRef(dm,point,array,&ptr); 145552f7358SJed Brown $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 146552f7358SJed Brown 147552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 148552f7358SJed Brown @*/ 149552f7358SJed Brown PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 150552f7358SJed Brown { 151552f7358SJed Brown PetscErrorCode ierr; 152552f7358SJed Brown PetscInt start; 153552f7358SJed Brown 154552f7358SJed Brown PetscFunctionBegin; 155552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 156552f7358SJed Brown PetscValidScalarPointer(array,3); 157552f7358SJed Brown PetscValidPointer(ptr,4); 158552f7358SJed Brown ierr = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr); 159552f7358SJed Brown *(PetscScalar**)ptr = array + start; 160552f7358SJed Brown PetscFunctionReturn(0); 161552f7358SJed Brown } 162552f7358SJed Brown 163552f7358SJed Brown #undef __FUNCT__ 164552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointGlobal" 165552f7358SJed Brown /*@ 166552f7358SJed Brown DMPlexGetPointGlobal - get location of point data in global Vec 167552f7358SJed Brown 168552f7358SJed Brown Not Collective 169552f7358SJed Brown 170552f7358SJed Brown Input Arguments: 171552f7358SJed Brown + dm - DM defining the topological space 172552f7358SJed Brown - point - topological point 173552f7358SJed Brown 174552f7358SJed Brown Output Arguments: 175552f7358SJed Brown + start - start of point data; returns -(global_start+1) if point is not owned 176552f7358SJed Brown - end - end of point data; returns -(global_end+1) if point is not owned 177552f7358SJed Brown 178552f7358SJed Brown Level: intermediate 179552f7358SJed Brown 180552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef() 181552f7358SJed Brown @*/ 182552f7358SJed Brown PetscErrorCode DMPlexGetPointGlobal(DM dm,PetscInt point,PetscInt *start,PetscInt *end) 183552f7358SJed Brown { 184552f7358SJed Brown PetscErrorCode ierr; 185552f7358SJed Brown PetscInt offset,dof; 186552f7358SJed Brown 187552f7358SJed Brown PetscFunctionBegin; 188552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 189552f7358SJed Brown ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,&offset);CHKERRQ(ierr); 190552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr); 191552f7358SJed Brown if (start) *start = offset; 192552f7358SJed Brown if (end) *end = offset + dof; 193552f7358SJed Brown PetscFunctionReturn(0); 194552f7358SJed Brown } 195552f7358SJed Brown 196552f7358SJed Brown #undef __FUNCT__ 197552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRead" 198552f7358SJed Brown /*@ 199552f7358SJed Brown DMPlexPointGlobalRead - return read access to a point in global array 200552f7358SJed Brown 201552f7358SJed Brown Not Collective 202552f7358SJed Brown 203552f7358SJed Brown Input Arguments: 204552f7358SJed Brown + dm - DM defining topological space 205552f7358SJed Brown . point - topological point 206552f7358SJed Brown - array - array to index into 207552f7358SJed Brown 208552f7358SJed Brown Output Arguments: 2090298fd71SBarry 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 210552f7358SJed Brown 211552f7358SJed Brown Level: intermediate 212552f7358SJed Brown 213552f7358SJed Brown Note: 214552f7358SJed Brown A common usage when data sizes are known statically: 215552f7358SJed Brown 216552f7358SJed Brown $ const struct { PetscScalar foo,bar,baz; } *ptr; 217552f7358SJed Brown $ DMPlexPointGlobalRead(dm,point,array,&ptr); 218552f7358SJed Brown $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 219552f7358SJed Brown 220552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef() 221552f7358SJed Brown @*/ 222552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 223552f7358SJed Brown { 224552f7358SJed Brown PetscErrorCode ierr; 225552f7358SJed Brown PetscInt start; 226552f7358SJed Brown 227552f7358SJed Brown PetscFunctionBegin; 228552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 229552f7358SJed Brown PetscValidScalarPointer(array,3); 230552f7358SJed Brown PetscValidPointer(ptr,4); 231552f7358SJed Brown ierr = DMPlexGetGlobalOffset_Private(dm,point,&start);CHKERRQ(ierr); 2320298fd71SBarry Smith *(const PetscScalar**)ptr = (start >= 0) ? array + start - dm->map->rstart : NULL; 233552f7358SJed Brown PetscFunctionReturn(0); 234552f7358SJed Brown } 235552f7358SJed Brown 236552f7358SJed Brown #undef __FUNCT__ 237552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRef" 238552f7358SJed Brown /*@ 239552f7358SJed Brown DMPlexPointGlobalRef - return read/write access to a point in global array 240552f7358SJed Brown 241552f7358SJed Brown Not Collective 242552f7358SJed Brown 243552f7358SJed Brown Input Arguments: 244552f7358SJed Brown + dm - DM defining topological space 245552f7358SJed Brown . point - topological point 246552f7358SJed Brown - array - array to index into 247552f7358SJed Brown 248552f7358SJed Brown Output Arguments: 2490298fd71SBarry Smith . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 250552f7358SJed Brown 251552f7358SJed Brown Level: intermediate 252552f7358SJed Brown 253552f7358SJed Brown Note: 254552f7358SJed Brown A common usage when data sizes are known statically: 255552f7358SJed Brown 256552f7358SJed Brown $ struct { PetscScalar foo,bar,baz; } *ptr; 257552f7358SJed Brown $ DMPlexPointGlobalRef(dm,point,array,&ptr); 258552f7358SJed Brown $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 259552f7358SJed Brown 260552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead() 261552f7358SJed Brown @*/ 262552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 263552f7358SJed Brown { 264552f7358SJed Brown PetscErrorCode ierr; 265552f7358SJed Brown PetscInt start; 266552f7358SJed Brown 267552f7358SJed Brown PetscFunctionBegin; 268552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 269552f7358SJed Brown PetscValidScalarPointer(array,3); 270552f7358SJed Brown PetscValidPointer(ptr,4); 271552f7358SJed Brown ierr = DMPlexGetGlobalOffset_Private(dm,point,&start);CHKERRQ(ierr); 2720298fd71SBarry Smith *(PetscScalar**)ptr = (start >= 0) ? array + start - dm->map->rstart : NULL; 273552f7358SJed Brown PetscFunctionReturn(0); 274552f7358SJed Brown } 275