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