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