xref: /petsc/src/dm/impls/composite/pack.c (revision f692024ea6ceda132bc9ff30ca7a31e85cfbbcb2)
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 = PetscMalloc((com->nDM)*sizeof(ISLocalToGlobalMapping),ltogs);CHKERRQ(ierr);
544   next = com->next;
545   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
546 
547   /* loop over packed objects, handling one at at time */
548   cnt = 0;
549   while (next) {
550     ISLocalToGlobalMapping ltog;
551     PetscMPIInt            size;
552     const PetscInt         *suboff,*indices;
553     Vec                    global;
554 
555     /* Get sub-DM global indices for each local dof */
556     ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
557     ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
558     ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
559     ierr = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
560 
561     /* Get the offsets for the sub-DM global vector */
562     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
563     ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
564     ierr = MPI_Comm_size(((PetscObject)global)->comm,&size);CHKERRQ(ierr);
565 
566     /* Shift the sub-DM definition of the global space to the composite global space */
567     for (i=0; i<n; i++) {
568       PetscInt subi = indices[i],lo = 0,hi = size,t;
569       /* Binary search to find which rank owns subi */
570       while (hi-lo > 1) {
571         t = lo + (hi-lo)/2;
572         if (suboff[t] > subi) hi = t;
573         else                  lo = t;
574       }
575       idx[i] = subi - suboff[lo] + next->grstarts[lo];
576     }
577     ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
578     ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
579     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
580     next = next->next;
581     cnt++;
582   }
583   PetscFunctionReturn(0);
584 }
585 
586 #undef __FUNCT__
587 #define __FUNCT__ "DMCompositeGetLocalISs"
588 /*@C
589    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
590 
591    Not Collective
592 
593    Input Arguments:
594 . dm - composite DM
595 
596    Output Arguments:
597 . is - array of serial index sets for each each component of the DMComposite
598 
599    Level: intermediate
600 
601    Notes:
602    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
603    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
604    scatter to a composite local vector.  The user should not typically need to know which is being done.
605 
606    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
607 
608    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
609 
610    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
611 
612 .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
613 @*/
614 PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
615 {
616   PetscErrorCode         ierr;
617   DM_Composite           *com = (DM_Composite*)dm->data;
618   struct DMCompositeLink *link;
619   PetscInt cnt,start;
620 
621   PetscFunctionBegin;
622   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
623   PetscValidPointer(is,2);
624   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
625   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
626     PetscInt bs;
627     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
628     ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
629     ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
630   }
631   PetscFunctionReturn(0);
632 }
633 
634 #undef __FUNCT__
635 #define __FUNCT__ "DMCompositeGetGlobalISs"
636 /*@C
637     DMCompositeGetGlobalISs - Gets the index sets for each composed object
638 
639     Collective on DMComposite
640 
641     Input Parameter:
642 .    dm - the packer object
643 
644     Output Parameters:
645 .    is - the array of index sets
646 
647     Level: advanced
648 
649     Notes:
650        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
651 
652        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
653 
654        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
655        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
656        indices.
657 
658 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
659          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
660          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
661 
662 @*/
663 
664 PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
665 {
666   PetscErrorCode         ierr;
667   PetscInt               cnt = 0,*idx,i;
668   struct DMCompositeLink *next;
669   PetscMPIInt            rank;
670   DM_Composite           *com = (DM_Composite*)dm->data;
671 
672   PetscFunctionBegin;
673   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
674   ierr = PetscMalloc((com->nDM)*sizeof(IS),is);CHKERRQ(ierr);
675   next = com->next;
676   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
677 
678   /* loop over packed objects, handling one at at time */
679   while (next) {
680     ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
681     for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
682     ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
683     cnt++;
684     next = next->next;
685   }
686   PetscFunctionReturn(0);
687 }
688 
689 /* -------------------------------------------------------------------------------------*/
690 #undef __FUNCT__
691 #define __FUNCT__ "DMCompositeGetLocalVectors"
692 /*@C
693     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
694        Use DMCompositeRestoreLocalVectors() to return them.
695 
696     Not Collective
697 
698     Input Parameter:
699 .    dm - the packer object
700 
701     Output Parameter:
702 .   Vec ... - the individual sequential Vecs
703 
704     Level: advanced
705 
706 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
707          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
708          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
709 
710 @*/
711 PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
712 {
713   va_list                Argp;
714   PetscErrorCode         ierr;
715   struct DMCompositeLink *next;
716   DM_Composite           *com = (DM_Composite*)dm->data;
717 
718   PetscFunctionBegin;
719   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
720   next = com->next;
721   /* loop over packed objects, handling one at at time */
722   va_start(Argp,dm);
723   while (next) {
724     Vec *vec;
725     vec = va_arg(Argp, Vec*);
726     ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);
727     next = next->next;
728   }
729   va_end(Argp);
730   PetscFunctionReturn(0);
731 }
732 
733 #undef __FUNCT__
734 #define __FUNCT__ "DMCompositeRestoreLocalVectors"
735 /*@C
736     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
737 
738     Not Collective
739 
740     Input Parameter:
741 .    dm - the packer object
742 
743     Output Parameter:
744 .   Vec ... - the individual sequential Vecs
745 
746     Level: advanced
747 
748 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
749          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
750          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
751 
752 @*/
753 PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
754 {
755   va_list                Argp;
756   PetscErrorCode         ierr;
757   struct DMCompositeLink *next;
758   DM_Composite           *com = (DM_Composite*)dm->data;
759 
760   PetscFunctionBegin;
761   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
762   next = com->next;
763   /* loop over packed objects, handling one at at time */
764   va_start(Argp,dm);
765   while (next) {
766     Vec *vec;
767     vec = va_arg(Argp, Vec*);
768     ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);
769     next = next->next;
770   }
771   va_end(Argp);
772   PetscFunctionReturn(0);
773 }
774 
775 /* -------------------------------------------------------------------------------------*/
776 #undef __FUNCT__
777 #define __FUNCT__ "DMCompositeGetEntries"
778 /*@C
779     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
780 
781     Not Collective
782 
783     Input Parameter:
784 .    dm - the packer object
785 
786     Output Parameter:
787 .   DM ... - the individual entries (DMs)
788 
789     Level: advanced
790 
791 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
792          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
793          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
794          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
795 
796 @*/
797 PetscErrorCode  DMCompositeGetEntries(DM dm,...)
798 {
799   va_list                Argp;
800   struct DMCompositeLink *next;
801   DM_Composite           *com = (DM_Composite*)dm->data;
802 
803   PetscFunctionBegin;
804   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
805   next = com->next;
806   /* loop over packed objects, handling one at at time */
807   va_start(Argp,dm);
808   while (next) {
809     DM *dmn;
810     dmn = va_arg(Argp, DM*);
811     if (dmn) *dmn = next->dm;
812     next = next->next;
813   }
814   va_end(Argp);
815   PetscFunctionReturn(0);
816 }
817 
818 #undef __FUNCT__
819 #define __FUNCT__ "DMCompositeGetEntriesArray"
820 /*@
821     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
822 
823     Not Collective
824 
825     Input Parameter:
826 +    dm - the packer object
827 -    dms - array of sufficient length (see DMCompositeGetNumberDM()), holds the DMs on output
828 
829     Level: advanced
830 
831 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
832          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
833          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
834          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
835 
836 @*/
837 PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
838 {
839   struct DMCompositeLink *next;
840   DM_Composite           *com = (DM_Composite*)dm->data;
841   PetscInt               i;
842 
843   PetscFunctionBegin;
844   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
845   /* loop over packed objects, handling one at at time */
846   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
847   PetscFunctionReturn(0);
848 }
849 
850 #undef __FUNCT__
851 #define __FUNCT__ "DMRefine_Composite"
852 PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
853 {
854   PetscErrorCode         ierr;
855   struct DMCompositeLink *next;
856   DM_Composite           *com = (DM_Composite*)dmi->data;
857   DM                     dm;
858 
859   PetscFunctionBegin;
860   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
861   ierr = DMSetUp(dm);CHKERRQ(ierr);
862   next = com->next;
863   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
864 
865   /* loop over packed objects, handling one at at time */
866   while (next) {
867     ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
868     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
869     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
870     next = next->next;
871   }
872   PetscFunctionReturn(0);
873 }
874 
875 #undef __FUNCT__
876 #define __FUNCT__ "DMCoarsen_Composite"
877 PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
878 {
879   PetscErrorCode         ierr;
880   struct DMCompositeLink *next;
881   DM_Composite           *com = (DM_Composite*)dmi->data;
882   DM                     dm;
883 
884   PetscFunctionBegin;
885   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
886   if (!comm) {
887     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
888   }
889   next = com->next;
890   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
891 
892   /* loop over packed objects, handling one at at time */
893   while (next) {
894     ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr);
895     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
896     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
897     next = next->next;
898   }
899   PetscFunctionReturn(0);
900 }
901 
902 #undef __FUNCT__
903 #define __FUNCT__ "DMCreateInterpolation_Composite"
904 PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
905 {
906   PetscErrorCode         ierr;
907   PetscInt               m,n,M,N,nDM,i;
908   struct DMCompositeLink *nextc;
909   struct DMCompositeLink *nextf;
910   Vec                    gcoarse,gfine,*vecs;
911   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
912   DM_Composite           *comfine = (DM_Composite*)fine->data;
913   Mat                    *mats;
914 
915   PetscFunctionBegin;
916   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
917   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
918   ierr = DMSetUp(coarse);CHKERRQ(ierr);
919   ierr = DMSetUp(fine);CHKERRQ(ierr);
920   /* use global vectors only for determining matrix layout */
921   ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
922   ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr);
923   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
924   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
925   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
926   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
927   ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
928   ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr);
929 
930   nDM = comfine->nDM;
931   if (nDM != comcoarse->nDM) SETERRQ2(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %D entries, but coarse has %D",nDM,comcoarse->nDM);
932   ierr = PetscMalloc(nDM*nDM*sizeof(Mat),&mats);CHKERRQ(ierr);
933   ierr = PetscMemzero(mats,nDM*nDM*sizeof(Mat));CHKERRQ(ierr);
934   if (v) {
935     ierr = PetscMalloc(nDM*sizeof(Vec),&vecs);CHKERRQ(ierr);
936     ierr = PetscMemzero(vecs,nDM*sizeof(Vec));CHKERRQ(ierr);
937   }
938 
939   /* loop over packed objects, handling one at at time */
940   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
941     if (!v) {
942       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],PETSC_NULL);CHKERRQ(ierr);
943     } else {
944       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr);
945     }
946   }
947   ierr = MatCreateNest(((PetscObject)fine)->comm,nDM,PETSC_NULL,nDM,PETSC_NULL,mats,A);CHKERRQ(ierr);
948   if (v) {
949     ierr = VecCreateNest(((PetscObject)fine)->comm,nDM,PETSC_NULL,vecs,v);CHKERRQ(ierr);
950   }
951   for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);}
952   ierr = PetscFree(mats);CHKERRQ(ierr);
953   if (v) {
954     for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);}
955     ierr = PetscFree(vecs);CHKERRQ(ierr);
956   }
957   PetscFunctionReturn(0);
958 }
959 
960 #undef __FUNCT__
961 #define __FUNCT__ "DMCreateLocalToGlobalMapping_Composite"
962 static PetscErrorCode DMCreateLocalToGlobalMapping_Composite(DM dm)
963 {
964   DM_Composite           *com = (DM_Composite*)dm->data;
965   ISLocalToGlobalMapping *ltogs;
966   PetscInt               i;
967   PetscErrorCode         ierr;
968 
969   PetscFunctionBegin;
970   /* Set the ISLocalToGlobalMapping on the new matrix */
971   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
972   ierr = ISLocalToGlobalMappingConcatenate(((PetscObject)dm)->comm,com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr);
973   for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
974   ierr = PetscFree(ltogs);CHKERRQ(ierr);
975   PetscFunctionReturn(0);
976 }
977 
978 
979 #undef __FUNCT__
980 #define __FUNCT__ "DMCreateColoring_Composite"
981 PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
982 {
983   PetscErrorCode         ierr;
984   PetscInt               n,i,cnt;
985   ISColoringValue        *colors;
986   PetscBool              dense = PETSC_FALSE;
987   ISColoringValue        maxcol = 0;
988   DM_Composite           *com = (DM_Composite*)dm->data;
989 
990   PetscFunctionBegin;
991   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
992   if (ctype == IS_COLORING_GHOSTED) {
993     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" );
994   } else if (ctype == IS_COLORING_GLOBAL) {
995     n = com->n;
996   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
997   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
998 
999   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
1000   if (dense) {
1001     for (i=0; i<n; i++) {
1002       colors[i] = (ISColoringValue)(com->rstart + i);
1003     }
1004     maxcol = com->N;
1005   } else {
1006     struct DMCompositeLink *next = com->next;
1007     PetscMPIInt            rank;
1008 
1009     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1010     cnt  = 0;
1011     while (next) {
1012       ISColoring     lcoloring;
1013 
1014       ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
1015       for (i=0; i<lcoloring->N; i++) {
1016         colors[cnt++] = maxcol + lcoloring->colors[i];
1017       }
1018       maxcol += lcoloring->n;
1019       ierr = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
1020       next = next->next;
1021     }
1022   }
1023   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
1024   PetscFunctionReturn(0);
1025 }
1026 
1027 #undef __FUNCT__
1028 #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
1029 PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1030 {
1031   PetscErrorCode         ierr;
1032   struct DMCompositeLink *next;
1033   PetscInt               cnt = 3;
1034   PetscMPIInt            rank;
1035   PetscScalar            *garray,*larray;
1036   DM_Composite           *com = (DM_Composite*)dm->data;
1037 
1038   PetscFunctionBegin;
1039   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1040   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
1041   next = com->next;
1042   if (!com->setup) {
1043     ierr = DMSetUp(dm);CHKERRQ(ierr);
1044   }
1045   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1046   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
1047   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
1048 
1049   /* loop over packed objects, handling one at at time */
1050   while (next) {
1051     Vec      local,global;
1052     PetscInt N;
1053 
1054     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
1055     ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
1056     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
1057     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
1058     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
1059     ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
1060     ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
1061     ierr = VecResetArray(global);CHKERRQ(ierr);
1062     ierr = VecResetArray(local);CHKERRQ(ierr);
1063     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
1064     ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
1065     cnt++;
1066     larray += next->nlocal;
1067     next    = next->next;
1068   }
1069 
1070   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
1071   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
1072   PetscFunctionReturn(0);
1073 }
1074 
1075 #undef __FUNCT__
1076 #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
1077 PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1078 {
1079   PetscFunctionBegin;
1080   PetscFunctionReturn(0);
1081 }
1082 
1083 EXTERN_C_BEGIN
1084 #undef __FUNCT__
1085 #define __FUNCT__ "DMCreate_Composite"
1086 PetscErrorCode  DMCreate_Composite(DM p)
1087 {
1088   PetscErrorCode ierr;
1089   DM_Composite   *com;
1090 
1091   PetscFunctionBegin;
1092   ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1093   p->data = com;
1094   ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1095   com->n            = 0;
1096   com->next         = PETSC_NULL;
1097   com->nDM          = 0;
1098 
1099   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1100   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1101   p->ops->createlocaltoglobalmapping      = DMCreateLocalToGlobalMapping_Composite;
1102   p->ops->createlocaltoglobalmappingblock = 0;
1103   p->ops->refine                          = DMRefine_Composite;
1104   p->ops->coarsen                         = DMCoarsen_Composite;
1105   p->ops->createinterpolation                = DMCreateInterpolation_Composite;
1106   p->ops->creatematrix                       = DMCreateMatrix_Composite;
1107   p->ops->getcoloring                     = DMCreateColoring_Composite;
1108   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1109   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1110   p->ops->destroy                         = DMDestroy_Composite;
1111   p->ops->view                            = DMView_Composite;
1112   p->ops->setup                           = DMSetUp_Composite;
1113   PetscFunctionReturn(0);
1114 }
1115 EXTERN_C_END
1116 
1117 #undef __FUNCT__
1118 #define __FUNCT__ "DMCompositeCreate"
1119 /*@C
1120     DMCompositeCreate - Creates a vector packer, used to generate "composite"
1121       vectors made up of several subvectors.
1122 
1123     Collective on MPI_Comm
1124 
1125     Input Parameter:
1126 .   comm - the processors that will share the global vector
1127 
1128     Output Parameters:
1129 .   packer - the packer object
1130 
1131     Level: advanced
1132 
1133 .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(),
1134          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
1135          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
1136 
1137 @*/
1138 PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
1139 {
1140   PetscErrorCode ierr;
1141 
1142   PetscFunctionBegin;
1143   PetscValidPointer(packer,2);
1144   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1145   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
1146   PetscFunctionReturn(0);
1147 }
1148