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