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