134541f0dSBarry 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; 174c0d72c2SMatthew G. Knepley *offset = s->atlasOff[point - s->pStart]; 18552f7358SJed Brown } 19552f7358SJed Brown #endif 20552f7358SJed Brown PetscFunctionReturn(0); 21552f7358SJed Brown } 22552f7358SJed Brown 23552f7358SJed Brown #undef __FUNCT__ 24*4824f456SMatthew G. Knepley #define __FUNCT__ "DMPlexGetLocalFieldOffset_Private" 25*4824f456SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMPlexGetLocalFieldOffset_Private(DM dm,PetscInt point,PetscInt field,PetscInt *offset) 26*4824f456SMatthew G. Knepley { 27*4824f456SMatthew G. Knepley PetscFunctionBegin; 28*4824f456SMatthew G. Knepley #if defined(PETSC_USE_DEBUG) 29*4824f456SMatthew G. Knepley { 30*4824f456SMatthew G. Knepley PetscErrorCode ierr; 31*4824f456SMatthew G. Knepley ierr = PetscSectionGetFieldOffset(dm->defaultSection,point,field,offset);CHKERRQ(ierr); 32*4824f456SMatthew G. Knepley } 33*4824f456SMatthew G. Knepley #else 34*4824f456SMatthew G. Knepley { 35*4824f456SMatthew G. Knepley PetscSection s = dm->defaultSection->field[field]; 36*4824f456SMatthew G. Knepley *offset = s->atlasOff[point - s->pStart]; 37*4824f456SMatthew G. Knepley } 38*4824f456SMatthew G. Knepley #endif 39*4824f456SMatthew G. Knepley PetscFunctionReturn(0); 40*4824f456SMatthew G. Knepley } 41*4824f456SMatthew G. Knepley 42*4824f456SMatthew G. Knepley #undef __FUNCT__ 43552f7358SJed Brown #define __FUNCT__ "DMPlexGetGlobalOffset_Private" 44552f7358SJed Brown PETSC_STATIC_INLINE PetscErrorCode DMPlexGetGlobalOffset_Private(DM dm,PetscInt point,PetscInt *offset) 45552f7358SJed Brown { 46552f7358SJed Brown PetscFunctionBegin; 47552f7358SJed Brown #if defined(PETSC_USE_DEBUG) 48552f7358SJed Brown { 49552f7358SJed Brown PetscErrorCode ierr; 50552f7358SJed Brown PetscInt dof,cdof; 51552f7358SJed Brown ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,offset);CHKERRQ(ierr); 52552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr); 53552f7358SJed Brown ierr = PetscSectionGetConstraintDof(dm->defaultGlobalSection,point,&cdof);CHKERRQ(ierr); 54552f7358SJed Brown if (dof-cdof <= 0) *offset = -1; /* Indicates no data */ 55552f7358SJed Brown } 56552f7358SJed Brown #else 57552f7358SJed Brown { 58552f7358SJed Brown PetscSection s = dm->defaultGlobalSection; 59552f7358SJed Brown PetscInt dof,cdof; 604c0d72c2SMatthew G. Knepley *offset = s->atlasOff[point - s->pStart]; 614c0d72c2SMatthew G. Knepley dof = s->atlasDof[point - s->pStart]; 624c0d72c2SMatthew G. Knepley cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0; 63552f7358SJed Brown if (dof-cdof <= 0) *offset = -1; 64552f7358SJed Brown } 65552f7358SJed Brown #endif 66552f7358SJed Brown PetscFunctionReturn(0); 67552f7358SJed Brown } 68552f7358SJed Brown 69552f7358SJed Brown #undef __FUNCT__ 70552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointLocal" 71552f7358SJed Brown /*@ 72552f7358SJed Brown DMPlexGetPointLocal - get location of point data in local Vec 73552f7358SJed Brown 74552f7358SJed Brown Not Collective 75552f7358SJed Brown 76552f7358SJed Brown Input Arguments: 77552f7358SJed Brown + dm - DM defining the topological space 78552f7358SJed Brown - point - topological point 79552f7358SJed Brown 80552f7358SJed Brown Output Arguments: 81552f7358SJed Brown + start - start of point data 82552f7358SJed Brown - end - end of point data 83552f7358SJed Brown 84552f7358SJed Brown Level: intermediate 85552f7358SJed Brown 86552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef() 87552f7358SJed Brown @*/ 88552f7358SJed Brown PetscErrorCode DMPlexGetPointLocal(DM dm,PetscInt point,PetscInt *start,PetscInt *end) 89552f7358SJed Brown { 90552f7358SJed Brown PetscErrorCode ierr; 91552f7358SJed Brown PetscInt offset,dof; 92552f7358SJed Brown 93552f7358SJed Brown PetscFunctionBegin; 94552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 95552f7358SJed Brown ierr = PetscSectionGetOffset(dm->defaultSection,point,&offset);CHKERRQ(ierr); 96552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultSection,point,&dof);CHKERRQ(ierr); 97552f7358SJed Brown if (start) *start = offset; 98552f7358SJed Brown if (end) *end = offset + dof; 99552f7358SJed Brown PetscFunctionReturn(0); 100552f7358SJed Brown } 101552f7358SJed Brown 102552f7358SJed Brown #undef __FUNCT__ 103552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRead" 104552f7358SJed Brown /*@ 105552f7358SJed Brown DMPlexPointLocalRead - return read access to a point in local array 106552f7358SJed Brown 107552f7358SJed Brown Not Collective 108552f7358SJed Brown 109552f7358SJed Brown Input Arguments: 110552f7358SJed Brown + dm - DM defining topological space 111552f7358SJed Brown . point - topological point 112552f7358SJed Brown - array - array to index into 113552f7358SJed Brown 114552f7358SJed Brown Output Arguments: 115552f7358SJed Brown . ptr - address of read reference to point data, type generic so user can place in structure 116552f7358SJed Brown 117552f7358SJed Brown Level: intermediate 118552f7358SJed Brown 119552f7358SJed Brown Note: 120552f7358SJed Brown A common usage when data sizes are known statically: 121552f7358SJed Brown 122552f7358SJed Brown $ const struct { PetscScalar foo,bar,baz; } *ptr; 123552f7358SJed Brown $ DMPlexPointLocalRead(dm,point,array,&ptr); 124552f7358SJed Brown $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 125552f7358SJed Brown 126552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead() 127552f7358SJed Brown @*/ 128552f7358SJed Brown PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 129552f7358SJed Brown { 130552f7358SJed Brown PetscErrorCode ierr; 131552f7358SJed Brown PetscInt start; 132552f7358SJed Brown 133552f7358SJed Brown PetscFunctionBegin; 134552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 135552f7358SJed Brown PetscValidScalarPointer(array,3); 136552f7358SJed Brown PetscValidPointer(ptr,4); 137552f7358SJed Brown ierr = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr); 138552f7358SJed Brown *(const PetscScalar**)ptr = array + start; 139552f7358SJed Brown PetscFunctionReturn(0); 140552f7358SJed Brown } 141552f7358SJed Brown 142552f7358SJed Brown #undef __FUNCT__ 143552f7358SJed Brown #define __FUNCT__ "DMPlexPointLocalRef" 144552f7358SJed Brown /*@ 145552f7358SJed Brown DMPlexPointLocalRef - return read/write access to a point in local array 146552f7358SJed Brown 147552f7358SJed Brown Not Collective 148552f7358SJed Brown 149552f7358SJed Brown Input Arguments: 150552f7358SJed Brown + dm - DM defining topological space 151552f7358SJed Brown . point - topological point 152552f7358SJed Brown - array - array to index into 153552f7358SJed Brown 154552f7358SJed Brown Output Arguments: 155552f7358SJed Brown . ptr - address of reference to point data, type generic so user can place in structure 156552f7358SJed Brown 157552f7358SJed Brown Level: intermediate 158552f7358SJed Brown 159552f7358SJed Brown Note: 160552f7358SJed Brown A common usage when data sizes are known statically: 161552f7358SJed Brown 162552f7358SJed Brown $ struct { PetscScalar foo,bar,baz; } *ptr; 163552f7358SJed Brown $ DMPlexPointLocalRef(dm,point,array,&ptr); 164552f7358SJed Brown $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 165552f7358SJed Brown 166552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 167552f7358SJed Brown @*/ 168552f7358SJed Brown PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 169552f7358SJed Brown { 170552f7358SJed Brown PetscErrorCode ierr; 171552f7358SJed Brown PetscInt start; 172552f7358SJed Brown 173552f7358SJed Brown PetscFunctionBegin; 174552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 175552f7358SJed Brown PetscValidScalarPointer(array,3); 176552f7358SJed Brown PetscValidPointer(ptr,4); 177552f7358SJed Brown ierr = DMPlexGetLocalOffset_Private(dm,point,&start);CHKERRQ(ierr); 178552f7358SJed Brown *(PetscScalar**)ptr = array + start; 179552f7358SJed Brown PetscFunctionReturn(0); 180552f7358SJed Brown } 181552f7358SJed Brown 182552f7358SJed Brown #undef __FUNCT__ 183*4824f456SMatthew G. Knepley #define __FUNCT__ "DMPlexPointLocalFieldRef" 184*4824f456SMatthew G. Knepley /*@ 185*4824f456SMatthew G. Knepley DMPlexPointLocalFieldRef - return read/write access to a field on a point in local array 186*4824f456SMatthew G. Knepley 187*4824f456SMatthew G. Knepley Not Collective 188*4824f456SMatthew G. Knepley 189*4824f456SMatthew G. Knepley Input Arguments: 190*4824f456SMatthew G. Knepley + dm - DM defining topological space 191*4824f456SMatthew G. Knepley . point - topological point 192*4824f456SMatthew G. Knepley . field - field number 193*4824f456SMatthew G. Knepley - array - array to index into 194*4824f456SMatthew G. Knepley 195*4824f456SMatthew G. Knepley Output Arguments: 196*4824f456SMatthew G. Knepley . ptr - address of reference to point data, type generic so user can place in structure 197*4824f456SMatthew G. Knepley 198*4824f456SMatthew G. Knepley Level: intermediate 199*4824f456SMatthew G. Knepley 200*4824f456SMatthew G. Knepley .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 201*4824f456SMatthew G. Knepley @*/ 202*4824f456SMatthew G. Knepley PetscErrorCode DMPlexPointLocalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr) 203*4824f456SMatthew G. Knepley { 204*4824f456SMatthew G. Knepley PetscErrorCode ierr; 205*4824f456SMatthew G. Knepley PetscInt start; 206*4824f456SMatthew G. Knepley 207*4824f456SMatthew G. Knepley PetscFunctionBegin; 208*4824f456SMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 209*4824f456SMatthew G. Knepley PetscValidScalarPointer(array,3); 210*4824f456SMatthew G. Knepley PetscValidPointer(ptr,4); 211*4824f456SMatthew G. Knepley ierr = DMPlexGetLocalFieldOffset_Private(dm,point,field,&start);CHKERRQ(ierr); 212*4824f456SMatthew G. Knepley *(PetscScalar**)ptr = array + start; 213*4824f456SMatthew G. Knepley PetscFunctionReturn(0); 214*4824f456SMatthew G. Knepley } 215*4824f456SMatthew G. Knepley 216*4824f456SMatthew G. Knepley #undef __FUNCT__ 217552f7358SJed Brown #define __FUNCT__ "DMPlexGetPointGlobal" 218552f7358SJed Brown /*@ 219552f7358SJed Brown DMPlexGetPointGlobal - get location of point data in global Vec 220552f7358SJed Brown 221552f7358SJed Brown Not Collective 222552f7358SJed Brown 223552f7358SJed Brown Input Arguments: 224552f7358SJed Brown + dm - DM defining the topological space 225552f7358SJed Brown - point - topological point 226552f7358SJed Brown 227552f7358SJed Brown Output Arguments: 228552f7358SJed Brown + start - start of point data; returns -(global_start+1) if point is not owned 229552f7358SJed Brown - end - end of point data; returns -(global_end+1) if point is not owned 230552f7358SJed Brown 231552f7358SJed Brown Level: intermediate 232552f7358SJed Brown 233552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef() 234552f7358SJed Brown @*/ 235552f7358SJed Brown PetscErrorCode DMPlexGetPointGlobal(DM dm,PetscInt point,PetscInt *start,PetscInt *end) 236552f7358SJed Brown { 237552f7358SJed Brown PetscErrorCode ierr; 238552f7358SJed Brown PetscInt offset,dof; 239552f7358SJed Brown 240552f7358SJed Brown PetscFunctionBegin; 241552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 242552f7358SJed Brown ierr = PetscSectionGetOffset(dm->defaultGlobalSection,point,&offset);CHKERRQ(ierr); 243552f7358SJed Brown ierr = PetscSectionGetDof(dm->defaultGlobalSection,point,&dof);CHKERRQ(ierr); 244552f7358SJed Brown if (start) *start = offset; 245552f7358SJed Brown if (end) *end = offset + dof; 246552f7358SJed Brown PetscFunctionReturn(0); 247552f7358SJed Brown } 248552f7358SJed Brown 249552f7358SJed Brown #undef __FUNCT__ 250552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRead" 251552f7358SJed Brown /*@ 252552f7358SJed Brown DMPlexPointGlobalRead - return read access to a point in global array 253552f7358SJed Brown 254552f7358SJed Brown Not Collective 255552f7358SJed Brown 256552f7358SJed Brown Input Arguments: 257552f7358SJed Brown + dm - DM defining topological space 258552f7358SJed Brown . point - topological point 259552f7358SJed Brown - array - array to index into 260552f7358SJed Brown 261552f7358SJed Brown Output Arguments: 2620298fd71SBarry 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 263552f7358SJed Brown 264552f7358SJed Brown Level: intermediate 265552f7358SJed Brown 266552f7358SJed Brown Note: 267552f7358SJed Brown A common usage when data sizes are known statically: 268552f7358SJed Brown 269552f7358SJed Brown $ const struct { PetscScalar foo,bar,baz; } *ptr; 270552f7358SJed Brown $ DMPlexPointGlobalRead(dm,point,array,&ptr); 271552f7358SJed Brown $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 272552f7358SJed Brown 273552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef() 274552f7358SJed Brown @*/ 275552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 276552f7358SJed Brown { 277552f7358SJed Brown PetscErrorCode ierr; 278552f7358SJed Brown PetscInt start; 279552f7358SJed Brown 280552f7358SJed Brown PetscFunctionBegin; 281552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 282552f7358SJed Brown PetscValidScalarPointer(array,3); 283552f7358SJed Brown PetscValidPointer(ptr,4); 284552f7358SJed Brown ierr = DMPlexGetGlobalOffset_Private(dm,point,&start);CHKERRQ(ierr); 2850298fd71SBarry Smith *(const PetscScalar**)ptr = (start >= 0) ? array + start - dm->map->rstart : NULL; 286552f7358SJed Brown PetscFunctionReturn(0); 287552f7358SJed Brown } 288552f7358SJed Brown 289552f7358SJed Brown #undef __FUNCT__ 290552f7358SJed Brown #define __FUNCT__ "DMPlexPointGlobalRef" 291552f7358SJed Brown /*@ 292552f7358SJed Brown DMPlexPointGlobalRef - return read/write access to a point in global array 293552f7358SJed Brown 294552f7358SJed Brown Not Collective 295552f7358SJed Brown 296552f7358SJed Brown Input Arguments: 297552f7358SJed Brown + dm - DM defining topological space 298552f7358SJed Brown . point - topological point 299552f7358SJed Brown - array - array to index into 300552f7358SJed Brown 301552f7358SJed Brown Output Arguments: 3020298fd71SBarry Smith . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 303552f7358SJed Brown 304552f7358SJed Brown Level: intermediate 305552f7358SJed Brown 306552f7358SJed Brown Note: 307552f7358SJed Brown A common usage when data sizes are known statically: 308552f7358SJed Brown 309552f7358SJed Brown $ struct { PetscScalar foo,bar,baz; } *ptr; 310552f7358SJed Brown $ DMPlexPointGlobalRef(dm,point,array,&ptr); 311552f7358SJed Brown $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 312552f7358SJed Brown 313552f7358SJed Brown .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead() 314552f7358SJed Brown @*/ 315552f7358SJed Brown PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 316552f7358SJed Brown { 317552f7358SJed Brown PetscErrorCode ierr; 318552f7358SJed Brown PetscInt start; 319552f7358SJed Brown 320552f7358SJed Brown PetscFunctionBegin; 321552f7358SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 322552f7358SJed Brown PetscValidScalarPointer(array,3); 323552f7358SJed Brown PetscValidPointer(ptr,4); 324552f7358SJed Brown ierr = DMPlexGetGlobalOffset_Private(dm,point,&start);CHKERRQ(ierr); 3250298fd71SBarry Smith *(PetscScalar**)ptr = (start >= 0) ? array + start - dm->map->rstart : NULL; 326552f7358SJed Brown PetscFunctionReturn(0); 327552f7358SJed Brown } 328