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