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