xref: /petsc/src/dm/impls/composite/pack.c (revision b41ce5d507ea9a58bfa83cf403107a702e77a67d)
1 
2 #include <../src/dm/impls/composite/packimpl.h>       /*I  "petscdmcomposite.h"  I*/
3 #include <petsc/private/isimpl.h>
4 #include <petscds.h>
5 
6 /*@C
7     DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
8       separate 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 PetscErrorCode  DMDestroy_Composite(DM dm)
33 {
34   PetscErrorCode         ierr;
35   struct DMCompositeLink *next, *prev;
36   DM_Composite           *com = (DM_Composite*)dm->data;
37 
38   PetscFunctionBegin;
39   next = com->next;
40   while (next) {
41     prev = next;
42     next = next->next;
43     ierr = DMDestroy(&prev->dm);CHKERRQ(ierr);
44     ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
45     ierr = PetscFree(prev);CHKERRQ(ierr);
46   }
47   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
48   ierr = PetscFree(com);CHKERRQ(ierr);
49   PetscFunctionReturn(0);
50 }
51 
52 PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
53 {
54   PetscErrorCode ierr;
55   PetscBool      iascii;
56   DM_Composite   *com = (DM_Composite*)dm->data;
57 
58   PetscFunctionBegin;
59   ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
60   if (iascii) {
61     struct DMCompositeLink *lnk = com->next;
62     PetscInt               i;
63 
64     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix");CHKERRQ(ierr);
65     ierr = PetscViewerASCIIPrintf(v,"  contains %D DMs\n",com->nDM);CHKERRQ(ierr);
66     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
67     for (i=0; lnk; lnk=lnk->next,i++) {
68       ierr = PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
69       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
70       ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
71       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
72     }
73     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
74   }
75   PetscFunctionReturn(0);
76 }
77 
78 /* --------------------------------------------------------------------------------------*/
79 PetscErrorCode  DMSetUp_Composite(DM dm)
80 {
81   PetscErrorCode         ierr;
82   PetscInt               nprev = 0;
83   PetscMPIInt            rank,size;
84   DM_Composite           *com  = (DM_Composite*)dm->data;
85   struct DMCompositeLink *next = com->next;
86   PetscLayout            map;
87 
88   PetscFunctionBegin;
89   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
90   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm),&map);CHKERRQ(ierr);
91   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
92   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
93   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
94   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
95   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
96   ierr = PetscLayoutGetRange(map,&com->rstart,NULL);CHKERRQ(ierr);
97   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
98 
99   /* now set the rstart for each linked vector */
100   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
101   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
102   while (next) {
103     next->rstart  = nprev;
104     nprev        += next->n;
105     next->grstart = com->rstart + next->rstart;
106     ierr          = PetscMalloc1(size,&next->grstarts);CHKERRQ(ierr);
107     ierr          = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
108     next          = next->next;
109   }
110   com->setup = PETSC_TRUE;
111   PetscFunctionReturn(0);
112 }
113 
114 /* ----------------------------------------------------------------------------------*/
115 
116 /*@
117     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
118        representation.
119 
120     Not Collective
121 
122     Input Parameter:
123 .    dm - the packer object
124 
125     Output Parameter:
126 .     nDM - the number of DMs
127 
128     Level: beginner
129 
130 @*/
131 PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
132 {
133   DM_Composite *com = (DM_Composite*)dm->data;
134 
135   PetscFunctionBegin;
136   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
137   *nDM = com->nDM;
138   PetscFunctionReturn(0);
139 }
140 
141 
142 /*@C
143     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
144        representation.
145 
146     Collective on DMComposite
147 
148     Input Parameters:
149 +    dm - the packer object
150 -    gvec - the global vector
151 
152     Output Parameters:
153 .    Vec* ... - the packed parallel vectors, NULL for those that are not needed
154 
155     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
156 
157     Fortran Notes:
158 
159     Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
160     or use the alternative interface DMCompositeGetAccessArray().
161 
162     Level: advanced
163 
164 .seealso: DMCompositeGetEntries(), DMCompositeScatter()
165 @*/
166 PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
167 {
168   va_list                Argp;
169   PetscErrorCode         ierr;
170   struct DMCompositeLink *next;
171   DM_Composite           *com = (DM_Composite*)dm->data;
172   PetscInt               readonly;
173 
174   PetscFunctionBegin;
175   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
176   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
177   next = com->next;
178   if (!com->setup) {
179     ierr = DMSetUp(dm);CHKERRQ(ierr);
180   }
181 
182   ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr);
183   /* loop over packed objects, handling one at at time */
184   va_start(Argp,gvec);
185   while (next) {
186     Vec *vec;
187     vec = va_arg(Argp, Vec*);
188     if (vec) {
189       ierr = DMGetGlobalVector(next->dm,vec);CHKERRQ(ierr);
190       if (readonly) {
191         const PetscScalar *array;
192         ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
193         ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
194         ierr = VecLockPush(*vec);CHKERRQ(ierr);
195         ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
196       } else {
197         PetscScalar *array;
198         ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
199         ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
200         ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
201       }
202     }
203     next = next->next;
204   }
205   va_end(Argp);
206   PetscFunctionReturn(0);
207 }
208 
209 /*@C
210     DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
211        representation.
212 
213     Collective on DMComposite
214 
215     Input Parameters:
216 +    dm - the packer object
217 .    pvec - packed vector
218 .    nwanted - number of vectors wanted
219 -    wanted - sorted array of vectors wanted, or NULL to get all vectors
220 
221     Output Parameters:
222 .    vecs - array of requested global vectors (must be allocated)
223 
224     Notes: Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them
225 
226     Level: advanced
227 
228 .seealso: DMCompositeGetAccess(), DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
229 @*/
230 PetscErrorCode  DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
231 {
232   PetscErrorCode         ierr;
233   struct DMCompositeLink *link;
234   PetscInt               i,wnum;
235   DM_Composite           *com = (DM_Composite*)dm->data;
236   PetscInt               readonly;
237 
238   PetscFunctionBegin;
239   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
240   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
241   if (!com->setup) {
242     ierr = DMSetUp(dm);CHKERRQ(ierr);
243   }
244 
245   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
246   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
247     if (!wanted || i == wanted[wnum]) {
248       Vec v;
249       ierr = DMGetGlobalVector(link->dm,&v);CHKERRQ(ierr);
250       if (readonly) {
251         const PetscScalar *array;
252         ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr);
253         ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
254         ierr = VecLockPush(v);CHKERRQ(ierr);
255         ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr);
256       } else {
257         PetscScalar *array;
258         ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
259         ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
260         ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
261       }
262       vecs[wnum++] = v;
263     }
264   }
265   PetscFunctionReturn(0);
266 }
267 
268 /*@C
269     DMCompositeGetLocalAccessArray - Allows one to access the individual
270     packed vectors in their local representation.
271 
272     Collective on DMComposite.
273 
274     Input Parameters:
275 +    dm - the packer object
276 .    pvec - packed vector
277 .    nwanted - number of vectors wanted
278 -    wanted - sorted array of vectors wanted, or NULL to get all vectors
279 
280     Output Parameters:
281 .    vecs - array of requested local vectors (must be allocated)
282 
283     Notes: Use DMCompositeRestoreLocalAccessArray() to return the vectors
284     when you no longer need them.
285 
286     Level: advanced
287 
288 .seealso: DMCompositeRestoreLocalAccessArray(), DMCompositeGetAccess(),
289 DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
290 @*/
291 PetscErrorCode  DMCompositeGetLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
292 {
293   PetscErrorCode         ierr;
294   struct DMCompositeLink *link;
295   PetscInt               i,wnum;
296   DM_Composite           *com = (DM_Composite*)dm->data;
297   PetscInt               readonly;
298   PetscInt               nlocal = 0;
299 
300   PetscFunctionBegin;
301   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
302   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
303   if (!com->setup) {
304     ierr = DMSetUp(dm);CHKERRQ(ierr);
305   }
306 
307   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
308   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
309     if (!wanted || i == wanted[wnum]) {
310       Vec v;
311       ierr = DMGetLocalVector(link->dm,&v);CHKERRQ(ierr);
312       if (readonly) {
313         const PetscScalar *array;
314         ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr);
315         ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr);
316         ierr = VecLockPush(v);CHKERRQ(ierr);
317         ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr);
318       } else {
319         PetscScalar *array;
320         ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
321         ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr);
322         ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
323       }
324       vecs[wnum++] = v;
325     }
326 
327     nlocal += link->nlocal;
328   }
329 
330   PetscFunctionReturn(0);
331 }
332 
333 /*@C
334     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
335        representation.
336 
337     Collective on DMComposite
338 
339     Input Parameters:
340 +    dm - the packer object
341 .    gvec - the global vector
342 -    Vec* ... - the individual parallel vectors, NULL for those that are not needed
343 
344     Level: advanced
345 
346 .seealso  DMCompositeAddDM(), DMCreateGlobalVector(),
347          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
348          DMCompositeRestoreAccess(), DMCompositeGetAccess()
349 
350 @*/
351 PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
352 {
353   va_list                Argp;
354   PetscErrorCode         ierr;
355   struct DMCompositeLink *next;
356   DM_Composite           *com = (DM_Composite*)dm->data;
357   PetscInt               readonly;
358 
359   PetscFunctionBegin;
360   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
361   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
362   next = com->next;
363   if (!com->setup) {
364     ierr = DMSetUp(dm);CHKERRQ(ierr);
365   }
366 
367   ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr);
368   /* loop over packed objects, handling one at at time */
369   va_start(Argp,gvec);
370   while (next) {
371     Vec *vec;
372     vec = va_arg(Argp, Vec*);
373     if (vec) {
374       ierr = VecResetArray(*vec);CHKERRQ(ierr);
375       if (readonly) {
376         ierr = VecLockPop(*vec);CHKERRQ(ierr);
377       }
378       ierr = DMRestoreGlobalVector(next->dm,vec);CHKERRQ(ierr);
379     }
380     next = next->next;
381   }
382   va_end(Argp);
383   PetscFunctionReturn(0);
384 }
385 
386 /*@C
387     DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray()
388 
389     Collective on DMComposite
390 
391     Input Parameters:
392 +    dm - the packer object
393 .    pvec - packed vector
394 .    nwanted - number of vectors wanted
395 .    wanted - sorted array of vectors wanted, or NULL to get all vectors
396 -    vecs - array of global vectors to return
397 
398     Level: advanced
399 
400 .seealso: DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), DMCompositeScatter(), DMCompositeGather()
401 @*/
402 PetscErrorCode  DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
403 {
404   PetscErrorCode         ierr;
405   struct DMCompositeLink *link;
406   PetscInt               i,wnum;
407   DM_Composite           *com = (DM_Composite*)dm->data;
408   PetscInt               readonly;
409 
410   PetscFunctionBegin;
411   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
412   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
413   if (!com->setup) {
414     ierr = DMSetUp(dm);CHKERRQ(ierr);
415   }
416 
417   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
418   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
419     if (!wanted || i == wanted[wnum]) {
420       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
421       if (readonly) {
422         ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr);
423       }
424       ierr = DMRestoreGlobalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
425       wnum++;
426     }
427   }
428   PetscFunctionReturn(0);
429 }
430 
431 /*@C
432     DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray().
433 
434     Collective on DMComposite.
435 
436     Input Parameters:
437 +    dm - the packer object
438 .    pvec - packed vector
439 .    nwanted - number of vectors wanted
440 .    wanted - sorted array of vectors wanted, or NULL to restore all vectors
441 -    vecs - array of local vectors to return
442 
443     Level: advanced
444 
445     Notes:
446     nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray()
447     otherwise the call will fail.
448 
449 .seealso: DMCompositeGetLocalAccessArray(), DMCompositeRestoreAccessArray(),
450 DMCompositeRestoreAccess(), DMCompositeRestoreEntries(),
451 DMCompositeScatter(), DMCompositeGather()
452 @*/
453 PetscErrorCode  DMCompositeRestoreLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
454 {
455   PetscErrorCode         ierr;
456   struct DMCompositeLink *link;
457   PetscInt               i,wnum;
458   DM_Composite           *com = (DM_Composite*)dm->data;
459   PetscInt               readonly;
460 
461   PetscFunctionBegin;
462   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
463   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
464   if (!com->setup) {
465     ierr = DMSetUp(dm);CHKERRQ(ierr);
466   }
467 
468   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
469   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
470     if (!wanted || i == wanted[wnum]) {
471       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
472       if (readonly) {
473         ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr);
474       }
475       ierr = DMRestoreLocalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
476       wnum++;
477     }
478   }
479   PetscFunctionReturn(0);
480 }
481 
482 /*@C
483     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
484 
485     Collective on DMComposite
486 
487     Input Parameters:
488 +    dm - the packer object
489 .    gvec - the global vector
490 -    Vec ... - the individual sequential vectors, NULL for those that are not needed
491 
492     Level: advanced
493 
494     Notes:
495     DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is
496     accessible from Fortran.
497 
498 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
499          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
500          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
501          DMCompositeScatterArray()
502 
503 @*/
504 PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
505 {
506   va_list                Argp;
507   PetscErrorCode         ierr;
508   struct DMCompositeLink *next;
509   PetscInt               cnt;
510   DM_Composite           *com = (DM_Composite*)dm->data;
511 
512   PetscFunctionBegin;
513   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
514   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
515   if (!com->setup) {
516     ierr = DMSetUp(dm);CHKERRQ(ierr);
517   }
518 
519   /* loop over packed objects, handling one at at time */
520   va_start(Argp,gvec);
521   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
522     Vec local;
523     local = va_arg(Argp, Vec);
524     if (local) {
525       Vec               global;
526       const PetscScalar *array;
527       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
528       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
529       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
530       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
531       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
532       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
533       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
534       ierr = VecResetArray(global);CHKERRQ(ierr);
535       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
536     }
537   }
538   va_end(Argp);
539   PetscFunctionReturn(0);
540 }
541 
542 /*@
543     DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
544 
545     Collective on DMComposite
546 
547     Input Parameters:
548 +    dm - the packer object
549 .    gvec - the global vector
550 .    lvecs - array of local vectors, NULL for any that are not needed
551 
552     Level: advanced
553 
554     Note:
555     This is a non-variadic alternative to DMCompositeScatter()
556 
557 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector()
558          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
559          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
560 
561 @*/
562 PetscErrorCode  DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs)
563 {
564   PetscErrorCode         ierr;
565   struct DMCompositeLink *next;
566   PetscInt               i;
567   DM_Composite           *com = (DM_Composite*)dm->data;
568 
569   PetscFunctionBegin;
570   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
571   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
572   if (!com->setup) {
573     ierr = DMSetUp(dm);CHKERRQ(ierr);
574   }
575 
576   /* loop over packed objects, handling one at at time */
577   for (i=0,next=com->next; next; next=next->next,i++) {
578     if (lvecs[i]) {
579       Vec         global;
580       const PetscScalar *array;
581       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
582       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
583       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
584       ierr = VecPlaceArray(global,(PetscScalar*)array+next->rstart);CHKERRQ(ierr);
585       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
586       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
587       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
588       ierr = VecResetArray(global);CHKERRQ(ierr);
589       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
590     }
591   }
592   PetscFunctionReturn(0);
593 }
594 
595 /*@C
596     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
597 
598     Collective on DMComposite
599 
600     Input Parameter:
601 +    dm - the packer object
602 .    gvec - the global vector
603 .    imode - INSERT_VALUES or ADD_VALUES
604 -    Vec ... - the individual sequential vectors, NULL for any that are not needed
605 
606     Level: advanced
607 
608 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
609          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
610          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
611 
612 @*/
613 PetscErrorCode  DMCompositeGather(DM dm,InsertMode imode,Vec gvec,...)
614 {
615   va_list                Argp;
616   PetscErrorCode         ierr;
617   struct DMCompositeLink *next;
618   DM_Composite           *com = (DM_Composite*)dm->data;
619   PetscInt               cnt;
620 
621   PetscFunctionBegin;
622   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
623   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
624   if (!com->setup) {
625     ierr = DMSetUp(dm);CHKERRQ(ierr);
626   }
627 
628   /* loop over packed objects, handling one at at time */
629   va_start(Argp,gvec);
630   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
631     Vec local;
632     local = va_arg(Argp, Vec);
633     if (local) {
634       PetscScalar *array;
635       Vec         global;
636       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
637       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
638       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
639       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
640       ierr = DMLocalToGlobalBegin(next->dm,local,imode,global);CHKERRQ(ierr);
641       ierr = DMLocalToGlobalEnd(next->dm,local,imode,global);CHKERRQ(ierr);
642       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
643       ierr = VecResetArray(global);CHKERRQ(ierr);
644       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
645     }
646   }
647   va_end(Argp);
648   PetscFunctionReturn(0);
649 }
650 
651 /*@
652     DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
653 
654     Collective on DMComposite
655 
656     Input Parameter:
657 +    dm - the packer object
658 .    gvec - the global vector
659 .    imode - INSERT_VALUES or ADD_VALUES
660 -    lvecs - the individual sequential vectors, NULL for any that are not needed
661 
662     Level: advanced
663 
664     Notes:
665     This is a non-variadic alternative to DMCompositeGather().
666 
667 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
668          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
669          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries(),
670 @*/
671 PetscErrorCode  DMCompositeGatherArray(DM dm,InsertMode imode,Vec gvec,Vec *lvecs)
672 {
673   PetscErrorCode         ierr;
674   struct DMCompositeLink *next;
675   DM_Composite           *com = (DM_Composite*)dm->data;
676   PetscInt               i;
677 
678   PetscFunctionBegin;
679   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
680   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
681   if (!com->setup) {
682     ierr = DMSetUp(dm);CHKERRQ(ierr);
683   }
684 
685   /* loop over packed objects, handling one at at time */
686   for (next=com->next,i=0; next; next=next->next,i++) {
687     if (lvecs[i]) {
688       PetscScalar *array;
689       Vec         global;
690       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
691       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
692       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
693       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
694       ierr = DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
695       ierr = DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
696       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
697       ierr = VecResetArray(global);CHKERRQ(ierr);
698       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
699     }
700   }
701   PetscFunctionReturn(0);
702 }
703 
704 /*@C
705     DMCompositeAddDM - adds a DM  vector to a DMComposite
706 
707     Collective on DMComposite
708 
709     Input Parameter:
710 +    dm - the packer object
711 -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
712 
713     Level: advanced
714 
715 .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
716          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
717          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
718 
719 @*/
720 PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
721 {
722   PetscErrorCode         ierr;
723   PetscInt               n,nlocal;
724   struct DMCompositeLink *mine,*next;
725   Vec                    global,local;
726   DM_Composite           *com = (DM_Composite*)dmc->data;
727 
728   PetscFunctionBegin;
729   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
730   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
731   next = com->next;
732   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
733 
734   /* create new link */
735   ierr = PetscNew(&mine);CHKERRQ(ierr);
736   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
737   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
738   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
739   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
740   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
741   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
742   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
743 
744   mine->n      = n;
745   mine->nlocal = nlocal;
746   mine->dm     = dm;
747   mine->next   = NULL;
748   com->n      += n;
749   com->nghost += nlocal;
750 
751   /* add to end of list */
752   if (!next) com->next = mine;
753   else {
754     while (next->next) next = next->next;
755     next->next = mine;
756   }
757   com->nDM++;
758   com->nmine++;
759   PetscFunctionReturn(0);
760 }
761 
762 #include <petscdraw.h>
763 PETSC_EXTERN PetscErrorCode  VecView_MPI(Vec,PetscViewer);
764 PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
765 {
766   DM                     dm;
767   PetscErrorCode         ierr;
768   struct DMCompositeLink *next;
769   PetscBool              isdraw;
770   DM_Composite           *com;
771 
772   PetscFunctionBegin;
773   ierr = VecGetDM(gvec, &dm);CHKERRQ(ierr);
774   if (!dm) SETERRQ(PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
775   com  = (DM_Composite*)dm->data;
776   next = com->next;
777 
778   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
779   if (!isdraw) {
780     /* do I really want to call this? */
781     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
782   } else {
783     PetscInt cnt = 0;
784 
785     /* loop over packed objects, handling one at at time */
786     while (next) {
787       Vec         vec;
788       PetscScalar *array;
789       PetscInt    bs;
790 
791       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
792       ierr = DMGetGlobalVector(next->dm,&vec);CHKERRQ(ierr);
793       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
794       ierr = VecPlaceArray(vec,array+next->rstart);CHKERRQ(ierr);
795       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
796       ierr = VecView(vec,viewer);CHKERRQ(ierr);
797       ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
798       ierr = VecResetArray(vec);CHKERRQ(ierr);
799       ierr = DMRestoreGlobalVector(next->dm,&vec);CHKERRQ(ierr);
800       ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
801       cnt += bs;
802       next = next->next;
803     }
804     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
805   }
806   PetscFunctionReturn(0);
807 }
808 
809 PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
810 {
811   PetscErrorCode ierr;
812   DM_Composite   *com = (DM_Composite*)dm->data;
813 
814   PetscFunctionBegin;
815   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
816   ierr = DMSetUp(dm);CHKERRQ(ierr);
817   ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),com->n,com->N,gvec);CHKERRQ(ierr);
818   ierr = VecSetDM(*gvec, dm);CHKERRQ(ierr);
819   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);CHKERRQ(ierr);
820   PetscFunctionReturn(0);
821 }
822 
823 PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
824 {
825   PetscErrorCode ierr;
826   DM_Composite   *com = (DM_Composite*)dm->data;
827 
828   PetscFunctionBegin;
829   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
830   if (!com->setup) {
831     ierr = DMSetUp(dm);CHKERRQ(ierr);
832   }
833   ierr = VecCreateSeq(PETSC_COMM_SELF,com->nghost,lvec);CHKERRQ(ierr);
834   ierr = VecSetDM(*lvec, dm);CHKERRQ(ierr);
835   PetscFunctionReturn(0);
836 }
837 
838 /*@C
839     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
840 
841     Collective on DM
842 
843     Input Parameter:
844 .    dm - the packer object
845 
846     Output Parameters:
847 .    ltogs - the individual mappings for each packed vector. Note that this includes
848            all the ghost points that individual ghosted DMDA's may have.
849 
850     Level: advanced
851 
852     Notes:
853        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
854 
855 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
856          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
857          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
858 
859 @*/
860 PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
861 {
862   PetscErrorCode         ierr;
863   PetscInt               i,*idx,n,cnt;
864   struct DMCompositeLink *next;
865   PetscMPIInt            rank;
866   DM_Composite           *com = (DM_Composite*)dm->data;
867 
868   PetscFunctionBegin;
869   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
870   ierr = DMSetUp(dm);CHKERRQ(ierr);
871   ierr = PetscMalloc1(com->nDM,ltogs);CHKERRQ(ierr);
872   next = com->next;
873   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
874 
875   /* loop over packed objects, handling one at at time */
876   cnt = 0;
877   while (next) {
878     ISLocalToGlobalMapping ltog;
879     PetscMPIInt            size;
880     const PetscInt         *suboff,*indices;
881     Vec                    global;
882 
883     /* Get sub-DM global indices for each local dof */
884     ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
885     ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
886     ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
887     ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr);
888 
889     /* Get the offsets for the sub-DM global vector */
890     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
891     ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
892     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)global),&size);CHKERRQ(ierr);
893 
894     /* Shift the sub-DM definition of the global space to the composite global space */
895     for (i=0; i<n; i++) {
896       PetscInt subi = indices[i],lo = 0,hi = size,t;
897       /* Binary search to find which rank owns subi */
898       while (hi-lo > 1) {
899         t = lo + (hi-lo)/2;
900         if (suboff[t] > subi) hi = t;
901         else                  lo = t;
902       }
903       idx[i] = subi - suboff[lo] + next->grstarts[lo];
904     }
905     ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
906     ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
907     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
908     next = next->next;
909     cnt++;
910   }
911   PetscFunctionReturn(0);
912 }
913 
914 /*@C
915    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
916 
917    Not Collective
918 
919    Input Arguments:
920 . dm - composite DM
921 
922    Output Arguments:
923 . is - array of serial index sets for each each component of the DMComposite
924 
925    Level: intermediate
926 
927    Notes:
928    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
929    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
930    scatter to a composite local vector.  The user should not typically need to know which is being done.
931 
932    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
933 
934    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
935 
936    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
937 
938 .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
939 @*/
940 PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
941 {
942   PetscErrorCode         ierr;
943   DM_Composite           *com = (DM_Composite*)dm->data;
944   struct DMCompositeLink *link;
945   PetscInt               cnt,start;
946 
947   PetscFunctionBegin;
948   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
949   PetscValidPointer(is,2);
950   ierr = PetscMalloc1(com->nmine,is);CHKERRQ(ierr);
951   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
952     PetscInt bs;
953     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
954     ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
955     ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
956   }
957   PetscFunctionReturn(0);
958 }
959 
960 /*@C
961     DMCompositeGetGlobalISs - Gets the index sets for each composed object
962 
963     Collective on DMComposite
964 
965     Input Parameter:
966 .    dm - the packer object
967 
968     Output Parameters:
969 .    is - the array of index sets
970 
971     Level: advanced
972 
973     Notes:
974        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
975 
976        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
977 
978        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
979        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
980        indices.
981 
982     Fortran Notes:
983 
984        The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM().
985 
986 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
987          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
988          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
989 
990 @*/
991 PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
992 {
993   PetscErrorCode         ierr;
994   PetscInt               cnt = 0;
995   struct DMCompositeLink *next;
996   PetscMPIInt            rank;
997   DM_Composite           *com = (DM_Composite*)dm->data;
998 
999   PetscFunctionBegin;
1000   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1001   ierr = PetscMalloc1(com->nDM,is);CHKERRQ(ierr);
1002   next = com->next;
1003   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1004 
1005   /* loop over packed objects, handling one at at time */
1006   while (next) {
1007     ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt]);CHKERRQ(ierr);
1008     if (dm->prob) {
1009       MatNullSpace space;
1010       Mat          pmat;
1011       PetscObject  disc;
1012       PetscInt     Nf;
1013 
1014       ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
1015       if (cnt < Nf) {
1016         ierr = PetscDSGetDiscretization(dm->prob, cnt, &disc);CHKERRQ(ierr);
1017         ierr = PetscObjectQuery(disc, "nullspace", (PetscObject*) &space);CHKERRQ(ierr);
1018         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space);CHKERRQ(ierr);}
1019         ierr = PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space);CHKERRQ(ierr);
1020         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space);CHKERRQ(ierr);}
1021         ierr = PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
1022         if (pmat) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat);CHKERRQ(ierr);}
1023       }
1024     }
1025     cnt++;
1026     next = next->next;
1027   }
1028   PetscFunctionReturn(0);
1029 }
1030 
1031 PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields)
1032 {
1033   PetscInt       nDM;
1034   DM             *dms;
1035   PetscInt       i;
1036   PetscErrorCode ierr;
1037 
1038   PetscFunctionBegin;
1039   ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
1040   if (numFields) *numFields = nDM;
1041   ierr = DMCompositeGetGlobalISs(dm, fields);CHKERRQ(ierr);
1042   if (fieldNames) {
1043     ierr = PetscMalloc1(nDM, &dms);CHKERRQ(ierr);
1044     ierr = PetscMalloc1(nDM, fieldNames);CHKERRQ(ierr);
1045     ierr = DMCompositeGetEntriesArray(dm, dms);CHKERRQ(ierr);
1046     for (i=0; i<nDM; i++) {
1047       char       buf[256];
1048       const char *splitname;
1049 
1050       /* Split naming precedence: object name, prefix, number */
1051       splitname = ((PetscObject) dm)->name;
1052       if (!splitname) {
1053         ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr);
1054         if (splitname) {
1055           size_t len;
1056           ierr                 = PetscStrncpy(buf,splitname,sizeof(buf));CHKERRQ(ierr);
1057           buf[sizeof(buf) - 1] = 0;
1058           ierr                 = PetscStrlen(buf,&len);CHKERRQ(ierr);
1059           if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
1060           splitname = buf;
1061         }
1062       }
1063       if (!splitname) {
1064         ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",i);CHKERRQ(ierr);
1065         splitname = buf;
1066       }
1067       ierr = PetscStrallocpy(splitname,&(*fieldNames)[i]);CHKERRQ(ierr);
1068     }
1069     ierr = PetscFree(dms);CHKERRQ(ierr);
1070   }
1071   PetscFunctionReturn(0);
1072 }
1073 
1074 /*
1075  This could take over from DMCreateFieldIS(), as it is more general,
1076  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1077  At this point it's probably best to be less intrusive, however.
1078  */
1079 PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
1080 {
1081   PetscInt       nDM;
1082   PetscInt       i;
1083   PetscErrorCode ierr;
1084 
1085   PetscFunctionBegin;
1086   ierr = DMCreateFieldIS_Composite(dm, len, namelist, islist);CHKERRQ(ierr);
1087   if (dmlist) {
1088     ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
1089     ierr = PetscMalloc1(nDM, dmlist);CHKERRQ(ierr);
1090     ierr = DMCompositeGetEntriesArray(dm, *dmlist);CHKERRQ(ierr);
1091     for (i=0; i<nDM; i++) {
1092       ierr = PetscObjectReference((PetscObject)((*dmlist)[i]));CHKERRQ(ierr);
1093     }
1094   }
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 
1099 
1100 /* -------------------------------------------------------------------------------------*/
1101 /*@C
1102     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
1103        Use DMCompositeRestoreLocalVectors() to return them.
1104 
1105     Not Collective
1106 
1107     Input Parameter:
1108 .    dm - the packer object
1109 
1110     Output Parameter:
1111 .   Vec ... - the individual sequential Vecs
1112 
1113     Level: advanced
1114 
1115 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
1116          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1117          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
1118 
1119 @*/
1120 PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
1121 {
1122   va_list                Argp;
1123   PetscErrorCode         ierr;
1124   struct DMCompositeLink *next;
1125   DM_Composite           *com = (DM_Composite*)dm->data;
1126 
1127   PetscFunctionBegin;
1128   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1129   next = com->next;
1130   /* loop over packed objects, handling one at at time */
1131   va_start(Argp,dm);
1132   while (next) {
1133     Vec *vec;
1134     vec = va_arg(Argp, Vec*);
1135     if (vec) {ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);}
1136     next = next->next;
1137   }
1138   va_end(Argp);
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 /*@C
1143     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
1144 
1145     Not Collective
1146 
1147     Input Parameter:
1148 .    dm - the packer object
1149 
1150     Output Parameter:
1151 .   Vec ... - the individual sequential Vecs
1152 
1153     Level: advanced
1154 
1155 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
1156          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1157          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
1158 
1159 @*/
1160 PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
1161 {
1162   va_list                Argp;
1163   PetscErrorCode         ierr;
1164   struct DMCompositeLink *next;
1165   DM_Composite           *com = (DM_Composite*)dm->data;
1166 
1167   PetscFunctionBegin;
1168   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1169   next = com->next;
1170   /* loop over packed objects, handling one at at time */
1171   va_start(Argp,dm);
1172   while (next) {
1173     Vec *vec;
1174     vec = va_arg(Argp, Vec*);
1175     if (vec) {ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);}
1176     next = next->next;
1177   }
1178   va_end(Argp);
1179   PetscFunctionReturn(0);
1180 }
1181 
1182 /* -------------------------------------------------------------------------------------*/
1183 /*@C
1184     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
1185 
1186     Not Collective
1187 
1188     Input Parameter:
1189 .    dm - the packer object
1190 
1191     Output Parameter:
1192 .   DM ... - the individual entries (DMs)
1193 
1194     Level: advanced
1195 
1196 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
1197          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1198          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
1199          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
1200 
1201 @*/
1202 PetscErrorCode  DMCompositeGetEntries(DM dm,...)
1203 {
1204   va_list                Argp;
1205   struct DMCompositeLink *next;
1206   DM_Composite           *com = (DM_Composite*)dm->data;
1207 
1208   PetscFunctionBegin;
1209   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1210   next = com->next;
1211   /* loop over packed objects, handling one at at time */
1212   va_start(Argp,dm);
1213   while (next) {
1214     DM *dmn;
1215     dmn = va_arg(Argp, DM*);
1216     if (dmn) *dmn = next->dm;
1217     next = next->next;
1218   }
1219   va_end(Argp);
1220   PetscFunctionReturn(0);
1221 }
1222 
1223 /*@C
1224     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
1225 
1226     Not Collective
1227 
1228     Input Parameter:
1229 .    dm - the packer object
1230 
1231     Output Parameter:
1232 .    dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs
1233 
1234     Level: advanced
1235 
1236 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
1237          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1238          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
1239          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
1240 
1241 @*/
1242 PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
1243 {
1244   struct DMCompositeLink *next;
1245   DM_Composite           *com = (DM_Composite*)dm->data;
1246   PetscInt               i;
1247 
1248   PetscFunctionBegin;
1249   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1250   /* loop over packed objects, handling one at at time */
1251   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
1252   PetscFunctionReturn(0);
1253 }
1254 
1255 PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
1256 {
1257   PetscErrorCode         ierr;
1258   struct DMCompositeLink *next;
1259   DM_Composite           *com = (DM_Composite*)dmi->data;
1260   DM                     dm;
1261 
1262   PetscFunctionBegin;
1263   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1264   if (comm == MPI_COMM_NULL) {
1265     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
1266   }
1267   ierr = DMSetUp(dmi);CHKERRQ(ierr);
1268   next = com->next;
1269   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
1270 
1271   /* loop over packed objects, handling one at at time */
1272   while (next) {
1273     ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
1274     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
1275     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1276     next = next->next;
1277   }
1278   PetscFunctionReturn(0);
1279 }
1280 
1281 PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
1282 {
1283   PetscErrorCode         ierr;
1284   struct DMCompositeLink *next;
1285   DM_Composite           *com = (DM_Composite*)dmi->data;
1286   DM                     dm;
1287 
1288   PetscFunctionBegin;
1289   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1290   ierr = DMSetUp(dmi);CHKERRQ(ierr);
1291   if (comm == MPI_COMM_NULL) {
1292     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
1293   }
1294   next = com->next;
1295   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
1296 
1297   /* loop over packed objects, handling one at at time */
1298   while (next) {
1299     ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr);
1300     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
1301     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1302     next = next->next;
1303   }
1304   PetscFunctionReturn(0);
1305 }
1306 
1307 PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
1308 {
1309   PetscErrorCode         ierr;
1310   PetscInt               m,n,M,N,nDM,i;
1311   struct DMCompositeLink *nextc;
1312   struct DMCompositeLink *nextf;
1313   Vec                    gcoarse,gfine,*vecs;
1314   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
1315   DM_Composite           *comfine   = (DM_Composite*)fine->data;
1316   Mat                    *mats;
1317 
1318   PetscFunctionBegin;
1319   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
1320   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
1321   ierr = DMSetUp(coarse);CHKERRQ(ierr);
1322   ierr = DMSetUp(fine);CHKERRQ(ierr);
1323   /* use global vectors only for determining matrix layout */
1324   ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
1325   ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr);
1326   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
1327   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
1328   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
1329   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
1330   ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
1331   ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr);
1332 
1333   nDM = comfine->nDM;
1334   if (nDM != comcoarse->nDM) SETERRQ2(PetscObjectComm((PetscObject)fine),PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %D entries, but coarse has %D",nDM,comcoarse->nDM);
1335   ierr = PetscCalloc1(nDM*nDM,&mats);CHKERRQ(ierr);
1336   if (v) {
1337     ierr = PetscCalloc1(nDM,&vecs);CHKERRQ(ierr);
1338   }
1339 
1340   /* loop over packed objects, handling one at at time */
1341   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
1342     if (!v) {
1343       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL);CHKERRQ(ierr);
1344     } else {
1345       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr);
1346     }
1347   }
1348   ierr = MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A);CHKERRQ(ierr);
1349   if (v) {
1350     ierr = VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v);CHKERRQ(ierr);
1351   }
1352   for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);}
1353   ierr = PetscFree(mats);CHKERRQ(ierr);
1354   if (v) {
1355     for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);}
1356     ierr = PetscFree(vecs);CHKERRQ(ierr);
1357   }
1358   PetscFunctionReturn(0);
1359 }
1360 
1361 static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
1362 {
1363   DM_Composite           *com = (DM_Composite*)dm->data;
1364   ISLocalToGlobalMapping *ltogs;
1365   PetscInt               i;
1366   PetscErrorCode         ierr;
1367 
1368   PetscFunctionBegin;
1369   /* Set the ISLocalToGlobalMapping on the new matrix */
1370   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
1371   ierr = ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr);
1372   for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
1373   ierr = PetscFree(ltogs);CHKERRQ(ierr);
1374   PetscFunctionReturn(0);
1375 }
1376 
1377 
1378 PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring)
1379 {
1380   PetscErrorCode  ierr;
1381   PetscInt        n,i,cnt;
1382   ISColoringValue *colors;
1383   PetscBool       dense  = PETSC_FALSE;
1384   ISColoringValue maxcol = 0;
1385   DM_Composite    *com   = (DM_Composite*)dm->data;
1386 
1387   PetscFunctionBegin;
1388   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1389   if (ctype == IS_COLORING_LOCAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported");
1390   else if (ctype == IS_COLORING_GLOBAL) {
1391     n = com->n;
1392   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1393   ierr = PetscMalloc1(n,&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
1394 
1395   ierr = PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL);CHKERRQ(ierr);
1396   if (dense) {
1397     for (i=0; i<n; i++) {
1398       colors[i] = (ISColoringValue)(com->rstart + i);
1399     }
1400     maxcol = com->N;
1401   } else {
1402     struct DMCompositeLink *next = com->next;
1403     PetscMPIInt            rank;
1404 
1405     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1406     cnt  = 0;
1407     while (next) {
1408       ISColoring lcoloring;
1409 
1410       ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring);CHKERRQ(ierr);
1411       for (i=0; i<lcoloring->N; i++) {
1412         colors[cnt++] = maxcol + lcoloring->colors[i];
1413       }
1414       maxcol += lcoloring->n;
1415       ierr    = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
1416       next    = next->next;
1417     }
1418   }
1419   ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,PETSC_OWN_POINTER,coloring);CHKERRQ(ierr);
1420   PetscFunctionReturn(0);
1421 }
1422 
1423 PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1424 {
1425   PetscErrorCode         ierr;
1426   struct DMCompositeLink *next;
1427   PetscScalar            *garray,*larray;
1428   DM_Composite           *com = (DM_Composite*)dm->data;
1429 
1430   PetscFunctionBegin;
1431   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1432   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
1433 
1434   if (!com->setup) {
1435     ierr = DMSetUp(dm);CHKERRQ(ierr);
1436   }
1437 
1438   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
1439   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
1440 
1441   /* loop over packed objects, handling one at at time */
1442   next = com->next;
1443   while (next) {
1444     Vec      local,global;
1445     PetscInt N;
1446 
1447     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
1448     ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
1449     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
1450     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
1451     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
1452     ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
1453     ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
1454     ierr = VecResetArray(global);CHKERRQ(ierr);
1455     ierr = VecResetArray(local);CHKERRQ(ierr);
1456     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
1457     ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr);
1458 
1459     larray += next->nlocal;
1460     garray += next->n;
1461     next    = next->next;
1462   }
1463 
1464   ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr);
1465   ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr);
1466   PetscFunctionReturn(0);
1467 }
1468 
1469 PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1470 {
1471   PetscFunctionBegin;
1472   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1473   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
1474   PetscValidHeaderSpecific(lvec,VEC_CLASSID,4);
1475   PetscFunctionReturn(0);
1476 }
1477 
1478 PetscErrorCode  DMLocalToGlobalBegin_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
1479 {
1480   PetscErrorCode         ierr;
1481   struct DMCompositeLink *next;
1482   PetscScalar            *larray,*garray;
1483   DM_Composite           *com = (DM_Composite*)dm->data;
1484 
1485   PetscFunctionBegin;
1486   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1487   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
1488   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
1489 
1490   if (!com->setup) {
1491     ierr = DMSetUp(dm);CHKERRQ(ierr);
1492   }
1493 
1494   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
1495   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
1496 
1497   /* loop over packed objects, handling one at at time */
1498   next = com->next;
1499   while (next) {
1500     Vec      global,local;
1501 
1502     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
1503     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
1504     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
1505     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
1506     ierr = DMLocalToGlobalBegin(next->dm,local,mode,global);CHKERRQ(ierr);
1507     ierr = DMLocalToGlobalEnd(next->dm,local,mode,global);CHKERRQ(ierr);
1508     ierr = VecResetArray(local);CHKERRQ(ierr);
1509     ierr = VecResetArray(global);CHKERRQ(ierr);
1510     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
1511     ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr);
1512 
1513     garray += next->n;
1514     larray += next->nlocal;
1515     next    = next->next;
1516   }
1517 
1518   ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr);
1519   ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr);
1520 
1521   PetscFunctionReturn(0);
1522 }
1523 
1524 PetscErrorCode  DMLocalToGlobalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
1525 {
1526   PetscFunctionBegin;
1527   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1528   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
1529   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
1530   PetscFunctionReturn(0);
1531 }
1532 
1533 PetscErrorCode  DMLocalToLocalBegin_Composite(DM dm,Vec vec1,InsertMode mode,Vec vec2)
1534 {
1535   PetscErrorCode         ierr;
1536   struct DMCompositeLink *next;
1537   PetscScalar            *array1,*array2;
1538   DM_Composite           *com = (DM_Composite*)dm->data;
1539 
1540   PetscFunctionBegin;
1541   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1542   PetscValidHeaderSpecific(vec1,VEC_CLASSID,2);
1543   PetscValidHeaderSpecific(vec2,VEC_CLASSID,4);
1544 
1545   if (!com->setup) {
1546     ierr = DMSetUp(dm);CHKERRQ(ierr);
1547   }
1548 
1549   ierr = VecGetArray(vec1,&array1);CHKERRQ(ierr);
1550   ierr = VecGetArray(vec2,&array2);CHKERRQ(ierr);
1551 
1552   /* loop over packed objects, handling one at at time */
1553   next = com->next;
1554   while (next) {
1555     Vec      local1,local2;
1556 
1557     ierr = DMGetLocalVector(next->dm,&local1);CHKERRQ(ierr);
1558     ierr = VecPlaceArray(local1,array1);CHKERRQ(ierr);
1559     ierr = DMGetLocalVector(next->dm,&local2);CHKERRQ(ierr);
1560     ierr = VecPlaceArray(local2,array2);CHKERRQ(ierr);
1561     ierr = DMLocalToLocalBegin(next->dm,local1,mode,local2);CHKERRQ(ierr);
1562     ierr = DMLocalToLocalEnd(next->dm,local1,mode,local2);CHKERRQ(ierr);
1563     ierr = VecResetArray(local2);CHKERRQ(ierr);
1564     ierr = DMRestoreLocalVector(next->dm,&local2);CHKERRQ(ierr);
1565     ierr = VecResetArray(local1);CHKERRQ(ierr);
1566     ierr = DMRestoreLocalVector(next->dm,&local1);CHKERRQ(ierr);
1567 
1568     array1 += next->nlocal;
1569     array2 += next->nlocal;
1570     next    = next->next;
1571   }
1572 
1573   ierr = VecRestoreArray(vec1,NULL);CHKERRQ(ierr);
1574   ierr = VecRestoreArray(vec2,NULL);CHKERRQ(ierr);
1575 
1576   PetscFunctionReturn(0);
1577 }
1578 
1579 PetscErrorCode  DMLocalToLocalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
1580 {
1581   PetscFunctionBegin;
1582   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1583   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
1584   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
1585   PetscFunctionReturn(0);
1586 }
1587 
1588 /*MC
1589    DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs
1590 
1591   Level: intermediate
1592 
1593 .seealso: DMType, DM, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate()
1594 M*/
1595 
1596 
1597 PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1598 {
1599   PetscErrorCode ierr;
1600   DM_Composite   *com;
1601 
1602   PetscFunctionBegin;
1603   ierr          = PetscNewLog(p,&com);CHKERRQ(ierr);
1604   p->data       = com;
1605   ierr          = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1606   com->n        = 0;
1607   com->nghost   = 0;
1608   com->next     = NULL;
1609   com->nDM      = 0;
1610 
1611   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1612   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1613   p->ops->getlocaltoglobalmapping         = DMGetLocalToGlobalMapping_Composite;
1614   p->ops->createfieldis                   = DMCreateFieldIS_Composite;
1615   p->ops->createfielddecomposition        = DMCreateFieldDecomposition_Composite;
1616   p->ops->refine                          = DMRefine_Composite;
1617   p->ops->coarsen                         = DMCoarsen_Composite;
1618   p->ops->createinterpolation             = DMCreateInterpolation_Composite;
1619   p->ops->creatematrix                    = DMCreateMatrix_Composite;
1620   p->ops->getcoloring                     = DMCreateColoring_Composite;
1621   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1622   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1623   p->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Composite;
1624   p->ops->localtoglobalend                = DMLocalToGlobalEnd_Composite;
1625   p->ops->localtolocalbegin               = DMLocalToLocalBegin_Composite;
1626   p->ops->localtolocalend                 = DMLocalToLocalEnd_Composite;
1627   p->ops->destroy                         = DMDestroy_Composite;
1628   p->ops->view                            = DMView_Composite;
1629   p->ops->setup                           = DMSetUp_Composite;
1630   PetscFunctionReturn(0);
1631 }
1632 
1633 /*@C
1634     DMCompositeCreate - Creates a vector packer, used to generate "composite"
1635       vectors made up of several subvectors.
1636 
1637     Collective on MPI_Comm
1638 
1639     Input Parameter:
1640 .   comm - the processors that will share the global vector
1641 
1642     Output Parameters:
1643 .   packer - the packer object
1644 
1645     Level: advanced
1646 
1647 .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(), DMCOMPOSITE,DMCreate()
1648          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
1649          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
1650 
1651 @*/
1652 PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
1653 {
1654   PetscErrorCode ierr;
1655 
1656   PetscFunctionBegin;
1657   PetscValidPointer(packer,2);
1658   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1659   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
1660   PetscFunctionReturn(0);
1661 }
1662