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