xref: /petsc/src/dm/impls/composite/pack.c (revision 8865f1ea4ea560cd84ab8db62e98b7095cdff96f)
1 
2 #include <../src/dm/impls/composite/packimpl.h>       /*I  "petscdmcomposite.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 (DMs) 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 DMSetApplicationContext() and DMGetApplicationContext() 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__ "DMDestroy_Composite"
34 PetscErrorCode  DMDestroy_Composite(DM dm)
35 {
36   PetscErrorCode         ierr;
37   struct DMCompositeLink *next, *prev;
38   DM_Composite           *com = (DM_Composite*)dm->data;
39 
40   PetscFunctionBegin;
41   next = com->next;
42   while (next) {
43     prev = next;
44     next = next->next;
45     ierr = DMDestroy(&prev->dm);CHKERRQ(ierr);
46     ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
47     ierr = PetscFree(prev);CHKERRQ(ierr);
48   }
49   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
50   ierr = PetscFree(com);CHKERRQ(ierr);
51   PetscFunctionReturn(0);
52 }
53 
54 #undef __FUNCT__
55 #define __FUNCT__ "DMView_Composite"
56 PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
57 {
58   PetscErrorCode ierr;
59   PetscBool      iascii;
60   DM_Composite   *com = (DM_Composite*)dm->data;
61 
62   PetscFunctionBegin;
63   ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
64   if (iascii) {
65     struct DMCompositeLink *lnk = com->next;
66     PetscInt               i;
67 
68     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix");CHKERRQ(ierr);
69     ierr = PetscViewerASCIIPrintf(v,"  contains %D DMs\n",com->nDM);CHKERRQ(ierr);
70     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
71     for (i=0; lnk; lnk=lnk->next,i++) {
72       ierr = PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
73       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
74       ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
75       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
76     }
77     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
78   }
79   PetscFunctionReturn(0);
80 }
81 
82 /* --------------------------------------------------------------------------------------*/
83 #undef __FUNCT__
84 #define __FUNCT__ "DMSetUp_Composite"
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;
93 
94   PetscFunctionBegin;
95   if (com->setup) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
96   ierr = PetscLayoutCreate(((PetscObject)dm)->comm,&map);CHKERRQ(ierr);
97   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
98   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
99   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
100   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
101   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
102   ierr = PetscLayoutGetRange(map,&com->rstart,PETSC_NULL);CHKERRQ(ierr);
103   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
104 
105   /* now set the rstart for each linked vector */
106   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
107   ierr = MPI_Comm_size(((PetscObject)dm)->comm,&size);CHKERRQ(ierr);
108   while (next) {
109     next->rstart  = nprev;
110     nprev        += next->n;
111     next->grstart = com->rstart + next->rstart;
112     ierr          = PetscMalloc(size*sizeof(PetscInt),&next->grstarts);CHKERRQ(ierr);
113     ierr          = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
114     next          = next->next;
115   }
116   com->setup = PETSC_TRUE;
117   PetscFunctionReturn(0);
118 }
119 
120 /* ----------------------------------------------------------------------------------*/
121 
122 #include <stdarg.h>
123 
124 #undef __FUNCT__
125 #define __FUNCT__ "DMCompositeGetNumberDM"
126 /*@C
127     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
128        representation.
129 
130     Not Collective
131 
132     Input Parameter:
133 .    dm - the packer object
134 
135     Output Parameter:
136 .     nDM - the number of DMs
137 
138     Level: beginner
139 
140 @*/
141 PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
142 {
143   DM_Composite *com = (DM_Composite*)dm->data;
144 
145   PetscFunctionBegin;
146   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
147   *nDM = com->nDM;
148   PetscFunctionReturn(0);
149 }
150 
151 
152 #undef __FUNCT__
153 #define __FUNCT__ "DMCompositeGetAccess"
154 /*@C
155     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
156        representation.
157 
158     Collective on DMComposite
159 
160     Input Parameters:
161 +    dm - the packer object
162 -    gvec - the global vector
163 
164     Output Parameters:
165 .    Vec* ... - the packed parallel vectors, PETSC_NULL for those that are not needed
166 
167     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
168 
169     Level: advanced
170 
171 @*/
172 PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
173 {
174   va_list                Argp;
175   PetscErrorCode         ierr;
176   struct DMCompositeLink *next;
177   DM_Composite           *com = (DM_Composite*)dm->data;
178 
179   PetscFunctionBegin;
180   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
181   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
182   next = com->next;
183   if (!com->setup) {
184     ierr = DMSetUp(dm);CHKERRQ(ierr);
185   }
186 
187   /* loop over packed objects, handling one at at time */
188   va_start(Argp,gvec);
189   while (next) {
190     Vec *vec;
191     vec = va_arg(Argp, Vec*);
192     if (vec) {
193       PetscScalar *array;
194       ierr = DMGetGlobalVector(next->dm,vec);CHKERRQ(ierr);
195       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
196       ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
197       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
198     }
199     next = next->next;
200   }
201   va_end(Argp);
202   PetscFunctionReturn(0);
203 }
204 
205 #undef __FUNCT__
206 #define __FUNCT__ "DMCompositeRestoreAccess"
207 /*@C
208     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
209        representation.
210 
211     Collective on DMComposite
212 
213     Input Parameters:
214 +    dm - the packer object
215 .    gvec - the global vector
216 -    Vec* ... - the individual parallel vectors, PETSC_NULL for those that are not needed
217 
218     Level: advanced
219 
220 .seealso  DMCompositeAddDM(), DMCreateGlobalVector(),
221          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
222          DMCompositeRestoreAccess(), DMCompositeGetAccess()
223 
224 @*/
225 PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
226 {
227   va_list                Argp;
228   PetscErrorCode         ierr;
229   struct DMCompositeLink *next;
230   DM_Composite           *com = (DM_Composite*)dm->data;
231 
232   PetscFunctionBegin;
233   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
234   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
235   next = com->next;
236   if (!com->setup) {
237     ierr = DMSetUp(dm);CHKERRQ(ierr);
238   }
239 
240   /* loop over packed objects, handling one at at time */
241   va_start(Argp,gvec);
242   while (next) {
243     Vec *vec;
244     vec = va_arg(Argp, Vec*);
245     if (vec) {
246       ierr = VecResetArray(*vec);CHKERRQ(ierr);
247       ierr = DMRestoreGlobalVector(next->dm,vec);CHKERRQ(ierr);
248     }
249     next = next->next;
250   }
251   va_end(Argp);
252   PetscFunctionReturn(0);
253 }
254 
255 #undef __FUNCT__
256 #define __FUNCT__ "DMCompositeScatter"
257 /*@C
258     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
259 
260     Collective on DMComposite
261 
262     Input Parameters:
263 +    dm - the packer object
264 .    gvec - the global vector
265 -    Vec ... - the individual sequential vectors, PETSC_NULL for those that are not needed
266 
267     Level: advanced
268 
269 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
270          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
271          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
272 
273 @*/
274 PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
275 {
276   va_list                Argp;
277   PetscErrorCode         ierr;
278   struct DMCompositeLink *next;
279   PetscInt               cnt;
280   DM_Composite           *com = (DM_Composite*)dm->data;
281 
282   PetscFunctionBegin;
283   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
284   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
285   if (!com->setup) {
286     ierr = DMSetUp(dm);CHKERRQ(ierr);
287   }
288 
289   /* loop over packed objects, handling one at at time */
290   va_start(Argp,gvec);
291   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
292     Vec local;
293     local = va_arg(Argp, Vec);
294     if (local) {
295       Vec         global;
296       PetscScalar *array;
297       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
298       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
299       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
300       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
301       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
302       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
303       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
304       ierr = VecResetArray(global);CHKERRQ(ierr);
305       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
306     }
307   }
308   va_end(Argp);
309   PetscFunctionReturn(0);
310 }
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "DMCompositeGather"
314 /*@C
315     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
316 
317     Collective on DMComposite
318 
319     Input Parameter:
320 +    dm - the packer object
321 .    gvec - the global vector
322 -    Vec ... - the individual sequential vectors, PETSC_NULL for any that are not needed
323 
324     Level: advanced
325 
326 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
327          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
328          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
329 
330 @*/
331 PetscErrorCode  DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
332 {
333   va_list                Argp;
334   PetscErrorCode         ierr;
335   struct DMCompositeLink *next;
336   DM_Composite           *com = (DM_Composite*)dm->data;
337   PetscInt               cnt;
338 
339   PetscFunctionBegin;
340   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
341   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
342   if (!com->setup) {
343     ierr = DMSetUp(dm);CHKERRQ(ierr);
344   }
345 
346   /* loop over packed objects, handling one at at time */
347   va_start(Argp,imode);
348   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
349     Vec local;
350     local = va_arg(Argp, Vec);
351     if (local) {
352       PetscScalar *array;
353       Vec         global;
354       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
355       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
356       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
357       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
358       ierr = DMLocalToGlobalBegin(next->dm,local,imode,global);CHKERRQ(ierr);
359       ierr = DMLocalToGlobalEnd(next->dm,local,imode,global);CHKERRQ(ierr);
360       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
361       ierr = VecResetArray(global);CHKERRQ(ierr);
362       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
363     }
364   }
365   va_end(Argp);
366   PetscFunctionReturn(0);
367 }
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "DMCompositeAddDM"
371 /*@C
372     DMCompositeAddDM - adds a DM  vector to a DMComposite
373 
374     Collective on DMComposite
375 
376     Input Parameter:
377 +    dm - the packer object
378 -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
379 
380     Level: advanced
381 
382 .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
383          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
384          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
385 
386 @*/
387 PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
388 {
389   PetscErrorCode         ierr;
390   PetscInt               n,nlocal;
391   struct DMCompositeLink *mine,*next;
392   Vec                    global,local;
393   DM_Composite           *com = (DM_Composite*)dmc->data;
394 
395   PetscFunctionBegin;
396   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
397   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
398   next = com->next;
399   if (com->setup) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
400 
401   /* create new link */
402   ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
403   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
404   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
405   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
406   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
407   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
408   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
409   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
410 
411   mine->n      = n;
412   mine->nlocal = nlocal;
413   mine->dm     = dm;
414   mine->next   = PETSC_NULL;
415   com->n      += n;
416 
417   /* add to end of list */
418   if (!next) com->next = mine;
419   else {
420     while (next->next) next = next->next;
421     next->next = mine;
422   }
423   com->nDM++;
424   com->nmine++;
425   PetscFunctionReturn(0);
426 }
427 
428 extern PetscErrorCode  VecView_MPI(Vec,PetscViewer);
429 EXTERN_C_BEGIN
430 #undef __FUNCT__
431 #define __FUNCT__ "VecView_DMComposite"
432 PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
433 {
434   DM                     dm;
435   PetscErrorCode         ierr;
436   struct DMCompositeLink *next;
437   PetscBool              isdraw;
438   DM_Composite           *com;
439 
440   PetscFunctionBegin;
441   ierr = VecGetDM(gvec, &dm);CHKERRQ(ierr);
442   if (!dm) SETERRQ(((PetscObject)gvec)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
443   com  = (DM_Composite*)dm->data;
444   next = com->next;
445 
446   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
447   if (!isdraw) {
448     /* do I really want to call this? */
449     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
450   } else {
451     PetscInt cnt = 0;
452 
453     /* loop over packed objects, handling one at at time */
454     while (next) {
455       Vec         vec;
456       PetscScalar *array;
457       PetscInt    bs;
458 
459       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
460       ierr = DMGetGlobalVector(next->dm,&vec);CHKERRQ(ierr);
461       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
462       ierr = VecPlaceArray(vec,array+next->rstart);CHKERRQ(ierr);
463       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
464       ierr = VecView(vec,viewer);CHKERRQ(ierr);
465       ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
466       ierr = VecResetArray(vec);CHKERRQ(ierr);
467       ierr = DMRestoreGlobalVector(next->dm,&vec);CHKERRQ(ierr);
468       ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
469       cnt += bs;
470       next = next->next;
471     }
472     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
473   }
474   PetscFunctionReturn(0);
475 }
476 EXTERN_C_END
477 
478 
479 #undef __FUNCT__
480 #define __FUNCT__ "DMCreateGlobalVector_Composite"
481 PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
482 {
483   PetscErrorCode ierr;
484   DM_Composite   *com = (DM_Composite*)dm->data;
485 
486   PetscFunctionBegin;
487   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
488   ierr = DMSetUp(dm);CHKERRQ(ierr);
489   ierr = VecCreateMPI(((PetscObject)dm)->comm,com->n,com->N,gvec);CHKERRQ(ierr);
490   ierr = VecSetDM(*gvec, dm);CHKERRQ(ierr);
491   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);CHKERRQ(ierr);
492   PetscFunctionReturn(0);
493 }
494 
495 #undef __FUNCT__
496 #define __FUNCT__ "DMCreateLocalVector_Composite"
497 PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
498 {
499   PetscErrorCode ierr;
500   DM_Composite   *com = (DM_Composite*)dm->data;
501 
502   PetscFunctionBegin;
503   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
504   if (!com->setup) {
505     ierr = DMSetUp(dm);CHKERRQ(ierr);
506   }
507   ierr = VecCreateSeq(((PetscObject)dm)->comm,com->nghost,lvec);CHKERRQ(ierr);
508   ierr = VecSetDM(*lvec, dm);CHKERRQ(ierr);
509   PetscFunctionReturn(0);
510 }
511 
512 #undef __FUNCT__
513 #define __FUNCT__ "DMCompositeGetISLocalToGlobalMappings"
514 /*@C
515     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
516 
517     Collective on DM
518 
519     Input Parameter:
520 .    dm - the packer object
521 
522     Output Parameters:
523 .    ltogs - the individual mappings for each packed vector. Note that this includes
524            all the ghost points that individual ghosted DMDA's may have.
525 
526     Level: advanced
527 
528     Notes:
529        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
530 
531 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
532          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
533          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
534 
535 @*/
536 PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
537 {
538   PetscErrorCode         ierr;
539   PetscInt               i,*idx,n,cnt;
540   struct DMCompositeLink *next;
541   PetscMPIInt            rank;
542   DM_Composite           *com = (DM_Composite*)dm->data;
543 
544   PetscFunctionBegin;
545   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
546   ierr = DMSetUp(dm);CHKERRQ(ierr);
547   ierr = PetscMalloc((com->nDM)*sizeof(ISLocalToGlobalMapping),ltogs);CHKERRQ(ierr);
548   next = com->next;
549   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
550 
551   /* loop over packed objects, handling one at at time */
552   cnt = 0;
553   while (next) {
554     ISLocalToGlobalMapping ltog;
555     PetscMPIInt            size;
556     const PetscInt         *suboff,*indices;
557     Vec                    global;
558 
559     /* Get sub-DM global indices for each local dof */
560     ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
561     ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
562     ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
563     ierr = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
564 
565     /* Get the offsets for the sub-DM global vector */
566     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
567     ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
568     ierr = MPI_Comm_size(((PetscObject)global)->comm,&size);CHKERRQ(ierr);
569 
570     /* Shift the sub-DM definition of the global space to the composite global space */
571     for (i=0; i<n; i++) {
572       PetscInt subi = indices[i],lo = 0,hi = size,t;
573       /* Binary search to find which rank owns subi */
574       while (hi-lo > 1) {
575         t = lo + (hi-lo)/2;
576         if (suboff[t] > subi) hi = t;
577         else                  lo = t;
578       }
579       idx[i] = subi - suboff[lo] + next->grstarts[lo];
580     }
581     ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
582     ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
583     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
584     next = next->next;
585     cnt++;
586   }
587   PetscFunctionReturn(0);
588 }
589 
590 #undef __FUNCT__
591 #define __FUNCT__ "DMCompositeGetLocalISs"
592 /*@C
593    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
594 
595    Not Collective
596 
597    Input Arguments:
598 . dm - composite DM
599 
600    Output Arguments:
601 . is - array of serial index sets for each each component of the DMComposite
602 
603    Level: intermediate
604 
605    Notes:
606    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
607    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
608    scatter to a composite local vector.  The user should not typically need to know which is being done.
609 
610    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
611 
612    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
613 
614    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
615 
616 .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
617 @*/
618 PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
619 {
620   PetscErrorCode         ierr;
621   DM_Composite           *com = (DM_Composite*)dm->data;
622   struct DMCompositeLink *link;
623   PetscInt               cnt,start;
624 
625   PetscFunctionBegin;
626   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
627   PetscValidPointer(is,2);
628   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
629   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
630     PetscInt bs;
631     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
632     ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
633     ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
634   }
635   PetscFunctionReturn(0);
636 }
637 
638 #undef __FUNCT__
639 #define __FUNCT__ "DMCompositeGetGlobalISs"
640 /*@C
641     DMCompositeGetGlobalISs - Gets the index sets for each composed object
642 
643     Collective on DMComposite
644 
645     Input Parameter:
646 .    dm - the packer object
647 
648     Output Parameters:
649 .    is - the array of index sets
650 
651     Level: advanced
652 
653     Notes:
654        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
655 
656        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
657 
658        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
659        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
660        indices.
661 
662 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
663          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
664          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
665 
666 @*/
667 
668 PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
669 {
670   PetscErrorCode         ierr;
671   PetscInt               cnt = 0,*idx,i;
672   struct DMCompositeLink *next;
673   PetscMPIInt            rank;
674   DM_Composite           *com = (DM_Composite*)dm->data;
675 
676   PetscFunctionBegin;
677   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
678   ierr = PetscMalloc((com->nDM)*sizeof(IS),is);CHKERRQ(ierr);
679   next = com->next;
680   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
681 
682   /* loop over packed objects, handling one at at time */
683   while (next) {
684     ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
685     for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
686     ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
687     cnt++;
688     next = next->next;
689   }
690   PetscFunctionReturn(0);
691 }
692 
693 #undef __FUNCT__
694 #define __FUNCT__ "DMCreateFieldIS_Composite"
695 PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields)
696 {
697   PetscInt       nDM;
698   DM             *dms;
699   PetscInt       i;
700   PetscErrorCode ierr;
701 
702   PetscFunctionBegin;
703   ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
704   if (numFields) *numFields = nDM;
705   ierr = DMCompositeGetGlobalISs(dm, fields);CHKERRQ(ierr);
706   if (fieldNames) {
707     ierr = PetscMalloc(nDM*sizeof(DM), &dms);CHKERRQ(ierr);
708     ierr = PetscMalloc(nDM*sizeof(const char*), fieldNames);CHKERRQ(ierr);
709     ierr = DMCompositeGetEntriesArray(dm, dms);CHKERRQ(ierr);
710     for (i=0; i<nDM; i++) {
711       char       buf[256];
712       const char *splitname;
713 
714       /* Split naming precedence: object name, prefix, number */
715       splitname = ((PetscObject) dm)->name;
716       if (!splitname) {
717         ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr);
718         if (splitname) {
719           size_t len;
720           ierr                 = PetscStrncpy(buf,splitname,sizeof(buf));CHKERRQ(ierr);
721           buf[sizeof(buf) - 1] = 0;
722           ierr                 = PetscStrlen(buf,&len);CHKERRQ(ierr);
723           if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
724           splitname = buf;
725         }
726       }
727       if (!splitname) {
728         ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",i);CHKERRQ(ierr);
729         splitname = buf;
730       }
731       ierr = PetscStrallocpy(splitname,&(*fieldNames)[i]);CHKERRQ(ierr);
732     }
733     ierr = PetscFree(dms);CHKERRQ(ierr);
734   }
735   PetscFunctionReturn(0);
736 }
737 
738 /*
739  This could take over from DMCreateFieldIS(), as it is more general,
740  making DMCreateFieldIS() a special case -- calling with dmlist == PETSC_NULL;
741  At this point it's probably best to be less intrusive, however.
742  */
743 #undef __FUNCT__
744 #define __FUNCT__ "DMCreateFieldDecomposition_Composite"
745 PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
746 {
747   PetscInt       nDM;
748   PetscInt       i;
749   PetscErrorCode ierr;
750 
751   PetscFunctionBegin;
752   ierr = DMCreateFieldIS_Composite(dm, len, namelist, islist);CHKERRQ(ierr);
753   if (dmlist) {
754     ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
755     ierr = PetscMalloc(nDM*sizeof(DM), dmlist);CHKERRQ(ierr);
756     ierr = DMCompositeGetEntriesArray(dm, *dmlist);CHKERRQ(ierr);
757     for (i=0; i<nDM; i++) {
758       ierr = PetscObjectReference((PetscObject)((*dmlist)[i]));CHKERRQ(ierr);
759     }
760   }
761   PetscFunctionReturn(0);
762 }
763 
764 
765 
766 /* -------------------------------------------------------------------------------------*/
767 #undef __FUNCT__
768 #define __FUNCT__ "DMCompositeGetLocalVectors"
769 /*@C
770     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
771        Use DMCompositeRestoreLocalVectors() to return them.
772 
773     Not Collective
774 
775     Input Parameter:
776 .    dm - the packer object
777 
778     Output Parameter:
779 .   Vec ... - the individual sequential Vecs
780 
781     Level: advanced
782 
783 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
784          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
785          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
786 
787 @*/
788 PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
789 {
790   va_list                Argp;
791   PetscErrorCode         ierr;
792   struct DMCompositeLink *next;
793   DM_Composite           *com = (DM_Composite*)dm->data;
794 
795   PetscFunctionBegin;
796   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
797   next = com->next;
798   /* loop over packed objects, handling one at at time */
799   va_start(Argp,dm);
800   while (next) {
801     Vec *vec;
802     vec = va_arg(Argp, Vec*);
803     if (vec) {ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);}
804     next = next->next;
805   }
806   va_end(Argp);
807   PetscFunctionReturn(0);
808 }
809 
810 #undef __FUNCT__
811 #define __FUNCT__ "DMCompositeRestoreLocalVectors"
812 /*@C
813     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
814 
815     Not Collective
816 
817     Input Parameter:
818 .    dm - the packer object
819 
820     Output Parameter:
821 .   Vec ... - the individual sequential Vecs
822 
823     Level: advanced
824 
825 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
826          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
827          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
828 
829 @*/
830 PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
831 {
832   va_list                Argp;
833   PetscErrorCode         ierr;
834   struct DMCompositeLink *next;
835   DM_Composite           *com = (DM_Composite*)dm->data;
836 
837   PetscFunctionBegin;
838   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
839   next = com->next;
840   /* loop over packed objects, handling one at at time */
841   va_start(Argp,dm);
842   while (next) {
843     Vec *vec;
844     vec = va_arg(Argp, Vec*);
845     if (vec) {ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);}
846     next = next->next;
847   }
848   va_end(Argp);
849   PetscFunctionReturn(0);
850 }
851 
852 /* -------------------------------------------------------------------------------------*/
853 #undef __FUNCT__
854 #define __FUNCT__ "DMCompositeGetEntries"
855 /*@C
856     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
857 
858     Not Collective
859 
860     Input Parameter:
861 .    dm - the packer object
862 
863     Output Parameter:
864 .   DM ... - the individual entries (DMs)
865 
866     Level: advanced
867 
868 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
869          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
870          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
871          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
872 
873 @*/
874 PetscErrorCode  DMCompositeGetEntries(DM dm,...)
875 {
876   va_list                Argp;
877   struct DMCompositeLink *next;
878   DM_Composite           *com = (DM_Composite*)dm->data;
879 
880   PetscFunctionBegin;
881   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
882   next = com->next;
883   /* loop over packed objects, handling one at at time */
884   va_start(Argp,dm);
885   while (next) {
886     DM *dmn;
887     dmn = va_arg(Argp, DM*);
888     if (dmn) *dmn = next->dm;
889     next = next->next;
890   }
891   va_end(Argp);
892   PetscFunctionReturn(0);
893 }
894 
895 #undef __FUNCT__
896 #define __FUNCT__ "DMCompositeGetEntriesArray"
897 /*@C
898     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
899 
900     Not Collective
901 
902     Input Parameter:
903 +    dm - the packer object
904 -    dms - array of sufficient length (see DMCompositeGetNumberDM()), holds the DMs on output
905 
906     Level: advanced
907 
908 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
909          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
910          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
911          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
912 
913 @*/
914 PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
915 {
916   struct DMCompositeLink *next;
917   DM_Composite           *com = (DM_Composite*)dm->data;
918   PetscInt               i;
919 
920   PetscFunctionBegin;
921   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
922   /* loop over packed objects, handling one at at time */
923   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
924   PetscFunctionReturn(0);
925 }
926 
927 #undef __FUNCT__
928 #define __FUNCT__ "DMRefine_Composite"
929 PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
930 {
931   PetscErrorCode         ierr;
932   struct DMCompositeLink *next;
933   DM_Composite           *com = (DM_Composite*)dmi->data;
934   DM                     dm;
935 
936   PetscFunctionBegin;
937   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
938   if (comm == MPI_COMM_NULL) comm = ((PetscObject)dmi)->comm;
939   ierr = DMSetUp(dmi);CHKERRQ(ierr);
940   next = com->next;
941   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
942 
943   /* loop over packed objects, handling one at at time */
944   while (next) {
945     ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
946     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
947     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
948     next = next->next;
949   }
950   PetscFunctionReturn(0);
951 }
952 
953 #undef __FUNCT__
954 #define __FUNCT__ "DMCoarsen_Composite"
955 PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
956 {
957   PetscErrorCode         ierr;
958   struct DMCompositeLink *next;
959   DM_Composite           *com = (DM_Composite*)dmi->data;
960   DM                     dm;
961 
962   PetscFunctionBegin;
963   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
964   ierr = DMSetUp(dmi);CHKERRQ(ierr);
965   if (comm == MPI_COMM_NULL) {
966     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
967   }
968   next = com->next;
969   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
970 
971   /* loop over packed objects, handling one at at time */
972   while (next) {
973     ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr);
974     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
975     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
976     next = next->next;
977   }
978   PetscFunctionReturn(0);
979 }
980 
981 #undef __FUNCT__
982 #define __FUNCT__ "DMCreateInterpolation_Composite"
983 PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
984 {
985   PetscErrorCode         ierr;
986   PetscInt               m,n,M,N,nDM,i;
987   struct DMCompositeLink *nextc;
988   struct DMCompositeLink *nextf;
989   Vec                    gcoarse,gfine,*vecs;
990   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
991   DM_Composite           *comfine   = (DM_Composite*)fine->data;
992   Mat                    *mats;
993 
994   PetscFunctionBegin;
995   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
996   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
997   ierr = DMSetUp(coarse);CHKERRQ(ierr);
998   ierr = DMSetUp(fine);CHKERRQ(ierr);
999   /* use global vectors only for determining matrix layout */
1000   ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
1001   ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr);
1002   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
1003   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
1004   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
1005   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
1006   ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
1007   ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr);
1008 
1009   nDM = comfine->nDM;
1010   if (nDM != comcoarse->nDM) SETERRQ2(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %D entries, but coarse has %D",nDM,comcoarse->nDM);
1011   ierr = PetscMalloc(nDM*nDM*sizeof(Mat),&mats);CHKERRQ(ierr);
1012   ierr = PetscMemzero(mats,nDM*nDM*sizeof(Mat));CHKERRQ(ierr);
1013   if (v) {
1014     ierr = PetscMalloc(nDM*sizeof(Vec),&vecs);CHKERRQ(ierr);
1015     ierr = PetscMemzero(vecs,nDM*sizeof(Vec));CHKERRQ(ierr);
1016   }
1017 
1018   /* loop over packed objects, handling one at at time */
1019   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
1020     if (!v) {
1021       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],PETSC_NULL);CHKERRQ(ierr);
1022     } else {
1023       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr);
1024     }
1025   }
1026   ierr = MatCreateNest(((PetscObject)fine)->comm,nDM,PETSC_NULL,nDM,PETSC_NULL,mats,A);CHKERRQ(ierr);
1027   if (v) {
1028     ierr = VecCreateNest(((PetscObject)fine)->comm,nDM,PETSC_NULL,vecs,v);CHKERRQ(ierr);
1029   }
1030   for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);}
1031   ierr = PetscFree(mats);CHKERRQ(ierr);
1032   if (v) {
1033     for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);}
1034     ierr = PetscFree(vecs);CHKERRQ(ierr);
1035   }
1036   PetscFunctionReturn(0);
1037 }
1038 
1039 #undef __FUNCT__
1040 #define __FUNCT__ "DMCreateLocalToGlobalMapping_Composite"
1041 static PetscErrorCode DMCreateLocalToGlobalMapping_Composite(DM dm)
1042 {
1043   DM_Composite           *com = (DM_Composite*)dm->data;
1044   ISLocalToGlobalMapping *ltogs;
1045   PetscInt               i;
1046   PetscErrorCode         ierr;
1047 
1048   PetscFunctionBegin;
1049   /* Set the ISLocalToGlobalMapping on the new matrix */
1050   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
1051   ierr = ISLocalToGlobalMappingConcatenate(((PetscObject)dm)->comm,com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr);
1052   for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
1053   ierr = PetscFree(ltogs);CHKERRQ(ierr);
1054   PetscFunctionReturn(0);
1055 }
1056 
1057 
1058 #undef __FUNCT__
1059 #define __FUNCT__ "DMCreateColoring_Composite"
1060 PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,MatType mtype,ISColoring *coloring)
1061 {
1062   PetscErrorCode  ierr;
1063   PetscInt        n,i,cnt;
1064   ISColoringValue *colors;
1065   PetscBool       dense  = PETSC_FALSE;
1066   ISColoringValue maxcol = 0;
1067   DM_Composite    *com   = (DM_Composite*)dm->data;
1068 
1069   PetscFunctionBegin;
1070   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1071   if (ctype == IS_COLORING_GHOSTED) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Only global coloring supported");
1072   else if (ctype == IS_COLORING_GLOBAL) {
1073     n = com->n;
1074   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1075   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
1076 
1077   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
1078   if (dense) {
1079     for (i=0; i<n; i++) {
1080       colors[i] = (ISColoringValue)(com->rstart + i);
1081     }
1082     maxcol = com->N;
1083   } else {
1084     struct DMCompositeLink *next = com->next;
1085     PetscMPIInt            rank;
1086 
1087     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1088     cnt  = 0;
1089     while (next) {
1090       ISColoring lcoloring;
1091 
1092       ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
1093       for (i=0; i<lcoloring->N; i++) {
1094         colors[cnt++] = maxcol + lcoloring->colors[i];
1095       }
1096       maxcol += lcoloring->n;
1097       ierr    = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
1098       next    = next->next;
1099     }
1100   }
1101   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
1102   PetscFunctionReturn(0);
1103 }
1104 
1105 #undef __FUNCT__
1106 #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
1107 PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1108 {
1109   PetscErrorCode         ierr;
1110   struct DMCompositeLink *next;
1111   PetscInt               cnt = 3;
1112   PetscMPIInt            rank;
1113   PetscScalar            *garray,*larray;
1114   DM_Composite           *com = (DM_Composite*)dm->data;
1115 
1116   PetscFunctionBegin;
1117   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1118   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
1119   next = com->next;
1120   if (!com->setup) {
1121     ierr = DMSetUp(dm);CHKERRQ(ierr);
1122   }
1123   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1124   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
1125   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
1126 
1127   /* loop over packed objects, handling one at at time */
1128   while (next) {
1129     Vec      local,global;
1130     PetscInt N;
1131 
1132     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
1133     ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
1134     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
1135     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
1136     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
1137     ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
1138     ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
1139     ierr = VecResetArray(global);CHKERRQ(ierr);
1140     ierr = VecResetArray(local);CHKERRQ(ierr);
1141     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
1142     ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
1143     cnt++;
1144     larray += next->nlocal;
1145     next    = next->next;
1146   }
1147 
1148   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
1149   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
1150   PetscFunctionReturn(0);
1151 }
1152 
1153 #undef __FUNCT__
1154 #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
1155 PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1156 {
1157   PetscFunctionBegin;
1158   PetscFunctionReturn(0);
1159 }
1160 
1161 /*MC
1162    DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs
1163 
1164 
1165 
1166   Level: intermediate
1167 
1168 .seealso: DMType, DMCOMPOSITE, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate()
1169 M*/
1170 
1171 
1172 EXTERN_C_BEGIN
1173 #undef __FUNCT__
1174 #define __FUNCT__ "DMCreate_Composite"
1175 PetscErrorCode  DMCreate_Composite(DM p)
1176 {
1177   PetscErrorCode ierr;
1178   DM_Composite   *com;
1179 
1180   PetscFunctionBegin;
1181   ierr      = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1182   p->data   = com;
1183   ierr      = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1184   com->n    = 0;
1185   com->next = PETSC_NULL;
1186   com->nDM  = 0;
1187 
1188   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1189   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1190   p->ops->createlocaltoglobalmapping      = DMCreateLocalToGlobalMapping_Composite;
1191   p->ops->createlocaltoglobalmappingblock = 0;
1192   p->ops->createfieldis                   = DMCreateFieldIS_Composite;
1193   p->ops->createfielddecomposition        = DMCreateFieldDecomposition_Composite;
1194   p->ops->refine                          = DMRefine_Composite;
1195   p->ops->coarsen                         = DMCoarsen_Composite;
1196   p->ops->createinterpolation             = DMCreateInterpolation_Composite;
1197   p->ops->creatematrix                    = DMCreateMatrix_Composite;
1198   p->ops->getcoloring                     = DMCreateColoring_Composite;
1199   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1200   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1201   p->ops->destroy                         = DMDestroy_Composite;
1202   p->ops->view                            = DMView_Composite;
1203   p->ops->setup                           = DMSetUp_Composite;
1204   PetscFunctionReturn(0);
1205 }
1206 EXTERN_C_END
1207 
1208 #undef __FUNCT__
1209 #define __FUNCT__ "DMCompositeCreate"
1210 /*@C
1211     DMCompositeCreate - Creates a vector packer, used to generate "composite"
1212       vectors made up of several subvectors.
1213 
1214     Collective on MPI_Comm
1215 
1216     Input Parameter:
1217 .   comm - the processors that will share the global vector
1218 
1219     Output Parameters:
1220 .   packer - the packer object
1221 
1222     Level: advanced
1223 
1224 .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(),
1225          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
1226          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
1227 
1228 @*/
1229 PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
1230 {
1231   PetscErrorCode ierr;
1232 
1233   PetscFunctionBegin;
1234   PetscValidPointer(packer,2);
1235   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1236   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
1237   PetscFunctionReturn(0);
1238 }
1239