xref: /petsc/src/dm/impls/composite/pack.c (revision 7ac2b803e3b92ba6d298600f281e6c25f3c4a7f6)
147c6ae99SBarry Smith 
2ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h>       /*I  "petscdmcomposite.h"  I*/
3af0996ceSBarry Smith #include <petsc/private/isimpl.h>
42764a2aaSMatthew G. Knepley #include <petscds.h>
547c6ae99SBarry Smith 
647c6ae99SBarry Smith #undef __FUNCT__
747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeSetCoupling"
847c6ae99SBarry Smith /*@C
947c6ae99SBarry Smith     DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
10bebe2cf6SSatish Balay       separate components (DMs) in a DMto build the correct matrix nonzero structure.
1147c6ae99SBarry Smith 
1247c6ae99SBarry Smith 
1347c6ae99SBarry Smith     Logically Collective on MPI_Comm
1447c6ae99SBarry Smith 
1547c6ae99SBarry Smith     Input Parameter:
1647c6ae99SBarry Smith +   dm - the composite object
1747c6ae99SBarry Smith -   formcouplelocations - routine to set the nonzero locations in the matrix
1847c6ae99SBarry Smith 
1947c6ae99SBarry Smith     Level: advanced
2047c6ae99SBarry Smith 
211b2093e4SBarry Smith     Notes: See DMSetApplicationContext() and DMGetApplicationContext() for how to get user information into
2247c6ae99SBarry Smith         this routine
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith @*/
257087cfbeSBarry Smith PetscErrorCode  DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt))
2647c6ae99SBarry Smith {
2747c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
2847c6ae99SBarry Smith 
2947c6ae99SBarry Smith   PetscFunctionBegin;
3047c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
3147c6ae99SBarry Smith   PetscFunctionReturn(0);
3247c6ae99SBarry Smith }
3347c6ae99SBarry Smith 
3447c6ae99SBarry Smith #undef __FUNCT__
350c010503SBarry Smith #define __FUNCT__ "DMDestroy_Composite"
366bf464f9SBarry Smith PetscErrorCode  DMDestroy_Composite(DM dm)
3747c6ae99SBarry Smith {
3847c6ae99SBarry Smith   PetscErrorCode         ierr;
3947c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
4047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
4147c6ae99SBarry Smith 
4247c6ae99SBarry Smith   PetscFunctionBegin;
4347c6ae99SBarry Smith   next = com->next;
4447c6ae99SBarry Smith   while (next) {
4547c6ae99SBarry Smith     prev = next;
4647c6ae99SBarry Smith     next = next->next;
47fcfd50ebSBarry Smith     ierr = DMDestroy(&prev->dm);CHKERRQ(ierr);
4847c6ae99SBarry Smith     ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
4947c6ae99SBarry Smith     ierr = PetscFree(prev);CHKERRQ(ierr);
5047c6ae99SBarry Smith   }
51435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
52435a35e8SMatthew G Knepley   ierr = PetscFree(com);CHKERRQ(ierr);
5347c6ae99SBarry Smith   PetscFunctionReturn(0);
5447c6ae99SBarry Smith }
5547c6ae99SBarry Smith 
5647c6ae99SBarry Smith #undef __FUNCT__
570c010503SBarry Smith #define __FUNCT__ "DMView_Composite"
587087cfbeSBarry Smith PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
5947c6ae99SBarry Smith {
6047c6ae99SBarry Smith   PetscErrorCode ierr;
6147c6ae99SBarry Smith   PetscBool      iascii;
6247c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
6347c6ae99SBarry Smith 
6447c6ae99SBarry Smith   PetscFunctionBegin;
65251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6647c6ae99SBarry Smith   if (iascii) {
6747c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
6847c6ae99SBarry Smith     PetscInt               i;
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix");CHKERRQ(ierr);
719ae5db72SJed Brown     ierr = PetscViewerASCIIPrintf(v,"  contains %D DMs\n",com->nDM);CHKERRQ(ierr);
7247c6ae99SBarry Smith     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
7347c6ae99SBarry Smith     for (i=0; lnk; lnk=lnk->next,i++) {
749ae5db72SJed Brown       ierr = PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
7547c6ae99SBarry Smith       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
7647c6ae99SBarry Smith       ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
7747c6ae99SBarry Smith       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
7847c6ae99SBarry Smith     }
7947c6ae99SBarry Smith     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
8047c6ae99SBarry Smith   }
8147c6ae99SBarry Smith   PetscFunctionReturn(0);
8247c6ae99SBarry Smith }
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
8547c6ae99SBarry Smith #undef __FUNCT__
86d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp_Composite"
877087cfbeSBarry Smith PetscErrorCode  DMSetUp_Composite(DM dm)
8847c6ae99SBarry Smith {
8947c6ae99SBarry Smith   PetscErrorCode         ierr;
9047c6ae99SBarry Smith   PetscInt               nprev = 0;
9147c6ae99SBarry Smith   PetscMPIInt            rank,size;
9247c6ae99SBarry Smith   DM_Composite           *com  = (DM_Composite*)dm->data;
9347c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
9447c6ae99SBarry Smith   PetscLayout            map;
9547c6ae99SBarry Smith 
9647c6ae99SBarry Smith   PetscFunctionBegin;
97ce94432eSBarry Smith   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
98ce94432eSBarry Smith   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm),&map);CHKERRQ(ierr);
9947c6ae99SBarry Smith   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
10047c6ae99SBarry Smith   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
10147c6ae99SBarry Smith   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
10247c6ae99SBarry Smith   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
10347c6ae99SBarry Smith   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
1040298fd71SBarry Smith   ierr = PetscLayoutGetRange(map,&com->rstart,NULL);CHKERRQ(ierr);
105fcfd50ebSBarry Smith   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
10647c6ae99SBarry Smith 
1079ae5db72SJed Brown   /* now set the rstart for each linked vector */
108ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
109ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
11047c6ae99SBarry Smith   while (next) {
11147c6ae99SBarry Smith     next->rstart  = nprev;
11206ebdd98SJed Brown     nprev        += next->n;
11347c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
114785e854fSJed Brown     ierr          = PetscMalloc1(size,&next->grstarts);CHKERRQ(ierr);
115ce94432eSBarry Smith     ierr          = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
11647c6ae99SBarry Smith     next          = next->next;
11747c6ae99SBarry Smith   }
11847c6ae99SBarry Smith   com->setup = PETSC_TRUE;
11947c6ae99SBarry Smith   PetscFunctionReturn(0);
12047c6ae99SBarry Smith }
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith #undef __FUNCT__
12547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetNumberDM"
12673e31fe2SJed Brown /*@
12747c6ae99SBarry Smith     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
12847c6ae99SBarry Smith        representation.
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith     Not Collective
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith     Input Parameter:
13347c6ae99SBarry Smith .    dm - the packer object
13447c6ae99SBarry Smith 
13547c6ae99SBarry Smith     Output Parameter:
13647c6ae99SBarry Smith .     nDM - the number of DMs
13747c6ae99SBarry Smith 
13847c6ae99SBarry Smith     Level: beginner
13947c6ae99SBarry Smith 
14047c6ae99SBarry Smith @*/
1417087cfbeSBarry Smith PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
14247c6ae99SBarry Smith {
14347c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
1445fd66863SKarl Rupp 
14547c6ae99SBarry Smith   PetscFunctionBegin;
14647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
14747c6ae99SBarry Smith   *nDM = com->nDM;
14847c6ae99SBarry Smith   PetscFunctionReturn(0);
14947c6ae99SBarry Smith }
15047c6ae99SBarry Smith 
15147c6ae99SBarry Smith 
15247c6ae99SBarry Smith #undef __FUNCT__
15347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess"
15447c6ae99SBarry Smith /*@C
15547c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
15647c6ae99SBarry Smith        representation.
15747c6ae99SBarry Smith 
15847c6ae99SBarry Smith     Collective on DMComposite
15947c6ae99SBarry Smith 
1609ae5db72SJed Brown     Input Parameters:
16147c6ae99SBarry Smith +    dm - the packer object
1629ae5db72SJed Brown -    gvec - the global vector
1639ae5db72SJed Brown 
1649ae5db72SJed Brown     Output Parameters:
1650298fd71SBarry Smith .    Vec* ... - the packed parallel vectors, NULL for those that are not needed
16647c6ae99SBarry Smith 
16747c6ae99SBarry Smith     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
16847c6ae99SBarry Smith 
169f73e5cebSJed Brown     Fortran Notes:
170f73e5cebSJed Brown 
171f73e5cebSJed Brown     Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
172f73e5cebSJed Brown     or use the alternative interface DMCompositeGetAccessArray().
173f73e5cebSJed Brown 
17447c6ae99SBarry Smith     Level: advanced
17547c6ae99SBarry Smith 
176f73e5cebSJed Brown .seealso: DMCompositeGetEntries(), DMCompositeScatter()
17747c6ae99SBarry Smith @*/
1787087cfbeSBarry Smith PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
17947c6ae99SBarry Smith {
18047c6ae99SBarry Smith   va_list                Argp;
18147c6ae99SBarry Smith   PetscErrorCode         ierr;
18247c6ae99SBarry Smith   struct DMCompositeLink *next;
18347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
1845edff71fSBarry Smith   PetscInt               readonly;
18547c6ae99SBarry Smith 
18647c6ae99SBarry Smith   PetscFunctionBegin;
18747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
18847c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
18947c6ae99SBarry Smith   next = com->next;
19047c6ae99SBarry Smith   if (!com->setup) {
191d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
19247c6ae99SBarry Smith   }
19347c6ae99SBarry Smith 
1945edff71fSBarry Smith   ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr);
19547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
19647c6ae99SBarry Smith   va_start(Argp,gvec);
19747c6ae99SBarry Smith   while (next) {
19847c6ae99SBarry Smith     Vec *vec;
19947c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
2009ae5db72SJed Brown     if (vec) {
2019ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,vec);CHKERRQ(ierr);
2025edff71fSBarry Smith       if (readonly) {
2035edff71fSBarry Smith         const PetscScalar *array;
2045edff71fSBarry Smith         ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
2055edff71fSBarry Smith         ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
2065edff71fSBarry Smith         ierr = VecLockPush(*vec);CHKERRQ(ierr);
2075edff71fSBarry Smith         ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
2085edff71fSBarry Smith       } else {
2095edff71fSBarry Smith         PetscScalar *array;
2109ae5db72SJed Brown         ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
2119ae5db72SJed Brown         ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
2129ae5db72SJed Brown         ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
21347c6ae99SBarry Smith       }
2145edff71fSBarry Smith     }
21547c6ae99SBarry Smith     next = next->next;
21647c6ae99SBarry Smith   }
21747c6ae99SBarry Smith   va_end(Argp);
21847c6ae99SBarry Smith   PetscFunctionReturn(0);
21947c6ae99SBarry Smith }
22047c6ae99SBarry Smith 
22147c6ae99SBarry Smith #undef __FUNCT__
222f73e5cebSJed Brown #define __FUNCT__ "DMCompositeGetAccessArray"
223f73e5cebSJed Brown /*@C
224f73e5cebSJed Brown     DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
225f73e5cebSJed Brown        representation.
226f73e5cebSJed Brown 
227f73e5cebSJed Brown     Collective on DMComposite
228f73e5cebSJed Brown 
229f73e5cebSJed Brown     Input Parameters:
230f73e5cebSJed Brown +    dm - the packer object
231f73e5cebSJed Brown .    pvec - packed vector
232f73e5cebSJed Brown .    nwanted - number of vectors wanted
2330298fd71SBarry Smith -    wanted - sorted array of vectors wanted, or NULL to get all vectors
234f73e5cebSJed Brown 
235f73e5cebSJed Brown     Output Parameters:
236f73e5cebSJed Brown .    vecs - array of requested global vectors (must be allocated)
237f73e5cebSJed Brown 
238f73e5cebSJed Brown     Notes: Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them
239f73e5cebSJed Brown 
240f73e5cebSJed Brown     Level: advanced
241f73e5cebSJed Brown 
242f73e5cebSJed Brown .seealso: DMCompositeGetAccess(), DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
243f73e5cebSJed Brown @*/
244f73e5cebSJed Brown PetscErrorCode  DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
245f73e5cebSJed Brown {
246f73e5cebSJed Brown   PetscErrorCode         ierr;
247f73e5cebSJed Brown   struct DMCompositeLink *link;
248f73e5cebSJed Brown   PetscInt               i,wnum;
249f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
250bee642f7SBarry Smith   PetscInt               readonly;
251f73e5cebSJed Brown 
252f73e5cebSJed Brown   PetscFunctionBegin;
253f73e5cebSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
254f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
255f73e5cebSJed Brown   if (!com->setup) {
256f73e5cebSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
257f73e5cebSJed Brown   }
258f73e5cebSJed Brown 
259bee642f7SBarry Smith   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
260f73e5cebSJed Brown   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
261f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
262f73e5cebSJed Brown       Vec v;
263f73e5cebSJed Brown       ierr = DMGetGlobalVector(link->dm,&v);CHKERRQ(ierr);
264bee642f7SBarry Smith       if (readonly) {
265bee642f7SBarry Smith         const PetscScalar *array;
266bee642f7SBarry Smith         ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr);
267bee642f7SBarry Smith         ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
268bee642f7SBarry Smith         ierr = VecLockPush(v);CHKERRQ(ierr);
269bee642f7SBarry Smith         ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr);
270bee642f7SBarry Smith       } else {
271bee642f7SBarry Smith         PetscScalar *array;
272f73e5cebSJed Brown         ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
273f73e5cebSJed Brown         ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
274f73e5cebSJed Brown         ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
275bee642f7SBarry Smith       }
276f73e5cebSJed Brown       vecs[wnum++] = v;
277f73e5cebSJed Brown     }
278f73e5cebSJed Brown   }
279f73e5cebSJed Brown   PetscFunctionReturn(0);
280f73e5cebSJed Brown }
281f73e5cebSJed Brown 
282f73e5cebSJed Brown #undef __FUNCT__
283*7ac2b803SAlex Fikl #define __FUNCT__ "DMCompositeGetLocalAccessArray"
284*7ac2b803SAlex Fikl /*@C
285*7ac2b803SAlex Fikl     DMCompositeGetLocalAccessArray - Allows one to access the individual
286*7ac2b803SAlex Fikl     packed vectors in their local representation.
287*7ac2b803SAlex Fikl 
288*7ac2b803SAlex Fikl     Collective on DMComposite.
289*7ac2b803SAlex Fikl 
290*7ac2b803SAlex Fikl     Input Parameters:
291*7ac2b803SAlex Fikl +    dm - the packer object
292*7ac2b803SAlex Fikl .    pvec - packed vector
293*7ac2b803SAlex Fikl .    nwanted - number of vectors wanted
294*7ac2b803SAlex Fikl -    wanted - sorted array of vectors wanted, or NULL to get all vectors
295*7ac2b803SAlex Fikl 
296*7ac2b803SAlex Fikl     Output Parameters:
297*7ac2b803SAlex Fikl .    vecs - array of requested local vectors (must be allocated)
298*7ac2b803SAlex Fikl 
299*7ac2b803SAlex Fikl     Notes: Use DMCompositeRestoreLocalAccessArray() to return the vectors
300*7ac2b803SAlex Fikl     when you no longer need them.
301*7ac2b803SAlex Fikl 
302*7ac2b803SAlex Fikl     Level: advanced
303*7ac2b803SAlex Fikl 
304*7ac2b803SAlex Fikl .seealso: DMCompositeRestoreLocalAccessArray(), DMCompositeGetAccess(),
305*7ac2b803SAlex Fikl DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
306*7ac2b803SAlex Fikl @*/
307*7ac2b803SAlex Fikl PetscErrorCode  DMCompositeGetLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
308*7ac2b803SAlex Fikl {
309*7ac2b803SAlex Fikl   PetscErrorCode         ierr;
310*7ac2b803SAlex Fikl   struct DMCompositeLink *link;
311*7ac2b803SAlex Fikl   PetscInt               i,wnum;
312*7ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite*)dm->data;
313*7ac2b803SAlex Fikl   PetscInt               readonly;
314*7ac2b803SAlex Fikl   PetscInt               nlocal = 0;
315*7ac2b803SAlex Fikl 
316*7ac2b803SAlex Fikl   PetscFunctionBegin;
317*7ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
318*7ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
319*7ac2b803SAlex Fikl   if (!com->setup) {
320*7ac2b803SAlex Fikl     ierr = DMSetUp(dm);CHKERRQ(ierr);
321*7ac2b803SAlex Fikl   }
322*7ac2b803SAlex Fikl 
323*7ac2b803SAlex Fikl   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
324*7ac2b803SAlex Fikl   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
325*7ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
326*7ac2b803SAlex Fikl       Vec v;
327*7ac2b803SAlex Fikl       ierr = DMGetLocalVector(link->dm,&v);CHKERRQ(ierr);
328*7ac2b803SAlex Fikl       if (readonly) {
329*7ac2b803SAlex Fikl         const PetscScalar *array;
330*7ac2b803SAlex Fikl         ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr);
331*7ac2b803SAlex Fikl         ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr);
332*7ac2b803SAlex Fikl         ierr = VecLockPush(v);CHKERRQ(ierr);
333*7ac2b803SAlex Fikl         ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr);
334*7ac2b803SAlex Fikl       } else {
335*7ac2b803SAlex Fikl         PetscScalar *array;
336*7ac2b803SAlex Fikl         ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
337*7ac2b803SAlex Fikl         ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr);
338*7ac2b803SAlex Fikl         ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
339*7ac2b803SAlex Fikl       }
340*7ac2b803SAlex Fikl       vecs[wnum++] = v;
341*7ac2b803SAlex Fikl     }
342*7ac2b803SAlex Fikl 
343*7ac2b803SAlex Fikl     nlocal += link->nlocal;
344*7ac2b803SAlex Fikl   }
345*7ac2b803SAlex Fikl 
346*7ac2b803SAlex Fikl   PetscFunctionReturn(0);
347*7ac2b803SAlex Fikl }
348*7ac2b803SAlex Fikl 
349*7ac2b803SAlex Fikl #undef __FUNCT__
35047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess"
35147c6ae99SBarry Smith /*@C
352aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
35347c6ae99SBarry Smith        representation.
35447c6ae99SBarry Smith 
35547c6ae99SBarry Smith     Collective on DMComposite
35647c6ae99SBarry Smith 
3579ae5db72SJed Brown     Input Parameters:
35847c6ae99SBarry Smith +    dm - the packer object
35947c6ae99SBarry Smith .    gvec - the global vector
3600298fd71SBarry Smith -    Vec* ... - the individual parallel vectors, NULL for those that are not needed
36147c6ae99SBarry Smith 
36247c6ae99SBarry Smith     Level: advanced
36347c6ae99SBarry Smith 
3649ae5db72SJed Brown .seealso  DMCompositeAddDM(), DMCreateGlobalVector(),
3656eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
366aa219208SBarry Smith          DMCompositeRestoreAccess(), DMCompositeGetAccess()
36747c6ae99SBarry Smith 
36847c6ae99SBarry Smith @*/
3697087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
37047c6ae99SBarry Smith {
37147c6ae99SBarry Smith   va_list                Argp;
37247c6ae99SBarry Smith   PetscErrorCode         ierr;
37347c6ae99SBarry Smith   struct DMCompositeLink *next;
37447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
3755edff71fSBarry Smith   PetscInt               readonly;
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith   PetscFunctionBegin;
37847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
37947c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
38047c6ae99SBarry Smith   next = com->next;
38147c6ae99SBarry Smith   if (!com->setup) {
382d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
38347c6ae99SBarry Smith   }
38447c6ae99SBarry Smith 
3855edff71fSBarry Smith   ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr);
38647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
38747c6ae99SBarry Smith   va_start(Argp,gvec);
38847c6ae99SBarry Smith   while (next) {
38947c6ae99SBarry Smith     Vec *vec;
39047c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
3919ae5db72SJed Brown     if (vec) {
3929ae5db72SJed Brown       ierr = VecResetArray(*vec);CHKERRQ(ierr);
3935edff71fSBarry Smith       if (readonly) {
3945edff71fSBarry Smith         ierr = VecLockPop(*vec);CHKERRQ(ierr);
3955edff71fSBarry Smith       }
396bee642f7SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,vec);CHKERRQ(ierr);
39747c6ae99SBarry Smith     }
39847c6ae99SBarry Smith     next = next->next;
39947c6ae99SBarry Smith   }
40047c6ae99SBarry Smith   va_end(Argp);
40147c6ae99SBarry Smith   PetscFunctionReturn(0);
40247c6ae99SBarry Smith }
40347c6ae99SBarry Smith 
40447c6ae99SBarry Smith #undef __FUNCT__
405f73e5cebSJed Brown #define __FUNCT__ "DMCompositeRestoreAccessArray"
406f73e5cebSJed Brown /*@C
407f73e5cebSJed Brown     DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray()
408f73e5cebSJed Brown 
409f73e5cebSJed Brown     Collective on DMComposite
410f73e5cebSJed Brown 
411f73e5cebSJed Brown     Input Parameters:
412f73e5cebSJed Brown +    dm - the packer object
413f73e5cebSJed Brown .    pvec - packed vector
414f73e5cebSJed Brown .    nwanted - number of vectors wanted
4150298fd71SBarry Smith .    wanted - sorted array of vectors wanted, or NULL to get all vectors
416f73e5cebSJed Brown -    vecs - array of global vectors to return
417f73e5cebSJed Brown 
418f73e5cebSJed Brown     Level: advanced
419f73e5cebSJed Brown 
420f73e5cebSJed Brown .seealso: DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), DMCompositeScatter(), DMCompositeGather()
421f73e5cebSJed Brown @*/
422f73e5cebSJed Brown PetscErrorCode  DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
423f73e5cebSJed Brown {
424f73e5cebSJed Brown   PetscErrorCode         ierr;
425f73e5cebSJed Brown   struct DMCompositeLink *link;
426f73e5cebSJed Brown   PetscInt               i,wnum;
427f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
428bee642f7SBarry Smith   PetscInt               readonly;
429f73e5cebSJed Brown 
430f73e5cebSJed Brown   PetscFunctionBegin;
431f73e5cebSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
432f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
433f73e5cebSJed Brown   if (!com->setup) {
434f73e5cebSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
435f73e5cebSJed Brown   }
436f73e5cebSJed Brown 
437bee642f7SBarry Smith   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
438f73e5cebSJed Brown   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
439f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
440f73e5cebSJed Brown       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
441bee642f7SBarry Smith       if (readonly) {
442bee642f7SBarry Smith         ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr);
443bee642f7SBarry Smith       }
444f73e5cebSJed Brown       ierr = DMRestoreGlobalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
445f73e5cebSJed Brown       wnum++;
446f73e5cebSJed Brown     }
447f73e5cebSJed Brown   }
448f73e5cebSJed Brown   PetscFunctionReturn(0);
449f73e5cebSJed Brown }
450f73e5cebSJed Brown 
451f73e5cebSJed Brown #undef __FUNCT__
452*7ac2b803SAlex Fikl #define __FUNCT__ "DMCompositeRestoreLocalAccessArray"
453*7ac2b803SAlex Fikl /*@C
454*7ac2b803SAlex Fikl     DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray().
455*7ac2b803SAlex Fikl 
456*7ac2b803SAlex Fikl     Collective on DMComposite.
457*7ac2b803SAlex Fikl 
458*7ac2b803SAlex Fikl     Input Parameters:
459*7ac2b803SAlex Fikl +    dm - the packer object
460*7ac2b803SAlex Fikl .    pvec - packed vector
461*7ac2b803SAlex Fikl .    nwanted - number of vectors wanted
462*7ac2b803SAlex Fikl .    wanted - sorted array of vectors wanted, or NULL to restore all vectors
463*7ac2b803SAlex Fikl -    vecs - array of local vectors to return
464*7ac2b803SAlex Fikl 
465*7ac2b803SAlex Fikl     Level: advanced
466*7ac2b803SAlex Fikl 
467*7ac2b803SAlex Fikl     Notes:
468*7ac2b803SAlex Fikl     nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray()
469*7ac2b803SAlex Fikl     otherwise the call will fail.
470*7ac2b803SAlex Fikl 
471*7ac2b803SAlex Fikl .seealso: DMCompositeGetLocalAccessArray(), DMCompositeRestoreAccessArray(),
472*7ac2b803SAlex Fikl DMCompositeRestoreAccess(), DMCompositeRestoreEntries(),
473*7ac2b803SAlex Fikl DMCompositeScatter(), DMCompositeGather()
474*7ac2b803SAlex Fikl @*/
475*7ac2b803SAlex Fikl PetscErrorCode  DMCompositeRestoreLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
476*7ac2b803SAlex Fikl {
477*7ac2b803SAlex Fikl   PetscErrorCode         ierr;
478*7ac2b803SAlex Fikl   struct DMCompositeLink *link;
479*7ac2b803SAlex Fikl   PetscInt               i,wnum;
480*7ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite*)dm->data;
481*7ac2b803SAlex Fikl   PetscInt               readonly;
482*7ac2b803SAlex Fikl 
483*7ac2b803SAlex Fikl   PetscFunctionBegin;
484*7ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
485*7ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
486*7ac2b803SAlex Fikl   if (!com->setup) {
487*7ac2b803SAlex Fikl     ierr = DMSetUp(dm);CHKERRQ(ierr);
488*7ac2b803SAlex Fikl   }
489*7ac2b803SAlex Fikl 
490*7ac2b803SAlex Fikl   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
491*7ac2b803SAlex Fikl   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
492*7ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
493*7ac2b803SAlex Fikl       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
494*7ac2b803SAlex Fikl       if (readonly) {
495*7ac2b803SAlex Fikl         ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr);
496*7ac2b803SAlex Fikl       }
497*7ac2b803SAlex Fikl       ierr = DMRestoreLocalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
498*7ac2b803SAlex Fikl       wnum++;
499*7ac2b803SAlex Fikl     }
500*7ac2b803SAlex Fikl   }
501*7ac2b803SAlex Fikl   PetscFunctionReturn(0);
502*7ac2b803SAlex Fikl }
503*7ac2b803SAlex Fikl 
504*7ac2b803SAlex Fikl #undef __FUNCT__
50547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter"
50647c6ae99SBarry Smith /*@C
50747c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith     Collective on DMComposite
51047c6ae99SBarry Smith 
5119ae5db72SJed Brown     Input Parameters:
51247c6ae99SBarry Smith +    dm - the packer object
51347c6ae99SBarry Smith .    gvec - the global vector
5140298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for those that are not needed
51547c6ae99SBarry Smith 
51647c6ae99SBarry Smith     Level: advanced
51747c6ae99SBarry Smith 
5186f3c3dcfSJed Brown     Notes:
5196f3c3dcfSJed Brown     DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is
5206f3c3dcfSJed Brown     accessible from Fortran.
5216f3c3dcfSJed Brown 
5229ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
5236eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
52447c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
5256f3c3dcfSJed Brown          DMCompositeScatterArray()
52647c6ae99SBarry Smith 
52747c6ae99SBarry Smith @*/
5287087cfbeSBarry Smith PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
52947c6ae99SBarry Smith {
53047c6ae99SBarry Smith   va_list                Argp;
53147c6ae99SBarry Smith   PetscErrorCode         ierr;
53247c6ae99SBarry Smith   struct DMCompositeLink *next;
5338fd8f222SJed Brown   PetscInt               cnt;
53447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
53547c6ae99SBarry Smith 
53647c6ae99SBarry Smith   PetscFunctionBegin;
53747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
53847c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
53947c6ae99SBarry Smith   if (!com->setup) {
540d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
54147c6ae99SBarry Smith   }
54247c6ae99SBarry Smith 
54347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
54447c6ae99SBarry Smith   va_start(Argp,gvec);
5458fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
5469ae5db72SJed Brown     Vec local;
5479ae5db72SJed Brown     local = va_arg(Argp, Vec);
5489ae5db72SJed Brown     if (local) {
5499ae5db72SJed Brown       Vec               global;
5505edff71fSBarry Smith       const PetscScalar *array;
5519ae5db72SJed Brown       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
5529ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
5535edff71fSBarry Smith       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
5549ae5db72SJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
5559ae5db72SJed Brown       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
5569ae5db72SJed Brown       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
5575edff71fSBarry Smith       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
5589ae5db72SJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
5599ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
56047c6ae99SBarry Smith     }
56147c6ae99SBarry Smith   }
56247c6ae99SBarry Smith   va_end(Argp);
56347c6ae99SBarry Smith   PetscFunctionReturn(0);
56447c6ae99SBarry Smith }
56547c6ae99SBarry Smith 
56647c6ae99SBarry Smith #undef __FUNCT__
5676f3c3dcfSJed Brown #define __FUNCT__ "DMCompositeScatterArray"
5686f3c3dcfSJed Brown /*@
5696f3c3dcfSJed Brown     DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5706f3c3dcfSJed Brown 
5716f3c3dcfSJed Brown     Collective on DMComposite
5726f3c3dcfSJed Brown 
5736f3c3dcfSJed Brown     Input Parameters:
5746f3c3dcfSJed Brown +    dm - the packer object
5756f3c3dcfSJed Brown .    gvec - the global vector
5766f3c3dcfSJed Brown .    lvecs - array of local vectors, NULL for any that are not needed
5776f3c3dcfSJed Brown 
5786f3c3dcfSJed Brown     Level: advanced
5796f3c3dcfSJed Brown 
5806f3c3dcfSJed Brown     Note:
581907376e6SBarry Smith     This is a non-variadic alternative to DMCompositeScatter()
5826f3c3dcfSJed Brown 
5836f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector()
5846f3c3dcfSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
5856f3c3dcfSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
5866f3c3dcfSJed Brown 
5876f3c3dcfSJed Brown @*/
5886f3c3dcfSJed Brown PetscErrorCode  DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs)
5896f3c3dcfSJed Brown {
5906f3c3dcfSJed Brown   PetscErrorCode         ierr;
5916f3c3dcfSJed Brown   struct DMCompositeLink *next;
5926f3c3dcfSJed Brown   PetscInt               i;
5936f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
5946f3c3dcfSJed Brown 
5956f3c3dcfSJed Brown   PetscFunctionBegin;
5966f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5976f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
5986f3c3dcfSJed Brown   if (!com->setup) {
5996f3c3dcfSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
6006f3c3dcfSJed Brown   }
6016f3c3dcfSJed Brown 
6026f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
6036f3c3dcfSJed Brown   for (i=0,next=com->next; next; next=next->next,i++) {
6046f3c3dcfSJed Brown     if (lvecs[i]) {
6056f3c3dcfSJed Brown       Vec         global;
606c5d31e75SLisandro Dalcin       const PetscScalar *array;
6076f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
6086f3c3dcfSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
609c5d31e75SLisandro Dalcin       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
610c5d31e75SLisandro Dalcin       ierr = VecPlaceArray(global,(PetscScalar*)array+next->rstart);CHKERRQ(ierr);
6116f3c3dcfSJed Brown       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
6126f3c3dcfSJed Brown       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
613c5d31e75SLisandro Dalcin       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
6146f3c3dcfSJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
6156f3c3dcfSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
6166f3c3dcfSJed Brown     }
6176f3c3dcfSJed Brown   }
6186f3c3dcfSJed Brown   PetscFunctionReturn(0);
6196f3c3dcfSJed Brown }
6206f3c3dcfSJed Brown 
6216f3c3dcfSJed Brown #undef __FUNCT__
62247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather"
62347c6ae99SBarry Smith /*@C
62447c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
62547c6ae99SBarry Smith 
62647c6ae99SBarry Smith     Collective on DMComposite
62747c6ae99SBarry Smith 
62847c6ae99SBarry Smith     Input Parameter:
62947c6ae99SBarry Smith +    dm - the packer object
63047c6ae99SBarry Smith .    gvec - the global vector
631907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6320298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for any that are not needed
63347c6ae99SBarry Smith 
63447c6ae99SBarry Smith     Level: advanced
63547c6ae99SBarry Smith 
6369ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
6376eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
63847c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
63947c6ae99SBarry Smith 
64047c6ae99SBarry Smith @*/
6417087cfbeSBarry Smith PetscErrorCode  DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
64247c6ae99SBarry Smith {
64347c6ae99SBarry Smith   va_list                Argp;
64447c6ae99SBarry Smith   PetscErrorCode         ierr;
64547c6ae99SBarry Smith   struct DMCompositeLink *next;
64647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
6478fd8f222SJed Brown   PetscInt               cnt;
64847c6ae99SBarry Smith 
64947c6ae99SBarry Smith   PetscFunctionBegin;
65047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
65147c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
65247c6ae99SBarry Smith   if (!com->setup) {
653d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
65447c6ae99SBarry Smith   }
65547c6ae99SBarry Smith 
65647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
657df0c820aSJed Brown   va_start(Argp,imode);
6588fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
6599ae5db72SJed Brown     Vec local;
6609ae5db72SJed Brown     local = va_arg(Argp, Vec);
6619ae5db72SJed Brown     if (local) {
66247c6ae99SBarry Smith       PetscScalar *array;
6639ae5db72SJed Brown       Vec         global;
6649ae5db72SJed Brown       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
6659ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
6669ae5db72SJed Brown       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
6679ae5db72SJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
6689ae5db72SJed Brown       ierr = DMLocalToGlobalBegin(next->dm,local,imode,global);CHKERRQ(ierr);
6699ae5db72SJed Brown       ierr = DMLocalToGlobalEnd(next->dm,local,imode,global);CHKERRQ(ierr);
6709ae5db72SJed Brown       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
6719ae5db72SJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
6729ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
67347c6ae99SBarry Smith     }
67447c6ae99SBarry Smith   }
67547c6ae99SBarry Smith   va_end(Argp);
67647c6ae99SBarry Smith   PetscFunctionReturn(0);
67747c6ae99SBarry Smith }
67847c6ae99SBarry Smith 
67947c6ae99SBarry Smith #undef __FUNCT__
6806f3c3dcfSJed Brown #define __FUNCT__ "DMCompositeGatherArray"
6816f3c3dcfSJed Brown /*@
6826f3c3dcfSJed Brown     DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6836f3c3dcfSJed Brown 
6846f3c3dcfSJed Brown     Collective on DMComposite
6856f3c3dcfSJed Brown 
6866f3c3dcfSJed Brown     Input Parameter:
6876f3c3dcfSJed Brown +    dm - the packer object
6886f3c3dcfSJed Brown .    gvec - the global vector
689907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6906f3c3dcfSJed Brown -    lvecs - the individual sequential vectors, NULL for any that are not needed
6916f3c3dcfSJed Brown 
6926f3c3dcfSJed Brown     Level: advanced
6936f3c3dcfSJed Brown 
6946f3c3dcfSJed Brown     Notes:
6956f3c3dcfSJed Brown     This is a non-variadic alternative to DMCompositeGather().
6966f3c3dcfSJed Brown 
6976f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
6986f3c3dcfSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
6996f3c3dcfSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries(),
7006f3c3dcfSJed Brown @*/
7016f3c3dcfSJed Brown PetscErrorCode  DMCompositeGatherArray(DM dm,Vec gvec,InsertMode imode,Vec *lvecs)
7026f3c3dcfSJed Brown {
7036f3c3dcfSJed Brown   PetscErrorCode         ierr;
7046f3c3dcfSJed Brown   struct DMCompositeLink *next;
7056f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
7066f3c3dcfSJed Brown   PetscInt               i;
7076f3c3dcfSJed Brown 
7086f3c3dcfSJed Brown   PetscFunctionBegin;
7096f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
7106f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
7116f3c3dcfSJed Brown   if (!com->setup) {
7126f3c3dcfSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
7136f3c3dcfSJed Brown   }
7146f3c3dcfSJed Brown 
7156f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
7166f3c3dcfSJed Brown   for (next=com->next,i=0; next; next=next->next,i++) {
7176f3c3dcfSJed Brown     if (lvecs[i]) {
7186f3c3dcfSJed Brown       PetscScalar *array;
7196f3c3dcfSJed Brown       Vec         global;
7206f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
7216f3c3dcfSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
7226f3c3dcfSJed Brown       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
7236f3c3dcfSJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
7246f3c3dcfSJed Brown       ierr = DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
7256f3c3dcfSJed Brown       ierr = DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
7266f3c3dcfSJed Brown       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
7276f3c3dcfSJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
7286f3c3dcfSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
7296f3c3dcfSJed Brown     }
7306f3c3dcfSJed Brown   }
7316f3c3dcfSJed Brown   PetscFunctionReturn(0);
7326f3c3dcfSJed Brown }
7336f3c3dcfSJed Brown 
7346f3c3dcfSJed Brown #undef __FUNCT__
73547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddDM"
73647c6ae99SBarry Smith /*@C
737aa219208SBarry Smith     DMCompositeAddDM - adds a DM  vector to a DMComposite
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith     Collective on DMComposite
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith     Input Parameter:
74247c6ae99SBarry Smith +    dm - the packer object
74347c6ae99SBarry Smith -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
74447c6ae99SBarry Smith 
74547c6ae99SBarry Smith     Level: advanced
74647c6ae99SBarry Smith 
7470c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
7486eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
74947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
75047c6ae99SBarry Smith 
75147c6ae99SBarry Smith @*/
7527087cfbeSBarry Smith PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
75347c6ae99SBarry Smith {
75447c6ae99SBarry Smith   PetscErrorCode         ierr;
75506ebdd98SJed Brown   PetscInt               n,nlocal;
75647c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
75706ebdd98SJed Brown   Vec                    global,local;
75847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmc->data;
75947c6ae99SBarry Smith 
76047c6ae99SBarry Smith   PetscFunctionBegin;
76147c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
76247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
76347c6ae99SBarry Smith   next = com->next;
764ce94432eSBarry Smith   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
76547c6ae99SBarry Smith 
76647c6ae99SBarry Smith   /* create new link */
767b00a9115SJed Brown   ierr = PetscNew(&mine);CHKERRQ(ierr);
76847c6ae99SBarry Smith   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
76947c6ae99SBarry Smith   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
77047c6ae99SBarry Smith   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
77147c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
77206ebdd98SJed Brown   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
77306ebdd98SJed Brown   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
77406ebdd98SJed Brown   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
7758865f1eaSKarl Rupp 
77647c6ae99SBarry Smith   mine->n      = n;
77706ebdd98SJed Brown   mine->nlocal = nlocal;
77847c6ae99SBarry Smith   mine->dm     = dm;
7790298fd71SBarry Smith   mine->next   = NULL;
78047c6ae99SBarry Smith   com->n      += n;
781*7ac2b803SAlex Fikl   com->nghost += nlocal;
78247c6ae99SBarry Smith 
78347c6ae99SBarry Smith   /* add to end of list */
7848865f1eaSKarl Rupp   if (!next) com->next = mine;
7858865f1eaSKarl Rupp   else {
78647c6ae99SBarry Smith     while (next->next) next = next->next;
78747c6ae99SBarry Smith     next->next = mine;
78847c6ae99SBarry Smith   }
78947c6ae99SBarry Smith   com->nDM++;
79047c6ae99SBarry Smith   com->nmine++;
79147c6ae99SBarry Smith   PetscFunctionReturn(0);
79247c6ae99SBarry Smith }
79347c6ae99SBarry Smith 
7949804daf3SBarry Smith #include <petscdraw.h>
79526887b52SJed Brown PETSC_EXTERN PetscErrorCode  VecView_MPI(Vec,PetscViewer);
79647c6ae99SBarry Smith #undef __FUNCT__
79747c6ae99SBarry Smith #define __FUNCT__ "VecView_DMComposite"
7987087cfbeSBarry Smith PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
79947c6ae99SBarry Smith {
80047c6ae99SBarry Smith   DM                     dm;
80147c6ae99SBarry Smith   PetscErrorCode         ierr;
80247c6ae99SBarry Smith   struct DMCompositeLink *next;
80347c6ae99SBarry Smith   PetscBool              isdraw;
804cef07954SSatish Balay   DM_Composite           *com;
80547c6ae99SBarry Smith 
80647c6ae99SBarry Smith   PetscFunctionBegin;
807c688c046SMatthew G Knepley   ierr = VecGetDM(gvec, &dm);CHKERRQ(ierr);
808ce94432eSBarry Smith   if (!dm) SETERRQ(PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
80947c6ae99SBarry Smith   com  = (DM_Composite*)dm->data;
81047c6ae99SBarry Smith   next = com->next;
81147c6ae99SBarry Smith 
812251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
81347c6ae99SBarry Smith   if (!isdraw) {
81447c6ae99SBarry Smith     /* do I really want to call this? */
81547c6ae99SBarry Smith     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
81647c6ae99SBarry Smith   } else {
81747c6ae99SBarry Smith     PetscInt cnt = 0;
81847c6ae99SBarry Smith 
81947c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
82047c6ae99SBarry Smith     while (next) {
82147c6ae99SBarry Smith       Vec         vec;
8229ae5db72SJed Brown       PetscScalar *array;
82347c6ae99SBarry Smith       PetscInt    bs;
82447c6ae99SBarry Smith 
8259ae5db72SJed Brown       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
8269ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&vec);CHKERRQ(ierr);
8279ae5db72SJed Brown       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
8289ae5db72SJed Brown       ierr = VecPlaceArray(vec,array+next->rstart);CHKERRQ(ierr);
8299ae5db72SJed Brown       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
83047c6ae99SBarry Smith       ierr = VecView(vec,viewer);CHKERRQ(ierr);
83147c6ae99SBarry Smith       ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
8329ae5db72SJed Brown       ierr = VecResetArray(vec);CHKERRQ(ierr);
8339ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&vec);CHKERRQ(ierr);
83447c6ae99SBarry Smith       ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
83547c6ae99SBarry Smith       cnt += bs;
83647c6ae99SBarry Smith       next = next->next;
83747c6ae99SBarry Smith     }
83847c6ae99SBarry Smith     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
83947c6ae99SBarry Smith   }
84047c6ae99SBarry Smith   PetscFunctionReturn(0);
84147c6ae99SBarry Smith }
84247c6ae99SBarry Smith 
84347c6ae99SBarry Smith #undef __FUNCT__
8440c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Composite"
8457087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
84647c6ae99SBarry Smith {
84747c6ae99SBarry Smith   PetscErrorCode ierr;
84847c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
84947c6ae99SBarry Smith 
85047c6ae99SBarry Smith   PetscFunctionBegin;
85147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
852d7bf68aeSBarry Smith   ierr = DMSetUp(dm);CHKERRQ(ierr);
853ce94432eSBarry Smith   ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),com->n,com->N,gvec);CHKERRQ(ierr);
854c688c046SMatthew G Knepley   ierr = VecSetDM(*gvec, dm);CHKERRQ(ierr);
85547c6ae99SBarry Smith   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);CHKERRQ(ierr);
85647c6ae99SBarry Smith   PetscFunctionReturn(0);
85747c6ae99SBarry Smith }
85847c6ae99SBarry Smith 
85947c6ae99SBarry Smith #undef __FUNCT__
8600c010503SBarry Smith #define __FUNCT__ "DMCreateLocalVector_Composite"
8617087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
86247c6ae99SBarry Smith {
86347c6ae99SBarry Smith   PetscErrorCode ierr;
86447c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
86547c6ae99SBarry Smith 
86647c6ae99SBarry Smith   PetscFunctionBegin;
86747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
86847c6ae99SBarry Smith   if (!com->setup) {
869d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
87047c6ae99SBarry Smith   }
871f0e01b1fSVincent Le Chenadec   ierr = VecCreateSeq(PETSC_COMM_SELF,com->nghost,lvec);CHKERRQ(ierr);
872c688c046SMatthew G Knepley   ierr = VecSetDM(*lvec, dm);CHKERRQ(ierr);
87347c6ae99SBarry Smith   PetscFunctionReturn(0);
87447c6ae99SBarry Smith }
87547c6ae99SBarry Smith 
87647c6ae99SBarry Smith #undef __FUNCT__
8776eb61c8cSJed Brown #define __FUNCT__ "DMCompositeGetISLocalToGlobalMappings"
87847c6ae99SBarry Smith /*@C
8799ae5db72SJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
88047c6ae99SBarry Smith 
88106ebdd98SJed Brown     Collective on DM
88247c6ae99SBarry Smith 
88347c6ae99SBarry Smith     Input Parameter:
88447c6ae99SBarry Smith .    dm - the packer object
88547c6ae99SBarry Smith 
88647c6ae99SBarry Smith     Output Parameters:
8879ae5db72SJed Brown .    ltogs - the individual mappings for each packed vector. Note that this includes
8889ae5db72SJed Brown            all the ghost points that individual ghosted DMDA's may have.
88947c6ae99SBarry Smith 
89047c6ae99SBarry Smith     Level: advanced
89147c6ae99SBarry Smith 
89247c6ae99SBarry Smith     Notes:
8936eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
89447c6ae99SBarry Smith 
8959ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
89647c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
89747c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
89847c6ae99SBarry Smith 
89947c6ae99SBarry Smith @*/
9007087cfbeSBarry Smith PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
90147c6ae99SBarry Smith {
90247c6ae99SBarry Smith   PetscErrorCode         ierr;
90347c6ae99SBarry Smith   PetscInt               i,*idx,n,cnt;
90447c6ae99SBarry Smith   struct DMCompositeLink *next;
90547c6ae99SBarry Smith   PetscMPIInt            rank;
90647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
90747c6ae99SBarry Smith 
90847c6ae99SBarry Smith   PetscFunctionBegin;
90947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
910728e99d6SJed Brown   ierr = DMSetUp(dm);CHKERRQ(ierr);
911854ce69bSBarry Smith   ierr = PetscMalloc1(com->nDM,ltogs);CHKERRQ(ierr);
91247c6ae99SBarry Smith   next = com->next;
913ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
91447c6ae99SBarry Smith 
91547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
91647c6ae99SBarry Smith   cnt = 0;
91747c6ae99SBarry Smith   while (next) {
9186eb61c8cSJed Brown     ISLocalToGlobalMapping ltog;
9196eb61c8cSJed Brown     PetscMPIInt            size;
92086994e45SJed Brown     const PetscInt         *suboff,*indices;
9216eb61c8cSJed Brown     Vec                    global;
92247c6ae99SBarry Smith 
9236eb61c8cSJed Brown     /* Get sub-DM global indices for each local dof */
9241411c6eeSJed Brown     ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
9256eb61c8cSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
92686994e45SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
927785e854fSJed Brown     ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr);
92847c6ae99SBarry Smith 
9296eb61c8cSJed Brown     /* Get the offsets for the sub-DM global vector */
9306eb61c8cSJed Brown     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
9316eb61c8cSJed Brown     ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
932ce94432eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)global),&size);CHKERRQ(ierr);
9336eb61c8cSJed Brown 
9346eb61c8cSJed Brown     /* Shift the sub-DM definition of the global space to the composite global space */
9356eb61c8cSJed Brown     for (i=0; i<n; i++) {
93686994e45SJed Brown       PetscInt subi = indices[i],lo = 0,hi = size,t;
9376eb61c8cSJed Brown       /* Binary search to find which rank owns subi */
9386eb61c8cSJed Brown       while (hi-lo > 1) {
9396eb61c8cSJed Brown         t = lo + (hi-lo)/2;
9406eb61c8cSJed Brown         if (suboff[t] > subi) hi = t;
9416eb61c8cSJed Brown         else                  lo = t;
9426eb61c8cSJed Brown       }
9436eb61c8cSJed Brown       idx[i] = subi - suboff[lo] + next->grstarts[lo];
9446eb61c8cSJed Brown     }
94586994e45SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
946f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
9476eb61c8cSJed Brown     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
94847c6ae99SBarry Smith     next = next->next;
94947c6ae99SBarry Smith     cnt++;
95047c6ae99SBarry Smith   }
95147c6ae99SBarry Smith   PetscFunctionReturn(0);
95247c6ae99SBarry Smith }
95347c6ae99SBarry Smith 
95447c6ae99SBarry Smith #undef __FUNCT__
95587c85e80SJed Brown #define __FUNCT__ "DMCompositeGetLocalISs"
95687c85e80SJed Brown /*@C
9579ae5db72SJed Brown    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
95887c85e80SJed Brown 
95987c85e80SJed Brown    Not Collective
96087c85e80SJed Brown 
96187c85e80SJed Brown    Input Arguments:
96287c85e80SJed Brown . dm - composite DM
96387c85e80SJed Brown 
96487c85e80SJed Brown    Output Arguments:
96587c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
96687c85e80SJed Brown 
96787c85e80SJed Brown    Level: intermediate
96887c85e80SJed Brown 
96987c85e80SJed Brown    Notes:
97087c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
97187c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
9729ae5db72SJed Brown    scatter to a composite local vector.  The user should not typically need to know which is being done.
97387c85e80SJed Brown 
97487c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
97587c85e80SJed Brown 
97687c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
97787c85e80SJed Brown 
97887c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
97987c85e80SJed Brown 
98087c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
98187c85e80SJed Brown @*/
9827087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
98387c85e80SJed Brown {
98487c85e80SJed Brown   PetscErrorCode         ierr;
98587c85e80SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
98687c85e80SJed Brown   struct DMCompositeLink *link;
98787c85e80SJed Brown   PetscInt               cnt,start;
98887c85e80SJed Brown 
98987c85e80SJed Brown   PetscFunctionBegin;
99087c85e80SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
99187c85e80SJed Brown   PetscValidPointer(is,2);
992785e854fSJed Brown   ierr = PetscMalloc1(com->nmine,is);CHKERRQ(ierr);
99306ebdd98SJed Brown   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
994520db06cSJed Brown     PetscInt bs;
9959ae5db72SJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
9961411c6eeSJed Brown     ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
997520db06cSJed Brown     ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
998520db06cSJed Brown   }
99987c85e80SJed Brown   PetscFunctionReturn(0);
100087c85e80SJed Brown }
100187c85e80SJed Brown 
100287c85e80SJed Brown #undef __FUNCT__
100347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetGlobalISs"
100447c6ae99SBarry Smith /*@C
100547c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
100647c6ae99SBarry Smith 
100747c6ae99SBarry Smith     Collective on DMComposite
100847c6ae99SBarry Smith 
100947c6ae99SBarry Smith     Input Parameter:
101047c6ae99SBarry Smith .    dm - the packer object
101147c6ae99SBarry Smith 
101247c6ae99SBarry Smith     Output Parameters:
101347c6ae99SBarry Smith .    is - the array of index sets
101447c6ae99SBarry Smith 
101547c6ae99SBarry Smith     Level: advanced
101647c6ae99SBarry Smith 
101747c6ae99SBarry Smith     Notes:
101847c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
101947c6ae99SBarry Smith 
102047c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
102147c6ae99SBarry Smith 
10226eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
10236eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
10246eb61c8cSJed Brown        indices.
102547c6ae99SBarry Smith 
1026f3cb0f7eSJed Brown     Fortran Notes:
1027f3cb0f7eSJed Brown 
1028f3cb0f7eSJed Brown        The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM().
1029f3cb0f7eSJed Brown 
10309ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
103147c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
103247c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
103347c6ae99SBarry Smith 
103447c6ae99SBarry Smith @*/
10357087cfbeSBarry Smith PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
103647c6ae99SBarry Smith {
103747c6ae99SBarry Smith   PetscErrorCode         ierr;
103866bb578eSMark Adams   PetscInt               cnt = 0;
103947c6ae99SBarry Smith   struct DMCompositeLink *next;
104047c6ae99SBarry Smith   PetscMPIInt            rank;
104147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
104247c6ae99SBarry Smith 
104347c6ae99SBarry Smith   PetscFunctionBegin;
104447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1045854ce69bSBarry Smith   ierr = PetscMalloc1(com->nDM,is);CHKERRQ(ierr);
104647c6ae99SBarry Smith   next = com->next;
1047ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
104847c6ae99SBarry Smith 
104947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
105047c6ae99SBarry Smith   while (next) {
105166bb578eSMark Adams     ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt]);CHKERRQ(ierr);
10520f21e855SMatthew G. Knepley     if (dm->prob) {
105365c226d8SMatthew G. Knepley       MatNullSpace space;
105465c226d8SMatthew G. Knepley       Mat          pmat;
10550f21e855SMatthew G. Knepley       PetscObject  disc;
10560f21e855SMatthew G. Knepley       PetscInt     Nf;
105765c226d8SMatthew G. Knepley 
10582764a2aaSMatthew G. Knepley       ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
1059f24dd8d2SMatthew G. Knepley       if (cnt < Nf) {
10602764a2aaSMatthew G. Knepley         ierr = PetscDSGetDiscretization(dm->prob, cnt, &disc);CHKERRQ(ierr);
10610f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nullspace", (PetscObject*) &space);CHKERRQ(ierr);
1062aac2dd2dSMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space);CHKERRQ(ierr);}
10630f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space);CHKERRQ(ierr);
1064aac2dd2dSMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space);CHKERRQ(ierr);}
10650f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
1066aac2dd2dSMatthew G. Knepley         if (pmat) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat);CHKERRQ(ierr);}
106765c226d8SMatthew G. Knepley       }
1068f24dd8d2SMatthew G. Knepley     }
106947c6ae99SBarry Smith     cnt++;
107047c6ae99SBarry Smith     next = next->next;
107147c6ae99SBarry Smith   }
107247c6ae99SBarry Smith   PetscFunctionReturn(0);
107347c6ae99SBarry Smith }
107447c6ae99SBarry Smith 
10754d343eeaSMatthew G Knepley #undef __FUNCT__
10764d343eeaSMatthew G Knepley #define __FUNCT__ "DMCreateFieldIS_Composite"
107721c9b008SJed Brown PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields)
10784d343eeaSMatthew G Knepley {
10794d343eeaSMatthew G Knepley   PetscInt       nDM;
10804d343eeaSMatthew G Knepley   DM             *dms;
10814d343eeaSMatthew G Knepley   PetscInt       i;
10824d343eeaSMatthew G Knepley   PetscErrorCode ierr;
10834d343eeaSMatthew G Knepley 
10844d343eeaSMatthew G Knepley   PetscFunctionBegin;
10854d343eeaSMatthew G Knepley   ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
10868865f1eaSKarl Rupp   if (numFields) *numFields = nDM;
10874d343eeaSMatthew G Knepley   ierr = DMCompositeGetGlobalISs(dm, fields);CHKERRQ(ierr);
10884d343eeaSMatthew G Knepley   if (fieldNames) {
1089785e854fSJed Brown     ierr = PetscMalloc1(nDM, &dms);CHKERRQ(ierr);
1090785e854fSJed Brown     ierr = PetscMalloc1(nDM, fieldNames);CHKERRQ(ierr);
10914d343eeaSMatthew G Knepley     ierr = DMCompositeGetEntriesArray(dm, dms);CHKERRQ(ierr);
10924d343eeaSMatthew G Knepley     for (i=0; i<nDM; i++) {
10934d343eeaSMatthew G Knepley       char       buf[256];
10944d343eeaSMatthew G Knepley       const char *splitname;
10954d343eeaSMatthew G Knepley 
10964d343eeaSMatthew G Knepley       /* Split naming precedence: object name, prefix, number */
10974d343eeaSMatthew G Knepley       splitname = ((PetscObject) dm)->name;
10984d343eeaSMatthew G Knepley       if (!splitname) {
10994d343eeaSMatthew G Knepley         ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr);
11004d343eeaSMatthew G Knepley         if (splitname) {
11014d343eeaSMatthew G Knepley           size_t len;
11028caf3d72SBarry Smith           ierr                 = PetscStrncpy(buf,splitname,sizeof(buf));CHKERRQ(ierr);
11038caf3d72SBarry Smith           buf[sizeof(buf) - 1] = 0;
11044d343eeaSMatthew G Knepley           ierr                 = PetscStrlen(buf,&len);CHKERRQ(ierr);
11054d343eeaSMatthew G Knepley           if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
11064d343eeaSMatthew G Knepley           splitname = buf;
11074d343eeaSMatthew G Knepley         }
11084d343eeaSMatthew G Knepley       }
11094d343eeaSMatthew G Knepley       if (!splitname) {
11108caf3d72SBarry Smith         ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",i);CHKERRQ(ierr);
11114d343eeaSMatthew G Knepley         splitname = buf;
11124d343eeaSMatthew G Knepley       }
111321c9b008SJed Brown       ierr = PetscStrallocpy(splitname,&(*fieldNames)[i]);CHKERRQ(ierr);
11144d343eeaSMatthew G Knepley     }
11154d343eeaSMatthew G Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
11164d343eeaSMatthew G Knepley   }
11174d343eeaSMatthew G Knepley   PetscFunctionReturn(0);
11184d343eeaSMatthew G Knepley }
11194d343eeaSMatthew G Knepley 
1120e7c4fc90SDmitry Karpeev /*
1121e7c4fc90SDmitry Karpeev  This could take over from DMCreateFieldIS(), as it is more general,
11220298fd71SBarry Smith  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1123e7c4fc90SDmitry Karpeev  At this point it's probably best to be less intrusive, however.
1124e7c4fc90SDmitry Karpeev  */
1125e7c4fc90SDmitry Karpeev #undef __FUNCT__
112616621825SDmitry Karpeev #define __FUNCT__ "DMCreateFieldDecomposition_Composite"
112716621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
1128e7c4fc90SDmitry Karpeev {
1129e7c4fc90SDmitry Karpeev   PetscInt       nDM;
1130e7c4fc90SDmitry Karpeev   PetscInt       i;
1131e7c4fc90SDmitry Karpeev   PetscErrorCode ierr;
1132e7c4fc90SDmitry Karpeev 
1133e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
1134e7c4fc90SDmitry Karpeev   ierr = DMCreateFieldIS_Composite(dm, len, namelist, islist);CHKERRQ(ierr);
1135e7c4fc90SDmitry Karpeev   if (dmlist) {
1136e7c4fc90SDmitry Karpeev     ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
1137785e854fSJed Brown     ierr = PetscMalloc1(nDM, dmlist);CHKERRQ(ierr);
1138e7c4fc90SDmitry Karpeev     ierr = DMCompositeGetEntriesArray(dm, *dmlist);CHKERRQ(ierr);
1139e7c4fc90SDmitry Karpeev     for (i=0; i<nDM; i++) {
1140e7c4fc90SDmitry Karpeev       ierr = PetscObjectReference((PetscObject)((*dmlist)[i]));CHKERRQ(ierr);
1141e7c4fc90SDmitry Karpeev     }
1142e7c4fc90SDmitry Karpeev   }
1143e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
1144e7c4fc90SDmitry Karpeev }
1145e7c4fc90SDmitry Karpeev 
1146e7c4fc90SDmitry Karpeev 
1147e7c4fc90SDmitry Karpeev 
114847c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
114947c6ae99SBarry Smith #undef __FUNCT__
115047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors"
115147c6ae99SBarry Smith /*@C
11529ae5db72SJed Brown     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
115347c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
115447c6ae99SBarry Smith 
115547c6ae99SBarry Smith     Not Collective
115647c6ae99SBarry Smith 
115747c6ae99SBarry Smith     Input Parameter:
115847c6ae99SBarry Smith .    dm - the packer object
115947c6ae99SBarry Smith 
116047c6ae99SBarry Smith     Output Parameter:
11619ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
116247c6ae99SBarry Smith 
116347c6ae99SBarry Smith     Level: advanced
116447c6ae99SBarry Smith 
11659ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
11666eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
116747c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
116847c6ae99SBarry Smith 
116947c6ae99SBarry Smith @*/
11707087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
117147c6ae99SBarry Smith {
117247c6ae99SBarry Smith   va_list                Argp;
117347c6ae99SBarry Smith   PetscErrorCode         ierr;
117447c6ae99SBarry Smith   struct DMCompositeLink *next;
117547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
117647c6ae99SBarry Smith 
117747c6ae99SBarry Smith   PetscFunctionBegin;
117847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
117947c6ae99SBarry Smith   next = com->next;
118047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
118147c6ae99SBarry Smith   va_start(Argp,dm);
118247c6ae99SBarry Smith   while (next) {
118347c6ae99SBarry Smith     Vec *vec;
118447c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
118506930112SJed Brown     if (vec) {ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);}
118647c6ae99SBarry Smith     next = next->next;
118747c6ae99SBarry Smith   }
118847c6ae99SBarry Smith   va_end(Argp);
118947c6ae99SBarry Smith   PetscFunctionReturn(0);
119047c6ae99SBarry Smith }
119147c6ae99SBarry Smith 
119247c6ae99SBarry Smith #undef __FUNCT__
119347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors"
119447c6ae99SBarry Smith /*@C
11959ae5db72SJed Brown     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
119647c6ae99SBarry Smith 
119747c6ae99SBarry Smith     Not Collective
119847c6ae99SBarry Smith 
119947c6ae99SBarry Smith     Input Parameter:
120047c6ae99SBarry Smith .    dm - the packer object
120147c6ae99SBarry Smith 
120247c6ae99SBarry Smith     Output Parameter:
12039ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
120447c6ae99SBarry Smith 
120547c6ae99SBarry Smith     Level: advanced
120647c6ae99SBarry Smith 
12079ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
12086eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
120947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
121047c6ae99SBarry Smith 
121147c6ae99SBarry Smith @*/
12127087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
121347c6ae99SBarry Smith {
121447c6ae99SBarry Smith   va_list                Argp;
121547c6ae99SBarry Smith   PetscErrorCode         ierr;
121647c6ae99SBarry Smith   struct DMCompositeLink *next;
121747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
121847c6ae99SBarry Smith 
121947c6ae99SBarry Smith   PetscFunctionBegin;
122047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
122147c6ae99SBarry Smith   next = com->next;
122247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
122347c6ae99SBarry Smith   va_start(Argp,dm);
122447c6ae99SBarry Smith   while (next) {
122547c6ae99SBarry Smith     Vec *vec;
122647c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
122706930112SJed Brown     if (vec) {ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);}
122847c6ae99SBarry Smith     next = next->next;
122947c6ae99SBarry Smith   }
123047c6ae99SBarry Smith   va_end(Argp);
123147c6ae99SBarry Smith   PetscFunctionReturn(0);
123247c6ae99SBarry Smith }
123347c6ae99SBarry Smith 
123447c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
123547c6ae99SBarry Smith #undef __FUNCT__
123647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries"
123747c6ae99SBarry Smith /*@C
12389ae5db72SJed Brown     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
123947c6ae99SBarry Smith 
124047c6ae99SBarry Smith     Not Collective
124147c6ae99SBarry Smith 
124247c6ae99SBarry Smith     Input Parameter:
124347c6ae99SBarry Smith .    dm - the packer object
124447c6ae99SBarry Smith 
124547c6ae99SBarry Smith     Output Parameter:
12469ae5db72SJed Brown .   DM ... - the individual entries (DMs)
124747c6ae99SBarry Smith 
124847c6ae99SBarry Smith     Level: advanced
124947c6ae99SBarry Smith 
12502fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
12516eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
125247c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
125347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
125447c6ae99SBarry Smith 
125547c6ae99SBarry Smith @*/
12567087cfbeSBarry Smith PetscErrorCode  DMCompositeGetEntries(DM dm,...)
125747c6ae99SBarry Smith {
125847c6ae99SBarry Smith   va_list                Argp;
125947c6ae99SBarry Smith   struct DMCompositeLink *next;
126047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
126147c6ae99SBarry Smith 
126247c6ae99SBarry Smith   PetscFunctionBegin;
126347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
126447c6ae99SBarry Smith   next = com->next;
126547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
126647c6ae99SBarry Smith   va_start(Argp,dm);
126747c6ae99SBarry Smith   while (next) {
126847c6ae99SBarry Smith     DM *dmn;
126947c6ae99SBarry Smith     dmn = va_arg(Argp, DM*);
12709ae5db72SJed Brown     if (dmn) *dmn = next->dm;
127147c6ae99SBarry Smith     next = next->next;
127247c6ae99SBarry Smith   }
127347c6ae99SBarry Smith   va_end(Argp);
127447c6ae99SBarry Smith   PetscFunctionReturn(0);
127547c6ae99SBarry Smith }
127647c6ae99SBarry Smith 
127747c6ae99SBarry Smith #undef __FUNCT__
12782fa5ba8aSJed Brown #define __FUNCT__ "DMCompositeGetEntriesArray"
1279dbab29e1SMark F. Adams /*@C
12802fa5ba8aSJed Brown     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
12812fa5ba8aSJed Brown 
12822fa5ba8aSJed Brown     Not Collective
12832fa5ba8aSJed Brown 
12842fa5ba8aSJed Brown     Input Parameter:
1285907376e6SBarry Smith .    dm - the packer object
1286907376e6SBarry Smith 
1287907376e6SBarry Smith     Output Parameter:
1288907376e6SBarry Smith .    dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs
12892fa5ba8aSJed Brown 
12902fa5ba8aSJed Brown     Level: advanced
12912fa5ba8aSJed Brown 
12922fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
12932fa5ba8aSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
12942fa5ba8aSJed Brown          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
12952fa5ba8aSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
12962fa5ba8aSJed Brown 
12972fa5ba8aSJed Brown @*/
12982fa5ba8aSJed Brown PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
12992fa5ba8aSJed Brown {
13002fa5ba8aSJed Brown   struct DMCompositeLink *next;
13012fa5ba8aSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
13022fa5ba8aSJed Brown   PetscInt               i;
13032fa5ba8aSJed Brown 
13042fa5ba8aSJed Brown   PetscFunctionBegin;
13052fa5ba8aSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
13062fa5ba8aSJed Brown   /* loop over packed objects, handling one at at time */
13072fa5ba8aSJed Brown   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
13082fa5ba8aSJed Brown   PetscFunctionReturn(0);
13092fa5ba8aSJed Brown }
13102fa5ba8aSJed Brown 
13112fa5ba8aSJed Brown #undef __FUNCT__
13120c010503SBarry Smith #define __FUNCT__ "DMRefine_Composite"
13137087cfbeSBarry Smith PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
131447c6ae99SBarry Smith {
131547c6ae99SBarry Smith   PetscErrorCode         ierr;
131647c6ae99SBarry Smith   struct DMCompositeLink *next;
131747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmi->data;
131847c6ae99SBarry Smith   DM                     dm;
131947c6ae99SBarry Smith 
132047c6ae99SBarry Smith   PetscFunctionBegin;
132147c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1322ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
1323ce94432eSBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
1324ce94432eSBarry Smith   }
13252ce3a92bSJed Brown   ierr = DMSetUp(dmi);CHKERRQ(ierr);
132647c6ae99SBarry Smith   next = com->next;
132747c6ae99SBarry Smith   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
132847c6ae99SBarry Smith 
132947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
133047c6ae99SBarry Smith   while (next) {
133147c6ae99SBarry Smith     ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
133247c6ae99SBarry Smith     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
133347c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
133447c6ae99SBarry Smith     next = next->next;
133547c6ae99SBarry Smith   }
133647c6ae99SBarry Smith   PetscFunctionReturn(0);
133747c6ae99SBarry Smith }
133847c6ae99SBarry Smith 
133914354c39SJed Brown #undef __FUNCT__
134014354c39SJed Brown #define __FUNCT__ "DMCoarsen_Composite"
134114354c39SJed Brown PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
134214354c39SJed Brown {
134314354c39SJed Brown   PetscErrorCode         ierr;
134414354c39SJed Brown   struct DMCompositeLink *next;
134514354c39SJed Brown   DM_Composite           *com = (DM_Composite*)dmi->data;
134614354c39SJed Brown   DM                     dm;
134714354c39SJed Brown 
134814354c39SJed Brown   PetscFunctionBegin;
134914354c39SJed Brown   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
13502ce3a92bSJed Brown   ierr = DMSetUp(dmi);CHKERRQ(ierr);
13512ee06e3bSJed Brown   if (comm == MPI_COMM_NULL) {
135225296bd5SBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
135325296bd5SBarry Smith   }
135414354c39SJed Brown   next = com->next;
135514354c39SJed Brown   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
135614354c39SJed Brown 
135714354c39SJed Brown   /* loop over packed objects, handling one at at time */
135814354c39SJed Brown   while (next) {
135914354c39SJed Brown     ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr);
136014354c39SJed Brown     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
136114354c39SJed Brown     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
136214354c39SJed Brown     next = next->next;
136314354c39SJed Brown   }
136414354c39SJed Brown   PetscFunctionReturn(0);
136514354c39SJed Brown }
136647c6ae99SBarry Smith 
136747c6ae99SBarry Smith #undef __FUNCT__
1368e727c939SJed Brown #define __FUNCT__ "DMCreateInterpolation_Composite"
1369e727c939SJed Brown PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
137047c6ae99SBarry Smith {
137147c6ae99SBarry Smith   PetscErrorCode         ierr;
13729ae5db72SJed Brown   PetscInt               m,n,M,N,nDM,i;
137347c6ae99SBarry Smith   struct DMCompositeLink *nextc;
137447c6ae99SBarry Smith   struct DMCompositeLink *nextf;
137525296bd5SBarry Smith   Vec                    gcoarse,gfine,*vecs;
137647c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
137747c6ae99SBarry Smith   DM_Composite           *comfine   = (DM_Composite*)fine->data;
13789ae5db72SJed Brown   Mat                    *mats;
137947c6ae99SBarry Smith 
138047c6ae99SBarry Smith   PetscFunctionBegin;
138147c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
138247c6ae99SBarry Smith   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
1383f692024eSJed Brown   ierr = DMSetUp(coarse);CHKERRQ(ierr);
1384f692024eSJed Brown   ierr = DMSetUp(fine);CHKERRQ(ierr);
138547c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
13869ae5db72SJed Brown   ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
13879ae5db72SJed Brown   ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr);
138847c6ae99SBarry Smith   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
138947c6ae99SBarry Smith   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
139047c6ae99SBarry Smith   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
139147c6ae99SBarry Smith   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
13929ae5db72SJed Brown   ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
13939ae5db72SJed Brown   ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr);
139447c6ae99SBarry Smith 
13959ae5db72SJed Brown   nDM = comfine->nDM;
1396ce94432eSBarry Smith   if (nDM != comcoarse->nDM) SETERRQ2(PetscObjectComm((PetscObject)fine),PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %D entries, but coarse has %D",nDM,comcoarse->nDM);
13971795a4d1SJed Brown   ierr = PetscCalloc1(nDM*nDM,&mats);CHKERRQ(ierr);
139825296bd5SBarry Smith   if (v) {
13991795a4d1SJed Brown     ierr = PetscCalloc1(nDM,&vecs);CHKERRQ(ierr);
140025296bd5SBarry Smith   }
140147c6ae99SBarry Smith 
140247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
14039ae5db72SJed Brown   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
140425296bd5SBarry Smith     if (!v) {
14050298fd71SBarry Smith       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL);CHKERRQ(ierr);
140625296bd5SBarry Smith     } else {
140725296bd5SBarry Smith       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr);
140825296bd5SBarry Smith     }
140947c6ae99SBarry Smith   }
1410ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A);CHKERRQ(ierr);
141125296bd5SBarry Smith   if (v) {
1412ce94432eSBarry Smith     ierr = VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v);CHKERRQ(ierr);
141325296bd5SBarry Smith   }
14149ae5db72SJed Brown   for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);}
14159ae5db72SJed Brown   ierr = PetscFree(mats);CHKERRQ(ierr);
141625296bd5SBarry Smith   if (v) {
141725296bd5SBarry Smith     for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);}
141825296bd5SBarry Smith     ierr = PetscFree(vecs);CHKERRQ(ierr);
141925296bd5SBarry Smith   }
142047c6ae99SBarry Smith   PetscFunctionReturn(0);
142147c6ae99SBarry Smith }
142247c6ae99SBarry Smith 
142347c6ae99SBarry Smith #undef __FUNCT__
1424184d77edSJed Brown #define __FUNCT__ "DMGetLocalToGlobalMapping_Composite"
1425184d77edSJed Brown static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
14261411c6eeSJed Brown {
14271411c6eeSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
14281411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1429f7efa3c7SJed Brown   PetscInt               i;
14301411c6eeSJed Brown   PetscErrorCode         ierr;
14311411c6eeSJed Brown 
14321411c6eeSJed Brown   PetscFunctionBegin;
14331411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
14341411c6eeSJed Brown   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
1435ce94432eSBarry Smith   ierr = ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr);
14369ae5db72SJed Brown   for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
14371411c6eeSJed Brown   ierr = PetscFree(ltogs);CHKERRQ(ierr);
14381411c6eeSJed Brown   PetscFunctionReturn(0);
14391411c6eeSJed Brown }
14401411c6eeSJed Brown 
14411411c6eeSJed Brown 
14421411c6eeSJed Brown #undef __FUNCT__
1443e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_Composite"
1444b412c318SBarry Smith PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring)
144547c6ae99SBarry Smith {
144647c6ae99SBarry Smith   PetscErrorCode  ierr;
144747c6ae99SBarry Smith   PetscInt        n,i,cnt;
144847c6ae99SBarry Smith   ISColoringValue *colors;
144947c6ae99SBarry Smith   PetscBool       dense  = PETSC_FALSE;
145047c6ae99SBarry Smith   ISColoringValue maxcol = 0;
145147c6ae99SBarry Smith   DM_Composite    *com   = (DM_Composite*)dm->data;
145247c6ae99SBarry Smith 
145347c6ae99SBarry Smith   PetscFunctionBegin;
145447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1455ce94432eSBarry Smith   if (ctype == IS_COLORING_GHOSTED) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported");
1456e3247f34SBarry Smith   else if (ctype == IS_COLORING_GLOBAL) {
145747c6ae99SBarry Smith     n = com->n;
1458ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1459785e854fSJed Brown   ierr = PetscMalloc1(n,&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
146047c6ae99SBarry Smith 
1461c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL);CHKERRQ(ierr);
146247c6ae99SBarry Smith   if (dense) {
146347c6ae99SBarry Smith     for (i=0; i<n; i++) {
146447c6ae99SBarry Smith       colors[i] = (ISColoringValue)(com->rstart + i);
146547c6ae99SBarry Smith     }
146647c6ae99SBarry Smith     maxcol = com->N;
146747c6ae99SBarry Smith   } else {
146847c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
146947c6ae99SBarry Smith     PetscMPIInt            rank;
147047c6ae99SBarry Smith 
1471ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
147247c6ae99SBarry Smith     cnt  = 0;
147347c6ae99SBarry Smith     while (next) {
147447c6ae99SBarry Smith       ISColoring lcoloring;
147547c6ae99SBarry Smith 
1476b412c318SBarry Smith       ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring);CHKERRQ(ierr);
147747c6ae99SBarry Smith       for (i=0; i<lcoloring->N; i++) {
147847c6ae99SBarry Smith         colors[cnt++] = maxcol + lcoloring->colors[i];
147947c6ae99SBarry Smith       }
148047c6ae99SBarry Smith       maxcol += lcoloring->n;
1481fcfd50ebSBarry Smith       ierr    = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
148247c6ae99SBarry Smith       next    = next->next;
148347c6ae99SBarry Smith     }
148447c6ae99SBarry Smith   }
1485aaf3ff59SMatthew G. Knepley   ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,PETSC_OWN_POINTER,coloring);CHKERRQ(ierr);
148647c6ae99SBarry Smith   PetscFunctionReturn(0);
148747c6ae99SBarry Smith }
148847c6ae99SBarry Smith 
148947c6ae99SBarry Smith #undef __FUNCT__
14900c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
14917087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
149247c6ae99SBarry Smith {
149347c6ae99SBarry Smith   PetscErrorCode         ierr;
149447c6ae99SBarry Smith   struct DMCompositeLink *next;
149547c6ae99SBarry Smith   PetscInt               cnt = 3;
149647c6ae99SBarry Smith   PetscMPIInt            rank;
149747c6ae99SBarry Smith   PetscScalar            *garray,*larray;
149847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
149947c6ae99SBarry Smith 
150047c6ae99SBarry Smith   PetscFunctionBegin;
150147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
150247c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
150347c6ae99SBarry Smith   next = com->next;
150447c6ae99SBarry Smith   if (!com->setup) {
1505d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
150647c6ae99SBarry Smith   }
1507ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
150847c6ae99SBarry Smith   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
150947c6ae99SBarry Smith   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
151047c6ae99SBarry Smith 
151147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
151247c6ae99SBarry Smith   while (next) {
151347c6ae99SBarry Smith     Vec      local,global;
151447c6ae99SBarry Smith     PetscInt N;
151547c6ae99SBarry Smith 
151647c6ae99SBarry Smith     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
151747c6ae99SBarry Smith     ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
151847c6ae99SBarry Smith     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
151947c6ae99SBarry Smith     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
152047c6ae99SBarry Smith     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
152147c6ae99SBarry Smith     ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
152247c6ae99SBarry Smith     ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
152347c6ae99SBarry Smith     ierr = VecResetArray(global);CHKERRQ(ierr);
152447c6ae99SBarry Smith     ierr = VecResetArray(local);CHKERRQ(ierr);
152547c6ae99SBarry Smith     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
152647c6ae99SBarry Smith     ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
152747c6ae99SBarry Smith     cnt++;
152806ebdd98SJed Brown     larray += next->nlocal;
152947c6ae99SBarry Smith     next    = next->next;
153047c6ae99SBarry Smith   }
153147c6ae99SBarry Smith 
15320298fd71SBarry Smith   ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr);
15330298fd71SBarry Smith   ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr);
153447c6ae99SBarry Smith   PetscFunctionReturn(0);
153547c6ae99SBarry Smith }
153647c6ae99SBarry Smith 
153747c6ae99SBarry Smith #undef __FUNCT__
15380c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
15397087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
15400c010503SBarry Smith {
15410c010503SBarry Smith   PetscFunctionBegin;
15420c010503SBarry Smith   PetscFunctionReturn(0);
15430c010503SBarry Smith }
154447c6ae99SBarry Smith 
15456ae3a549SBarry Smith /*MC
15466ae3a549SBarry Smith    DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs
15476ae3a549SBarry Smith 
15486ae3a549SBarry Smith   Level: intermediate
15496ae3a549SBarry Smith 
15501abcec8cSBarry Smith .seealso: DMType, DM, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate()
15516ae3a549SBarry Smith M*/
15526ae3a549SBarry Smith 
15536ae3a549SBarry Smith 
1554a4121054SBarry Smith #undef __FUNCT__
1555a4121054SBarry Smith #define __FUNCT__ "DMCreate_Composite"
15568cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1557a4121054SBarry Smith {
1558a4121054SBarry Smith   PetscErrorCode ierr;
1559a4121054SBarry Smith   DM_Composite   *com;
1560a4121054SBarry Smith 
1561a4121054SBarry Smith   PetscFunctionBegin;
1562b00a9115SJed Brown   ierr          = PetscNewLog(p,&com);CHKERRQ(ierr);
1563a4121054SBarry Smith   p->data       = com;
1564a4121054SBarry Smith   ierr          = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1565a4121054SBarry Smith   com->n        = 0;
1566*7ac2b803SAlex Fikl   com->nghost   = 0;
15670298fd71SBarry Smith   com->next     = NULL;
1568a4121054SBarry Smith   com->nDM      = 0;
1569a4121054SBarry Smith 
1570a4121054SBarry Smith   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1571a4121054SBarry Smith   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1572184d77edSJed Brown   p->ops->getlocaltoglobalmapping         = DMGetLocalToGlobalMapping_Composite;
15734d343eeaSMatthew G Knepley   p->ops->createfieldis                   = DMCreateFieldIS_Composite;
157416621825SDmitry Karpeev   p->ops->createfielddecomposition        = DMCreateFieldDecomposition_Composite;
1575a4121054SBarry Smith   p->ops->refine                          = DMRefine_Composite;
157614354c39SJed Brown   p->ops->coarsen                         = DMCoarsen_Composite;
157725296bd5SBarry Smith   p->ops->createinterpolation             = DMCreateInterpolation_Composite;
157825296bd5SBarry Smith   p->ops->creatematrix                    = DMCreateMatrix_Composite;
1579e727c939SJed Brown   p->ops->getcoloring                     = DMCreateColoring_Composite;
1580a4121054SBarry Smith   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1581a4121054SBarry Smith   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1582a4121054SBarry Smith   p->ops->destroy                         = DMDestroy_Composite;
1583a4121054SBarry Smith   p->ops->view                            = DMView_Composite;
1584a4121054SBarry Smith   p->ops->setup                           = DMSetUp_Composite;
1585a4121054SBarry Smith   PetscFunctionReturn(0);
1586a4121054SBarry Smith }
1587a4121054SBarry Smith 
15880c010503SBarry Smith #undef __FUNCT__
15890c010503SBarry Smith #define __FUNCT__ "DMCompositeCreate"
15900c010503SBarry Smith /*@C
15910c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
15920c010503SBarry Smith       vectors made up of several subvectors.
15930c010503SBarry Smith 
15940c010503SBarry Smith     Collective on MPI_Comm
159547c6ae99SBarry Smith 
159647c6ae99SBarry Smith     Input Parameter:
15970c010503SBarry Smith .   comm - the processors that will share the global vector
15980c010503SBarry Smith 
15990c010503SBarry Smith     Output Parameters:
16000c010503SBarry Smith .   packer - the packer object
160147c6ae99SBarry Smith 
160247c6ae99SBarry Smith     Level: advanced
160347c6ae99SBarry Smith 
16041abcec8cSBarry Smith .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(), DMCOMPOSITE,DMCreate()
16056eb61c8cSJed Brown          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
160647c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
160747c6ae99SBarry Smith 
160847c6ae99SBarry Smith @*/
16097087cfbeSBarry Smith PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
161047c6ae99SBarry Smith {
16110c010503SBarry Smith   PetscErrorCode ierr;
16120c010503SBarry Smith 
161347c6ae99SBarry Smith   PetscFunctionBegin;
16140c010503SBarry Smith   PetscValidPointer(packer,2);
1615a4121054SBarry Smith   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1616a4121054SBarry Smith   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
161747c6ae99SBarry Smith   PetscFunctionReturn(0);
161847c6ae99SBarry Smith }
1619