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