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