xref: /petsc/src/dm/impls/composite/pack.c (revision efe48dd8cf8b935accbbb9f4bcb20bc83865fa4d)
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(), DMCompositeGetGlobalIndices(), 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(), DMCompositeGetGlobalIndices(), 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(), DMCompositeGetGlobalIndices(), 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(), DMCompositeGetGlobalIndices(), 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(), DMCompositeGetGlobalIndices(), 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__ "DMCompositeGetLocalISs"
836 /*@C
837     DMCompositeGetLocalISs - gets an IS for each DM/array in the DMComposite, include ghost points
838 
839     Collective on DMComposite
840 
841     Input Parameter:
842 .    dm - the packer object
843 
844     Output Parameters:
845 .    is - the individual indices 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            is for EACH redundant array (not just the local redundant arrays).
848 
849     Level: advanced
850 
851     Notes:
852        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
853 
854        Use DMCompositeGetGlobalISs() for non-ghosted ISs.
855 
856 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
857          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
858          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
859 
860 @*/
861 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetGlobalIndices(DM dm,IS *is[])
862 {
863   PetscErrorCode         ierr;
864   PetscInt               i,*idx,n,cnt;
865   struct DMCompositeLink *next;
866   Vec                    global,dglobal;
867   PF                     pf;
868   PetscScalar            *array;
869   PetscMPIInt            rank;
870   DM_Composite           *com = (DM_Composite*)dm->data;
871 
872   PetscFunctionBegin;
873   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
874   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
875   next = com->next;
876   ierr = DMCreateGlobalVector(dm,&global);CHKERRQ(ierr);
877   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
878 
879   /* put 0 to N-1 into the global vector */
880   ierr = PFCreate(PETSC_COMM_WORLD,1,1,&pf);CHKERRQ(ierr);
881   ierr = PFSetType(pf,PFIDENTITY,PETSC_NULL);CHKERRQ(ierr);
882   ierr = PFApplyVec(pf,PETSC_NULL,global);CHKERRQ(ierr);
883   ierr = PFDestroy(pf);CHKERRQ(ierr);
884 
885   /* loop over packed objects, handling one at at time */
886   cnt = 0;
887   while (next) {
888 
889     if (next->type == DMCOMPOSITE_ARRAY) {
890 
891       ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
892       if (rank == next->rank) {
893         ierr   = VecGetArray(global,&array);CHKERRQ(ierr);
894         array += next->rstart;
895         for (i=0; i<next->n; i++) idx[i] = (PetscInt)PetscRealPart(array[i]);
896         array -= next->rstart;
897         ierr   = VecRestoreArray(global,&array);CHKERRQ(ierr);
898       }
899       ierr = MPI_Bcast(idx,next->n,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
900       ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
901     } else if (next->type == DMCOMPOSITE_DM) {
902       Vec local;
903 
904       ierr   = DMCreateLocalVector(next->dm,&local);CHKERRQ(ierr);
905       ierr   = VecGetArray(global,&array);CHKERRQ(ierr);
906       array += next->rstart;
907       ierr   = DMGetGlobalVector(next->dm,&dglobal);CHKERRQ(ierr);
908       ierr   = VecPlaceArray(dglobal,array);CHKERRQ(ierr);
909       ierr   = DMGlobalToLocalBegin(next->dm,dglobal,INSERT_VALUES,local);CHKERRQ(ierr);
910       ierr   = DMGlobalToLocalEnd(next->dm,dglobal,INSERT_VALUES,local);CHKERRQ(ierr);
911       array -= next->rstart;
912       ierr   = VecRestoreArray(global,&array);CHKERRQ(ierr);
913       ierr   = VecResetArray(dglobal);CHKERRQ(ierr);
914       ierr   = DMRestoreGlobalVector(next->dm,&dglobal);CHKERRQ(ierr);
915 
916       ierr   = VecGetArray(local,&array);CHKERRQ(ierr);
917       ierr   = VecGetSize(local,&n);CHKERRQ(ierr);
918       ierr   = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
919       for (i=0; i<n; i++) idx[i] = (PetscInt)PetscRealPart(array[i]);
920       ierr    = VecRestoreArray(local,&array);CHKERRQ(ierr);
921       ierr    = VecDestroy(local);CHKERRQ(ierr);
922       ierr    = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
923 
924     } else SETERRQ(((PetscObject)global)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
925     next = next->next;
926     cnt++;
927   }
928   ierr = VecDestroy(global);CHKERRQ(ierr);
929   PetscFunctionReturn(0);
930 }
931 
932 #undef __FUNCT__
933 #define __FUNCT__ "DMCompositeGetGlobalISs"
934 /*@C
935     DMCompositeGetGlobalISs - Gets the index sets for each composed object
936 
937     Collective on DMComposite
938 
939     Input Parameter:
940 .    dm - the packer object
941 
942     Output Parameters:
943 .    is - the array of index sets
944 
945     Level: advanced
946 
947     Notes:
948        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
949 
950        The number of IS on each process will/may be different when redundant arrays are used
951 
952        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
953 
954        Use DMCompositeGetLocalISs() for index sets that include ghost points
955 
956 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
957          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
958          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
959 
960 @*/
961 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetGlobalISs(DM dm,IS *is[])
962 {
963   PetscErrorCode         ierr;
964   PetscInt               cnt = 0,*idx,i;
965   struct DMCompositeLink *next;
966   PetscMPIInt            rank;
967   DM_Composite           *com = (DM_Composite*)dm->data;
968 
969   PetscFunctionBegin;
970   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
971   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
972   next = com->next;
973   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
974 
975   /* loop over packed objects, handling one at at time */
976   while (next) {
977 
978     if (next->type == DMCOMPOSITE_ARRAY) {
979 
980       if (rank == next->rank) {
981         ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
982         for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
983         ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
984         cnt++;
985       }
986 
987     } else if (next->type == DMCOMPOSITE_DM) {
988       ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
989       for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
990       ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
991       cnt++;
992     } else {
993       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
994     }
995     next = next->next;
996   }
997   PetscFunctionReturn(0);
998 }
999 
1000 /* -------------------------------------------------------------------------------------*/
1001 #undef __FUNCT__
1002 #define __FUNCT__ "DMCompositeGetLocalVectors_Array"
1003 PetscErrorCode DMCompositeGetLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
1004 {
1005   PetscErrorCode ierr;
1006   PetscFunctionBegin;
1007   if (array) {
1008     ierr = PetscMalloc(mine->n*sizeof(PetscScalar),array);CHKERRQ(ierr);
1009   }
1010   PetscFunctionReturn(0);
1011 }
1012 
1013 #undef __FUNCT__
1014 #define __FUNCT__ "DMCompositeGetLocalVectors_DM"
1015 PetscErrorCode DMCompositeGetLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
1016 {
1017   PetscErrorCode ierr;
1018   PetscFunctionBegin;
1019   if (local) {
1020     ierr = DMGetLocalVector(mine->dm,local);CHKERRQ(ierr);
1021   }
1022   PetscFunctionReturn(0);
1023 }
1024 
1025 #undef __FUNCT__
1026 #define __FUNCT__ "DMCompositeRestoreLocalVectors_Array"
1027 PetscErrorCode DMCompositeRestoreLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
1028 {
1029   PetscErrorCode ierr;
1030   PetscFunctionBegin;
1031   if (array) {
1032     ierr = PetscFree(*array);CHKERRQ(ierr);
1033   }
1034   PetscFunctionReturn(0);
1035 }
1036 
1037 #undef __FUNCT__
1038 #define __FUNCT__ "DMCompositeRestoreLocalVectors_DM"
1039 PetscErrorCode DMCompositeRestoreLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
1040 {
1041   PetscErrorCode ierr;
1042   PetscFunctionBegin;
1043   if (local) {
1044     ierr = DMRestoreLocalVector(mine->dm,local);CHKERRQ(ierr);
1045   }
1046   PetscFunctionReturn(0);
1047 }
1048 
1049 #undef __FUNCT__
1050 #define __FUNCT__ "DMCompositeGetLocalVectors"
1051 /*@C
1052     DMCompositeGetLocalVectors - Gets local vectors and arrays for each part of a DMComposite.
1053        Use DMCompositeRestoreLocalVectors() to return them.
1054 
1055     Not Collective
1056 
1057     Input Parameter:
1058 .    dm - the packer object
1059 
1060     Output Parameter:
1061 .    ... - the individual sequential objects (arrays or vectors)
1062 
1063     Level: advanced
1064 
1065 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
1066          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(),
1067          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
1068 
1069 @*/
1070 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetLocalVectors(DM dm,...)
1071 {
1072   va_list                Argp;
1073   PetscErrorCode         ierr;
1074   struct DMCompositeLink *next;
1075   DM_Composite           *com = (DM_Composite*)dm->data;
1076 
1077   PetscFunctionBegin;
1078   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1079   next = com->next;
1080   /* loop over packed objects, handling one at at time */
1081   va_start(Argp,dm);
1082   while (next) {
1083     if (next->type == DMCOMPOSITE_ARRAY) {
1084       PetscScalar **array;
1085       array = va_arg(Argp, PetscScalar**);
1086       ierr = DMCompositeGetLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
1087     } else if (next->type == DMCOMPOSITE_DM) {
1088       Vec *vec;
1089       vec = va_arg(Argp, Vec*);
1090       ierr = DMCompositeGetLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
1091     } else {
1092       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1093     }
1094     next = next->next;
1095   }
1096   va_end(Argp);
1097   PetscFunctionReturn(0);
1098 }
1099 
1100 #undef __FUNCT__
1101 #define __FUNCT__ "DMCompositeRestoreLocalVectors"
1102 /*@C
1103     DMCompositeRestoreLocalVectors - Restores local vectors and arrays for each part of a DMComposite.
1104 
1105     Not Collective
1106 
1107     Input Parameter:
1108 .    dm - the packer object
1109 
1110     Output Parameter:
1111 .    ... - the individual sequential objects (arrays or vectors)
1112 
1113     Level: advanced
1114 
1115 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
1116          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(),
1117          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
1118 
1119 @*/
1120 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreLocalVectors(DM dm,...)
1121 {
1122   va_list                Argp;
1123   PetscErrorCode         ierr;
1124   struct DMCompositeLink *next;
1125   DM_Composite           *com = (DM_Composite*)dm->data;
1126 
1127   PetscFunctionBegin;
1128   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1129   next = com->next;
1130   /* loop over packed objects, handling one at at time */
1131   va_start(Argp,dm);
1132   while (next) {
1133     if (next->type == DMCOMPOSITE_ARRAY) {
1134       PetscScalar **array;
1135       array = va_arg(Argp, PetscScalar**);
1136       ierr = DMCompositeRestoreLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
1137     } else if (next->type == DMCOMPOSITE_DM) {
1138       Vec *vec;
1139       vec = va_arg(Argp, Vec*);
1140       ierr = DMCompositeRestoreLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
1141     } else {
1142       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1143     }
1144     next = next->next;
1145   }
1146   va_end(Argp);
1147   PetscFunctionReturn(0);
1148 }
1149 
1150 /* -------------------------------------------------------------------------------------*/
1151 #undef __FUNCT__
1152 #define __FUNCT__ "DMCompositeGetEntries_Array"
1153 PetscErrorCode DMCompositeGetEntries_Array(DM dm,struct DMCompositeLink *mine,PetscInt *n)
1154 {
1155   PetscFunctionBegin;
1156   if (n) *n = mine->n;
1157   PetscFunctionReturn(0);
1158 }
1159 
1160 #undef __FUNCT__
1161 #define __FUNCT__ "DMCompositeGetEntries_DM"
1162 PetscErrorCode DMCompositeGetEntries_DM(DM dmi,struct DMCompositeLink *mine,DM *dm)
1163 {
1164   PetscFunctionBegin;
1165   if (dm) *dm = mine->dm;
1166   PetscFunctionReturn(0);
1167 }
1168 
1169 #undef __FUNCT__
1170 #define __FUNCT__ "DMCompositeGetEntries"
1171 /*@C
1172     DMCompositeGetEntries - Gets the DM, redundant size, etc for each entry in a DMComposite.
1173 
1174     Not Collective
1175 
1176     Input Parameter:
1177 .    dm - the packer object
1178 
1179     Output Parameter:
1180 .    ... - the individual entries, DMs or integer sizes)
1181 
1182     Level: advanced
1183 
1184 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
1185          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(),
1186          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
1187          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
1188 
1189 @*/
1190 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetEntries(DM dm,...)
1191 {
1192   va_list                Argp;
1193   PetscErrorCode         ierr;
1194   struct DMCompositeLink *next;
1195   DM_Composite           *com = (DM_Composite*)dm->data;
1196 
1197   PetscFunctionBegin;
1198   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1199   next = com->next;
1200   /* loop over packed objects, handling one at at time */
1201   va_start(Argp,dm);
1202   while (next) {
1203     if (next->type == DMCOMPOSITE_ARRAY) {
1204       PetscInt *n;
1205       n = va_arg(Argp, PetscInt*);
1206       ierr = DMCompositeGetEntries_Array(dm,next,n);CHKERRQ(ierr);
1207     } else if (next->type == DMCOMPOSITE_DM) {
1208       DM *dmn;
1209       dmn = va_arg(Argp, DM*);
1210       ierr = DMCompositeGetEntries_DM(dm,next,dmn);CHKERRQ(ierr);
1211     } else {
1212       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1213     }
1214     next = next->next;
1215   }
1216   va_end(Argp);
1217   PetscFunctionReturn(0);
1218 }
1219 
1220 #undef __FUNCT__
1221 #define __FUNCT__ "DMRefine_Composite"
1222 PetscErrorCode PETSCDM_DLLEXPORT DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
1223 {
1224   PetscErrorCode         ierr;
1225   struct DMCompositeLink *next;
1226   DM_Composite           *com = (DM_Composite*)dmi->data;
1227   DM                     dm;
1228 
1229   PetscFunctionBegin;
1230   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1231   next = com->next;
1232   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
1233 
1234   /* loop over packed objects, handling one at at time */
1235   while (next) {
1236     if (next->type == DMCOMPOSITE_ARRAY) {
1237       ierr = DMCompositeAddArray(*fine,next->rank,next->n);CHKERRQ(ierr);
1238     } else if (next->type == DMCOMPOSITE_DM) {
1239       ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
1240       ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
1241       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1242     } else {
1243       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1244     }
1245     next = next->next;
1246   }
1247   PetscFunctionReturn(0);
1248 }
1249 
1250 #include "petscmat.h"
1251 
1252 struct MatPackLink {
1253   Mat                A;
1254   struct MatPackLink *next;
1255 };
1256 
1257 struct MatPack {
1258   DM                 right,left;
1259   struct MatPackLink *next;
1260 };
1261 
1262 #undef __FUNCT__
1263 #define __FUNCT__ "MatMultBoth_Shell_Pack"
1264 PetscErrorCode MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscBool  add)
1265 {
1266   struct MatPack         *mpack;
1267   struct DMCompositeLink *xnext,*ynext;
1268   struct MatPackLink     *anext;
1269   PetscScalar            *xarray,*yarray;
1270   PetscErrorCode         ierr;
1271   PetscInt               i;
1272   Vec                    xglobal,yglobal;
1273   PetscMPIInt            rank;
1274   DM_Composite           *comright;
1275   DM_Composite           *comleft;
1276 
1277   PetscFunctionBegin;
1278   ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
1279   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
1280   comright = (DM_Composite*)mpack->right->data;
1281   comleft = (DM_Composite*)mpack->left->data;
1282   xnext = comright->next;
1283   ynext = comleft->next;
1284   anext = mpack->next;
1285 
1286   while (xnext) {
1287     if (xnext->type == DMCOMPOSITE_ARRAY) {
1288       if (rank == xnext->rank) {
1289         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
1290         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
1291         if (add) {
1292           for (i=0; i<xnext->n; i++) {
1293             yarray[ynext->rstart+i] += xarray[xnext->rstart+i];
1294           }
1295         } else {
1296           ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
1297         }
1298         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1299         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1300       }
1301     } else if (xnext->type == DMCOMPOSITE_DM) {
1302       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
1303       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
1304       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
1305       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
1306       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
1307       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
1308       if (add) {
1309         ierr  = MatMultAdd(anext->A,xglobal,yglobal,yglobal);CHKERRQ(ierr);
1310       } else {
1311         ierr  = MatMult(anext->A,xglobal,yglobal);CHKERRQ(ierr);
1312       }
1313       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1314       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1315       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
1316       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
1317       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
1318       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
1319       anext = anext->next;
1320     } else {
1321       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1322     }
1323     xnext = xnext->next;
1324     ynext = ynext->next;
1325   }
1326   PetscFunctionReturn(0);
1327 }
1328 
1329 #undef __FUNCT__
1330 #define __FUNCT__ "MatMultAdd_Shell_Pack"
1331 PetscErrorCode MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z)
1332 {
1333   PetscErrorCode ierr;
1334   PetscFunctionBegin;
1335   if (z != y) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Handles y == z only");
1336   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
1337   PetscFunctionReturn(0);
1338 }
1339 
1340 #undef __FUNCT__
1341 #define __FUNCT__ "MatMult_Shell_Pack"
1342 PetscErrorCode MatMult_Shell_Pack(Mat A,Vec x,Vec y)
1343 {
1344   PetscErrorCode ierr;
1345   PetscFunctionBegin;
1346   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
1347   PetscFunctionReturn(0);
1348 }
1349 
1350 #undef __FUNCT__
1351 #define __FUNCT__ "MatMultTranspose_Shell_Pack"
1352 PetscErrorCode MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y)
1353 {
1354   struct MatPack         *mpack;
1355   struct DMCompositeLink *xnext,*ynext;
1356   struct MatPackLink     *anext;
1357   PetscScalar            *xarray,*yarray;
1358   PetscErrorCode         ierr;
1359   Vec                    xglobal,yglobal;
1360   PetscMPIInt            rank;
1361   DM_Composite           *comright;
1362   DM_Composite           *comleft;
1363 
1364   PetscFunctionBegin;
1365   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
1366   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
1367   comright = (DM_Composite*)mpack->right->data;
1368   comleft = (DM_Composite*)mpack->left->data;
1369   ynext = comright->next;
1370   xnext = comleft->next;
1371   anext = mpack->next;
1372 
1373   while (xnext) {
1374     if (xnext->type == DMCOMPOSITE_ARRAY) {
1375       if (rank == ynext->rank) {
1376         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
1377         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
1378         ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
1379         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1380         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1381       }
1382     } else if (xnext->type == DMCOMPOSITE_DM) {
1383       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
1384       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
1385       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
1386       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
1387       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
1388       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
1389       ierr  = MatMultTranspose(anext->A,xglobal,yglobal);CHKERRQ(ierr);
1390       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
1391       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
1392       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
1393       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
1394       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
1395       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
1396       anext = anext->next;
1397     } else {
1398       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1399     }
1400     xnext = xnext->next;
1401     ynext = ynext->next;
1402   }
1403   PetscFunctionReturn(0);
1404 }
1405 
1406 #undef __FUNCT__
1407 #define __FUNCT__ "MatDestroy_Shell_Pack"
1408 PetscErrorCode MatDestroy_Shell_Pack(Mat A)
1409 {
1410   struct MatPack     *mpack;
1411   struct MatPackLink *anext,*oldanext;
1412   PetscErrorCode     ierr;
1413 
1414   PetscFunctionBegin;
1415   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
1416   anext = mpack->next;
1417 
1418   while (anext) {
1419     ierr     = MatDestroy(anext->A);CHKERRQ(ierr);
1420     oldanext = anext;
1421     anext    = anext->next;
1422     ierr     = PetscFree(oldanext);CHKERRQ(ierr);
1423   }
1424   ierr = PetscFree(mpack);CHKERRQ(ierr);
1425   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
1426   PetscFunctionReturn(0);
1427 }
1428 
1429 #undef __FUNCT__
1430 #define __FUNCT__ "DMGetInterpolation_Composite"
1431 PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
1432 {
1433   PetscErrorCode         ierr;
1434   PetscInt               m,n,M,N;
1435   struct DMCompositeLink *nextc;
1436   struct DMCompositeLink *nextf;
1437   struct MatPackLink     *nextmat,*pnextmat = 0;
1438   struct MatPack         *mpack;
1439   Vec                    gcoarse,gfine;
1440   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
1441   DM_Composite           *comfine = (DM_Composite*)fine->data;
1442 
1443   PetscFunctionBegin;
1444   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
1445   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
1446   nextc = comcoarse->next;
1447   nextf = comfine->next;
1448   /* use global vectors only for determining matrix layout */
1449   ierr = DMCreateGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
1450   ierr = DMCreateGlobalVector(fine,&gfine);CHKERRQ(ierr);
1451   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
1452   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
1453   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
1454   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
1455   ierr = VecDestroy(gcoarse);CHKERRQ(ierr);
1456   ierr = VecDestroy(gfine);CHKERRQ(ierr);
1457 
1458   ierr         = PetscNew(struct MatPack,&mpack);CHKERRQ(ierr);
1459   mpack->right = coarse;
1460   mpack->left  = fine;
1461   ierr  = MatCreate(((PetscObject)fine)->comm,A);CHKERRQ(ierr);
1462   ierr  = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1463   ierr  = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
1464   ierr  = MatShellSetContext(*A,mpack);CHKERRQ(ierr);
1465   ierr  = MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);CHKERRQ(ierr);
1466   ierr  = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);CHKERRQ(ierr);
1467   ierr  = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);CHKERRQ(ierr);
1468   ierr  = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);CHKERRQ(ierr);
1469 
1470   /* loop over packed objects, handling one at at time */
1471   while (nextc) {
1472     if (nextc->type != nextf->type) SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Two DM have different layout");
1473 
1474     if (nextc->type == DMCOMPOSITE_ARRAY) {
1475       ;
1476     } else if (nextc->type == DMCOMPOSITE_DM) {
1477       ierr          = PetscNew(struct MatPackLink,&nextmat);CHKERRQ(ierr);
1478       nextmat->next = 0;
1479       if (pnextmat) {
1480         pnextmat->next = nextmat;
1481         pnextmat       = nextmat;
1482       } else {
1483         pnextmat    = nextmat;
1484         mpack->next = nextmat;
1485       }
1486       ierr = DMGetInterpolation(nextc->dm,nextf->dm,&nextmat->A,PETSC_NULL);CHKERRQ(ierr);
1487     } else {
1488       SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1489     }
1490     nextc = nextc->next;
1491     nextf = nextf->next;
1492   }
1493   PetscFunctionReturn(0);
1494 }
1495 
1496 #undef __FUNCT__
1497 #define __FUNCT__ "DMGetMatrix_Composite"
1498 PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix_Composite(DM dm, const MatType mtype,Mat *J)
1499 {
1500   PetscErrorCode         ierr;
1501   DM_Composite           *com = (DM_Composite*)dm->data;
1502   struct DMCompositeLink *next = com->next;
1503   PetscInt               m,*dnz,*onz,i,j,mA;
1504   Mat                    Atmp;
1505   PetscMPIInt            rank;
1506   PetscScalar            zero = 0.0;
1507   PetscBool              dense = PETSC_FALSE;
1508 
1509   PetscFunctionBegin;
1510   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1511 
1512   /* use global vector to determine layout needed for matrix */
1513   m = com->n;
1514   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1515   ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr);
1516   ierr = MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1517   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
1518 
1519   /*
1520      Extremely inefficient but will compute entire Jacobian for testing
1521   */
1522   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
1523   if (dense) {
1524     PetscInt    rstart,rend,*indices;
1525     PetscScalar *values;
1526 
1527     mA    = com->N;
1528     ierr = MatMPIAIJSetPreallocation(*J,mA,PETSC_NULL,mA-m,PETSC_NULL);CHKERRQ(ierr);
1529     ierr = MatSeqAIJSetPreallocation(*J,mA,PETSC_NULL);CHKERRQ(ierr);
1530 
1531     ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr);
1532     ierr = PetscMalloc2(mA,PetscScalar,&values,mA,PetscInt,&indices);CHKERRQ(ierr);
1533     ierr = PetscMemzero(values,mA*sizeof(PetscScalar));CHKERRQ(ierr);
1534     for (i=0; i<mA; i++) indices[i] = i;
1535     for (i=rstart; i<rend; i++) {
1536       ierr = MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES);CHKERRQ(ierr);
1537     }
1538     ierr = PetscFree2(values,indices);CHKERRQ(ierr);
1539     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1540     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1541     PetscFunctionReturn(0);
1542   }
1543 
1544   ierr = MatPreallocateInitialize(((PetscObject)dm)->comm,m,m,dnz,onz);CHKERRQ(ierr);
1545   /* loop over packed objects, handling one at at time */
1546   next = com->next;
1547   while (next) {
1548     if (next->type == DMCOMPOSITE_ARRAY) {
1549       if (rank == next->rank) {  /* zero the "little" block */
1550         for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) {
1551           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
1552             ierr = MatPreallocateSet(j,1,&i,dnz,onz);CHKERRQ(ierr);
1553           }
1554         }
1555       }
1556     } else if (next->type == DMCOMPOSITE_DM) {
1557       PetscInt       nc,rstart,*ccols,maxnc;
1558       const PetscInt *cols,*rstarts;
1559       PetscMPIInt    proc;
1560 
1561       ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr);
1562       ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr);
1563       ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
1564       ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr);
1565 
1566       maxnc = 0;
1567       for (i=0; i<mA; i++) {
1568         ierr  = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1569         ierr  = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1570         maxnc = PetscMax(nc,maxnc);
1571       }
1572       ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
1573       for (i=0; i<mA; i++) {
1574         ierr = MatGetRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr);
1575         /* remap the columns taking into how much they are shifted on each process */
1576         for (j=0; j<nc; j++) {
1577           proc = 0;
1578           while (cols[j] >= rstarts[proc+1]) proc++;
1579           ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
1580         }
1581         ierr = MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz);CHKERRQ(ierr);
1582         ierr = MatRestoreRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr);
1583       }
1584       ierr = PetscFree(ccols);CHKERRQ(ierr);
1585       ierr = MatDestroy(Atmp);CHKERRQ(ierr);
1586     } else {
1587       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1588     }
1589     next = next->next;
1590   }
1591   if (com->FormCoupleLocations) {
1592     ierr = (*com->FormCoupleLocations)(dm,PETSC_NULL,dnz,onz,__rstart,__nrows,__start,__end);CHKERRQ(ierr);
1593   }
1594   ierr = MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);CHKERRQ(ierr);
1595   ierr = MatSeqAIJSetPreallocation(*J,0,dnz);CHKERRQ(ierr);
1596   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1597 
1598   next = com->next;
1599   while (next) {
1600     if (next->type == DMCOMPOSITE_ARRAY) {
1601       if (rank == next->rank) {
1602         for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) {
1603           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
1604             ierr = MatSetValues(*J,1,&j,1,&i,&zero,INSERT_VALUES);CHKERRQ(ierr);
1605           }
1606         }
1607       }
1608     } else if (next->type == DMCOMPOSITE_DM) {
1609       PetscInt          nc,rstart,row,maxnc,*ccols;
1610       const PetscInt    *cols,*rstarts;
1611       const PetscScalar *values;
1612       PetscMPIInt       proc;
1613 
1614       ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr);
1615       ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr);
1616       ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
1617       ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr);
1618       maxnc = 0;
1619       for (i=0; i<mA; i++) {
1620         ierr  = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1621         ierr  = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
1622         maxnc = PetscMax(nc,maxnc);
1623       }
1624       ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
1625       for (i=0; i<mA; i++) {
1626         ierr = MatGetRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr);
1627         for (j=0; j<nc; j++) {
1628           proc = 0;
1629           while (cols[j] >= rstarts[proc+1]) proc++;
1630           ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
1631         }
1632         row  = com->rstart+next->rstart+i;
1633         ierr = MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES);CHKERRQ(ierr);
1634         ierr = MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr);
1635       }
1636       ierr = PetscFree(ccols);CHKERRQ(ierr);
1637       ierr = MatDestroy(Atmp);CHKERRQ(ierr);
1638     } else {
1639       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1640     }
1641     next = next->next;
1642   }
1643   if (com->FormCoupleLocations) {
1644     PetscInt __rstart;
1645     ierr = MatGetOwnershipRange(*J,&__rstart,PETSC_NULL);CHKERRQ(ierr);
1646     ierr = (*com->FormCoupleLocations)(dm,*J,PETSC_NULL,PETSC_NULL,__rstart,0,0,0);CHKERRQ(ierr);
1647   }
1648   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1649   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1650   PetscFunctionReturn(0);
1651 }
1652 
1653 #undef __FUNCT__
1654 #define __FUNCT__ "DMGetColoring_Composite"
1655 PetscErrorCode PETSCDM_DLLEXPORT DMGetColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
1656 {
1657   PetscErrorCode         ierr;
1658   PetscInt               n,i,cnt;
1659   ISColoringValue        *colors;
1660   PetscBool              dense = PETSC_FALSE;
1661   ISColoringValue        maxcol = 0;
1662   DM_Composite           *com = (DM_Composite*)dm->data;
1663 
1664   PetscFunctionBegin;
1665   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1666   if (ctype == IS_COLORING_GHOSTED) {
1667     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" );
1668   } else if (ctype == IS_COLORING_GLOBAL) {
1669     n = com->n;
1670   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1671   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
1672 
1673   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
1674   if (dense) {
1675     for (i=0; i<n; i++) {
1676       colors[i] = (ISColoringValue)(com->rstart + i);
1677     }
1678     maxcol = com->N;
1679   } else {
1680     struct DMCompositeLink *next = com->next;
1681     PetscMPIInt            rank;
1682 
1683     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1684     cnt  = 0;
1685     while (next) {
1686       if (next->type == DMCOMPOSITE_ARRAY) {
1687         if (rank == next->rank) {  /* each column gets is own color */
1688           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
1689             colors[cnt++] = maxcol++;
1690           }
1691         }
1692         ierr = MPI_Bcast(&maxcol,1,MPIU_COLORING_VALUE,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
1693       } else if (next->type == DMCOMPOSITE_DM) {
1694         ISColoring     lcoloring;
1695 
1696         ierr = DMGetColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
1697         for (i=0; i<lcoloring->N; i++) {
1698           colors[cnt++] = maxcol + lcoloring->colors[i];
1699         }
1700         maxcol += lcoloring->n;
1701         ierr = ISColoringDestroy(lcoloring);CHKERRQ(ierr);
1702       } else {
1703         SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1704       }
1705       next = next->next;
1706     }
1707   }
1708   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
1709   PetscFunctionReturn(0);
1710 }
1711 
1712 #undef __FUNCT__
1713 #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
1714 PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1715 {
1716   PetscErrorCode         ierr;
1717   struct DMCompositeLink *next;
1718   PetscInt               cnt = 3;
1719   PetscMPIInt            rank;
1720   PetscScalar            *garray,*larray;
1721   DM_Composite           *com = (DM_Composite*)dm->data;
1722 
1723   PetscFunctionBegin;
1724   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1725   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
1726   next = com->next;
1727   if (!com->setup) {
1728     ierr = DMSetUp(dm);CHKERRQ(ierr);
1729   }
1730   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1731   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
1732   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
1733 
1734   /* loop over packed objects, handling one at at time */
1735   while (next) {
1736     if (next->type == DMCOMPOSITE_ARRAY) {
1737       if (rank == next->rank) {
1738         ierr    = PetscMemcpy(larray,garray,next->n*sizeof(PetscScalar));CHKERRQ(ierr);
1739         garray += next->n;
1740       }
1741       /* does not handle ADD_VALUES */
1742       ierr = MPI_Bcast(larray,next->n,MPIU_SCALAR,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
1743       larray += next->n;
1744     } else if (next->type == DMCOMPOSITE_DM) {
1745       Vec      local,global;
1746       PetscInt N;
1747 
1748       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
1749       ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
1750       ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
1751       ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
1752       ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
1753       ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
1754       ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
1755       ierr = VecResetArray(global);CHKERRQ(ierr);
1756       ierr = VecResetArray(local);CHKERRQ(ierr);
1757       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
1758       ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
1759       larray += next->n;
1760     } else {
1761       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
1762     }
1763     cnt++;
1764     next = next->next;
1765   }
1766 
1767   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
1768   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
1769   PetscFunctionReturn(0);
1770 }
1771 
1772 #undef __FUNCT__
1773 #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
1774 PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1775 {
1776   PetscFunctionBegin;
1777   PetscFunctionReturn(0);
1778 }
1779 
1780 EXTERN_C_BEGIN
1781 #undef __FUNCT__
1782 #define __FUNCT__ "DMCreate_Composite"
1783 PetscErrorCode PETSCDM_DLLEXPORT DMCreate_Composite(DM p)
1784 {
1785   PetscErrorCode ierr;
1786   DM_Composite   *com;
1787 
1788   PetscFunctionBegin;
1789   ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1790   p->data = com;
1791   ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1792   com->n            = 0;
1793   com->next         = PETSC_NULL;
1794   com->nredundant   = 0;
1795   com->nDM          = 0;
1796 
1797   p->ops->createglobalvector = DMCreateGlobalVector_Composite;
1798   p->ops->createlocalvector  = DMCreateLocalVector_Composite;
1799   p->ops->refine             = DMRefine_Composite;
1800   p->ops->getinterpolation   = DMGetInterpolation_Composite;
1801   p->ops->getmatrix          = DMGetMatrix_Composite;
1802   p->ops->getcoloring        = DMGetColoring_Composite;
1803   p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite;
1804   p->ops->globaltolocalend   = DMGlobalToLocalEnd_Composite;
1805   p->ops->destroy            = DMDestroy_Composite;
1806   p->ops->view               = DMView_Composite;
1807   p->ops->setup              = DMSetUp_Composite;
1808   PetscFunctionReturn(0);
1809 }
1810 EXTERN_C_END
1811 
1812 #undef __FUNCT__
1813 #define __FUNCT__ "DMCompositeCreate"
1814 /*@C
1815     DMCompositeCreate - Creates a vector packer, used to generate "composite"
1816       vectors made up of several subvectors.
1817 
1818     Collective on MPI_Comm
1819 
1820     Input Parameter:
1821 .   comm - the processors that will share the global vector
1822 
1823     Output Parameters:
1824 .   packer - the packer object
1825 
1826     Level: advanced
1827 
1828 .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCompositeScatter(),
1829          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess()
1830          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
1831 
1832 @*/
1833 PetscErrorCode PETSCDM_DLLEXPORT DMCompositeCreate(MPI_Comm comm,DM *packer)
1834 {
1835   PetscErrorCode ierr;
1836 
1837   PetscFunctionBegin;
1838   PetscValidPointer(packer,2);
1839   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1840   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
1841   PetscFunctionReturn(0);
1842 }
1843