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