xref: /petsc/src/dm/impls/composite/pack.c (revision 1fd49c2520df80b4b23f7bac1a7aae848a1a6ec1)
147c6ae99SBarry Smith 
2ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h>       /*I  "petscdmcomposite.h"  I*/
3af0996ceSBarry Smith #include <petsc/private/isimpl.h>
4e10fd815SStefano Zampini #include <petsc/private/glvisviewerimpl.h>
52764a2aaSMatthew G. Knepley #include <petscds.h>
647c6ae99SBarry Smith 
747c6ae99SBarry Smith /*@C
847c6ae99SBarry Smith     DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
9bebe2cf6SSatish Balay       separate components (DMs) in a DMto build the correct matrix nonzero structure.
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith 
1247c6ae99SBarry Smith     Logically Collective on MPI_Comm
1347c6ae99SBarry Smith 
1447c6ae99SBarry Smith     Input Parameter:
1547c6ae99SBarry Smith +   dm - the composite object
1647c6ae99SBarry Smith -   formcouplelocations - routine to set the nonzero locations in the matrix
1747c6ae99SBarry Smith 
1847c6ae99SBarry Smith     Level: advanced
1947c6ae99SBarry Smith 
201b2093e4SBarry Smith     Notes: See DMSetApplicationContext() and DMGetApplicationContext() for how to get user information into
2147c6ae99SBarry Smith         this routine
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith @*/
247087cfbeSBarry Smith PetscErrorCode  DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt))
2547c6ae99SBarry Smith {
2647c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith   PetscFunctionBegin;
2947c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
3047c6ae99SBarry Smith   PetscFunctionReturn(0);
3147c6ae99SBarry Smith }
3247c6ae99SBarry Smith 
336bf464f9SBarry Smith PetscErrorCode  DMDestroy_Composite(DM dm)
3447c6ae99SBarry Smith {
3547c6ae99SBarry Smith   PetscErrorCode         ierr;
3647c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
3747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
3847c6ae99SBarry Smith 
3947c6ae99SBarry Smith   PetscFunctionBegin;
4047c6ae99SBarry Smith   next = com->next;
4147c6ae99SBarry Smith   while (next) {
4247c6ae99SBarry Smith     prev = next;
4347c6ae99SBarry Smith     next = next->next;
44fcfd50ebSBarry Smith     ierr = DMDestroy(&prev->dm);CHKERRQ(ierr);
4547c6ae99SBarry Smith     ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
4647c6ae99SBarry Smith     ierr = PetscFree(prev);CHKERRQ(ierr);
4747c6ae99SBarry Smith   }
48e10fd815SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL);CHKERRQ(ierr);
49435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
50435a35e8SMatthew G Knepley   ierr = PetscFree(com);CHKERRQ(ierr);
5147c6ae99SBarry Smith   PetscFunctionReturn(0);
5247c6ae99SBarry Smith }
5347c6ae99SBarry Smith 
547087cfbeSBarry Smith PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
5547c6ae99SBarry Smith {
5647c6ae99SBarry Smith   PetscErrorCode ierr;
5747c6ae99SBarry Smith   PetscBool      iascii;
5847c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
5947c6ae99SBarry Smith 
6047c6ae99SBarry Smith   PetscFunctionBegin;
61251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6247c6ae99SBarry Smith   if (iascii) {
6347c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
6447c6ae99SBarry Smith     PetscInt               i;
6547c6ae99SBarry Smith 
6647c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix");CHKERRQ(ierr);
679ae5db72SJed Brown     ierr = PetscViewerASCIIPrintf(v,"  contains %D DMs\n",com->nDM);CHKERRQ(ierr);
6847c6ae99SBarry Smith     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
6947c6ae99SBarry Smith     for (i=0; lnk; lnk=lnk->next,i++) {
709ae5db72SJed Brown       ierr = PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
7147c6ae99SBarry Smith       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
7247c6ae99SBarry Smith       ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
7347c6ae99SBarry Smith       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
7447c6ae99SBarry Smith     }
7547c6ae99SBarry Smith     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
7647c6ae99SBarry Smith   }
7747c6ae99SBarry Smith   PetscFunctionReturn(0);
7847c6ae99SBarry Smith }
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
817087cfbeSBarry Smith PetscErrorCode  DMSetUp_Composite(DM dm)
8247c6ae99SBarry Smith {
8347c6ae99SBarry Smith   PetscErrorCode         ierr;
8447c6ae99SBarry Smith   PetscInt               nprev = 0;
8547c6ae99SBarry Smith   PetscMPIInt            rank,size;
8647c6ae99SBarry Smith   DM_Composite           *com  = (DM_Composite*)dm->data;
8747c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
8847c6ae99SBarry Smith   PetscLayout            map;
8947c6ae99SBarry Smith 
9047c6ae99SBarry Smith   PetscFunctionBegin;
91ce94432eSBarry Smith   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
92ce94432eSBarry Smith   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm),&map);CHKERRQ(ierr);
9347c6ae99SBarry Smith   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
9447c6ae99SBarry Smith   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
9547c6ae99SBarry Smith   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
9647c6ae99SBarry Smith   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
9747c6ae99SBarry Smith   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
980298fd71SBarry Smith   ierr = PetscLayoutGetRange(map,&com->rstart,NULL);CHKERRQ(ierr);
99fcfd50ebSBarry Smith   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
10047c6ae99SBarry Smith 
1019ae5db72SJed Brown   /* now set the rstart for each linked vector */
102ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
103ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
10447c6ae99SBarry Smith   while (next) {
10547c6ae99SBarry Smith     next->rstart  = nprev;
10606ebdd98SJed Brown     nprev        += next->n;
10747c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
108785e854fSJed Brown     ierr          = PetscMalloc1(size,&next->grstarts);CHKERRQ(ierr);
109ce94432eSBarry Smith     ierr          = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
11047c6ae99SBarry Smith     next          = next->next;
11147c6ae99SBarry Smith   }
11247c6ae99SBarry Smith   com->setup = PETSC_TRUE;
11347c6ae99SBarry Smith   PetscFunctionReturn(0);
11447c6ae99SBarry Smith }
11547c6ae99SBarry Smith 
11647c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
11747c6ae99SBarry Smith 
11873e31fe2SJed Brown /*@
11947c6ae99SBarry Smith     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
12047c6ae99SBarry Smith        representation.
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith     Not Collective
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith     Input Parameter:
12547c6ae99SBarry Smith .    dm - the packer object
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith     Output Parameter:
12847c6ae99SBarry Smith .     nDM - the number of DMs
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith     Level: beginner
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith @*/
1337087cfbeSBarry Smith PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
13447c6ae99SBarry Smith {
13547c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
1365fd66863SKarl Rupp 
13747c6ae99SBarry Smith   PetscFunctionBegin;
13847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
13947c6ae99SBarry Smith   *nDM = com->nDM;
14047c6ae99SBarry Smith   PetscFunctionReturn(0);
14147c6ae99SBarry Smith }
14247c6ae99SBarry Smith 
14347c6ae99SBarry Smith 
14447c6ae99SBarry Smith /*@C
14547c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
14647c6ae99SBarry Smith        representation.
14747c6ae99SBarry Smith 
14847c6ae99SBarry Smith     Collective on DMComposite
14947c6ae99SBarry Smith 
1509ae5db72SJed Brown     Input Parameters:
15147c6ae99SBarry Smith +    dm - the packer object
1529ae5db72SJed Brown -    gvec - the global vector
1539ae5db72SJed Brown 
1549ae5db72SJed Brown     Output Parameters:
1550298fd71SBarry Smith .    Vec* ... - the packed parallel vectors, NULL for those that are not needed
15647c6ae99SBarry Smith 
15747c6ae99SBarry Smith     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
15847c6ae99SBarry Smith 
159f73e5cebSJed Brown     Fortran Notes:
160f73e5cebSJed Brown 
161f73e5cebSJed Brown     Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
162f73e5cebSJed Brown     or use the alternative interface DMCompositeGetAccessArray().
163f73e5cebSJed Brown 
16447c6ae99SBarry Smith     Level: advanced
16547c6ae99SBarry Smith 
166f73e5cebSJed Brown .seealso: DMCompositeGetEntries(), DMCompositeScatter()
16747c6ae99SBarry Smith @*/
1687087cfbeSBarry Smith PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
16947c6ae99SBarry Smith {
17047c6ae99SBarry Smith   va_list                Argp;
17147c6ae99SBarry Smith   PetscErrorCode         ierr;
17247c6ae99SBarry Smith   struct DMCompositeLink *next;
17347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
1745edff71fSBarry Smith   PetscInt               readonly;
17547c6ae99SBarry Smith 
17647c6ae99SBarry Smith   PetscFunctionBegin;
17747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
17847c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
17947c6ae99SBarry Smith   next = com->next;
18047c6ae99SBarry Smith   if (!com->setup) {
181d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
18247c6ae99SBarry Smith   }
18347c6ae99SBarry Smith 
1845edff71fSBarry Smith   ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr);
18547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
18647c6ae99SBarry Smith   va_start(Argp,gvec);
18747c6ae99SBarry Smith   while (next) {
18847c6ae99SBarry Smith     Vec *vec;
18947c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
1909ae5db72SJed Brown     if (vec) {
1919ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,vec);CHKERRQ(ierr);
1925edff71fSBarry Smith       if (readonly) {
1935edff71fSBarry Smith         const PetscScalar *array;
1945edff71fSBarry Smith         ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
1955edff71fSBarry Smith         ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
1965edff71fSBarry Smith         ierr = VecLockPush(*vec);CHKERRQ(ierr);
1975edff71fSBarry Smith         ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
1985edff71fSBarry Smith       } else {
1995edff71fSBarry Smith         PetscScalar *array;
2009ae5db72SJed Brown         ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
2019ae5db72SJed Brown         ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
2029ae5db72SJed Brown         ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
20347c6ae99SBarry Smith       }
2045edff71fSBarry Smith     }
20547c6ae99SBarry Smith     next = next->next;
20647c6ae99SBarry Smith   }
20747c6ae99SBarry Smith   va_end(Argp);
20847c6ae99SBarry Smith   PetscFunctionReturn(0);
20947c6ae99SBarry Smith }
21047c6ae99SBarry Smith 
211f73e5cebSJed Brown /*@C
212f73e5cebSJed Brown     DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
213f73e5cebSJed Brown        representation.
214f73e5cebSJed Brown 
215f73e5cebSJed Brown     Collective on DMComposite
216f73e5cebSJed Brown 
217f73e5cebSJed Brown     Input Parameters:
218f73e5cebSJed Brown +    dm - the packer object
219f73e5cebSJed Brown .    pvec - packed vector
220f73e5cebSJed Brown .    nwanted - number of vectors wanted
2210298fd71SBarry Smith -    wanted - sorted array of vectors wanted, or NULL to get all vectors
222f73e5cebSJed Brown 
223f73e5cebSJed Brown     Output Parameters:
224f73e5cebSJed Brown .    vecs - array of requested global vectors (must be allocated)
225f73e5cebSJed Brown 
226f73e5cebSJed Brown     Notes: Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them
227f73e5cebSJed Brown 
228f73e5cebSJed Brown     Level: advanced
229f73e5cebSJed Brown 
230f73e5cebSJed Brown .seealso: DMCompositeGetAccess(), DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
231f73e5cebSJed Brown @*/
232f73e5cebSJed Brown PetscErrorCode  DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
233f73e5cebSJed Brown {
234f73e5cebSJed Brown   PetscErrorCode         ierr;
235f73e5cebSJed Brown   struct DMCompositeLink *link;
236f73e5cebSJed Brown   PetscInt               i,wnum;
237f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
238bee642f7SBarry Smith   PetscInt               readonly;
239f73e5cebSJed Brown 
240f73e5cebSJed Brown   PetscFunctionBegin;
241f73e5cebSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
242f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
243f73e5cebSJed Brown   if (!com->setup) {
244f73e5cebSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
245f73e5cebSJed Brown   }
246f73e5cebSJed Brown 
247bee642f7SBarry Smith   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
248f73e5cebSJed Brown   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
249f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
250f73e5cebSJed Brown       Vec v;
251f73e5cebSJed Brown       ierr = DMGetGlobalVector(link->dm,&v);CHKERRQ(ierr);
252bee642f7SBarry Smith       if (readonly) {
253bee642f7SBarry Smith         const PetscScalar *array;
254bee642f7SBarry Smith         ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr);
255bee642f7SBarry Smith         ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
256bee642f7SBarry Smith         ierr = VecLockPush(v);CHKERRQ(ierr);
257bee642f7SBarry Smith         ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr);
258bee642f7SBarry Smith       } else {
259bee642f7SBarry Smith         PetscScalar *array;
260f73e5cebSJed Brown         ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
261f73e5cebSJed Brown         ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
262f73e5cebSJed Brown         ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
263bee642f7SBarry Smith       }
264f73e5cebSJed Brown       vecs[wnum++] = v;
265f73e5cebSJed Brown     }
266f73e5cebSJed Brown   }
267f73e5cebSJed Brown   PetscFunctionReturn(0);
268f73e5cebSJed Brown }
269f73e5cebSJed Brown 
2707ac2b803SAlex Fikl /*@C
2717ac2b803SAlex Fikl     DMCompositeGetLocalAccessArray - Allows one to access the individual
2727ac2b803SAlex Fikl     packed vectors in their local representation.
2737ac2b803SAlex Fikl 
2747ac2b803SAlex Fikl     Collective on DMComposite.
2757ac2b803SAlex Fikl 
2767ac2b803SAlex Fikl     Input Parameters:
2777ac2b803SAlex Fikl +    dm - the packer object
2787ac2b803SAlex Fikl .    pvec - packed vector
2797ac2b803SAlex Fikl .    nwanted - number of vectors wanted
2807ac2b803SAlex Fikl -    wanted - sorted array of vectors wanted, or NULL to get all vectors
2817ac2b803SAlex Fikl 
2827ac2b803SAlex Fikl     Output Parameters:
2837ac2b803SAlex Fikl .    vecs - array of requested local vectors (must be allocated)
2847ac2b803SAlex Fikl 
2857ac2b803SAlex Fikl     Notes: Use DMCompositeRestoreLocalAccessArray() to return the vectors
2867ac2b803SAlex Fikl     when you no longer need them.
2877ac2b803SAlex Fikl 
2887ac2b803SAlex Fikl     Level: advanced
2897ac2b803SAlex Fikl 
2907ac2b803SAlex Fikl .seealso: DMCompositeRestoreLocalAccessArray(), DMCompositeGetAccess(),
2917ac2b803SAlex Fikl DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
2927ac2b803SAlex Fikl @*/
2937ac2b803SAlex Fikl PetscErrorCode  DMCompositeGetLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
2947ac2b803SAlex Fikl {
2957ac2b803SAlex Fikl   PetscErrorCode         ierr;
2967ac2b803SAlex Fikl   struct DMCompositeLink *link;
2977ac2b803SAlex Fikl   PetscInt               i,wnum;
2987ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite*)dm->data;
2997ac2b803SAlex Fikl   PetscInt               readonly;
3007ac2b803SAlex Fikl   PetscInt               nlocal = 0;
3017ac2b803SAlex Fikl 
3027ac2b803SAlex Fikl   PetscFunctionBegin;
3037ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3047ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
3057ac2b803SAlex Fikl   if (!com->setup) {
3067ac2b803SAlex Fikl     ierr = DMSetUp(dm);CHKERRQ(ierr);
3077ac2b803SAlex Fikl   }
3087ac2b803SAlex Fikl 
3097ac2b803SAlex Fikl   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
3107ac2b803SAlex Fikl   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
3117ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
3127ac2b803SAlex Fikl       Vec v;
3137ac2b803SAlex Fikl       ierr = DMGetLocalVector(link->dm,&v);CHKERRQ(ierr);
3147ac2b803SAlex Fikl       if (readonly) {
3157ac2b803SAlex Fikl         const PetscScalar *array;
3167ac2b803SAlex Fikl         ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr);
3177ac2b803SAlex Fikl         ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr);
3187ac2b803SAlex Fikl         ierr = VecLockPush(v);CHKERRQ(ierr);
3197ac2b803SAlex Fikl         ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr);
3207ac2b803SAlex Fikl       } else {
3217ac2b803SAlex Fikl         PetscScalar *array;
3227ac2b803SAlex Fikl         ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
3237ac2b803SAlex Fikl         ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr);
3247ac2b803SAlex Fikl         ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
3257ac2b803SAlex Fikl       }
3267ac2b803SAlex Fikl       vecs[wnum++] = v;
3277ac2b803SAlex Fikl     }
3287ac2b803SAlex Fikl 
3297ac2b803SAlex Fikl     nlocal += link->nlocal;
3307ac2b803SAlex Fikl   }
3317ac2b803SAlex Fikl 
3327ac2b803SAlex Fikl   PetscFunctionReturn(0);
3337ac2b803SAlex Fikl }
3347ac2b803SAlex Fikl 
33547c6ae99SBarry Smith /*@C
336aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
33747c6ae99SBarry Smith        representation.
33847c6ae99SBarry Smith 
33947c6ae99SBarry Smith     Collective on DMComposite
34047c6ae99SBarry Smith 
3419ae5db72SJed Brown     Input Parameters:
34247c6ae99SBarry Smith +    dm - the packer object
34347c6ae99SBarry Smith .    gvec - the global vector
3440298fd71SBarry Smith -    Vec* ... - the individual parallel vectors, NULL for those that are not needed
34547c6ae99SBarry Smith 
34647c6ae99SBarry Smith     Level: advanced
34747c6ae99SBarry Smith 
3489ae5db72SJed Brown .seealso  DMCompositeAddDM(), DMCreateGlobalVector(),
3496eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
350aa219208SBarry Smith          DMCompositeRestoreAccess(), DMCompositeGetAccess()
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith @*/
3537087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
35447c6ae99SBarry Smith {
35547c6ae99SBarry Smith   va_list                Argp;
35647c6ae99SBarry Smith   PetscErrorCode         ierr;
35747c6ae99SBarry Smith   struct DMCompositeLink *next;
35847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
3595edff71fSBarry Smith   PetscInt               readonly;
36047c6ae99SBarry Smith 
36147c6ae99SBarry Smith   PetscFunctionBegin;
36247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36347c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
36447c6ae99SBarry Smith   next = com->next;
36547c6ae99SBarry Smith   if (!com->setup) {
366d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
36747c6ae99SBarry Smith   }
36847c6ae99SBarry Smith 
3695edff71fSBarry Smith   ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr);
37047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
37147c6ae99SBarry Smith   va_start(Argp,gvec);
37247c6ae99SBarry Smith   while (next) {
37347c6ae99SBarry Smith     Vec *vec;
37447c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
3759ae5db72SJed Brown     if (vec) {
3769ae5db72SJed Brown       ierr = VecResetArray(*vec);CHKERRQ(ierr);
3775edff71fSBarry Smith       if (readonly) {
3785edff71fSBarry Smith         ierr = VecLockPop(*vec);CHKERRQ(ierr);
3795edff71fSBarry Smith       }
380bee642f7SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,vec);CHKERRQ(ierr);
38147c6ae99SBarry Smith     }
38247c6ae99SBarry Smith     next = next->next;
38347c6ae99SBarry Smith   }
38447c6ae99SBarry Smith   va_end(Argp);
38547c6ae99SBarry Smith   PetscFunctionReturn(0);
38647c6ae99SBarry Smith }
38747c6ae99SBarry Smith 
388f73e5cebSJed Brown /*@C
389f73e5cebSJed Brown     DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray()
390f73e5cebSJed Brown 
391f73e5cebSJed Brown     Collective on DMComposite
392f73e5cebSJed Brown 
393f73e5cebSJed Brown     Input Parameters:
394f73e5cebSJed Brown +    dm - the packer object
395f73e5cebSJed Brown .    pvec - packed vector
396f73e5cebSJed Brown .    nwanted - number of vectors wanted
3970298fd71SBarry Smith .    wanted - sorted array of vectors wanted, or NULL to get all vectors
398f73e5cebSJed Brown -    vecs - array of global vectors to return
399f73e5cebSJed Brown 
400f73e5cebSJed Brown     Level: advanced
401f73e5cebSJed Brown 
402f73e5cebSJed Brown .seealso: DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), DMCompositeScatter(), DMCompositeGather()
403f73e5cebSJed Brown @*/
404f73e5cebSJed Brown PetscErrorCode  DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
405f73e5cebSJed Brown {
406f73e5cebSJed Brown   PetscErrorCode         ierr;
407f73e5cebSJed Brown   struct DMCompositeLink *link;
408f73e5cebSJed Brown   PetscInt               i,wnum;
409f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
410bee642f7SBarry Smith   PetscInt               readonly;
411f73e5cebSJed Brown 
412f73e5cebSJed Brown   PetscFunctionBegin;
413f73e5cebSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
414f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
415f73e5cebSJed Brown   if (!com->setup) {
416f73e5cebSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
417f73e5cebSJed Brown   }
418f73e5cebSJed Brown 
419bee642f7SBarry Smith   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
420f73e5cebSJed Brown   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
421f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
422f73e5cebSJed Brown       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
423bee642f7SBarry Smith       if (readonly) {
424bee642f7SBarry Smith         ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr);
425bee642f7SBarry Smith       }
426f73e5cebSJed Brown       ierr = DMRestoreGlobalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
427f73e5cebSJed Brown       wnum++;
428f73e5cebSJed Brown     }
429f73e5cebSJed Brown   }
430f73e5cebSJed Brown   PetscFunctionReturn(0);
431f73e5cebSJed Brown }
432f73e5cebSJed Brown 
4337ac2b803SAlex Fikl /*@C
4347ac2b803SAlex Fikl     DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray().
4357ac2b803SAlex Fikl 
4367ac2b803SAlex Fikl     Collective on DMComposite.
4377ac2b803SAlex Fikl 
4387ac2b803SAlex Fikl     Input Parameters:
4397ac2b803SAlex Fikl +    dm - the packer object
4407ac2b803SAlex Fikl .    pvec - packed vector
4417ac2b803SAlex Fikl .    nwanted - number of vectors wanted
4427ac2b803SAlex Fikl .    wanted - sorted array of vectors wanted, or NULL to restore all vectors
4437ac2b803SAlex Fikl -    vecs - array of local vectors to return
4447ac2b803SAlex Fikl 
4457ac2b803SAlex Fikl     Level: advanced
4467ac2b803SAlex Fikl 
4477ac2b803SAlex Fikl     Notes:
4487ac2b803SAlex Fikl     nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray()
4497ac2b803SAlex Fikl     otherwise the call will fail.
4507ac2b803SAlex Fikl 
4517ac2b803SAlex Fikl .seealso: DMCompositeGetLocalAccessArray(), DMCompositeRestoreAccessArray(),
4527ac2b803SAlex Fikl DMCompositeRestoreAccess(), DMCompositeRestoreEntries(),
4537ac2b803SAlex Fikl DMCompositeScatter(), DMCompositeGather()
4547ac2b803SAlex Fikl @*/
4557ac2b803SAlex Fikl PetscErrorCode  DMCompositeRestoreLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
4567ac2b803SAlex Fikl {
4577ac2b803SAlex Fikl   PetscErrorCode         ierr;
4587ac2b803SAlex Fikl   struct DMCompositeLink *link;
4597ac2b803SAlex Fikl   PetscInt               i,wnum;
4607ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite*)dm->data;
4617ac2b803SAlex Fikl   PetscInt               readonly;
4627ac2b803SAlex Fikl 
4637ac2b803SAlex Fikl   PetscFunctionBegin;
4647ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4657ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
4667ac2b803SAlex Fikl   if (!com->setup) {
4677ac2b803SAlex Fikl     ierr = DMSetUp(dm);CHKERRQ(ierr);
4687ac2b803SAlex Fikl   }
4697ac2b803SAlex Fikl 
4707ac2b803SAlex Fikl   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
4717ac2b803SAlex Fikl   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
4727ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
4737ac2b803SAlex Fikl       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
4747ac2b803SAlex Fikl       if (readonly) {
4757ac2b803SAlex Fikl         ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr);
4767ac2b803SAlex Fikl       }
4777ac2b803SAlex Fikl       ierr = DMRestoreLocalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
4787ac2b803SAlex Fikl       wnum++;
4797ac2b803SAlex Fikl     }
4807ac2b803SAlex Fikl   }
4817ac2b803SAlex Fikl   PetscFunctionReturn(0);
4827ac2b803SAlex Fikl }
4837ac2b803SAlex Fikl 
48447c6ae99SBarry Smith /*@C
48547c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
48647c6ae99SBarry Smith 
48747c6ae99SBarry Smith     Collective on DMComposite
48847c6ae99SBarry Smith 
4899ae5db72SJed Brown     Input Parameters:
49047c6ae99SBarry Smith +    dm - the packer object
49147c6ae99SBarry Smith .    gvec - the global vector
4920298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for those that are not needed
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith     Level: advanced
49547c6ae99SBarry Smith 
4966f3c3dcfSJed Brown     Notes:
4976f3c3dcfSJed Brown     DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is
4986f3c3dcfSJed Brown     accessible from Fortran.
4996f3c3dcfSJed Brown 
5009ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
5016eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
50247c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
5036f3c3dcfSJed Brown          DMCompositeScatterArray()
50447c6ae99SBarry Smith 
50547c6ae99SBarry Smith @*/
5067087cfbeSBarry Smith PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
50747c6ae99SBarry Smith {
50847c6ae99SBarry Smith   va_list                Argp;
50947c6ae99SBarry Smith   PetscErrorCode         ierr;
51047c6ae99SBarry Smith   struct DMCompositeLink *next;
5118fd8f222SJed Brown   PetscInt               cnt;
51247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
51347c6ae99SBarry Smith 
51447c6ae99SBarry Smith   PetscFunctionBegin;
51547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
51647c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
51747c6ae99SBarry Smith   if (!com->setup) {
518d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
51947c6ae99SBarry Smith   }
52047c6ae99SBarry Smith 
52147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
52247c6ae99SBarry Smith   va_start(Argp,gvec);
5238fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
5249ae5db72SJed Brown     Vec local;
5259ae5db72SJed Brown     local = va_arg(Argp, Vec);
5269ae5db72SJed Brown     if (local) {
5279ae5db72SJed Brown       Vec               global;
5285edff71fSBarry Smith       const PetscScalar *array;
5299ae5db72SJed Brown       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
5309ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
5315edff71fSBarry Smith       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
5329ae5db72SJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
5339ae5db72SJed Brown       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
5349ae5db72SJed Brown       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
5355edff71fSBarry Smith       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
5369ae5db72SJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
5379ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
53847c6ae99SBarry Smith     }
53947c6ae99SBarry Smith   }
54047c6ae99SBarry Smith   va_end(Argp);
54147c6ae99SBarry Smith   PetscFunctionReturn(0);
54247c6ae99SBarry Smith }
54347c6ae99SBarry Smith 
5446f3c3dcfSJed Brown /*@
5456f3c3dcfSJed Brown     DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5466f3c3dcfSJed Brown 
5476f3c3dcfSJed Brown     Collective on DMComposite
5486f3c3dcfSJed Brown 
5496f3c3dcfSJed Brown     Input Parameters:
5506f3c3dcfSJed Brown +    dm - the packer object
5516f3c3dcfSJed Brown .    gvec - the global vector
5526f3c3dcfSJed Brown .    lvecs - array of local vectors, NULL for any that are not needed
5536f3c3dcfSJed Brown 
5546f3c3dcfSJed Brown     Level: advanced
5556f3c3dcfSJed Brown 
5566f3c3dcfSJed Brown     Note:
557907376e6SBarry Smith     This is a non-variadic alternative to DMCompositeScatter()
5586f3c3dcfSJed Brown 
5596f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector()
5606f3c3dcfSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
5616f3c3dcfSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
5626f3c3dcfSJed Brown 
5636f3c3dcfSJed Brown @*/
5646f3c3dcfSJed Brown PetscErrorCode  DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs)
5656f3c3dcfSJed Brown {
5666f3c3dcfSJed Brown   PetscErrorCode         ierr;
5676f3c3dcfSJed Brown   struct DMCompositeLink *next;
5686f3c3dcfSJed Brown   PetscInt               i;
5696f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
5706f3c3dcfSJed Brown 
5716f3c3dcfSJed Brown   PetscFunctionBegin;
5726f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5736f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
5746f3c3dcfSJed Brown   if (!com->setup) {
5756f3c3dcfSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
5766f3c3dcfSJed Brown   }
5776f3c3dcfSJed Brown 
5786f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
5796f3c3dcfSJed Brown   for (i=0,next=com->next; next; next=next->next,i++) {
5806f3c3dcfSJed Brown     if (lvecs[i]) {
5816f3c3dcfSJed Brown       Vec         global;
582c5d31e75SLisandro Dalcin       const PetscScalar *array;
5836f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
5846f3c3dcfSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
585c5d31e75SLisandro Dalcin       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
586c5d31e75SLisandro Dalcin       ierr = VecPlaceArray(global,(PetscScalar*)array+next->rstart);CHKERRQ(ierr);
5876f3c3dcfSJed Brown       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
5886f3c3dcfSJed Brown       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
589c5d31e75SLisandro Dalcin       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
5906f3c3dcfSJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
5916f3c3dcfSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
5926f3c3dcfSJed Brown     }
5936f3c3dcfSJed Brown   }
5946f3c3dcfSJed Brown   PetscFunctionReturn(0);
5956f3c3dcfSJed Brown }
5966f3c3dcfSJed Brown 
59747c6ae99SBarry Smith /*@C
59847c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
59947c6ae99SBarry Smith 
60047c6ae99SBarry Smith     Collective on DMComposite
60147c6ae99SBarry Smith 
60247c6ae99SBarry Smith     Input Parameter:
60347c6ae99SBarry Smith +    dm - the packer object
60447c6ae99SBarry Smith .    gvec - the global vector
605907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6060298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for any that are not needed
60747c6ae99SBarry Smith 
60847c6ae99SBarry Smith     Level: advanced
60947c6ae99SBarry Smith 
6109ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
6116eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
61247c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
61347c6ae99SBarry Smith 
61447c6ae99SBarry Smith @*/
6151dac896bSSatish Balay PetscErrorCode  DMCompositeGather(DM dm,InsertMode imode,Vec gvec,...)
61647c6ae99SBarry Smith {
61747c6ae99SBarry Smith   va_list                Argp;
61847c6ae99SBarry Smith   PetscErrorCode         ierr;
61947c6ae99SBarry Smith   struct DMCompositeLink *next;
62047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
6218fd8f222SJed Brown   PetscInt               cnt;
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith   PetscFunctionBegin;
62447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
62547c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
62647c6ae99SBarry Smith   if (!com->setup) {
627d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
62847c6ae99SBarry Smith   }
62947c6ae99SBarry Smith 
63047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
6311dac896bSSatish Balay   va_start(Argp,gvec);
6328fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
6339ae5db72SJed Brown     Vec local;
6349ae5db72SJed Brown     local = va_arg(Argp, Vec);
6359ae5db72SJed Brown     if (local) {
63647c6ae99SBarry Smith       PetscScalar *array;
6379ae5db72SJed Brown       Vec         global;
6389ae5db72SJed Brown       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
6399ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
6409ae5db72SJed Brown       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
6419ae5db72SJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
6429ae5db72SJed Brown       ierr = DMLocalToGlobalBegin(next->dm,local,imode,global);CHKERRQ(ierr);
6439ae5db72SJed Brown       ierr = DMLocalToGlobalEnd(next->dm,local,imode,global);CHKERRQ(ierr);
6449ae5db72SJed Brown       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
6459ae5db72SJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
6469ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
64747c6ae99SBarry Smith     }
64847c6ae99SBarry Smith   }
64947c6ae99SBarry Smith   va_end(Argp);
65047c6ae99SBarry Smith   PetscFunctionReturn(0);
65147c6ae99SBarry Smith }
65247c6ae99SBarry Smith 
6536f3c3dcfSJed Brown /*@
6546f3c3dcfSJed Brown     DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6556f3c3dcfSJed Brown 
6566f3c3dcfSJed Brown     Collective on DMComposite
6576f3c3dcfSJed Brown 
6586f3c3dcfSJed Brown     Input Parameter:
6596f3c3dcfSJed Brown +    dm - the packer object
6606f3c3dcfSJed Brown .    gvec - the global vector
661907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6626f3c3dcfSJed Brown -    lvecs - the individual sequential vectors, NULL for any that are not needed
6636f3c3dcfSJed Brown 
6646f3c3dcfSJed Brown     Level: advanced
6656f3c3dcfSJed Brown 
6666f3c3dcfSJed Brown     Notes:
6676f3c3dcfSJed Brown     This is a non-variadic alternative to DMCompositeGather().
6686f3c3dcfSJed Brown 
6696f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
6706f3c3dcfSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
6716f3c3dcfSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries(),
6726f3c3dcfSJed Brown @*/
6731dac896bSSatish Balay PetscErrorCode  DMCompositeGatherArray(DM dm,InsertMode imode,Vec gvec,Vec *lvecs)
6746f3c3dcfSJed Brown {
6756f3c3dcfSJed Brown   PetscErrorCode         ierr;
6766f3c3dcfSJed Brown   struct DMCompositeLink *next;
6776f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
6786f3c3dcfSJed Brown   PetscInt               i;
6796f3c3dcfSJed Brown 
6806f3c3dcfSJed Brown   PetscFunctionBegin;
6816f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6826f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
6836f3c3dcfSJed Brown   if (!com->setup) {
6846f3c3dcfSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
6856f3c3dcfSJed Brown   }
6866f3c3dcfSJed Brown 
6876f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
6886f3c3dcfSJed Brown   for (next=com->next,i=0; next; next=next->next,i++) {
6896f3c3dcfSJed Brown     if (lvecs[i]) {
6906f3c3dcfSJed Brown       PetscScalar *array;
6916f3c3dcfSJed Brown       Vec         global;
6926f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
6936f3c3dcfSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
6946f3c3dcfSJed Brown       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
6956f3c3dcfSJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
6966f3c3dcfSJed Brown       ierr = DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
6976f3c3dcfSJed Brown       ierr = DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
6986f3c3dcfSJed Brown       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
6996f3c3dcfSJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
7006f3c3dcfSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
7016f3c3dcfSJed Brown     }
7026f3c3dcfSJed Brown   }
7036f3c3dcfSJed Brown   PetscFunctionReturn(0);
7046f3c3dcfSJed Brown }
7056f3c3dcfSJed Brown 
70647c6ae99SBarry Smith /*@C
707aa219208SBarry Smith     DMCompositeAddDM - adds a DM vector to a DMComposite
70847c6ae99SBarry Smith 
70947c6ae99SBarry Smith     Collective on DMComposite
71047c6ae99SBarry Smith 
71147c6ae99SBarry Smith     Input Parameter:
7129b52a9b5SPatrick Sanan +    dmc - the DMComposite (packer) object
7139b52a9b5SPatrick Sanan -    dm - the DM object; if the DM is a DMDA you will need to cast it with a (DM)
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith     Level: advanced
71647c6ae99SBarry Smith 
7170c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
7186eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
71947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
72047c6ae99SBarry Smith 
72147c6ae99SBarry Smith @*/
7227087cfbeSBarry Smith PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
72347c6ae99SBarry Smith {
72447c6ae99SBarry Smith   PetscErrorCode         ierr;
72506ebdd98SJed Brown   PetscInt               n,nlocal;
72647c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
72706ebdd98SJed Brown   Vec                    global,local;
72847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmc->data;
72947c6ae99SBarry Smith 
73047c6ae99SBarry Smith   PetscFunctionBegin;
73147c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
73247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
73347c6ae99SBarry Smith   next = com->next;
734ce94432eSBarry Smith   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
73547c6ae99SBarry Smith 
73647c6ae99SBarry Smith   /* create new link */
737b00a9115SJed Brown   ierr = PetscNew(&mine);CHKERRQ(ierr);
73847c6ae99SBarry Smith   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
73947c6ae99SBarry Smith   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
74047c6ae99SBarry Smith   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
74147c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
74206ebdd98SJed Brown   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
74306ebdd98SJed Brown   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
74406ebdd98SJed Brown   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
7458865f1eaSKarl Rupp 
74647c6ae99SBarry Smith   mine->n      = n;
74706ebdd98SJed Brown   mine->nlocal = nlocal;
74847c6ae99SBarry Smith   mine->dm     = dm;
7490298fd71SBarry Smith   mine->next   = NULL;
75047c6ae99SBarry Smith   com->n      += n;
7517ac2b803SAlex Fikl   com->nghost += nlocal;
75247c6ae99SBarry Smith 
75347c6ae99SBarry Smith   /* add to end of list */
7548865f1eaSKarl Rupp   if (!next) com->next = mine;
7558865f1eaSKarl Rupp   else {
75647c6ae99SBarry Smith     while (next->next) next = next->next;
75747c6ae99SBarry Smith     next->next = mine;
75847c6ae99SBarry Smith   }
75947c6ae99SBarry Smith   com->nDM++;
76047c6ae99SBarry Smith   com->nmine++;
76147c6ae99SBarry Smith   PetscFunctionReturn(0);
76247c6ae99SBarry Smith }
76347c6ae99SBarry Smith 
7649804daf3SBarry Smith #include <petscdraw.h>
76526887b52SJed Brown PETSC_EXTERN PetscErrorCode  VecView_MPI(Vec,PetscViewer);
7667087cfbeSBarry Smith PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
76747c6ae99SBarry Smith {
76847c6ae99SBarry Smith   DM                     dm;
76947c6ae99SBarry Smith   PetscErrorCode         ierr;
77047c6ae99SBarry Smith   struct DMCompositeLink *next;
77147c6ae99SBarry Smith   PetscBool              isdraw;
772cef07954SSatish Balay   DM_Composite           *com;
77347c6ae99SBarry Smith 
77447c6ae99SBarry Smith   PetscFunctionBegin;
775c688c046SMatthew G Knepley   ierr = VecGetDM(gvec, &dm);CHKERRQ(ierr);
776ce94432eSBarry Smith   if (!dm) SETERRQ(PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
77747c6ae99SBarry Smith   com  = (DM_Composite*)dm->data;
77847c6ae99SBarry Smith   next = com->next;
77947c6ae99SBarry Smith 
780251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
78147c6ae99SBarry Smith   if (!isdraw) {
78247c6ae99SBarry Smith     /* do I really want to call this? */
78347c6ae99SBarry Smith     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
78447c6ae99SBarry Smith   } else {
78547c6ae99SBarry Smith     PetscInt cnt = 0;
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
78847c6ae99SBarry Smith     while (next) {
78947c6ae99SBarry Smith       Vec               vec;
790ca4278abSLisandro Dalcin       const PetscScalar *array;
79147c6ae99SBarry Smith       PetscInt          bs;
79247c6ae99SBarry Smith 
7939ae5db72SJed Brown       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
7949ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&vec);CHKERRQ(ierr);
795ca4278abSLisandro Dalcin       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
796ca4278abSLisandro Dalcin       ierr = VecPlaceArray(vec,(PetscScalar*)array+next->rstart);CHKERRQ(ierr);
797ca4278abSLisandro Dalcin       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
79847c6ae99SBarry Smith       ierr = VecView(vec,viewer);CHKERRQ(ierr);
7999ae5db72SJed Brown       ierr = VecResetArray(vec);CHKERRQ(ierr);
800ca4278abSLisandro Dalcin       ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
8019ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&vec);CHKERRQ(ierr);
80247c6ae99SBarry Smith       ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
80347c6ae99SBarry Smith       cnt += bs;
80447c6ae99SBarry Smith       next = next->next;
80547c6ae99SBarry Smith     }
80647c6ae99SBarry Smith     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
80747c6ae99SBarry Smith   }
80847c6ae99SBarry Smith   PetscFunctionReturn(0);
80947c6ae99SBarry Smith }
81047c6ae99SBarry Smith 
8117087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
81247c6ae99SBarry Smith {
81347c6ae99SBarry Smith   PetscErrorCode ierr;
81447c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
81547c6ae99SBarry Smith 
81647c6ae99SBarry Smith   PetscFunctionBegin;
81747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
818d7bf68aeSBarry Smith   ierr = DMSetUp(dm);CHKERRQ(ierr);
819ce94432eSBarry Smith   ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),com->n,com->N,gvec);CHKERRQ(ierr);
820c688c046SMatthew G Knepley   ierr = VecSetDM(*gvec, dm);CHKERRQ(ierr);
82147c6ae99SBarry Smith   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);CHKERRQ(ierr);
82247c6ae99SBarry Smith   PetscFunctionReturn(0);
82347c6ae99SBarry Smith }
82447c6ae99SBarry Smith 
8257087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
82647c6ae99SBarry Smith {
82747c6ae99SBarry Smith   PetscErrorCode ierr;
82847c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith   PetscFunctionBegin;
83147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
83247c6ae99SBarry Smith   if (!com->setup) {
833d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
83447c6ae99SBarry Smith   }
835f0e01b1fSVincent Le Chenadec   ierr = VecCreateSeq(PETSC_COMM_SELF,com->nghost,lvec);CHKERRQ(ierr);
836c688c046SMatthew G Knepley   ierr = VecSetDM(*lvec, dm);CHKERRQ(ierr);
83747c6ae99SBarry Smith   PetscFunctionReturn(0);
83847c6ae99SBarry Smith }
83947c6ae99SBarry Smith 
84047c6ae99SBarry Smith /*@C
8419ae5db72SJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
84247c6ae99SBarry Smith 
84306ebdd98SJed Brown     Collective on DM
84447c6ae99SBarry Smith 
84547c6ae99SBarry Smith     Input Parameter:
84647c6ae99SBarry Smith .    dm - the packer object
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith     Output Parameters:
8499ae5db72SJed Brown .    ltogs - the individual mappings for each packed vector. Note that this includes
8509ae5db72SJed Brown            all the ghost points that individual ghosted DMDA's may have.
85147c6ae99SBarry Smith 
85247c6ae99SBarry Smith     Level: advanced
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith     Notes:
8556eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
85647c6ae99SBarry Smith 
8579ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
85847c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
85947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith @*/
8627087cfbeSBarry Smith PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
86347c6ae99SBarry Smith {
86447c6ae99SBarry Smith   PetscErrorCode         ierr;
86547c6ae99SBarry Smith   PetscInt               i,*idx,n,cnt;
86647c6ae99SBarry Smith   struct DMCompositeLink *next;
86747c6ae99SBarry Smith   PetscMPIInt            rank;
86847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
86947c6ae99SBarry Smith 
87047c6ae99SBarry Smith   PetscFunctionBegin;
87147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
872728e99d6SJed Brown   ierr = DMSetUp(dm);CHKERRQ(ierr);
873854ce69bSBarry Smith   ierr = PetscMalloc1(com->nDM,ltogs);CHKERRQ(ierr);
87447c6ae99SBarry Smith   next = com->next;
875ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
87647c6ae99SBarry Smith 
87747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
87847c6ae99SBarry Smith   cnt = 0;
87947c6ae99SBarry Smith   while (next) {
8806eb61c8cSJed Brown     ISLocalToGlobalMapping ltog;
8816eb61c8cSJed Brown     PetscMPIInt            size;
88286994e45SJed Brown     const PetscInt         *suboff,*indices;
8836eb61c8cSJed Brown     Vec                    global;
88447c6ae99SBarry Smith 
8856eb61c8cSJed Brown     /* Get sub-DM global indices for each local dof */
8861411c6eeSJed Brown     ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
8876eb61c8cSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
88886994e45SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
889785e854fSJed Brown     ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr);
89047c6ae99SBarry Smith 
8916eb61c8cSJed Brown     /* Get the offsets for the sub-DM global vector */
8926eb61c8cSJed Brown     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
8936eb61c8cSJed Brown     ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
894ce94432eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)global),&size);CHKERRQ(ierr);
8956eb61c8cSJed Brown 
8966eb61c8cSJed Brown     /* Shift the sub-DM definition of the global space to the composite global space */
8976eb61c8cSJed Brown     for (i=0; i<n; i++) {
89886994e45SJed Brown       PetscInt subi = indices[i],lo = 0,hi = size,t;
8996eb61c8cSJed Brown       /* Binary search to find which rank owns subi */
9006eb61c8cSJed Brown       while (hi-lo > 1) {
9016eb61c8cSJed Brown         t = lo + (hi-lo)/2;
9026eb61c8cSJed Brown         if (suboff[t] > subi) hi = t;
9036eb61c8cSJed Brown         else                  lo = t;
9046eb61c8cSJed Brown       }
9056eb61c8cSJed Brown       idx[i] = subi - suboff[lo] + next->grstarts[lo];
9066eb61c8cSJed Brown     }
90786994e45SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
908f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
9096eb61c8cSJed Brown     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
91047c6ae99SBarry Smith     next = next->next;
91147c6ae99SBarry Smith     cnt++;
91247c6ae99SBarry Smith   }
91347c6ae99SBarry Smith   PetscFunctionReturn(0);
91447c6ae99SBarry Smith }
91547c6ae99SBarry Smith 
91687c85e80SJed Brown /*@C
9179ae5db72SJed Brown    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
91887c85e80SJed Brown 
91987c85e80SJed Brown    Not Collective
92087c85e80SJed Brown 
92187c85e80SJed Brown    Input Arguments:
92287c85e80SJed Brown . dm - composite DM
92387c85e80SJed Brown 
92487c85e80SJed Brown    Output Arguments:
92587c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
92687c85e80SJed Brown 
92787c85e80SJed Brown    Level: intermediate
92887c85e80SJed Brown 
92987c85e80SJed Brown    Notes:
93087c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
93187c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
9329ae5db72SJed Brown    scatter to a composite local vector.  The user should not typically need to know which is being done.
93387c85e80SJed Brown 
93487c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
93587c85e80SJed Brown 
93687c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
93787c85e80SJed Brown 
93887c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
93987c85e80SJed Brown 
94087c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
94187c85e80SJed Brown @*/
9427087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
94387c85e80SJed Brown {
94487c85e80SJed Brown   PetscErrorCode         ierr;
94587c85e80SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
94687c85e80SJed Brown   struct DMCompositeLink *link;
94787c85e80SJed Brown   PetscInt               cnt,start;
94887c85e80SJed Brown 
94987c85e80SJed Brown   PetscFunctionBegin;
95087c85e80SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
95187c85e80SJed Brown   PetscValidPointer(is,2);
952785e854fSJed Brown   ierr = PetscMalloc1(com->nmine,is);CHKERRQ(ierr);
95306ebdd98SJed Brown   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
954520db06cSJed Brown     PetscInt bs;
9559ae5db72SJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
9561411c6eeSJed Brown     ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
957520db06cSJed Brown     ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
958520db06cSJed Brown   }
95987c85e80SJed Brown   PetscFunctionReturn(0);
96087c85e80SJed Brown }
96187c85e80SJed Brown 
96247c6ae99SBarry Smith /*@C
96347c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
96447c6ae99SBarry Smith 
96547c6ae99SBarry Smith     Collective on DMComposite
96647c6ae99SBarry Smith 
96747c6ae99SBarry Smith     Input Parameter:
96847c6ae99SBarry Smith .    dm - the packer object
96947c6ae99SBarry Smith 
97047c6ae99SBarry Smith     Output Parameters:
97147c6ae99SBarry Smith .    is - the array of index sets
97247c6ae99SBarry Smith 
97347c6ae99SBarry Smith     Level: advanced
97447c6ae99SBarry Smith 
97547c6ae99SBarry Smith     Notes:
97647c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
97747c6ae99SBarry Smith 
97847c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
97947c6ae99SBarry Smith 
9806eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
9816eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
9826eb61c8cSJed Brown        indices.
98347c6ae99SBarry Smith 
984f3cb0f7eSJed Brown     Fortran Notes:
985f3cb0f7eSJed Brown 
986f3cb0f7eSJed Brown        The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM().
987f3cb0f7eSJed Brown 
9889ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
98947c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
99047c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
99147c6ae99SBarry Smith 
99247c6ae99SBarry Smith @*/
9937087cfbeSBarry Smith PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
99447c6ae99SBarry Smith {
99547c6ae99SBarry Smith   PetscErrorCode         ierr;
99666bb578eSMark Adams   PetscInt               cnt = 0;
99747c6ae99SBarry Smith   struct DMCompositeLink *next;
99847c6ae99SBarry Smith   PetscMPIInt            rank;
99947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
100047c6ae99SBarry Smith 
100147c6ae99SBarry Smith   PetscFunctionBegin;
100247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1003854ce69bSBarry Smith   ierr = PetscMalloc1(com->nDM,is);CHKERRQ(ierr);
100447c6ae99SBarry Smith   next = com->next;
1005ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
100647c6ae99SBarry Smith 
100747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
100847c6ae99SBarry Smith   while (next) {
100966bb578eSMark Adams     ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt]);CHKERRQ(ierr);
10100f21e855SMatthew G. Knepley     if (dm->prob) {
101165c226d8SMatthew G. Knepley       MatNullSpace space;
101265c226d8SMatthew G. Knepley       Mat          pmat;
10130f21e855SMatthew G. Knepley       PetscObject  disc;
10140f21e855SMatthew G. Knepley       PetscInt     Nf;
101565c226d8SMatthew G. Knepley 
10162764a2aaSMatthew G. Knepley       ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
1017f24dd8d2SMatthew G. Knepley       if (cnt < Nf) {
10182764a2aaSMatthew G. Knepley         ierr = PetscDSGetDiscretization(dm->prob, cnt, &disc);CHKERRQ(ierr);
10190f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nullspace", (PetscObject*) &space);CHKERRQ(ierr);
1020aac2dd2dSMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space);CHKERRQ(ierr);}
10210f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space);CHKERRQ(ierr);
1022aac2dd2dSMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space);CHKERRQ(ierr);}
10230f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
1024aac2dd2dSMatthew G. Knepley         if (pmat) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat);CHKERRQ(ierr);}
102565c226d8SMatthew G. Knepley       }
1026f24dd8d2SMatthew G. Knepley     }
102747c6ae99SBarry Smith     cnt++;
102847c6ae99SBarry Smith     next = next->next;
102947c6ae99SBarry Smith   }
103047c6ae99SBarry Smith   PetscFunctionReturn(0);
103147c6ae99SBarry Smith }
103247c6ae99SBarry Smith 
103321c9b008SJed Brown PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields)
10344d343eeaSMatthew G Knepley {
10354d343eeaSMatthew G Knepley   PetscInt       nDM;
10364d343eeaSMatthew G Knepley   DM             *dms;
10374d343eeaSMatthew G Knepley   PetscInt       i;
10384d343eeaSMatthew G Knepley   PetscErrorCode ierr;
10394d343eeaSMatthew G Knepley 
10404d343eeaSMatthew G Knepley   PetscFunctionBegin;
10414d343eeaSMatthew G Knepley   ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
10428865f1eaSKarl Rupp   if (numFields) *numFields = nDM;
10434d343eeaSMatthew G Knepley   ierr = DMCompositeGetGlobalISs(dm, fields);CHKERRQ(ierr);
10444d343eeaSMatthew G Knepley   if (fieldNames) {
1045785e854fSJed Brown     ierr = PetscMalloc1(nDM, &dms);CHKERRQ(ierr);
1046785e854fSJed Brown     ierr = PetscMalloc1(nDM, fieldNames);CHKERRQ(ierr);
10474d343eeaSMatthew G Knepley     ierr = DMCompositeGetEntriesArray(dm, dms);CHKERRQ(ierr);
10484d343eeaSMatthew G Knepley     for (i=0; i<nDM; i++) {
10494d343eeaSMatthew G Knepley       char       buf[256];
10504d343eeaSMatthew G Knepley       const char *splitname;
10514d343eeaSMatthew G Knepley 
10524d343eeaSMatthew G Knepley       /* Split naming precedence: object name, prefix, number */
10534d343eeaSMatthew G Knepley       splitname = ((PetscObject) dm)->name;
10544d343eeaSMatthew G Knepley       if (!splitname) {
10554d343eeaSMatthew G Knepley         ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr);
10564d343eeaSMatthew G Knepley         if (splitname) {
10574d343eeaSMatthew G Knepley           size_t len;
10588caf3d72SBarry Smith           ierr                 = PetscStrncpy(buf,splitname,sizeof(buf));CHKERRQ(ierr);
10598caf3d72SBarry Smith           buf[sizeof(buf) - 1] = 0;
10604d343eeaSMatthew G Knepley           ierr                 = PetscStrlen(buf,&len);CHKERRQ(ierr);
10614d343eeaSMatthew G Knepley           if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
10624d343eeaSMatthew G Knepley           splitname = buf;
10634d343eeaSMatthew G Knepley         }
10644d343eeaSMatthew G Knepley       }
10654d343eeaSMatthew G Knepley       if (!splitname) {
10668caf3d72SBarry Smith         ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",i);CHKERRQ(ierr);
10674d343eeaSMatthew G Knepley         splitname = buf;
10684d343eeaSMatthew G Knepley       }
106921c9b008SJed Brown       ierr = PetscStrallocpy(splitname,&(*fieldNames)[i]);CHKERRQ(ierr);
10704d343eeaSMatthew G Knepley     }
10714d343eeaSMatthew G Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
10724d343eeaSMatthew G Knepley   }
10734d343eeaSMatthew G Knepley   PetscFunctionReturn(0);
10744d343eeaSMatthew G Knepley }
10754d343eeaSMatthew G Knepley 
1076e7c4fc90SDmitry Karpeev /*
1077e7c4fc90SDmitry Karpeev  This could take over from DMCreateFieldIS(), as it is more general,
10780298fd71SBarry Smith  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1079e7c4fc90SDmitry Karpeev  At this point it's probably best to be less intrusive, however.
1080e7c4fc90SDmitry Karpeev  */
108116621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
1082e7c4fc90SDmitry Karpeev {
1083e7c4fc90SDmitry Karpeev   PetscInt       nDM;
1084e7c4fc90SDmitry Karpeev   PetscInt       i;
1085e7c4fc90SDmitry Karpeev   PetscErrorCode ierr;
1086e7c4fc90SDmitry Karpeev 
1087e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
1088e7c4fc90SDmitry Karpeev   ierr = DMCreateFieldIS_Composite(dm, len, namelist, islist);CHKERRQ(ierr);
1089e7c4fc90SDmitry Karpeev   if (dmlist) {
1090e7c4fc90SDmitry Karpeev     ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
1091785e854fSJed Brown     ierr = PetscMalloc1(nDM, dmlist);CHKERRQ(ierr);
1092e7c4fc90SDmitry Karpeev     ierr = DMCompositeGetEntriesArray(dm, *dmlist);CHKERRQ(ierr);
1093e7c4fc90SDmitry Karpeev     for (i=0; i<nDM; i++) {
1094e7c4fc90SDmitry Karpeev       ierr = PetscObjectReference((PetscObject)((*dmlist)[i]));CHKERRQ(ierr);
1095e7c4fc90SDmitry Karpeev     }
1096e7c4fc90SDmitry Karpeev   }
1097e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
1098e7c4fc90SDmitry Karpeev }
1099e7c4fc90SDmitry Karpeev 
1100e7c4fc90SDmitry Karpeev 
1101e7c4fc90SDmitry Karpeev 
110247c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
110347c6ae99SBarry Smith /*@C
11049ae5db72SJed Brown     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
110547c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
110647c6ae99SBarry Smith 
110747c6ae99SBarry Smith     Not Collective
110847c6ae99SBarry Smith 
110947c6ae99SBarry Smith     Input Parameter:
111047c6ae99SBarry Smith .    dm - the packer object
111147c6ae99SBarry Smith 
111247c6ae99SBarry Smith     Output Parameter:
11139ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
111447c6ae99SBarry Smith 
111547c6ae99SBarry Smith     Level: advanced
111647c6ae99SBarry Smith 
11179ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
11186eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
111947c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
112047c6ae99SBarry Smith 
112147c6ae99SBarry Smith @*/
11227087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
112347c6ae99SBarry Smith {
112447c6ae99SBarry Smith   va_list                Argp;
112547c6ae99SBarry Smith   PetscErrorCode         ierr;
112647c6ae99SBarry Smith   struct DMCompositeLink *next;
112747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
112847c6ae99SBarry Smith 
112947c6ae99SBarry Smith   PetscFunctionBegin;
113047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
113147c6ae99SBarry Smith   next = com->next;
113247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
113347c6ae99SBarry Smith   va_start(Argp,dm);
113447c6ae99SBarry Smith   while (next) {
113547c6ae99SBarry Smith     Vec *vec;
113647c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
113706930112SJed Brown     if (vec) {ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);}
113847c6ae99SBarry Smith     next = next->next;
113947c6ae99SBarry Smith   }
114047c6ae99SBarry Smith   va_end(Argp);
114147c6ae99SBarry Smith   PetscFunctionReturn(0);
114247c6ae99SBarry Smith }
114347c6ae99SBarry Smith 
114447c6ae99SBarry Smith /*@C
11459ae5db72SJed Brown     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
114647c6ae99SBarry Smith 
114747c6ae99SBarry Smith     Not Collective
114847c6ae99SBarry Smith 
114947c6ae99SBarry Smith     Input Parameter:
115047c6ae99SBarry Smith .    dm - the packer object
115147c6ae99SBarry Smith 
115247c6ae99SBarry Smith     Output Parameter:
11539ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
115447c6ae99SBarry Smith 
115547c6ae99SBarry Smith     Level: advanced
115647c6ae99SBarry Smith 
11579ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
11586eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
115947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
116047c6ae99SBarry Smith 
116147c6ae99SBarry Smith @*/
11627087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
116347c6ae99SBarry Smith {
116447c6ae99SBarry Smith   va_list                Argp;
116547c6ae99SBarry Smith   PetscErrorCode         ierr;
116647c6ae99SBarry Smith   struct DMCompositeLink *next;
116747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
116847c6ae99SBarry Smith 
116947c6ae99SBarry Smith   PetscFunctionBegin;
117047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
117147c6ae99SBarry Smith   next = com->next;
117247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
117347c6ae99SBarry Smith   va_start(Argp,dm);
117447c6ae99SBarry Smith   while (next) {
117547c6ae99SBarry Smith     Vec *vec;
117647c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
117706930112SJed Brown     if (vec) {ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);}
117847c6ae99SBarry Smith     next = next->next;
117947c6ae99SBarry Smith   }
118047c6ae99SBarry Smith   va_end(Argp);
118147c6ae99SBarry Smith   PetscFunctionReturn(0);
118247c6ae99SBarry Smith }
118347c6ae99SBarry Smith 
118447c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
118547c6ae99SBarry Smith /*@C
11869ae5db72SJed Brown     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
118747c6ae99SBarry Smith 
118847c6ae99SBarry Smith     Not Collective
118947c6ae99SBarry Smith 
119047c6ae99SBarry Smith     Input Parameter:
119147c6ae99SBarry Smith .    dm - the packer object
119247c6ae99SBarry Smith 
119347c6ae99SBarry Smith     Output Parameter:
11949ae5db72SJed Brown .   DM ... - the individual entries (DMs)
119547c6ae99SBarry Smith 
119647c6ae99SBarry Smith     Level: advanced
119747c6ae99SBarry Smith 
11982fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
11996eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
120047c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
120147c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
120247c6ae99SBarry Smith 
120347c6ae99SBarry Smith @*/
12047087cfbeSBarry Smith PetscErrorCode  DMCompositeGetEntries(DM dm,...)
120547c6ae99SBarry Smith {
120647c6ae99SBarry Smith   va_list                Argp;
120747c6ae99SBarry Smith   struct DMCompositeLink *next;
120847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
120947c6ae99SBarry Smith 
121047c6ae99SBarry Smith   PetscFunctionBegin;
121147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
121247c6ae99SBarry Smith   next = com->next;
121347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
121447c6ae99SBarry Smith   va_start(Argp,dm);
121547c6ae99SBarry Smith   while (next) {
121647c6ae99SBarry Smith     DM *dmn;
121747c6ae99SBarry Smith     dmn = va_arg(Argp, DM*);
12189ae5db72SJed Brown     if (dmn) *dmn = next->dm;
121947c6ae99SBarry Smith     next = next->next;
122047c6ae99SBarry Smith   }
122147c6ae99SBarry Smith   va_end(Argp);
122247c6ae99SBarry Smith   PetscFunctionReturn(0);
122347c6ae99SBarry Smith }
122447c6ae99SBarry Smith 
1225dbab29e1SMark F. Adams /*@C
12262fa5ba8aSJed Brown     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
12272fa5ba8aSJed Brown 
12282fa5ba8aSJed Brown     Not Collective
12292fa5ba8aSJed Brown 
12302fa5ba8aSJed Brown     Input Parameter:
1231907376e6SBarry Smith .    dm - the packer object
1232907376e6SBarry Smith 
1233907376e6SBarry Smith     Output Parameter:
1234907376e6SBarry Smith .    dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs
12352fa5ba8aSJed Brown 
12362fa5ba8aSJed Brown     Level: advanced
12372fa5ba8aSJed Brown 
12382fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
12392fa5ba8aSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
12402fa5ba8aSJed Brown          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
12412fa5ba8aSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
12422fa5ba8aSJed Brown 
12432fa5ba8aSJed Brown @*/
12442fa5ba8aSJed Brown PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
12452fa5ba8aSJed Brown {
12462fa5ba8aSJed Brown   struct DMCompositeLink *next;
12472fa5ba8aSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
12482fa5ba8aSJed Brown   PetscInt               i;
12492fa5ba8aSJed Brown 
12502fa5ba8aSJed Brown   PetscFunctionBegin;
12512fa5ba8aSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12522fa5ba8aSJed Brown   /* loop over packed objects, handling one at at time */
12532fa5ba8aSJed Brown   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
12542fa5ba8aSJed Brown   PetscFunctionReturn(0);
12552fa5ba8aSJed Brown }
12562fa5ba8aSJed Brown 
1257e10fd815SStefano Zampini typedef struct {
1258e10fd815SStefano Zampini   DM          dm;
1259e10fd815SStefano Zampini   PetscViewer *subv;
1260e10fd815SStefano Zampini   Vec         *vecs;
1261e10fd815SStefano Zampini } GLVisViewerCtx;
1262e10fd815SStefano Zampini 
1263e10fd815SStefano Zampini static PetscErrorCode  DestroyGLVisViewerCtx_Private(void *vctx)
1264e10fd815SStefano Zampini {
1265e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
1266e10fd815SStefano Zampini   PetscInt       i,n;
1267e10fd815SStefano Zampini   PetscErrorCode ierr;
1268e10fd815SStefano Zampini 
1269e10fd815SStefano Zampini   PetscFunctionBegin;
1270e10fd815SStefano Zampini   ierr = DMCompositeGetNumberDM(ctx->dm,&n);CHKERRQ(ierr);
1271e10fd815SStefano Zampini   for (i = 0; i < n; i++) {
1272e10fd815SStefano Zampini     ierr = PetscViewerDestroy(&ctx->subv[i]);CHKERRQ(ierr);
1273e10fd815SStefano Zampini   }
1274e10fd815SStefano Zampini   ierr = PetscFree2(ctx->subv,ctx->vecs);CHKERRQ(ierr);
1275e10fd815SStefano Zampini   ierr = DMDestroy(&ctx->dm);CHKERRQ(ierr);
1276e10fd815SStefano Zampini   ierr = PetscFree(ctx);CHKERRQ(ierr);
1277e10fd815SStefano Zampini   PetscFunctionReturn(0);
1278e10fd815SStefano Zampini }
1279e10fd815SStefano Zampini 
1280e10fd815SStefano Zampini static PetscErrorCode  DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
1281e10fd815SStefano Zampini {
1282e10fd815SStefano Zampini   Vec            X = (Vec)oX;
1283e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
1284e10fd815SStefano Zampini   PetscInt       i,n,cumf;
1285e10fd815SStefano Zampini   PetscErrorCode ierr;
1286e10fd815SStefano Zampini 
1287e10fd815SStefano Zampini   PetscFunctionBegin;
1288e10fd815SStefano Zampini   ierr = DMCompositeGetNumberDM(ctx->dm,&n);CHKERRQ(ierr);
1289e10fd815SStefano Zampini   ierr = DMCompositeGetAccessArray(ctx->dm,X,n,NULL,ctx->vecs);CHKERRQ(ierr);
1290e10fd815SStefano Zampini   for (i = 0, cumf = 0; i < n; i++) {
1291e10fd815SStefano Zampini     PetscErrorCode (*g2l)(PetscObject,PetscInt,PetscObject[],void*);
1292e10fd815SStefano Zampini     void           *fctx;
1293e10fd815SStefano Zampini     PetscInt       nfi;
1294e10fd815SStefano Zampini 
1295e10fd815SStefano Zampini     ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nfi,NULL,NULL,&g2l,NULL,&fctx);CHKERRQ(ierr);
1296e10fd815SStefano Zampini     if (!nfi) continue;
1297e10fd815SStefano Zampini     if (g2l) {
1298e10fd815SStefano Zampini       ierr = (*g2l)((PetscObject)ctx->vecs[i],nfi,oXfield+cumf,fctx);CHKERRQ(ierr);
1299e10fd815SStefano Zampini     } else {
1300e10fd815SStefano Zampini       ierr = VecCopy(ctx->vecs[i],(Vec)(oXfield[cumf]));CHKERRQ(ierr);
1301e10fd815SStefano Zampini     }
1302e10fd815SStefano Zampini     cumf += nfi;
1303e10fd815SStefano Zampini   }
1304e10fd815SStefano Zampini   ierr = DMCompositeRestoreAccessArray(ctx->dm,X,n,NULL,ctx->vecs);CHKERRQ(ierr);
1305e10fd815SStefano Zampini   PetscFunctionReturn(0);
1306e10fd815SStefano Zampini }
1307e10fd815SStefano Zampini 
1308e10fd815SStefano Zampini static PetscErrorCode  DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer)
1309e10fd815SStefano Zampini {
1310e10fd815SStefano Zampini   DM             dm = (DM)odm, *dms;
1311e10fd815SStefano Zampini   Vec            *Ufds;
1312e10fd815SStefano Zampini   GLVisViewerCtx *ctx;
1313e10fd815SStefano Zampini   PetscInt       i,n,tnf,*sdim;
1314e10fd815SStefano Zampini   char           **fecs;
1315e10fd815SStefano Zampini   PetscErrorCode ierr;
1316e10fd815SStefano Zampini 
1317e10fd815SStefano Zampini   PetscFunctionBegin;
1318e10fd815SStefano Zampini   ierr = PetscNew(&ctx);CHKERRQ(ierr);
1319e10fd815SStefano Zampini   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
1320e10fd815SStefano Zampini   ctx->dm = dm;
1321e10fd815SStefano Zampini   ierr = DMCompositeGetNumberDM(dm,&n);CHKERRQ(ierr);
1322e10fd815SStefano Zampini   ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr);
1323e10fd815SStefano Zampini   ierr = DMCompositeGetEntriesArray(dm,dms);CHKERRQ(ierr);
1324e10fd815SStefano Zampini   ierr = PetscMalloc2(n,&ctx->subv,n,&ctx->vecs);CHKERRQ(ierr);
1325e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1326e10fd815SStefano Zampini     PetscInt nf;
1327e10fd815SStefano Zampini 
1328e10fd815SStefano Zampini     ierr = PetscViewerCreate(PetscObjectComm(odm),&ctx->subv[i]);CHKERRQ(ierr);
1329e10fd815SStefano Zampini     ierr = PetscViewerSetType(ctx->subv[i],PETSCVIEWERGLVIS);CHKERRQ(ierr);
1330e10fd815SStefano Zampini     ierr = PetscViewerGLVisSetDM_Private(ctx->subv[i],(PetscObject)dms[i]);CHKERRQ(ierr);
1331e10fd815SStefano Zampini     ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1332e10fd815SStefano Zampini     tnf += nf;
1333e10fd815SStefano Zampini   }
1334e10fd815SStefano Zampini   ierr = PetscFree(dms);CHKERRQ(ierr);
1335e10fd815SStefano Zampini   ierr = PetscMalloc3(tnf,&fecs,tnf,&sdim,tnf,&Ufds);CHKERRQ(ierr);
1336e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1337e10fd815SStefano Zampini     PetscInt   *sd,nf,f;
1338e10fd815SStefano Zampini     const char **fec;
1339e10fd815SStefano Zampini     Vec        *Uf;
1340e10fd815SStefano Zampini 
1341e10fd815SStefano Zampini     ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,&fec,&sd,NULL,(PetscObject**)&Uf,NULL);CHKERRQ(ierr);
1342e10fd815SStefano Zampini     for (f = 0; f < nf; f++) {
1343e10fd815SStefano Zampini       ierr = PetscStrallocpy(fec[f],&fecs[tnf+f]);CHKERRQ(ierr);
1344e10fd815SStefano Zampini       Ufds[tnf+f] = Uf[f];
1345e10fd815SStefano Zampini       sdim[tnf+f] = sd[f];
1346e10fd815SStefano Zampini     }
1347e10fd815SStefano Zampini     tnf += nf;
1348e10fd815SStefano Zampini   }
1349e10fd815SStefano Zampini   ierr = PetscViewerGLVisSetFields(viewer,tnf,(const char**)fecs,sdim,DMCompositeSampleGLVisFields_Private,(PetscObject*)Ufds,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
1350e10fd815SStefano Zampini   for (i = 0; i < tnf; i++) {
1351e10fd815SStefano Zampini     ierr = PetscFree(fecs[i]);CHKERRQ(ierr);
1352e10fd815SStefano Zampini   }
1353e10fd815SStefano Zampini   ierr = PetscFree3(fecs,sdim,Ufds);CHKERRQ(ierr);
1354e10fd815SStefano Zampini   PetscFunctionReturn(0);
1355e10fd815SStefano Zampini }
1356e10fd815SStefano Zampini 
13577087cfbeSBarry Smith PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
135847c6ae99SBarry Smith {
135947c6ae99SBarry Smith   PetscErrorCode         ierr;
136047c6ae99SBarry Smith   struct DMCompositeLink *next;
136147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmi->data;
136247c6ae99SBarry Smith   DM                     dm;
136347c6ae99SBarry Smith 
136447c6ae99SBarry Smith   PetscFunctionBegin;
136547c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1366ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
1367ce94432eSBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
1368ce94432eSBarry Smith   }
13692ce3a92bSJed Brown   ierr = DMSetUp(dmi);CHKERRQ(ierr);
137047c6ae99SBarry Smith   next = com->next;
137147c6ae99SBarry Smith   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
137247c6ae99SBarry Smith 
137347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
137447c6ae99SBarry Smith   while (next) {
137547c6ae99SBarry Smith     ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
137647c6ae99SBarry Smith     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
137747c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
137847c6ae99SBarry Smith     next = next->next;
137947c6ae99SBarry Smith   }
138047c6ae99SBarry Smith   PetscFunctionReturn(0);
138147c6ae99SBarry Smith }
138247c6ae99SBarry Smith 
138314354c39SJed Brown PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
138414354c39SJed Brown {
138514354c39SJed Brown   PetscErrorCode         ierr;
138614354c39SJed Brown   struct DMCompositeLink *next;
138714354c39SJed Brown   DM_Composite           *com = (DM_Composite*)dmi->data;
138814354c39SJed Brown   DM                     dm;
138914354c39SJed Brown 
139014354c39SJed Brown   PetscFunctionBegin;
139114354c39SJed Brown   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
13922ce3a92bSJed Brown   ierr = DMSetUp(dmi);CHKERRQ(ierr);
13932ee06e3bSJed Brown   if (comm == MPI_COMM_NULL) {
139425296bd5SBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
139525296bd5SBarry Smith   }
139614354c39SJed Brown   next = com->next;
139714354c39SJed Brown   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
139814354c39SJed Brown 
139914354c39SJed Brown   /* loop over packed objects, handling one at at time */
140014354c39SJed Brown   while (next) {
140114354c39SJed Brown     ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr);
140214354c39SJed Brown     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
140314354c39SJed Brown     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
140414354c39SJed Brown     next = next->next;
140514354c39SJed Brown   }
140614354c39SJed Brown   PetscFunctionReturn(0);
140714354c39SJed Brown }
140847c6ae99SBarry Smith 
1409e727c939SJed Brown PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
141047c6ae99SBarry Smith {
141147c6ae99SBarry Smith   PetscErrorCode         ierr;
14129ae5db72SJed Brown   PetscInt               m,n,M,N,nDM,i;
141347c6ae99SBarry Smith   struct DMCompositeLink *nextc;
141447c6ae99SBarry Smith   struct DMCompositeLink *nextf;
141525296bd5SBarry Smith   Vec                    gcoarse,gfine,*vecs;
141647c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
141747c6ae99SBarry Smith   DM_Composite           *comfine   = (DM_Composite*)fine->data;
14189ae5db72SJed Brown   Mat                    *mats;
141947c6ae99SBarry Smith 
142047c6ae99SBarry Smith   PetscFunctionBegin;
142147c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
142247c6ae99SBarry Smith   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
1423f692024eSJed Brown   ierr = DMSetUp(coarse);CHKERRQ(ierr);
1424f692024eSJed Brown   ierr = DMSetUp(fine);CHKERRQ(ierr);
142547c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14269ae5db72SJed Brown   ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
14279ae5db72SJed Brown   ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr);
142847c6ae99SBarry Smith   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
142947c6ae99SBarry Smith   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
143047c6ae99SBarry Smith   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
143147c6ae99SBarry Smith   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
14329ae5db72SJed Brown   ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
14339ae5db72SJed Brown   ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr);
143447c6ae99SBarry Smith 
14359ae5db72SJed Brown   nDM = comfine->nDM;
1436ce94432eSBarry 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);
14371795a4d1SJed Brown   ierr = PetscCalloc1(nDM*nDM,&mats);CHKERRQ(ierr);
143825296bd5SBarry Smith   if (v) {
14391795a4d1SJed Brown     ierr = PetscCalloc1(nDM,&vecs);CHKERRQ(ierr);
144025296bd5SBarry Smith   }
144147c6ae99SBarry Smith 
144247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
14439ae5db72SJed Brown   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
144425296bd5SBarry Smith     if (!v) {
14450298fd71SBarry Smith       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL);CHKERRQ(ierr);
144625296bd5SBarry Smith     } else {
144725296bd5SBarry Smith       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr);
144825296bd5SBarry Smith     }
144947c6ae99SBarry Smith   }
1450ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A);CHKERRQ(ierr);
145125296bd5SBarry Smith   if (v) {
1452ce94432eSBarry Smith     ierr = VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v);CHKERRQ(ierr);
145325296bd5SBarry Smith   }
14549ae5db72SJed Brown   for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);}
14559ae5db72SJed Brown   ierr = PetscFree(mats);CHKERRQ(ierr);
145625296bd5SBarry Smith   if (v) {
145725296bd5SBarry Smith     for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);}
145825296bd5SBarry Smith     ierr = PetscFree(vecs);CHKERRQ(ierr);
145925296bd5SBarry Smith   }
146047c6ae99SBarry Smith   PetscFunctionReturn(0);
146147c6ae99SBarry Smith }
146247c6ae99SBarry Smith 
1463184d77edSJed Brown static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
14641411c6eeSJed Brown {
14651411c6eeSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
14661411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1467f7efa3c7SJed Brown   PetscInt               i;
14681411c6eeSJed Brown   PetscErrorCode         ierr;
14691411c6eeSJed Brown 
14701411c6eeSJed Brown   PetscFunctionBegin;
14711411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
14721411c6eeSJed Brown   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
1473ce94432eSBarry Smith   ierr = ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr);
14749ae5db72SJed Brown   for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
14751411c6eeSJed Brown   ierr = PetscFree(ltogs);CHKERRQ(ierr);
14761411c6eeSJed Brown   PetscFunctionReturn(0);
14771411c6eeSJed Brown }
14781411c6eeSJed Brown 
14791411c6eeSJed Brown 
1480b412c318SBarry Smith PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring)
148147c6ae99SBarry Smith {
148247c6ae99SBarry Smith   PetscErrorCode  ierr;
148347c6ae99SBarry Smith   PetscInt        n,i,cnt;
148447c6ae99SBarry Smith   ISColoringValue *colors;
148547c6ae99SBarry Smith   PetscBool       dense  = PETSC_FALSE;
148647c6ae99SBarry Smith   ISColoringValue maxcol = 0;
148747c6ae99SBarry Smith   DM_Composite    *com   = (DM_Composite*)dm->data;
148847c6ae99SBarry Smith 
148947c6ae99SBarry Smith   PetscFunctionBegin;
149047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
14915bdb020cSBarry Smith   if (ctype == IS_COLORING_LOCAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported");
1492e3247f34SBarry Smith   else if (ctype == IS_COLORING_GLOBAL) {
149347c6ae99SBarry Smith     n = com->n;
1494ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1495785e854fSJed Brown   ierr = PetscMalloc1(n,&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
149647c6ae99SBarry Smith 
1497c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL);CHKERRQ(ierr);
149847c6ae99SBarry Smith   if (dense) {
149947c6ae99SBarry Smith     for (i=0; i<n; i++) {
150047c6ae99SBarry Smith       colors[i] = (ISColoringValue)(com->rstart + i);
150147c6ae99SBarry Smith     }
150247c6ae99SBarry Smith     maxcol = com->N;
150347c6ae99SBarry Smith   } else {
150447c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
150547c6ae99SBarry Smith     PetscMPIInt            rank;
150647c6ae99SBarry Smith 
1507ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
150847c6ae99SBarry Smith     cnt  = 0;
150947c6ae99SBarry Smith     while (next) {
151047c6ae99SBarry Smith       ISColoring lcoloring;
151147c6ae99SBarry Smith 
1512b412c318SBarry Smith       ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring);CHKERRQ(ierr);
151347c6ae99SBarry Smith       for (i=0; i<lcoloring->N; i++) {
151447c6ae99SBarry Smith         colors[cnt++] = maxcol + lcoloring->colors[i];
151547c6ae99SBarry Smith       }
151647c6ae99SBarry Smith       maxcol += lcoloring->n;
1517fcfd50ebSBarry Smith       ierr    = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
151847c6ae99SBarry Smith       next    = next->next;
151947c6ae99SBarry Smith     }
152047c6ae99SBarry Smith   }
1521aaf3ff59SMatthew G. Knepley   ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,PETSC_OWN_POINTER,coloring);CHKERRQ(ierr);
152247c6ae99SBarry Smith   PetscFunctionReturn(0);
152347c6ae99SBarry Smith }
152447c6ae99SBarry Smith 
15257087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
152647c6ae99SBarry Smith {
152747c6ae99SBarry Smith   PetscErrorCode         ierr;
152847c6ae99SBarry Smith   struct DMCompositeLink *next;
152947c6ae99SBarry Smith   PetscScalar            *garray,*larray;
153047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
153147c6ae99SBarry Smith 
153247c6ae99SBarry Smith   PetscFunctionBegin;
153347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
153447c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
153539d35262SVincent Le Chenadec 
153647c6ae99SBarry Smith   if (!com->setup) {
1537d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
153847c6ae99SBarry Smith   }
153939d35262SVincent Le Chenadec 
154047c6ae99SBarry Smith   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
154147c6ae99SBarry Smith   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
154247c6ae99SBarry Smith 
154347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
154439d35262SVincent Le Chenadec   next = com->next;
154547c6ae99SBarry Smith   while (next) {
154647c6ae99SBarry Smith     Vec      local,global;
154747c6ae99SBarry Smith     PetscInt N;
154847c6ae99SBarry Smith 
154947c6ae99SBarry Smith     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
155047c6ae99SBarry Smith     ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
155147c6ae99SBarry Smith     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
155247c6ae99SBarry Smith     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
155347c6ae99SBarry Smith     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
155447c6ae99SBarry Smith     ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
155547c6ae99SBarry Smith     ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
155647c6ae99SBarry Smith     ierr = VecResetArray(global);CHKERRQ(ierr);
155747c6ae99SBarry Smith     ierr = VecResetArray(local);CHKERRQ(ierr);
155847c6ae99SBarry Smith     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
155939d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr);
156039d35262SVincent Le Chenadec 
156106ebdd98SJed Brown     larray += next->nlocal;
156239d35262SVincent Le Chenadec     garray += next->n;
156347c6ae99SBarry Smith     next    = next->next;
156447c6ae99SBarry Smith   }
156547c6ae99SBarry Smith 
15660298fd71SBarry Smith   ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr);
15670298fd71SBarry Smith   ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr);
156847c6ae99SBarry Smith   PetscFunctionReturn(0);
156947c6ae99SBarry Smith }
157047c6ae99SBarry Smith 
15717087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
15720c010503SBarry Smith {
15730c010503SBarry Smith   PetscFunctionBegin;
157439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
157539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
157639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,4);
157739d35262SVincent Le Chenadec   PetscFunctionReturn(0);
157839d35262SVincent Le Chenadec }
157939d35262SVincent Le Chenadec 
158039d35262SVincent Le Chenadec PetscErrorCode  DMLocalToGlobalBegin_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
158139d35262SVincent Le Chenadec {
158239d35262SVincent Le Chenadec   PetscErrorCode         ierr;
158339d35262SVincent Le Chenadec   struct DMCompositeLink *next;
158439d35262SVincent Le Chenadec   PetscScalar            *larray,*garray;
158539d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite*)dm->data;
158639d35262SVincent Le Chenadec 
158739d35262SVincent Le Chenadec   PetscFunctionBegin;
158839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
158939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
159039d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
159139d35262SVincent Le Chenadec 
159239d35262SVincent Le Chenadec   if (!com->setup) {
159339d35262SVincent Le Chenadec     ierr = DMSetUp(dm);CHKERRQ(ierr);
159439d35262SVincent Le Chenadec   }
159539d35262SVincent Le Chenadec 
159639d35262SVincent Le Chenadec   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
159739d35262SVincent Le Chenadec   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
159839d35262SVincent Le Chenadec 
159939d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
160039d35262SVincent Le Chenadec   next = com->next;
160139d35262SVincent Le Chenadec   while (next) {
160239d35262SVincent Le Chenadec     Vec      global,local;
160339d35262SVincent Le Chenadec 
160439d35262SVincent Le Chenadec     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
160539d35262SVincent Le Chenadec     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
160639d35262SVincent Le Chenadec     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
160739d35262SVincent Le Chenadec     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
160839d35262SVincent Le Chenadec     ierr = DMLocalToGlobalBegin(next->dm,local,mode,global);CHKERRQ(ierr);
160939d35262SVincent Le Chenadec     ierr = DMLocalToGlobalEnd(next->dm,local,mode,global);CHKERRQ(ierr);
161039d35262SVincent Le Chenadec     ierr = VecResetArray(local);CHKERRQ(ierr);
161139d35262SVincent Le Chenadec     ierr = VecResetArray(global);CHKERRQ(ierr);
161239d35262SVincent Le Chenadec     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
161339d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr);
161439d35262SVincent Le Chenadec 
161539d35262SVincent Le Chenadec     garray += next->n;
161639d35262SVincent Le Chenadec     larray += next->nlocal;
161739d35262SVincent Le Chenadec     next    = next->next;
161839d35262SVincent Le Chenadec   }
161939d35262SVincent Le Chenadec 
162039d35262SVincent Le Chenadec   ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr);
162139d35262SVincent Le Chenadec   ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr);
162239d35262SVincent Le Chenadec 
162339d35262SVincent Le Chenadec   PetscFunctionReturn(0);
162439d35262SVincent Le Chenadec }
162539d35262SVincent Le Chenadec 
162639d35262SVincent Le Chenadec PetscErrorCode  DMLocalToGlobalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
162739d35262SVincent Le Chenadec {
162839d35262SVincent Le Chenadec   PetscFunctionBegin;
162939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
163039d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
163139d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
163239d35262SVincent Le Chenadec   PetscFunctionReturn(0);
163339d35262SVincent Le Chenadec }
163439d35262SVincent Le Chenadec 
163539d35262SVincent Le Chenadec PetscErrorCode  DMLocalToLocalBegin_Composite(DM dm,Vec vec1,InsertMode mode,Vec vec2)
163639d35262SVincent Le Chenadec {
163739d35262SVincent Le Chenadec   PetscErrorCode         ierr;
163839d35262SVincent Le Chenadec   struct DMCompositeLink *next;
163939d35262SVincent Le Chenadec   PetscScalar            *array1,*array2;
164039d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite*)dm->data;
164139d35262SVincent Le Chenadec 
164239d35262SVincent Le Chenadec   PetscFunctionBegin;
164339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
164439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec1,VEC_CLASSID,2);
164539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec2,VEC_CLASSID,4);
164639d35262SVincent Le Chenadec 
164739d35262SVincent Le Chenadec   if (!com->setup) {
164839d35262SVincent Le Chenadec     ierr = DMSetUp(dm);CHKERRQ(ierr);
164939d35262SVincent Le Chenadec   }
165039d35262SVincent Le Chenadec 
165139d35262SVincent Le Chenadec   ierr = VecGetArray(vec1,&array1);CHKERRQ(ierr);
165239d35262SVincent Le Chenadec   ierr = VecGetArray(vec2,&array2);CHKERRQ(ierr);
165339d35262SVincent Le Chenadec 
165439d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
165539d35262SVincent Le Chenadec   next = com->next;
165639d35262SVincent Le Chenadec   while (next) {
165739d35262SVincent Le Chenadec     Vec      local1,local2;
165839d35262SVincent Le Chenadec 
165939d35262SVincent Le Chenadec     ierr = DMGetLocalVector(next->dm,&local1);CHKERRQ(ierr);
166039d35262SVincent Le Chenadec     ierr = VecPlaceArray(local1,array1);CHKERRQ(ierr);
166139d35262SVincent Le Chenadec     ierr = DMGetLocalVector(next->dm,&local2);CHKERRQ(ierr);
166239d35262SVincent Le Chenadec     ierr = VecPlaceArray(local2,array2);CHKERRQ(ierr);
166339d35262SVincent Le Chenadec     ierr = DMLocalToLocalBegin(next->dm,local1,mode,local2);CHKERRQ(ierr);
166439d35262SVincent Le Chenadec     ierr = DMLocalToLocalEnd(next->dm,local1,mode,local2);CHKERRQ(ierr);
166539d35262SVincent Le Chenadec     ierr = VecResetArray(local2);CHKERRQ(ierr);
166639d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local2);CHKERRQ(ierr);
166739d35262SVincent Le Chenadec     ierr = VecResetArray(local1);CHKERRQ(ierr);
166839d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local1);CHKERRQ(ierr);
166939d35262SVincent Le Chenadec 
167039d35262SVincent Le Chenadec     array1 += next->nlocal;
167139d35262SVincent Le Chenadec     array2 += next->nlocal;
167239d35262SVincent Le Chenadec     next    = next->next;
167339d35262SVincent Le Chenadec   }
167439d35262SVincent Le Chenadec 
167539d35262SVincent Le Chenadec   ierr = VecRestoreArray(vec1,NULL);CHKERRQ(ierr);
167639d35262SVincent Le Chenadec   ierr = VecRestoreArray(vec2,NULL);CHKERRQ(ierr);
167739d35262SVincent Le Chenadec 
167839d35262SVincent Le Chenadec   PetscFunctionReturn(0);
167939d35262SVincent Le Chenadec }
168039d35262SVincent Le Chenadec 
168139d35262SVincent Le Chenadec PetscErrorCode  DMLocalToLocalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
168239d35262SVincent Le Chenadec {
168339d35262SVincent Le Chenadec   PetscFunctionBegin;
168439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
168539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
168639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
16870c010503SBarry Smith   PetscFunctionReturn(0);
16880c010503SBarry Smith }
168947c6ae99SBarry Smith 
16906ae3a549SBarry Smith /*MC
16916ae3a549SBarry Smith    DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs
16926ae3a549SBarry Smith 
16936ae3a549SBarry Smith   Level: intermediate
16946ae3a549SBarry Smith 
16951abcec8cSBarry Smith .seealso: DMType, DM, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate()
16966ae3a549SBarry Smith M*/
16976ae3a549SBarry Smith 
16986ae3a549SBarry Smith 
16998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1700a4121054SBarry Smith {
1701a4121054SBarry Smith   PetscErrorCode ierr;
1702a4121054SBarry Smith   DM_Composite   *com;
1703a4121054SBarry Smith 
1704a4121054SBarry Smith   PetscFunctionBegin;
1705b00a9115SJed Brown   ierr          = PetscNewLog(p,&com);CHKERRQ(ierr);
1706a4121054SBarry Smith   p->data       = com;
1707a4121054SBarry Smith   ierr          = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1708a4121054SBarry Smith   com->n        = 0;
17097ac2b803SAlex Fikl   com->nghost   = 0;
17100298fd71SBarry Smith   com->next     = NULL;
1711a4121054SBarry Smith   com->nDM      = 0;
1712a4121054SBarry Smith 
1713a4121054SBarry Smith   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1714a4121054SBarry Smith   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1715184d77edSJed Brown   p->ops->getlocaltoglobalmapping         = DMGetLocalToGlobalMapping_Composite;
17164d343eeaSMatthew G Knepley   p->ops->createfieldis                   = DMCreateFieldIS_Composite;
171716621825SDmitry Karpeev   p->ops->createfielddecomposition        = DMCreateFieldDecomposition_Composite;
1718a4121054SBarry Smith   p->ops->refine                          = DMRefine_Composite;
171914354c39SJed Brown   p->ops->coarsen                         = DMCoarsen_Composite;
172025296bd5SBarry Smith   p->ops->createinterpolation             = DMCreateInterpolation_Composite;
172125296bd5SBarry Smith   p->ops->creatematrix                    = DMCreateMatrix_Composite;
1722e727c939SJed Brown   p->ops->getcoloring                     = DMCreateColoring_Composite;
1723a4121054SBarry Smith   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1724a4121054SBarry Smith   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
172539d35262SVincent Le Chenadec   p->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Composite;
172639d35262SVincent Le Chenadec   p->ops->localtoglobalend                = DMLocalToGlobalEnd_Composite;
172739d35262SVincent Le Chenadec   p->ops->localtolocalbegin               = DMLocalToLocalBegin_Composite;
172839d35262SVincent Le Chenadec   p->ops->localtolocalend                 = DMLocalToLocalEnd_Composite;
1729a4121054SBarry Smith   p->ops->destroy                         = DMDestroy_Composite;
1730a4121054SBarry Smith   p->ops->view                            = DMView_Composite;
1731a4121054SBarry Smith   p->ops->setup                           = DMSetUp_Composite;
1732e10fd815SStefano Zampini 
1733e10fd815SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)p,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Composite);CHKERRQ(ierr);
1734a4121054SBarry Smith   PetscFunctionReturn(0);
1735a4121054SBarry Smith }
1736a4121054SBarry Smith 
1737*1fd49c25SBarry Smith /*@
17380c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
17390c010503SBarry Smith       vectors made up of several subvectors.
17400c010503SBarry Smith 
17410c010503SBarry Smith     Collective on MPI_Comm
174247c6ae99SBarry Smith 
174347c6ae99SBarry Smith     Input Parameter:
17440c010503SBarry Smith .   comm - the processors that will share the global vector
17450c010503SBarry Smith 
17460c010503SBarry Smith     Output Parameters:
17470c010503SBarry Smith .   packer - the packer object
174847c6ae99SBarry Smith 
174947c6ae99SBarry Smith     Level: advanced
175047c6ae99SBarry Smith 
17511abcec8cSBarry Smith .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(), DMCOMPOSITE,DMCreate()
17526eb61c8cSJed Brown          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
175347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
175447c6ae99SBarry Smith 
175547c6ae99SBarry Smith @*/
17567087cfbeSBarry Smith PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
175747c6ae99SBarry Smith {
17580c010503SBarry Smith   PetscErrorCode ierr;
17590c010503SBarry Smith 
176047c6ae99SBarry Smith   PetscFunctionBegin;
17610c010503SBarry Smith   PetscValidPointer(packer,2);
1762a4121054SBarry Smith   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1763a4121054SBarry Smith   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
176447c6ae99SBarry Smith   PetscFunctionReturn(0);
176547c6ae99SBarry Smith }
1766