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