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