Actual source code: pack.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <../src/dm/impls/composite/packimpl.h>       /*I  "petscdmcomposite.h"  I*/
  3: #include <petscds.h>

  7: /*@C
  8:     DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
  9:       seperate components (DMs) in a DMto build the correct matrix nonzero structure.


 12:     Logically Collective on MPI_Comm

 14:     Input Parameter:
 15: +   dm - the composite object
 16: -   formcouplelocations - routine to set the nonzero locations in the matrix

 18:     Level: advanced

 20:     Notes: See DMSetApplicationContext() and DMGetApplicationContext() for how to get user information into
 21:         this routine

 23: @*/
 24: PetscErrorCode  DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt))
 25: {
 26:   DM_Composite *com = (DM_Composite*)dm->data;

 29:   com->FormCoupleLocations = FormCoupleLocations;
 30:   return(0);
 31: }

 35: PetscErrorCode  DMDestroy_Composite(DM dm)
 36: {
 37:   PetscErrorCode         ierr;
 38:   struct DMCompositeLink *next, *prev;
 39:   DM_Composite           *com = (DM_Composite*)dm->data;

 42:   next = com->next;
 43:   while (next) {
 44:     prev = next;
 45:     next = next->next;
 46:     DMDestroy(&prev->dm);
 47:     PetscFree(prev->grstarts);
 48:     PetscFree(prev);
 49:   }
 50:   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
 51:   PetscFree(com);
 52:   return(0);
 53: }

 57: PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
 58: {
 60:   PetscBool      iascii;
 61:   DM_Composite   *com = (DM_Composite*)dm->data;

 64:   PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);
 65:   if (iascii) {
 66:     struct DMCompositeLink *lnk = com->next;
 67:     PetscInt               i;

 69:     PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix");
 70:     PetscViewerASCIIPrintf(v,"  contains %D DMs\n",com->nDM);
 71:     PetscViewerASCIIPushTab(v);
 72:     for (i=0; lnk; lnk=lnk->next,i++) {
 73:       PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);
 74:       PetscViewerASCIIPushTab(v);
 75:       DMView(lnk->dm,v);
 76:       PetscViewerASCIIPopTab(v);
 77:     }
 78:     PetscViewerASCIIPopTab(v);
 79:   }
 80:   return(0);
 81: }

 83: /* --------------------------------------------------------------------------------------*/
 86: PetscErrorCode  DMSetUp_Composite(DM dm)
 87: {
 88:   PetscErrorCode         ierr;
 89:   PetscInt               nprev = 0;
 90:   PetscMPIInt            rank,size;
 91:   DM_Composite           *com  = (DM_Composite*)dm->data;
 92:   struct DMCompositeLink *next = com->next;
 93:   PetscLayout            map;

 96:   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
 97:   PetscLayoutCreate(PetscObjectComm((PetscObject)dm),&map);
 98:   PetscLayoutSetLocalSize(map,com->n);
 99:   PetscLayoutSetSize(map,PETSC_DETERMINE);
100:   PetscLayoutSetBlockSize(map,1);
101:   PetscLayoutSetUp(map);
102:   PetscLayoutGetSize(map,&com->N);
103:   PetscLayoutGetRange(map,&com->rstart,NULL);
104:   PetscLayoutDestroy(&map);

106:   /* now set the rstart for each linked vector */
107:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
108:   MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);
109:   while (next) {
110:     next->rstart  = nprev;
111:     nprev        += next->n;
112:     next->grstart = com->rstart + next->rstart;
113:     PetscMalloc1(size,&next->grstarts);
114:     MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,PetscObjectComm((PetscObject)dm));
115:     next          = next->next;
116:   }
117:   com->setup = PETSC_TRUE;
118:   return(0);
119: }

121: /* ----------------------------------------------------------------------------------*/

125: /*@
126:     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
127:        representation.

129:     Not Collective

131:     Input Parameter:
132: .    dm - the packer object

134:     Output Parameter:
135: .     nDM - the number of DMs

137:     Level: beginner

139: @*/
140: PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
141: {
142:   DM_Composite *com = (DM_Composite*)dm->data;

146:   *nDM = com->nDM;
147:   return(0);
148: }


153: /*@C
154:     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
155:        representation.

157:     Collective on DMComposite

159:     Input Parameters:
160: +    dm - the packer object
161: -    gvec - the global vector

163:     Output Parameters:
164: .    Vec* ... - the packed parallel vectors, NULL for those that are not needed

166:     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them

168:     Fortran Notes:

170:     Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
171:     or use the alternative interface DMCompositeGetAccessArray().

173:     Level: advanced

175: .seealso: DMCompositeGetEntries(), DMCompositeScatter()
176: @*/
177: PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
178: {
179:   va_list                Argp;
180:   PetscErrorCode         ierr;
181:   struct DMCompositeLink *next;
182:   DM_Composite           *com = (DM_Composite*)dm->data;

187:   next = com->next;
188:   if (!com->setup) {
189:     DMSetUp(dm);
190:   }

192:   /* loop over packed objects, handling one at at time */
193:   va_start(Argp,gvec);
194:   while (next) {
195:     Vec *vec;
196:     vec = va_arg(Argp, Vec*);
197:     if (vec) {
198:       PetscScalar *array;
199:       DMGetGlobalVector(next->dm,vec);
200:       VecGetArray(gvec,&array);
201:       VecPlaceArray(*vec,array+next->rstart);
202:       VecRestoreArray(gvec,&array);
203:     }
204:     next = next->next;
205:   }
206:   va_end(Argp);
207:   return(0);
208: }

212: /*@C
213:     DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
214:        representation.

216:     Collective on DMComposite

218:     Input Parameters:
219: +    dm - the packer object
220: .    pvec - packed vector
221: .    nwanted - number of vectors wanted
222: -    wanted - sorted array of vectors wanted, or NULL to get all vectors

224:     Output Parameters:
225: .    vecs - array of requested global vectors (must be allocated)

227:     Notes: Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them

229:     Level: advanced

231: .seealso: DMCompositeGetAccess(), DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
232: @*/
233: PetscErrorCode  DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
234: {
235:   PetscErrorCode         ierr;
236:   struct DMCompositeLink *link;
237:   PetscInt               i,wnum;
238:   DM_Composite           *com = (DM_Composite*)dm->data;

243:   if (!com->setup) {
244:     DMSetUp(dm);
245:   }

247:   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
248:     if (!wanted || i == wanted[wnum]) {
249:       PetscScalar *array;
250:       Vec v;
251:       DMGetGlobalVector(link->dm,&v);
252:       VecGetArray(pvec,&array);
253:       VecPlaceArray(v,array+link->rstart);
254:       VecRestoreArray(pvec,&array);
255:       vecs[wnum++] = v;
256:     }
257:   }
258:   return(0);
259: }

263: /*@C
264:     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
265:        representation.

267:     Collective on DMComposite

269:     Input Parameters:
270: +    dm - the packer object
271: .    gvec - the global vector
272: -    Vec* ... - the individual parallel vectors, NULL for those that are not needed

274:     Level: advanced

276: .seealso  DMCompositeAddDM(), DMCreateGlobalVector(),
277:          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
278:          DMCompositeRestoreAccess(), DMCompositeGetAccess()

280: @*/
281: PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
282: {
283:   va_list                Argp;
284:   PetscErrorCode         ierr;
285:   struct DMCompositeLink *next;
286:   DM_Composite           *com = (DM_Composite*)dm->data;

291:   next = com->next;
292:   if (!com->setup) {
293:     DMSetUp(dm);
294:   }

296:   /* loop over packed objects, handling one at at time */
297:   va_start(Argp,gvec);
298:   while (next) {
299:     Vec *vec;
300:     vec = va_arg(Argp, Vec*);
301:     if (vec) {
302:       VecResetArray(*vec);
303:       DMRestoreGlobalVector(next->dm,vec);
304:     }
305:     next = next->next;
306:   }
307:   va_end(Argp);
308:   return(0);
309: }

313: /*@C
314:     DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray()

316:     Collective on DMComposite

318:     Input Parameters:
319: +    dm - the packer object
320: .    pvec - packed vector
321: .    nwanted - number of vectors wanted
322: .    wanted - sorted array of vectors wanted, or NULL to get all vectors
323: -    vecs - array of global vectors to return

325:     Level: advanced

327: .seealso: DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), DMCompositeScatter(), DMCompositeGather()
328: @*/
329: PetscErrorCode  DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
330: {
331:   PetscErrorCode         ierr;
332:   struct DMCompositeLink *link;
333:   PetscInt               i,wnum;
334:   DM_Composite           *com = (DM_Composite*)dm->data;

339:   if (!com->setup) {
340:     DMSetUp(dm);
341:   }

343:   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
344:     if (!wanted || i == wanted[wnum]) {
345:       VecResetArray(vecs[wnum]);
346:       DMRestoreGlobalVector(link->dm,&vecs[wnum]);
347:       wnum++;
348:     }
349:   }
350:   return(0);
351: }

355: /*@C
356:     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors

358:     Collective on DMComposite

360:     Input Parameters:
361: +    dm - the packer object
362: .    gvec - the global vector
363: -    Vec ... - the individual sequential vectors, NULL for those that are not needed

365:     Level: advanced

367:     Notes:
368:     DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is
369:     accessible from Fortran.

371: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
372:          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
373:          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
374:          DMCompositeScatterArray()

376: @*/
377: PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
378: {
379:   va_list                Argp;
380:   PetscErrorCode         ierr;
381:   struct DMCompositeLink *next;
382:   PetscInt               cnt;
383:   DM_Composite           *com = (DM_Composite*)dm->data;

388:   if (!com->setup) {
389:     DMSetUp(dm);
390:   }

392:   /* loop over packed objects, handling one at at time */
393:   va_start(Argp,gvec);
394:   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
395:     Vec local;
396:     local = va_arg(Argp, Vec);
397:     if (local) {
398:       Vec         global;
399:       PetscScalar *array;
401:       DMGetGlobalVector(next->dm,&global);
402:       VecGetArray(gvec,&array);
403:       VecPlaceArray(global,array+next->rstart);
404:       DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);
405:       DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);
406:       VecRestoreArray(gvec,&array);
407:       VecResetArray(global);
408:       DMRestoreGlobalVector(next->dm,&global);
409:     }
410:   }
411:   va_end(Argp);
412:   return(0);
413: }

417: /*@
418:     DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors

420:     Collective on DMComposite

422:     Input Parameters:
423: +    dm - the packer object
424: .    gvec - the global vector
425: .    lvecs - array of local vectors, NULL for any that are not needed

427:     Level: advanced

429:     Note:
430:     This is a non-variadic alternative to DMCompositeScatterArray()

432: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector()
433:          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
434:          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()

436: @*/
437: PetscErrorCode  DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs)
438: {
439:   PetscErrorCode         ierr;
440:   struct DMCompositeLink *next;
441:   PetscInt               i;
442:   DM_Composite           *com = (DM_Composite*)dm->data;

447:   if (!com->setup) {
448:     DMSetUp(dm);
449:   }

451:   /* loop over packed objects, handling one at at time */
452:   for (i=0,next=com->next; next; next=next->next,i++) {
453:     if (lvecs[i]) {
454:       Vec         global;
455:       PetscScalar *array;
457:       DMGetGlobalVector(next->dm,&global);
458:       VecGetArray(gvec,&array);
459:       VecPlaceArray(global,array+next->rstart);
460:       DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]);
461:       DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]);
462:       VecRestoreArray(gvec,&array);
463:       VecResetArray(global);
464:       DMRestoreGlobalVector(next->dm,&global);
465:     }
466:   }
467:   return(0);
468: }

472: /*@C
473:     DMCompositeGather - Gathers into a global packed vector from its individual local vectors

475:     Collective on DMComposite

477:     Input Parameter:
478: +    dm - the packer object
479: .    gvec - the global vector
480: -    Vec ... - the individual sequential vectors, NULL for any that are not needed

482:     Level: advanced

484: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
485:          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
486:          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()

488: @*/
489: PetscErrorCode  DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
490: {
491:   va_list                Argp;
492:   PetscErrorCode         ierr;
493:   struct DMCompositeLink *next;
494:   DM_Composite           *com = (DM_Composite*)dm->data;
495:   PetscInt               cnt;

500:   if (!com->setup) {
501:     DMSetUp(dm);
502:   }

504:   /* loop over packed objects, handling one at at time */
505:   va_start(Argp,imode);
506:   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
507:     Vec local;
508:     local = va_arg(Argp, Vec);
509:     if (local) {
510:       PetscScalar *array;
511:       Vec         global;
513:       DMGetGlobalVector(next->dm,&global);
514:       VecGetArray(gvec,&array);
515:       VecPlaceArray(global,array+next->rstart);
516:       DMLocalToGlobalBegin(next->dm,local,imode,global);
517:       DMLocalToGlobalEnd(next->dm,local,imode,global);
518:       VecRestoreArray(gvec,&array);
519:       VecResetArray(global);
520:       DMRestoreGlobalVector(next->dm,&global);
521:     }
522:   }
523:   va_end(Argp);
524:   return(0);
525: }

529: /*@
530:     DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors

532:     Collective on DMComposite

534:     Input Parameter:
535: +    dm - the packer object
536: .    gvec - the global vector
537: -    lvecs - the individual sequential vectors, NULL for any that are not needed

539:     Level: advanced

541:     Notes:
542:     This is a non-variadic alternative to DMCompositeGather().

544: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
545:          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
546:          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries(),
547: @*/
548: PetscErrorCode  DMCompositeGatherArray(DM dm,Vec gvec,InsertMode imode,Vec *lvecs)
549: {
550:   PetscErrorCode         ierr;
551:   struct DMCompositeLink *next;
552:   DM_Composite           *com = (DM_Composite*)dm->data;
553:   PetscInt               i;

558:   if (!com->setup) {
559:     DMSetUp(dm);
560:   }

562:   /* loop over packed objects, handling one at at time */
563:   for (next=com->next,i=0; next; next=next->next,i++) {
564:     if (lvecs[i]) {
565:       PetscScalar *array;
566:       Vec         global;
568:       DMGetGlobalVector(next->dm,&global);
569:       VecGetArray(gvec,&array);
570:       VecPlaceArray(global,array+next->rstart);
571:       DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global);
572:       DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global);
573:       VecRestoreArray(gvec,&array);
574:       VecResetArray(global);
575:       DMRestoreGlobalVector(next->dm,&global);
576:     }
577:   }
578:   return(0);
579: }

583: /*@C
584:     DMCompositeAddDM - adds a DM  vector to a DMComposite

586:     Collective on DMComposite

588:     Input Parameter:
589: +    dm - the packer object
590: -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)

592:     Level: advanced

594: .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
595:          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
596:          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()

598: @*/
599: PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
600: {
601:   PetscErrorCode         ierr;
602:   PetscInt               n,nlocal;
603:   struct DMCompositeLink *mine,*next;
604:   Vec                    global,local;
605:   DM_Composite           *com = (DM_Composite*)dmc->data;

610:   next = com->next;
611:   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");

613:   /* create new link */
614:   PetscNew(&mine);
615:   PetscObjectReference((PetscObject)dm);
616:   DMGetGlobalVector(dm,&global);
617:   VecGetLocalSize(global,&n);
618:   DMRestoreGlobalVector(dm,&global);
619:   DMGetLocalVector(dm,&local);
620:   VecGetSize(local,&nlocal);
621:   DMRestoreLocalVector(dm,&local);

623:   mine->n      = n;
624:   mine->nlocal = nlocal;
625:   mine->dm     = dm;
626:   mine->next   = NULL;
627:   com->n      += n;

629:   /* add to end of list */
630:   if (!next) com->next = mine;
631:   else {
632:     while (next->next) next = next->next;
633:     next->next = mine;
634:   }
635:   com->nDM++;
636:   com->nmine++;
637:   return(0);
638: }

640: #include <petscdraw.h>
641: PETSC_EXTERN PetscErrorCode  VecView_MPI(Vec,PetscViewer);
644: PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
645: {
646:   DM                     dm;
647:   PetscErrorCode         ierr;
648:   struct DMCompositeLink *next;
649:   PetscBool              isdraw;
650:   DM_Composite           *com;

653:   VecGetDM(gvec, &dm);
654:   if (!dm) SETERRQ(PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
655:   com  = (DM_Composite*)dm->data;
656:   next = com->next;

658:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
659:   if (!isdraw) {
660:     /* do I really want to call this? */
661:     VecView_MPI(gvec,viewer);
662:   } else {
663:     PetscInt cnt = 0;

665:     /* loop over packed objects, handling one at at time */
666:     while (next) {
667:       Vec         vec;
668:       PetscScalar *array;
669:       PetscInt    bs;

671:       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
672:       DMGetGlobalVector(next->dm,&vec);
673:       VecGetArray(gvec,&array);
674:       VecPlaceArray(vec,array+next->rstart);
675:       VecRestoreArray(gvec,&array);
676:       VecView(vec,viewer);
677:       VecGetBlockSize(vec,&bs);
678:       VecResetArray(vec);
679:       DMRestoreGlobalVector(next->dm,&vec);
680:       PetscViewerDrawBaseAdd(viewer,bs);
681:       cnt += bs;
682:       next = next->next;
683:     }
684:     PetscViewerDrawBaseAdd(viewer,-cnt);
685:   }
686:   return(0);
687: }

691: PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
692: {
694:   DM_Composite   *com = (DM_Composite*)dm->data;

698:   DMSetUp(dm);
699:   VecCreateMPI(PetscObjectComm((PetscObject)dm),com->n,com->N,gvec);
700:   VecSetDM(*gvec, dm);
701:   VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);
702:   return(0);
703: }

707: PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
708: {
710:   DM_Composite   *com = (DM_Composite*)dm->data;

714:   if (!com->setup) {
715:     DMSetUp(dm);
716:   }
717:   VecCreateSeq(PetscObjectComm((PetscObject)dm),com->nghost,lvec);
718:   VecSetDM(*lvec, dm);
719:   return(0);
720: }

724: /*@C
725:     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space

727:     Collective on DM

729:     Input Parameter:
730: .    dm - the packer object

732:     Output Parameters:
733: .    ltogs - the individual mappings for each packed vector. Note that this includes
734:            all the ghost points that individual ghosted DMDA's may have.

736:     Level: advanced

738:     Notes:
739:        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().

741: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
742:          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
743:          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()

745: @*/
746: PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
747: {
748:   PetscErrorCode         ierr;
749:   PetscInt               i,*idx,n,cnt;
750:   struct DMCompositeLink *next;
751:   PetscMPIInt            rank;
752:   DM_Composite           *com = (DM_Composite*)dm->data;

756:   DMSetUp(dm);
757:   PetscMalloc1((com->nDM),ltogs);
758:   next = com->next;
759:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);

761:   /* loop over packed objects, handling one at at time */
762:   cnt = 0;
763:   while (next) {
764:     ISLocalToGlobalMapping ltog;
765:     PetscMPIInt            size;
766:     const PetscInt         *suboff,*indices;
767:     Vec                    global;

769:     /* Get sub-DM global indices for each local dof */
770:     DMGetLocalToGlobalMapping(next->dm,&ltog);
771:     ISLocalToGlobalMappingGetSize(ltog,&n);
772:     ISLocalToGlobalMappingGetIndices(ltog,&indices);
773:     PetscMalloc1(n,&idx);

775:     /* Get the offsets for the sub-DM global vector */
776:     DMGetGlobalVector(next->dm,&global);
777:     VecGetOwnershipRanges(global,&suboff);
778:     MPI_Comm_size(PetscObjectComm((PetscObject)global),&size);

780:     /* Shift the sub-DM definition of the global space to the composite global space */
781:     for (i=0; i<n; i++) {
782:       PetscInt subi = indices[i],lo = 0,hi = size,t;
783:       /* Binary search to find which rank owns subi */
784:       while (hi-lo > 1) {
785:         t = lo + (hi-lo)/2;
786:         if (suboff[t] > subi) hi = t;
787:         else                  lo = t;
788:       }
789:       idx[i] = subi - suboff[lo] + next->grstarts[lo];
790:     }
791:     ISLocalToGlobalMappingRestoreIndices(ltog,&indices);
792:     ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);
793:     DMRestoreGlobalVector(next->dm,&global);
794:     next = next->next;
795:     cnt++;
796:   }
797:   return(0);
798: }

802: /*@C
803:    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector

805:    Not Collective

807:    Input Arguments:
808: . dm - composite DM

810:    Output Arguments:
811: . is - array of serial index sets for each each component of the DMComposite

813:    Level: intermediate

815:    Notes:
816:    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
817:    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
818:    scatter to a composite local vector.  The user should not typically need to know which is being done.

820:    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().

822:    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().

824:    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().

826: .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
827: @*/
828: PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
829: {
830:   PetscErrorCode         ierr;
831:   DM_Composite           *com = (DM_Composite*)dm->data;
832:   struct DMCompositeLink *link;
833:   PetscInt               cnt,start;

838:   PetscMalloc1(com->nmine,is);
839:   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
840:     PetscInt bs;
841:     ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);
842:     DMGetBlockSize(link->dm,&bs);
843:     ISSetBlockSize((*is)[cnt],bs);
844:   }
845:   return(0);
846: }

850: /*@C
851:     DMCompositeGetGlobalISs - Gets the index sets for each composed object

853:     Collective on DMComposite

855:     Input Parameter:
856: .    dm - the packer object

858:     Output Parameters:
859: .    is - the array of index sets

861:     Level: advanced

863:     Notes:
864:        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()

866:        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner

868:        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
869:        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
870:        indices.

872: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
873:          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
874:          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()

876: @*/

878: PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
879: {
880:   PetscErrorCode         ierr;
881:   PetscInt               cnt = 0;
882:   struct DMCompositeLink *next;
883:   PetscMPIInt            rank;
884:   DM_Composite           *com = (DM_Composite*)dm->data;

888:   PetscMalloc1((com->nDM),is);
889:   next = com->next;
890:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);

892:   /* loop over packed objects, handling one at at time */
893:   while (next) {
894:     ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt]);
895:     if (dm->prob) {
896:       MatNullSpace space;
897:       Mat          pmat;
898:       PetscObject  disc;
899:       PetscInt     Nf;

901:       PetscDSGetNumFields(dm->prob, &Nf);
902:       if (cnt < Nf) {
903:         PetscDSGetDiscretization(dm->prob, cnt, &disc);
904:         PetscObjectQuery(disc, "nullspace", (PetscObject*) &space);
905:         if (space) {PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space);}
906:         PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space);
907:         if (space) {PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space);}
908:         PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat);
909:         if (pmat) {PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat);}
910:       }
911:     }
912:     cnt++;
913:     next = next->next;
914:   }
915:   return(0);
916: }

920: PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields)
921: {
922:   PetscInt       nDM;
923:   DM             *dms;
924:   PetscInt       i;

928:   DMCompositeGetNumberDM(dm, &nDM);
929:   if (numFields) *numFields = nDM;
930:   DMCompositeGetGlobalISs(dm, fields);
931:   if (fieldNames) {
932:     PetscMalloc1(nDM, &dms);
933:     PetscMalloc1(nDM, fieldNames);
934:     DMCompositeGetEntriesArray(dm, dms);
935:     for (i=0; i<nDM; i++) {
936:       char       buf[256];
937:       const char *splitname;

939:       /* Split naming precedence: object name, prefix, number */
940:       splitname = ((PetscObject) dm)->name;
941:       if (!splitname) {
942:         PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);
943:         if (splitname) {
944:           size_t len;
945:           PetscStrncpy(buf,splitname,sizeof(buf));
946:           buf[sizeof(buf) - 1] = 0;
947:           PetscStrlen(buf,&len);
948:           if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
949:           splitname = buf;
950:         }
951:       }
952:       if (!splitname) {
953:         PetscSNPrintf(buf,sizeof(buf),"%D",i);
954:         splitname = buf;
955:       }
956:       PetscStrallocpy(splitname,&(*fieldNames)[i]);
957:     }
958:     PetscFree(dms);
959:   }
960:   return(0);
961: }

963: /*
964:  This could take over from DMCreateFieldIS(), as it is more general,
965:  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
966:  At this point it's probably best to be less intrusive, however.
967:  */
970: PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
971: {
972:   PetscInt       nDM;
973:   PetscInt       i;

977:   DMCreateFieldIS_Composite(dm, len, namelist, islist);
978:   if (dmlist) {
979:     DMCompositeGetNumberDM(dm, &nDM);
980:     PetscMalloc1(nDM, dmlist);
981:     DMCompositeGetEntriesArray(dm, *dmlist);
982:     for (i=0; i<nDM; i++) {
983:       PetscObjectReference((PetscObject)((*dmlist)[i]));
984:     }
985:   }
986:   return(0);
987: }



991: /* -------------------------------------------------------------------------------------*/
994: /*@C
995:     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
996:        Use DMCompositeRestoreLocalVectors() to return them.

998:     Not Collective

1000:     Input Parameter:
1001: .    dm - the packer object

1003:     Output Parameter:
1004: .   Vec ... - the individual sequential Vecs

1006:     Level: advanced

1008: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
1009:          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1010:          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()

1012: @*/
1013: PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
1014: {
1015:   va_list                Argp;
1016:   PetscErrorCode         ierr;
1017:   struct DMCompositeLink *next;
1018:   DM_Composite           *com = (DM_Composite*)dm->data;

1022:   next = com->next;
1023:   /* loop over packed objects, handling one at at time */
1024:   va_start(Argp,dm);
1025:   while (next) {
1026:     Vec *vec;
1027:     vec = va_arg(Argp, Vec*);
1028:     if (vec) {DMGetLocalVector(next->dm,vec);}
1029:     next = next->next;
1030:   }
1031:   va_end(Argp);
1032:   return(0);
1033: }

1037: /*@C
1038:     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.

1040:     Not Collective

1042:     Input Parameter:
1043: .    dm - the packer object

1045:     Output Parameter:
1046: .   Vec ... - the individual sequential Vecs

1048:     Level: advanced

1050: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
1051:          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1052:          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()

1054: @*/
1055: PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
1056: {
1057:   va_list                Argp;
1058:   PetscErrorCode         ierr;
1059:   struct DMCompositeLink *next;
1060:   DM_Composite           *com = (DM_Composite*)dm->data;

1064:   next = com->next;
1065:   /* loop over packed objects, handling one at at time */
1066:   va_start(Argp,dm);
1067:   while (next) {
1068:     Vec *vec;
1069:     vec = va_arg(Argp, Vec*);
1070:     if (vec) {DMRestoreLocalVector(next->dm,vec);}
1071:     next = next->next;
1072:   }
1073:   va_end(Argp);
1074:   return(0);
1075: }

1077: /* -------------------------------------------------------------------------------------*/
1080: /*@C
1081:     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.

1083:     Not Collective

1085:     Input Parameter:
1086: .    dm - the packer object

1088:     Output Parameter:
1089: .   DM ... - the individual entries (DMs)

1091:     Level: advanced

1093: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
1094:          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1095:          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
1096:          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()

1098: @*/
1099: PetscErrorCode  DMCompositeGetEntries(DM dm,...)
1100: {
1101:   va_list                Argp;
1102:   struct DMCompositeLink *next;
1103:   DM_Composite           *com = (DM_Composite*)dm->data;

1107:   next = com->next;
1108:   /* loop over packed objects, handling one at at time */
1109:   va_start(Argp,dm);
1110:   while (next) {
1111:     DM *dmn;
1112:     dmn = va_arg(Argp, DM*);
1113:     if (dmn) *dmn = next->dm;
1114:     next = next->next;
1115:   }
1116:   va_end(Argp);
1117:   return(0);
1118: }

1122: /*@C
1123:     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.

1125:     Not Collective

1127:     Input Parameter:
1128: +    dm - the packer object
1129: -    dms - array of sufficient length (see DMCompositeGetNumberDM()), holds the DMs on output

1131:     Level: advanced

1133: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
1134:          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1135:          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
1136:          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()

1138: @*/
1139: PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
1140: {
1141:   struct DMCompositeLink *next;
1142:   DM_Composite           *com = (DM_Composite*)dm->data;
1143:   PetscInt               i;

1147:   /* loop over packed objects, handling one at at time */
1148:   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
1149:   return(0);
1150: }

1154: PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
1155: {
1156:   PetscErrorCode         ierr;
1157:   struct DMCompositeLink *next;
1158:   DM_Composite           *com = (DM_Composite*)dmi->data;
1159:   DM                     dm;

1163:   if (comm == MPI_COMM_NULL) {
1164:     PetscObjectGetComm((PetscObject)dmi,&comm);
1165:   }
1166:   DMSetUp(dmi);
1167:   next = com->next;
1168:   DMCompositeCreate(comm,fine);

1170:   /* loop over packed objects, handling one at at time */
1171:   while (next) {
1172:     DMRefine(next->dm,comm,&dm);
1173:     DMCompositeAddDM(*fine,dm);
1174:     PetscObjectDereference((PetscObject)dm);
1175:     next = next->next;
1176:   }
1177:   return(0);
1178: }

1182: PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
1183: {
1184:   PetscErrorCode         ierr;
1185:   struct DMCompositeLink *next;
1186:   DM_Composite           *com = (DM_Composite*)dmi->data;
1187:   DM                     dm;

1191:   DMSetUp(dmi);
1192:   if (comm == MPI_COMM_NULL) {
1193:     PetscObjectGetComm((PetscObject)dmi,&comm);
1194:   }
1195:   next = com->next;
1196:   DMCompositeCreate(comm,fine);

1198:   /* loop over packed objects, handling one at at time */
1199:   while (next) {
1200:     DMCoarsen(next->dm,comm,&dm);
1201:     DMCompositeAddDM(*fine,dm);
1202:     PetscObjectDereference((PetscObject)dm);
1203:     next = next->next;
1204:   }
1205:   return(0);
1206: }

1210: PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
1211: {
1212:   PetscErrorCode         ierr;
1213:   PetscInt               m,n,M,N,nDM,i;
1214:   struct DMCompositeLink *nextc;
1215:   struct DMCompositeLink *nextf;
1216:   Vec                    gcoarse,gfine,*vecs;
1217:   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
1218:   DM_Composite           *comfine   = (DM_Composite*)fine->data;
1219:   Mat                    *mats;

1224:   DMSetUp(coarse);
1225:   DMSetUp(fine);
1226:   /* use global vectors only for determining matrix layout */
1227:   DMGetGlobalVector(coarse,&gcoarse);
1228:   DMGetGlobalVector(fine,&gfine);
1229:   VecGetLocalSize(gcoarse,&n);
1230:   VecGetLocalSize(gfine,&m);
1231:   VecGetSize(gcoarse,&N);
1232:   VecGetSize(gfine,&M);
1233:   DMRestoreGlobalVector(coarse,&gcoarse);
1234:   DMRestoreGlobalVector(fine,&gfine);

1236:   nDM = comfine->nDM;
1237:   if (nDM != comcoarse->nDM) SETERRQ2(PetscObjectComm((PetscObject)fine),PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %D entries, but coarse has %D",nDM,comcoarse->nDM);
1238:   PetscCalloc1(nDM*nDM,&mats);
1239:   if (v) {
1240:     PetscCalloc1(nDM,&vecs);
1241:   }

1243:   /* loop over packed objects, handling one at at time */
1244:   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
1245:     if (!v) {
1246:       DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL);
1247:     } else {
1248:       DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);
1249:     }
1250:   }
1251:   MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A);
1252:   if (v) {
1253:     VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v);
1254:   }
1255:   for (i=0; i<nDM*nDM; i++) {MatDestroy(&mats[i]);}
1256:   PetscFree(mats);
1257:   if (v) {
1258:     for (i=0; i<nDM; i++) {VecDestroy(&vecs[i]);}
1259:     PetscFree(vecs);
1260:   }
1261:   return(0);
1262: }

1266: static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
1267: {
1268:   DM_Composite           *com = (DM_Composite*)dm->data;
1269:   ISLocalToGlobalMapping *ltogs;
1270:   PetscInt               i;
1271:   PetscErrorCode         ierr;

1274:   /* Set the ISLocalToGlobalMapping on the new matrix */
1275:   DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);
1276:   ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap);
1277:   for (i=0; i<com->nDM; i++) {ISLocalToGlobalMappingDestroy(&ltogs[i]);}
1278:   PetscFree(ltogs);
1279:   return(0);
1280: }


1285: PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring)
1286: {
1287:   PetscErrorCode  ierr;
1288:   PetscInt        n,i,cnt;
1289:   ISColoringValue *colors;
1290:   PetscBool       dense  = PETSC_FALSE;
1291:   ISColoringValue maxcol = 0;
1292:   DM_Composite    *com   = (DM_Composite*)dm->data;

1296:   if (ctype == IS_COLORING_GHOSTED) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported");
1297:   else if (ctype == IS_COLORING_GLOBAL) {
1298:     n = com->n;
1299:   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1300:   PetscMalloc1(n,&colors); /* freed in ISColoringDestroy() */

1302:   PetscOptionsGetBool(NULL,"-dmcomposite_dense_jacobian",&dense,NULL);
1303:   if (dense) {
1304:     for (i=0; i<n; i++) {
1305:       colors[i] = (ISColoringValue)(com->rstart + i);
1306:     }
1307:     maxcol = com->N;
1308:   } else {
1309:     struct DMCompositeLink *next = com->next;
1310:     PetscMPIInt            rank;

1312:     MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
1313:     cnt  = 0;
1314:     while (next) {
1315:       ISColoring lcoloring;

1317:       DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring);
1318:       for (i=0; i<lcoloring->N; i++) {
1319:         colors[cnt++] = maxcol + lcoloring->colors[i];
1320:       }
1321:       maxcol += lcoloring->n;
1322:       ISColoringDestroy(&lcoloring);
1323:       next    = next->next;
1324:     }
1325:   }
1326:   ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,coloring);
1327:   return(0);
1328: }

1332: PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1333: {
1334:   PetscErrorCode         ierr;
1335:   struct DMCompositeLink *next;
1336:   PetscInt               cnt = 3;
1337:   PetscMPIInt            rank;
1338:   PetscScalar            *garray,*larray;
1339:   DM_Composite           *com = (DM_Composite*)dm->data;

1344:   next = com->next;
1345:   if (!com->setup) {
1346:     DMSetUp(dm);
1347:   }
1348:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
1349:   VecGetArray(gvec,&garray);
1350:   VecGetArray(lvec,&larray);

1352:   /* loop over packed objects, handling one at at time */
1353:   while (next) {
1354:     Vec      local,global;
1355:     PetscInt N;

1357:     DMGetGlobalVector(next->dm,&global);
1358:     VecGetLocalSize(global,&N);
1359:     VecPlaceArray(global,garray);
1360:     DMGetLocalVector(next->dm,&local);
1361:     VecPlaceArray(local,larray);
1362:     DMGlobalToLocalBegin(next->dm,global,mode,local);
1363:     DMGlobalToLocalEnd(next->dm,global,mode,local);
1364:     VecResetArray(global);
1365:     VecResetArray(local);
1366:     DMRestoreGlobalVector(next->dm,&global);
1367:     DMRestoreGlobalVector(next->dm,&local);
1368:     cnt++;
1369:     larray += next->nlocal;
1370:     next    = next->next;
1371:   }

1373:   VecRestoreArray(gvec,NULL);
1374:   VecRestoreArray(lvec,NULL);
1375:   return(0);
1376: }

1380: PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1381: {
1383:   return(0);
1384: }

1386: /*MC
1387:    DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs

1389:   Level: intermediate

1391: .seealso: DMType, DMCOMPOSITE, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate()
1392: M*/


1397: PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1398: {
1400:   DM_Composite   *com;

1403:   PetscNewLog(p,&com);
1404:   p->data   = com;
1405:   PetscObjectChangeTypeName((PetscObject)p,"DMComposite");
1406:   com->n    = 0;
1407:   com->next = NULL;
1408:   com->nDM  = 0;

1410:   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1411:   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1412:   p->ops->getlocaltoglobalmapping         = DMGetLocalToGlobalMapping_Composite;
1413:   p->ops->createfieldis                   = DMCreateFieldIS_Composite;
1414:   p->ops->createfielddecomposition        = DMCreateFieldDecomposition_Composite;
1415:   p->ops->refine                          = DMRefine_Composite;
1416:   p->ops->coarsen                         = DMCoarsen_Composite;
1417:   p->ops->createinterpolation             = DMCreateInterpolation_Composite;
1418:   p->ops->creatematrix                    = DMCreateMatrix_Composite;
1419:   p->ops->getcoloring                     = DMCreateColoring_Composite;
1420:   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1421:   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1422:   p->ops->destroy                         = DMDestroy_Composite;
1423:   p->ops->view                            = DMView_Composite;
1424:   p->ops->setup                           = DMSetUp_Composite;
1425:   return(0);
1426: }

1430: /*@C
1431:     DMCompositeCreate - Creates a vector packer, used to generate "composite"
1432:       vectors made up of several subvectors.

1434:     Collective on MPI_Comm

1436:     Input Parameter:
1437: .   comm - the processors that will share the global vector

1439:     Output Parameters:
1440: .   packer - the packer object

1442:     Level: advanced

1444: .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(),
1445:          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
1446:          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()

1448: @*/
1449: PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
1450: {

1455:   DMCreate(comm,packer);
1456:   DMSetType(*packer,DMCOMPOSITE);
1457:   return(0);
1458: }