Actual source code: pack.c
petsc-3.4.5 2014-06-29
2: #include <../src/dm/impls/composite/packimpl.h> /*I "petscdmcomposite.h" I*/
6: /*@C
7: DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
8: seperate components (DMs) in a DMto build the correct matrix nonzero structure.
11: Logically Collective on MPI_Comm
13: Input Parameter:
14: + dm - the composite object
15: - formcouplelocations - routine to set the nonzero locations in the matrix
17: Level: advanced
19: Notes: See DMSetApplicationContext() and DMGetApplicationContext() for how to get user information into
20: this routine
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;
28: com->FormCoupleLocations = FormCoupleLocations;
29: return(0);
30: }
34: PetscErrorCode DMDestroy_Composite(DM dm)
35: {
36: PetscErrorCode ierr;
37: struct DMCompositeLink *next, *prev;
38: DM_Composite *com = (DM_Composite*)dm->data;
41: next = com->next;
42: while (next) {
43: prev = next;
44: next = next->next;
45: DMDestroy(&prev->dm);
46: PetscFree(prev->grstarts);
47: PetscFree(prev);
48: }
49: /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
50: PetscFree(com);
51: return(0);
52: }
56: PetscErrorCode DMView_Composite(DM dm,PetscViewer v)
57: {
59: PetscBool iascii;
60: DM_Composite *com = (DM_Composite*)dm->data;
63: PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);
64: if (iascii) {
65: struct DMCompositeLink *lnk = com->next;
66: PetscInt i;
68: PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix");
69: PetscViewerASCIIPrintf(v," contains %D DMs\n",com->nDM);
70: PetscViewerASCIIPushTab(v);
71: for (i=0; lnk; lnk=lnk->next,i++) {
72: PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);
73: PetscViewerASCIIPushTab(v);
74: DMView(lnk->dm,v);
75: PetscViewerASCIIPopTab(v);
76: }
77: PetscViewerASCIIPopTab(v);
78: }
79: return(0);
80: }
82: /* --------------------------------------------------------------------------------------*/
85: PetscErrorCode DMSetUp_Composite(DM dm)
86: {
87: PetscErrorCode ierr;
88: PetscInt nprev = 0;
89: PetscMPIInt rank,size;
90: DM_Composite *com = (DM_Composite*)dm->data;
91: struct DMCompositeLink *next = com->next;
92: PetscLayout map;
95: if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
96: PetscLayoutCreate(PetscObjectComm((PetscObject)dm),&map);
97: PetscLayoutSetLocalSize(map,com->n);
98: PetscLayoutSetSize(map,PETSC_DETERMINE);
99: PetscLayoutSetBlockSize(map,1);
100: PetscLayoutSetUp(map);
101: PetscLayoutGetSize(map,&com->N);
102: PetscLayoutGetRange(map,&com->rstart,NULL);
103: PetscLayoutDestroy(&map);
105: /* now set the rstart for each linked vector */
106: MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
107: MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);
108: while (next) {
109: next->rstart = nprev;
110: nprev += next->n;
111: next->grstart = com->rstart + next->rstart;
112: PetscMalloc(size*sizeof(PetscInt),&next->grstarts);
113: MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,PetscObjectComm((PetscObject)dm));
114: next = next->next;
115: }
116: com->setup = PETSC_TRUE;
117: return(0);
118: }
120: /* ----------------------------------------------------------------------------------*/
124: /*@
125: DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
126: representation.
128: Not Collective
130: Input Parameter:
131: . dm - the packer object
133: Output Parameter:
134: . nDM - the number of DMs
136: Level: beginner
138: @*/
139: PetscErrorCode DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
140: {
141: DM_Composite *com = (DM_Composite*)dm->data;
145: *nDM = com->nDM;
146: return(0);
147: }
152: /*@C
153: DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
154: representation.
156: Collective on DMComposite
158: Input Parameters:
159: + dm - the packer object
160: - gvec - the global vector
162: Output Parameters:
163: . Vec* ... - the packed parallel vectors, NULL for those that are not needed
165: Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
167: Fortran Notes:
169: Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
170: or use the alternative interface DMCompositeGetAccessArray().
172: Level: advanced
174: .seealso: DMCompositeGetEntries(), DMCompositeScatter()
175: @*/
176: PetscErrorCode DMCompositeGetAccess(DM dm,Vec gvec,...)
177: {
178: va_list Argp;
179: PetscErrorCode ierr;
180: struct DMCompositeLink *next;
181: DM_Composite *com = (DM_Composite*)dm->data;
186: next = com->next;
187: if (!com->setup) {
188: DMSetUp(dm);
189: }
191: /* loop over packed objects, handling one at at time */
192: va_start(Argp,gvec);
193: while (next) {
194: Vec *vec;
195: vec = va_arg(Argp, Vec*);
196: if (vec) {
197: PetscScalar *array;
198: DMGetGlobalVector(next->dm,vec);
199: VecGetArray(gvec,&array);
200: VecPlaceArray(*vec,array+next->rstart);
201: VecRestoreArray(gvec,&array);
202: }
203: next = next->next;
204: }
205: va_end(Argp);
206: return(0);
207: }
211: /*@C
212: DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
213: representation.
215: Collective on DMComposite
217: Input Parameters:
218: + dm - the packer object
219: . pvec - packed vector
220: . nwanted - number of vectors wanted
221: - wanted - sorted array of vectors wanted, or NULL to get all vectors
223: Output Parameters:
224: . vecs - array of requested global vectors (must be allocated)
226: Notes: Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them
228: Level: advanced
230: .seealso: DMCompositeGetAccess(), DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
231: @*/
232: PetscErrorCode DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
233: {
234: PetscErrorCode ierr;
235: struct DMCompositeLink *link;
236: PetscInt i,wnum;
237: DM_Composite *com = (DM_Composite*)dm->data;
242: if (!com->setup) {
243: DMSetUp(dm);
244: }
246: for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
247: if (!wanted || i == wanted[wnum]) {
248: PetscScalar *array;
249: Vec v;
250: DMGetGlobalVector(link->dm,&v);
251: VecGetArray(pvec,&array);
252: VecPlaceArray(v,array+link->rstart);
253: VecRestoreArray(pvec,&array);
254: vecs[wnum++] = v;
255: }
256: }
257: return(0);
258: }
262: /*@C
263: DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
264: representation.
266: Collective on DMComposite
268: Input Parameters:
269: + dm - the packer object
270: . gvec - the global vector
271: - Vec* ... - the individual parallel vectors, NULL for those that are not needed
273: Level: advanced
275: .seealso DMCompositeAddDM(), DMCreateGlobalVector(),
276: DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
277: DMCompositeRestoreAccess(), DMCompositeGetAccess()
279: @*/
280: PetscErrorCode DMCompositeRestoreAccess(DM dm,Vec gvec,...)
281: {
282: va_list Argp;
283: PetscErrorCode ierr;
284: struct DMCompositeLink *next;
285: DM_Composite *com = (DM_Composite*)dm->data;
290: next = com->next;
291: if (!com->setup) {
292: DMSetUp(dm);
293: }
295: /* loop over packed objects, handling one at at time */
296: va_start(Argp,gvec);
297: while (next) {
298: Vec *vec;
299: vec = va_arg(Argp, Vec*);
300: if (vec) {
301: VecResetArray(*vec);
302: DMRestoreGlobalVector(next->dm,vec);
303: }
304: next = next->next;
305: }
306: va_end(Argp);
307: return(0);
308: }
312: /*@C
313: DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray()
315: Collective on DMComposite
317: Input Parameters:
318: + dm - the packer object
319: . pvec - packed vector
320: . nwanted - number of vectors wanted
321: . wanted - sorted array of vectors wanted, or NULL to get all vectors
322: - vecs - array of global vectors to return
324: Level: advanced
326: .seealso: DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), DMCompositeScatter(), DMCompositeGather()
327: @*/
328: PetscErrorCode DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
329: {
330: PetscErrorCode ierr;
331: struct DMCompositeLink *link;
332: PetscInt i,wnum;
333: DM_Composite *com = (DM_Composite*)dm->data;
338: if (!com->setup) {
339: DMSetUp(dm);
340: }
342: for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
343: if (!wanted || i == wanted[wnum]) {
344: VecResetArray(vecs[wnum]);
345: DMRestoreGlobalVector(link->dm,&vecs[wnum]);
346: wnum++;
347: }
348: }
349: return(0);
350: }
354: /*@C
355: DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
357: Collective on DMComposite
359: Input Parameters:
360: + dm - the packer object
361: . gvec - the global vector
362: - Vec ... - the individual sequential vectors, NULL for those that are not needed
364: Level: advanced
366: Notes:
367: DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is
368: accessible from Fortran.
370: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
371: DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
372: DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
373: DMCompositeScatterArray()
375: @*/
376: PetscErrorCode DMCompositeScatter(DM dm,Vec gvec,...)
377: {
378: va_list Argp;
379: PetscErrorCode ierr;
380: struct DMCompositeLink *next;
381: PetscInt cnt;
382: DM_Composite *com = (DM_Composite*)dm->data;
387: if (!com->setup) {
388: DMSetUp(dm);
389: }
391: /* loop over packed objects, handling one at at time */
392: va_start(Argp,gvec);
393: for (cnt=3,next=com->next; next; cnt++,next=next->next) {
394: Vec local;
395: local = va_arg(Argp, Vec);
396: if (local) {
397: Vec global;
398: PetscScalar *array;
400: DMGetGlobalVector(next->dm,&global);
401: VecGetArray(gvec,&array);
402: VecPlaceArray(global,array+next->rstart);
403: DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);
404: DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);
405: VecRestoreArray(gvec,&array);
406: VecResetArray(global);
407: DMRestoreGlobalVector(next->dm,&global);
408: }
409: }
410: va_end(Argp);
411: return(0);
412: }
416: /*@
417: DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
419: Collective on DMComposite
421: Input Parameters:
422: + dm - the packer object
423: . gvec - the global vector
424: . lvecs - array of local vectors, NULL for any that are not needed
426: Level: advanced
428: Note:
429: This is a non-variadic alternative to DMCompositeScatterArray()
431: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector()
432: DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
433: DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
435: @*/
436: PetscErrorCode DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs)
437: {
438: PetscErrorCode ierr;
439: struct DMCompositeLink *next;
440: PetscInt i;
441: DM_Composite *com = (DM_Composite*)dm->data;
446: if (!com->setup) {
447: DMSetUp(dm);
448: }
450: /* loop over packed objects, handling one at at time */
451: for (i=0,next=com->next; next; next=next->next,i++) {
452: if (lvecs[i]) {
453: Vec global;
454: PetscScalar *array;
456: DMGetGlobalVector(next->dm,&global);
457: VecGetArray(gvec,&array);
458: VecPlaceArray(global,array+next->rstart);
459: DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]);
460: DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]);
461: VecRestoreArray(gvec,&array);
462: VecResetArray(global);
463: DMRestoreGlobalVector(next->dm,&global);
464: }
465: }
466: return(0);
467: }
471: /*@C
472: DMCompositeGather - Gathers into a global packed vector from its individual local vectors
474: Collective on DMComposite
476: Input Parameter:
477: + dm - the packer object
478: . gvec - the global vector
479: - Vec ... - the individual sequential vectors, NULL for any that are not needed
481: Level: advanced
483: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
484: DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
485: DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
487: @*/
488: PetscErrorCode DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
489: {
490: va_list Argp;
491: PetscErrorCode ierr;
492: struct DMCompositeLink *next;
493: DM_Composite *com = (DM_Composite*)dm->data;
494: PetscInt cnt;
499: if (!com->setup) {
500: DMSetUp(dm);
501: }
503: /* loop over packed objects, handling one at at time */
504: va_start(Argp,imode);
505: for (cnt=3,next=com->next; next; cnt++,next=next->next) {
506: Vec local;
507: local = va_arg(Argp, Vec);
508: if (local) {
509: PetscScalar *array;
510: Vec global;
512: DMGetGlobalVector(next->dm,&global);
513: VecGetArray(gvec,&array);
514: VecPlaceArray(global,array+next->rstart);
515: DMLocalToGlobalBegin(next->dm,local,imode,global);
516: DMLocalToGlobalEnd(next->dm,local,imode,global);
517: VecRestoreArray(gvec,&array);
518: VecResetArray(global);
519: DMRestoreGlobalVector(next->dm,&global);
520: }
521: }
522: va_end(Argp);
523: return(0);
524: }
528: /*@
529: DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
531: Collective on DMComposite
533: Input Parameter:
534: + dm - the packer object
535: . gvec - the global vector
536: - lvecs - the individual sequential vectors, NULL for any that are not needed
538: Level: advanced
540: Notes:
541: This is a non-variadic alternative to DMCompositeGather().
543: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
544: DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
545: DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries(),
546: @*/
547: PetscErrorCode DMCompositeGatherArray(DM dm,Vec gvec,InsertMode imode,Vec *lvecs)
548: {
549: PetscErrorCode ierr;
550: struct DMCompositeLink *next;
551: DM_Composite *com = (DM_Composite*)dm->data;
552: PetscInt i;
557: if (!com->setup) {
558: DMSetUp(dm);
559: }
561: /* loop over packed objects, handling one at at time */
562: for (next=com->next,i=0; next; next=next->next,i++) {
563: if (lvecs[i]) {
564: PetscScalar *array;
565: Vec global;
567: DMGetGlobalVector(next->dm,&global);
568: VecGetArray(gvec,&array);
569: VecPlaceArray(global,array+next->rstart);
570: DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global);
571: DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global);
572: VecRestoreArray(gvec,&array);
573: VecResetArray(global);
574: DMRestoreGlobalVector(next->dm,&global);
575: }
576: }
577: return(0);
578: }
582: /*@C
583: DMCompositeAddDM - adds a DM vector to a DMComposite
585: Collective on DMComposite
587: Input Parameter:
588: + dm - the packer object
589: - dm - the DM object, if the DM is a da you will need to caste it with a (DM)
591: Level: advanced
593: .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
594: DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
595: DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
597: @*/
598: PetscErrorCode DMCompositeAddDM(DM dmc,DM dm)
599: {
600: PetscErrorCode ierr;
601: PetscInt n,nlocal;
602: struct DMCompositeLink *mine,*next;
603: Vec global,local;
604: DM_Composite *com = (DM_Composite*)dmc->data;
609: next = com->next;
610: if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
612: /* create new link */
613: PetscNew(struct DMCompositeLink,&mine);
614: PetscObjectReference((PetscObject)dm);
615: DMGetGlobalVector(dm,&global);
616: VecGetLocalSize(global,&n);
617: DMRestoreGlobalVector(dm,&global);
618: DMGetLocalVector(dm,&local);
619: VecGetSize(local,&nlocal);
620: DMRestoreLocalVector(dm,&local);
622: mine->n = n;
623: mine->nlocal = nlocal;
624: mine->dm = dm;
625: mine->next = NULL;
626: com->n += n;
628: /* add to end of list */
629: if (!next) com->next = mine;
630: else {
631: while (next->next) next = next->next;
632: next->next = mine;
633: }
634: com->nDM++;
635: com->nmine++;
636: return(0);
637: }
639: #include <petscdraw.h>
640: PETSC_EXTERN PetscErrorCode VecView_MPI(Vec,PetscViewer);
643: PetscErrorCode VecView_DMComposite(Vec gvec,PetscViewer viewer)
644: {
645: DM dm;
646: PetscErrorCode ierr;
647: struct DMCompositeLink *next;
648: PetscBool isdraw;
649: DM_Composite *com;
652: VecGetDM(gvec, &dm);
653: if (!dm) SETERRQ(PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
654: com = (DM_Composite*)dm->data;
655: next = com->next;
657: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
658: if (!isdraw) {
659: /* do I really want to call this? */
660: VecView_MPI(gvec,viewer);
661: } else {
662: PetscInt cnt = 0;
664: /* loop over packed objects, handling one at at time */
665: while (next) {
666: Vec vec;
667: PetscScalar *array;
668: PetscInt bs;
670: /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
671: DMGetGlobalVector(next->dm,&vec);
672: VecGetArray(gvec,&array);
673: VecPlaceArray(vec,array+next->rstart);
674: VecRestoreArray(gvec,&array);
675: VecView(vec,viewer);
676: VecGetBlockSize(vec,&bs);
677: VecResetArray(vec);
678: DMRestoreGlobalVector(next->dm,&vec);
679: PetscViewerDrawBaseAdd(viewer,bs);
680: cnt += bs;
681: next = next->next;
682: }
683: PetscViewerDrawBaseAdd(viewer,-cnt);
684: }
685: return(0);
686: }
690: PetscErrorCode DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
691: {
693: DM_Composite *com = (DM_Composite*)dm->data;
697: DMSetUp(dm);
698: VecCreateMPI(PetscObjectComm((PetscObject)dm),com->n,com->N,gvec);
699: VecSetDM(*gvec, dm);
700: VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);
701: return(0);
702: }
706: PetscErrorCode DMCreateLocalVector_Composite(DM dm,Vec *lvec)
707: {
709: DM_Composite *com = (DM_Composite*)dm->data;
713: if (!com->setup) {
714: DMSetUp(dm);
715: }
716: VecCreateSeq(PetscObjectComm((PetscObject)dm),com->nghost,lvec);
717: VecSetDM(*lvec, dm);
718: return(0);
719: }
723: /*@C
724: DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
726: Collective on DM
728: Input Parameter:
729: . dm - the packer object
731: Output Parameters:
732: . ltogs - the individual mappings for each packed vector. Note that this includes
733: all the ghost points that individual ghosted DMDA's may have.
735: Level: advanced
737: Notes:
738: Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
740: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
741: DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
742: DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
744: @*/
745: PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
746: {
747: PetscErrorCode ierr;
748: PetscInt i,*idx,n,cnt;
749: struct DMCompositeLink *next;
750: PetscMPIInt rank;
751: DM_Composite *com = (DM_Composite*)dm->data;
755: DMSetUp(dm);
756: PetscMalloc((com->nDM)*sizeof(ISLocalToGlobalMapping),ltogs);
757: next = com->next;
758: MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
760: /* loop over packed objects, handling one at at time */
761: cnt = 0;
762: while (next) {
763: ISLocalToGlobalMapping ltog;
764: PetscMPIInt size;
765: const PetscInt *suboff,*indices;
766: Vec global;
768: /* Get sub-DM global indices for each local dof */
769: DMGetLocalToGlobalMapping(next->dm,<og);
770: ISLocalToGlobalMappingGetSize(ltog,&n);
771: ISLocalToGlobalMappingGetIndices(ltog,&indices);
772: PetscMalloc(n*sizeof(PetscInt),&idx);
774: /* Get the offsets for the sub-DM global vector */
775: DMGetGlobalVector(next->dm,&global);
776: VecGetOwnershipRanges(global,&suboff);
777: MPI_Comm_size(PetscObjectComm((PetscObject)global),&size);
779: /* Shift the sub-DM definition of the global space to the composite global space */
780: for (i=0; i<n; i++) {
781: PetscInt subi = indices[i],lo = 0,hi = size,t;
782: /* Binary search to find which rank owns subi */
783: while (hi-lo > 1) {
784: t = lo + (hi-lo)/2;
785: if (suboff[t] > subi) hi = t;
786: else lo = t;
787: }
788: idx[i] = subi - suboff[lo] + next->grstarts[lo];
789: }
790: ISLocalToGlobalMappingRestoreIndices(ltog,&indices);
791: ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);
792: DMRestoreGlobalVector(next->dm,&global);
793: next = next->next;
794: cnt++;
795: }
796: return(0);
797: }
801: /*@C
802: DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
804: Not Collective
806: Input Arguments:
807: . dm - composite DM
809: Output Arguments:
810: . is - array of serial index sets for each each component of the DMComposite
812: Level: intermediate
814: Notes:
815: At present, a composite local vector does not normally exist. This function is used to provide index sets for
816: MatGetLocalSubMatrix(). In the future, the scatters for each entry in the DMComposite may be be merged into a single
817: scatter to a composite local vector. The user should not typically need to know which is being done.
819: To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
821: To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
823: Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
825: .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
826: @*/
827: PetscErrorCode DMCompositeGetLocalISs(DM dm,IS **is)
828: {
829: PetscErrorCode ierr;
830: DM_Composite *com = (DM_Composite*)dm->data;
831: struct DMCompositeLink *link;
832: PetscInt cnt,start;
837: PetscMalloc(com->nmine*sizeof(IS),is);
838: for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
839: PetscInt bs;
840: ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);
841: DMGetBlockSize(link->dm,&bs);
842: ISSetBlockSize((*is)[cnt],bs);
843: }
844: return(0);
845: }
849: /*@C
850: DMCompositeGetGlobalISs - Gets the index sets for each composed object
852: Collective on DMComposite
854: Input Parameter:
855: . dm - the packer object
857: Output Parameters:
858: . is - the array of index sets
860: Level: advanced
862: Notes:
863: The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
865: These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
867: Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
868: DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
869: indices.
871: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
872: DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
873: DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
875: @*/
877: PetscErrorCode DMCompositeGetGlobalISs(DM dm,IS *is[])
878: {
879: PetscErrorCode ierr;
880: PetscInt cnt = 0,*idx,i;
881: struct DMCompositeLink *next;
882: PetscMPIInt rank;
883: DM_Composite *com = (DM_Composite*)dm->data;
887: PetscMalloc((com->nDM)*sizeof(IS),is);
888: next = com->next;
889: MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
891: /* loop over packed objects, handling one at at time */
892: while (next) {
893: PetscMalloc(next->n*sizeof(PetscInt),&idx);
894: for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
895: ISCreateGeneral(PetscObjectComm((PetscObject)dm),next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);
896: cnt++;
897: next = next->next;
898: }
899: return(0);
900: }
904: PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields)
905: {
906: PetscInt nDM;
907: DM *dms;
908: PetscInt i;
912: DMCompositeGetNumberDM(dm, &nDM);
913: if (numFields) *numFields = nDM;
914: DMCompositeGetGlobalISs(dm, fields);
915: if (fieldNames) {
916: PetscMalloc(nDM*sizeof(DM), &dms);
917: PetscMalloc(nDM*sizeof(const char*), fieldNames);
918: DMCompositeGetEntriesArray(dm, dms);
919: for (i=0; i<nDM; i++) {
920: char buf[256];
921: const char *splitname;
923: /* Split naming precedence: object name, prefix, number */
924: splitname = ((PetscObject) dm)->name;
925: if (!splitname) {
926: PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);
927: if (splitname) {
928: size_t len;
929: PetscStrncpy(buf,splitname,sizeof(buf));
930: buf[sizeof(buf) - 1] = 0;
931: PetscStrlen(buf,&len);
932: if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
933: splitname = buf;
934: }
935: }
936: if (!splitname) {
937: PetscSNPrintf(buf,sizeof(buf),"%D",i);
938: splitname = buf;
939: }
940: PetscStrallocpy(splitname,&(*fieldNames)[i]);
941: }
942: PetscFree(dms);
943: }
944: return(0);
945: }
947: /*
948: This could take over from DMCreateFieldIS(), as it is more general,
949: making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
950: At this point it's probably best to be less intrusive, however.
951: */
954: PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
955: {
956: PetscInt nDM;
957: PetscInt i;
961: DMCreateFieldIS_Composite(dm, len, namelist, islist);
962: if (dmlist) {
963: DMCompositeGetNumberDM(dm, &nDM);
964: PetscMalloc(nDM*sizeof(DM), dmlist);
965: DMCompositeGetEntriesArray(dm, *dmlist);
966: for (i=0; i<nDM; i++) {
967: PetscObjectReference((PetscObject)((*dmlist)[i]));
968: }
969: }
970: return(0);
971: }
975: /* -------------------------------------------------------------------------------------*/
978: /*@C
979: DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
980: Use DMCompositeRestoreLocalVectors() to return them.
982: Not Collective
984: Input Parameter:
985: . dm - the packer object
987: Output Parameter:
988: . Vec ... - the individual sequential Vecs
990: Level: advanced
992: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
993: DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
994: DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
996: @*/
997: PetscErrorCode DMCompositeGetLocalVectors(DM dm,...)
998: {
999: va_list Argp;
1000: PetscErrorCode ierr;
1001: struct DMCompositeLink *next;
1002: DM_Composite *com = (DM_Composite*)dm->data;
1006: next = com->next;
1007: /* loop over packed objects, handling one at at time */
1008: va_start(Argp,dm);
1009: while (next) {
1010: Vec *vec;
1011: vec = va_arg(Argp, Vec*);
1012: if (vec) {DMGetLocalVector(next->dm,vec);}
1013: next = next->next;
1014: }
1015: va_end(Argp);
1016: return(0);
1017: }
1021: /*@C
1022: DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
1024: Not Collective
1026: Input Parameter:
1027: . dm - the packer object
1029: Output Parameter:
1030: . Vec ... - the individual sequential Vecs
1032: Level: advanced
1034: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
1035: DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1036: DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
1038: @*/
1039: PetscErrorCode DMCompositeRestoreLocalVectors(DM dm,...)
1040: {
1041: va_list Argp;
1042: PetscErrorCode ierr;
1043: struct DMCompositeLink *next;
1044: DM_Composite *com = (DM_Composite*)dm->data;
1048: next = com->next;
1049: /* loop over packed objects, handling one at at time */
1050: va_start(Argp,dm);
1051: while (next) {
1052: Vec *vec;
1053: vec = va_arg(Argp, Vec*);
1054: if (vec) {DMRestoreLocalVector(next->dm,vec);}
1055: next = next->next;
1056: }
1057: va_end(Argp);
1058: return(0);
1059: }
1061: /* -------------------------------------------------------------------------------------*/
1064: /*@C
1065: DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
1067: Not Collective
1069: Input Parameter:
1070: . dm - the packer object
1072: Output Parameter:
1073: . DM ... - the individual entries (DMs)
1075: Level: advanced
1077: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
1078: DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1079: DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(), DMCompositeScatter(),
1080: DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
1082: @*/
1083: PetscErrorCode DMCompositeGetEntries(DM dm,...)
1084: {
1085: va_list Argp;
1086: struct DMCompositeLink *next;
1087: DM_Composite *com = (DM_Composite*)dm->data;
1091: next = com->next;
1092: /* loop over packed objects, handling one at at time */
1093: va_start(Argp,dm);
1094: while (next) {
1095: DM *dmn;
1096: dmn = va_arg(Argp, DM*);
1097: if (dmn) *dmn = next->dm;
1098: next = next->next;
1099: }
1100: va_end(Argp);
1101: return(0);
1102: }
1106: /*@C
1107: DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
1109: Not Collective
1111: Input Parameter:
1112: + dm - the packer object
1113: - dms - array of sufficient length (see DMCompositeGetNumberDM()), holds the DMs on output
1115: Level: advanced
1117: .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
1118: DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1119: DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(), DMCompositeScatter(),
1120: DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
1122: @*/
1123: PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
1124: {
1125: struct DMCompositeLink *next;
1126: DM_Composite *com = (DM_Composite*)dm->data;
1127: PetscInt i;
1131: /* loop over packed objects, handling one at at time */
1132: for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
1133: return(0);
1134: }
1138: PetscErrorCode DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
1139: {
1140: PetscErrorCode ierr;
1141: struct DMCompositeLink *next;
1142: DM_Composite *com = (DM_Composite*)dmi->data;
1143: DM dm;
1147: if (comm == MPI_COMM_NULL) {
1148: PetscObjectGetComm((PetscObject)dmi,&comm);
1149: }
1150: DMSetUp(dmi);
1151: next = com->next;
1152: DMCompositeCreate(comm,fine);
1154: /* loop over packed objects, handling one at at time */
1155: while (next) {
1156: DMRefine(next->dm,comm,&dm);
1157: DMCompositeAddDM(*fine,dm);
1158: PetscObjectDereference((PetscObject)dm);
1159: next = next->next;
1160: }
1161: return(0);
1162: }
1166: PetscErrorCode DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
1167: {
1168: PetscErrorCode ierr;
1169: struct DMCompositeLink *next;
1170: DM_Composite *com = (DM_Composite*)dmi->data;
1171: DM dm;
1175: DMSetUp(dmi);
1176: if (comm == MPI_COMM_NULL) {
1177: PetscObjectGetComm((PetscObject)dmi,&comm);
1178: }
1179: next = com->next;
1180: DMCompositeCreate(comm,fine);
1182: /* loop over packed objects, handling one at at time */
1183: while (next) {
1184: DMCoarsen(next->dm,comm,&dm);
1185: DMCompositeAddDM(*fine,dm);
1186: PetscObjectDereference((PetscObject)dm);
1187: next = next->next;
1188: }
1189: return(0);
1190: }
1194: PetscErrorCode DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
1195: {
1196: PetscErrorCode ierr;
1197: PetscInt m,n,M,N,nDM,i;
1198: struct DMCompositeLink *nextc;
1199: struct DMCompositeLink *nextf;
1200: Vec gcoarse,gfine,*vecs;
1201: DM_Composite *comcoarse = (DM_Composite*)coarse->data;
1202: DM_Composite *comfine = (DM_Composite*)fine->data;
1203: Mat *mats;
1208: DMSetUp(coarse);
1209: DMSetUp(fine);
1210: /* use global vectors only for determining matrix layout */
1211: DMGetGlobalVector(coarse,&gcoarse);
1212: DMGetGlobalVector(fine,&gfine);
1213: VecGetLocalSize(gcoarse,&n);
1214: VecGetLocalSize(gfine,&m);
1215: VecGetSize(gcoarse,&N);
1216: VecGetSize(gfine,&M);
1217: DMRestoreGlobalVector(coarse,&gcoarse);
1218: DMRestoreGlobalVector(fine,&gfine);
1220: nDM = comfine->nDM;
1221: if (nDM != comcoarse->nDM) SETERRQ2(PetscObjectComm((PetscObject)fine),PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %D entries, but coarse has %D",nDM,comcoarse->nDM);
1222: PetscMalloc(nDM*nDM*sizeof(Mat),&mats);
1223: PetscMemzero(mats,nDM*nDM*sizeof(Mat));
1224: if (v) {
1225: PetscMalloc(nDM*sizeof(Vec),&vecs);
1226: PetscMemzero(vecs,nDM*sizeof(Vec));
1227: }
1229: /* loop over packed objects, handling one at at time */
1230: for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
1231: if (!v) {
1232: DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL);
1233: } else {
1234: DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);
1235: }
1236: }
1237: MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A);
1238: if (v) {
1239: VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v);
1240: }
1241: for (i=0; i<nDM*nDM; i++) {MatDestroy(&mats[i]);}
1242: PetscFree(mats);
1243: if (v) {
1244: for (i=0; i<nDM; i++) {VecDestroy(&vecs[i]);}
1245: PetscFree(vecs);
1246: }
1247: return(0);
1248: }
1252: static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
1253: {
1254: DM_Composite *com = (DM_Composite*)dm->data;
1255: ISLocalToGlobalMapping *ltogs;
1256: PetscInt i;
1257: PetscErrorCode ierr;
1260: /* Set the ISLocalToGlobalMapping on the new matrix */
1261: DMCompositeGetISLocalToGlobalMappings(dm,<ogs);
1262: ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap);
1263: for (i=0; i<com->nDM; i++) {ISLocalToGlobalMappingDestroy(<ogs[i]);}
1264: PetscFree(ltogs);
1265: return(0);
1266: }
1271: PetscErrorCode DMCreateColoring_Composite(DM dm,ISColoringType ctype,MatType mtype,ISColoring *coloring)
1272: {
1273: PetscErrorCode ierr;
1274: PetscInt n,i,cnt;
1275: ISColoringValue *colors;
1276: PetscBool dense = PETSC_FALSE;
1277: ISColoringValue maxcol = 0;
1278: DM_Composite *com = (DM_Composite*)dm->data;
1282: if (ctype == IS_COLORING_GHOSTED) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported");
1283: else if (ctype == IS_COLORING_GLOBAL) {
1284: n = com->n;
1285: } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1286: PetscMalloc(n*sizeof(ISColoringValue),&colors); /* freed in ISColoringDestroy() */
1288: PetscOptionsGetBool(NULL,"-dmcomposite_dense_jacobian",&dense,NULL);
1289: if (dense) {
1290: for (i=0; i<n; i++) {
1291: colors[i] = (ISColoringValue)(com->rstart + i);
1292: }
1293: maxcol = com->N;
1294: } else {
1295: struct DMCompositeLink *next = com->next;
1296: PetscMPIInt rank;
1298: MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
1299: cnt = 0;
1300: while (next) {
1301: ISColoring lcoloring;
1303: DMCreateColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);
1304: for (i=0; i<lcoloring->N; i++) {
1305: colors[cnt++] = maxcol + lcoloring->colors[i];
1306: }
1307: maxcol += lcoloring->n;
1308: ISColoringDestroy(&lcoloring);
1309: next = next->next;
1310: }
1311: }
1312: ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,coloring);
1313: return(0);
1314: }
1318: PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1319: {
1320: PetscErrorCode ierr;
1321: struct DMCompositeLink *next;
1322: PetscInt cnt = 3;
1323: PetscMPIInt rank;
1324: PetscScalar *garray,*larray;
1325: DM_Composite *com = (DM_Composite*)dm->data;
1330: next = com->next;
1331: if (!com->setup) {
1332: DMSetUp(dm);
1333: }
1334: MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
1335: VecGetArray(gvec,&garray);
1336: VecGetArray(lvec,&larray);
1338: /* loop over packed objects, handling one at at time */
1339: while (next) {
1340: Vec local,global;
1341: PetscInt N;
1343: DMGetGlobalVector(next->dm,&global);
1344: VecGetLocalSize(global,&N);
1345: VecPlaceArray(global,garray);
1346: DMGetLocalVector(next->dm,&local);
1347: VecPlaceArray(local,larray);
1348: DMGlobalToLocalBegin(next->dm,global,mode,local);
1349: DMGlobalToLocalEnd(next->dm,global,mode,local);
1350: VecResetArray(global);
1351: VecResetArray(local);
1352: DMRestoreGlobalVector(next->dm,&global);
1353: DMRestoreGlobalVector(next->dm,&local);
1354: cnt++;
1355: larray += next->nlocal;
1356: next = next->next;
1357: }
1359: VecRestoreArray(gvec,NULL);
1360: VecRestoreArray(lvec,NULL);
1361: return(0);
1362: }
1366: PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1367: {
1369: return(0);
1370: }
1372: /*MC
1373: DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs
1377: Level: intermediate
1379: .seealso: DMType, DMCOMPOSITE, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate()
1380: M*/
1385: PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1386: {
1388: DM_Composite *com;
1391: PetscNewLog(p,DM_Composite,&com);
1392: p->data = com;
1393: PetscObjectChangeTypeName((PetscObject)p,"DMComposite");
1394: com->n = 0;
1395: com->next = NULL;
1396: com->nDM = 0;
1398: p->ops->createglobalvector = DMCreateGlobalVector_Composite;
1399: p->ops->createlocalvector = DMCreateLocalVector_Composite;
1400: p->ops->getlocaltoglobalmapping = DMGetLocalToGlobalMapping_Composite;
1401: p->ops->getlocaltoglobalmappingblock = 0;
1402: p->ops->createfieldis = DMCreateFieldIS_Composite;
1403: p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite;
1404: p->ops->refine = DMRefine_Composite;
1405: p->ops->coarsen = DMCoarsen_Composite;
1406: p->ops->createinterpolation = DMCreateInterpolation_Composite;
1407: p->ops->creatematrix = DMCreateMatrix_Composite;
1408: p->ops->getcoloring = DMCreateColoring_Composite;
1409: p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite;
1410: p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite;
1411: p->ops->destroy = DMDestroy_Composite;
1412: p->ops->view = DMView_Composite;
1413: p->ops->setup = DMSetUp_Composite;
1414: return(0);
1415: }
1419: /*@C
1420: DMCompositeCreate - Creates a vector packer, used to generate "composite"
1421: vectors made up of several subvectors.
1423: Collective on MPI_Comm
1425: Input Parameter:
1426: . comm - the processors that will share the global vector
1428: Output Parameters:
1429: . packer - the packer object
1431: Level: advanced
1433: .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(),
1434: DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
1435: DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
1437: @*/
1438: PetscErrorCode DMCompositeCreate(MPI_Comm comm,DM *packer)
1439: {
1444: DMCreate(comm,packer);
1445: DMSetType(*packer,DMCOMPOSITE);
1446: return(0);
1447: }