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