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