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