xref: /petsc/src/dm/impls/composite/pack.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
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 
20f5f57ec0SBarry Smith     Not available from Fortran
21f5f57ec0SBarry Smith 
22*95452b02SPatrick Sanan     Notes:
23*95452b02SPatrick Sanan     See DMSetApplicationContext() and DMGetApplicationContext() for how to get user information into
2447c6ae99SBarry Smith         this routine
2547c6ae99SBarry Smith 
2647c6ae99SBarry Smith @*/
277087cfbeSBarry Smith PetscErrorCode  DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt))
2847c6ae99SBarry Smith {
2947c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
3047c6ae99SBarry Smith 
3147c6ae99SBarry Smith   PetscFunctionBegin;
3247c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
3347c6ae99SBarry Smith   PetscFunctionReturn(0);
3447c6ae99SBarry Smith }
3547c6ae99SBarry Smith 
366bf464f9SBarry Smith PetscErrorCode  DMDestroy_Composite(DM dm)
3747c6ae99SBarry Smith {
3847c6ae99SBarry Smith   PetscErrorCode         ierr;
3947c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
4047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
4147c6ae99SBarry Smith 
4247c6ae99SBarry Smith   PetscFunctionBegin;
4347c6ae99SBarry Smith   next = com->next;
4447c6ae99SBarry Smith   while (next) {
4547c6ae99SBarry Smith     prev = next;
4647c6ae99SBarry Smith     next = next->next;
47fcfd50ebSBarry Smith     ierr = DMDestroy(&prev->dm);CHKERRQ(ierr);
4847c6ae99SBarry Smith     ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
4947c6ae99SBarry Smith     ierr = PetscFree(prev);CHKERRQ(ierr);
5047c6ae99SBarry Smith   }
51e10fd815SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL);CHKERRQ(ierr);
52435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
53435a35e8SMatthew G Knepley   ierr = PetscFree(com);CHKERRQ(ierr);
5447c6ae99SBarry Smith   PetscFunctionReturn(0);
5547c6ae99SBarry Smith }
5647c6ae99SBarry Smith 
577087cfbeSBarry Smith PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
5847c6ae99SBarry Smith {
5947c6ae99SBarry Smith   PetscErrorCode ierr;
6047c6ae99SBarry Smith   PetscBool      iascii;
6147c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
6247c6ae99SBarry Smith 
6347c6ae99SBarry Smith   PetscFunctionBegin;
64251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6547c6ae99SBarry Smith   if (iascii) {
6647c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
6747c6ae99SBarry Smith     PetscInt               i;
6847c6ae99SBarry Smith 
6947c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix");CHKERRQ(ierr);
709ae5db72SJed Brown     ierr = PetscViewerASCIIPrintf(v,"  contains %D DMs\n",com->nDM);CHKERRQ(ierr);
7147c6ae99SBarry Smith     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
7247c6ae99SBarry Smith     for (i=0; lnk; lnk=lnk->next,i++) {
739ae5db72SJed Brown       ierr = PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
7447c6ae99SBarry Smith       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
7547c6ae99SBarry Smith       ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
7647c6ae99SBarry Smith       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
7747c6ae99SBarry Smith     }
7847c6ae99SBarry Smith     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
7947c6ae99SBarry Smith   }
8047c6ae99SBarry Smith   PetscFunctionReturn(0);
8147c6ae99SBarry Smith }
8247c6ae99SBarry Smith 
8347c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
847087cfbeSBarry Smith PetscErrorCode  DMSetUp_Composite(DM dm)
8547c6ae99SBarry Smith {
8647c6ae99SBarry Smith   PetscErrorCode         ierr;
8747c6ae99SBarry Smith   PetscInt               nprev = 0;
8847c6ae99SBarry Smith   PetscMPIInt            rank,size;
8947c6ae99SBarry Smith   DM_Composite           *com  = (DM_Composite*)dm->data;
9047c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
9147c6ae99SBarry Smith   PetscLayout            map;
9247c6ae99SBarry Smith 
9347c6ae99SBarry Smith   PetscFunctionBegin;
94ce94432eSBarry Smith   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
95ce94432eSBarry Smith   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm),&map);CHKERRQ(ierr);
9647c6ae99SBarry Smith   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
9747c6ae99SBarry Smith   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
9847c6ae99SBarry Smith   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
9947c6ae99SBarry Smith   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
10047c6ae99SBarry Smith   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
1010298fd71SBarry Smith   ierr = PetscLayoutGetRange(map,&com->rstart,NULL);CHKERRQ(ierr);
102fcfd50ebSBarry Smith   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
10347c6ae99SBarry Smith 
1049ae5db72SJed Brown   /* now set the rstart for each linked vector */
105ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
106ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
10747c6ae99SBarry Smith   while (next) {
10847c6ae99SBarry Smith     next->rstart  = nprev;
10906ebdd98SJed Brown     nprev        += next->n;
11047c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
111785e854fSJed Brown     ierr          = PetscMalloc1(size,&next->grstarts);CHKERRQ(ierr);
112ce94432eSBarry Smith     ierr          = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
11347c6ae99SBarry Smith     next          = next->next;
11447c6ae99SBarry Smith   }
11547c6ae99SBarry Smith   com->setup = PETSC_TRUE;
11647c6ae99SBarry Smith   PetscFunctionReturn(0);
11747c6ae99SBarry Smith }
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
12047c6ae99SBarry Smith 
12173e31fe2SJed Brown /*@
12247c6ae99SBarry Smith     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
12347c6ae99SBarry Smith        representation.
12447c6ae99SBarry Smith 
12547c6ae99SBarry Smith     Not Collective
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith     Input Parameter:
12847c6ae99SBarry Smith .    dm - the packer object
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith     Output Parameter:
13147c6ae99SBarry Smith .     nDM - the number of DMs
13247c6ae99SBarry Smith 
13347c6ae99SBarry Smith     Level: beginner
13447c6ae99SBarry Smith 
13547c6ae99SBarry Smith @*/
1367087cfbeSBarry Smith PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
13747c6ae99SBarry Smith {
13847c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
1395fd66863SKarl Rupp 
14047c6ae99SBarry Smith   PetscFunctionBegin;
14147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
14247c6ae99SBarry Smith   *nDM = com->nDM;
14347c6ae99SBarry Smith   PetscFunctionReturn(0);
14447c6ae99SBarry Smith }
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith 
14747c6ae99SBarry Smith /*@C
14847c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
14947c6ae99SBarry Smith        representation.
15047c6ae99SBarry Smith 
15147c6ae99SBarry Smith     Collective on DMComposite
15247c6ae99SBarry Smith 
1539ae5db72SJed Brown     Input Parameters:
15447c6ae99SBarry Smith +    dm - the packer object
1559ae5db72SJed Brown -    gvec - the global vector
1569ae5db72SJed Brown 
1579ae5db72SJed Brown     Output Parameters:
1580298fd71SBarry Smith .    Vec* ... - the packed parallel vectors, NULL for those that are not needed
15947c6ae99SBarry Smith 
160*95452b02SPatrick Sanan     Notes:
161*95452b02SPatrick Sanan     Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
16247c6ae99SBarry Smith 
163f73e5cebSJed Brown     Fortran Notes:
164f73e5cebSJed Brown 
165f73e5cebSJed Brown     Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
166f73e5cebSJed Brown     or use the alternative interface DMCompositeGetAccessArray().
167f73e5cebSJed Brown 
16847c6ae99SBarry Smith     Level: advanced
16947c6ae99SBarry Smith 
170f73e5cebSJed Brown .seealso: DMCompositeGetEntries(), DMCompositeScatter()
17147c6ae99SBarry Smith @*/
1727087cfbeSBarry Smith PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
17347c6ae99SBarry Smith {
17447c6ae99SBarry Smith   va_list                Argp;
17547c6ae99SBarry Smith   PetscErrorCode         ierr;
17647c6ae99SBarry Smith   struct DMCompositeLink *next;
17747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
1785edff71fSBarry Smith   PetscInt               readonly;
17947c6ae99SBarry Smith 
18047c6ae99SBarry Smith   PetscFunctionBegin;
18147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
18247c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
18347c6ae99SBarry Smith   next = com->next;
18447c6ae99SBarry Smith   if (!com->setup) {
185d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
18647c6ae99SBarry Smith   }
18747c6ae99SBarry Smith 
1885edff71fSBarry Smith   ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr);
18947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
19047c6ae99SBarry Smith   va_start(Argp,gvec);
19147c6ae99SBarry Smith   while (next) {
19247c6ae99SBarry Smith     Vec *vec;
19347c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
1949ae5db72SJed Brown     if (vec) {
1959ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,vec);CHKERRQ(ierr);
1965edff71fSBarry Smith       if (readonly) {
1975edff71fSBarry Smith         const PetscScalar *array;
1985edff71fSBarry Smith         ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
1995edff71fSBarry Smith         ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
2005edff71fSBarry Smith         ierr = VecLockPush(*vec);CHKERRQ(ierr);
2015edff71fSBarry Smith         ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
2025edff71fSBarry Smith       } else {
2035edff71fSBarry Smith         PetscScalar *array;
2049ae5db72SJed Brown         ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
2059ae5db72SJed Brown         ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
2069ae5db72SJed Brown         ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
20747c6ae99SBarry Smith       }
2085edff71fSBarry Smith     }
20947c6ae99SBarry Smith     next = next->next;
21047c6ae99SBarry Smith   }
21147c6ae99SBarry Smith   va_end(Argp);
21247c6ae99SBarry Smith   PetscFunctionReturn(0);
21347c6ae99SBarry Smith }
21447c6ae99SBarry Smith 
215f73e5cebSJed Brown /*@C
216f73e5cebSJed Brown     DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
217f73e5cebSJed Brown        representation.
218f73e5cebSJed Brown 
219f73e5cebSJed Brown     Collective on DMComposite
220f73e5cebSJed Brown 
221f73e5cebSJed Brown     Input Parameters:
222f73e5cebSJed Brown +    dm - the packer object
223f73e5cebSJed Brown .    pvec - packed vector
224f73e5cebSJed Brown .    nwanted - number of vectors wanted
2250298fd71SBarry Smith -    wanted - sorted array of vectors wanted, or NULL to get all vectors
226f73e5cebSJed Brown 
227f73e5cebSJed Brown     Output Parameters:
228f73e5cebSJed Brown .    vecs - array of requested global vectors (must be allocated)
229f73e5cebSJed Brown 
230*95452b02SPatrick Sanan     Notes:
231*95452b02SPatrick Sanan     Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them
232f73e5cebSJed Brown 
233f73e5cebSJed Brown     Level: advanced
234f73e5cebSJed Brown 
235f73e5cebSJed Brown .seealso: DMCompositeGetAccess(), DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
236f73e5cebSJed Brown @*/
237f73e5cebSJed Brown PetscErrorCode  DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
238f73e5cebSJed Brown {
239f73e5cebSJed Brown   PetscErrorCode         ierr;
240f73e5cebSJed Brown   struct DMCompositeLink *link;
241f73e5cebSJed Brown   PetscInt               i,wnum;
242f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
243bee642f7SBarry Smith   PetscInt               readonly;
244f73e5cebSJed Brown 
245f73e5cebSJed Brown   PetscFunctionBegin;
246f73e5cebSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
247f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
248f73e5cebSJed Brown   if (!com->setup) {
249f73e5cebSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
250f73e5cebSJed Brown   }
251f73e5cebSJed Brown 
252bee642f7SBarry Smith   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
253f73e5cebSJed Brown   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
254f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
255f73e5cebSJed Brown       Vec v;
256f73e5cebSJed Brown       ierr = DMGetGlobalVector(link->dm,&v);CHKERRQ(ierr);
257bee642f7SBarry Smith       if (readonly) {
258bee642f7SBarry Smith         const PetscScalar *array;
259bee642f7SBarry Smith         ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr);
260bee642f7SBarry Smith         ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
261bee642f7SBarry Smith         ierr = VecLockPush(v);CHKERRQ(ierr);
262bee642f7SBarry Smith         ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr);
263bee642f7SBarry Smith       } else {
264bee642f7SBarry Smith         PetscScalar *array;
265f73e5cebSJed Brown         ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
266f73e5cebSJed Brown         ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
267f73e5cebSJed Brown         ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
268bee642f7SBarry Smith       }
269f73e5cebSJed Brown       vecs[wnum++] = v;
270f73e5cebSJed Brown     }
271f73e5cebSJed Brown   }
272f73e5cebSJed Brown   PetscFunctionReturn(0);
273f73e5cebSJed Brown }
274f73e5cebSJed Brown 
2757ac2b803SAlex Fikl /*@C
2767ac2b803SAlex Fikl     DMCompositeGetLocalAccessArray - Allows one to access the individual
2777ac2b803SAlex Fikl     packed vectors in their local representation.
2787ac2b803SAlex Fikl 
2797ac2b803SAlex Fikl     Collective on DMComposite.
2807ac2b803SAlex Fikl 
2817ac2b803SAlex Fikl     Input Parameters:
2827ac2b803SAlex Fikl +    dm - the packer object
2837ac2b803SAlex Fikl .    pvec - packed vector
2847ac2b803SAlex Fikl .    nwanted - number of vectors wanted
2857ac2b803SAlex Fikl -    wanted - sorted array of vectors wanted, or NULL to get all vectors
2867ac2b803SAlex Fikl 
2877ac2b803SAlex Fikl     Output Parameters:
2887ac2b803SAlex Fikl .    vecs - array of requested local vectors (must be allocated)
2897ac2b803SAlex Fikl 
290*95452b02SPatrick Sanan     Notes:
291*95452b02SPatrick Sanan     Use DMCompositeRestoreLocalAccessArray() to return the vectors
2927ac2b803SAlex Fikl     when you no longer need them.
2937ac2b803SAlex Fikl 
2947ac2b803SAlex Fikl     Level: advanced
2957ac2b803SAlex Fikl 
2967ac2b803SAlex Fikl .seealso: DMCompositeRestoreLocalAccessArray(), DMCompositeGetAccess(),
2977ac2b803SAlex Fikl DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
2987ac2b803SAlex Fikl @*/
2997ac2b803SAlex Fikl PetscErrorCode  DMCompositeGetLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
3007ac2b803SAlex Fikl {
3017ac2b803SAlex Fikl   PetscErrorCode         ierr;
3027ac2b803SAlex Fikl   struct DMCompositeLink *link;
3037ac2b803SAlex Fikl   PetscInt               i,wnum;
3047ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite*)dm->data;
3057ac2b803SAlex Fikl   PetscInt               readonly;
3067ac2b803SAlex Fikl   PetscInt               nlocal = 0;
3077ac2b803SAlex Fikl 
3087ac2b803SAlex Fikl   PetscFunctionBegin;
3097ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3107ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
3117ac2b803SAlex Fikl   if (!com->setup) {
3127ac2b803SAlex Fikl     ierr = DMSetUp(dm);CHKERRQ(ierr);
3137ac2b803SAlex Fikl   }
3147ac2b803SAlex Fikl 
3157ac2b803SAlex Fikl   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
3167ac2b803SAlex Fikl   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
3177ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
3187ac2b803SAlex Fikl       Vec v;
3197ac2b803SAlex Fikl       ierr = DMGetLocalVector(link->dm,&v);CHKERRQ(ierr);
3207ac2b803SAlex Fikl       if (readonly) {
3217ac2b803SAlex Fikl         const PetscScalar *array;
3227ac2b803SAlex Fikl         ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr);
3237ac2b803SAlex Fikl         ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr);
3247ac2b803SAlex Fikl         ierr = VecLockPush(v);CHKERRQ(ierr);
3257ac2b803SAlex Fikl         ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr);
3267ac2b803SAlex Fikl       } else {
3277ac2b803SAlex Fikl         PetscScalar *array;
3287ac2b803SAlex Fikl         ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
3297ac2b803SAlex Fikl         ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr);
3307ac2b803SAlex Fikl         ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
3317ac2b803SAlex Fikl       }
3327ac2b803SAlex Fikl       vecs[wnum++] = v;
3337ac2b803SAlex Fikl     }
3347ac2b803SAlex Fikl 
3357ac2b803SAlex Fikl     nlocal += link->nlocal;
3367ac2b803SAlex Fikl   }
3377ac2b803SAlex Fikl 
3387ac2b803SAlex Fikl   PetscFunctionReturn(0);
3397ac2b803SAlex Fikl }
3407ac2b803SAlex Fikl 
34147c6ae99SBarry Smith /*@C
342aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
34347c6ae99SBarry Smith        representation.
34447c6ae99SBarry Smith 
34547c6ae99SBarry Smith     Collective on DMComposite
34647c6ae99SBarry Smith 
3479ae5db72SJed Brown     Input Parameters:
34847c6ae99SBarry Smith +    dm - the packer object
34947c6ae99SBarry Smith .    gvec - the global vector
3500298fd71SBarry Smith -    Vec* ... - the individual parallel vectors, NULL for those that are not needed
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith     Level: advanced
35347c6ae99SBarry Smith 
3549ae5db72SJed Brown .seealso  DMCompositeAddDM(), DMCreateGlobalVector(),
3556eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
356aa219208SBarry Smith          DMCompositeRestoreAccess(), DMCompositeGetAccess()
35747c6ae99SBarry Smith 
35847c6ae99SBarry Smith @*/
3597087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
36047c6ae99SBarry Smith {
36147c6ae99SBarry Smith   va_list                Argp;
36247c6ae99SBarry Smith   PetscErrorCode         ierr;
36347c6ae99SBarry Smith   struct DMCompositeLink *next;
36447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
3655edff71fSBarry Smith   PetscInt               readonly;
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith   PetscFunctionBegin;
36847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36947c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
37047c6ae99SBarry Smith   next = com->next;
37147c6ae99SBarry Smith   if (!com->setup) {
372d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
37347c6ae99SBarry Smith   }
37447c6ae99SBarry Smith 
3755edff71fSBarry Smith   ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr);
37647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
37747c6ae99SBarry Smith   va_start(Argp,gvec);
37847c6ae99SBarry Smith   while (next) {
37947c6ae99SBarry Smith     Vec *vec;
38047c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
3819ae5db72SJed Brown     if (vec) {
3829ae5db72SJed Brown       ierr = VecResetArray(*vec);CHKERRQ(ierr);
3835edff71fSBarry Smith       if (readonly) {
3845edff71fSBarry Smith         ierr = VecLockPop(*vec);CHKERRQ(ierr);
3855edff71fSBarry Smith       }
386bee642f7SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,vec);CHKERRQ(ierr);
38747c6ae99SBarry Smith     }
38847c6ae99SBarry Smith     next = next->next;
38947c6ae99SBarry Smith   }
39047c6ae99SBarry Smith   va_end(Argp);
39147c6ae99SBarry Smith   PetscFunctionReturn(0);
39247c6ae99SBarry Smith }
39347c6ae99SBarry Smith 
394f73e5cebSJed Brown /*@C
395f73e5cebSJed Brown     DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray()
396f73e5cebSJed Brown 
397f73e5cebSJed Brown     Collective on DMComposite
398f73e5cebSJed Brown 
399f73e5cebSJed Brown     Input Parameters:
400f73e5cebSJed Brown +    dm - the packer object
401f73e5cebSJed Brown .    pvec - packed vector
402f73e5cebSJed Brown .    nwanted - number of vectors wanted
4030298fd71SBarry Smith .    wanted - sorted array of vectors wanted, or NULL to get all vectors
404f73e5cebSJed Brown -    vecs - array of global vectors to return
405f73e5cebSJed Brown 
406f73e5cebSJed Brown     Level: advanced
407f73e5cebSJed Brown 
408f73e5cebSJed Brown .seealso: DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), DMCompositeScatter(), DMCompositeGather()
409f73e5cebSJed Brown @*/
410f73e5cebSJed Brown PetscErrorCode  DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
411f73e5cebSJed Brown {
412f73e5cebSJed Brown   PetscErrorCode         ierr;
413f73e5cebSJed Brown   struct DMCompositeLink *link;
414f73e5cebSJed Brown   PetscInt               i,wnum;
415f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
416bee642f7SBarry Smith   PetscInt               readonly;
417f73e5cebSJed Brown 
418f73e5cebSJed Brown   PetscFunctionBegin;
419f73e5cebSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
420f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
421f73e5cebSJed Brown   if (!com->setup) {
422f73e5cebSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
423f73e5cebSJed Brown   }
424f73e5cebSJed Brown 
425bee642f7SBarry Smith   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
426f73e5cebSJed Brown   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
427f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
428f73e5cebSJed Brown       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
429bee642f7SBarry Smith       if (readonly) {
430bee642f7SBarry Smith         ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr);
431bee642f7SBarry Smith       }
432f73e5cebSJed Brown       ierr = DMRestoreGlobalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
433f73e5cebSJed Brown       wnum++;
434f73e5cebSJed Brown     }
435f73e5cebSJed Brown   }
436f73e5cebSJed Brown   PetscFunctionReturn(0);
437f73e5cebSJed Brown }
438f73e5cebSJed Brown 
4397ac2b803SAlex Fikl /*@C
4407ac2b803SAlex Fikl     DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray().
4417ac2b803SAlex Fikl 
4427ac2b803SAlex Fikl     Collective on DMComposite.
4437ac2b803SAlex Fikl 
4447ac2b803SAlex Fikl     Input Parameters:
4457ac2b803SAlex Fikl +    dm - the packer object
4467ac2b803SAlex Fikl .    pvec - packed vector
4477ac2b803SAlex Fikl .    nwanted - number of vectors wanted
4487ac2b803SAlex Fikl .    wanted - sorted array of vectors wanted, or NULL to restore all vectors
4497ac2b803SAlex Fikl -    vecs - array of local vectors to return
4507ac2b803SAlex Fikl 
4517ac2b803SAlex Fikl     Level: advanced
4527ac2b803SAlex Fikl 
4537ac2b803SAlex Fikl     Notes:
4547ac2b803SAlex Fikl     nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray()
4557ac2b803SAlex Fikl     otherwise the call will fail.
4567ac2b803SAlex Fikl 
4577ac2b803SAlex Fikl .seealso: DMCompositeGetLocalAccessArray(), DMCompositeRestoreAccessArray(),
4587ac2b803SAlex Fikl DMCompositeRestoreAccess(), DMCompositeRestoreEntries(),
4597ac2b803SAlex Fikl DMCompositeScatter(), DMCompositeGather()
4607ac2b803SAlex Fikl @*/
4617ac2b803SAlex Fikl PetscErrorCode  DMCompositeRestoreLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
4627ac2b803SAlex Fikl {
4637ac2b803SAlex Fikl   PetscErrorCode         ierr;
4647ac2b803SAlex Fikl   struct DMCompositeLink *link;
4657ac2b803SAlex Fikl   PetscInt               i,wnum;
4667ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite*)dm->data;
4677ac2b803SAlex Fikl   PetscInt               readonly;
4687ac2b803SAlex Fikl 
4697ac2b803SAlex Fikl   PetscFunctionBegin;
4707ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4717ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
4727ac2b803SAlex Fikl   if (!com->setup) {
4737ac2b803SAlex Fikl     ierr = DMSetUp(dm);CHKERRQ(ierr);
4747ac2b803SAlex Fikl   }
4757ac2b803SAlex Fikl 
4767ac2b803SAlex Fikl   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
4777ac2b803SAlex Fikl   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
4787ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
4797ac2b803SAlex Fikl       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
4807ac2b803SAlex Fikl       if (readonly) {
4817ac2b803SAlex Fikl         ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr);
4827ac2b803SAlex Fikl       }
4837ac2b803SAlex Fikl       ierr = DMRestoreLocalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
4847ac2b803SAlex Fikl       wnum++;
4857ac2b803SAlex Fikl     }
4867ac2b803SAlex Fikl   }
4877ac2b803SAlex Fikl   PetscFunctionReturn(0);
4887ac2b803SAlex Fikl }
4897ac2b803SAlex Fikl 
49047c6ae99SBarry Smith /*@C
49147c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith     Collective on DMComposite
49447c6ae99SBarry Smith 
4959ae5db72SJed Brown     Input Parameters:
49647c6ae99SBarry Smith +    dm - the packer object
49747c6ae99SBarry Smith .    gvec - the global vector
4980298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for those that are not needed
49947c6ae99SBarry Smith 
50047c6ae99SBarry Smith     Level: advanced
50147c6ae99SBarry Smith 
5026f3c3dcfSJed Brown     Notes:
5036f3c3dcfSJed Brown     DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is
5046f3c3dcfSJed Brown     accessible from Fortran.
5056f3c3dcfSJed Brown 
5069ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
5076eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
50847c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
5096f3c3dcfSJed Brown          DMCompositeScatterArray()
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith @*/
5127087cfbeSBarry Smith PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
51347c6ae99SBarry Smith {
51447c6ae99SBarry Smith   va_list                Argp;
51547c6ae99SBarry Smith   PetscErrorCode         ierr;
51647c6ae99SBarry Smith   struct DMCompositeLink *next;
5178fd8f222SJed Brown   PetscInt               cnt;
51847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith   PetscFunctionBegin;
52147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
52247c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
52347c6ae99SBarry Smith   if (!com->setup) {
524d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
52547c6ae99SBarry Smith   }
52647c6ae99SBarry Smith 
52747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
52847c6ae99SBarry Smith   va_start(Argp,gvec);
5298fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
5309ae5db72SJed Brown     Vec local;
5319ae5db72SJed Brown     local = va_arg(Argp, Vec);
5329ae5db72SJed Brown     if (local) {
5339ae5db72SJed Brown       Vec               global;
5345edff71fSBarry Smith       const PetscScalar *array;
5359ae5db72SJed Brown       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
5369ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
5375edff71fSBarry Smith       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
5389ae5db72SJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
5399ae5db72SJed Brown       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
5409ae5db72SJed Brown       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
5415edff71fSBarry Smith       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
5429ae5db72SJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
5439ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
54447c6ae99SBarry Smith     }
54547c6ae99SBarry Smith   }
54647c6ae99SBarry Smith   va_end(Argp);
54747c6ae99SBarry Smith   PetscFunctionReturn(0);
54847c6ae99SBarry Smith }
54947c6ae99SBarry Smith 
5506f3c3dcfSJed Brown /*@
5516f3c3dcfSJed Brown     DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5526f3c3dcfSJed Brown 
5536f3c3dcfSJed Brown     Collective on DMComposite
5546f3c3dcfSJed Brown 
5556f3c3dcfSJed Brown     Input Parameters:
5566f3c3dcfSJed Brown +    dm - the packer object
5576f3c3dcfSJed Brown .    gvec - the global vector
5586f3c3dcfSJed Brown .    lvecs - array of local vectors, NULL for any that are not needed
5596f3c3dcfSJed Brown 
5606f3c3dcfSJed Brown     Level: advanced
5616f3c3dcfSJed Brown 
5626f3c3dcfSJed Brown     Note:
563907376e6SBarry Smith     This is a non-variadic alternative to DMCompositeScatter()
5646f3c3dcfSJed Brown 
5656f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector()
5666f3c3dcfSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
5676f3c3dcfSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
5686f3c3dcfSJed Brown 
5696f3c3dcfSJed Brown @*/
5706f3c3dcfSJed Brown PetscErrorCode  DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs)
5716f3c3dcfSJed Brown {
5726f3c3dcfSJed Brown   PetscErrorCode         ierr;
5736f3c3dcfSJed Brown   struct DMCompositeLink *next;
5746f3c3dcfSJed Brown   PetscInt               i;
5756f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
5766f3c3dcfSJed Brown 
5776f3c3dcfSJed Brown   PetscFunctionBegin;
5786f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5796f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
5806f3c3dcfSJed Brown   if (!com->setup) {
5816f3c3dcfSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
5826f3c3dcfSJed Brown   }
5836f3c3dcfSJed Brown 
5846f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
5856f3c3dcfSJed Brown   for (i=0,next=com->next; next; next=next->next,i++) {
5866f3c3dcfSJed Brown     if (lvecs[i]) {
5876f3c3dcfSJed Brown       Vec         global;
588c5d31e75SLisandro Dalcin       const PetscScalar *array;
5896f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
5906f3c3dcfSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
591c5d31e75SLisandro Dalcin       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
592c5d31e75SLisandro Dalcin       ierr = VecPlaceArray(global,(PetscScalar*)array+next->rstart);CHKERRQ(ierr);
5936f3c3dcfSJed Brown       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
5946f3c3dcfSJed Brown       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
595c5d31e75SLisandro Dalcin       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
5966f3c3dcfSJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
5976f3c3dcfSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
5986f3c3dcfSJed Brown     }
5996f3c3dcfSJed Brown   }
6006f3c3dcfSJed Brown   PetscFunctionReturn(0);
6016f3c3dcfSJed Brown }
6026f3c3dcfSJed Brown 
60347c6ae99SBarry Smith /*@C
60447c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
60547c6ae99SBarry Smith 
60647c6ae99SBarry Smith     Collective on DMComposite
60747c6ae99SBarry Smith 
60847c6ae99SBarry Smith     Input Parameter:
60947c6ae99SBarry Smith +    dm - the packer object
61047c6ae99SBarry Smith .    gvec - the global vector
611907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6120298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for any that are not needed
61347c6ae99SBarry Smith 
61447c6ae99SBarry Smith     Level: advanced
61547c6ae99SBarry Smith 
616f5f57ec0SBarry Smith     Not available from Fortran, Fortran users can use DMCompositeGatherArray()
617f5f57ec0SBarry Smith 
6189ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
6196eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
62047c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
62147c6ae99SBarry Smith 
62247c6ae99SBarry Smith @*/
6231dac896bSSatish Balay PetscErrorCode  DMCompositeGather(DM dm,InsertMode imode,Vec gvec,...)
62447c6ae99SBarry Smith {
62547c6ae99SBarry Smith   va_list                Argp;
62647c6ae99SBarry Smith   PetscErrorCode         ierr;
62747c6ae99SBarry Smith   struct DMCompositeLink *next;
62847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
6298fd8f222SJed Brown   PetscInt               cnt;
63047c6ae99SBarry Smith 
63147c6ae99SBarry Smith   PetscFunctionBegin;
63247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
63347c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
63447c6ae99SBarry Smith   if (!com->setup) {
635d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
63647c6ae99SBarry Smith   }
63747c6ae99SBarry Smith 
63847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
6391dac896bSSatish Balay   va_start(Argp,gvec);
6408fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
6419ae5db72SJed Brown     Vec local;
6429ae5db72SJed Brown     local = va_arg(Argp, Vec);
6439ae5db72SJed Brown     if (local) {
64447c6ae99SBarry Smith       PetscScalar *array;
6459ae5db72SJed Brown       Vec         global;
6469ae5db72SJed Brown       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
6479ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
6489ae5db72SJed Brown       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
6499ae5db72SJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
6509ae5db72SJed Brown       ierr = DMLocalToGlobalBegin(next->dm,local,imode,global);CHKERRQ(ierr);
6519ae5db72SJed Brown       ierr = DMLocalToGlobalEnd(next->dm,local,imode,global);CHKERRQ(ierr);
6529ae5db72SJed Brown       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
6539ae5db72SJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
6549ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
65547c6ae99SBarry Smith     }
65647c6ae99SBarry Smith   }
65747c6ae99SBarry Smith   va_end(Argp);
65847c6ae99SBarry Smith   PetscFunctionReturn(0);
65947c6ae99SBarry Smith }
66047c6ae99SBarry Smith 
6616f3c3dcfSJed Brown /*@
6626f3c3dcfSJed Brown     DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6636f3c3dcfSJed Brown 
6646f3c3dcfSJed Brown     Collective on DMComposite
6656f3c3dcfSJed Brown 
6666f3c3dcfSJed Brown     Input Parameter:
6676f3c3dcfSJed Brown +    dm - the packer object
6686f3c3dcfSJed Brown .    gvec - the global vector
669907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6706f3c3dcfSJed Brown -    lvecs - the individual sequential vectors, NULL for any that are not needed
6716f3c3dcfSJed Brown 
6726f3c3dcfSJed Brown     Level: advanced
6736f3c3dcfSJed Brown 
6746f3c3dcfSJed Brown     Notes:
6756f3c3dcfSJed Brown     This is a non-variadic alternative to DMCompositeGather().
6766f3c3dcfSJed Brown 
6776f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
6786f3c3dcfSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
6796f3c3dcfSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries(),
6806f3c3dcfSJed Brown @*/
6811dac896bSSatish Balay PetscErrorCode  DMCompositeGatherArray(DM dm,InsertMode imode,Vec gvec,Vec *lvecs)
6826f3c3dcfSJed Brown {
6836f3c3dcfSJed Brown   PetscErrorCode         ierr;
6846f3c3dcfSJed Brown   struct DMCompositeLink *next;
6856f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
6866f3c3dcfSJed Brown   PetscInt               i;
6876f3c3dcfSJed Brown 
6886f3c3dcfSJed Brown   PetscFunctionBegin;
6896f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6906f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
6916f3c3dcfSJed Brown   if (!com->setup) {
6926f3c3dcfSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
6936f3c3dcfSJed Brown   }
6946f3c3dcfSJed Brown 
6956f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
6966f3c3dcfSJed Brown   for (next=com->next,i=0; next; next=next->next,i++) {
6976f3c3dcfSJed Brown     if (lvecs[i]) {
6986f3c3dcfSJed Brown       PetscScalar *array;
6996f3c3dcfSJed Brown       Vec         global;
7006f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
7016f3c3dcfSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
7026f3c3dcfSJed Brown       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
7036f3c3dcfSJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
7046f3c3dcfSJed Brown       ierr = DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
7056f3c3dcfSJed Brown       ierr = DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
7066f3c3dcfSJed Brown       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
7076f3c3dcfSJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
7086f3c3dcfSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
7096f3c3dcfSJed Brown     }
7106f3c3dcfSJed Brown   }
7116f3c3dcfSJed Brown   PetscFunctionReturn(0);
7126f3c3dcfSJed Brown }
7136f3c3dcfSJed Brown 
714f5f57ec0SBarry Smith /*@
715aa219208SBarry Smith     DMCompositeAddDM - adds a DM vector to a DMComposite
71647c6ae99SBarry Smith 
71747c6ae99SBarry Smith     Collective on DMComposite
71847c6ae99SBarry Smith 
71947c6ae99SBarry Smith     Input Parameter:
7209b52a9b5SPatrick Sanan +    dmc - the DMComposite (packer) object
721f5f57ec0SBarry Smith -    dm - the DM object
72247c6ae99SBarry Smith 
72347c6ae99SBarry Smith     Level: advanced
72447c6ae99SBarry Smith 
7250c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
7266eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
72747c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith @*/
7307087cfbeSBarry Smith PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
73147c6ae99SBarry Smith {
73247c6ae99SBarry Smith   PetscErrorCode         ierr;
73306ebdd98SJed Brown   PetscInt               n,nlocal;
73447c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
73506ebdd98SJed Brown   Vec                    global,local;
73647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmc->data;
73747c6ae99SBarry Smith 
73847c6ae99SBarry Smith   PetscFunctionBegin;
73947c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
74047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
74147c6ae99SBarry Smith   next = com->next;
742ce94432eSBarry Smith   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
74347c6ae99SBarry Smith 
74447c6ae99SBarry Smith   /* create new link */
745b00a9115SJed Brown   ierr = PetscNew(&mine);CHKERRQ(ierr);
74647c6ae99SBarry Smith   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
74747c6ae99SBarry Smith   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
74847c6ae99SBarry Smith   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
74947c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
75006ebdd98SJed Brown   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
75106ebdd98SJed Brown   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
75206ebdd98SJed Brown   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
7538865f1eaSKarl Rupp 
75447c6ae99SBarry Smith   mine->n      = n;
75506ebdd98SJed Brown   mine->nlocal = nlocal;
75647c6ae99SBarry Smith   mine->dm     = dm;
7570298fd71SBarry Smith   mine->next   = NULL;
75847c6ae99SBarry Smith   com->n      += n;
7597ac2b803SAlex Fikl   com->nghost += nlocal;
76047c6ae99SBarry Smith 
76147c6ae99SBarry Smith   /* add to end of list */
7628865f1eaSKarl Rupp   if (!next) com->next = mine;
7638865f1eaSKarl Rupp   else {
76447c6ae99SBarry Smith     while (next->next) next = next->next;
76547c6ae99SBarry Smith     next->next = mine;
76647c6ae99SBarry Smith   }
76747c6ae99SBarry Smith   com->nDM++;
76847c6ae99SBarry Smith   com->nmine++;
76947c6ae99SBarry Smith   PetscFunctionReturn(0);
77047c6ae99SBarry Smith }
77147c6ae99SBarry Smith 
7729804daf3SBarry Smith #include <petscdraw.h>
77326887b52SJed Brown PETSC_EXTERN PetscErrorCode  VecView_MPI(Vec,PetscViewer);
7747087cfbeSBarry Smith PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
77547c6ae99SBarry Smith {
77647c6ae99SBarry Smith   DM                     dm;
77747c6ae99SBarry Smith   PetscErrorCode         ierr;
77847c6ae99SBarry Smith   struct DMCompositeLink *next;
77947c6ae99SBarry Smith   PetscBool              isdraw;
780cef07954SSatish Balay   DM_Composite           *com;
78147c6ae99SBarry Smith 
78247c6ae99SBarry Smith   PetscFunctionBegin;
783c688c046SMatthew G Knepley   ierr = VecGetDM(gvec, &dm);CHKERRQ(ierr);
784ce94432eSBarry Smith   if (!dm) SETERRQ(PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
78547c6ae99SBarry Smith   com  = (DM_Composite*)dm->data;
78647c6ae99SBarry Smith   next = com->next;
78747c6ae99SBarry Smith 
788251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
78947c6ae99SBarry Smith   if (!isdraw) {
79047c6ae99SBarry Smith     /* do I really want to call this? */
79147c6ae99SBarry Smith     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
79247c6ae99SBarry Smith   } else {
79347c6ae99SBarry Smith     PetscInt cnt = 0;
79447c6ae99SBarry Smith 
79547c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
79647c6ae99SBarry Smith     while (next) {
79747c6ae99SBarry Smith       Vec               vec;
798ca4278abSLisandro Dalcin       const PetscScalar *array;
79947c6ae99SBarry Smith       PetscInt          bs;
80047c6ae99SBarry Smith 
8019ae5db72SJed Brown       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
8029ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&vec);CHKERRQ(ierr);
803ca4278abSLisandro Dalcin       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
804ca4278abSLisandro Dalcin       ierr = VecPlaceArray(vec,(PetscScalar*)array+next->rstart);CHKERRQ(ierr);
805ca4278abSLisandro Dalcin       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
80647c6ae99SBarry Smith       ierr = VecView(vec,viewer);CHKERRQ(ierr);
8079ae5db72SJed Brown       ierr = VecResetArray(vec);CHKERRQ(ierr);
808ca4278abSLisandro Dalcin       ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
8099ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&vec);CHKERRQ(ierr);
81047c6ae99SBarry Smith       ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
81147c6ae99SBarry Smith       cnt += bs;
81247c6ae99SBarry Smith       next = next->next;
81347c6ae99SBarry Smith     }
81447c6ae99SBarry Smith     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
81547c6ae99SBarry Smith   }
81647c6ae99SBarry Smith   PetscFunctionReturn(0);
81747c6ae99SBarry Smith }
81847c6ae99SBarry Smith 
8197087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
82047c6ae99SBarry Smith {
82147c6ae99SBarry Smith   PetscErrorCode ierr;
82247c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
82347c6ae99SBarry Smith 
82447c6ae99SBarry Smith   PetscFunctionBegin;
82547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
826d7bf68aeSBarry Smith   ierr = DMSetUp(dm);CHKERRQ(ierr);
827ce94432eSBarry Smith   ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),com->n,com->N,gvec);CHKERRQ(ierr);
828c688c046SMatthew G Knepley   ierr = VecSetDM(*gvec, dm);CHKERRQ(ierr);
82947c6ae99SBarry Smith   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);CHKERRQ(ierr);
83047c6ae99SBarry Smith   PetscFunctionReturn(0);
83147c6ae99SBarry Smith }
83247c6ae99SBarry Smith 
8337087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
83447c6ae99SBarry Smith {
83547c6ae99SBarry Smith   PetscErrorCode ierr;
83647c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
83747c6ae99SBarry Smith 
83847c6ae99SBarry Smith   PetscFunctionBegin;
83947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
84047c6ae99SBarry Smith   if (!com->setup) {
841d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
84247c6ae99SBarry Smith   }
843f0e01b1fSVincent Le Chenadec   ierr = VecCreateSeq(PETSC_COMM_SELF,com->nghost,lvec);CHKERRQ(ierr);
844c688c046SMatthew G Knepley   ierr = VecSetDM(*lvec, dm);CHKERRQ(ierr);
84547c6ae99SBarry Smith   PetscFunctionReturn(0);
84647c6ae99SBarry Smith }
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith /*@C
8499ae5db72SJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
85047c6ae99SBarry Smith 
85106ebdd98SJed Brown     Collective on DM
85247c6ae99SBarry Smith 
85347c6ae99SBarry Smith     Input Parameter:
85447c6ae99SBarry Smith .    dm - the packer object
85547c6ae99SBarry Smith 
85647c6ae99SBarry Smith     Output Parameters:
8579ae5db72SJed Brown .    ltogs - the individual mappings for each packed vector. Note that this includes
8589ae5db72SJed Brown            all the ghost points that individual ghosted DMDA's may have.
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith     Level: advanced
86147c6ae99SBarry Smith 
86247c6ae99SBarry Smith     Notes:
8636eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
86447c6ae99SBarry Smith 
865f5f57ec0SBarry Smith     Not available from Fortran
866f5f57ec0SBarry Smith 
8679ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
86847c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
86947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
87047c6ae99SBarry Smith 
87147c6ae99SBarry Smith @*/
8727087cfbeSBarry Smith PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
87347c6ae99SBarry Smith {
87447c6ae99SBarry Smith   PetscErrorCode         ierr;
87547c6ae99SBarry Smith   PetscInt               i,*idx,n,cnt;
87647c6ae99SBarry Smith   struct DMCompositeLink *next;
87747c6ae99SBarry Smith   PetscMPIInt            rank;
87847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
87947c6ae99SBarry Smith 
88047c6ae99SBarry Smith   PetscFunctionBegin;
88147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
882728e99d6SJed Brown   ierr = DMSetUp(dm);CHKERRQ(ierr);
883854ce69bSBarry Smith   ierr = PetscMalloc1(com->nDM,ltogs);CHKERRQ(ierr);
88447c6ae99SBarry Smith   next = com->next;
885ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
88647c6ae99SBarry Smith 
88747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
88847c6ae99SBarry Smith   cnt = 0;
88947c6ae99SBarry Smith   while (next) {
8906eb61c8cSJed Brown     ISLocalToGlobalMapping ltog;
8916eb61c8cSJed Brown     PetscMPIInt            size;
89286994e45SJed Brown     const PetscInt         *suboff,*indices;
8936eb61c8cSJed Brown     Vec                    global;
89447c6ae99SBarry Smith 
8956eb61c8cSJed Brown     /* Get sub-DM global indices for each local dof */
8961411c6eeSJed Brown     ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
8976eb61c8cSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
89886994e45SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
899785e854fSJed Brown     ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr);
90047c6ae99SBarry Smith 
9016eb61c8cSJed Brown     /* Get the offsets for the sub-DM global vector */
9026eb61c8cSJed Brown     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
9036eb61c8cSJed Brown     ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
904ce94432eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)global),&size);CHKERRQ(ierr);
9056eb61c8cSJed Brown 
9066eb61c8cSJed Brown     /* Shift the sub-DM definition of the global space to the composite global space */
9076eb61c8cSJed Brown     for (i=0; i<n; i++) {
90886994e45SJed Brown       PetscInt subi = indices[i],lo = 0,hi = size,t;
9096eb61c8cSJed Brown       /* Binary search to find which rank owns subi */
9106eb61c8cSJed Brown       while (hi-lo > 1) {
9116eb61c8cSJed Brown         t = lo + (hi-lo)/2;
9126eb61c8cSJed Brown         if (suboff[t] > subi) hi = t;
9136eb61c8cSJed Brown         else                  lo = t;
9146eb61c8cSJed Brown       }
9156eb61c8cSJed Brown       idx[i] = subi - suboff[lo] + next->grstarts[lo];
9166eb61c8cSJed Brown     }
91786994e45SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
918f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
9196eb61c8cSJed Brown     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
92047c6ae99SBarry Smith     next = next->next;
92147c6ae99SBarry Smith     cnt++;
92247c6ae99SBarry Smith   }
92347c6ae99SBarry Smith   PetscFunctionReturn(0);
92447c6ae99SBarry Smith }
92547c6ae99SBarry Smith 
92687c85e80SJed Brown /*@C
9279ae5db72SJed Brown    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
92887c85e80SJed Brown 
92987c85e80SJed Brown    Not Collective
93087c85e80SJed Brown 
93187c85e80SJed Brown    Input Arguments:
93287c85e80SJed Brown . dm - composite DM
93387c85e80SJed Brown 
93487c85e80SJed Brown    Output Arguments:
93587c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
93687c85e80SJed Brown 
93787c85e80SJed Brown    Level: intermediate
93887c85e80SJed Brown 
93987c85e80SJed Brown    Notes:
94087c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
94187c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
9429ae5db72SJed Brown    scatter to a composite local vector.  The user should not typically need to know which is being done.
94387c85e80SJed Brown 
94487c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
94587c85e80SJed Brown 
94687c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
94787c85e80SJed Brown 
94887c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
94987c85e80SJed Brown 
950f5f57ec0SBarry Smith    Not available from Fortran
951f5f57ec0SBarry Smith 
95287c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
95387c85e80SJed Brown @*/
9547087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
95587c85e80SJed Brown {
95687c85e80SJed Brown   PetscErrorCode         ierr;
95787c85e80SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
95887c85e80SJed Brown   struct DMCompositeLink *link;
95987c85e80SJed Brown   PetscInt               cnt,start;
96087c85e80SJed Brown 
96187c85e80SJed Brown   PetscFunctionBegin;
96287c85e80SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
96387c85e80SJed Brown   PetscValidPointer(is,2);
964785e854fSJed Brown   ierr = PetscMalloc1(com->nmine,is);CHKERRQ(ierr);
96506ebdd98SJed Brown   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
966520db06cSJed Brown     PetscInt bs;
9679ae5db72SJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
9681411c6eeSJed Brown     ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
969520db06cSJed Brown     ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
970520db06cSJed Brown   }
97187c85e80SJed Brown   PetscFunctionReturn(0);
97287c85e80SJed Brown }
97387c85e80SJed Brown 
97447c6ae99SBarry Smith /*@C
97547c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith     Collective on DMComposite
97847c6ae99SBarry Smith 
97947c6ae99SBarry Smith     Input Parameter:
98047c6ae99SBarry Smith .    dm - the packer object
98147c6ae99SBarry Smith 
98247c6ae99SBarry Smith     Output Parameters:
98347c6ae99SBarry Smith .    is - the array of index sets
98447c6ae99SBarry Smith 
98547c6ae99SBarry Smith     Level: advanced
98647c6ae99SBarry Smith 
98747c6ae99SBarry Smith     Notes:
98847c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
98947c6ae99SBarry Smith 
99047c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
99147c6ae99SBarry Smith 
9926eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
9936eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
9946eb61c8cSJed Brown        indices.
99547c6ae99SBarry Smith 
996f3cb0f7eSJed Brown     Fortran Notes:
997f3cb0f7eSJed Brown 
998f3cb0f7eSJed Brown        The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM().
999f3cb0f7eSJed Brown 
10009ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
100147c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
100247c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
100347c6ae99SBarry Smith 
100447c6ae99SBarry Smith @*/
10057087cfbeSBarry Smith PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
100647c6ae99SBarry Smith {
100747c6ae99SBarry Smith   PetscErrorCode         ierr;
100866bb578eSMark Adams   PetscInt               cnt = 0;
100947c6ae99SBarry Smith   struct DMCompositeLink *next;
101047c6ae99SBarry Smith   PetscMPIInt            rank;
101147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
101247c6ae99SBarry Smith 
101347c6ae99SBarry Smith   PetscFunctionBegin;
101447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1015854ce69bSBarry Smith   ierr = PetscMalloc1(com->nDM,is);CHKERRQ(ierr);
101647c6ae99SBarry Smith   next = com->next;
1017ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
101847c6ae99SBarry Smith 
101947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
102047c6ae99SBarry Smith   while (next) {
102166bb578eSMark Adams     ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt]);CHKERRQ(ierr);
10220f21e855SMatthew G. Knepley     if (dm->prob) {
102365c226d8SMatthew G. Knepley       MatNullSpace space;
102465c226d8SMatthew G. Knepley       Mat          pmat;
10250f21e855SMatthew G. Knepley       PetscObject  disc;
10260f21e855SMatthew G. Knepley       PetscInt     Nf;
102765c226d8SMatthew G. Knepley 
10282764a2aaSMatthew G. Knepley       ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
1029f24dd8d2SMatthew G. Knepley       if (cnt < Nf) {
10302764a2aaSMatthew G. Knepley         ierr = PetscDSGetDiscretization(dm->prob, cnt, &disc);CHKERRQ(ierr);
10310f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nullspace", (PetscObject*) &space);CHKERRQ(ierr);
1032aac2dd2dSMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space);CHKERRQ(ierr);}
10330f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space);CHKERRQ(ierr);
1034aac2dd2dSMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space);CHKERRQ(ierr);}
10350f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
1036aac2dd2dSMatthew G. Knepley         if (pmat) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat);CHKERRQ(ierr);}
103765c226d8SMatthew G. Knepley       }
1038f24dd8d2SMatthew G. Knepley     }
103947c6ae99SBarry Smith     cnt++;
104047c6ae99SBarry Smith     next = next->next;
104147c6ae99SBarry Smith   }
104247c6ae99SBarry Smith   PetscFunctionReturn(0);
104347c6ae99SBarry Smith }
104447c6ae99SBarry Smith 
104521c9b008SJed Brown PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields)
10464d343eeaSMatthew G Knepley {
10474d343eeaSMatthew G Knepley   PetscInt       nDM;
10484d343eeaSMatthew G Knepley   DM             *dms;
10494d343eeaSMatthew G Knepley   PetscInt       i;
10504d343eeaSMatthew G Knepley   PetscErrorCode ierr;
10514d343eeaSMatthew G Knepley 
10524d343eeaSMatthew G Knepley   PetscFunctionBegin;
10534d343eeaSMatthew G Knepley   ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
10548865f1eaSKarl Rupp   if (numFields) *numFields = nDM;
10554d343eeaSMatthew G Knepley   ierr = DMCompositeGetGlobalISs(dm, fields);CHKERRQ(ierr);
10564d343eeaSMatthew G Knepley   if (fieldNames) {
1057785e854fSJed Brown     ierr = PetscMalloc1(nDM, &dms);CHKERRQ(ierr);
1058785e854fSJed Brown     ierr = PetscMalloc1(nDM, fieldNames);CHKERRQ(ierr);
10594d343eeaSMatthew G Knepley     ierr = DMCompositeGetEntriesArray(dm, dms);CHKERRQ(ierr);
10604d343eeaSMatthew G Knepley     for (i=0; i<nDM; i++) {
10614d343eeaSMatthew G Knepley       char       buf[256];
10624d343eeaSMatthew G Knepley       const char *splitname;
10634d343eeaSMatthew G Knepley 
10644d343eeaSMatthew G Knepley       /* Split naming precedence: object name, prefix, number */
10654d343eeaSMatthew G Knepley       splitname = ((PetscObject) dm)->name;
10664d343eeaSMatthew G Knepley       if (!splitname) {
10674d343eeaSMatthew G Knepley         ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr);
10684d343eeaSMatthew G Knepley         if (splitname) {
10694d343eeaSMatthew G Knepley           size_t len;
10708caf3d72SBarry Smith           ierr                 = PetscStrncpy(buf,splitname,sizeof(buf));CHKERRQ(ierr);
10718caf3d72SBarry Smith           buf[sizeof(buf) - 1] = 0;
10724d343eeaSMatthew G Knepley           ierr                 = PetscStrlen(buf,&len);CHKERRQ(ierr);
10734d343eeaSMatthew G Knepley           if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
10744d343eeaSMatthew G Knepley           splitname = buf;
10754d343eeaSMatthew G Knepley         }
10764d343eeaSMatthew G Knepley       }
10774d343eeaSMatthew G Knepley       if (!splitname) {
10788caf3d72SBarry Smith         ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",i);CHKERRQ(ierr);
10794d343eeaSMatthew G Knepley         splitname = buf;
10804d343eeaSMatthew G Knepley       }
108121c9b008SJed Brown       ierr = PetscStrallocpy(splitname,&(*fieldNames)[i]);CHKERRQ(ierr);
10824d343eeaSMatthew G Knepley     }
10834d343eeaSMatthew G Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
10844d343eeaSMatthew G Knepley   }
10854d343eeaSMatthew G Knepley   PetscFunctionReturn(0);
10864d343eeaSMatthew G Knepley }
10874d343eeaSMatthew G Knepley 
1088e7c4fc90SDmitry Karpeev /*
1089e7c4fc90SDmitry Karpeev  This could take over from DMCreateFieldIS(), as it is more general,
10900298fd71SBarry Smith  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1091e7c4fc90SDmitry Karpeev  At this point it's probably best to be less intrusive, however.
1092e7c4fc90SDmitry Karpeev  */
109316621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
1094e7c4fc90SDmitry Karpeev {
1095e7c4fc90SDmitry Karpeev   PetscInt       nDM;
1096e7c4fc90SDmitry Karpeev   PetscInt       i;
1097e7c4fc90SDmitry Karpeev   PetscErrorCode ierr;
1098e7c4fc90SDmitry Karpeev 
1099e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
1100e7c4fc90SDmitry Karpeev   ierr = DMCreateFieldIS_Composite(dm, len, namelist, islist);CHKERRQ(ierr);
1101e7c4fc90SDmitry Karpeev   if (dmlist) {
1102e7c4fc90SDmitry Karpeev     ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
1103785e854fSJed Brown     ierr = PetscMalloc1(nDM, dmlist);CHKERRQ(ierr);
1104e7c4fc90SDmitry Karpeev     ierr = DMCompositeGetEntriesArray(dm, *dmlist);CHKERRQ(ierr);
1105e7c4fc90SDmitry Karpeev     for (i=0; i<nDM; i++) {
1106e7c4fc90SDmitry Karpeev       ierr = PetscObjectReference((PetscObject)((*dmlist)[i]));CHKERRQ(ierr);
1107e7c4fc90SDmitry Karpeev     }
1108e7c4fc90SDmitry Karpeev   }
1109e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
1110e7c4fc90SDmitry Karpeev }
1111e7c4fc90SDmitry Karpeev 
1112e7c4fc90SDmitry Karpeev 
1113e7c4fc90SDmitry Karpeev 
111447c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
111547c6ae99SBarry Smith /*@C
11169ae5db72SJed Brown     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
111747c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
111847c6ae99SBarry Smith 
111947c6ae99SBarry Smith     Not Collective
112047c6ae99SBarry Smith 
112147c6ae99SBarry Smith     Input Parameter:
112247c6ae99SBarry Smith .    dm - the packer object
112347c6ae99SBarry Smith 
112447c6ae99SBarry Smith     Output Parameter:
11259ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
112647c6ae99SBarry Smith 
112747c6ae99SBarry Smith     Level: advanced
112847c6ae99SBarry Smith 
1129f5f57ec0SBarry Smith     Not available from Fortran
1130f5f57ec0SBarry Smith 
11319ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
11326eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
113347c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
113447c6ae99SBarry Smith 
113547c6ae99SBarry Smith @*/
11367087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
113747c6ae99SBarry Smith {
113847c6ae99SBarry Smith   va_list                Argp;
113947c6ae99SBarry Smith   PetscErrorCode         ierr;
114047c6ae99SBarry Smith   struct DMCompositeLink *next;
114147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
114247c6ae99SBarry Smith 
114347c6ae99SBarry Smith   PetscFunctionBegin;
114447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
114547c6ae99SBarry Smith   next = com->next;
114647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
114747c6ae99SBarry Smith   va_start(Argp,dm);
114847c6ae99SBarry Smith   while (next) {
114947c6ae99SBarry Smith     Vec *vec;
115047c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
115106930112SJed Brown     if (vec) {ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);}
115247c6ae99SBarry Smith     next = next->next;
115347c6ae99SBarry Smith   }
115447c6ae99SBarry Smith   va_end(Argp);
115547c6ae99SBarry Smith   PetscFunctionReturn(0);
115647c6ae99SBarry Smith }
115747c6ae99SBarry Smith 
115847c6ae99SBarry Smith /*@C
11599ae5db72SJed Brown     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
116047c6ae99SBarry Smith 
116147c6ae99SBarry Smith     Not Collective
116247c6ae99SBarry Smith 
116347c6ae99SBarry Smith     Input Parameter:
116447c6ae99SBarry Smith .    dm - the packer object
116547c6ae99SBarry Smith 
116647c6ae99SBarry Smith     Output Parameter:
11679ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
116847c6ae99SBarry Smith 
116947c6ae99SBarry Smith     Level: advanced
117047c6ae99SBarry Smith 
1171f5f57ec0SBarry Smith     Not available from Fortran
1172f5f57ec0SBarry Smith 
11739ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
11746eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
117547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
117647c6ae99SBarry Smith 
117747c6ae99SBarry Smith @*/
11787087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
117947c6ae99SBarry Smith {
118047c6ae99SBarry Smith   va_list                Argp;
118147c6ae99SBarry Smith   PetscErrorCode         ierr;
118247c6ae99SBarry Smith   struct DMCompositeLink *next;
118347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
118447c6ae99SBarry Smith 
118547c6ae99SBarry Smith   PetscFunctionBegin;
118647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
118747c6ae99SBarry Smith   next = com->next;
118847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
118947c6ae99SBarry Smith   va_start(Argp,dm);
119047c6ae99SBarry Smith   while (next) {
119147c6ae99SBarry Smith     Vec *vec;
119247c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
119306930112SJed Brown     if (vec) {ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);}
119447c6ae99SBarry Smith     next = next->next;
119547c6ae99SBarry Smith   }
119647c6ae99SBarry Smith   va_end(Argp);
119747c6ae99SBarry Smith   PetscFunctionReturn(0);
119847c6ae99SBarry Smith }
119947c6ae99SBarry Smith 
120047c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
120147c6ae99SBarry Smith /*@C
12029ae5db72SJed Brown     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
120347c6ae99SBarry Smith 
120447c6ae99SBarry Smith     Not Collective
120547c6ae99SBarry Smith 
120647c6ae99SBarry Smith     Input Parameter:
120747c6ae99SBarry Smith .    dm - the packer object
120847c6ae99SBarry Smith 
120947c6ae99SBarry Smith     Output Parameter:
12109ae5db72SJed Brown .   DM ... - the individual entries (DMs)
121147c6ae99SBarry Smith 
121247c6ae99SBarry Smith     Level: advanced
121347c6ae99SBarry Smith 
1214*95452b02SPatrick Sanan     Fortran Notes:
1215*95452b02SPatrick Sanan     Available as DMCompositeGetEntries() for one output DM, DMCompositeGetEntries2() for 2, etc
1216f5f57ec0SBarry Smith 
12172fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
12186eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
121947c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
122047c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
122147c6ae99SBarry Smith 
122247c6ae99SBarry Smith @*/
12237087cfbeSBarry Smith PetscErrorCode  DMCompositeGetEntries(DM dm,...)
122447c6ae99SBarry Smith {
122547c6ae99SBarry Smith   va_list                Argp;
122647c6ae99SBarry Smith   struct DMCompositeLink *next;
122747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
122847c6ae99SBarry Smith 
122947c6ae99SBarry Smith   PetscFunctionBegin;
123047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
123147c6ae99SBarry Smith   next = com->next;
123247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
123347c6ae99SBarry Smith   va_start(Argp,dm);
123447c6ae99SBarry Smith   while (next) {
123547c6ae99SBarry Smith     DM *dmn;
123647c6ae99SBarry Smith     dmn = va_arg(Argp, DM*);
12379ae5db72SJed Brown     if (dmn) *dmn = next->dm;
123847c6ae99SBarry Smith     next = next->next;
123947c6ae99SBarry Smith   }
124047c6ae99SBarry Smith   va_end(Argp);
124147c6ae99SBarry Smith   PetscFunctionReturn(0);
124247c6ae99SBarry Smith }
124347c6ae99SBarry Smith 
1244dbab29e1SMark F. Adams /*@C
12452fa5ba8aSJed Brown     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
12462fa5ba8aSJed Brown 
12472fa5ba8aSJed Brown     Not Collective
12482fa5ba8aSJed Brown 
12492fa5ba8aSJed Brown     Input Parameter:
1250907376e6SBarry Smith .    dm - the packer object
1251907376e6SBarry Smith 
1252907376e6SBarry Smith     Output Parameter:
1253907376e6SBarry Smith .    dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs
12542fa5ba8aSJed Brown 
12552fa5ba8aSJed Brown     Level: advanced
12562fa5ba8aSJed Brown 
12572fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
12582fa5ba8aSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
12592fa5ba8aSJed Brown          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
12602fa5ba8aSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
12612fa5ba8aSJed Brown 
12622fa5ba8aSJed Brown @*/
12632fa5ba8aSJed Brown PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
12642fa5ba8aSJed Brown {
12652fa5ba8aSJed Brown   struct DMCompositeLink *next;
12662fa5ba8aSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
12672fa5ba8aSJed Brown   PetscInt               i;
12682fa5ba8aSJed Brown 
12692fa5ba8aSJed Brown   PetscFunctionBegin;
12702fa5ba8aSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12712fa5ba8aSJed Brown   /* loop over packed objects, handling one at at time */
12722fa5ba8aSJed Brown   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
12732fa5ba8aSJed Brown   PetscFunctionReturn(0);
12742fa5ba8aSJed Brown }
12752fa5ba8aSJed Brown 
1276e10fd815SStefano Zampini typedef struct {
1277e10fd815SStefano Zampini   DM          dm;
1278e10fd815SStefano Zampini   PetscViewer *subv;
1279e10fd815SStefano Zampini   Vec         *vecs;
1280e10fd815SStefano Zampini } GLVisViewerCtx;
1281e10fd815SStefano Zampini 
1282e10fd815SStefano Zampini static PetscErrorCode  DestroyGLVisViewerCtx_Private(void *vctx)
1283e10fd815SStefano Zampini {
1284e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
1285e10fd815SStefano Zampini   PetscInt       i,n;
1286e10fd815SStefano Zampini   PetscErrorCode ierr;
1287e10fd815SStefano Zampini 
1288e10fd815SStefano Zampini   PetscFunctionBegin;
1289e10fd815SStefano Zampini   ierr = DMCompositeGetNumberDM(ctx->dm,&n);CHKERRQ(ierr);
1290e10fd815SStefano Zampini   for (i = 0; i < n; i++) {
1291e10fd815SStefano Zampini     ierr = PetscViewerDestroy(&ctx->subv[i]);CHKERRQ(ierr);
1292e10fd815SStefano Zampini   }
1293e10fd815SStefano Zampini   ierr = PetscFree2(ctx->subv,ctx->vecs);CHKERRQ(ierr);
1294e10fd815SStefano Zampini   ierr = DMDestroy(&ctx->dm);CHKERRQ(ierr);
1295e10fd815SStefano Zampini   ierr = PetscFree(ctx);CHKERRQ(ierr);
1296e10fd815SStefano Zampini   PetscFunctionReturn(0);
1297e10fd815SStefano Zampini }
1298e10fd815SStefano Zampini 
1299e10fd815SStefano Zampini static PetscErrorCode  DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
1300e10fd815SStefano Zampini {
1301e10fd815SStefano Zampini   Vec            X = (Vec)oX;
1302e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
1303e10fd815SStefano Zampini   PetscInt       i,n,cumf;
1304e10fd815SStefano Zampini   PetscErrorCode ierr;
1305e10fd815SStefano Zampini 
1306e10fd815SStefano Zampini   PetscFunctionBegin;
1307e10fd815SStefano Zampini   ierr = DMCompositeGetNumberDM(ctx->dm,&n);CHKERRQ(ierr);
1308e10fd815SStefano Zampini   ierr = DMCompositeGetAccessArray(ctx->dm,X,n,NULL,ctx->vecs);CHKERRQ(ierr);
1309e10fd815SStefano Zampini   for (i = 0, cumf = 0; i < n; i++) {
1310e10fd815SStefano Zampini     PetscErrorCode (*g2l)(PetscObject,PetscInt,PetscObject[],void*);
1311e10fd815SStefano Zampini     void           *fctx;
1312e10fd815SStefano Zampini     PetscInt       nfi;
1313e10fd815SStefano Zampini 
1314e10fd815SStefano Zampini     ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nfi,NULL,NULL,&g2l,NULL,&fctx);CHKERRQ(ierr);
1315e10fd815SStefano Zampini     if (!nfi) continue;
1316e10fd815SStefano Zampini     if (g2l) {
1317e10fd815SStefano Zampini       ierr = (*g2l)((PetscObject)ctx->vecs[i],nfi,oXfield+cumf,fctx);CHKERRQ(ierr);
1318e10fd815SStefano Zampini     } else {
1319e10fd815SStefano Zampini       ierr = VecCopy(ctx->vecs[i],(Vec)(oXfield[cumf]));CHKERRQ(ierr);
1320e10fd815SStefano Zampini     }
1321e10fd815SStefano Zampini     cumf += nfi;
1322e10fd815SStefano Zampini   }
1323e10fd815SStefano Zampini   ierr = DMCompositeRestoreAccessArray(ctx->dm,X,n,NULL,ctx->vecs);CHKERRQ(ierr);
1324e10fd815SStefano Zampini   PetscFunctionReturn(0);
1325e10fd815SStefano Zampini }
1326e10fd815SStefano Zampini 
1327e10fd815SStefano Zampini static PetscErrorCode  DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer)
1328e10fd815SStefano Zampini {
1329e10fd815SStefano Zampini   DM             dm = (DM)odm, *dms;
1330e10fd815SStefano Zampini   Vec            *Ufds;
1331e10fd815SStefano Zampini   GLVisViewerCtx *ctx;
1332e10fd815SStefano Zampini   PetscInt       i,n,tnf,*sdim;
1333e10fd815SStefano Zampini   char           **fecs;
1334e10fd815SStefano Zampini   PetscErrorCode ierr;
1335e10fd815SStefano Zampini 
1336e10fd815SStefano Zampini   PetscFunctionBegin;
1337e10fd815SStefano Zampini   ierr = PetscNew(&ctx);CHKERRQ(ierr);
1338e10fd815SStefano Zampini   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
1339e10fd815SStefano Zampini   ctx->dm = dm;
1340e10fd815SStefano Zampini   ierr = DMCompositeGetNumberDM(dm,&n);CHKERRQ(ierr);
1341e10fd815SStefano Zampini   ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr);
1342e10fd815SStefano Zampini   ierr = DMCompositeGetEntriesArray(dm,dms);CHKERRQ(ierr);
1343e10fd815SStefano Zampini   ierr = PetscMalloc2(n,&ctx->subv,n,&ctx->vecs);CHKERRQ(ierr);
1344e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1345e10fd815SStefano Zampini     PetscInt nf;
1346e10fd815SStefano Zampini 
1347e10fd815SStefano Zampini     ierr = PetscViewerCreate(PetscObjectComm(odm),&ctx->subv[i]);CHKERRQ(ierr);
1348e10fd815SStefano Zampini     ierr = PetscViewerSetType(ctx->subv[i],PETSCVIEWERGLVIS);CHKERRQ(ierr);
1349e10fd815SStefano Zampini     ierr = PetscViewerGLVisSetDM_Private(ctx->subv[i],(PetscObject)dms[i]);CHKERRQ(ierr);
1350e10fd815SStefano Zampini     ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1351e10fd815SStefano Zampini     tnf += nf;
1352e10fd815SStefano Zampini   }
1353e10fd815SStefano Zampini   ierr = PetscFree(dms);CHKERRQ(ierr);
1354e10fd815SStefano Zampini   ierr = PetscMalloc3(tnf,&fecs,tnf,&sdim,tnf,&Ufds);CHKERRQ(ierr);
1355e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1356e10fd815SStefano Zampini     PetscInt   *sd,nf,f;
1357e10fd815SStefano Zampini     const char **fec;
1358e10fd815SStefano Zampini     Vec        *Uf;
1359e10fd815SStefano Zampini 
1360e10fd815SStefano Zampini     ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,&fec,&sd,NULL,(PetscObject**)&Uf,NULL);CHKERRQ(ierr);
1361e10fd815SStefano Zampini     for (f = 0; f < nf; f++) {
1362e10fd815SStefano Zampini       ierr = PetscStrallocpy(fec[f],&fecs[tnf+f]);CHKERRQ(ierr);
1363e10fd815SStefano Zampini       Ufds[tnf+f] = Uf[f];
1364e10fd815SStefano Zampini       sdim[tnf+f] = sd[f];
1365e10fd815SStefano Zampini     }
1366e10fd815SStefano Zampini     tnf += nf;
1367e10fd815SStefano Zampini   }
1368e10fd815SStefano Zampini   ierr = PetscViewerGLVisSetFields(viewer,tnf,(const char**)fecs,sdim,DMCompositeSampleGLVisFields_Private,(PetscObject*)Ufds,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
1369e10fd815SStefano Zampini   for (i = 0; i < tnf; i++) {
1370e10fd815SStefano Zampini     ierr = PetscFree(fecs[i]);CHKERRQ(ierr);
1371e10fd815SStefano Zampini   }
1372e10fd815SStefano Zampini   ierr = PetscFree3(fecs,sdim,Ufds);CHKERRQ(ierr);
1373e10fd815SStefano Zampini   PetscFunctionReturn(0);
1374e10fd815SStefano Zampini }
1375e10fd815SStefano Zampini 
13767087cfbeSBarry Smith PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
137747c6ae99SBarry Smith {
137847c6ae99SBarry Smith   PetscErrorCode         ierr;
137947c6ae99SBarry Smith   struct DMCompositeLink *next;
138047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmi->data;
138147c6ae99SBarry Smith   DM                     dm;
138247c6ae99SBarry Smith 
138347c6ae99SBarry Smith   PetscFunctionBegin;
138447c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1385ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
1386ce94432eSBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
1387ce94432eSBarry Smith   }
13882ce3a92bSJed Brown   ierr = DMSetUp(dmi);CHKERRQ(ierr);
138947c6ae99SBarry Smith   next = com->next;
139047c6ae99SBarry Smith   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
139147c6ae99SBarry Smith 
139247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
139347c6ae99SBarry Smith   while (next) {
139447c6ae99SBarry Smith     ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
139547c6ae99SBarry Smith     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
139647c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
139747c6ae99SBarry Smith     next = next->next;
139847c6ae99SBarry Smith   }
139947c6ae99SBarry Smith   PetscFunctionReturn(0);
140047c6ae99SBarry Smith }
140147c6ae99SBarry Smith 
140214354c39SJed Brown PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
140314354c39SJed Brown {
140414354c39SJed Brown   PetscErrorCode         ierr;
140514354c39SJed Brown   struct DMCompositeLink *next;
140614354c39SJed Brown   DM_Composite           *com = (DM_Composite*)dmi->data;
140714354c39SJed Brown   DM                     dm;
140814354c39SJed Brown 
140914354c39SJed Brown   PetscFunctionBegin;
141014354c39SJed Brown   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
14112ce3a92bSJed Brown   ierr = DMSetUp(dmi);CHKERRQ(ierr);
14122ee06e3bSJed Brown   if (comm == MPI_COMM_NULL) {
141325296bd5SBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
141425296bd5SBarry Smith   }
141514354c39SJed Brown   next = com->next;
141614354c39SJed Brown   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
141714354c39SJed Brown 
141814354c39SJed Brown   /* loop over packed objects, handling one at at time */
141914354c39SJed Brown   while (next) {
142014354c39SJed Brown     ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr);
142114354c39SJed Brown     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
142214354c39SJed Brown     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
142314354c39SJed Brown     next = next->next;
142414354c39SJed Brown   }
142514354c39SJed Brown   PetscFunctionReturn(0);
142614354c39SJed Brown }
142747c6ae99SBarry Smith 
1428e727c939SJed Brown PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
142947c6ae99SBarry Smith {
143047c6ae99SBarry Smith   PetscErrorCode         ierr;
14319ae5db72SJed Brown   PetscInt               m,n,M,N,nDM,i;
143247c6ae99SBarry Smith   struct DMCompositeLink *nextc;
143347c6ae99SBarry Smith   struct DMCompositeLink *nextf;
143425296bd5SBarry Smith   Vec                    gcoarse,gfine,*vecs;
143547c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
143647c6ae99SBarry Smith   DM_Composite           *comfine   = (DM_Composite*)fine->data;
14379ae5db72SJed Brown   Mat                    *mats;
143847c6ae99SBarry Smith 
143947c6ae99SBarry Smith   PetscFunctionBegin;
144047c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
144147c6ae99SBarry Smith   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
1442f692024eSJed Brown   ierr = DMSetUp(coarse);CHKERRQ(ierr);
1443f692024eSJed Brown   ierr = DMSetUp(fine);CHKERRQ(ierr);
144447c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14459ae5db72SJed Brown   ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
14469ae5db72SJed Brown   ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr);
144747c6ae99SBarry Smith   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
144847c6ae99SBarry Smith   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
144947c6ae99SBarry Smith   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
145047c6ae99SBarry Smith   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
14519ae5db72SJed Brown   ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
14529ae5db72SJed Brown   ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr);
145347c6ae99SBarry Smith 
14549ae5db72SJed Brown   nDM = comfine->nDM;
1455ce94432eSBarry 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);
14561795a4d1SJed Brown   ierr = PetscCalloc1(nDM*nDM,&mats);CHKERRQ(ierr);
145725296bd5SBarry Smith   if (v) {
14581795a4d1SJed Brown     ierr = PetscCalloc1(nDM,&vecs);CHKERRQ(ierr);
145925296bd5SBarry Smith   }
146047c6ae99SBarry Smith 
146147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
14629ae5db72SJed Brown   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
146325296bd5SBarry Smith     if (!v) {
14640298fd71SBarry Smith       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL);CHKERRQ(ierr);
146525296bd5SBarry Smith     } else {
146625296bd5SBarry Smith       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr);
146725296bd5SBarry Smith     }
146847c6ae99SBarry Smith   }
1469ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A);CHKERRQ(ierr);
147025296bd5SBarry Smith   if (v) {
1471ce94432eSBarry Smith     ierr = VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v);CHKERRQ(ierr);
147225296bd5SBarry Smith   }
14739ae5db72SJed Brown   for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);}
14749ae5db72SJed Brown   ierr = PetscFree(mats);CHKERRQ(ierr);
147525296bd5SBarry Smith   if (v) {
147625296bd5SBarry Smith     for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);}
147725296bd5SBarry Smith     ierr = PetscFree(vecs);CHKERRQ(ierr);
147825296bd5SBarry Smith   }
147947c6ae99SBarry Smith   PetscFunctionReturn(0);
148047c6ae99SBarry Smith }
148147c6ae99SBarry Smith 
1482184d77edSJed Brown static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
14831411c6eeSJed Brown {
14841411c6eeSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
14851411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1486f7efa3c7SJed Brown   PetscInt               i;
14871411c6eeSJed Brown   PetscErrorCode         ierr;
14881411c6eeSJed Brown 
14891411c6eeSJed Brown   PetscFunctionBegin;
14901411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
14911411c6eeSJed Brown   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
1492ce94432eSBarry Smith   ierr = ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr);
14939ae5db72SJed Brown   for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
14941411c6eeSJed Brown   ierr = PetscFree(ltogs);CHKERRQ(ierr);
14951411c6eeSJed Brown   PetscFunctionReturn(0);
14961411c6eeSJed Brown }
14971411c6eeSJed Brown 
14981411c6eeSJed Brown 
1499b412c318SBarry Smith PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring)
150047c6ae99SBarry Smith {
150147c6ae99SBarry Smith   PetscErrorCode  ierr;
150247c6ae99SBarry Smith   PetscInt        n,i,cnt;
150347c6ae99SBarry Smith   ISColoringValue *colors;
150447c6ae99SBarry Smith   PetscBool       dense  = PETSC_FALSE;
150547c6ae99SBarry Smith   ISColoringValue maxcol = 0;
150647c6ae99SBarry Smith   DM_Composite    *com   = (DM_Composite*)dm->data;
150747c6ae99SBarry Smith 
150847c6ae99SBarry Smith   PetscFunctionBegin;
150947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
15105bdb020cSBarry Smith   if (ctype == IS_COLORING_LOCAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported");
1511e3247f34SBarry Smith   else if (ctype == IS_COLORING_GLOBAL) {
151247c6ae99SBarry Smith     n = com->n;
1513ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1514785e854fSJed Brown   ierr = PetscMalloc1(n,&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
151547c6ae99SBarry Smith 
1516c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL);CHKERRQ(ierr);
151747c6ae99SBarry Smith   if (dense) {
151847c6ae99SBarry Smith     for (i=0; i<n; i++) {
151947c6ae99SBarry Smith       colors[i] = (ISColoringValue)(com->rstart + i);
152047c6ae99SBarry Smith     }
152147c6ae99SBarry Smith     maxcol = com->N;
152247c6ae99SBarry Smith   } else {
152347c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
152447c6ae99SBarry Smith     PetscMPIInt            rank;
152547c6ae99SBarry Smith 
1526ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
152747c6ae99SBarry Smith     cnt  = 0;
152847c6ae99SBarry Smith     while (next) {
152947c6ae99SBarry Smith       ISColoring lcoloring;
153047c6ae99SBarry Smith 
1531b412c318SBarry Smith       ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring);CHKERRQ(ierr);
153247c6ae99SBarry Smith       for (i=0; i<lcoloring->N; i++) {
153347c6ae99SBarry Smith         colors[cnt++] = maxcol + lcoloring->colors[i];
153447c6ae99SBarry Smith       }
153547c6ae99SBarry Smith       maxcol += lcoloring->n;
1536fcfd50ebSBarry Smith       ierr    = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
153747c6ae99SBarry Smith       next    = next->next;
153847c6ae99SBarry Smith     }
153947c6ae99SBarry Smith   }
1540aaf3ff59SMatthew G. Knepley   ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,PETSC_OWN_POINTER,coloring);CHKERRQ(ierr);
154147c6ae99SBarry Smith   PetscFunctionReturn(0);
154247c6ae99SBarry Smith }
154347c6ae99SBarry Smith 
15447087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
154547c6ae99SBarry Smith {
154647c6ae99SBarry Smith   PetscErrorCode         ierr;
154747c6ae99SBarry Smith   struct DMCompositeLink *next;
154847c6ae99SBarry Smith   PetscScalar            *garray,*larray;
154947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
155047c6ae99SBarry Smith 
155147c6ae99SBarry Smith   PetscFunctionBegin;
155247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
155347c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
155439d35262SVincent Le Chenadec 
155547c6ae99SBarry Smith   if (!com->setup) {
1556d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
155747c6ae99SBarry Smith   }
155839d35262SVincent Le Chenadec 
155947c6ae99SBarry Smith   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
156047c6ae99SBarry Smith   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
156147c6ae99SBarry Smith 
156247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
156339d35262SVincent Le Chenadec   next = com->next;
156447c6ae99SBarry Smith   while (next) {
156547c6ae99SBarry Smith     Vec      local,global;
156647c6ae99SBarry Smith     PetscInt N;
156747c6ae99SBarry Smith 
156847c6ae99SBarry Smith     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
156947c6ae99SBarry Smith     ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
157047c6ae99SBarry Smith     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
157147c6ae99SBarry Smith     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
157247c6ae99SBarry Smith     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
157347c6ae99SBarry Smith     ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
157447c6ae99SBarry Smith     ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
157547c6ae99SBarry Smith     ierr = VecResetArray(global);CHKERRQ(ierr);
157647c6ae99SBarry Smith     ierr = VecResetArray(local);CHKERRQ(ierr);
157747c6ae99SBarry Smith     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
157839d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr);
157939d35262SVincent Le Chenadec 
158006ebdd98SJed Brown     larray += next->nlocal;
158139d35262SVincent Le Chenadec     garray += next->n;
158247c6ae99SBarry Smith     next    = next->next;
158347c6ae99SBarry Smith   }
158447c6ae99SBarry Smith 
15850298fd71SBarry Smith   ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr);
15860298fd71SBarry Smith   ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr);
158747c6ae99SBarry Smith   PetscFunctionReturn(0);
158847c6ae99SBarry Smith }
158947c6ae99SBarry Smith 
15907087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
15910c010503SBarry Smith {
15920c010503SBarry Smith   PetscFunctionBegin;
159339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
159439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
159539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,4);
159639d35262SVincent Le Chenadec   PetscFunctionReturn(0);
159739d35262SVincent Le Chenadec }
159839d35262SVincent Le Chenadec 
159939d35262SVincent Le Chenadec PetscErrorCode  DMLocalToGlobalBegin_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
160039d35262SVincent Le Chenadec {
160139d35262SVincent Le Chenadec   PetscErrorCode         ierr;
160239d35262SVincent Le Chenadec   struct DMCompositeLink *next;
160339d35262SVincent Le Chenadec   PetscScalar            *larray,*garray;
160439d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite*)dm->data;
160539d35262SVincent Le Chenadec 
160639d35262SVincent Le Chenadec   PetscFunctionBegin;
160739d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
160839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
160939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
161039d35262SVincent Le Chenadec 
161139d35262SVincent Le Chenadec   if (!com->setup) {
161239d35262SVincent Le Chenadec     ierr = DMSetUp(dm);CHKERRQ(ierr);
161339d35262SVincent Le Chenadec   }
161439d35262SVincent Le Chenadec 
161539d35262SVincent Le Chenadec   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
161639d35262SVincent Le Chenadec   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
161739d35262SVincent Le Chenadec 
161839d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
161939d35262SVincent Le Chenadec   next = com->next;
162039d35262SVincent Le Chenadec   while (next) {
162139d35262SVincent Le Chenadec     Vec      global,local;
162239d35262SVincent Le Chenadec 
162339d35262SVincent Le Chenadec     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
162439d35262SVincent Le Chenadec     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
162539d35262SVincent Le Chenadec     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
162639d35262SVincent Le Chenadec     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
162739d35262SVincent Le Chenadec     ierr = DMLocalToGlobalBegin(next->dm,local,mode,global);CHKERRQ(ierr);
162839d35262SVincent Le Chenadec     ierr = DMLocalToGlobalEnd(next->dm,local,mode,global);CHKERRQ(ierr);
162939d35262SVincent Le Chenadec     ierr = VecResetArray(local);CHKERRQ(ierr);
163039d35262SVincent Le Chenadec     ierr = VecResetArray(global);CHKERRQ(ierr);
163139d35262SVincent Le Chenadec     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
163239d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr);
163339d35262SVincent Le Chenadec 
163439d35262SVincent Le Chenadec     garray += next->n;
163539d35262SVincent Le Chenadec     larray += next->nlocal;
163639d35262SVincent Le Chenadec     next    = next->next;
163739d35262SVincent Le Chenadec   }
163839d35262SVincent Le Chenadec 
163939d35262SVincent Le Chenadec   ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr);
164039d35262SVincent Le Chenadec   ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr);
164139d35262SVincent Le Chenadec 
164239d35262SVincent Le Chenadec   PetscFunctionReturn(0);
164339d35262SVincent Le Chenadec }
164439d35262SVincent Le Chenadec 
164539d35262SVincent Le Chenadec PetscErrorCode  DMLocalToGlobalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
164639d35262SVincent Le Chenadec {
164739d35262SVincent Le Chenadec   PetscFunctionBegin;
164839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
164939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
165039d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
165139d35262SVincent Le Chenadec   PetscFunctionReturn(0);
165239d35262SVincent Le Chenadec }
165339d35262SVincent Le Chenadec 
165439d35262SVincent Le Chenadec PetscErrorCode  DMLocalToLocalBegin_Composite(DM dm,Vec vec1,InsertMode mode,Vec vec2)
165539d35262SVincent Le Chenadec {
165639d35262SVincent Le Chenadec   PetscErrorCode         ierr;
165739d35262SVincent Le Chenadec   struct DMCompositeLink *next;
165839d35262SVincent Le Chenadec   PetscScalar            *array1,*array2;
165939d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite*)dm->data;
166039d35262SVincent Le Chenadec 
166139d35262SVincent Le Chenadec   PetscFunctionBegin;
166239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
166339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec1,VEC_CLASSID,2);
166439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec2,VEC_CLASSID,4);
166539d35262SVincent Le Chenadec 
166639d35262SVincent Le Chenadec   if (!com->setup) {
166739d35262SVincent Le Chenadec     ierr = DMSetUp(dm);CHKERRQ(ierr);
166839d35262SVincent Le Chenadec   }
166939d35262SVincent Le Chenadec 
167039d35262SVincent Le Chenadec   ierr = VecGetArray(vec1,&array1);CHKERRQ(ierr);
167139d35262SVincent Le Chenadec   ierr = VecGetArray(vec2,&array2);CHKERRQ(ierr);
167239d35262SVincent Le Chenadec 
167339d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
167439d35262SVincent Le Chenadec   next = com->next;
167539d35262SVincent Le Chenadec   while (next) {
167639d35262SVincent Le Chenadec     Vec      local1,local2;
167739d35262SVincent Le Chenadec 
167839d35262SVincent Le Chenadec     ierr = DMGetLocalVector(next->dm,&local1);CHKERRQ(ierr);
167939d35262SVincent Le Chenadec     ierr = VecPlaceArray(local1,array1);CHKERRQ(ierr);
168039d35262SVincent Le Chenadec     ierr = DMGetLocalVector(next->dm,&local2);CHKERRQ(ierr);
168139d35262SVincent Le Chenadec     ierr = VecPlaceArray(local2,array2);CHKERRQ(ierr);
168239d35262SVincent Le Chenadec     ierr = DMLocalToLocalBegin(next->dm,local1,mode,local2);CHKERRQ(ierr);
168339d35262SVincent Le Chenadec     ierr = DMLocalToLocalEnd(next->dm,local1,mode,local2);CHKERRQ(ierr);
168439d35262SVincent Le Chenadec     ierr = VecResetArray(local2);CHKERRQ(ierr);
168539d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local2);CHKERRQ(ierr);
168639d35262SVincent Le Chenadec     ierr = VecResetArray(local1);CHKERRQ(ierr);
168739d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local1);CHKERRQ(ierr);
168839d35262SVincent Le Chenadec 
168939d35262SVincent Le Chenadec     array1 += next->nlocal;
169039d35262SVincent Le Chenadec     array2 += next->nlocal;
169139d35262SVincent Le Chenadec     next    = next->next;
169239d35262SVincent Le Chenadec   }
169339d35262SVincent Le Chenadec 
169439d35262SVincent Le Chenadec   ierr = VecRestoreArray(vec1,NULL);CHKERRQ(ierr);
169539d35262SVincent Le Chenadec   ierr = VecRestoreArray(vec2,NULL);CHKERRQ(ierr);
169639d35262SVincent Le Chenadec 
169739d35262SVincent Le Chenadec   PetscFunctionReturn(0);
169839d35262SVincent Le Chenadec }
169939d35262SVincent Le Chenadec 
170039d35262SVincent Le Chenadec PetscErrorCode  DMLocalToLocalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
170139d35262SVincent Le Chenadec {
170239d35262SVincent Le Chenadec   PetscFunctionBegin;
170339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
170439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
170539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
17060c010503SBarry Smith   PetscFunctionReturn(0);
17070c010503SBarry Smith }
170847c6ae99SBarry Smith 
17096ae3a549SBarry Smith /*MC
17106ae3a549SBarry Smith    DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs
17116ae3a549SBarry Smith 
17126ae3a549SBarry Smith   Level: intermediate
17136ae3a549SBarry Smith 
17141abcec8cSBarry Smith .seealso: DMType, DM, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate()
17156ae3a549SBarry Smith M*/
17166ae3a549SBarry Smith 
17176ae3a549SBarry Smith 
17188cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1719a4121054SBarry Smith {
1720a4121054SBarry Smith   PetscErrorCode ierr;
1721a4121054SBarry Smith   DM_Composite   *com;
1722a4121054SBarry Smith 
1723a4121054SBarry Smith   PetscFunctionBegin;
1724b00a9115SJed Brown   ierr          = PetscNewLog(p,&com);CHKERRQ(ierr);
1725a4121054SBarry Smith   p->data       = com;
1726a4121054SBarry Smith   ierr          = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1727a4121054SBarry Smith   com->n        = 0;
17287ac2b803SAlex Fikl   com->nghost   = 0;
17290298fd71SBarry Smith   com->next     = NULL;
1730a4121054SBarry Smith   com->nDM      = 0;
1731a4121054SBarry Smith 
1732a4121054SBarry Smith   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1733a4121054SBarry Smith   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1734184d77edSJed Brown   p->ops->getlocaltoglobalmapping         = DMGetLocalToGlobalMapping_Composite;
17354d343eeaSMatthew G Knepley   p->ops->createfieldis                   = DMCreateFieldIS_Composite;
173616621825SDmitry Karpeev   p->ops->createfielddecomposition        = DMCreateFieldDecomposition_Composite;
1737a4121054SBarry Smith   p->ops->refine                          = DMRefine_Composite;
173814354c39SJed Brown   p->ops->coarsen                         = DMCoarsen_Composite;
173925296bd5SBarry Smith   p->ops->createinterpolation             = DMCreateInterpolation_Composite;
174025296bd5SBarry Smith   p->ops->creatematrix                    = DMCreateMatrix_Composite;
1741e727c939SJed Brown   p->ops->getcoloring                     = DMCreateColoring_Composite;
1742a4121054SBarry Smith   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1743a4121054SBarry Smith   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
174439d35262SVincent Le Chenadec   p->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Composite;
174539d35262SVincent Le Chenadec   p->ops->localtoglobalend                = DMLocalToGlobalEnd_Composite;
174639d35262SVincent Le Chenadec   p->ops->localtolocalbegin               = DMLocalToLocalBegin_Composite;
174739d35262SVincent Le Chenadec   p->ops->localtolocalend                 = DMLocalToLocalEnd_Composite;
1748a4121054SBarry Smith   p->ops->destroy                         = DMDestroy_Composite;
1749a4121054SBarry Smith   p->ops->view                            = DMView_Composite;
1750a4121054SBarry Smith   p->ops->setup                           = DMSetUp_Composite;
1751e10fd815SStefano Zampini 
1752e10fd815SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)p,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Composite);CHKERRQ(ierr);
1753a4121054SBarry Smith   PetscFunctionReturn(0);
1754a4121054SBarry Smith }
1755a4121054SBarry Smith 
17561fd49c25SBarry Smith /*@
17570c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
17580c010503SBarry Smith       vectors made up of several subvectors.
17590c010503SBarry Smith 
17600c010503SBarry Smith     Collective on MPI_Comm
176147c6ae99SBarry Smith 
176247c6ae99SBarry Smith     Input Parameter:
17630c010503SBarry Smith .   comm - the processors that will share the global vector
17640c010503SBarry Smith 
17650c010503SBarry Smith     Output Parameters:
17660c010503SBarry Smith .   packer - the packer object
176747c6ae99SBarry Smith 
176847c6ae99SBarry Smith     Level: advanced
176947c6ae99SBarry Smith 
17701abcec8cSBarry Smith .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(), DMCOMPOSITE,DMCreate()
17716eb61c8cSJed Brown          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
177247c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
177347c6ae99SBarry Smith 
177447c6ae99SBarry Smith @*/
17757087cfbeSBarry Smith PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
177647c6ae99SBarry Smith {
17770c010503SBarry Smith   PetscErrorCode ierr;
17780c010503SBarry Smith 
177947c6ae99SBarry Smith   PetscFunctionBegin;
17800c010503SBarry Smith   PetscValidPointer(packer,2);
1781a4121054SBarry Smith   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1782a4121054SBarry Smith   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
178347c6ae99SBarry Smith   PetscFunctionReturn(0);
178447c6ae99SBarry Smith }
1785