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