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