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