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