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