xref: /petsc/src/dm/impls/composite/pack.c (revision 907376e66e56134b42a3e3ff22609655dbbb89d0)
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       seperate 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 
185   PetscFunctionBegin;
186   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
187   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
188   next = com->next;
189   if (!com->setup) {
190     ierr = DMSetUp(dm);CHKERRQ(ierr);
191   }
192 
193   /* loop over packed objects, handling one at at time */
194   va_start(Argp,gvec);
195   while (next) {
196     Vec *vec;
197     vec = va_arg(Argp, Vec*);
198     if (vec) {
199       PetscScalar *array;
200       ierr = DMGetGlobalVector(next->dm,vec);CHKERRQ(ierr);
201       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
202       ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
203       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
204     }
205     next = next->next;
206   }
207   va_end(Argp);
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "DMCompositeGetAccessArray"
213 /*@C
214     DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
215        representation.
216 
217     Collective on DMComposite
218 
219     Input Parameters:
220 +    dm - the packer object
221 .    pvec - packed vector
222 .    nwanted - number of vectors wanted
223 -    wanted - sorted array of vectors wanted, or NULL to get all vectors
224 
225     Output Parameters:
226 .    vecs - array of requested global vectors (must be allocated)
227 
228     Notes: Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them
229 
230     Level: advanced
231 
232 .seealso: DMCompositeGetAccess(), DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
233 @*/
234 PetscErrorCode  DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
235 {
236   PetscErrorCode         ierr;
237   struct DMCompositeLink *link;
238   PetscInt               i,wnum;
239   DM_Composite           *com = (DM_Composite*)dm->data;
240 
241   PetscFunctionBegin;
242   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
243   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
244   if (!com->setup) {
245     ierr = DMSetUp(dm);CHKERRQ(ierr);
246   }
247 
248   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
249     if (!wanted || i == wanted[wnum]) {
250       PetscScalar *array;
251       Vec v;
252       ierr = DMGetGlobalVector(link->dm,&v);CHKERRQ(ierr);
253       ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
254       ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
255       ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
256       vecs[wnum++] = v;
257     }
258   }
259   PetscFunctionReturn(0);
260 }
261 
262 #undef __FUNCT__
263 #define __FUNCT__ "DMCompositeRestoreAccess"
264 /*@C
265     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
266        representation.
267 
268     Collective on DMComposite
269 
270     Input Parameters:
271 +    dm - the packer object
272 .    gvec - the global vector
273 -    Vec* ... - the individual parallel vectors, NULL for those that are not needed
274 
275     Level: advanced
276 
277 .seealso  DMCompositeAddDM(), DMCreateGlobalVector(),
278          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
279          DMCompositeRestoreAccess(), DMCompositeGetAccess()
280 
281 @*/
282 PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
283 {
284   va_list                Argp;
285   PetscErrorCode         ierr;
286   struct DMCompositeLink *next;
287   DM_Composite           *com = (DM_Composite*)dm->data;
288 
289   PetscFunctionBegin;
290   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
291   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
292   next = com->next;
293   if (!com->setup) {
294     ierr = DMSetUp(dm);CHKERRQ(ierr);
295   }
296 
297   /* loop over packed objects, handling one at at time */
298   va_start(Argp,gvec);
299   while (next) {
300     Vec *vec;
301     vec = va_arg(Argp, Vec*);
302     if (vec) {
303       ierr = VecResetArray(*vec);CHKERRQ(ierr);
304       ierr = DMRestoreGlobalVector(next->dm,vec);CHKERRQ(ierr);
305     }
306     next = next->next;
307   }
308   va_end(Argp);
309   PetscFunctionReturn(0);
310 }
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "DMCompositeRestoreAccessArray"
314 /*@C
315     DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray()
316 
317     Collective on DMComposite
318 
319     Input Parameters:
320 +    dm - the packer object
321 .    pvec - packed vector
322 .    nwanted - number of vectors wanted
323 .    wanted - sorted array of vectors wanted, or NULL to get all vectors
324 -    vecs - array of global vectors to return
325 
326     Level: advanced
327 
328 .seealso: DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), DMCompositeScatter(), DMCompositeGather()
329 @*/
330 PetscErrorCode  DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
331 {
332   PetscErrorCode         ierr;
333   struct DMCompositeLink *link;
334   PetscInt               i,wnum;
335   DM_Composite           *com = (DM_Composite*)dm->data;
336 
337   PetscFunctionBegin;
338   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
339   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
340   if (!com->setup) {
341     ierr = DMSetUp(dm);CHKERRQ(ierr);
342   }
343 
344   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
345     if (!wanted || i == wanted[wnum]) {
346       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
347       ierr = DMRestoreGlobalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
348       wnum++;
349     }
350   }
351   PetscFunctionReturn(0);
352 }
353 
354 #undef __FUNCT__
355 #define __FUNCT__ "DMCompositeScatter"
356 /*@C
357     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
358 
359     Collective on DMComposite
360 
361     Input Parameters:
362 +    dm - the packer object
363 .    gvec - the global vector
364 -    Vec ... - the individual sequential vectors, NULL for those that are not needed
365 
366     Level: advanced
367 
368     Notes:
369     DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is
370     accessible from Fortran.
371 
372 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
373          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
374          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
375          DMCompositeScatterArray()
376 
377 @*/
378 PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
379 {
380   va_list                Argp;
381   PetscErrorCode         ierr;
382   struct DMCompositeLink *next;
383   PetscInt               cnt;
384   DM_Composite           *com = (DM_Composite*)dm->data;
385 
386   PetscFunctionBegin;
387   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
388   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
389   if (!com->setup) {
390     ierr = DMSetUp(dm);CHKERRQ(ierr);
391   }
392 
393   /* loop over packed objects, handling one at at time */
394   va_start(Argp,gvec);
395   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
396     Vec local;
397     local = va_arg(Argp, Vec);
398     if (local) {
399       Vec         global;
400       PetscScalar *array;
401       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
402       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
403       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
404       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
405       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
406       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
407       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
408       ierr = VecResetArray(global);CHKERRQ(ierr);
409       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
410     }
411   }
412   va_end(Argp);
413   PetscFunctionReturn(0);
414 }
415 
416 #undef __FUNCT__
417 #define __FUNCT__ "DMCompositeScatterArray"
418 /*@
419     DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
420 
421     Collective on DMComposite
422 
423     Input Parameters:
424 +    dm - the packer object
425 .    gvec - the global vector
426 .    lvecs - array of local vectors, NULL for any that are not needed
427 
428     Level: advanced
429 
430     Note:
431     This is a non-variadic alternative to DMCompositeScatter()
432 
433 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector()
434          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
435          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
436 
437 @*/
438 PetscErrorCode  DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs)
439 {
440   PetscErrorCode         ierr;
441   struct DMCompositeLink *next;
442   PetscInt               i;
443   DM_Composite           *com = (DM_Composite*)dm->data;
444 
445   PetscFunctionBegin;
446   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
447   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
448   if (!com->setup) {
449     ierr = DMSetUp(dm);CHKERRQ(ierr);
450   }
451 
452   /* loop over packed objects, handling one at at time */
453   for (i=0,next=com->next; next; next=next->next,i++) {
454     if (lvecs[i]) {
455       Vec         global;
456       PetscScalar *array;
457       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
458       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
459       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
460       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
461       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
462       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
463       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
464       ierr = VecResetArray(global);CHKERRQ(ierr);
465       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
466     }
467   }
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "DMCompositeGather"
473 /*@C
474     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
475 
476     Collective on DMComposite
477 
478     Input Parameter:
479 +    dm - the packer object
480 .    gvec - the global vector
481 .    imode - INSERT_VALUES or ADD_VALUES
482 -    Vec ... - the individual sequential vectors, NULL for any that are not needed
483 
484     Level: advanced
485 
486 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
487          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
488          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
489 
490 @*/
491 PetscErrorCode  DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
492 {
493   va_list                Argp;
494   PetscErrorCode         ierr;
495   struct DMCompositeLink *next;
496   DM_Composite           *com = (DM_Composite*)dm->data;
497   PetscInt               cnt;
498 
499   PetscFunctionBegin;
500   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
501   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
502   if (!com->setup) {
503     ierr = DMSetUp(dm);CHKERRQ(ierr);
504   }
505 
506   /* loop over packed objects, handling one at at time */
507   va_start(Argp,imode);
508   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
509     Vec local;
510     local = va_arg(Argp, Vec);
511     if (local) {
512       PetscScalar *array;
513       Vec         global;
514       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
515       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
516       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
517       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
518       ierr = DMLocalToGlobalBegin(next->dm,local,imode,global);CHKERRQ(ierr);
519       ierr = DMLocalToGlobalEnd(next->dm,local,imode,global);CHKERRQ(ierr);
520       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
521       ierr = VecResetArray(global);CHKERRQ(ierr);
522       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
523     }
524   }
525   va_end(Argp);
526   PetscFunctionReturn(0);
527 }
528 
529 #undef __FUNCT__
530 #define __FUNCT__ "DMCompositeGatherArray"
531 /*@
532     DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
533 
534     Collective on DMComposite
535 
536     Input Parameter:
537 +    dm - the packer object
538 .    gvec - the global vector
539 .    imode - INSERT_VALUES or ADD_VALUES
540 -    lvecs - the individual sequential vectors, NULL for any that are not needed
541 
542     Level: advanced
543 
544     Notes:
545     This is a non-variadic alternative to DMCompositeGather().
546 
547 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
548          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
549          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries(),
550 @*/
551 PetscErrorCode  DMCompositeGatherArray(DM dm,Vec gvec,InsertMode imode,Vec *lvecs)
552 {
553   PetscErrorCode         ierr;
554   struct DMCompositeLink *next;
555   DM_Composite           *com = (DM_Composite*)dm->data;
556   PetscInt               i;
557 
558   PetscFunctionBegin;
559   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
560   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
561   if (!com->setup) {
562     ierr = DMSetUp(dm);CHKERRQ(ierr);
563   }
564 
565   /* loop over packed objects, handling one at at time */
566   for (next=com->next,i=0; next; next=next->next,i++) {
567     if (lvecs[i]) {
568       PetscScalar *array;
569       Vec         global;
570       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
571       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
572       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
573       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
574       ierr = DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
575       ierr = DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
576       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
577       ierr = VecResetArray(global);CHKERRQ(ierr);
578       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
579     }
580   }
581   PetscFunctionReturn(0);
582 }
583 
584 #undef __FUNCT__
585 #define __FUNCT__ "DMCompositeAddDM"
586 /*@C
587     DMCompositeAddDM - adds a DM  vector to a DMComposite
588 
589     Collective on DMComposite
590 
591     Input Parameter:
592 +    dm - the packer object
593 -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
594 
595     Level: advanced
596 
597 .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
598          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
599          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
600 
601 @*/
602 PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
603 {
604   PetscErrorCode         ierr;
605   PetscInt               n,nlocal;
606   struct DMCompositeLink *mine,*next;
607   Vec                    global,local;
608   DM_Composite           *com = (DM_Composite*)dmc->data;
609 
610   PetscFunctionBegin;
611   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
612   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
613   next = com->next;
614   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
615 
616   /* create new link */
617   ierr = PetscNew(&mine);CHKERRQ(ierr);
618   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
619   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
620   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
621   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
622   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
623   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
624   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
625 
626   mine->n      = n;
627   mine->nlocal = nlocal;
628   mine->dm     = dm;
629   mine->next   = NULL;
630   com->n      += n;
631 
632   /* add to end of list */
633   if (!next) com->next = mine;
634   else {
635     while (next->next) next = next->next;
636     next->next = mine;
637   }
638   com->nDM++;
639   com->nmine++;
640   PetscFunctionReturn(0);
641 }
642 
643 #include <petscdraw.h>
644 PETSC_EXTERN PetscErrorCode  VecView_MPI(Vec,PetscViewer);
645 #undef __FUNCT__
646 #define __FUNCT__ "VecView_DMComposite"
647 PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
648 {
649   DM                     dm;
650   PetscErrorCode         ierr;
651   struct DMCompositeLink *next;
652   PetscBool              isdraw;
653   DM_Composite           *com;
654 
655   PetscFunctionBegin;
656   ierr = VecGetDM(gvec, &dm);CHKERRQ(ierr);
657   if (!dm) SETERRQ(PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
658   com  = (DM_Composite*)dm->data;
659   next = com->next;
660 
661   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
662   if (!isdraw) {
663     /* do I really want to call this? */
664     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
665   } else {
666     PetscInt cnt = 0;
667 
668     /* loop over packed objects, handling one at at time */
669     while (next) {
670       Vec         vec;
671       PetscScalar *array;
672       PetscInt    bs;
673 
674       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
675       ierr = DMGetGlobalVector(next->dm,&vec);CHKERRQ(ierr);
676       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
677       ierr = VecPlaceArray(vec,array+next->rstart);CHKERRQ(ierr);
678       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
679       ierr = VecView(vec,viewer);CHKERRQ(ierr);
680       ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
681       ierr = VecResetArray(vec);CHKERRQ(ierr);
682       ierr = DMRestoreGlobalVector(next->dm,&vec);CHKERRQ(ierr);
683       ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
684       cnt += bs;
685       next = next->next;
686     }
687     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
688   }
689   PetscFunctionReturn(0);
690 }
691 
692 #undef __FUNCT__
693 #define __FUNCT__ "DMCreateGlobalVector_Composite"
694 PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
695 {
696   PetscErrorCode ierr;
697   DM_Composite   *com = (DM_Composite*)dm->data;
698 
699   PetscFunctionBegin;
700   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
701   ierr = DMSetUp(dm);CHKERRQ(ierr);
702   ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),com->n,com->N,gvec);CHKERRQ(ierr);
703   ierr = VecSetDM(*gvec, dm);CHKERRQ(ierr);
704   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);CHKERRQ(ierr);
705   PetscFunctionReturn(0);
706 }
707 
708 #undef __FUNCT__
709 #define __FUNCT__ "DMCreateLocalVector_Composite"
710 PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
711 {
712   PetscErrorCode ierr;
713   DM_Composite   *com = (DM_Composite*)dm->data;
714 
715   PetscFunctionBegin;
716   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
717   if (!com->setup) {
718     ierr = DMSetUp(dm);CHKERRQ(ierr);
719   }
720   ierr = VecCreateSeq(PetscObjectComm((PetscObject)dm),com->nghost,lvec);CHKERRQ(ierr);
721   ierr = VecSetDM(*lvec, dm);CHKERRQ(ierr);
722   PetscFunctionReturn(0);
723 }
724 
725 #undef __FUNCT__
726 #define __FUNCT__ "DMCompositeGetISLocalToGlobalMappings"
727 /*@C
728     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
729 
730     Collective on DM
731 
732     Input Parameter:
733 .    dm - the packer object
734 
735     Output Parameters:
736 .    ltogs - the individual mappings for each packed vector. Note that this includes
737            all the ghost points that individual ghosted DMDA's may have.
738 
739     Level: advanced
740 
741     Notes:
742        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
743 
744 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
745          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
746          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
747 
748 @*/
749 PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
750 {
751   PetscErrorCode         ierr;
752   PetscInt               i,*idx,n,cnt;
753   struct DMCompositeLink *next;
754   PetscMPIInt            rank;
755   DM_Composite           *com = (DM_Composite*)dm->data;
756 
757   PetscFunctionBegin;
758   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
759   ierr = DMSetUp(dm);CHKERRQ(ierr);
760   ierr = PetscMalloc1(com->nDM,ltogs);CHKERRQ(ierr);
761   next = com->next;
762   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
763 
764   /* loop over packed objects, handling one at at time */
765   cnt = 0;
766   while (next) {
767     ISLocalToGlobalMapping ltog;
768     PetscMPIInt            size;
769     const PetscInt         *suboff,*indices;
770     Vec                    global;
771 
772     /* Get sub-DM global indices for each local dof */
773     ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
774     ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
775     ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
776     ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr);
777 
778     /* Get the offsets for the sub-DM global vector */
779     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
780     ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
781     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)global),&size);CHKERRQ(ierr);
782 
783     /* Shift the sub-DM definition of the global space to the composite global space */
784     for (i=0; i<n; i++) {
785       PetscInt subi = indices[i],lo = 0,hi = size,t;
786       /* Binary search to find which rank owns subi */
787       while (hi-lo > 1) {
788         t = lo + (hi-lo)/2;
789         if (suboff[t] > subi) hi = t;
790         else                  lo = t;
791       }
792       idx[i] = subi - suboff[lo] + next->grstarts[lo];
793     }
794     ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
795     ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
796     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
797     next = next->next;
798     cnt++;
799   }
800   PetscFunctionReturn(0);
801 }
802 
803 #undef __FUNCT__
804 #define __FUNCT__ "DMCompositeGetLocalISs"
805 /*@C
806    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
807 
808    Not Collective
809 
810    Input Arguments:
811 . dm - composite DM
812 
813    Output Arguments:
814 . is - array of serial index sets for each each component of the DMComposite
815 
816    Level: intermediate
817 
818    Notes:
819    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
820    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
821    scatter to a composite local vector.  The user should not typically need to know which is being done.
822 
823    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
824 
825    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
826 
827    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
828 
829 .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
830 @*/
831 PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
832 {
833   PetscErrorCode         ierr;
834   DM_Composite           *com = (DM_Composite*)dm->data;
835   struct DMCompositeLink *link;
836   PetscInt               cnt,start;
837 
838   PetscFunctionBegin;
839   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
840   PetscValidPointer(is,2);
841   ierr = PetscMalloc1(com->nmine,is);CHKERRQ(ierr);
842   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
843     PetscInt bs;
844     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
845     ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
846     ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
847   }
848   PetscFunctionReturn(0);
849 }
850 
851 #undef __FUNCT__
852 #define __FUNCT__ "DMCompositeGetGlobalISs"
853 /*@C
854     DMCompositeGetGlobalISs - Gets the index sets for each composed object
855 
856     Collective on DMComposite
857 
858     Input Parameter:
859 .    dm - the packer object
860 
861     Output Parameters:
862 .    is - the array of index sets
863 
864     Level: advanced
865 
866     Notes:
867        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
868 
869        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
870 
871        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
872        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
873        indices.
874 
875 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
876          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
877          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
878 
879 @*/
880 
881 PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
882 {
883   PetscErrorCode         ierr;
884   PetscInt               cnt = 0;
885   struct DMCompositeLink *next;
886   PetscMPIInt            rank;
887   DM_Composite           *com = (DM_Composite*)dm->data;
888 
889   PetscFunctionBegin;
890   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
891   ierr = PetscMalloc1(com->nDM,is);CHKERRQ(ierr);
892   next = com->next;
893   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
894 
895   /* loop over packed objects, handling one at at time */
896   while (next) {
897     ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt]);CHKERRQ(ierr);
898     if (dm->prob) {
899       MatNullSpace space;
900       Mat          pmat;
901       PetscObject  disc;
902       PetscInt     Nf;
903 
904       ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
905       if (cnt < Nf) {
906         ierr = PetscDSGetDiscretization(dm->prob, cnt, &disc);CHKERRQ(ierr);
907         ierr = PetscObjectQuery(disc, "nullspace", (PetscObject*) &space);CHKERRQ(ierr);
908         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space);CHKERRQ(ierr);}
909         ierr = PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space);CHKERRQ(ierr);
910         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space);CHKERRQ(ierr);}
911         ierr = PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
912         if (pmat) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat);CHKERRQ(ierr);}
913       }
914     }
915     cnt++;
916     next = next->next;
917   }
918   PetscFunctionReturn(0);
919 }
920 
921 #undef __FUNCT__
922 #define __FUNCT__ "DMCreateFieldIS_Composite"
923 PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields)
924 {
925   PetscInt       nDM;
926   DM             *dms;
927   PetscInt       i;
928   PetscErrorCode ierr;
929 
930   PetscFunctionBegin;
931   ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
932   if (numFields) *numFields = nDM;
933   ierr = DMCompositeGetGlobalISs(dm, fields);CHKERRQ(ierr);
934   if (fieldNames) {
935     ierr = PetscMalloc1(nDM, &dms);CHKERRQ(ierr);
936     ierr = PetscMalloc1(nDM, fieldNames);CHKERRQ(ierr);
937     ierr = DMCompositeGetEntriesArray(dm, dms);CHKERRQ(ierr);
938     for (i=0; i<nDM; i++) {
939       char       buf[256];
940       const char *splitname;
941 
942       /* Split naming precedence: object name, prefix, number */
943       splitname = ((PetscObject) dm)->name;
944       if (!splitname) {
945         ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr);
946         if (splitname) {
947           size_t len;
948           ierr                 = PetscStrncpy(buf,splitname,sizeof(buf));CHKERRQ(ierr);
949           buf[sizeof(buf) - 1] = 0;
950           ierr                 = PetscStrlen(buf,&len);CHKERRQ(ierr);
951           if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
952           splitname = buf;
953         }
954       }
955       if (!splitname) {
956         ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",i);CHKERRQ(ierr);
957         splitname = buf;
958       }
959       ierr = PetscStrallocpy(splitname,&(*fieldNames)[i]);CHKERRQ(ierr);
960     }
961     ierr = PetscFree(dms);CHKERRQ(ierr);
962   }
963   PetscFunctionReturn(0);
964 }
965 
966 /*
967  This could take over from DMCreateFieldIS(), as it is more general,
968  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
969  At this point it's probably best to be less intrusive, however.
970  */
971 #undef __FUNCT__
972 #define __FUNCT__ "DMCreateFieldDecomposition_Composite"
973 PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
974 {
975   PetscInt       nDM;
976   PetscInt       i;
977   PetscErrorCode ierr;
978 
979   PetscFunctionBegin;
980   ierr = DMCreateFieldIS_Composite(dm, len, namelist, islist);CHKERRQ(ierr);
981   if (dmlist) {
982     ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
983     ierr = PetscMalloc1(nDM, dmlist);CHKERRQ(ierr);
984     ierr = DMCompositeGetEntriesArray(dm, *dmlist);CHKERRQ(ierr);
985     for (i=0; i<nDM; i++) {
986       ierr = PetscObjectReference((PetscObject)((*dmlist)[i]));CHKERRQ(ierr);
987     }
988   }
989   PetscFunctionReturn(0);
990 }
991 
992 
993 
994 /* -------------------------------------------------------------------------------------*/
995 #undef __FUNCT__
996 #define __FUNCT__ "DMCompositeGetLocalVectors"
997 /*@C
998     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
999        Use DMCompositeRestoreLocalVectors() to return them.
1000 
1001     Not Collective
1002 
1003     Input Parameter:
1004 .    dm - the packer object
1005 
1006     Output Parameter:
1007 .   Vec ... - the individual sequential Vecs
1008 
1009     Level: advanced
1010 
1011 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
1012          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1013          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
1014 
1015 @*/
1016 PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
1017 {
1018   va_list                Argp;
1019   PetscErrorCode         ierr;
1020   struct DMCompositeLink *next;
1021   DM_Composite           *com = (DM_Composite*)dm->data;
1022 
1023   PetscFunctionBegin;
1024   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1025   next = com->next;
1026   /* loop over packed objects, handling one at at time */
1027   va_start(Argp,dm);
1028   while (next) {
1029     Vec *vec;
1030     vec = va_arg(Argp, Vec*);
1031     if (vec) {ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);}
1032     next = next->next;
1033   }
1034   va_end(Argp);
1035   PetscFunctionReturn(0);
1036 }
1037 
1038 #undef __FUNCT__
1039 #define __FUNCT__ "DMCompositeRestoreLocalVectors"
1040 /*@C
1041     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
1042 
1043     Not Collective
1044 
1045     Input Parameter:
1046 .    dm - the packer object
1047 
1048     Output Parameter:
1049 .   Vec ... - the individual sequential Vecs
1050 
1051     Level: advanced
1052 
1053 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
1054          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1055          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
1056 
1057 @*/
1058 PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
1059 {
1060   va_list                Argp;
1061   PetscErrorCode         ierr;
1062   struct DMCompositeLink *next;
1063   DM_Composite           *com = (DM_Composite*)dm->data;
1064 
1065   PetscFunctionBegin;
1066   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1067   next = com->next;
1068   /* loop over packed objects, handling one at at time */
1069   va_start(Argp,dm);
1070   while (next) {
1071     Vec *vec;
1072     vec = va_arg(Argp, Vec*);
1073     if (vec) {ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);}
1074     next = next->next;
1075   }
1076   va_end(Argp);
1077   PetscFunctionReturn(0);
1078 }
1079 
1080 /* -------------------------------------------------------------------------------------*/
1081 #undef __FUNCT__
1082 #define __FUNCT__ "DMCompositeGetEntries"
1083 /*@C
1084     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
1085 
1086     Not Collective
1087 
1088     Input Parameter:
1089 .    dm - the packer object
1090 
1091     Output Parameter:
1092 .   DM ... - the individual entries (DMs)
1093 
1094     Level: advanced
1095 
1096 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
1097          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1098          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
1099          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
1100 
1101 @*/
1102 PetscErrorCode  DMCompositeGetEntries(DM dm,...)
1103 {
1104   va_list                Argp;
1105   struct DMCompositeLink *next;
1106   DM_Composite           *com = (DM_Composite*)dm->data;
1107 
1108   PetscFunctionBegin;
1109   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1110   next = com->next;
1111   /* loop over packed objects, handling one at at time */
1112   va_start(Argp,dm);
1113   while (next) {
1114     DM *dmn;
1115     dmn = va_arg(Argp, DM*);
1116     if (dmn) *dmn = next->dm;
1117     next = next->next;
1118   }
1119   va_end(Argp);
1120   PetscFunctionReturn(0);
1121 }
1122 
1123 #undef __FUNCT__
1124 #define __FUNCT__ "DMCompositeGetEntriesArray"
1125 /*@C
1126     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
1127 
1128     Not Collective
1129 
1130     Input Parameter:
1131 .    dm - the packer object
1132 
1133     Output Parameter:
1134 .    dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs
1135 
1136     Level: advanced
1137 
1138 .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
1139          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
1140          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
1141          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
1142 
1143 @*/
1144 PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
1145 {
1146   struct DMCompositeLink *next;
1147   DM_Composite           *com = (DM_Composite*)dm->data;
1148   PetscInt               i;
1149 
1150   PetscFunctionBegin;
1151   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1152   /* loop over packed objects, handling one at at time */
1153   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
1154   PetscFunctionReturn(0);
1155 }
1156 
1157 #undef __FUNCT__
1158 #define __FUNCT__ "DMRefine_Composite"
1159 PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
1160 {
1161   PetscErrorCode         ierr;
1162   struct DMCompositeLink *next;
1163   DM_Composite           *com = (DM_Composite*)dmi->data;
1164   DM                     dm;
1165 
1166   PetscFunctionBegin;
1167   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1168   if (comm == MPI_COMM_NULL) {
1169     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
1170   }
1171   ierr = DMSetUp(dmi);CHKERRQ(ierr);
1172   next = com->next;
1173   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
1174 
1175   /* loop over packed objects, handling one at at time */
1176   while (next) {
1177     ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
1178     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
1179     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1180     next = next->next;
1181   }
1182   PetscFunctionReturn(0);
1183 }
1184 
1185 #undef __FUNCT__
1186 #define __FUNCT__ "DMCoarsen_Composite"
1187 PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
1188 {
1189   PetscErrorCode         ierr;
1190   struct DMCompositeLink *next;
1191   DM_Composite           *com = (DM_Composite*)dmi->data;
1192   DM                     dm;
1193 
1194   PetscFunctionBegin;
1195   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1196   ierr = DMSetUp(dmi);CHKERRQ(ierr);
1197   if (comm == MPI_COMM_NULL) {
1198     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
1199   }
1200   next = com->next;
1201   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
1202 
1203   /* loop over packed objects, handling one at at time */
1204   while (next) {
1205     ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr);
1206     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
1207     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
1208     next = next->next;
1209   }
1210   PetscFunctionReturn(0);
1211 }
1212 
1213 #undef __FUNCT__
1214 #define __FUNCT__ "DMCreateInterpolation_Composite"
1215 PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
1216 {
1217   PetscErrorCode         ierr;
1218   PetscInt               m,n,M,N,nDM,i;
1219   struct DMCompositeLink *nextc;
1220   struct DMCompositeLink *nextf;
1221   Vec                    gcoarse,gfine,*vecs;
1222   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
1223   DM_Composite           *comfine   = (DM_Composite*)fine->data;
1224   Mat                    *mats;
1225 
1226   PetscFunctionBegin;
1227   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
1228   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
1229   ierr = DMSetUp(coarse);CHKERRQ(ierr);
1230   ierr = DMSetUp(fine);CHKERRQ(ierr);
1231   /* use global vectors only for determining matrix layout */
1232   ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
1233   ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr);
1234   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
1235   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
1236   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
1237   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
1238   ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
1239   ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr);
1240 
1241   nDM = comfine->nDM;
1242   if (nDM != comcoarse->nDM) SETERRQ2(PetscObjectComm((PetscObject)fine),PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %D entries, but coarse has %D",nDM,comcoarse->nDM);
1243   ierr = PetscCalloc1(nDM*nDM,&mats);CHKERRQ(ierr);
1244   if (v) {
1245     ierr = PetscCalloc1(nDM,&vecs);CHKERRQ(ierr);
1246   }
1247 
1248   /* loop over packed objects, handling one at at time */
1249   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
1250     if (!v) {
1251       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL);CHKERRQ(ierr);
1252     } else {
1253       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr);
1254     }
1255   }
1256   ierr = MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A);CHKERRQ(ierr);
1257   if (v) {
1258     ierr = VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v);CHKERRQ(ierr);
1259   }
1260   for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);}
1261   ierr = PetscFree(mats);CHKERRQ(ierr);
1262   if (v) {
1263     for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);}
1264     ierr = PetscFree(vecs);CHKERRQ(ierr);
1265   }
1266   PetscFunctionReturn(0);
1267 }
1268 
1269 #undef __FUNCT__
1270 #define __FUNCT__ "DMGetLocalToGlobalMapping_Composite"
1271 static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
1272 {
1273   DM_Composite           *com = (DM_Composite*)dm->data;
1274   ISLocalToGlobalMapping *ltogs;
1275   PetscInt               i;
1276   PetscErrorCode         ierr;
1277 
1278   PetscFunctionBegin;
1279   /* Set the ISLocalToGlobalMapping on the new matrix */
1280   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
1281   ierr = ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr);
1282   for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
1283   ierr = PetscFree(ltogs);CHKERRQ(ierr);
1284   PetscFunctionReturn(0);
1285 }
1286 
1287 
1288 #undef __FUNCT__
1289 #define __FUNCT__ "DMCreateColoring_Composite"
1290 PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring)
1291 {
1292   PetscErrorCode  ierr;
1293   PetscInt        n,i,cnt;
1294   ISColoringValue *colors;
1295   PetscBool       dense  = PETSC_FALSE;
1296   ISColoringValue maxcol = 0;
1297   DM_Composite    *com   = (DM_Composite*)dm->data;
1298 
1299   PetscFunctionBegin;
1300   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1301   if (ctype == IS_COLORING_GHOSTED) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported");
1302   else if (ctype == IS_COLORING_GLOBAL) {
1303     n = com->n;
1304   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1305   ierr = PetscMalloc1(n,&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
1306 
1307   ierr = PetscOptionsGetBool(NULL,"-dmcomposite_dense_jacobian",&dense,NULL);CHKERRQ(ierr);
1308   if (dense) {
1309     for (i=0; i<n; i++) {
1310       colors[i] = (ISColoringValue)(com->rstart + i);
1311     }
1312     maxcol = com->N;
1313   } else {
1314     struct DMCompositeLink *next = com->next;
1315     PetscMPIInt            rank;
1316 
1317     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1318     cnt  = 0;
1319     while (next) {
1320       ISColoring lcoloring;
1321 
1322       ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring);CHKERRQ(ierr);
1323       for (i=0; i<lcoloring->N; i++) {
1324         colors[cnt++] = maxcol + lcoloring->colors[i];
1325       }
1326       maxcol += lcoloring->n;
1327       ierr    = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
1328       next    = next->next;
1329     }
1330   }
1331   ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,PETSC_OWN_POINTER,coloring);CHKERRQ(ierr);
1332   PetscFunctionReturn(0);
1333 }
1334 
1335 #undef __FUNCT__
1336 #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
1337 PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1338 {
1339   PetscErrorCode         ierr;
1340   struct DMCompositeLink *next;
1341   PetscInt               cnt = 3;
1342   PetscMPIInt            rank;
1343   PetscScalar            *garray,*larray;
1344   DM_Composite           *com = (DM_Composite*)dm->data;
1345 
1346   PetscFunctionBegin;
1347   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1348   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
1349   next = com->next;
1350   if (!com->setup) {
1351     ierr = DMSetUp(dm);CHKERRQ(ierr);
1352   }
1353   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1354   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
1355   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
1356 
1357   /* loop over packed objects, handling one at at time */
1358   while (next) {
1359     Vec      local,global;
1360     PetscInt N;
1361 
1362     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
1363     ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
1364     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
1365     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
1366     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
1367     ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
1368     ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
1369     ierr = VecResetArray(global);CHKERRQ(ierr);
1370     ierr = VecResetArray(local);CHKERRQ(ierr);
1371     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
1372     ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
1373     cnt++;
1374     larray += next->nlocal;
1375     next    = next->next;
1376   }
1377 
1378   ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr);
1379   ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr);
1380   PetscFunctionReturn(0);
1381 }
1382 
1383 #undef __FUNCT__
1384 #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
1385 PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
1386 {
1387   PetscFunctionBegin;
1388   PetscFunctionReturn(0);
1389 }
1390 
1391 /*MC
1392    DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs
1393 
1394   Level: intermediate
1395 
1396 .seealso: DMType, DM, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate()
1397 M*/
1398 
1399 
1400 #undef __FUNCT__
1401 #define __FUNCT__ "DMCreate_Composite"
1402 PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1403 {
1404   PetscErrorCode ierr;
1405   DM_Composite   *com;
1406 
1407   PetscFunctionBegin;
1408   ierr      = PetscNewLog(p,&com);CHKERRQ(ierr);
1409   p->data   = com;
1410   ierr      = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1411   com->n    = 0;
1412   com->next = NULL;
1413   com->nDM  = 0;
1414 
1415   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1416   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1417   p->ops->getlocaltoglobalmapping         = DMGetLocalToGlobalMapping_Composite;
1418   p->ops->createfieldis                   = DMCreateFieldIS_Composite;
1419   p->ops->createfielddecomposition        = DMCreateFieldDecomposition_Composite;
1420   p->ops->refine                          = DMRefine_Composite;
1421   p->ops->coarsen                         = DMCoarsen_Composite;
1422   p->ops->createinterpolation             = DMCreateInterpolation_Composite;
1423   p->ops->creatematrix                    = DMCreateMatrix_Composite;
1424   p->ops->getcoloring                     = DMCreateColoring_Composite;
1425   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1426   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1427   p->ops->destroy                         = DMDestroy_Composite;
1428   p->ops->view                            = DMView_Composite;
1429   p->ops->setup                           = DMSetUp_Composite;
1430   PetscFunctionReturn(0);
1431 }
1432 
1433 #undef __FUNCT__
1434 #define __FUNCT__ "DMCompositeCreate"
1435 /*@C
1436     DMCompositeCreate - Creates a vector packer, used to generate "composite"
1437       vectors made up of several subvectors.
1438 
1439     Collective on MPI_Comm
1440 
1441     Input Parameter:
1442 .   comm - the processors that will share the global vector
1443 
1444     Output Parameters:
1445 .   packer - the packer object
1446 
1447     Level: advanced
1448 
1449 .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(), DMCOMPOSITE,DMCreate()
1450          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
1451          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
1452 
1453 @*/
1454 PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
1455 {
1456   PetscErrorCode ierr;
1457 
1458   PetscFunctionBegin;
1459   PetscValidPointer(packer,2);
1460   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1461   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
1462   PetscFunctionReturn(0);
1463 }
1464