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