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