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