xref: /petsc/src/dm/impls/composite/pack.c (revision f08646a875d6c8bb0dc8fb3d7c6d74e67952873a)
1 
2 #include "packimpl.h" /*I   "petscdmcomposite.h"   I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMCompositeSetCoupling"
6 /*@C
7     DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
8       seperate components (DMs) in a DMto build the correct matrix nonzero structure.
9 
10 
11     Logically Collective on MPI_Comm
12 
13     Input Parameter:
14 +   dm - the composite object
15 -   formcouplelocations - routine to set the nonzero locations in the matrix
16 
17     Level: advanced
18 
19     Notes: See DMSetApplicationContext() and DMGetApplicationContext() for how to get user information into
20         this routine
21 
22 @*/
23 PetscErrorCode  DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt))
24 {
25   DM_Composite *com = (DM_Composite*)dm->data;
26 
27   PetscFunctionBegin;
28   com->FormCoupleLocations = FormCoupleLocations;
29   PetscFunctionReturn(0);
30 }
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "DMDestroy_Composite"
34 PetscErrorCode  DMDestroy_Composite(DM dm)
35 {
36   PetscErrorCode         ierr;
37   struct DMCompositeLink *next, *prev;
38   DM_Composite           *com = (DM_Composite *)dm->data;
39 
40   PetscFunctionBegin;
41   next = com->next;
42   while (next) {
43     prev = next;
44     next = next->next;
45     ierr = DMDestroy(&prev->dm);CHKERRQ(ierr);
46     ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
47     ierr = PetscFree(prev);CHKERRQ(ierr);
48   }
49   PetscFunctionReturn(0);
50 }
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "DMView_Composite"
54 PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
55 {
56   PetscErrorCode ierr;
57   PetscBool      iascii;
58   DM_Composite   *com = (DM_Composite *)dm->data;
59 
60   PetscFunctionBegin;
61   ierr = PetscTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
62   if (iascii) {
63     struct DMCompositeLink *lnk = com->next;
64     PetscInt               i;
65 
66     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix?((PetscObject)dm)->prefix:"no prefix");CHKERRQ(ierr);
67     ierr = PetscViewerASCIIPrintf(v,"  contains %D DMs\n",com->nDM);CHKERRQ(ierr);
68     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
69     for (i=0; lnk; lnk=lnk->next,i++) {
70       ierr = PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
71       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
72       ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
73       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
74     }
75     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
76   }
77   PetscFunctionReturn(0);
78 }
79 
80 /* --------------------------------------------------------------------------------------*/
81 #undef __FUNCT__
82 #define __FUNCT__ "DMSetUp_Composite"
83 PetscErrorCode  DMSetUp_Composite(DM dm)
84 {
85   PetscErrorCode         ierr;
86   PetscInt               nprev = 0;
87   PetscMPIInt            rank,size;
88   DM_Composite           *com = (DM_Composite*)dm->data;
89   struct DMCompositeLink *next = com->next;
90   PetscLayout            map;
91 
92   PetscFunctionBegin;
93   if (com->setup) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
94   ierr = PetscLayoutCreate(((PetscObject)dm)->comm,&map);CHKERRQ(ierr);
95   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
96   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
97   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
98   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
99   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
100   ierr = PetscLayoutGetRange(map,&com->rstart,PETSC_NULL);CHKERRQ(ierr);
101   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
102 
103   /* now set the rstart for each linked vector */
104   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
105   ierr = MPI_Comm_size(((PetscObject)dm)->comm,&size);CHKERRQ(ierr);
106   while (next) {
107     next->rstart = nprev;
108     nprev += next->n;
109     next->grstart = com->rstart + next->rstart;
110     ierr = PetscMalloc(size*sizeof(PetscInt),&next->grstarts);CHKERRQ(ierr);
111     ierr = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
112     next = next->next;
113   }
114   com->setup = PETSC_TRUE;
115   PetscFunctionReturn(0);
116 }
117 
118 /* ----------------------------------------------------------------------------------*/
119 
120 #include <stdarg.h>
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "DMCompositeGetNumberDM"
124 /*@C
125     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
126        representation.
127 
128     Not Collective
129 
130     Input Parameter:
131 .    dm - the packer object
132 
133     Output Parameter:
134 .     nDM - the number of DMs
135 
136     Level: beginner
137 
138 @*/
139 PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
140 {
141   DM_Composite *com = (DM_Composite*)dm->data;
142   PetscFunctionBegin;
143   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
144   *nDM = com->nDM;
145   PetscFunctionReturn(0);
146 }
147 
148 
149 #undef __FUNCT__
150 #define __FUNCT__ "DMCompositeGetAccess"
151 /*@C
152     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
153        representation.
154 
155     Collective on DMComposite
156 
157     Input Parameters:
158 +    dm - the packer object
159 -    gvec - the global vector
160 
161     Output Parameters:
162 .    Vec* ... - the packed parallel vectors, PETSC_NULL for those that are not needed
163 
164     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
165 
166     Level: advanced
167 
168 @*/
169 PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
170 {
171   va_list                Argp;
172   PetscErrorCode         ierr;
173   struct DMCompositeLink *next;
174   DM_Composite           *com = (DM_Composite*)dm->data;
175 
176   PetscFunctionBegin;
177   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
178   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
179   next = com->next;
180   if (!com->setup) {
181     ierr = DMSetUp(dm);CHKERRQ(ierr);
182   }
183 
184   /* loop over packed objects, handling one at at time */
185   va_start(Argp,gvec);
186   while (next) {
187     Vec *vec;
188     vec = va_arg(Argp, Vec*);
189     if (vec) {
190       PetscScalar *array;
191       ierr = DMGetGlobalVector(next->dm,vec);CHKERRQ(ierr);
192       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
193       ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
194       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
195     }
196     next = next->next;
197   }
198   va_end(Argp);
199   PetscFunctionReturn(0);
200 }
201 
202 #undef __FUNCT__
203 #define __FUNCT__ "DMCompositeRestoreAccess"
204 /*@C
205     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
206        representation.
207 
208     Collective on DMComposite
209 
210     Input Parameters:
211 +    dm - the packer object
212 .    gvec - the global vector
213 -    Vec* ... - the individual parallel vectors, PETSC_NULL for those that are not needed
214 
215     Level: advanced
216 
217 .seealso  DMCompositeAddDM(), DMCreateGlobalVector(),
218          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
219          DMCompositeRestoreAccess(), DMCompositeGetAccess()
220 
221 @*/
222 PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
223 {
224   va_list                Argp;
225   PetscErrorCode         ierr;
226   struct DMCompositeLink *next;
227   DM_Composite           *com = (DM_Composite*)dm->data;
228 
229   PetscFunctionBegin;
230   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
231   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
232   next = com->next;
233   if (!com->setup) {
234     ierr = DMSetUp(dm);CHKERRQ(ierr);
235   }
236 
237   /* loop over packed objects, handling one at at time */
238   va_start(Argp,gvec);
239   while (next) {
240     Vec *vec;
241     vec = va_arg(Argp, Vec*);
242     if (vec) {
243       ierr = VecResetArray(*vec);CHKERRQ(ierr);
244       ierr = DMRestoreGlobalVector(next->dm,vec);CHKERRQ(ierr);
245     }
246     next = next->next;
247   }
248   va_end(Argp);
249   PetscFunctionReturn(0);
250 }
251 
252 #undef __FUNCT__
253 #define __FUNCT__ "DMCompositeScatter"
254 /*@C
255     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
256 
257     Collective on DMComposite
258 
259     Input Parameters:
260 +    dm - the packer object
261 .    gvec - the global vector
262 -    Vec ... - the individual sequential vectors, PETSC_NULL for those that are not needed
263 
264     Level: advanced
265 
266 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
267          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
268          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
269 
270 @*/
271 PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
272 {
273   va_list                Argp;
274   PetscErrorCode         ierr;
275   struct DMCompositeLink *next;
276   PetscInt               cnt;
277   DM_Composite           *com = (DM_Composite*)dm->data;
278 
279   PetscFunctionBegin;
280   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
281   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
282   if (!com->setup) {
283     ierr = DMSetUp(dm);CHKERRQ(ierr);
284   }
285 
286   /* loop over packed objects, handling one at at time */
287   va_start(Argp,gvec);
288   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
289     Vec local;
290     local = va_arg(Argp, Vec);
291     if (local) {
292       Vec global;
293       PetscScalar *array;
294       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
295       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
296       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
297       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
298       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
299       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
300       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
301       ierr = VecResetArray(global);CHKERRQ(ierr);
302       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
303     }
304   }
305   va_end(Argp);
306   PetscFunctionReturn(0);
307 }
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "DMCompositeGather"
311 /*@C
312     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
313 
314     Collective on DMComposite
315 
316     Input Parameter:
317 +    dm - the packer object
318 .    gvec - the global vector
319 -    Vec ... - the individual sequential vectors, PETSC_NULL for any that are not needed
320 
321     Level: advanced
322 
323 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
324          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
325          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
326 
327 @*/
328 PetscErrorCode  DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
329 {
330   va_list                Argp;
331   PetscErrorCode         ierr;
332   struct DMCompositeLink *next;
333   DM_Composite           *com = (DM_Composite*)dm->data;
334   PetscInt               cnt;
335 
336   PetscFunctionBegin;
337   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
338   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
339   if (!com->setup) {
340     ierr = DMSetUp(dm);CHKERRQ(ierr);
341   }
342 
343   /* loop over packed objects, handling one at at time */
344   va_start(Argp,imode);
345   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
346     Vec local;
347     local = va_arg(Argp, Vec);
348     if (local) {
349       PetscScalar    *array;
350       Vec            global;
351       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
352       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
353       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
354       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
355       ierr = DMLocalToGlobalBegin(next->dm,local,imode,global);CHKERRQ(ierr);
356       ierr = DMLocalToGlobalEnd(next->dm,local,imode,global);CHKERRQ(ierr);
357       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
358       ierr = VecResetArray(global);CHKERRQ(ierr);
359       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
360     }
361   }
362   va_end(Argp);
363   PetscFunctionReturn(0);
364 }
365 
366 #undef __FUNCT__
367 #define __FUNCT__ "DMCompositeAddDM"
368 /*@C
369     DMCompositeAddDM - adds a DM  vector to a DMComposite
370 
371     Collective on DMComposite
372 
373     Input Parameter:
374 +    dm - the packer object
375 -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
376 
377     Level: advanced
378 
379 .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
380          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
381          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
382 
383 @*/
384 PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
385 {
386   PetscErrorCode         ierr;
387   PetscInt               n,nlocal;
388   struct DMCompositeLink *mine,*next;
389   Vec                    global,local;
390   DM_Composite           *com = (DM_Composite*)dmc->data;
391 
392   PetscFunctionBegin;
393   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
394   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
395   next = com->next;
396   if (com->setup) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
397 
398   /* create new link */
399   ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
400   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
401   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
402   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
403   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
404   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
405   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
406   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
407   mine->n      = n;
408   mine->nlocal = nlocal;
409   mine->dm     = dm;
410   mine->next   = PETSC_NULL;
411   com->n       += n;
412 
413   /* add to end of list */
414   if (!next) {
415     com->next = mine;
416   } else {
417     while (next->next) next = next->next;
418     next->next = mine;
419   }
420   com->nDM++;
421   com->nmine++;
422   PetscFunctionReturn(0);
423 }
424 
425 extern PetscErrorCode  VecView_MPI(Vec,PetscViewer);
426 EXTERN_C_BEGIN
427 #undef __FUNCT__
428 #define __FUNCT__ "VecView_DMComposite"
429 PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
430 {
431   DM                     dm;
432   PetscErrorCode         ierr;
433   struct DMCompositeLink *next;
434   PetscBool              isdraw;
435   DM_Composite           *com;
436 
437   PetscFunctionBegin;
438   ierr = PetscObjectQuery((PetscObject)gvec,"DM",(PetscObject*)&dm);CHKERRQ(ierr);
439   if (!dm) SETERRQ(((PetscObject)gvec)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
440   com = (DM_Composite*)dm->data;
441   next = com->next;
442 
443   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
444   if (!isdraw) {
445     /* do I really want to call this? */
446     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
447   } else {
448     PetscInt cnt = 0;
449 
450     /* loop over packed objects, handling one at at time */
451     while (next) {
452       Vec         vec;
453       PetscScalar *array;
454       PetscInt    bs;
455 
456       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
457       ierr = DMGetGlobalVector(next->dm,&vec);CHKERRQ(ierr);
458       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
459       ierr = VecPlaceArray(vec,array+next->rstart);CHKERRQ(ierr);
460       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
461       ierr = VecView(vec,viewer);CHKERRQ(ierr);
462       ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
463       ierr = VecResetArray(vec);CHKERRQ(ierr);
464       ierr = DMRestoreGlobalVector(next->dm,&vec);CHKERRQ(ierr);
465       ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
466       cnt += bs;
467       next = next->next;
468     }
469     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
470   }
471   PetscFunctionReturn(0);
472 }
473 EXTERN_C_END
474 
475 
476 #undef __FUNCT__
477 #define __FUNCT__ "DMCreateGlobalVector_Composite"
478 PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
479 {
480   PetscErrorCode         ierr;
481   DM_Composite           *com = (DM_Composite*)dm->data;
482 
483   PetscFunctionBegin;
484   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
485   if (!com->setup) {
486     ierr = DMSetUp(dm);CHKERRQ(ierr);
487   }
488   ierr = VecCreateMPI(((PetscObject)dm)->comm,com->n,com->N,gvec);CHKERRQ(ierr);
489   ierr = PetscObjectCompose((PetscObject)*gvec,"DM",(PetscObject)dm);CHKERRQ(ierr);
490   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void(*)(void))VecView_DMComposite);CHKERRQ(ierr);
491   PetscFunctionReturn(0);
492 }
493 
494 #undef __FUNCT__
495 #define __FUNCT__ "DMCreateLocalVector_Composite"
496 PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
497 {
498   PetscErrorCode         ierr;
499   DM_Composite           *com = (DM_Composite*)dm->data;
500 
501   PetscFunctionBegin;
502   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
503   if (!com->setup) {
504     ierr = DMSetUp(dm);CHKERRQ(ierr);
505   }
506   ierr = VecCreateSeq(((PetscObject)dm)->comm,com->nghost,lvec);CHKERRQ(ierr);
507   ierr = PetscObjectCompose((PetscObject)*lvec,"DM",(PetscObject)dm);CHKERRQ(ierr);
508   PetscFunctionReturn(0);
509 }
510 
511 #undef __FUNCT__
512 #define __FUNCT__ "DMCompositeGetISLocalToGlobalMappings"
513 /*@C
514     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
515 
516     Collective on DM
517 
518     Input Parameter:
519 .    dm - the packer object
520 
521     Output Parameters:
522 .    ltogs - the individual mappings for each packed vector. Note that this includes
523            all the ghost points that individual ghosted DMDA's may have.
524 
525     Level: advanced
526 
527     Notes:
528        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
529 
530 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
531          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
532          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
533 
534 @*/
535 PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
536 {
537   PetscErrorCode         ierr;
538   PetscInt               i,*idx,n,cnt;
539   struct DMCompositeLink *next;
540   PetscMPIInt            rank;
541   DM_Composite           *com = (DM_Composite*)dm->data;
542 
543   PetscFunctionBegin;
544   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
545   ierr = PetscMalloc((com->nDM)*sizeof(ISLocalToGlobalMapping),ltogs);CHKERRQ(ierr);
546   next = com->next;
547   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
548 
549   /* loop over packed objects, handling one at at time */
550   cnt = 0;
551   while (next) {
552     ISLocalToGlobalMapping ltog;
553     PetscMPIInt            size;
554     const PetscInt         *suboff,*indices;
555     Vec                    global;
556 
557     /* Get sub-DM global indices for each local dof */
558     ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
559     ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
560     ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
561     ierr = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
562 
563     /* Get the offsets for the sub-DM global vector */
564     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
565     ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
566     ierr = MPI_Comm_size(((PetscObject)global)->comm,&size);CHKERRQ(ierr);
567 
568     /* Shift the sub-DM definition of the global space to the composite global space */
569     for (i=0; i<n; i++) {
570       PetscInt subi = indices[i],lo = 0,hi = size,t;
571       /* Binary search to find which rank owns subi */
572       while (hi-lo > 1) {
573         t = lo + (hi-lo)/2;
574         if (suboff[t] > subi) hi = t;
575         else                  lo = t;
576       }
577       idx[i] = subi - suboff[lo] + next->grstarts[lo];
578     }
579     ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
580     ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
581     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
582     next = next->next;
583     cnt++;
584   }
585   PetscFunctionReturn(0);
586 }
587 
588 #undef __FUNCT__
589 #define __FUNCT__ "DMCompositeGetLocalISs"
590 /*@C
591    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
592 
593    Not Collective
594 
595    Input Arguments:
596 . dm - composite DM
597 
598    Output Arguments:
599 . is - array of serial index sets for each each component of the DMComposite
600 
601    Level: intermediate
602 
603    Notes:
604    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
605    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
606    scatter to a composite local vector.  The user should not typically need to know which is being done.
607 
608    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
609 
610    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
611 
612    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
613 
614 .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
615 @*/
616 PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
617 {
618   PetscErrorCode         ierr;
619   DM_Composite           *com = (DM_Composite*)dm->data;
620   struct DMCompositeLink *link;
621   PetscInt cnt,start;
622 
623   PetscFunctionBegin;
624   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
625   PetscValidPointer(is,2);
626   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
627   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
628     PetscInt bs;
629     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
630     ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
631     ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
632   }
633   PetscFunctionReturn(0);
634 }
635 
636 #undef __FUNCT__
637 #define __FUNCT__ "DMCompositeGetGlobalISs"
638 /*@C
639     DMCompositeGetGlobalISs - Gets the index sets for each composed object
640 
641     Collective on DMComposite
642 
643     Input Parameter:
644 .    dm - the packer object
645 
646     Output Parameters:
647 .    is - the array of index sets
648 
649     Level: advanced
650 
651     Notes:
652        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
653 
654        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
655 
656        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
657        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
658        indices.
659 
660 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
661          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
662          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
663 
664 @*/
665 
666 PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
667 {
668   PetscErrorCode         ierr;
669   PetscInt               cnt = 0,*idx,i;
670   struct DMCompositeLink *next;
671   PetscMPIInt            rank;
672   DM_Composite           *com = (DM_Composite*)dm->data;
673 
674   PetscFunctionBegin;
675   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
676   ierr = PetscMalloc((com->nDM)*sizeof(IS),is);CHKERRQ(ierr);
677   next = com->next;
678   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
679 
680   /* loop over packed objects, handling one at at time */
681   while (next) {
682     ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
683     for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
684     ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
685     cnt++;
686     next = next->next;
687   }
688   PetscFunctionReturn(0);
689 }
690 
691 /* -------------------------------------------------------------------------------------*/
692 #undef __FUNCT__
693 #define __FUNCT__ "DMCompositeGetLocalVectors"
694 /*@C
695     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
696        Use DMCompositeRestoreLocalVectors() to return them.
697 
698     Not Collective
699 
700     Input Parameter:
701 .    dm - the packer object
702 
703     Output Parameter:
704 .   Vec ... - the individual sequential Vecs
705 
706     Level: advanced
707 
708 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
709          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
710          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
711 
712 @*/
713 PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
714 {
715   va_list                Argp;
716   PetscErrorCode         ierr;
717   struct DMCompositeLink *next;
718   DM_Composite           *com = (DM_Composite*)dm->data;
719 
720   PetscFunctionBegin;
721   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
722   next = com->next;
723   /* loop over packed objects, handling one at at time */
724   va_start(Argp,dm);
725   while (next) {
726     Vec *vec;
727     vec = va_arg(Argp, Vec*);
728     ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);
729     next = next->next;
730   }
731   va_end(Argp);
732   PetscFunctionReturn(0);
733 }
734 
735 #undef __FUNCT__
736 #define __FUNCT__ "DMCompositeRestoreLocalVectors"
737 /*@C
738     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
739 
740     Not Collective
741 
742     Input Parameter:
743 .    dm - the packer object
744 
745     Output Parameter:
746 .   Vec ... - the individual sequential Vecs
747 
748     Level: advanced
749 
750 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
751          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
752          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
753 
754 @*/
755 PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
756 {
757   va_list                Argp;
758   PetscErrorCode         ierr;
759   struct DMCompositeLink *next;
760   DM_Composite           *com = (DM_Composite*)dm->data;
761 
762   PetscFunctionBegin;
763   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
764   next = com->next;
765   /* loop over packed objects, handling one at at time */
766   va_start(Argp,dm);
767   while (next) {
768     Vec *vec;
769     vec = va_arg(Argp, Vec*);
770     ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);
771     next = next->next;
772   }
773   va_end(Argp);
774   PetscFunctionReturn(0);
775 }
776 
777 /* -------------------------------------------------------------------------------------*/
778 #undef __FUNCT__
779 #define __FUNCT__ "DMCompositeGetEntries"
780 /*@C
781     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
782 
783     Not Collective
784 
785     Input Parameter:
786 .    dm - the packer object
787 
788     Output Parameter:
789 .   DM ... - the individual entries (DMs)
790 
791     Level: advanced
792 
793 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
794          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
795          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
796          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
797 
798 @*/
799 PetscErrorCode  DMCompositeGetEntries(DM dm,...)
800 {
801   va_list                Argp;
802   struct DMCompositeLink *next;
803   DM_Composite           *com = (DM_Composite*)dm->data;
804 
805   PetscFunctionBegin;
806   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
807   next = com->next;
808   /* loop over packed objects, handling one at at time */
809   va_start(Argp,dm);
810   while (next) {
811     DM *dmn;
812     dmn = va_arg(Argp, DM*);
813     if (dmn) *dmn = next->dm;
814     next = next->next;
815   }
816   va_end(Argp);
817   PetscFunctionReturn(0);
818 }
819 
820 #undef __FUNCT__
821 #define __FUNCT__ "DMCompositeGetEntriesArray"
822 /*@
823     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
824 
825     Not Collective
826 
827     Input Parameter:
828 +    dm - the packer object
829 -    dms - array of sufficient length (see DMCompositeGetNumberDM()), holds the DMs on output
830 
831     Level: advanced
832 
833 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
834          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
835          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
836          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
837 
838 @*/
839 PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
840 {
841   struct DMCompositeLink *next;
842   DM_Composite           *com = (DM_Composite*)dm->data;
843   PetscInt               i;
844 
845   PetscFunctionBegin;
846   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
847   /* loop over packed objects, handling one at at time */
848   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
849   PetscFunctionReturn(0);
850 }
851 
852 #undef __FUNCT__
853 #define __FUNCT__ "DMRefine_Composite"
854 PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
855 {
856   PetscErrorCode         ierr;
857   struct DMCompositeLink *next;
858   DM_Composite           *com = (DM_Composite*)dmi->data;
859   DM                     dm;
860 
861   PetscFunctionBegin;
862   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
863   next = com->next;
864   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
865 
866   /* loop over packed objects, handling one at at time */
867   while (next) {
868     ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
869     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
870     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
871     next = next->next;
872   }
873   PetscFunctionReturn(0);
874 }
875 
876 #undef __FUNCT__
877 #define __FUNCT__ "DMCoarsen_Composite"
878 PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
879 {
880   PetscErrorCode         ierr;
881   struct DMCompositeLink *next;
882   DM_Composite           *com = (DM_Composite*)dmi->data;
883   DM                     dm;
884 
885   PetscFunctionBegin;
886   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
887   next = com->next;
888   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
889 
890   /* loop over packed objects, handling one at at time */
891   while (next) {
892     ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr);
893     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
894     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
895     next = next->next;
896   }
897   PetscFunctionReturn(0);
898 }
899 
900 #undef __FUNCT__
901 #define __FUNCT__ "DMCreateInterpolation_Composite"
902 PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
903 {
904   PetscErrorCode         ierr;
905   PetscInt               m,n,M,N,nDM,i;
906   struct DMCompositeLink *nextc;
907   struct DMCompositeLink *nextf;
908   Vec                    gcoarse,gfine;
909   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
910   DM_Composite           *comfine = (DM_Composite*)fine->data;
911   Mat                    *mats;
912 
913   PetscFunctionBegin;
914   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
915   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
916   /* use global vectors only for determining matrix layout */
917   ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
918   ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr);
919   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
920   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
921   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
922   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
923   ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
924   ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr);
925 
926   nDM = comfine->nDM;
927   if (nDM != comcoarse->nDM) SETERRQ2(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %D entries, but coarse has %D",nDM,comcoarse->nDM);
928   ierr = PetscMalloc(nDM*nDM*sizeof(Mat),&mats);CHKERRQ(ierr);
929   ierr = PetscMemzero(mats,nDM*nDM*sizeof(Mat));CHKERRQ(ierr);
930 
931   /* loop over packed objects, handling one at at time */
932   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
933     ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],PETSC_NULL);CHKERRQ(ierr);
934   }
935   ierr = MatCreateNest(((PetscObject)fine)->comm,nDM,PETSC_NULL,nDM,PETSC_NULL,mats,A);CHKERRQ(ierr);
936   for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);}
937   ierr = PetscFree(mats);CHKERRQ(ierr);
938   PetscFunctionReturn(0);
939 }
940 
941 #undef __FUNCT__
942 #define __FUNCT__ "DMCreateLocalToGlobalMapping_Composite"
943 static PetscErrorCode DMCreateLocalToGlobalMapping_Composite(DM dm)
944 {
945   DM_Composite           *com = (DM_Composite*)dm->data;
946   ISLocalToGlobalMapping *ltogs;
947   PetscInt               i;
948   PetscErrorCode         ierr;
949 
950   PetscFunctionBegin;
951   /* Set the ISLocalToGlobalMapping on the new matrix */
952   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
953   ierr = ISLocalToGlobalMappingConcatenate(((PetscObject)dm)->comm,com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr);
954   for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
955   ierr = PetscFree(ltogs);CHKERRQ(ierr);
956   PetscFunctionReturn(0);
957 }
958 
959 
960 #undef __FUNCT__
961 #define __FUNCT__ "DMCreateColoring_Composite"
962 PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
963 {
964   PetscErrorCode         ierr;
965   PetscInt               n,i,cnt;
966   ISColoringValue        *colors;
967   PetscBool              dense = PETSC_FALSE;
968   ISColoringValue        maxcol = 0;
969   DM_Composite           *com = (DM_Composite*)dm->data;
970 
971   PetscFunctionBegin;
972   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
973   if (ctype == IS_COLORING_GHOSTED) {
974     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" );
975   } else if (ctype == IS_COLORING_GLOBAL) {
976     n = com->n;
977   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
978   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
979 
980   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
981   if (dense) {
982     for (i=0; i<n; i++) {
983       colors[i] = (ISColoringValue)(com->rstart + i);
984     }
985     maxcol = com->N;
986   } else {
987     struct DMCompositeLink *next = com->next;
988     PetscMPIInt            rank;
989 
990     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
991     cnt  = 0;
992     while (next) {
993       ISColoring     lcoloring;
994 
995       ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
996       for (i=0; i<lcoloring->N; i++) {
997         colors[cnt++] = maxcol + lcoloring->colors[i];
998       }
999       maxcol += lcoloring->n;
1000       ierr = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
1001       next = next->next;
1002     }
1003   }
1004   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
1005   PetscFunctionReturn(0);
1006 }
1007 
1008 #undef __FUNCT__
1009 #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
1010 PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1011 {
1012   PetscErrorCode         ierr;
1013   struct DMCompositeLink *next;
1014   PetscInt               cnt = 3;
1015   PetscMPIInt            rank;
1016   PetscScalar            *garray,*larray;
1017   DM_Composite           *com = (DM_Composite*)dm->data;
1018 
1019   PetscFunctionBegin;
1020   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1021   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
1022   next = com->next;
1023   if (!com->setup) {
1024     ierr = DMSetUp(dm);CHKERRQ(ierr);
1025   }
1026   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
1027   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
1028   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
1029 
1030   /* loop over packed objects, handling one at at time */
1031   while (next) {
1032     Vec      local,global;
1033     PetscInt N;
1034 
1035     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
1036     ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
1037     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
1038     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
1039     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
1040     ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
1041     ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
1042     ierr = VecResetArray(global);CHKERRQ(ierr);
1043     ierr = VecResetArray(local);CHKERRQ(ierr);
1044     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
1045     ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
1046     cnt++;
1047     larray += next->nlocal;
1048     next    = next->next;
1049   }
1050 
1051   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
1052   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
1053   PetscFunctionReturn(0);
1054 }
1055 
1056 #undef __FUNCT__
1057 #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
1058 PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1059 {
1060   PetscFunctionBegin;
1061   PetscFunctionReturn(0);
1062 }
1063 
1064 EXTERN_C_BEGIN
1065 #undef __FUNCT__
1066 #define __FUNCT__ "DMCreate_Composite"
1067 PetscErrorCode  DMCreate_Composite(DM p)
1068 {
1069   PetscErrorCode ierr;
1070   DM_Composite   *com;
1071 
1072   PetscFunctionBegin;
1073   ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1074   p->data = com;
1075   ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1076   com->n            = 0;
1077   com->next         = PETSC_NULL;
1078   com->nDM          = 0;
1079 
1080   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1081   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1082   p->ops->createlocaltoglobalmapping      = DMCreateLocalToGlobalMapping_Composite;
1083   p->ops->createlocaltoglobalmappingblock = 0;
1084   p->ops->refine                          = DMRefine_Composite;
1085   p->ops->coarsen                         = DMCoarsen_Composite;
1086   p->ops->getinterpolation                = DMCreateInterpolation_Composite;
1087   p->ops->getmatrix                       = DMCreateMatrix_Composite;
1088   p->ops->getcoloring                     = DMCreateColoring_Composite;
1089   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1090   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1091   p->ops->destroy                         = DMDestroy_Composite;
1092   p->ops->view                            = DMView_Composite;
1093   p->ops->setup                           = DMSetUp_Composite;
1094   PetscFunctionReturn(0);
1095 }
1096 EXTERN_C_END
1097 
1098 #undef __FUNCT__
1099 #define __FUNCT__ "DMCompositeCreate"
1100 /*@C
1101     DMCompositeCreate - Creates a vector packer, used to generate "composite"
1102       vectors made up of several subvectors.
1103 
1104     Collective on MPI_Comm
1105 
1106     Input Parameter:
1107 .   comm - the processors that will share the global vector
1108 
1109     Output Parameters:
1110 .   packer - the packer object
1111 
1112     Level: advanced
1113 
1114 .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(),
1115          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
1116          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
1117 
1118 @*/
1119 PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
1120 {
1121   PetscErrorCode ierr;
1122 
1123   PetscFunctionBegin;
1124   PetscValidPointer(packer,2);
1125   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1126   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
1127   PetscFunctionReturn(0);
1128 }
1129