xref: /petsc/src/dm/impls/composite/pack.c (revision 0c010503e280889cfcd2455cff2d09d8471c4db0)
1 #define PETSCDM_DLL
2 
3 #include "petscda.h"             /*I      "petscda.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 DA 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 (DA'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__ "DMCompositeSetUp"
187 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeSetUp(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,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     ierr    = PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(PetscScalar));CHKERRQ(ierr);
341     ierr    = VecRestoreArray(vec,&varray);CHKERRQ(ierr);
342   }
343   PetscFunctionReturn(0);
344 }
345 
346 #undef __FUNCT__
347 #define __FUNCT__ "DMCompositeGather_DM"
348 PetscErrorCode DMCompositeGather_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec local)
349 {
350   PetscErrorCode ierr;
351   PetscScalar    *array;
352   Vec            global;
353 
354   PetscFunctionBegin;
355   ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr);
356   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
357   ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr);
358   ierr = DMLocalToGlobal(mine->dm,local,INSERT_VALUES,global);CHKERRQ(ierr);
359   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
360   ierr = VecResetArray(global);CHKERRQ(ierr);
361   ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr);
362   PetscFunctionReturn(0);
363 }
364 
365 /* ----------------------------------------------------------------------------------*/
366 
367 #include <stdarg.h>
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "DMCompositeGetNumberDM"
371 /*@C
372     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
373        representation.
374 
375     Not Collective
376 
377     Input Parameter:
378 .    dm - the packer object
379 
380     Output Parameter:
381 .     nDM - the number of DMs
382 
383     Level: beginner
384 
385 @*/
386 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
387 {
388   DM_Composite *com = (DM_Composite*)dm->data;
389   PetscFunctionBegin;
390   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
391   *nDM = com->nDM;
392   PetscFunctionReturn(0);
393 }
394 
395 
396 #undef __FUNCT__
397 #define __FUNCT__ "DMCompositeGetAccess"
398 /*@C
399     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
400        representation.
401 
402     Collective on DMComposite
403 
404     Input Parameter:
405 +    dm - the packer object
406 .    gvec - the global vector
407 -    ... - the individual sequential or parallel objects (arrays or vectors)
408 
409     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
410 
411     Level: advanced
412 
413 @*/
414 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetAccess(DM dm,Vec gvec,...)
415 {
416   va_list                Argp;
417   PetscErrorCode         ierr;
418   struct DMCompositeLink *next;
419   DM_Composite           *com = (DM_Composite*)dm->data;
420 
421   PetscFunctionBegin;
422   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
423   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
424   next = com->next;
425   if (!com->setup) {
426     ierr = DMCompositeSetUp(dm);CHKERRQ(ierr);
427   }
428 
429   /* loop over packed objects, handling one at at time */
430   va_start(Argp,gvec);
431   while (next) {
432     if (next->type == DMCOMPOSITE_ARRAY) {
433       PetscScalar **array;
434       array = va_arg(Argp, PetscScalar**);
435       ierr  = DMCompositeGetAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
436     } else if (next->type == DMCOMPOSITE_DM) {
437       Vec *vec;
438       vec  = va_arg(Argp, Vec*);
439       ierr = DMCompositeGetAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
440     } else {
441       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
442     }
443     next = next->next;
444   }
445   va_end(Argp);
446   PetscFunctionReturn(0);
447 }
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "DMCompositeRestoreAccess"
451 /*@C
452     DMCompositeRestoreAccess - Returns the vectors obtained with DACompositeGetAccess()
453        representation.
454 
455     Collective on DMComposite
456 
457     Input Parameter:
458 +    dm - the packer object
459 .    gvec - the global vector
460 -    ... - the individual sequential or parallel objects (arrays or vectors)
461 
462     Level: advanced
463 
464 .seealso  DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
465          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeScatter(),
466          DMCompositeRestoreAccess(), DACompositeGetAccess()
467 
468 @*/
469 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreAccess(DM dm,Vec gvec,...)
470 {
471   va_list                Argp;
472   PetscErrorCode         ierr;
473   struct DMCompositeLink *next;
474   DM_Composite           *com = (DM_Composite*)dm->data;
475 
476   PetscFunctionBegin;
477   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
478   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
479   next = com->next;
480   if (!com->setup) {
481     ierr = DMCompositeSetUp(dm);CHKERRQ(ierr);
482   }
483 
484   /* loop over packed objects, handling one at at time */
485   va_start(Argp,gvec);
486   while (next) {
487     if (next->type == DMCOMPOSITE_ARRAY) {
488       PetscScalar **array;
489       array = va_arg(Argp, PetscScalar**);
490       ierr  = DMCompositeRestoreAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
491     } else if (next->type == DMCOMPOSITE_DM) {
492       Vec *vec;
493       vec  = va_arg(Argp, Vec*);
494       ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
495     } else {
496       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
497     }
498     next = next->next;
499   }
500   va_end(Argp);
501   PetscFunctionReturn(0);
502 }
503 
504 #undef __FUNCT__
505 #define __FUNCT__ "DMCompositeScatter"
506 /*@C
507     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
508 
509     Collective on DMComposite
510 
511     Input Parameter:
512 +    dm - the packer object
513 .    gvec - the global vector
514 -    ... - the individual sequential objects (arrays or vectors)
515 
516     Level: advanced
517 
518 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
519          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(),
520          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
521 
522 @*/
523 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeScatter(DM dm,Vec gvec,...)
524 {
525   va_list                Argp;
526   PetscErrorCode         ierr;
527   struct DMCompositeLink *next;
528   PetscInt               cnt = 3;
529   DM_Composite           *com = (DM_Composite*)dm->data;
530 
531   PetscFunctionBegin;
532   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
533   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
534   next = com->next;
535   if (!com->setup) {
536     ierr = DMCompositeSetUp(dm);CHKERRQ(ierr);
537   }
538 
539   /* loop over packed objects, handling one at at time */
540   va_start(Argp,gvec);
541   while (next) {
542     if (next->type == DMCOMPOSITE_ARRAY) {
543       PetscScalar *array;
544       array = va_arg(Argp, PetscScalar*);
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       PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt);
550       ierr = DMCompositeScatter_DM(dm,next,gvec,vec);CHKERRQ(ierr);
551     } else {
552       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
553     }
554     cnt++;
555     next = next->next;
556   }
557   va_end(Argp);
558   PetscFunctionReturn(0);
559 }
560 
561 #undef __FUNCT__
562 #define __FUNCT__ "DMCompositeGather"
563 /*@C
564     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
565 
566     Collective on DMComposite
567 
568     Input Parameter:
569 +    dm - the packer object
570 .    gvec - the global vector
571 -    ... - the individual sequential objects (arrays or vectors)
572 
573     Level: advanced
574 
575 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
576          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(),
577          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
578 
579 @*/
580 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGather(DM dm,Vec gvec,...)
581 {
582   va_list                Argp;
583   PetscErrorCode         ierr;
584   struct DMCompositeLink *next;
585   DM_Composite           *com = (DM_Composite*)dm->data;
586 
587   PetscFunctionBegin;
588   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
589   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
590   next = com->next;
591   if (!com->setup) {
592     ierr = DMCompositeSetUp(dm);CHKERRQ(ierr);
593   }
594 
595   /* loop over packed objects, handling one at at time */
596   va_start(Argp,gvec);
597   while (next) {
598     if (next->type == DMCOMPOSITE_ARRAY) {
599       PetscScalar *array;
600       array = va_arg(Argp, PetscScalar*);
601       ierr  = DMCompositeGather_Array(dm,next,gvec,array);CHKERRQ(ierr);
602     } else if (next->type == DMCOMPOSITE_DM) {
603       Vec vec;
604       vec = va_arg(Argp, Vec);
605       PetscValidHeaderSpecific(vec,VEC_CLASSID,3);
606       ierr = DMCompositeGather_DM(dm,next,gvec,vec);CHKERRQ(ierr);
607     } else {
608       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
609     }
610     next = next->next;
611   }
612   va_end(Argp);
613   PetscFunctionReturn(0);
614 }
615 
616 #undef __FUNCT__
617 #define __FUNCT__ "DMCompositeAddArray"
618 /*@C
619     DMCompositeAddArray - adds an "redundant" array to a DMComposite. The array values will
620        be stored in part of the array on process orank.
621 
622     Collective on DMComposite
623 
624     Input Parameter:
625 +    dm - the packer object
626 .    orank - the process on which the array entries officially live, this number must be
627              the same on all processes.
628 -    n - the length of the array
629 
630     Level: advanced
631 
632 .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
633          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(),
634          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
635 
636 @*/
637 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeAddArray(DM dm,PetscMPIInt orank,PetscInt n)
638 {
639   struct DMCompositeLink *mine,*next;
640   PetscErrorCode         ierr;
641   PetscMPIInt            rank;
642   DM_Composite           *com = (DM_Composite*)dm->data;
643 
644   PetscFunctionBegin;
645   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
646   next = com->next;
647   if (com->setup) {
648     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add an array once you have used the DMComposite");
649   }
650 #if defined(PETSC_USE_DEBUG)
651   {
652     PetscMPIInt orankmax;
653     ierr = MPI_Allreduce(&orank,&orankmax,1,MPI_INT,MPI_MAX,((PetscObject)dm)->comm);CHKERRQ(ierr);
654     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);
655   }
656 #endif
657 
658   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
659   /* create new link */
660   ierr                = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
661   mine->n             = n;
662   mine->rank          = orank;
663   mine->dm            = PETSC_NULL;
664   mine->type          = DMCOMPOSITE_ARRAY;
665   mine->next          = PETSC_NULL;
666   if (rank == mine->rank) {com->n += n;com->nmine++;}
667 
668   /* add to end of list */
669   if (!next) {
670     com->next = mine;
671   } else {
672     while (next->next) next = next->next;
673     next->next = mine;
674   }
675   com->nredundant++;
676   PetscFunctionReturn(0);
677 }
678 
679 #undef __FUNCT__
680 #define __FUNCT__ "DMCompositeAddDM"
681 /*@C
682     DMCompositeAddDM - adds a DM (includes DA) vector to a DMComposite
683 
684     Collective on DMComposite
685 
686     Input Parameter:
687 +    dm - the packer object
688 -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
689 
690     Level: advanced
691 
692 .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
693          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(),
694          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
695 
696 @*/
697 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeAddDM(DM dmc,DM dm)
698 {
699   PetscErrorCode         ierr;
700   PetscInt               n;
701   struct DMCompositeLink *mine,*next;
702   Vec                    global;
703   DM_Composite           *com = (DM_Composite*)dmc->data;
704 
705   PetscFunctionBegin;
706   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
707   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
708   next = com->next;
709   if (com->setup) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DA once you have used the DMComposite");
710 
711   /* create new link */
712   ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
713   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
714   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
715   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
716   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
717   mine->n      = n;
718   mine->dm     = dm;
719   mine->type   = DMCOMPOSITE_DM;
720   mine->next   = PETSC_NULL;
721   com->n       += n;
722 
723   /* add to end of list */
724   if (!next) {
725     com->next = mine;
726   } else {
727     while (next->next) next = next->next;
728     next->next = mine;
729   }
730   com->nDM++;
731   com->nmine++;
732   PetscFunctionReturn(0);
733 }
734 
735 extern PetscErrorCode PETSCDM_DLLEXPORT VecView_MPI(Vec,PetscViewer);
736 EXTERN_C_BEGIN
737 #undef __FUNCT__
738 #define __FUNCT__ "VecView_DMComposite"
739 PetscErrorCode PETSCDM_DLLEXPORT VecView_DMComposite(Vec gvec,PetscViewer viewer)
740 {
741   DM                     dm;
742   PetscErrorCode         ierr;
743   struct DMCompositeLink *next;
744   PetscBool              isdraw;
745   DM_Composite           *com = (DM_Composite*)dm->data;
746 
747   PetscFunctionBegin;
748   ierr = PetscObjectQuery((PetscObject)gvec,"DMComposite",(PetscObject*)&dm);CHKERRQ(ierr);
749   if (!dm) SETERRQ(((PetscObject)gvec)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
750   com = (DM_Composite*)dm->data;
751   next = com->next;
752 
753   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
754   if (!isdraw) {
755     /* do I really want to call this? */
756     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
757   } else {
758     PetscInt cnt = 0;
759 
760     /* loop over packed objects, handling one at at time */
761     while (next) {
762       if (next->type == DMCOMPOSITE_ARRAY) {
763 	PetscScalar *array;
764 	ierr  = DMCompositeGetAccess_Array(dm,next,gvec,&array);CHKERRQ(ierr);
765 
766 	/*skip it for now */
767       } else if (next->type == DMCOMPOSITE_DM) {
768 	Vec      vec;
769         PetscInt bs;
770 
771 	ierr = DMCompositeGetAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
772 	ierr = VecView(vec,viewer);CHKERRQ(ierr);
773         ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
774 	ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
775         ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
776         cnt += bs;
777       } else {
778 	SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
779       }
780       next = next->next;
781     }
782     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
783   }
784   PetscFunctionReturn(0);
785 }
786 EXTERN_C_END
787 
788 
789 #undef __FUNCT__
790 #define __FUNCT__ "DMCreateGlobalVector_Composite"
791 PetscErrorCode PETSCDM_DLLEXPORT DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
792 {
793   PetscErrorCode         ierr;
794   DM_Composite           *com = (DM_Composite*)dm->data;
795 
796   PetscFunctionBegin;
797   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
798   if (!com->setup) {
799     ierr = DMCompositeSetUp(dm);CHKERRQ(ierr);
800   }
801   ierr = VecCreateMPI(((PetscObject)dm)->comm,com->n,com->N,gvec);CHKERRQ(ierr);
802   ierr = PetscObjectCompose((PetscObject)*gvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr);
803   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void(*)(void))VecView_DMComposite);CHKERRQ(ierr);
804   PetscFunctionReturn(0);
805 }
806 
807 #undef __FUNCT__
808 #define __FUNCT__ "DMCreateLocalVector_Composite"
809 PetscErrorCode PETSCDM_DLLEXPORT DMCreateLocalVector_Composite(DM dm,Vec *lvec)
810 {
811   PetscErrorCode         ierr;
812   DM_Composite           *com = (DM_Composite*)dm->data;
813 
814   PetscFunctionBegin;
815   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
816   if (!com->setup) {
817     ierr = DMCompositeSetUp(dm);CHKERRQ(ierr);
818   }
819   ierr = VecCreateSeq(((PetscObject)dm)->comm,com->nghost,lvec);CHKERRQ(ierr);
820   ierr = PetscObjectCompose((PetscObject)*lvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr);
821   PetscFunctionReturn(0);
822 }
823 
824 #undef __FUNCT__
825 #define __FUNCT__ "DMCompositeGetLocalISs"
826 /*@C
827     DMCompositeGetLocalISs - gets an IS for each DM/array in the DMComposite, include ghost points
828 
829     Collective on DMComposite
830 
831     Input Parameter:
832 .    dm - the packer object
833 
834     Output Parameters:
835 .    is - the individual indices for each packed vector/array. Note that this includes
836            all the ghost points that individual ghosted DA's may have. Also each process has an
837            is for EACH redundant array (not just the local redundant arrays).
838 
839     Level: advanced
840 
841     Notes:
842        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
843 
844        Use DMCompositeGetGlobalISs() for non-ghosted ISs.
845 
846 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
847          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
848          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
849 
850 @*/
851 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetGlobalIndices(DM dm,IS *is[])
852 {
853   PetscErrorCode         ierr;
854   PetscInt               i,*idx,n,cnt;
855   struct DMCompositeLink *next;
856   Vec                    global,dglobal;
857   PF                     pf;
858   PetscScalar            *array;
859   PetscMPIInt            rank;
860   DM_Composite           *com = (DM_Composite*)dm->data;
861 
862   PetscFunctionBegin;
863   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
864   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
865   next = com->next;
866   ierr = DMCreateGlobalVector(dm,&global);CHKERRQ(ierr);
867   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
868 
869   /* put 0 to N-1 into the global vector */
870   ierr = PFCreate(PETSC_COMM_WORLD,1,1,&pf);CHKERRQ(ierr);
871   ierr = PFSetType(pf,PFIDENTITY,PETSC_NULL);CHKERRQ(ierr);
872   ierr = PFApplyVec(pf,PETSC_NULL,global);CHKERRQ(ierr);
873   ierr = PFDestroy(pf);CHKERRQ(ierr);
874 
875   /* loop over packed objects, handling one at at time */
876   cnt = 0;
877   while (next) {
878 
879     if (next->type == DMCOMPOSITE_ARRAY) {
880 
881       ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
882       if (rank == next->rank) {
883         ierr   = VecGetArray(global,&array);CHKERRQ(ierr);
884         array += next->rstart;
885         for (i=0; i<next->n; i++) idx[i] = (PetscInt)PetscRealPart(array[i]);
886         array -= next->rstart;
887         ierr   = VecRestoreArray(global,&array);CHKERRQ(ierr);
888       }
889       ierr = MPI_Bcast(idx,next->n,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
890       ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
891     } else if (next->type == DMCOMPOSITE_DM) {
892       Vec local;
893 
894       ierr   = DMCreateLocalVector(next->dm,&local);CHKERRQ(ierr);
895       ierr   = VecGetArray(global,&array);CHKERRQ(ierr);
896       array += next->rstart;
897       ierr   = DMGetGlobalVector(next->dm,&dglobal);CHKERRQ(ierr);
898       ierr   = VecPlaceArray(dglobal,array);CHKERRQ(ierr);
899       ierr   = DMGlobalToLocalBegin(next->dm,dglobal,INSERT_VALUES,local);CHKERRQ(ierr);
900       ierr   = DMGlobalToLocalEnd(next->dm,dglobal,INSERT_VALUES,local);CHKERRQ(ierr);
901       array -= next->rstart;
902       ierr   = VecRestoreArray(global,&array);CHKERRQ(ierr);
903       ierr   = VecResetArray(dglobal);CHKERRQ(ierr);
904       ierr   = DMRestoreGlobalVector(next->dm,&dglobal);CHKERRQ(ierr);
905 
906       ierr   = VecGetArray(local,&array);CHKERRQ(ierr);
907       ierr   = VecGetSize(local,&n);CHKERRQ(ierr);
908       ierr   = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
909       for (i=0; i<n; i++) idx[i] = (PetscInt)PetscRealPart(array[i]);
910       ierr    = VecRestoreArray(local,&array);CHKERRQ(ierr);
911       ierr    = VecDestroy(local);CHKERRQ(ierr);
912       ierr    = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
913 
914     } else SETERRQ(((PetscObject)global)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
915     next = next->next;
916     cnt++;
917   }
918   ierr = VecDestroy(global);CHKERRQ(ierr);
919   PetscFunctionReturn(0);
920 }
921 
922 #undef __FUNCT__
923 #define __FUNCT__ "DMCompositeGetGlobalISs"
924 /*@C
925     DMCompositeGetGlobalISs - Gets the index sets for each composed object
926 
927     Collective on DMComposite
928 
929     Input Parameter:
930 .    dm - the packer object
931 
932     Output Parameters:
933 .    is - the array of index sets
934 
935     Level: advanced
936 
937     Notes:
938        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
939 
940        The number of IS on each process will/may be different when redundant arrays are used
941 
942        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
943 
944        Use DMCompositeGetLocalISs() for index sets that include ghost points
945 
946 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
947          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
948          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
949 
950 @*/
951 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetGlobalISs(DM dm,IS *is[])
952 {
953   PetscErrorCode         ierr;
954   PetscInt               cnt = 0,*idx,i;
955   struct DMCompositeLink *next;
956   PetscMPIInt            rank;
957   DM_Composite           *com = (DM_Composite*)dm->data;
958 
959   PetscFunctionBegin;
960   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
961   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
962   next = com->next;
963   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
964 
965   /* loop over packed objects, handling one at at time */
966   while (next) {
967 
968     if (next->type == DMCOMPOSITE_ARRAY) {
969 
970       if (rank == next->rank) {
971         ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
972         for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
973         ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
974         cnt++;
975       }
976 
977     } else if (next->type == DMCOMPOSITE_DM) {
978       ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
979       for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
980       ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
981       cnt++;
982     } else {
983       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
984     }
985     next = next->next;
986   }
987   PetscFunctionReturn(0);
988 }
989 
990 /* -------------------------------------------------------------------------------------*/
991 #undef __FUNCT__
992 #define __FUNCT__ "DMCompositeGetLocalVectors_Array"
993 PetscErrorCode DMCompositeGetLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
994 {
995   PetscErrorCode ierr;
996   PetscFunctionBegin;
997   if (array) {
998     ierr = PetscMalloc(mine->n*sizeof(PetscScalar),array);CHKERRQ(ierr);
999   }
1000   PetscFunctionReturn(0);
1001 }
1002 
1003 #undef __FUNCT__
1004 #define __FUNCT__ "DMCompositeGetLocalVectors_DM"
1005 PetscErrorCode DMCompositeGetLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
1006 {
1007   PetscErrorCode ierr;
1008   PetscFunctionBegin;
1009   if (local) {
1010     ierr = DMGetLocalVector(mine->dm,local);CHKERRQ(ierr);
1011   }
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 #undef __FUNCT__
1016 #define __FUNCT__ "DMCompositeRestoreLocalVectors_Array"
1017 PetscErrorCode DMCompositeRestoreLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
1018 {
1019   PetscErrorCode ierr;
1020   PetscFunctionBegin;
1021   if (array) {
1022     ierr = PetscFree(*array);CHKERRQ(ierr);
1023   }
1024   PetscFunctionReturn(0);
1025 }
1026 
1027 #undef __FUNCT__
1028 #define __FUNCT__ "DMCompositeRestoreLocalVectors_DM"
1029 PetscErrorCode DMCompositeRestoreLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
1030 {
1031   PetscErrorCode ierr;
1032   PetscFunctionBegin;
1033   if (local) {
1034     ierr = DMRestoreLocalVector(mine->dm,local);CHKERRQ(ierr);
1035   }
1036   PetscFunctionReturn(0);
1037 }
1038 
1039 #undef __FUNCT__
1040 #define __FUNCT__ "DMCompositeGetLocalVectors"
1041 /*@C
1042     DMCompositeGetLocalVectors - Gets local vectors and arrays for each part of a DMComposite.
1043        Use DMCompositeRestoreLocalVectors() to return them.
1044 
1045     Not Collective
1046 
1047     Input Parameter:
1048 .    dm - the packer object
1049 
1050     Output Parameter:
1051 .    ... - the individual sequential objects (arrays or vectors)
1052 
1053     Level: advanced
1054 
1055 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
1056          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(),
1057          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
1058 
1059 @*/
1060 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetLocalVectors(DM dm,...)
1061 {
1062   va_list                Argp;
1063   PetscErrorCode         ierr;
1064   struct DMCompositeLink *next;
1065   DM_Composite           *com = (DM_Composite*)dm->data;
1066 
1067   PetscFunctionBegin;
1068   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1069   next = com->next;
1070   /* loop over packed objects, handling one at at time */
1071   va_start(Argp,dm);
1072   while (next) {
1073     if (next->type == DMCOMPOSITE_ARRAY) {
1074       PetscScalar **array;
1075       array = va_arg(Argp, PetscScalar**);
1076       ierr = DMCompositeGetLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
1077     } else if (next->type == DMCOMPOSITE_DM) {
1078       Vec *vec;
1079       vec = va_arg(Argp, Vec*);
1080       ierr = DMCompositeGetLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
1081     } else {
1082       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1083     }
1084     next = next->next;
1085   }
1086   va_end(Argp);
1087   PetscFunctionReturn(0);
1088 }
1089 
1090 #undef __FUNCT__
1091 #define __FUNCT__ "DMCompositeRestoreLocalVectors"
1092 /*@C
1093     DMCompositeRestoreLocalVectors - Restores local vectors and arrays for each part of a DMComposite.
1094 
1095     Not Collective
1096 
1097     Input Parameter:
1098 .    dm - the packer object
1099 
1100     Output Parameter:
1101 .    ... - the individual sequential objects (arrays or vectors)
1102 
1103     Level: advanced
1104 
1105 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
1106          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(),
1107          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
1108 
1109 @*/
1110 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreLocalVectors(DM dm,...)
1111 {
1112   va_list                Argp;
1113   PetscErrorCode         ierr;
1114   struct DMCompositeLink *next;
1115   DM_Composite           *com = (DM_Composite*)dm->data;
1116 
1117   PetscFunctionBegin;
1118   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1119   next = com->next;
1120   /* loop over packed objects, handling one at at time */
1121   va_start(Argp,dm);
1122   while (next) {
1123     if (next->type == DMCOMPOSITE_ARRAY) {
1124       PetscScalar **array;
1125       array = va_arg(Argp, PetscScalar**);
1126       ierr = DMCompositeRestoreLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
1127     } else if (next->type == DMCOMPOSITE_DM) {
1128       Vec *vec;
1129       vec = va_arg(Argp, Vec*);
1130       ierr = DMCompositeRestoreLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
1131     } else {
1132       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1133     }
1134     next = next->next;
1135   }
1136   va_end(Argp);
1137   PetscFunctionReturn(0);
1138 }
1139 
1140 /* -------------------------------------------------------------------------------------*/
1141 #undef __FUNCT__
1142 #define __FUNCT__ "DMCompositeGetEntries_Array"
1143 PetscErrorCode DMCompositeGetEntries_Array(DM dm,struct DMCompositeLink *mine,PetscInt *n)
1144 {
1145   PetscFunctionBegin;
1146   if (n) *n = mine->n;
1147   PetscFunctionReturn(0);
1148 }
1149 
1150 #undef __FUNCT__
1151 #define __FUNCT__ "DMCompositeGetEntries_DM"
1152 PetscErrorCode DMCompositeGetEntries_DM(DM dmi,struct DMCompositeLink *mine,DM *dm)
1153 {
1154   PetscFunctionBegin;
1155   if (dm) *dm = mine->dm;
1156   PetscFunctionReturn(0);
1157 }
1158 
1159 #undef __FUNCT__
1160 #define __FUNCT__ "DMCompositeGetEntries"
1161 /*@C
1162     DMCompositeGetEntries - Gets the DA, redundant size, etc for each entry in a DMComposite.
1163 
1164     Not Collective
1165 
1166     Input Parameter:
1167 .    dm - the packer object
1168 
1169     Output Parameter:
1170 .    ... - the individual entries, DAs or integer sizes)
1171 
1172     Level: advanced
1173 
1174 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
1175          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(),
1176          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
1177          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
1178 
1179 @*/
1180 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetEntries(DM dm,...)
1181 {
1182   va_list                Argp;
1183   PetscErrorCode         ierr;
1184   struct DMCompositeLink *next;
1185   DM_Composite           *com = (DM_Composite*)dm->data;
1186 
1187   PetscFunctionBegin;
1188   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1189   next = com->next;
1190   /* loop over packed objects, handling one at at time */
1191   va_start(Argp,dm);
1192   while (next) {
1193     if (next->type == DMCOMPOSITE_ARRAY) {
1194       PetscInt *n;
1195       n = va_arg(Argp, PetscInt*);
1196       ierr = DMCompositeGetEntries_Array(dm,next,n);CHKERRQ(ierr);
1197     } else if (next->type == DMCOMPOSITE_DM) {
1198       DM *dmn;
1199       dmn = va_arg(Argp, DM*);
1200       ierr = DMCompositeGetEntries_DM(dm,next,dmn);CHKERRQ(ierr);
1201     } else {
1202       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1203     }
1204     next = next->next;
1205   }
1206   va_end(Argp);
1207   PetscFunctionReturn(0);
1208 }
1209 
1210 #undef __FUNCT__
1211 #define __FUNCT__ "DMRefine_Composite"
1212 PetscErrorCode PETSCDM_DLLEXPORT DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
1213 {
1214   PetscErrorCode         ierr;
1215   struct DMCompositeLink *next;
1216   DM_Composite           *com = (DM_Composite*)dmi->data;
1217   DM                     dm;
1218 
1219   PetscFunctionBegin;
1220   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1221   next = com->next;
1222   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
1223 
1224   /* loop over packed objects, handling one at at time */
1225   while (next) {
1226     if (next->type == DMCOMPOSITE_ARRAY) {
1227       ierr = DMCompositeAddArray(*fine,next->rank,next->n);CHKERRQ(ierr);
1228     } else if (next->type == DMCOMPOSITE_DM) {
1229       ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
1230       ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
1231       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1232     } else {
1233       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1234     }
1235     next = next->next;
1236   }
1237   PetscFunctionReturn(0);
1238 }
1239 
1240 #include "petscmat.h"
1241 
1242 struct MatPackLink {
1243   Mat                A;
1244   struct MatPackLink *next;
1245 };
1246 
1247 struct MatPack {
1248   DM                 right,left;
1249   struct MatPackLink *next;
1250 };
1251 
1252 #undef __FUNCT__
1253 #define __FUNCT__ "MatMultBoth_Shell_Pack"
1254 PetscErrorCode MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscBool  add)
1255 {
1256   struct MatPack         *mpack;
1257   struct DMCompositeLink *xnext,*ynext;
1258   struct MatPackLink     *anext;
1259   PetscScalar            *xarray,*yarray;
1260   PetscErrorCode         ierr;
1261   PetscInt               i;
1262   Vec                    xglobal,yglobal;
1263   PetscMPIInt            rank;
1264   DM_Composite           *comright;
1265   DM_Composite           *comleft;
1266 
1267   PetscFunctionBegin;
1268   ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
1269   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
1270   comright = (DM_Composite*)mpack->right->data;
1271   comleft = (DM_Composite*)mpack->left->data;
1272   xnext = comright->next;
1273   ynext = comleft->next;
1274   anext = mpack->next;
1275 
1276   while (xnext) {
1277     if (xnext->type == DMCOMPOSITE_ARRAY) {
1278       if (rank == xnext->rank) {
1279         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
1280         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
1281         if (add) {
1282           for (i=0; i<xnext->n; i++) {
1283             yarray[ynext->rstart+i] += xarray[xnext->rstart+i];
1284           }
1285         } else {
1286           ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
1287         }
1288         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1289         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1290       }
1291     } else if (xnext->type == DMCOMPOSITE_DM) {
1292       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
1293       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
1294       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
1295       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
1296       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
1297       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
1298       if (add) {
1299         ierr  = MatMultAdd(anext->A,xglobal,yglobal,yglobal);CHKERRQ(ierr);
1300       } else {
1301         ierr  = MatMult(anext->A,xglobal,yglobal);CHKERRQ(ierr);
1302       }
1303       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1304       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1305       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
1306       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
1307       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
1308       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
1309       anext = anext->next;
1310     } else {
1311       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1312     }
1313     xnext = xnext->next;
1314     ynext = ynext->next;
1315   }
1316   PetscFunctionReturn(0);
1317 }
1318 
1319 #undef __FUNCT__
1320 #define __FUNCT__ "MatMultAdd_Shell_Pack"
1321 PetscErrorCode MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z)
1322 {
1323   PetscErrorCode ierr;
1324   PetscFunctionBegin;
1325   if (z != y) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Handles y == z only");
1326   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
1327   PetscFunctionReturn(0);
1328 }
1329 
1330 #undef __FUNCT__
1331 #define __FUNCT__ "MatMult_Shell_Pack"
1332 PetscErrorCode MatMult_Shell_Pack(Mat A,Vec x,Vec y)
1333 {
1334   PetscErrorCode ierr;
1335   PetscFunctionBegin;
1336   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
1337   PetscFunctionReturn(0);
1338 }
1339 
1340 #undef __FUNCT__
1341 #define __FUNCT__ "MatMultTranspose_Shell_Pack"
1342 PetscErrorCode MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y)
1343 {
1344   struct MatPack         *mpack;
1345   struct DMCompositeLink *xnext,*ynext;
1346   struct MatPackLink     *anext;
1347   PetscScalar            *xarray,*yarray;
1348   PetscErrorCode         ierr;
1349   Vec                    xglobal,yglobal;
1350   PetscMPIInt            rank;
1351   DM_Composite           *comright;
1352   DM_Composite           *comleft;
1353 
1354   PetscFunctionBegin;
1355   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
1356   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
1357   comright = (DM_Composite*)mpack->right->data;
1358   comleft = (DM_Composite*)mpack->left->data;
1359   ynext = comright->next;
1360   xnext = comleft->next;
1361   anext = mpack->next;
1362 
1363   while (xnext) {
1364     if (xnext->type == DMCOMPOSITE_ARRAY) {
1365       if (rank == ynext->rank) {
1366         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
1367         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
1368         ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
1369         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1370         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1371       }
1372     } else if (xnext->type == DMCOMPOSITE_DM) {
1373       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
1374       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
1375       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
1376       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
1377       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
1378       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
1379       ierr  = MatMultTranspose(anext->A,xglobal,yglobal);CHKERRQ(ierr);
1380       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1381       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1382       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
1383       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
1384       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
1385       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
1386       anext = anext->next;
1387     } else {
1388       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1389     }
1390     xnext = xnext->next;
1391     ynext = ynext->next;
1392   }
1393   PetscFunctionReturn(0);
1394 }
1395 
1396 #undef __FUNCT__
1397 #define __FUNCT__ "MatDestroy_Shell_Pack"
1398 PetscErrorCode MatDestroy_Shell_Pack(Mat A)
1399 {
1400   struct MatPack     *mpack;
1401   struct MatPackLink *anext,*oldanext;
1402   PetscErrorCode     ierr;
1403 
1404   PetscFunctionBegin;
1405   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
1406   anext = mpack->next;
1407 
1408   while (anext) {
1409     ierr     = MatDestroy(anext->A);CHKERRQ(ierr);
1410     oldanext = anext;
1411     anext    = anext->next;
1412     ierr     = PetscFree(oldanext);CHKERRQ(ierr);
1413   }
1414   ierr = PetscFree(mpack);CHKERRQ(ierr);
1415   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
1416   PetscFunctionReturn(0);
1417 }
1418 
1419 #undef __FUNCT__
1420 #define __FUNCT__ "DMGetInterpolation_Composite"
1421 PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
1422 {
1423   PetscErrorCode         ierr;
1424   PetscInt               m,n,M,N;
1425   struct DMCompositeLink *nextc;
1426   struct DMCompositeLink *nextf;
1427   struct MatPackLink     *nextmat,*pnextmat = 0;
1428   struct MatPack         *mpack;
1429   Vec                    gcoarse,gfine;
1430   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
1431   DM_Composite           *comfine = (DM_Composite*)fine->data;
1432 
1433   PetscFunctionBegin;
1434   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
1435   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
1436   nextc = comcoarse->next;
1437   nextf = comfine->next;
1438   /* use global vectors only for determining matrix layout */
1439   ierr = DMCreateGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
1440   ierr = DMCreateGlobalVector(fine,&gfine);CHKERRQ(ierr);
1441   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
1442   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
1443   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
1444   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
1445   ierr = VecDestroy(gcoarse);CHKERRQ(ierr);
1446   ierr = VecDestroy(gfine);CHKERRQ(ierr);
1447 
1448   ierr         = PetscNew(struct MatPack,&mpack);CHKERRQ(ierr);
1449   mpack->right = coarse;
1450   mpack->left  = fine;
1451   ierr  = MatCreate(((PetscObject)fine)->comm,A);CHKERRQ(ierr);
1452   ierr  = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1453   ierr  = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
1454   ierr  = MatShellSetContext(*A,mpack);CHKERRQ(ierr);
1455   ierr  = MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);CHKERRQ(ierr);
1456   ierr  = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);CHKERRQ(ierr);
1457   ierr  = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);CHKERRQ(ierr);
1458   ierr  = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);CHKERRQ(ierr);
1459 
1460   /* loop over packed objects, handling one at at time */
1461   while (nextc) {
1462     if (nextc->type != nextf->type) SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Two DM have different layout");
1463 
1464     if (nextc->type == DMCOMPOSITE_ARRAY) {
1465       ;
1466     } else if (nextc->type == DMCOMPOSITE_DM) {
1467       ierr          = PetscNew(struct MatPackLink,&nextmat);CHKERRQ(ierr);
1468       nextmat->next = 0;
1469       if (pnextmat) {
1470         pnextmat->next = nextmat;
1471         pnextmat       = nextmat;
1472       } else {
1473         pnextmat    = nextmat;
1474         mpack->next = nextmat;
1475       }
1476       ierr = DMGetInterpolation(nextc->dm,nextf->dm,&nextmat->A,PETSC_NULL);CHKERRQ(ierr);
1477     } else {
1478       SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1479     }
1480     nextc = nextc->next;
1481     nextf = nextf->next;
1482   }
1483   PetscFunctionReturn(0);
1484 }
1485 
1486 #undef __FUNCT__
1487 #define __FUNCT__ "DMGetMatrix_Composite"
1488 PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix_Composite(DM dm, const MatType mtype,Mat *J)
1489 {
1490   PetscErrorCode         ierr;
1491   DM_Composite           *com = (DM_Composite*)dm->data;
1492   struct DMCompositeLink *next = com->next;
1493   PetscInt               m,*dnz,*onz,i,j,mA;
1494   Mat                    Atmp;
1495   PetscMPIInt            rank;
1496   PetscScalar            zero = 0.0;
1497   PetscBool              dense = PETSC_FALSE;
1498 
1499   PetscFunctionBegin;
1500   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1501 
1502   /* use global vector to determine layout needed for matrix */
1503   m = com->n;
1504   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1505   ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr);
1506   ierr = MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1507   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
1508 
1509   /*
1510      Extremely inefficient but will compute entire Jacobian for testing
1511   */
1512   ierr = PetscOptionsGetTruth(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
1513   if (dense) {
1514     PetscInt    rstart,rend,*indices;
1515     PetscScalar *values;
1516 
1517     mA    = com->N;
1518     ierr = MatMPIAIJSetPreallocation(*J,mA,PETSC_NULL,mA-m,PETSC_NULL);CHKERRQ(ierr);
1519     ierr = MatSeqAIJSetPreallocation(*J,mA,PETSC_NULL);CHKERRQ(ierr);
1520 
1521     ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr);
1522     ierr = PetscMalloc2(mA,PetscScalar,&values,mA,PetscInt,&indices);CHKERRQ(ierr);
1523     ierr = PetscMemzero(values,mA*sizeof(PetscScalar));CHKERRQ(ierr);
1524     for (i=0; i<mA; i++) indices[i] = i;
1525     for (i=rstart; i<rend; i++) {
1526       ierr = MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES);CHKERRQ(ierr);
1527     }
1528     ierr = PetscFree2(values,indices);CHKERRQ(ierr);
1529     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1530     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1531     PetscFunctionReturn(0);
1532   }
1533 
1534   ierr = MatPreallocateInitialize(((PetscObject)dm)->comm,m,m,dnz,onz);CHKERRQ(ierr);
1535   /* loop over packed objects, handling one at at time */
1536   next = com->next;
1537   while (next) {
1538     if (next->type == DMCOMPOSITE_ARRAY) {
1539       if (rank == next->rank) {  /* zero the "little" block */
1540         for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) {
1541           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
1542             ierr = MatPreallocateSet(j,1,&i,dnz,onz);CHKERRQ(ierr);
1543           }
1544         }
1545       }
1546     } else if (next->type == DMCOMPOSITE_DM) {
1547       PetscInt       nc,rstart,*ccols,maxnc;
1548       const PetscInt *cols,*rstarts;
1549       PetscMPIInt    proc;
1550 
1551       ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr);
1552       ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr);
1553       ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
1554       ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr);
1555 
1556       maxnc = 0;
1557       for (i=0; i<mA; i++) {
1558         ierr  = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1559         ierr  = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1560         maxnc = PetscMax(nc,maxnc);
1561       }
1562       ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
1563       for (i=0; i<mA; i++) {
1564         ierr = MatGetRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr);
1565         /* remap the columns taking into how much they are shifted on each process */
1566         for (j=0; j<nc; j++) {
1567           proc = 0;
1568           while (cols[j] >= rstarts[proc+1]) proc++;
1569           ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
1570         }
1571         ierr = MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz);CHKERRQ(ierr);
1572         ierr = MatRestoreRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr);
1573       }
1574       ierr = PetscFree(ccols);CHKERRQ(ierr);
1575       ierr = MatDestroy(Atmp);CHKERRQ(ierr);
1576     } else {
1577       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1578     }
1579     next = next->next;
1580   }
1581   if (com->FormCoupleLocations) {
1582     ierr = (*com->FormCoupleLocations)(dm,PETSC_NULL,dnz,onz,__rstart,__nrows,__start,__end);CHKERRQ(ierr);
1583   }
1584   ierr = MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);CHKERRQ(ierr);
1585   ierr = MatSeqAIJSetPreallocation(*J,0,dnz);CHKERRQ(ierr);
1586   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1587 
1588   next = com->next;
1589   while (next) {
1590     if (next->type == DMCOMPOSITE_ARRAY) {
1591       if (rank == next->rank) {
1592         for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) {
1593           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
1594             ierr = MatSetValues(*J,1,&j,1,&i,&zero,INSERT_VALUES);CHKERRQ(ierr);
1595           }
1596         }
1597       }
1598     } else if (next->type == DMCOMPOSITE_DM) {
1599       PetscInt          nc,rstart,row,maxnc,*ccols;
1600       const PetscInt    *cols,*rstarts;
1601       const PetscScalar *values;
1602       PetscMPIInt       proc;
1603 
1604       ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr);
1605       ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr);
1606       ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
1607       ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr);
1608       maxnc = 0;
1609       for (i=0; i<mA; i++) {
1610         ierr  = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1611         ierr  = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1612         maxnc = PetscMax(nc,maxnc);
1613       }
1614       ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
1615       for (i=0; i<mA; i++) {
1616         ierr = MatGetRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr);
1617         for (j=0; j<nc; j++) {
1618           proc = 0;
1619           while (cols[j] >= rstarts[proc+1]) proc++;
1620           ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
1621         }
1622         row  = com->rstart+next->rstart+i;
1623         ierr = MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES);CHKERRQ(ierr);
1624         ierr = MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr);
1625       }
1626       ierr = PetscFree(ccols);CHKERRQ(ierr);
1627       ierr = MatDestroy(Atmp);CHKERRQ(ierr);
1628     } else {
1629       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1630     }
1631     next = next->next;
1632   }
1633   if (com->FormCoupleLocations) {
1634     PetscInt __rstart;
1635     ierr = MatGetOwnershipRange(*J,&__rstart,PETSC_NULL);CHKERRQ(ierr);
1636     ierr = (*com->FormCoupleLocations)(dm,*J,PETSC_NULL,PETSC_NULL,__rstart,0,0,0);CHKERRQ(ierr);
1637   }
1638   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1639   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1640   PetscFunctionReturn(0);
1641 }
1642 
1643 #undef __FUNCT__
1644 #define __FUNCT__ "DMGetColoring_Composite"
1645 PetscErrorCode PETSCDM_DLLEXPORT DMGetColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
1646 {
1647   PetscErrorCode         ierr;
1648   PetscInt               n,i,cnt;
1649   ISColoringValue        *colors;
1650   PetscBool              dense = PETSC_FALSE;
1651   ISColoringValue        maxcol = 0;
1652   DM_Composite           *com = (DM_Composite*)dm->data;
1653 
1654   PetscFunctionBegin;
1655   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1656   if (ctype == IS_COLORING_GHOSTED) {
1657     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" );
1658   } else if (ctype == IS_COLORING_GLOBAL) {
1659     n = com->n;
1660   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1661   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
1662 
1663   ierr = PetscOptionsGetTruth(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
1664   if (dense) {
1665     for (i=0; i<n; i++) {
1666       colors[i] = (ISColoringValue)(com->rstart + i);
1667     }
1668     maxcol = com->N;
1669   } else {
1670     struct DMCompositeLink *next = com->next;
1671     PetscMPIInt            rank;
1672 
1673     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1674     cnt  = 0;
1675     while (next) {
1676       if (next->type == DMCOMPOSITE_ARRAY) {
1677         if (rank == next->rank) {  /* each column gets is own color */
1678           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
1679             colors[cnt++] = maxcol++;
1680           }
1681         }
1682         ierr = MPI_Bcast(&maxcol,1,MPIU_COLORING_VALUE,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
1683       } else if (next->type == DMCOMPOSITE_DM) {
1684         ISColoring     lcoloring;
1685 
1686         ierr = DMGetColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
1687         for (i=0; i<lcoloring->N; i++) {
1688           colors[cnt++] = maxcol + lcoloring->colors[i];
1689         }
1690         maxcol += lcoloring->n;
1691         ierr = ISColoringDestroy(lcoloring);CHKERRQ(ierr);
1692       } else {
1693         SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1694       }
1695       next = next->next;
1696     }
1697   }
1698   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
1699   PetscFunctionReturn(0);
1700 }
1701 
1702 #undef __FUNCT__
1703 #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
1704 PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1705 {
1706   PetscErrorCode         ierr;
1707   struct DMCompositeLink *next;
1708   PetscInt               cnt = 3;
1709   PetscMPIInt            rank;
1710   PetscScalar            *garray,*larray;
1711   DM_Composite           *com = (DM_Composite*)dm->data;
1712 
1713   PetscFunctionBegin;
1714   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1715   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
1716   next = com->next;
1717   if (!com->setup) {
1718     ierr = DMCompositeSetUp(dm);CHKERRQ(ierr);
1719   }
1720   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1721   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
1722   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
1723 
1724   /* loop over packed objects, handling one at at time */
1725   while (next) {
1726     if (next->type == DMCOMPOSITE_ARRAY) {
1727       if (rank == next->rank) {
1728         ierr    = PetscMemcpy(larray,garray,next->n*sizeof(PetscScalar));CHKERRQ(ierr);
1729         garray += next->n;
1730       }
1731       /* does not handle ADD_VALUES */
1732       ierr = MPI_Bcast(larray,next->n,MPIU_SCALAR,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
1733       larray += next->n;
1734     } else if (next->type == DMCOMPOSITE_DM) {
1735       Vec      local,global;
1736       PetscInt N;
1737 
1738       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
1739       ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
1740       ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
1741       ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
1742       ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
1743       ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
1744       ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
1745       ierr = VecResetArray(global);CHKERRQ(ierr);
1746       ierr = VecResetArray(local);CHKERRQ(ierr);
1747       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
1748       ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
1749       larray += next->n;
1750     } else {
1751       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1752     }
1753     cnt++;
1754     next = next->next;
1755   }
1756 
1757   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
1758   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
1759   PetscFunctionReturn(0);
1760 }
1761 
1762 #undef __FUNCT__
1763 #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
1764 PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1765 {
1766   PetscFunctionBegin;
1767   PetscFunctionReturn(0);
1768 }
1769 
1770 #undef __FUNCT__
1771 #define __FUNCT__ "DMCompositeCreate"
1772 /*@C
1773     DMCompositeCreate - Creates a vector packer, used to generate "composite"
1774       vectors made up of several subvectors.
1775 
1776     Collective on MPI_Comm
1777 
1778     Input Parameter:
1779 .   comm - the processors that will share the global vector
1780 
1781     Output Parameters:
1782 .   packer - the packer object
1783 
1784     Level: advanced
1785 
1786 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCompositeScatter(),
1787          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess()
1788          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
1789 
1790 @*/
1791 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeCreate(MPI_Comm comm,DM *packer)
1792 {
1793   PetscErrorCode ierr;
1794   DM             p;
1795   DM_Composite   *com;
1796 
1797   PetscFunctionBegin;
1798   PetscValidPointer(packer,2);
1799   *packer = PETSC_NULL;
1800 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
1801   ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1802 #endif
1803 
1804   ierr = PetscHeaderCreate(p,_p_DM,struct _DMOps,DM_CLASSID,0,"DM",comm,DMDestroy,DMView);CHKERRQ(ierr);
1805   ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1806   p->data = com;
1807   ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1808   com->n            = 0;
1809   com->next         = PETSC_NULL;
1810   com->nredundant   = 0;
1811   com->nDM          = 0;
1812 
1813   p->ops->createglobalvector = DMCreateGlobalVector_Composite;
1814   p->ops->createlocalvector  = DMCreateLocalVector_Composite;
1815   p->ops->refine             = DMRefine_Composite;
1816   p->ops->getinterpolation   = DMGetInterpolation_Composite;
1817   p->ops->getmatrix          = DMGetMatrix_Composite;
1818   p->ops->getcoloring        = DMGetColoring_Composite;
1819   p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite;
1820   p->ops->globaltolocalend   = DMGlobalToLocalEnd_Composite;
1821   p->ops->destroy            = DMDestroy_Composite;
1822   p->ops->view               = DMView_Composite;
1823 
1824   *packer = p;
1825   PetscFunctionReturn(0);
1826 }
1827