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