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