xref: /petsc/src/dm/impls/composite/pack.c (revision b17ce1af89b448f824bb09659a9de6779534cd8e)
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 = PetscTypeCompare((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 = PetscTypeCompare((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 /* -------------------------------------------------------------------------------------*/
691 #undef __FUNCT__
692 #define __FUNCT__ "DMCompositeGetLocalVectors"
693 /*@C
694     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
695        Use DMCompositeRestoreLocalVectors() to return them.
696 
697     Not Collective
698 
699     Input Parameter:
700 .    dm - the packer object
701 
702     Output Parameter:
703 .   Vec ... - the individual sequential Vecs
704 
705     Level: advanced
706 
707 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
708          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
709          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
710 
711 @*/
712 PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
713 {
714   va_list                Argp;
715   PetscErrorCode         ierr;
716   struct DMCompositeLink *next;
717   DM_Composite           *com = (DM_Composite*)dm->data;
718 
719   PetscFunctionBegin;
720   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
721   next = com->next;
722   /* loop over packed objects, handling one at at time */
723   va_start(Argp,dm);
724   while (next) {
725     Vec *vec;
726     vec = va_arg(Argp, Vec*);
727     ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);
728     next = next->next;
729   }
730   va_end(Argp);
731   PetscFunctionReturn(0);
732 }
733 
734 #undef __FUNCT__
735 #define __FUNCT__ "DMCompositeRestoreLocalVectors"
736 /*@C
737     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
738 
739     Not Collective
740 
741     Input Parameter:
742 .    dm - the packer object
743 
744     Output Parameter:
745 .   Vec ... - the individual sequential Vecs
746 
747     Level: advanced
748 
749 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
750          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
751          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
752 
753 @*/
754 PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
755 {
756   va_list                Argp;
757   PetscErrorCode         ierr;
758   struct DMCompositeLink *next;
759   DM_Composite           *com = (DM_Composite*)dm->data;
760 
761   PetscFunctionBegin;
762   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
763   next = com->next;
764   /* loop over packed objects, handling one at at time */
765   va_start(Argp,dm);
766   while (next) {
767     Vec *vec;
768     vec = va_arg(Argp, Vec*);
769     ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);
770     next = next->next;
771   }
772   va_end(Argp);
773   PetscFunctionReturn(0);
774 }
775 
776 /* -------------------------------------------------------------------------------------*/
777 #undef __FUNCT__
778 #define __FUNCT__ "DMCompositeGetEntries"
779 /*@C
780     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
781 
782     Not Collective
783 
784     Input Parameter:
785 .    dm - the packer object
786 
787     Output Parameter:
788 .   DM ... - the individual entries (DMs)
789 
790     Level: advanced
791 
792 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
793          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
794          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
795          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
796 
797 @*/
798 PetscErrorCode  DMCompositeGetEntries(DM dm,...)
799 {
800   va_list                Argp;
801   struct DMCompositeLink *next;
802   DM_Composite           *com = (DM_Composite*)dm->data;
803 
804   PetscFunctionBegin;
805   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
806   next = com->next;
807   /* loop over packed objects, handling one at at time */
808   va_start(Argp,dm);
809   while (next) {
810     DM *dmn;
811     dmn = va_arg(Argp, DM*);
812     if (dmn) *dmn = next->dm;
813     next = next->next;
814   }
815   va_end(Argp);
816   PetscFunctionReturn(0);
817 }
818 
819 #undef __FUNCT__
820 #define __FUNCT__ "DMCompositeGetEntriesArray"
821 /*@
822     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
823 
824     Not Collective
825 
826     Input Parameter:
827 +    dm - the packer object
828 -    dms - array of sufficient length (see DMCompositeGetNumberDM()), holds the DMs on output
829 
830     Level: advanced
831 
832 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
833          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
834          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
835          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
836 
837 @*/
838 PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
839 {
840   struct DMCompositeLink *next;
841   DM_Composite           *com = (DM_Composite*)dm->data;
842   PetscInt               i;
843 
844   PetscFunctionBegin;
845   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
846   /* loop over packed objects, handling one at at time */
847   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
848   PetscFunctionReturn(0);
849 }
850 
851 #undef __FUNCT__
852 #define __FUNCT__ "DMRefine_Composite"
853 PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
854 {
855   PetscErrorCode         ierr;
856   struct DMCompositeLink *next;
857   DM_Composite           *com = (DM_Composite*)dmi->data;
858   DM                     dm;
859 
860   PetscFunctionBegin;
861   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
862   ierr = DMSetUp(dmi);CHKERRQ(ierr);
863   next = com->next;
864   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
865 
866   /* loop over packed objects, handling one at at time */
867   while (next) {
868     ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
869     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
870     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
871     next = next->next;
872   }
873   PetscFunctionReturn(0);
874 }
875 
876 #undef __FUNCT__
877 #define __FUNCT__ "DMCoarsen_Composite"
878 PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
879 {
880   PetscErrorCode         ierr;
881   struct DMCompositeLink *next;
882   DM_Composite           *com = (DM_Composite*)dmi->data;
883   DM                     dm;
884 
885   PetscFunctionBegin;
886   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
887   ierr = DMSetUp(dmi);CHKERRQ(ierr);
888   if (!comm) {
889     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
890   }
891   next = com->next;
892   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
893 
894   /* loop over packed objects, handling one at at time */
895   while (next) {
896     ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr);
897     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
898     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
899     next = next->next;
900   }
901   PetscFunctionReturn(0);
902 }
903 
904 #undef __FUNCT__
905 #define __FUNCT__ "DMCreateInterpolation_Composite"
906 PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
907 {
908   PetscErrorCode         ierr;
909   PetscInt               m,n,M,N,nDM,i;
910   struct DMCompositeLink *nextc;
911   struct DMCompositeLink *nextf;
912   Vec                    gcoarse,gfine,*vecs;
913   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
914   DM_Composite           *comfine = (DM_Composite*)fine->data;
915   Mat                    *mats;
916 
917   PetscFunctionBegin;
918   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
919   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
920   ierr = DMSetUp(coarse);CHKERRQ(ierr);
921   ierr = DMSetUp(fine);CHKERRQ(ierr);
922   /* use global vectors only for determining matrix layout */
923   ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
924   ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr);
925   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
926   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
927   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
928   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
929   ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
930   ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr);
931 
932   nDM = comfine->nDM;
933   if (nDM != comcoarse->nDM) SETERRQ2(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %D entries, but coarse has %D",nDM,comcoarse->nDM);
934   ierr = PetscMalloc(nDM*nDM*sizeof(Mat),&mats);CHKERRQ(ierr);
935   ierr = PetscMemzero(mats,nDM*nDM*sizeof(Mat));CHKERRQ(ierr);
936   if (v) {
937     ierr = PetscMalloc(nDM*sizeof(Vec),&vecs);CHKERRQ(ierr);
938     ierr = PetscMemzero(vecs,nDM*sizeof(Vec));CHKERRQ(ierr);
939   }
940 
941   /* loop over packed objects, handling one at at time */
942   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
943     if (!v) {
944       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],PETSC_NULL);CHKERRQ(ierr);
945     } else {
946       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr);
947     }
948   }
949   ierr = MatCreateNest(((PetscObject)fine)->comm,nDM,PETSC_NULL,nDM,PETSC_NULL,mats,A);CHKERRQ(ierr);
950   if (v) {
951     ierr = VecCreateNest(((PetscObject)fine)->comm,nDM,PETSC_NULL,vecs,v);CHKERRQ(ierr);
952   }
953   for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);}
954   ierr = PetscFree(mats);CHKERRQ(ierr);
955   if (v) {
956     for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);}
957     ierr = PetscFree(vecs);CHKERRQ(ierr);
958   }
959   PetscFunctionReturn(0);
960 }
961 
962 #undef __FUNCT__
963 #define __FUNCT__ "DMCreateLocalToGlobalMapping_Composite"
964 static PetscErrorCode DMCreateLocalToGlobalMapping_Composite(DM dm)
965 {
966   DM_Composite           *com = (DM_Composite*)dm->data;
967   ISLocalToGlobalMapping *ltogs;
968   PetscInt               i;
969   PetscErrorCode         ierr;
970 
971   PetscFunctionBegin;
972   /* Set the ISLocalToGlobalMapping on the new matrix */
973   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
974   ierr = ISLocalToGlobalMappingConcatenate(((PetscObject)dm)->comm,com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr);
975   for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
976   ierr = PetscFree(ltogs);CHKERRQ(ierr);
977   PetscFunctionReturn(0);
978 }
979 
980 
981 #undef __FUNCT__
982 #define __FUNCT__ "DMCreateColoring_Composite"
983 PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
984 {
985   PetscErrorCode         ierr;
986   PetscInt               n,i,cnt;
987   ISColoringValue        *colors;
988   PetscBool              dense = PETSC_FALSE;
989   ISColoringValue        maxcol = 0;
990   DM_Composite           *com = (DM_Composite*)dm->data;
991 
992   PetscFunctionBegin;
993   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
994   if (ctype == IS_COLORING_GHOSTED) {
995     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" );
996   } else if (ctype == IS_COLORING_GLOBAL) {
997     n = com->n;
998   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
999   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
1000 
1001   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
1002   if (dense) {
1003     for (i=0; i<n; i++) {
1004       colors[i] = (ISColoringValue)(com->rstart + i);
1005     }
1006     maxcol = com->N;
1007   } else {
1008     struct DMCompositeLink *next = com->next;
1009     PetscMPIInt            rank;
1010 
1011     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1012     cnt  = 0;
1013     while (next) {
1014       ISColoring     lcoloring;
1015 
1016       ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
1017       for (i=0; i<lcoloring->N; i++) {
1018         colors[cnt++] = maxcol + lcoloring->colors[i];
1019       }
1020       maxcol += lcoloring->n;
1021       ierr = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
1022       next = next->next;
1023     }
1024   }
1025   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
1026   PetscFunctionReturn(0);
1027 }
1028 
1029 #undef __FUNCT__
1030 #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
1031 PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1032 {
1033   PetscErrorCode         ierr;
1034   struct DMCompositeLink *next;
1035   PetscInt               cnt = 3;
1036   PetscMPIInt            rank;
1037   PetscScalar            *garray,*larray;
1038   DM_Composite           *com = (DM_Composite*)dm->data;
1039 
1040   PetscFunctionBegin;
1041   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1042   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
1043   next = com->next;
1044   if (!com->setup) {
1045     ierr = DMSetUp(dm);CHKERRQ(ierr);
1046   }
1047   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1048   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
1049   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
1050 
1051   /* loop over packed objects, handling one at at time */
1052   while (next) {
1053     Vec      local,global;
1054     PetscInt N;
1055 
1056     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
1057     ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
1058     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
1059     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
1060     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
1061     ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
1062     ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
1063     ierr = VecResetArray(global);CHKERRQ(ierr);
1064     ierr = VecResetArray(local);CHKERRQ(ierr);
1065     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
1066     ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
1067     cnt++;
1068     larray += next->nlocal;
1069     next    = next->next;
1070   }
1071 
1072   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
1073   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
1074   PetscFunctionReturn(0);
1075 }
1076 
1077 #undef __FUNCT__
1078 #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
1079 PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1080 {
1081   PetscFunctionBegin;
1082   PetscFunctionReturn(0);
1083 }
1084 
1085 EXTERN_C_BEGIN
1086 #undef __FUNCT__
1087 #define __FUNCT__ "DMCreate_Composite"
1088 PetscErrorCode  DMCreate_Composite(DM p)
1089 {
1090   PetscErrorCode ierr;
1091   DM_Composite   *com;
1092 
1093   PetscFunctionBegin;
1094   ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1095   p->data = com;
1096   ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1097   com->n            = 0;
1098   com->next         = PETSC_NULL;
1099   com->nDM          = 0;
1100 
1101   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1102   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1103   p->ops->createlocaltoglobalmapping      = DMCreateLocalToGlobalMapping_Composite;
1104   p->ops->createlocaltoglobalmappingblock = 0;
1105   p->ops->refine                          = DMRefine_Composite;
1106   p->ops->coarsen                         = DMCoarsen_Composite;
1107   p->ops->createinterpolation                = DMCreateInterpolation_Composite;
1108   p->ops->creatematrix                       = DMCreateMatrix_Composite;
1109   p->ops->getcoloring                     = DMCreateColoring_Composite;
1110   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1111   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1112   p->ops->destroy                         = DMDestroy_Composite;
1113   p->ops->view                            = DMView_Composite;
1114   p->ops->setup                           = DMSetUp_Composite;
1115   PetscFunctionReturn(0);
1116 }
1117 EXTERN_C_END
1118 
1119 #undef __FUNCT__
1120 #define __FUNCT__ "DMCompositeCreate"
1121 /*@C
1122     DMCompositeCreate - Creates a vector packer, used to generate "composite"
1123       vectors made up of several subvectors.
1124 
1125     Collective on MPI_Comm
1126 
1127     Input Parameter:
1128 .   comm - the processors that will share the global vector
1129 
1130     Output Parameters:
1131 .   packer - the packer object
1132 
1133     Level: advanced
1134 
1135 .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(),
1136          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
1137          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
1138 
1139 @*/
1140 PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
1141 {
1142   PetscErrorCode ierr;
1143 
1144   PetscFunctionBegin;
1145   PetscValidPointer(packer,2);
1146   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1147   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
1148   PetscFunctionReturn(0);
1149 }
1150