1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "DMPlexGetPointLocal" 5 /*@ 6 DMPlexGetPointLocal - get location of point data in local Vec 7 8 Not Collective 9 10 Input Arguments: 11 + dm - DM defining the topological space 12 - point - topological point 13 14 Output Arguments: 15 + start - start of point data 16 - end - end of point data 17 18 Note: This is a half open interval [start, end) 19 20 Level: intermediate 21 22 .seealso: DMPlexGetPointLocalField(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef() 23 @*/ 24 PetscErrorCode DMPlexGetPointLocal(DM dm, PetscInt point, PetscInt *start, PetscInt *end) 25 { 26 PetscInt s, e; 27 PetscErrorCode ierr; 28 29 PetscFunctionBegin; 30 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 31 if (start) PetscValidPointer(start, 3); 32 if (end) PetscValidPointer(end, 4); 33 ierr = DMPlexGetLocalOffset_Private(dm, point, &s, &e);CHKERRQ(ierr); 34 if (start) *start = s; 35 if (end) *end = e; 36 PetscFunctionReturn(0); 37 } 38 39 #undef __FUNCT__ 40 #define __FUNCT__ "DMPlexPointLocalRead" 41 /*@ 42 DMPlexPointLocalRead - return read access to a point in local array 43 44 Not Collective 45 46 Input Arguments: 47 + dm - DM defining topological space 48 . point - topological point 49 - array - array to index into 50 51 Output Arguments: 52 . ptr - address of read reference to point data, type generic so user can place in structure 53 54 Level: intermediate 55 56 Note: 57 A common usage when data sizes are known statically: 58 59 $ const struct { PetscScalar foo,bar,baz; } *ptr; 60 $ DMPlexPointLocalRead(dm,point,array,&ptr); 61 $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 62 63 .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRead() 64 @*/ 65 PetscErrorCode DMPlexPointLocalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 66 { 67 PetscErrorCode ierr; 68 PetscInt start, end; 69 70 PetscFunctionBegin; 71 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 72 PetscValidScalarPointer(array,3); 73 PetscValidPointer(ptr,4); 74 ierr = DMPlexGetLocalOffset_Private(dm,point,&start,&end);CHKERRQ(ierr); 75 *(const PetscScalar**)ptr = (start < end) ? array + start : NULL; 76 PetscFunctionReturn(0); 77 } 78 79 #undef __FUNCT__ 80 #define __FUNCT__ "DMPlexPointLocalRef" 81 /*@ 82 DMPlexPointLocalRef - return read/write access to a point in local array 83 84 Not Collective 85 86 Input Arguments: 87 + dm - DM defining topological space 88 . point - topological point 89 - array - array to index into 90 91 Output Arguments: 92 . ptr - address of reference to point data, type generic so user can place in structure 93 94 Level: intermediate 95 96 Note: 97 A common usage when data sizes are known statically: 98 99 $ struct { PetscScalar foo,bar,baz; } *ptr; 100 $ DMPlexPointLocalRef(dm,point,array,&ptr); 101 $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 102 103 .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 104 @*/ 105 PetscErrorCode DMPlexPointLocalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 106 { 107 PetscErrorCode ierr; 108 PetscInt start, end; 109 110 PetscFunctionBegin; 111 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 112 PetscValidScalarPointer(array,3); 113 PetscValidPointer(ptr,4); 114 ierr = DMPlexGetLocalOffset_Private(dm,point,&start,&end);CHKERRQ(ierr); 115 *(PetscScalar**)ptr = (start < end) ? array + start : NULL; 116 PetscFunctionReturn(0); 117 } 118 119 #undef __FUNCT__ 120 #define __FUNCT__ "DMPlexGetPointLocalField" 121 /*@ 122 DMPlexGetPointLocalField - get location of point field data in local Vec 123 124 Not Collective 125 126 Input Arguments: 127 + dm - DM defining the topological space 128 . point - topological point 129 - field - the field number 130 131 Output Arguments: 132 + start - start of point data 133 - end - end of point data 134 135 Note: This is a half open interval [start, end) 136 137 Level: intermediate 138 139 .seealso: DMPlexGetPointLocal(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointLocalRead(), DMPlexPointLocalRead(), DMPlexPointLocalRef() 140 @*/ 141 PetscErrorCode DMPlexGetPointLocalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end) 142 { 143 PetscInt s, e; 144 PetscErrorCode ierr; 145 146 PetscFunctionBegin; 147 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 148 if (start) PetscValidPointer(start, 4); 149 if (end) PetscValidPointer(end, 5); 150 ierr = DMPlexGetLocalFieldOffset_Private(dm, point, field, &s, &e);CHKERRQ(ierr); 151 if (start) *start = s; 152 if (end) *end = e; 153 PetscFunctionReturn(0); 154 } 155 156 #undef __FUNCT__ 157 #define __FUNCT__ "DMPlexPointLocalFieldRead" 158 /*@ 159 DMPlexPointLocalFieldRead - return read access to a field on a point in local array 160 161 Not Collective 162 163 Input Arguments: 164 + dm - DM defining topological space 165 . point - topological point 166 . field - field number 167 - array - array to index into 168 169 Output Arguments: 170 . ptr - address of read reference to point data, type generic so user can place in structure 171 172 Level: intermediate 173 174 .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 175 @*/ 176 PetscErrorCode DMPlexPointLocalFieldRead(DM dm, PetscInt point,PetscInt field,const PetscScalar *array,const void *ptr) 177 { 178 PetscErrorCode ierr; 179 PetscInt start, end; 180 181 PetscFunctionBegin; 182 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 183 PetscValidScalarPointer(array,3); 184 PetscValidPointer(ptr,4); 185 ierr = DMPlexGetLocalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 186 *(const PetscScalar**)ptr = array + start; 187 PetscFunctionReturn(0); 188 } 189 190 #undef __FUNCT__ 191 #define __FUNCT__ "DMPlexPointLocalFieldRef" 192 /*@ 193 DMPlexPointLocalFieldRef - return read/write access to a field on a point in local array 194 195 Not Collective 196 197 Input Arguments: 198 + dm - DM defining topological space 199 . point - topological point 200 . field - field number 201 - array - array to index into 202 203 Output Arguments: 204 . ptr - address of reference to point data, type generic so user can place in structure 205 206 Level: intermediate 207 208 .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointLocal(), DMPlexPointGlobalRef() 209 @*/ 210 PetscErrorCode DMPlexPointLocalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr) 211 { 212 PetscErrorCode ierr; 213 PetscInt start, end; 214 215 PetscFunctionBegin; 216 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 217 PetscValidScalarPointer(array,3); 218 PetscValidPointer(ptr,4); 219 ierr = DMPlexGetLocalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 220 *(PetscScalar**)ptr = array + start; 221 PetscFunctionReturn(0); 222 } 223 224 #undef __FUNCT__ 225 #define __FUNCT__ "DMPlexGetPointGlobal" 226 /*@ 227 DMPlexGetPointGlobal - get location of point data in global Vec 228 229 Not Collective 230 231 Input Arguments: 232 + dm - DM defining the topological space 233 - point - topological point 234 235 Output Arguments: 236 + start - start of point data; returns -(globalStart+1) if point is not owned 237 - end - end of point data; returns -(globalEnd+1) if point is not owned 238 239 Note: This is a half open interval [start, end) 240 241 Level: intermediate 242 243 .seealso: DMPlexGetPointGlobalField(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef() 244 @*/ 245 PetscErrorCode DMPlexGetPointGlobal(DM dm, PetscInt point, PetscInt *start, PetscInt *end) 246 { 247 PetscInt s, e; 248 PetscErrorCode ierr; 249 250 PetscFunctionBegin; 251 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 252 if (start) PetscValidPointer(start, 3); 253 if (end) PetscValidPointer(end, 4); 254 ierr = DMPlexGetGlobalOffset_Private(dm, point, &s, &e);CHKERRQ(ierr); 255 if (start) *start = s; 256 if (end) *end = e; 257 PetscFunctionReturn(0); 258 } 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "DMPlexPointGlobalRead" 262 /*@ 263 DMPlexPointGlobalRead - return read access to a point in global array 264 265 Not Collective 266 267 Input Arguments: 268 + dm - DM defining topological space 269 . point - topological point 270 - array - array to index into 271 272 Output Arguments: 273 . ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 274 275 Level: intermediate 276 277 Note: 278 A common usage when data sizes are known statically: 279 280 $ const struct { PetscScalar foo,bar,baz; } *ptr; 281 $ DMPlexPointGlobalRead(dm,point,array,&ptr); 282 $ x = 2*ptr->foo + 3*ptr->bar + 5*ptr->baz; 283 284 .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef() 285 @*/ 286 PetscErrorCode DMPlexPointGlobalRead(DM dm,PetscInt point,const PetscScalar *array,const void *ptr) 287 { 288 PetscInt start, end; 289 PetscErrorCode ierr; 290 291 PetscFunctionBegin; 292 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 293 PetscValidScalarPointer(array, 3); 294 PetscValidPointer(ptr, 4); 295 ierr = DMPlexGetGlobalOffset_Private(dm, point, &start, &end);CHKERRQ(ierr); 296 *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 297 PetscFunctionReturn(0); 298 } 299 300 #undef __FUNCT__ 301 #define __FUNCT__ "DMPlexPointGlobalRef" 302 /*@ 303 DMPlexPointGlobalRef - return read/write access to a point in global array 304 305 Not Collective 306 307 Input Arguments: 308 + dm - DM defining topological space 309 . point - topological point 310 - array - array to index into 311 312 Output Arguments: 313 . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 314 315 Level: intermediate 316 317 Note: 318 A common usage when data sizes are known statically: 319 320 $ struct { PetscScalar foo,bar,baz; } *ptr; 321 $ DMPlexPointGlobalRef(dm,point,array,&ptr); 322 $ ptr->foo = 2; ptr->bar = 3; ptr->baz = 5; 323 324 .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead() 325 @*/ 326 PetscErrorCode DMPlexPointGlobalRef(DM dm,PetscInt point,PetscScalar *array,void *ptr) 327 { 328 PetscInt start, end; 329 PetscErrorCode ierr; 330 331 PetscFunctionBegin; 332 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 333 PetscValidScalarPointer(array, 3); 334 PetscValidPointer(ptr, 4); 335 ierr = DMPlexGetGlobalOffset_Private(dm, point, &start, &end);CHKERRQ(ierr); 336 *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 337 PetscFunctionReturn(0); 338 } 339 340 #undef __FUNCT__ 341 #define __FUNCT__ "DMPlexGetPointGlobalField" 342 /*@ 343 DMPlexGetPointGlobalField - get location of point field data in global Vec 344 345 Not Collective 346 347 Input Arguments: 348 + dm - DM defining the topological space 349 . point - topological point 350 - field - the field number 351 352 Output Arguments: 353 + start - start of point data; returns -(globalStart+1) if point is not owned 354 - end - end of point data; returns -(globalEnd+1) if point is not owned 355 356 Note: This is a half open interval [start, end) 357 358 Level: intermediate 359 360 .seealso: DMPlexGetPointGlobal(), DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexPointGlobalRead(), DMPlexGetPointLocal(), DMPlexPointGlobalRead(), DMPlexPointGlobalRef() 361 @*/ 362 PetscErrorCode DMPlexGetPointGlobalField(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end) 363 { 364 PetscInt s, e; 365 PetscErrorCode ierr; 366 367 PetscFunctionBegin; 368 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 369 if (start) PetscValidPointer(start, 4); 370 if (end) PetscValidPointer(end, 5); 371 ierr = DMPlexGetGlobalFieldOffset_Private(dm, point, field, &s, &e);CHKERRQ(ierr); 372 if (start) *start = s; 373 if (end) *end = e; 374 PetscFunctionReturn(0); 375 } 376 377 #undef __FUNCT__ 378 #define __FUNCT__ "DMPlexPointGlobalFieldRead" 379 /*@ 380 DMPlexPointGlobalFieldRead - return read access to a field on a point in global array 381 382 Not Collective 383 384 Input Arguments: 385 + dm - DM defining topological space 386 . point - topological point 387 . field - field number 388 - array - array to index into 389 390 Output Arguments: 391 . ptr - address of read reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 392 393 Level: intermediate 394 395 .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRead(), DMPlexPointGlobalRef() 396 @*/ 397 PetscErrorCode DMPlexPointGlobalFieldRead(DM dm,PetscInt point,PetscInt field,const PetscScalar *array,const void *ptr) 398 { 399 PetscInt start, end; 400 PetscErrorCode ierr; 401 402 PetscFunctionBegin; 403 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 404 PetscValidScalarPointer(array, 3); 405 PetscValidPointer(ptr, 4); 406 ierr = DMPlexGetGlobalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 407 *(const PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 408 PetscFunctionReturn(0); 409 } 410 411 #undef __FUNCT__ 412 #define __FUNCT__ "DMPlexPointGlobalFieldRef" 413 /*@ 414 DMPlexPointGlobalFieldRef - return read/write access to a field on a point in global array 415 416 Not Collective 417 418 Input Arguments: 419 + dm - DM defining topological space 420 . point - topological point 421 . field - field number 422 - array - array to index into 423 424 Output Arguments: 425 . ptr - address of reference to point data, type generic so user can place in structure; returns NULL if global point is not owned 426 427 Level: intermediate 428 429 .seealso: DMGetDefaultSection(), PetscSectionGetOffset(), PetscSectionGetDof(), DMPlexGetPointGlobal(), DMPlexPointLocalRef(), DMPlexPointGlobalRead() 430 @*/ 431 PetscErrorCode DMPlexPointGlobalFieldRef(DM dm,PetscInt point,PetscInt field,PetscScalar *array,void *ptr) 432 { 433 PetscInt start, end; 434 PetscErrorCode ierr; 435 436 PetscFunctionBegin; 437 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 438 PetscValidScalarPointer(array, 3); 439 PetscValidPointer(ptr, 4); 440 ierr = DMPlexGetGlobalFieldOffset_Private(dm, point, field, &start, &end);CHKERRQ(ierr); 441 *(PetscScalar**) ptr = (start < end) ? array + start - dm->map->rstart : NULL; 442 PetscFunctionReturn(0); 443 } 444