1 2 #include <../src/dm/impls/composite/packimpl.h> /*I "petscdmcomposite.h" I*/ 3 #include <petsc/private/isimpl.h> 4 #include <petsc/private/glvisviewerimpl.h> 5 #include <petscds.h> 6 7 /*@C 8 DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the 9 separate components `DM` in a `DMCOMPOSITE` to build the correct matrix nonzero structure. 10 11 Logically Collective 12 13 Input Parameters: 14 + dm - the composite object 15 - formcouplelocations - routine to set the nonzero locations in the matrix 16 17 Level: advanced 18 19 Note: 20 See `DMSetApplicationContext()` and `DMGetApplicationContext()` for how to get user information into 21 this routine 22 23 Fortran Note: 24 Not available from Fortran 25 26 .seealso: `DMCOMPOSITE`, `DM` 27 @*/ 28 PetscErrorCode DMCompositeSetCoupling(DM dm, PetscErrorCode (*FormCoupleLocations)(DM, Mat, PetscInt *, PetscInt *, PetscInt, PetscInt, PetscInt, PetscInt)) 29 { 30 DM_Composite *com = (DM_Composite *)dm->data; 31 PetscBool flg; 32 33 PetscFunctionBegin; 34 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 35 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 36 com->FormCoupleLocations = FormCoupleLocations; 37 PetscFunctionReturn(PETSC_SUCCESS); 38 } 39 40 PetscErrorCode DMDestroy_Composite(DM dm) 41 { 42 struct DMCompositeLink *next, *prev; 43 DM_Composite *com = (DM_Composite *)dm->data; 44 45 PetscFunctionBegin; 46 next = com->next; 47 while (next) { 48 prev = next; 49 next = next->next; 50 PetscCall(DMDestroy(&prev->dm)); 51 PetscCall(PetscFree(prev->grstarts)); 52 PetscCall(PetscFree(prev)); 53 } 54 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL)); 55 /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 56 PetscCall(PetscFree(com)); 57 PetscFunctionReturn(PETSC_SUCCESS); 58 } 59 60 PetscErrorCode DMView_Composite(DM dm, PetscViewer v) 61 { 62 PetscBool iascii; 63 DM_Composite *com = (DM_Composite *)dm->data; 64 65 PetscFunctionBegin; 66 PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii)); 67 if (iascii) { 68 struct DMCompositeLink *lnk = com->next; 69 PetscInt i; 70 71 PetscCall(PetscViewerASCIIPrintf(v, "DM (%s)\n", ((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix")); 72 PetscCall(PetscViewerASCIIPrintf(v, " contains %" PetscInt_FMT " DMs\n", com->nDM)); 73 PetscCall(PetscViewerASCIIPushTab(v)); 74 for (i = 0; lnk; lnk = lnk->next, i++) { 75 PetscCall(PetscViewerASCIIPrintf(v, "Link %" PetscInt_FMT ": DM of type %s\n", i, ((PetscObject)lnk->dm)->type_name)); 76 PetscCall(PetscViewerASCIIPushTab(v)); 77 PetscCall(DMView(lnk->dm, v)); 78 PetscCall(PetscViewerASCIIPopTab(v)); 79 } 80 PetscCall(PetscViewerASCIIPopTab(v)); 81 } 82 PetscFunctionReturn(PETSC_SUCCESS); 83 } 84 85 /* --------------------------------------------------------------------------------------*/ 86 PetscErrorCode DMSetUp_Composite(DM dm) 87 { 88 PetscInt nprev = 0; 89 PetscMPIInt rank, size; 90 DM_Composite *com = (DM_Composite *)dm->data; 91 struct DMCompositeLink *next = com->next; 92 PetscLayout map; 93 94 PetscFunctionBegin; 95 PetscCheck(!com->setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Packer has already been setup"); 96 PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &map)); 97 PetscCall(PetscLayoutSetLocalSize(map, com->n)); 98 PetscCall(PetscLayoutSetSize(map, PETSC_DETERMINE)); 99 PetscCall(PetscLayoutSetBlockSize(map, 1)); 100 PetscCall(PetscLayoutSetUp(map)); 101 PetscCall(PetscLayoutGetSize(map, &com->N)); 102 PetscCall(PetscLayoutGetRange(map, &com->rstart, NULL)); 103 PetscCall(PetscLayoutDestroy(&map)); 104 105 /* now set the rstart for each linked vector */ 106 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 107 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 108 while (next) { 109 next->rstart = nprev; 110 nprev += next->n; 111 next->grstart = com->rstart + next->rstart; 112 PetscCall(PetscMalloc1(size, &next->grstarts)); 113 PetscCallMPI(MPI_Allgather(&next->grstart, 1, MPIU_INT, next->grstarts, 1, MPIU_INT, PetscObjectComm((PetscObject)dm))); 114 next = next->next; 115 } 116 com->setup = PETSC_TRUE; 117 PetscFunctionReturn(PETSC_SUCCESS); 118 } 119 120 /* ----------------------------------------------------------------------------------*/ 121 122 /*@ 123 DMCompositeGetNumberDM - Get's the number of `DM` objects in the `DMCOMPOSITE` 124 representation. 125 126 Not Collective 127 128 Input Parameter: 129 . dm - the `DMCOMPOSITE` object 130 131 Output Parameter: 132 . nDM - the number of `DM` 133 134 Level: beginner 135 136 .seealso: `DMCOMPOSITE`, `DM` 137 @*/ 138 PetscErrorCode DMCompositeGetNumberDM(DM dm, PetscInt *nDM) 139 { 140 DM_Composite *com = (DM_Composite *)dm->data; 141 PetscBool flg; 142 143 PetscFunctionBegin; 144 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 145 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 146 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 147 *nDM = com->nDM; 148 PetscFunctionReturn(PETSC_SUCCESS); 149 } 150 151 /*@C 152 DMCompositeGetAccess - Allows one to access the individual packed vectors in their global 153 representation. 154 155 Collective on dm 156 157 Input Parameters: 158 + dm - the `DMCOMPOSITE` object 159 - gvec - the global vector 160 161 Output Parameters: 162 . Vec* ... - the packed parallel vectors, NULL for those that are not needed 163 164 Level: advanced 165 166 Note: 167 Use `DMCompositeRestoreAccess()` to return the vectors when you no longer need them 168 169 Fortran Note: 170 Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4) 171 or use the alternative interface `DMCompositeGetAccessArray()`. 172 173 .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetEntries()`, `DMCompositeScatter()` 174 @*/ 175 PetscErrorCode DMCompositeGetAccess(DM dm, Vec gvec, ...) 176 { 177 va_list Argp; 178 struct DMCompositeLink *next; 179 DM_Composite *com = (DM_Composite *)dm->data; 180 PetscInt readonly; 181 PetscBool flg; 182 183 PetscFunctionBegin; 184 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 185 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 186 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 187 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 188 next = com->next; 189 if (!com->setup) PetscCall(DMSetUp(dm)); 190 191 PetscCall(VecLockGet(gvec, &readonly)); 192 /* loop over packed objects, handling one at at time */ 193 va_start(Argp, gvec); 194 while (next) { 195 Vec *vec; 196 vec = va_arg(Argp, Vec *); 197 if (vec) { 198 PetscCall(DMGetGlobalVector(next->dm, vec)); 199 if (readonly) { 200 const PetscScalar *array; 201 PetscCall(VecGetArrayRead(gvec, &array)); 202 PetscCall(VecPlaceArray(*vec, array + next->rstart)); 203 PetscCall(VecLockReadPush(*vec)); 204 PetscCall(VecRestoreArrayRead(gvec, &array)); 205 } else { 206 PetscScalar *array; 207 PetscCall(VecGetArray(gvec, &array)); 208 PetscCall(VecPlaceArray(*vec, array + next->rstart)); 209 PetscCall(VecRestoreArray(gvec, &array)); 210 } 211 } 212 next = next->next; 213 } 214 va_end(Argp); 215 PetscFunctionReturn(PETSC_SUCCESS); 216 } 217 218 /*@C 219 DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global 220 representation. 221 222 Collective on dm 223 224 Input Parameters: 225 + dm - the `DMCOMPOSITE` 226 . pvec - packed vector 227 . nwanted - number of vectors wanted 228 - wanted - sorted array of vectors wanted, or NULL to get all vectors 229 230 Output Parameters: 231 . vecs - array of requested global vectors (must be allocated) 232 233 Level: advanced 234 235 Note: 236 Use `DMCompositeRestoreAccessArray()` to return the vectors when you no longer need them 237 238 .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 239 @*/ 240 PetscErrorCode DMCompositeGetAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 241 { 242 struct DMCompositeLink *link; 243 PetscInt i, wnum; 244 DM_Composite *com = (DM_Composite *)dm->data; 245 PetscInt readonly; 246 PetscBool flg; 247 248 PetscFunctionBegin; 249 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 250 PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 251 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 252 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 253 if (!com->setup) PetscCall(DMSetUp(dm)); 254 255 PetscCall(VecLockGet(pvec, &readonly)); 256 for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 257 if (!wanted || i == wanted[wnum]) { 258 Vec v; 259 PetscCall(DMGetGlobalVector(link->dm, &v)); 260 if (readonly) { 261 const PetscScalar *array; 262 PetscCall(VecGetArrayRead(pvec, &array)); 263 PetscCall(VecPlaceArray(v, array + link->rstart)); 264 PetscCall(VecLockReadPush(v)); 265 PetscCall(VecRestoreArrayRead(pvec, &array)); 266 } else { 267 PetscScalar *array; 268 PetscCall(VecGetArray(pvec, &array)); 269 PetscCall(VecPlaceArray(v, array + link->rstart)); 270 PetscCall(VecRestoreArray(pvec, &array)); 271 } 272 vecs[wnum++] = v; 273 } 274 } 275 PetscFunctionReturn(PETSC_SUCCESS); 276 } 277 278 /*@C 279 DMCompositeGetLocalAccessArray - Allows one to access the individual 280 packed vectors in their local representation. 281 282 Collective on dm. 283 284 Input Parameters: 285 + dm - the `DMCOMPOSITE` 286 . pvec - packed vector 287 . nwanted - number of vectors wanted 288 - wanted - sorted array of vectors wanted, or NULL to get all vectors 289 290 Output Parameters: 291 . vecs - array of requested local vectors (must be allocated) 292 293 Level: advanced 294 295 Note: 296 Use `DMCompositeRestoreLocalAccessArray()` to return the vectors 297 when you no longer need them. 298 299 .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`, 300 `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 301 @*/ 302 PetscErrorCode DMCompositeGetLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 303 { 304 struct DMCompositeLink *link; 305 PetscInt i, wnum; 306 DM_Composite *com = (DM_Composite *)dm->data; 307 PetscInt readonly; 308 PetscInt nlocal = 0; 309 PetscBool flg; 310 311 PetscFunctionBegin; 312 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 313 PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 314 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 315 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 316 if (!com->setup) PetscCall(DMSetUp(dm)); 317 318 PetscCall(VecLockGet(pvec, &readonly)); 319 for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 320 if (!wanted || i == wanted[wnum]) { 321 Vec v; 322 PetscCall(DMGetLocalVector(link->dm, &v)); 323 if (readonly) { 324 const PetscScalar *array; 325 PetscCall(VecGetArrayRead(pvec, &array)); 326 PetscCall(VecPlaceArray(v, array + nlocal)); 327 // this method does not make sense. The local vectors are not updated with a global-to-local and the user can not do it because it is locked 328 PetscCall(VecLockReadPush(v)); 329 PetscCall(VecRestoreArrayRead(pvec, &array)); 330 } else { 331 PetscScalar *array; 332 PetscCall(VecGetArray(pvec, &array)); 333 PetscCall(VecPlaceArray(v, array + nlocal)); 334 PetscCall(VecRestoreArray(pvec, &array)); 335 } 336 vecs[wnum++] = v; 337 } 338 339 nlocal += link->nlocal; 340 } 341 342 PetscFunctionReturn(PETSC_SUCCESS); 343 } 344 345 /*@C 346 DMCompositeRestoreAccess - Returns the vectors obtained with `DMCompositeGetAccess()` 347 representation. 348 349 Collective on dm 350 351 Input Parameters: 352 + dm - the `DMCOMPOSITE` object 353 . gvec - the global vector 354 - Vec* ... - the individual parallel vectors, NULL for those that are not needed 355 356 Level: advanced 357 358 .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 359 `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`, 360 `DMCompositeRestoreAccess()`, `DMCompositeGetAccess()` 361 @*/ 362 PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...) 363 { 364 va_list Argp; 365 struct DMCompositeLink *next; 366 DM_Composite *com = (DM_Composite *)dm->data; 367 PetscInt readonly; 368 PetscBool flg; 369 370 PetscFunctionBegin; 371 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 372 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 373 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 374 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 375 next = com->next; 376 if (!com->setup) PetscCall(DMSetUp(dm)); 377 378 PetscCall(VecLockGet(gvec, &readonly)); 379 /* loop over packed objects, handling one at at time */ 380 va_start(Argp, gvec); 381 while (next) { 382 Vec *vec; 383 vec = va_arg(Argp, Vec *); 384 if (vec) { 385 PetscCall(VecResetArray(*vec)); 386 if (readonly) PetscCall(VecLockReadPop(*vec)); 387 PetscCall(DMRestoreGlobalVector(next->dm, vec)); 388 } 389 next = next->next; 390 } 391 va_end(Argp); 392 PetscFunctionReturn(PETSC_SUCCESS); 393 } 394 395 /*@C 396 DMCompositeRestoreAccessArray - Returns the vectors obtained with `DMCompositeGetAccessArray()` 397 398 Collective on dm 399 400 Input Parameters: 401 + dm - the `DMCOMPOSITE` object 402 . pvec - packed vector 403 . nwanted - number of vectors wanted 404 . wanted - sorted array of vectors wanted, or NULL to get all vectors 405 - vecs - array of global vectors to return 406 407 Level: advanced 408 409 .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 410 @*/ 411 PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 412 { 413 struct DMCompositeLink *link; 414 PetscInt i, wnum; 415 DM_Composite *com = (DM_Composite *)dm->data; 416 PetscInt readonly; 417 PetscBool flg; 418 419 PetscFunctionBegin; 420 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 421 PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 422 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 423 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 424 if (!com->setup) PetscCall(DMSetUp(dm)); 425 426 PetscCall(VecLockGet(pvec, &readonly)); 427 for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 428 if (!wanted || i == wanted[wnum]) { 429 PetscCall(VecResetArray(vecs[wnum])); 430 if (readonly) PetscCall(VecLockReadPop(vecs[wnum])); 431 PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum])); 432 wnum++; 433 } 434 } 435 PetscFunctionReturn(PETSC_SUCCESS); 436 } 437 438 /*@C 439 DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with `DMCompositeGetLocalAccessArray()`. 440 441 Collective on dm. 442 443 Input Parameters: 444 + dm - the `DMCOMPOSITE` object 445 . pvec - packed vector 446 . nwanted - number of vectors wanted 447 . wanted - sorted array of vectors wanted, or NULL to restore all vectors 448 - vecs - array of local vectors to return 449 450 Level: advanced 451 452 Note: 453 nwanted and wanted must match the values given to `DMCompositeGetLocalAccessArray()` 454 otherwise the call will fail. 455 456 .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`, 457 `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, 458 `DMCompositeScatter()`, `DMCompositeGather()` 459 @*/ 460 PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 461 { 462 struct DMCompositeLink *link; 463 PetscInt i, wnum; 464 DM_Composite *com = (DM_Composite *)dm->data; 465 PetscInt readonly; 466 PetscBool flg; 467 468 PetscFunctionBegin; 469 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 470 PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 471 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 472 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 473 if (!com->setup) PetscCall(DMSetUp(dm)); 474 475 PetscCall(VecLockGet(pvec, &readonly)); 476 for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 477 if (!wanted || i == wanted[wnum]) { 478 PetscCall(VecResetArray(vecs[wnum])); 479 if (readonly) PetscCall(VecLockReadPop(vecs[wnum])); 480 PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum])); 481 wnum++; 482 } 483 } 484 PetscFunctionReturn(PETSC_SUCCESS); 485 } 486 487 /*@C 488 DMCompositeScatter - Scatters from a global packed vector into its individual local vectors 489 490 Collective on dm 491 492 Input Parameters: 493 + dm - the `DMCOMPOSITE` object 494 . gvec - the global vector 495 - Vec ... - the individual sequential vectors, NULL for those that are not needed 496 497 Level: advanced 498 499 Note: 500 `DMCompositeScatterArray()` is a non-variadic alternative that is often more convenient for library callers and is 501 accessible from Fortran. 502 503 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 504 `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 505 `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 506 `DMCompositeScatterArray()` 507 @*/ 508 PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...) 509 { 510 va_list Argp; 511 struct DMCompositeLink *next; 512 PETSC_UNUSED PetscInt cnt; 513 DM_Composite *com = (DM_Composite *)dm->data; 514 PetscBool flg; 515 516 PetscFunctionBegin; 517 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 518 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 519 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 520 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 521 if (!com->setup) PetscCall(DMSetUp(dm)); 522 523 /* loop over packed objects, handling one at at time */ 524 va_start(Argp, gvec); 525 for (cnt = 3, next = com->next; next; cnt++, next = next->next) { 526 Vec local; 527 local = va_arg(Argp, Vec); 528 if (local) { 529 Vec global; 530 const PetscScalar *array; 531 PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt)); 532 PetscCall(DMGetGlobalVector(next->dm, &global)); 533 PetscCall(VecGetArrayRead(gvec, &array)); 534 PetscCall(VecPlaceArray(global, array + next->rstart)); 535 PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local)); 536 PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local)); 537 PetscCall(VecRestoreArrayRead(gvec, &array)); 538 PetscCall(VecResetArray(global)); 539 PetscCall(DMRestoreGlobalVector(next->dm, &global)); 540 } 541 } 542 va_end(Argp); 543 PetscFunctionReturn(PETSC_SUCCESS); 544 } 545 546 /*@ 547 DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors 548 549 Collective on dm 550 551 Input Parameters: 552 + dm - the `DMCOMPOSITE` object 553 . gvec - the global vector 554 - lvecs - array of local vectors, NULL for any that are not needed 555 556 Level: advanced 557 558 Note: 559 This is a non-variadic alternative to `DMCompositeScatter()` 560 561 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()` 562 `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 563 `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 564 @*/ 565 PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs) 566 { 567 struct DMCompositeLink *next; 568 PetscInt i; 569 DM_Composite *com = (DM_Composite *)dm->data; 570 PetscBool flg; 571 572 PetscFunctionBegin; 573 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 574 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 575 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 576 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 577 if (!com->setup) PetscCall(DMSetUp(dm)); 578 579 /* loop over packed objects, handling one at at time */ 580 for (i = 0, next = com->next; next; next = next->next, i++) { 581 if (lvecs[i]) { 582 Vec global; 583 const PetscScalar *array; 584 PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3); 585 PetscCall(DMGetGlobalVector(next->dm, &global)); 586 PetscCall(VecGetArrayRead(gvec, &array)); 587 PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart)); 588 PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i])); 589 PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i])); 590 PetscCall(VecRestoreArrayRead(gvec, &array)); 591 PetscCall(VecResetArray(global)); 592 PetscCall(DMRestoreGlobalVector(next->dm, &global)); 593 } 594 } 595 PetscFunctionReturn(PETSC_SUCCESS); 596 } 597 598 /*@C 599 DMCompositeGather - Gathers into a global packed vector from its individual local vectors 600 601 Collective on dm 602 603 Input Parameters: 604 + dm - the `DMCOMPOSITE` object 605 . gvec - the global vector 606 . imode - `INSERT_VALUES` or `ADD_VALUES` 607 - Vec ... - the individual sequential vectors, NULL for any that are not needed 608 609 Level: advanced 610 611 Fortran Note: 612 Fortran users should use `DMCompositeGatherArray()` 613 614 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 615 `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 616 `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 617 @*/ 618 PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...) 619 { 620 va_list Argp; 621 struct DMCompositeLink *next; 622 DM_Composite *com = (DM_Composite *)dm->data; 623 PETSC_UNUSED PetscInt cnt; 624 PetscBool flg; 625 626 PetscFunctionBegin; 627 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 628 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3); 629 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 630 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 631 if (!com->setup) PetscCall(DMSetUp(dm)); 632 633 /* loop over packed objects, handling one at at time */ 634 va_start(Argp, gvec); 635 for (cnt = 3, next = com->next; next; cnt++, next = next->next) { 636 Vec local; 637 local = va_arg(Argp, Vec); 638 if (local) { 639 PetscScalar *array; 640 Vec global; 641 PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt)); 642 PetscCall(DMGetGlobalVector(next->dm, &global)); 643 PetscCall(VecGetArray(gvec, &array)); 644 PetscCall(VecPlaceArray(global, array + next->rstart)); 645 PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global)); 646 PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global)); 647 PetscCall(VecRestoreArray(gvec, &array)); 648 PetscCall(VecResetArray(global)); 649 PetscCall(DMRestoreGlobalVector(next->dm, &global)); 650 } 651 } 652 va_end(Argp); 653 PetscFunctionReturn(PETSC_SUCCESS); 654 } 655 656 /*@ 657 DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors 658 659 Collective on dm 660 661 Input Parameters: 662 + dm - the `DMCOMPOSITE` object 663 . gvec - the global vector 664 . imode - `INSERT_VALUES` or `ADD_VALUES` 665 - lvecs - the individual sequential vectors, NULL for any that are not needed 666 667 Level: advanced 668 669 Note: 670 This is a non-variadic alternative to `DMCompositeGather()`. 671 672 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 673 `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 674 `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`, 675 @*/ 676 PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs) 677 { 678 struct DMCompositeLink *next; 679 DM_Composite *com = (DM_Composite *)dm->data; 680 PetscInt i; 681 PetscBool flg; 682 683 PetscFunctionBegin; 684 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 685 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3); 686 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 687 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 688 if (!com->setup) PetscCall(DMSetUp(dm)); 689 690 /* loop over packed objects, handling one at at time */ 691 for (next = com->next, i = 0; next; next = next->next, i++) { 692 if (lvecs[i]) { 693 PetscScalar *array; 694 Vec global; 695 PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4); 696 PetscCall(DMGetGlobalVector(next->dm, &global)); 697 PetscCall(VecGetArray(gvec, &array)); 698 PetscCall(VecPlaceArray(global, array + next->rstart)); 699 PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global)); 700 PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global)); 701 PetscCall(VecRestoreArray(gvec, &array)); 702 PetscCall(VecResetArray(global)); 703 PetscCall(DMRestoreGlobalVector(next->dm, &global)); 704 } 705 } 706 PetscFunctionReturn(PETSC_SUCCESS); 707 } 708 709 /*@ 710 DMCompositeAddDM - adds a `DM` vector to a `DMCOMPOSITE` 711 712 Collective on dm 713 714 Input Parameters: 715 + dmc - the `DMCOMPOSITE` object 716 - dm - the `DM` object 717 718 Level: advanced 719 720 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeGather()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 721 `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 722 `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 723 @*/ 724 PetscErrorCode DMCompositeAddDM(DM dmc, DM dm) 725 { 726 PetscInt n, nlocal; 727 struct DMCompositeLink *mine, *next; 728 Vec global, local; 729 DM_Composite *com = (DM_Composite *)dmc->data; 730 PetscBool flg; 731 732 PetscFunctionBegin; 733 PetscValidHeaderSpecific(dmc, DM_CLASSID, 1); 734 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 735 PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg)); 736 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 737 next = com->next; 738 PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite"); 739 740 /* create new link */ 741 PetscCall(PetscNew(&mine)); 742 PetscCall(PetscObjectReference((PetscObject)dm)); 743 PetscCall(DMGetGlobalVector(dm, &global)); 744 PetscCall(VecGetLocalSize(global, &n)); 745 PetscCall(DMRestoreGlobalVector(dm, &global)); 746 PetscCall(DMGetLocalVector(dm, &local)); 747 PetscCall(VecGetSize(local, &nlocal)); 748 PetscCall(DMRestoreLocalVector(dm, &local)); 749 750 mine->n = n; 751 mine->nlocal = nlocal; 752 mine->dm = dm; 753 mine->next = NULL; 754 com->n += n; 755 com->nghost += nlocal; 756 757 /* add to end of list */ 758 if (!next) com->next = mine; 759 else { 760 while (next->next) next = next->next; 761 next->next = mine; 762 } 763 com->nDM++; 764 com->nmine++; 765 PetscFunctionReturn(PETSC_SUCCESS); 766 } 767 768 #include <petscdraw.h> 769 PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer); 770 PetscErrorCode VecView_DMComposite(Vec gvec, PetscViewer viewer) 771 { 772 DM dm; 773 struct DMCompositeLink *next; 774 PetscBool isdraw; 775 DM_Composite *com; 776 777 PetscFunctionBegin; 778 PetscCall(VecGetDM(gvec, &dm)); 779 PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite"); 780 com = (DM_Composite *)dm->data; 781 next = com->next; 782 783 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 784 if (!isdraw) { 785 /* do I really want to call this? */ 786 PetscCall(VecView_MPI(gvec, viewer)); 787 } else { 788 PetscInt cnt = 0; 789 790 /* loop over packed objects, handling one at at time */ 791 while (next) { 792 Vec vec; 793 const PetscScalar *array; 794 PetscInt bs; 795 796 /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */ 797 PetscCall(DMGetGlobalVector(next->dm, &vec)); 798 PetscCall(VecGetArrayRead(gvec, &array)); 799 PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart)); 800 PetscCall(VecRestoreArrayRead(gvec, &array)); 801 PetscCall(VecView(vec, viewer)); 802 PetscCall(VecResetArray(vec)); 803 PetscCall(VecGetBlockSize(vec, &bs)); 804 PetscCall(DMRestoreGlobalVector(next->dm, &vec)); 805 PetscCall(PetscViewerDrawBaseAdd(viewer, bs)); 806 cnt += bs; 807 next = next->next; 808 } 809 PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt)); 810 } 811 PetscFunctionReturn(PETSC_SUCCESS); 812 } 813 814 PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec) 815 { 816 DM_Composite *com = (DM_Composite *)dm->data; 817 818 PetscFunctionBegin; 819 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 820 PetscCall(DMSetFromOptions(dm)); 821 PetscCall(DMSetUp(dm)); 822 PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec)); 823 PetscCall(VecSetType(*gvec, dm->vectype)); 824 PetscCall(VecSetSizes(*gvec, com->n, com->N)); 825 PetscCall(VecSetDM(*gvec, dm)); 826 PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (void (*)(void))VecView_DMComposite)); 827 PetscFunctionReturn(PETSC_SUCCESS); 828 } 829 830 PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec) 831 { 832 DM_Composite *com = (DM_Composite *)dm->data; 833 834 PetscFunctionBegin; 835 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 836 if (!com->setup) { 837 PetscCall(DMSetFromOptions(dm)); 838 PetscCall(DMSetUp(dm)); 839 } 840 PetscCall(VecCreate(PETSC_COMM_SELF, lvec)); 841 PetscCall(VecSetType(*lvec, dm->vectype)); 842 PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE)); 843 PetscCall(VecSetDM(*lvec, dm)); 844 PetscFunctionReturn(PETSC_SUCCESS); 845 } 846 847 /*@C 848 DMCompositeGetISLocalToGlobalMappings - gets an `ISLocalToGlobalMapping` for each `DM` in the `DMCOMPOSITE`, maps to the composite global space 849 850 Collective on dm 851 852 Input Parameter: 853 . dm - the `DMCOMPOSITE` object 854 855 Output Parameters: 856 . ltogs - the individual mappings for each packed vector. Note that this includes 857 all the ghost points that individual ghosted `DMDA` may have. 858 859 Level: advanced 860 861 Note: 862 Each entry of ltogs should be destroyed with `ISLocalToGlobalMappingDestroy()`, the ltogs array should be freed with `PetscFree()`. 863 864 Fortran Note: 865 Not available from Fortran 866 867 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 868 `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 869 `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 870 @*/ 871 PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping **ltogs) 872 { 873 PetscInt i, *idx, n, cnt; 874 struct DMCompositeLink *next; 875 PetscMPIInt rank; 876 DM_Composite *com = (DM_Composite *)dm->data; 877 PetscBool flg; 878 879 PetscFunctionBegin; 880 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 881 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 882 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 883 PetscCall(DMSetUp(dm)); 884 PetscCall(PetscMalloc1(com->nDM, ltogs)); 885 next = com->next; 886 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 887 888 /* loop over packed objects, handling one at at time */ 889 cnt = 0; 890 while (next) { 891 ISLocalToGlobalMapping ltog; 892 PetscMPIInt size; 893 const PetscInt *suboff, *indices; 894 Vec global; 895 896 /* Get sub-DM global indices for each local dof */ 897 PetscCall(DMGetLocalToGlobalMapping(next->dm, <og)); 898 PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n)); 899 PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices)); 900 PetscCall(PetscMalloc1(n, &idx)); 901 902 /* Get the offsets for the sub-DM global vector */ 903 PetscCall(DMGetGlobalVector(next->dm, &global)); 904 PetscCall(VecGetOwnershipRanges(global, &suboff)); 905 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size)); 906 907 /* Shift the sub-DM definition of the global space to the composite global space */ 908 for (i = 0; i < n; i++) { 909 PetscInt subi = indices[i], lo = 0, hi = size, t; 910 /* There's no consensus on what a negative index means, 911 except for skipping when setting the values in vectors and matrices */ 912 if (subi < 0) { 913 idx[i] = subi - next->grstarts[rank]; 914 continue; 915 } 916 /* Binary search to find which rank owns subi */ 917 while (hi - lo > 1) { 918 t = lo + (hi - lo) / 2; 919 if (suboff[t] > subi) hi = t; 920 else lo = t; 921 } 922 idx[i] = subi - suboff[lo] + next->grstarts[lo]; 923 } 924 PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices)); 925 PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt])); 926 PetscCall(DMRestoreGlobalVector(next->dm, &global)); 927 next = next->next; 928 cnt++; 929 } 930 PetscFunctionReturn(PETSC_SUCCESS); 931 } 932 933 /*@C 934 DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector 935 936 Not Collective 937 938 Input Parameter: 939 . dm - the `DMCOMPOSITE` 940 941 Output Parameter: 942 . is - array of serial index sets for each each component of the `DMCOMPOSITE` 943 944 Level: intermediate 945 946 Notes: 947 At present, a composite local vector does not normally exist. This function is used to provide index sets for 948 `MatGetLocalSubMatrix()`. In the future, the scatters for each entry in the `DMCOMPOSITE` may be be merged into a single 949 scatter to a composite local vector. The user should not typically need to know which is being done. 950 951 To get the composite global indices at all local points (including ghosts), use `DMCompositeGetISLocalToGlobalMappings()`. 952 953 To get index sets for pieces of the composite global vector, use `DMCompositeGetGlobalISs()`. 954 955 Each returned `IS` should be destroyed with `ISDestroy()`, the array should be freed with `PetscFree()`. 956 957 Fortran Note: 958 Not available from Fortran 959 960 .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`, `MatCreateLocalRef()` 961 @*/ 962 PetscErrorCode DMCompositeGetLocalISs(DM dm, IS **is) 963 { 964 DM_Composite *com = (DM_Composite *)dm->data; 965 struct DMCompositeLink *link; 966 PetscInt cnt, start; 967 PetscBool flg; 968 969 PetscFunctionBegin; 970 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 971 PetscValidPointer(is, 2); 972 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 973 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 974 PetscCall(PetscMalloc1(com->nmine, is)); 975 for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) { 976 PetscInt bs; 977 PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt])); 978 PetscCall(DMGetBlockSize(link->dm, &bs)); 979 PetscCall(ISSetBlockSize((*is)[cnt], bs)); 980 } 981 PetscFunctionReturn(PETSC_SUCCESS); 982 } 983 984 /*@C 985 DMCompositeGetGlobalISs - Gets the index sets for each composed object in a `DMCOMPOSITE` 986 987 Collective on dm 988 989 Input Parameter: 990 . dm - the `DMCOMPOSITE` object 991 992 Output Parameters: 993 . is - the array of index sets 994 995 Level: advanced 996 997 Notes: 998 The is entries should be destroyed with `ISDestroy()`, the is array should be freed with `PetscFree()` 999 1000 These could be used to extract a subset of vector entries for a "multi-physics" preconditioner 1001 1002 Use `DMCompositeGetLocalISs()` for index sets in the packed local numbering, and 1003 `DMCompositeGetISLocalToGlobalMappings()` for to map local sub-`DM` (including ghost) indices to packed global 1004 indices. 1005 1006 Fortran Note: 1007 The output argument 'is' must be an allocated array of sufficient length, which can be learned using `DMCompositeGetNumberDM()`. 1008 1009 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1010 `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 1011 `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 1012 @*/ 1013 PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[]) 1014 { 1015 PetscInt cnt = 0; 1016 struct DMCompositeLink *next; 1017 PetscMPIInt rank; 1018 DM_Composite *com = (DM_Composite *)dm->data; 1019 PetscBool flg; 1020 1021 PetscFunctionBegin; 1022 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1023 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1024 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 1025 PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before"); 1026 PetscCall(PetscMalloc1(com->nDM, is)); 1027 next = com->next; 1028 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1029 1030 /* loop over packed objects, handling one at at time */ 1031 while (next) { 1032 PetscDS prob; 1033 1034 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt])); 1035 PetscCall(DMGetDS(dm, &prob)); 1036 if (prob) { 1037 MatNullSpace space; 1038 Mat pmat; 1039 PetscObject disc; 1040 PetscInt Nf; 1041 1042 PetscCall(PetscDSGetNumFields(prob, &Nf)); 1043 if (cnt < Nf) { 1044 PetscCall(PetscDSGetDiscretization(prob, cnt, &disc)); 1045 PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space)); 1046 if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space)); 1047 PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space)); 1048 if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space)); 1049 PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat)); 1050 if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat)); 1051 } 1052 } 1053 cnt++; 1054 next = next->next; 1055 } 1056 PetscFunctionReturn(PETSC_SUCCESS); 1057 } 1058 1059 PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) 1060 { 1061 PetscInt nDM; 1062 DM *dms; 1063 PetscInt i; 1064 1065 PetscFunctionBegin; 1066 PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 1067 if (numFields) *numFields = nDM; 1068 PetscCall(DMCompositeGetGlobalISs(dm, fields)); 1069 if (fieldNames) { 1070 PetscCall(PetscMalloc1(nDM, &dms)); 1071 PetscCall(PetscMalloc1(nDM, fieldNames)); 1072 PetscCall(DMCompositeGetEntriesArray(dm, dms)); 1073 for (i = 0; i < nDM; i++) { 1074 char buf[256]; 1075 const char *splitname; 1076 1077 /* Split naming precedence: object name, prefix, number */ 1078 splitname = ((PetscObject)dm)->name; 1079 if (!splitname) { 1080 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname)); 1081 if (splitname) { 1082 size_t len; 1083 PetscCall(PetscStrncpy(buf, splitname, sizeof(buf))); 1084 buf[sizeof(buf) - 1] = 0; 1085 PetscCall(PetscStrlen(buf, &len)); 1086 if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */ 1087 splitname = buf; 1088 } 1089 } 1090 if (!splitname) { 1091 PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i)); 1092 splitname = buf; 1093 } 1094 PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i])); 1095 } 1096 PetscCall(PetscFree(dms)); 1097 } 1098 PetscFunctionReturn(PETSC_SUCCESS); 1099 } 1100 1101 /* 1102 This could take over from DMCreateFieldIS(), as it is more general, 1103 making DMCreateFieldIS() a special case -- calling with dmlist == NULL; 1104 At this point it's probably best to be less intrusive, however. 1105 */ 1106 PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 1107 { 1108 PetscInt nDM; 1109 PetscInt i; 1110 1111 PetscFunctionBegin; 1112 PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist)); 1113 if (dmlist) { 1114 PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 1115 PetscCall(PetscMalloc1(nDM, dmlist)); 1116 PetscCall(DMCompositeGetEntriesArray(dm, *dmlist)); 1117 for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i]))); 1118 } 1119 PetscFunctionReturn(PETSC_SUCCESS); 1120 } 1121 1122 /* -------------------------------------------------------------------------------------*/ 1123 /*@C 1124 DMCompositeGetLocalVectors - Gets local vectors for each part of a `DMCOMPOSITE` 1125 Use `DMCompositeRestoreLocalVectors()` to return them. 1126 1127 Not Collective 1128 1129 Input Parameter: 1130 . dm - the `DMCOMPOSITE` object 1131 1132 Output Parameter: 1133 . Vec ... - the individual sequential Vecs 1134 1135 Level: advanced 1136 1137 Fortran Note: 1138 Not available from Fortran 1139 1140 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1141 `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1142 `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 1143 @*/ 1144 PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...) 1145 { 1146 va_list Argp; 1147 struct DMCompositeLink *next; 1148 DM_Composite *com = (DM_Composite *)dm->data; 1149 PetscBool flg; 1150 1151 PetscFunctionBegin; 1152 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1153 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1154 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 1155 next = com->next; 1156 /* loop over packed objects, handling one at at time */ 1157 va_start(Argp, dm); 1158 while (next) { 1159 Vec *vec; 1160 vec = va_arg(Argp, Vec *); 1161 if (vec) PetscCall(DMGetLocalVector(next->dm, vec)); 1162 next = next->next; 1163 } 1164 va_end(Argp); 1165 PetscFunctionReturn(PETSC_SUCCESS); 1166 } 1167 1168 /*@C 1169 DMCompositeRestoreLocalVectors - Restores local vectors for each part of a `DMCOMPOSITE` 1170 1171 Not Collective 1172 1173 Input Parameter: 1174 . dm - the `DMCOMPOSITE` object 1175 1176 Output Parameter: 1177 . Vec ... - the individual sequential `Vec` 1178 1179 Level: advanced 1180 1181 Fortran Note: 1182 Not available from Fortran 1183 1184 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1185 `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1186 `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 1187 @*/ 1188 PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...) 1189 { 1190 va_list Argp; 1191 struct DMCompositeLink *next; 1192 DM_Composite *com = (DM_Composite *)dm->data; 1193 PetscBool flg; 1194 1195 PetscFunctionBegin; 1196 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1197 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1198 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 1199 next = com->next; 1200 /* loop over packed objects, handling one at at time */ 1201 va_start(Argp, dm); 1202 while (next) { 1203 Vec *vec; 1204 vec = va_arg(Argp, Vec *); 1205 if (vec) PetscCall(DMRestoreLocalVector(next->dm, vec)); 1206 next = next->next; 1207 } 1208 va_end(Argp); 1209 PetscFunctionReturn(PETSC_SUCCESS); 1210 } 1211 1212 /* -------------------------------------------------------------------------------------*/ 1213 /*@C 1214 DMCompositeGetEntries - Gets the `DM` for each entry in a `DMCOMPOSITE`. 1215 1216 Not Collective 1217 1218 Input Parameter: 1219 . dm - the `DMCOMPOSITE` object 1220 1221 Output Parameter: 1222 . DM ... - the individual entries `DM` 1223 1224 Level: advanced 1225 1226 Fortran Note: 1227 Available as `DMCompositeGetEntries()` for one output `DM`, DMCompositeGetEntries2() for 2, etc 1228 1229 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()` 1230 `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1231 `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, 1232 `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()` 1233 @*/ 1234 PetscErrorCode DMCompositeGetEntries(DM dm, ...) 1235 { 1236 va_list Argp; 1237 struct DMCompositeLink *next; 1238 DM_Composite *com = (DM_Composite *)dm->data; 1239 PetscBool flg; 1240 1241 PetscFunctionBegin; 1242 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1243 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1244 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 1245 next = com->next; 1246 /* loop over packed objects, handling one at at time */ 1247 va_start(Argp, dm); 1248 while (next) { 1249 DM *dmn; 1250 dmn = va_arg(Argp, DM *); 1251 if (dmn) *dmn = next->dm; 1252 next = next->next; 1253 } 1254 va_end(Argp); 1255 PetscFunctionReturn(PETSC_SUCCESS); 1256 } 1257 1258 /*@C 1259 DMCompositeGetEntriesArray - Gets the DM for each entry in a `DMCOMPOSITE` 1260 1261 Not Collective 1262 1263 Input Parameter: 1264 . dm - the `DMCOMPOSITE` object 1265 1266 Output Parameter: 1267 . dms - array of sufficient length (see `DMCompositeGetNumberDM()`) to hold the individual `DM` 1268 1269 Level: advanced 1270 1271 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()` 1272 `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1273 `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, 1274 `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()` 1275 @*/ 1276 PetscErrorCode DMCompositeGetEntriesArray(DM dm, DM dms[]) 1277 { 1278 struct DMCompositeLink *next; 1279 DM_Composite *com = (DM_Composite *)dm->data; 1280 PetscInt i; 1281 PetscBool flg; 1282 1283 PetscFunctionBegin; 1284 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1285 PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1286 PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 1287 /* loop over packed objects, handling one at at time */ 1288 for (next = com->next, i = 0; next; next = next->next, i++) dms[i] = next->dm; 1289 PetscFunctionReturn(PETSC_SUCCESS); 1290 } 1291 1292 typedef struct { 1293 DM dm; 1294 PetscViewer *subv; 1295 Vec *vecs; 1296 } GLVisViewerCtx; 1297 1298 static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 1299 { 1300 GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 1301 PetscInt i, n; 1302 1303 PetscFunctionBegin; 1304 PetscCall(DMCompositeGetNumberDM(ctx->dm, &n)); 1305 for (i = 0; i < n; i++) PetscCall(PetscViewerDestroy(&ctx->subv[i])); 1306 PetscCall(PetscFree2(ctx->subv, ctx->vecs)); 1307 PetscCall(DMDestroy(&ctx->dm)); 1308 PetscCall(PetscFree(ctx)); 1309 PetscFunctionReturn(PETSC_SUCCESS); 1310 } 1311 1312 static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 1313 { 1314 Vec X = (Vec)oX; 1315 GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 1316 PetscInt i, n, cumf; 1317 1318 PetscFunctionBegin; 1319 PetscCall(DMCompositeGetNumberDM(ctx->dm, &n)); 1320 PetscCall(DMCompositeGetAccessArray(ctx->dm, X, n, NULL, ctx->vecs)); 1321 for (i = 0, cumf = 0; i < n; i++) { 1322 PetscErrorCode (*g2l)(PetscObject, PetscInt, PetscObject[], void *); 1323 void *fctx; 1324 PetscInt nfi; 1325 1326 PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nfi, NULL, NULL, &g2l, NULL, &fctx)); 1327 if (!nfi) continue; 1328 if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i], nfi, oXfield + cumf, fctx)); 1329 else PetscCall(VecCopy(ctx->vecs[i], (Vec)(oXfield[cumf]))); 1330 cumf += nfi; 1331 } 1332 PetscCall(DMCompositeRestoreAccessArray(ctx->dm, X, n, NULL, ctx->vecs)); 1333 PetscFunctionReturn(PETSC_SUCCESS); 1334 } 1335 1336 static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer) 1337 { 1338 DM dm = (DM)odm, *dms; 1339 Vec *Ufds; 1340 GLVisViewerCtx *ctx; 1341 PetscInt i, n, tnf, *sdim; 1342 char **fecs; 1343 1344 PetscFunctionBegin; 1345 PetscCall(PetscNew(&ctx)); 1346 PetscCall(PetscObjectReference((PetscObject)dm)); 1347 ctx->dm = dm; 1348 PetscCall(DMCompositeGetNumberDM(dm, &n)); 1349 PetscCall(PetscMalloc1(n, &dms)); 1350 PetscCall(DMCompositeGetEntriesArray(dm, dms)); 1351 PetscCall(PetscMalloc2(n, &ctx->subv, n, &ctx->vecs)); 1352 for (i = 0, tnf = 0; i < n; i++) { 1353 PetscInt nf; 1354 1355 PetscCall(PetscViewerCreate(PetscObjectComm(odm), &ctx->subv[i])); 1356 PetscCall(PetscViewerSetType(ctx->subv[i], PETSCVIEWERGLVIS)); 1357 PetscCall(PetscViewerGLVisSetDM_Private(ctx->subv[i], (PetscObject)dms[i])); 1358 PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nf, NULL, NULL, NULL, NULL, NULL)); 1359 tnf += nf; 1360 } 1361 PetscCall(PetscFree(dms)); 1362 PetscCall(PetscMalloc3(tnf, &fecs, tnf, &sdim, tnf, &Ufds)); 1363 for (i = 0, tnf = 0; i < n; i++) { 1364 PetscInt *sd, nf, f; 1365 const char **fec; 1366 Vec *Uf; 1367 1368 PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nf, &fec, &sd, NULL, (PetscObject **)&Uf, NULL)); 1369 for (f = 0; f < nf; f++) { 1370 PetscCall(PetscStrallocpy(fec[f], &fecs[tnf + f])); 1371 Ufds[tnf + f] = Uf[f]; 1372 sdim[tnf + f] = sd[f]; 1373 } 1374 tnf += nf; 1375 } 1376 PetscCall(PetscViewerGLVisSetFields(viewer, tnf, (const char **)fecs, sdim, DMCompositeSampleGLVisFields_Private, (PetscObject *)Ufds, ctx, DestroyGLVisViewerCtx_Private)); 1377 for (i = 0; i < tnf; i++) PetscCall(PetscFree(fecs[i])); 1378 PetscCall(PetscFree3(fecs, sdim, Ufds)); 1379 PetscFunctionReturn(PETSC_SUCCESS); 1380 } 1381 1382 PetscErrorCode DMRefine_Composite(DM dmi, MPI_Comm comm, DM *fine) 1383 { 1384 struct DMCompositeLink *next; 1385 DM_Composite *com = (DM_Composite *)dmi->data; 1386 DM dm; 1387 1388 PetscFunctionBegin; 1389 PetscValidHeaderSpecific(dmi, DM_CLASSID, 1); 1390 if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm)); 1391 PetscCall(DMSetUp(dmi)); 1392 next = com->next; 1393 PetscCall(DMCompositeCreate(comm, fine)); 1394 1395 /* loop over packed objects, handling one at at time */ 1396 while (next) { 1397 PetscCall(DMRefine(next->dm, comm, &dm)); 1398 PetscCall(DMCompositeAddDM(*fine, dm)); 1399 PetscCall(PetscObjectDereference((PetscObject)dm)); 1400 next = next->next; 1401 } 1402 PetscFunctionReturn(PETSC_SUCCESS); 1403 } 1404 1405 PetscErrorCode DMCoarsen_Composite(DM dmi, MPI_Comm comm, DM *fine) 1406 { 1407 struct DMCompositeLink *next; 1408 DM_Composite *com = (DM_Composite *)dmi->data; 1409 DM dm; 1410 1411 PetscFunctionBegin; 1412 PetscValidHeaderSpecific(dmi, DM_CLASSID, 1); 1413 PetscCall(DMSetUp(dmi)); 1414 if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm)); 1415 next = com->next; 1416 PetscCall(DMCompositeCreate(comm, fine)); 1417 1418 /* loop over packed objects, handling one at at time */ 1419 while (next) { 1420 PetscCall(DMCoarsen(next->dm, comm, &dm)); 1421 PetscCall(DMCompositeAddDM(*fine, dm)); 1422 PetscCall(PetscObjectDereference((PetscObject)dm)); 1423 next = next->next; 1424 } 1425 PetscFunctionReturn(PETSC_SUCCESS); 1426 } 1427 1428 PetscErrorCode DMCreateInterpolation_Composite(DM coarse, DM fine, Mat *A, Vec *v) 1429 { 1430 PetscInt m, n, M, N, nDM, i; 1431 struct DMCompositeLink *nextc; 1432 struct DMCompositeLink *nextf; 1433 Vec gcoarse, gfine, *vecs; 1434 DM_Composite *comcoarse = (DM_Composite *)coarse->data; 1435 DM_Composite *comfine = (DM_Composite *)fine->data; 1436 Mat *mats; 1437 1438 PetscFunctionBegin; 1439 PetscValidHeaderSpecific(coarse, DM_CLASSID, 1); 1440 PetscValidHeaderSpecific(fine, DM_CLASSID, 2); 1441 PetscCall(DMSetUp(coarse)); 1442 PetscCall(DMSetUp(fine)); 1443 /* use global vectors only for determining matrix layout */ 1444 PetscCall(DMGetGlobalVector(coarse, &gcoarse)); 1445 PetscCall(DMGetGlobalVector(fine, &gfine)); 1446 PetscCall(VecGetLocalSize(gcoarse, &n)); 1447 PetscCall(VecGetLocalSize(gfine, &m)); 1448 PetscCall(VecGetSize(gcoarse, &N)); 1449 PetscCall(VecGetSize(gfine, &M)); 1450 PetscCall(DMRestoreGlobalVector(coarse, &gcoarse)); 1451 PetscCall(DMRestoreGlobalVector(fine, &gfine)); 1452 1453 nDM = comfine->nDM; 1454 PetscCheck(nDM == comcoarse->nDM, PetscObjectComm((PetscObject)fine), PETSC_ERR_ARG_INCOMP, "Fine DMComposite has %" PetscInt_FMT " entries, but coarse has %" PetscInt_FMT, nDM, comcoarse->nDM); 1455 PetscCall(PetscCalloc1(nDM * nDM, &mats)); 1456 if (v) PetscCall(PetscCalloc1(nDM, &vecs)); 1457 1458 /* loop over packed objects, handling one at at time */ 1459 for (nextc = comcoarse->next, nextf = comfine->next, i = 0; nextc; nextc = nextc->next, nextf = nextf->next, i++) { 1460 if (!v) PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], NULL)); 1461 else PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], &vecs[i])); 1462 } 1463 PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, nDM, NULL, mats, A)); 1464 if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, vecs, v)); 1465 for (i = 0; i < nDM * nDM; i++) PetscCall(MatDestroy(&mats[i])); 1466 PetscCall(PetscFree(mats)); 1467 if (v) { 1468 for (i = 0; i < nDM; i++) PetscCall(VecDestroy(&vecs[i])); 1469 PetscCall(PetscFree(vecs)); 1470 } 1471 PetscFunctionReturn(PETSC_SUCCESS); 1472 } 1473 1474 static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm) 1475 { 1476 DM_Composite *com = (DM_Composite *)dm->data; 1477 ISLocalToGlobalMapping *ltogs; 1478 PetscInt i; 1479 1480 PetscFunctionBegin; 1481 /* Set the ISLocalToGlobalMapping on the new matrix */ 1482 PetscCall(DMCompositeGetISLocalToGlobalMappings(dm, <ogs)); 1483 PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm), com->nDM, ltogs, &dm->ltogmap)); 1484 for (i = 0; i < com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(<ogs[i])); 1485 PetscCall(PetscFree(ltogs)); 1486 PetscFunctionReturn(PETSC_SUCCESS); 1487 } 1488 1489 PetscErrorCode DMCreateColoring_Composite(DM dm, ISColoringType ctype, ISColoring *coloring) 1490 { 1491 PetscInt n, i, cnt; 1492 ISColoringValue *colors; 1493 PetscBool dense = PETSC_FALSE; 1494 ISColoringValue maxcol = 0; 1495 DM_Composite *com = (DM_Composite *)dm->data; 1496 1497 PetscFunctionBegin; 1498 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1499 PetscCheck(ctype != IS_COLORING_LOCAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only global coloring supported"); 1500 if (ctype == IS_COLORING_GLOBAL) { 1501 n = com->n; 1502 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unknown ISColoringType"); 1503 PetscCall(PetscMalloc1(n, &colors)); /* freed in ISColoringDestroy() */ 1504 1505 PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL)); 1506 if (dense) { 1507 for (i = 0; i < n; i++) colors[i] = (ISColoringValue)(com->rstart + i); 1508 maxcol = com->N; 1509 } else { 1510 struct DMCompositeLink *next = com->next; 1511 PetscMPIInt rank; 1512 1513 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1514 cnt = 0; 1515 while (next) { 1516 ISColoring lcoloring; 1517 1518 PetscCall(DMCreateColoring(next->dm, IS_COLORING_GLOBAL, &lcoloring)); 1519 for (i = 0; i < lcoloring->N; i++) colors[cnt++] = maxcol + lcoloring->colors[i]; 1520 maxcol += lcoloring->n; 1521 PetscCall(ISColoringDestroy(&lcoloring)); 1522 next = next->next; 1523 } 1524 } 1525 PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), maxcol, n, colors, PETSC_OWN_POINTER, coloring)); 1526 PetscFunctionReturn(PETSC_SUCCESS); 1527 } 1528 1529 PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) 1530 { 1531 struct DMCompositeLink *next; 1532 PetscScalar *garray, *larray; 1533 DM_Composite *com = (DM_Composite *)dm->data; 1534 1535 PetscFunctionBegin; 1536 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1537 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 1538 1539 if (!com->setup) PetscCall(DMSetUp(dm)); 1540 1541 PetscCall(VecGetArray(gvec, &garray)); 1542 PetscCall(VecGetArray(lvec, &larray)); 1543 1544 /* loop over packed objects, handling one at at time */ 1545 next = com->next; 1546 while (next) { 1547 Vec local, global; 1548 PetscInt N; 1549 1550 PetscCall(DMGetGlobalVector(next->dm, &global)); 1551 PetscCall(VecGetLocalSize(global, &N)); 1552 PetscCall(VecPlaceArray(global, garray)); 1553 PetscCall(DMGetLocalVector(next->dm, &local)); 1554 PetscCall(VecPlaceArray(local, larray)); 1555 PetscCall(DMGlobalToLocalBegin(next->dm, global, mode, local)); 1556 PetscCall(DMGlobalToLocalEnd(next->dm, global, mode, local)); 1557 PetscCall(VecResetArray(global)); 1558 PetscCall(VecResetArray(local)); 1559 PetscCall(DMRestoreGlobalVector(next->dm, &global)); 1560 PetscCall(DMRestoreLocalVector(next->dm, &local)); 1561 1562 larray += next->nlocal; 1563 garray += next->n; 1564 next = next->next; 1565 } 1566 1567 PetscCall(VecRestoreArray(gvec, NULL)); 1568 PetscCall(VecRestoreArray(lvec, NULL)); 1569 PetscFunctionReturn(PETSC_SUCCESS); 1570 } 1571 1572 PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) 1573 { 1574 PetscFunctionBegin; 1575 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1576 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 1577 PetscValidHeaderSpecific(lvec, VEC_CLASSID, 4); 1578 PetscFunctionReturn(PETSC_SUCCESS); 1579 } 1580 1581 PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1582 { 1583 struct DMCompositeLink *next; 1584 PetscScalar *larray, *garray; 1585 DM_Composite *com = (DM_Composite *)dm->data; 1586 1587 PetscFunctionBegin; 1588 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1589 PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 1590 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 1591 1592 if (!com->setup) PetscCall(DMSetUp(dm)); 1593 1594 PetscCall(VecGetArray(lvec, &larray)); 1595 PetscCall(VecGetArray(gvec, &garray)); 1596 1597 /* loop over packed objects, handling one at at time */ 1598 next = com->next; 1599 while (next) { 1600 Vec global, local; 1601 1602 PetscCall(DMGetLocalVector(next->dm, &local)); 1603 PetscCall(VecPlaceArray(local, larray)); 1604 PetscCall(DMGetGlobalVector(next->dm, &global)); 1605 PetscCall(VecPlaceArray(global, garray)); 1606 PetscCall(DMLocalToGlobalBegin(next->dm, local, mode, global)); 1607 PetscCall(DMLocalToGlobalEnd(next->dm, local, mode, global)); 1608 PetscCall(VecResetArray(local)); 1609 PetscCall(VecResetArray(global)); 1610 PetscCall(DMRestoreGlobalVector(next->dm, &global)); 1611 PetscCall(DMRestoreLocalVector(next->dm, &local)); 1612 1613 garray += next->n; 1614 larray += next->nlocal; 1615 next = next->next; 1616 } 1617 1618 PetscCall(VecRestoreArray(gvec, NULL)); 1619 PetscCall(VecRestoreArray(lvec, NULL)); 1620 1621 PetscFunctionReturn(PETSC_SUCCESS); 1622 } 1623 1624 PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1625 { 1626 PetscFunctionBegin; 1627 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1628 PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 1629 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 1630 PetscFunctionReturn(PETSC_SUCCESS); 1631 } 1632 1633 PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2) 1634 { 1635 struct DMCompositeLink *next; 1636 PetscScalar *array1, *array2; 1637 DM_Composite *com = (DM_Composite *)dm->data; 1638 1639 PetscFunctionBegin; 1640 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1641 PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2); 1642 PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4); 1643 1644 if (!com->setup) PetscCall(DMSetUp(dm)); 1645 1646 PetscCall(VecGetArray(vec1, &array1)); 1647 PetscCall(VecGetArray(vec2, &array2)); 1648 1649 /* loop over packed objects, handling one at at time */ 1650 next = com->next; 1651 while (next) { 1652 Vec local1, local2; 1653 1654 PetscCall(DMGetLocalVector(next->dm, &local1)); 1655 PetscCall(VecPlaceArray(local1, array1)); 1656 PetscCall(DMGetLocalVector(next->dm, &local2)); 1657 PetscCall(VecPlaceArray(local2, array2)); 1658 PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2)); 1659 PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2)); 1660 PetscCall(VecResetArray(local2)); 1661 PetscCall(DMRestoreLocalVector(next->dm, &local2)); 1662 PetscCall(VecResetArray(local1)); 1663 PetscCall(DMRestoreLocalVector(next->dm, &local1)); 1664 1665 array1 += next->nlocal; 1666 array2 += next->nlocal; 1667 next = next->next; 1668 } 1669 1670 PetscCall(VecRestoreArray(vec1, NULL)); 1671 PetscCall(VecRestoreArray(vec2, NULL)); 1672 1673 PetscFunctionReturn(PETSC_SUCCESS); 1674 } 1675 1676 PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1677 { 1678 PetscFunctionBegin; 1679 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1680 PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 1681 PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 1682 PetscFunctionReturn(PETSC_SUCCESS); 1683 } 1684 1685 /*MC 1686 DMCOMPOSITE = "composite" - A `DM` object that is used to manage data for a collection of `DM` 1687 1688 Level: intermediate 1689 1690 .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()` 1691 M*/ 1692 1693 PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p) 1694 { 1695 DM_Composite *com; 1696 1697 PetscFunctionBegin; 1698 PetscCall(PetscNew(&com)); 1699 p->data = com; 1700 com->n = 0; 1701 com->nghost = 0; 1702 com->next = NULL; 1703 com->nDM = 0; 1704 1705 p->ops->createglobalvector = DMCreateGlobalVector_Composite; 1706 p->ops->createlocalvector = DMCreateLocalVector_Composite; 1707 p->ops->getlocaltoglobalmapping = DMGetLocalToGlobalMapping_Composite; 1708 p->ops->createfieldis = DMCreateFieldIS_Composite; 1709 p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite; 1710 p->ops->refine = DMRefine_Composite; 1711 p->ops->coarsen = DMCoarsen_Composite; 1712 p->ops->createinterpolation = DMCreateInterpolation_Composite; 1713 p->ops->creatematrix = DMCreateMatrix_Composite; 1714 p->ops->getcoloring = DMCreateColoring_Composite; 1715 p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite; 1716 p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite; 1717 p->ops->localtoglobalbegin = DMLocalToGlobalBegin_Composite; 1718 p->ops->localtoglobalend = DMLocalToGlobalEnd_Composite; 1719 p->ops->localtolocalbegin = DMLocalToLocalBegin_Composite; 1720 p->ops->localtolocalend = DMLocalToLocalEnd_Composite; 1721 p->ops->destroy = DMDestroy_Composite; 1722 p->ops->view = DMView_Composite; 1723 p->ops->setup = DMSetUp_Composite; 1724 1725 PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite)); 1726 PetscFunctionReturn(PETSC_SUCCESS); 1727 } 1728 1729 /*@ 1730 DMCompositeCreate - Creates a `DMCOMPOSITE`, used to generate "composite" 1731 vectors made up of several subvectors. 1732 1733 Collective 1734 1735 Input Parameter: 1736 . comm - the processors that will share the global vector 1737 1738 Output Parameters: 1739 . packer - the `DMCOMPOSITE` object 1740 1741 Level: advanced 1742 1743 .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCOMPOSITE`, `DMCreate()` 1744 `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()` 1745 `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 1746 @*/ 1747 PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer) 1748 { 1749 PetscFunctionBegin; 1750 PetscValidPointer(packer, 2); 1751 PetscCall(DMCreate(comm, packer)); 1752 PetscCall(DMSetType(*packer, DMCOMPOSITE)); 1753 PetscFunctionReturn(PETSC_SUCCESS); 1754 } 1755