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