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