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