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