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