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