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