xref: /petsc/src/dm/impls/composite/pack.c (revision f5f57ec0aaa1e5be6a6f3c2ef714a1a38ff0b64a)
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 
20*f5f57ec0SBarry Smith     Not available from Fortran
21*f5f57ec0SBarry Smith 
221b2093e4SBarry Smith     Notes: See DMSetApplicationContext() and DMGetApplicationContext() for how to get user information into
2347c6ae99SBarry Smith         this routine
2447c6ae99SBarry Smith 
2547c6ae99SBarry Smith @*/
267087cfbeSBarry Smith PetscErrorCode  DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt))
2747c6ae99SBarry Smith {
2847c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
2947c6ae99SBarry Smith 
3047c6ae99SBarry Smith   PetscFunctionBegin;
3147c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
3247c6ae99SBarry Smith   PetscFunctionReturn(0);
3347c6ae99SBarry Smith }
3447c6ae99SBarry Smith 
356bf464f9SBarry Smith PetscErrorCode  DMDestroy_Composite(DM dm)
3647c6ae99SBarry Smith {
3747c6ae99SBarry Smith   PetscErrorCode         ierr;
3847c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
3947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
4047c6ae99SBarry Smith 
4147c6ae99SBarry Smith   PetscFunctionBegin;
4247c6ae99SBarry Smith   next = com->next;
4347c6ae99SBarry Smith   while (next) {
4447c6ae99SBarry Smith     prev = next;
4547c6ae99SBarry Smith     next = next->next;
46fcfd50ebSBarry Smith     ierr = DMDestroy(&prev->dm);CHKERRQ(ierr);
4747c6ae99SBarry Smith     ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
4847c6ae99SBarry Smith     ierr = PetscFree(prev);CHKERRQ(ierr);
4947c6ae99SBarry Smith   }
50e10fd815SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL);CHKERRQ(ierr);
51435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
52435a35e8SMatthew G Knepley   ierr = PetscFree(com);CHKERRQ(ierr);
5347c6ae99SBarry Smith   PetscFunctionReturn(0);
5447c6ae99SBarry Smith }
5547c6ae99SBarry Smith 
567087cfbeSBarry Smith PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
5747c6ae99SBarry Smith {
5847c6ae99SBarry Smith   PetscErrorCode ierr;
5947c6ae99SBarry Smith   PetscBool      iascii;
6047c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
6147c6ae99SBarry Smith 
6247c6ae99SBarry Smith   PetscFunctionBegin;
63251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6447c6ae99SBarry Smith   if (iascii) {
6547c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
6647c6ae99SBarry Smith     PetscInt               i;
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix");CHKERRQ(ierr);
699ae5db72SJed Brown     ierr = PetscViewerASCIIPrintf(v,"  contains %D DMs\n",com->nDM);CHKERRQ(ierr);
7047c6ae99SBarry Smith     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
7147c6ae99SBarry Smith     for (i=0; lnk; lnk=lnk->next,i++) {
729ae5db72SJed Brown       ierr = PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
7347c6ae99SBarry Smith       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
7447c6ae99SBarry Smith       ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
7547c6ae99SBarry Smith       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
7647c6ae99SBarry Smith     }
7747c6ae99SBarry Smith     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
7847c6ae99SBarry Smith   }
7947c6ae99SBarry Smith   PetscFunctionReturn(0);
8047c6ae99SBarry Smith }
8147c6ae99SBarry Smith 
8247c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
837087cfbeSBarry Smith PetscErrorCode  DMSetUp_Composite(DM dm)
8447c6ae99SBarry Smith {
8547c6ae99SBarry Smith   PetscErrorCode         ierr;
8647c6ae99SBarry Smith   PetscInt               nprev = 0;
8747c6ae99SBarry Smith   PetscMPIInt            rank,size;
8847c6ae99SBarry Smith   DM_Composite           *com  = (DM_Composite*)dm->data;
8947c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
9047c6ae99SBarry Smith   PetscLayout            map;
9147c6ae99SBarry Smith 
9247c6ae99SBarry Smith   PetscFunctionBegin;
93ce94432eSBarry Smith   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
94ce94432eSBarry Smith   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm),&map);CHKERRQ(ierr);
9547c6ae99SBarry Smith   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
9647c6ae99SBarry Smith   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
9747c6ae99SBarry Smith   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
9847c6ae99SBarry Smith   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
9947c6ae99SBarry Smith   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
1000298fd71SBarry Smith   ierr = PetscLayoutGetRange(map,&com->rstart,NULL);CHKERRQ(ierr);
101fcfd50ebSBarry Smith   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
10247c6ae99SBarry Smith 
1039ae5db72SJed Brown   /* now set the rstart for each linked vector */
104ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
105ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
10647c6ae99SBarry Smith   while (next) {
10747c6ae99SBarry Smith     next->rstart  = nprev;
10806ebdd98SJed Brown     nprev        += next->n;
10947c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
110785e854fSJed Brown     ierr          = PetscMalloc1(size,&next->grstarts);CHKERRQ(ierr);
111ce94432eSBarry Smith     ierr          = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
11247c6ae99SBarry Smith     next          = next->next;
11347c6ae99SBarry Smith   }
11447c6ae99SBarry Smith   com->setup = PETSC_TRUE;
11547c6ae99SBarry Smith   PetscFunctionReturn(0);
11647c6ae99SBarry Smith }
11747c6ae99SBarry Smith 
11847c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
11947c6ae99SBarry Smith 
12073e31fe2SJed Brown /*@
12147c6ae99SBarry Smith     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
12247c6ae99SBarry Smith        representation.
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith     Not Collective
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith     Input Parameter:
12747c6ae99SBarry Smith .    dm - the packer object
12847c6ae99SBarry Smith 
12947c6ae99SBarry Smith     Output Parameter:
13047c6ae99SBarry Smith .     nDM - the number of DMs
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith     Level: beginner
13347c6ae99SBarry Smith 
13447c6ae99SBarry Smith @*/
1357087cfbeSBarry Smith PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
13647c6ae99SBarry Smith {
13747c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
1385fd66863SKarl Rupp 
13947c6ae99SBarry Smith   PetscFunctionBegin;
14047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
14147c6ae99SBarry Smith   *nDM = com->nDM;
14247c6ae99SBarry Smith   PetscFunctionReturn(0);
14347c6ae99SBarry Smith }
14447c6ae99SBarry Smith 
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith /*@C
14747c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
14847c6ae99SBarry Smith        representation.
14947c6ae99SBarry Smith 
15047c6ae99SBarry Smith     Collective on DMComposite
15147c6ae99SBarry Smith 
1529ae5db72SJed Brown     Input Parameters:
15347c6ae99SBarry Smith +    dm - the packer object
1549ae5db72SJed Brown -    gvec - the global vector
1559ae5db72SJed Brown 
1569ae5db72SJed Brown     Output Parameters:
1570298fd71SBarry Smith .    Vec* ... - the packed parallel vectors, NULL for those that are not needed
15847c6ae99SBarry Smith 
15947c6ae99SBarry Smith     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
16047c6ae99SBarry Smith 
161f73e5cebSJed Brown     Fortran Notes:
162f73e5cebSJed Brown 
163f73e5cebSJed Brown     Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
164f73e5cebSJed Brown     or use the alternative interface DMCompositeGetAccessArray().
165f73e5cebSJed Brown 
16647c6ae99SBarry Smith     Level: advanced
16747c6ae99SBarry Smith 
168f73e5cebSJed Brown .seealso: DMCompositeGetEntries(), DMCompositeScatter()
16947c6ae99SBarry Smith @*/
1707087cfbeSBarry Smith PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
17147c6ae99SBarry Smith {
17247c6ae99SBarry Smith   va_list                Argp;
17347c6ae99SBarry Smith   PetscErrorCode         ierr;
17447c6ae99SBarry Smith   struct DMCompositeLink *next;
17547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
1765edff71fSBarry Smith   PetscInt               readonly;
17747c6ae99SBarry Smith 
17847c6ae99SBarry Smith   PetscFunctionBegin;
17947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
18047c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
18147c6ae99SBarry Smith   next = com->next;
18247c6ae99SBarry Smith   if (!com->setup) {
183d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
18447c6ae99SBarry Smith   }
18547c6ae99SBarry Smith 
1865edff71fSBarry Smith   ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr);
18747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
18847c6ae99SBarry Smith   va_start(Argp,gvec);
18947c6ae99SBarry Smith   while (next) {
19047c6ae99SBarry Smith     Vec *vec;
19147c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
1929ae5db72SJed Brown     if (vec) {
1939ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,vec);CHKERRQ(ierr);
1945edff71fSBarry Smith       if (readonly) {
1955edff71fSBarry Smith         const PetscScalar *array;
1965edff71fSBarry Smith         ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
1975edff71fSBarry Smith         ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
1985edff71fSBarry Smith         ierr = VecLockPush(*vec);CHKERRQ(ierr);
1995edff71fSBarry Smith         ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
2005edff71fSBarry Smith       } else {
2015edff71fSBarry Smith         PetscScalar *array;
2029ae5db72SJed Brown         ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
2039ae5db72SJed Brown         ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr);
2049ae5db72SJed Brown         ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
20547c6ae99SBarry Smith       }
2065edff71fSBarry Smith     }
20747c6ae99SBarry Smith     next = next->next;
20847c6ae99SBarry Smith   }
20947c6ae99SBarry Smith   va_end(Argp);
21047c6ae99SBarry Smith   PetscFunctionReturn(0);
21147c6ae99SBarry Smith }
21247c6ae99SBarry Smith 
213f73e5cebSJed Brown /*@C
214f73e5cebSJed Brown     DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
215f73e5cebSJed Brown        representation.
216f73e5cebSJed Brown 
217f73e5cebSJed Brown     Collective on DMComposite
218f73e5cebSJed Brown 
219f73e5cebSJed Brown     Input Parameters:
220f73e5cebSJed Brown +    dm - the packer object
221f73e5cebSJed Brown .    pvec - packed vector
222f73e5cebSJed Brown .    nwanted - number of vectors wanted
2230298fd71SBarry Smith -    wanted - sorted array of vectors wanted, or NULL to get all vectors
224f73e5cebSJed Brown 
225f73e5cebSJed Brown     Output Parameters:
226f73e5cebSJed Brown .    vecs - array of requested global vectors (must be allocated)
227f73e5cebSJed Brown 
228f73e5cebSJed Brown     Notes: Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them
229f73e5cebSJed Brown 
230f73e5cebSJed Brown     Level: advanced
231f73e5cebSJed Brown 
232f73e5cebSJed Brown .seealso: DMCompositeGetAccess(), DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
233f73e5cebSJed Brown @*/
234f73e5cebSJed Brown PetscErrorCode  DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
235f73e5cebSJed Brown {
236f73e5cebSJed Brown   PetscErrorCode         ierr;
237f73e5cebSJed Brown   struct DMCompositeLink *link;
238f73e5cebSJed Brown   PetscInt               i,wnum;
239f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
240bee642f7SBarry Smith   PetscInt               readonly;
241f73e5cebSJed Brown 
242f73e5cebSJed Brown   PetscFunctionBegin;
243f73e5cebSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
244f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
245f73e5cebSJed Brown   if (!com->setup) {
246f73e5cebSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
247f73e5cebSJed Brown   }
248f73e5cebSJed Brown 
249bee642f7SBarry Smith   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
250f73e5cebSJed Brown   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
251f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
252f73e5cebSJed Brown       Vec v;
253f73e5cebSJed Brown       ierr = DMGetGlobalVector(link->dm,&v);CHKERRQ(ierr);
254bee642f7SBarry Smith       if (readonly) {
255bee642f7SBarry Smith         const PetscScalar *array;
256bee642f7SBarry Smith         ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr);
257bee642f7SBarry Smith         ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
258bee642f7SBarry Smith         ierr = VecLockPush(v);CHKERRQ(ierr);
259bee642f7SBarry Smith         ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr);
260bee642f7SBarry Smith       } else {
261bee642f7SBarry Smith         PetscScalar *array;
262f73e5cebSJed Brown         ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
263f73e5cebSJed Brown         ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr);
264f73e5cebSJed Brown         ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
265bee642f7SBarry Smith       }
266f73e5cebSJed Brown       vecs[wnum++] = v;
267f73e5cebSJed Brown     }
268f73e5cebSJed Brown   }
269f73e5cebSJed Brown   PetscFunctionReturn(0);
270f73e5cebSJed Brown }
271f73e5cebSJed Brown 
2727ac2b803SAlex Fikl /*@C
2737ac2b803SAlex Fikl     DMCompositeGetLocalAccessArray - Allows one to access the individual
2747ac2b803SAlex Fikl     packed vectors in their local representation.
2757ac2b803SAlex Fikl 
2767ac2b803SAlex Fikl     Collective on DMComposite.
2777ac2b803SAlex Fikl 
2787ac2b803SAlex Fikl     Input Parameters:
2797ac2b803SAlex Fikl +    dm - the packer object
2807ac2b803SAlex Fikl .    pvec - packed vector
2817ac2b803SAlex Fikl .    nwanted - number of vectors wanted
2827ac2b803SAlex Fikl -    wanted - sorted array of vectors wanted, or NULL to get all vectors
2837ac2b803SAlex Fikl 
2847ac2b803SAlex Fikl     Output Parameters:
2857ac2b803SAlex Fikl .    vecs - array of requested local vectors (must be allocated)
2867ac2b803SAlex Fikl 
2877ac2b803SAlex Fikl     Notes: Use DMCompositeRestoreLocalAccessArray() to return the vectors
2887ac2b803SAlex Fikl     when you no longer need them.
2897ac2b803SAlex Fikl 
2907ac2b803SAlex Fikl     Level: advanced
2917ac2b803SAlex Fikl 
2927ac2b803SAlex Fikl .seealso: DMCompositeRestoreLocalAccessArray(), DMCompositeGetAccess(),
2937ac2b803SAlex Fikl DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather()
2947ac2b803SAlex Fikl @*/
2957ac2b803SAlex Fikl PetscErrorCode  DMCompositeGetLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
2967ac2b803SAlex Fikl {
2977ac2b803SAlex Fikl   PetscErrorCode         ierr;
2987ac2b803SAlex Fikl   struct DMCompositeLink *link;
2997ac2b803SAlex Fikl   PetscInt               i,wnum;
3007ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite*)dm->data;
3017ac2b803SAlex Fikl   PetscInt               readonly;
3027ac2b803SAlex Fikl   PetscInt               nlocal = 0;
3037ac2b803SAlex Fikl 
3047ac2b803SAlex Fikl   PetscFunctionBegin;
3057ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3067ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
3077ac2b803SAlex Fikl   if (!com->setup) {
3087ac2b803SAlex Fikl     ierr = DMSetUp(dm);CHKERRQ(ierr);
3097ac2b803SAlex Fikl   }
3107ac2b803SAlex Fikl 
3117ac2b803SAlex Fikl   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
3127ac2b803SAlex Fikl   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
3137ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
3147ac2b803SAlex Fikl       Vec v;
3157ac2b803SAlex Fikl       ierr = DMGetLocalVector(link->dm,&v);CHKERRQ(ierr);
3167ac2b803SAlex Fikl       if (readonly) {
3177ac2b803SAlex Fikl         const PetscScalar *array;
3187ac2b803SAlex Fikl         ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr);
3197ac2b803SAlex Fikl         ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr);
3207ac2b803SAlex Fikl         ierr = VecLockPush(v);CHKERRQ(ierr);
3217ac2b803SAlex Fikl         ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr);
3227ac2b803SAlex Fikl       } else {
3237ac2b803SAlex Fikl         PetscScalar *array;
3247ac2b803SAlex Fikl         ierr = VecGetArray(pvec,&array);CHKERRQ(ierr);
3257ac2b803SAlex Fikl         ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr);
3267ac2b803SAlex Fikl         ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr);
3277ac2b803SAlex Fikl       }
3287ac2b803SAlex Fikl       vecs[wnum++] = v;
3297ac2b803SAlex Fikl     }
3307ac2b803SAlex Fikl 
3317ac2b803SAlex Fikl     nlocal += link->nlocal;
3327ac2b803SAlex Fikl   }
3337ac2b803SAlex Fikl 
3347ac2b803SAlex Fikl   PetscFunctionReturn(0);
3357ac2b803SAlex Fikl }
3367ac2b803SAlex Fikl 
33747c6ae99SBarry Smith /*@C
338aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
33947c6ae99SBarry Smith        representation.
34047c6ae99SBarry Smith 
34147c6ae99SBarry Smith     Collective on DMComposite
34247c6ae99SBarry Smith 
3439ae5db72SJed Brown     Input Parameters:
34447c6ae99SBarry Smith +    dm - the packer object
34547c6ae99SBarry Smith .    gvec - the global vector
3460298fd71SBarry Smith -    Vec* ... - the individual parallel vectors, NULL for those that are not needed
34747c6ae99SBarry Smith 
34847c6ae99SBarry Smith     Level: advanced
34947c6ae99SBarry Smith 
3509ae5db72SJed Brown .seealso  DMCompositeAddDM(), DMCreateGlobalVector(),
3516eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
352aa219208SBarry Smith          DMCompositeRestoreAccess(), DMCompositeGetAccess()
35347c6ae99SBarry Smith 
35447c6ae99SBarry Smith @*/
3557087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
35647c6ae99SBarry Smith {
35747c6ae99SBarry Smith   va_list                Argp;
35847c6ae99SBarry Smith   PetscErrorCode         ierr;
35947c6ae99SBarry Smith   struct DMCompositeLink *next;
36047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
3615edff71fSBarry Smith   PetscInt               readonly;
36247c6ae99SBarry Smith 
36347c6ae99SBarry Smith   PetscFunctionBegin;
36447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36547c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
36647c6ae99SBarry Smith   next = com->next;
36747c6ae99SBarry Smith   if (!com->setup) {
368d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
36947c6ae99SBarry Smith   }
37047c6ae99SBarry Smith 
3715edff71fSBarry Smith   ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr);
37247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
37347c6ae99SBarry Smith   va_start(Argp,gvec);
37447c6ae99SBarry Smith   while (next) {
37547c6ae99SBarry Smith     Vec *vec;
37647c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
3779ae5db72SJed Brown     if (vec) {
3789ae5db72SJed Brown       ierr = VecResetArray(*vec);CHKERRQ(ierr);
3795edff71fSBarry Smith       if (readonly) {
3805edff71fSBarry Smith         ierr = VecLockPop(*vec);CHKERRQ(ierr);
3815edff71fSBarry Smith       }
382bee642f7SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,vec);CHKERRQ(ierr);
38347c6ae99SBarry Smith     }
38447c6ae99SBarry Smith     next = next->next;
38547c6ae99SBarry Smith   }
38647c6ae99SBarry Smith   va_end(Argp);
38747c6ae99SBarry Smith   PetscFunctionReturn(0);
38847c6ae99SBarry Smith }
38947c6ae99SBarry Smith 
390f73e5cebSJed Brown /*@C
391f73e5cebSJed Brown     DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray()
392f73e5cebSJed Brown 
393f73e5cebSJed Brown     Collective on DMComposite
394f73e5cebSJed Brown 
395f73e5cebSJed Brown     Input Parameters:
396f73e5cebSJed Brown +    dm - the packer object
397f73e5cebSJed Brown .    pvec - packed vector
398f73e5cebSJed Brown .    nwanted - number of vectors wanted
3990298fd71SBarry Smith .    wanted - sorted array of vectors wanted, or NULL to get all vectors
400f73e5cebSJed Brown -    vecs - array of global vectors to return
401f73e5cebSJed Brown 
402f73e5cebSJed Brown     Level: advanced
403f73e5cebSJed Brown 
404f73e5cebSJed Brown .seealso: DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), DMCompositeScatter(), DMCompositeGather()
405f73e5cebSJed Brown @*/
406f73e5cebSJed Brown PetscErrorCode  DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
407f73e5cebSJed Brown {
408f73e5cebSJed Brown   PetscErrorCode         ierr;
409f73e5cebSJed Brown   struct DMCompositeLink *link;
410f73e5cebSJed Brown   PetscInt               i,wnum;
411f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
412bee642f7SBarry Smith   PetscInt               readonly;
413f73e5cebSJed Brown 
414f73e5cebSJed Brown   PetscFunctionBegin;
415f73e5cebSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
416f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
417f73e5cebSJed Brown   if (!com->setup) {
418f73e5cebSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
419f73e5cebSJed Brown   }
420f73e5cebSJed Brown 
421bee642f7SBarry Smith   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
422f73e5cebSJed Brown   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
423f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
424f73e5cebSJed Brown       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
425bee642f7SBarry Smith       if (readonly) {
426bee642f7SBarry Smith         ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr);
427bee642f7SBarry Smith       }
428f73e5cebSJed Brown       ierr = DMRestoreGlobalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
429f73e5cebSJed Brown       wnum++;
430f73e5cebSJed Brown     }
431f73e5cebSJed Brown   }
432f73e5cebSJed Brown   PetscFunctionReturn(0);
433f73e5cebSJed Brown }
434f73e5cebSJed Brown 
4357ac2b803SAlex Fikl /*@C
4367ac2b803SAlex Fikl     DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray().
4377ac2b803SAlex Fikl 
4387ac2b803SAlex Fikl     Collective on DMComposite.
4397ac2b803SAlex Fikl 
4407ac2b803SAlex Fikl     Input Parameters:
4417ac2b803SAlex Fikl +    dm - the packer object
4427ac2b803SAlex Fikl .    pvec - packed vector
4437ac2b803SAlex Fikl .    nwanted - number of vectors wanted
4447ac2b803SAlex Fikl .    wanted - sorted array of vectors wanted, or NULL to restore all vectors
4457ac2b803SAlex Fikl -    vecs - array of local vectors to return
4467ac2b803SAlex Fikl 
4477ac2b803SAlex Fikl     Level: advanced
4487ac2b803SAlex Fikl 
4497ac2b803SAlex Fikl     Notes:
4507ac2b803SAlex Fikl     nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray()
4517ac2b803SAlex Fikl     otherwise the call will fail.
4527ac2b803SAlex Fikl 
4537ac2b803SAlex Fikl .seealso: DMCompositeGetLocalAccessArray(), DMCompositeRestoreAccessArray(),
4547ac2b803SAlex Fikl DMCompositeRestoreAccess(), DMCompositeRestoreEntries(),
4557ac2b803SAlex Fikl DMCompositeScatter(), DMCompositeGather()
4567ac2b803SAlex Fikl @*/
4577ac2b803SAlex Fikl PetscErrorCode  DMCompositeRestoreLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
4587ac2b803SAlex Fikl {
4597ac2b803SAlex Fikl   PetscErrorCode         ierr;
4607ac2b803SAlex Fikl   struct DMCompositeLink *link;
4617ac2b803SAlex Fikl   PetscInt               i,wnum;
4627ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite*)dm->data;
4637ac2b803SAlex Fikl   PetscInt               readonly;
4647ac2b803SAlex Fikl 
4657ac2b803SAlex Fikl   PetscFunctionBegin;
4667ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4677ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
4687ac2b803SAlex Fikl   if (!com->setup) {
4697ac2b803SAlex Fikl     ierr = DMSetUp(dm);CHKERRQ(ierr);
4707ac2b803SAlex Fikl   }
4717ac2b803SAlex Fikl 
4727ac2b803SAlex Fikl   ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr);
4737ac2b803SAlex Fikl   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
4747ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
4757ac2b803SAlex Fikl       ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr);
4767ac2b803SAlex Fikl       if (readonly) {
4777ac2b803SAlex Fikl         ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr);
4787ac2b803SAlex Fikl       }
4797ac2b803SAlex Fikl       ierr = DMRestoreLocalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr);
4807ac2b803SAlex Fikl       wnum++;
4817ac2b803SAlex Fikl     }
4827ac2b803SAlex Fikl   }
4837ac2b803SAlex Fikl   PetscFunctionReturn(0);
4847ac2b803SAlex Fikl }
4857ac2b803SAlex Fikl 
48647c6ae99SBarry Smith /*@C
48747c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
48847c6ae99SBarry Smith 
48947c6ae99SBarry Smith     Collective on DMComposite
49047c6ae99SBarry Smith 
4919ae5db72SJed Brown     Input Parameters:
49247c6ae99SBarry Smith +    dm - the packer object
49347c6ae99SBarry Smith .    gvec - the global vector
4940298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for those that are not needed
49547c6ae99SBarry Smith 
49647c6ae99SBarry Smith     Level: advanced
49747c6ae99SBarry Smith 
4986f3c3dcfSJed Brown     Notes:
4996f3c3dcfSJed Brown     DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is
5006f3c3dcfSJed Brown     accessible from Fortran.
5016f3c3dcfSJed Brown 
5029ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
5036eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
50447c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
5056f3c3dcfSJed Brown          DMCompositeScatterArray()
50647c6ae99SBarry Smith 
50747c6ae99SBarry Smith @*/
5087087cfbeSBarry Smith PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
50947c6ae99SBarry Smith {
51047c6ae99SBarry Smith   va_list                Argp;
51147c6ae99SBarry Smith   PetscErrorCode         ierr;
51247c6ae99SBarry Smith   struct DMCompositeLink *next;
5138fd8f222SJed Brown   PetscInt               cnt;
51447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
51547c6ae99SBarry Smith 
51647c6ae99SBarry Smith   PetscFunctionBegin;
51747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
51847c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
51947c6ae99SBarry Smith   if (!com->setup) {
520d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
52147c6ae99SBarry Smith   }
52247c6ae99SBarry Smith 
52347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
52447c6ae99SBarry Smith   va_start(Argp,gvec);
5258fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
5269ae5db72SJed Brown     Vec local;
5279ae5db72SJed Brown     local = va_arg(Argp, Vec);
5289ae5db72SJed Brown     if (local) {
5299ae5db72SJed Brown       Vec               global;
5305edff71fSBarry Smith       const PetscScalar *array;
5319ae5db72SJed Brown       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
5329ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
5335edff71fSBarry Smith       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
5349ae5db72SJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
5359ae5db72SJed Brown       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
5369ae5db72SJed Brown       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
5375edff71fSBarry Smith       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
5389ae5db72SJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
5399ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
54047c6ae99SBarry Smith     }
54147c6ae99SBarry Smith   }
54247c6ae99SBarry Smith   va_end(Argp);
54347c6ae99SBarry Smith   PetscFunctionReturn(0);
54447c6ae99SBarry Smith }
54547c6ae99SBarry Smith 
5466f3c3dcfSJed Brown /*@
5476f3c3dcfSJed Brown     DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5486f3c3dcfSJed Brown 
5496f3c3dcfSJed Brown     Collective on DMComposite
5506f3c3dcfSJed Brown 
5516f3c3dcfSJed Brown     Input Parameters:
5526f3c3dcfSJed Brown +    dm - the packer object
5536f3c3dcfSJed Brown .    gvec - the global vector
5546f3c3dcfSJed Brown .    lvecs - array of local vectors, NULL for any that are not needed
5556f3c3dcfSJed Brown 
5566f3c3dcfSJed Brown     Level: advanced
5576f3c3dcfSJed Brown 
5586f3c3dcfSJed Brown     Note:
559907376e6SBarry Smith     This is a non-variadic alternative to DMCompositeScatter()
5606f3c3dcfSJed Brown 
5616f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector()
5626f3c3dcfSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
5636f3c3dcfSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
5646f3c3dcfSJed Brown 
5656f3c3dcfSJed Brown @*/
5666f3c3dcfSJed Brown PetscErrorCode  DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs)
5676f3c3dcfSJed Brown {
5686f3c3dcfSJed Brown   PetscErrorCode         ierr;
5696f3c3dcfSJed Brown   struct DMCompositeLink *next;
5706f3c3dcfSJed Brown   PetscInt               i;
5716f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
5726f3c3dcfSJed Brown 
5736f3c3dcfSJed Brown   PetscFunctionBegin;
5746f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5756f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
5766f3c3dcfSJed Brown   if (!com->setup) {
5776f3c3dcfSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
5786f3c3dcfSJed Brown   }
5796f3c3dcfSJed Brown 
5806f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
5816f3c3dcfSJed Brown   for (i=0,next=com->next; next; next=next->next,i++) {
5826f3c3dcfSJed Brown     if (lvecs[i]) {
5836f3c3dcfSJed Brown       Vec         global;
584c5d31e75SLisandro Dalcin       const PetscScalar *array;
5856f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
5866f3c3dcfSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
587c5d31e75SLisandro Dalcin       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
588c5d31e75SLisandro Dalcin       ierr = VecPlaceArray(global,(PetscScalar*)array+next->rstart);CHKERRQ(ierr);
5896f3c3dcfSJed Brown       ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
5906f3c3dcfSJed Brown       ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr);
591c5d31e75SLisandro Dalcin       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
5926f3c3dcfSJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
5936f3c3dcfSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
5946f3c3dcfSJed Brown     }
5956f3c3dcfSJed Brown   }
5966f3c3dcfSJed Brown   PetscFunctionReturn(0);
5976f3c3dcfSJed Brown }
5986f3c3dcfSJed Brown 
59947c6ae99SBarry Smith /*@C
60047c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
60147c6ae99SBarry Smith 
60247c6ae99SBarry Smith     Collective on DMComposite
60347c6ae99SBarry Smith 
60447c6ae99SBarry Smith     Input Parameter:
60547c6ae99SBarry Smith +    dm - the packer object
60647c6ae99SBarry Smith .    gvec - the global vector
607907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6080298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for any that are not needed
60947c6ae99SBarry Smith 
61047c6ae99SBarry Smith     Level: advanced
61147c6ae99SBarry Smith 
612*f5f57ec0SBarry Smith     Not available from Fortran, Fortran users can use DMCompositeGatherArray()
613*f5f57ec0SBarry Smith 
6149ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
6156eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
61647c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
61747c6ae99SBarry Smith 
61847c6ae99SBarry Smith @*/
6191dac896bSSatish Balay PetscErrorCode  DMCompositeGather(DM dm,InsertMode imode,Vec gvec,...)
62047c6ae99SBarry Smith {
62147c6ae99SBarry Smith   va_list                Argp;
62247c6ae99SBarry Smith   PetscErrorCode         ierr;
62347c6ae99SBarry Smith   struct DMCompositeLink *next;
62447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
6258fd8f222SJed Brown   PetscInt               cnt;
62647c6ae99SBarry Smith 
62747c6ae99SBarry Smith   PetscFunctionBegin;
62847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
62947c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
63047c6ae99SBarry Smith   if (!com->setup) {
631d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
63247c6ae99SBarry Smith   }
63347c6ae99SBarry Smith 
63447c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
6351dac896bSSatish Balay   va_start(Argp,gvec);
6368fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
6379ae5db72SJed Brown     Vec local;
6389ae5db72SJed Brown     local = va_arg(Argp, Vec);
6399ae5db72SJed Brown     if (local) {
64047c6ae99SBarry Smith       PetscScalar *array;
6419ae5db72SJed Brown       Vec         global;
6429ae5db72SJed Brown       PetscValidHeaderSpecific(local,VEC_CLASSID,cnt);
6439ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
6449ae5db72SJed Brown       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
6459ae5db72SJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
6469ae5db72SJed Brown       ierr = DMLocalToGlobalBegin(next->dm,local,imode,global);CHKERRQ(ierr);
6479ae5db72SJed Brown       ierr = DMLocalToGlobalEnd(next->dm,local,imode,global);CHKERRQ(ierr);
6489ae5db72SJed Brown       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
6499ae5db72SJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
6509ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
65147c6ae99SBarry Smith     }
65247c6ae99SBarry Smith   }
65347c6ae99SBarry Smith   va_end(Argp);
65447c6ae99SBarry Smith   PetscFunctionReturn(0);
65547c6ae99SBarry Smith }
65647c6ae99SBarry Smith 
6576f3c3dcfSJed Brown /*@
6586f3c3dcfSJed Brown     DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6596f3c3dcfSJed Brown 
6606f3c3dcfSJed Brown     Collective on DMComposite
6616f3c3dcfSJed Brown 
6626f3c3dcfSJed Brown     Input Parameter:
6636f3c3dcfSJed Brown +    dm - the packer object
6646f3c3dcfSJed Brown .    gvec - the global vector
665907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6666f3c3dcfSJed Brown -    lvecs - the individual sequential vectors, NULL for any that are not needed
6676f3c3dcfSJed Brown 
6686f3c3dcfSJed Brown     Level: advanced
6696f3c3dcfSJed Brown 
6706f3c3dcfSJed Brown     Notes:
6716f3c3dcfSJed Brown     This is a non-variadic alternative to DMCompositeGather().
6726f3c3dcfSJed Brown 
6736f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
6746f3c3dcfSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
6756f3c3dcfSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries(),
6766f3c3dcfSJed Brown @*/
6771dac896bSSatish Balay PetscErrorCode  DMCompositeGatherArray(DM dm,InsertMode imode,Vec gvec,Vec *lvecs)
6786f3c3dcfSJed Brown {
6796f3c3dcfSJed Brown   PetscErrorCode         ierr;
6806f3c3dcfSJed Brown   struct DMCompositeLink *next;
6816f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
6826f3c3dcfSJed Brown   PetscInt               i;
6836f3c3dcfSJed Brown 
6846f3c3dcfSJed Brown   PetscFunctionBegin;
6856f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6866f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
6876f3c3dcfSJed Brown   if (!com->setup) {
6886f3c3dcfSJed Brown     ierr = DMSetUp(dm);CHKERRQ(ierr);
6896f3c3dcfSJed Brown   }
6906f3c3dcfSJed Brown 
6916f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
6926f3c3dcfSJed Brown   for (next=com->next,i=0; next; next=next->next,i++) {
6936f3c3dcfSJed Brown     if (lvecs[i]) {
6946f3c3dcfSJed Brown       PetscScalar *array;
6956f3c3dcfSJed Brown       Vec         global;
6966f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
6976f3c3dcfSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
6986f3c3dcfSJed Brown       ierr = VecGetArray(gvec,&array);CHKERRQ(ierr);
6996f3c3dcfSJed Brown       ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr);
7006f3c3dcfSJed Brown       ierr = DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
7016f3c3dcfSJed Brown       ierr = DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global);CHKERRQ(ierr);
7026f3c3dcfSJed Brown       ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr);
7036f3c3dcfSJed Brown       ierr = VecResetArray(global);CHKERRQ(ierr);
7046f3c3dcfSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
7056f3c3dcfSJed Brown     }
7066f3c3dcfSJed Brown   }
7076f3c3dcfSJed Brown   PetscFunctionReturn(0);
7086f3c3dcfSJed Brown }
7096f3c3dcfSJed Brown 
710*f5f57ec0SBarry Smith /*@
711aa219208SBarry Smith     DMCompositeAddDM - adds a DM vector to a DMComposite
71247c6ae99SBarry Smith 
71347c6ae99SBarry Smith     Collective on DMComposite
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith     Input Parameter:
7169b52a9b5SPatrick Sanan +    dmc - the DMComposite (packer) object
717*f5f57ec0SBarry Smith -    dm - the DM object
71847c6ae99SBarry Smith 
71947c6ae99SBarry Smith     Level: advanced
72047c6ae99SBarry Smith 
7210c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
7226eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
72347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
72447c6ae99SBarry Smith 
72547c6ae99SBarry Smith @*/
7267087cfbeSBarry Smith PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
72747c6ae99SBarry Smith {
72847c6ae99SBarry Smith   PetscErrorCode         ierr;
72906ebdd98SJed Brown   PetscInt               n,nlocal;
73047c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
73106ebdd98SJed Brown   Vec                    global,local;
73247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmc->data;
73347c6ae99SBarry Smith 
73447c6ae99SBarry Smith   PetscFunctionBegin;
73547c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
73647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
73747c6ae99SBarry Smith   next = com->next;
738ce94432eSBarry Smith   if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
73947c6ae99SBarry Smith 
74047c6ae99SBarry Smith   /* create new link */
741b00a9115SJed Brown   ierr = PetscNew(&mine);CHKERRQ(ierr);
74247c6ae99SBarry Smith   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
74347c6ae99SBarry Smith   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
74447c6ae99SBarry Smith   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
74547c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
74606ebdd98SJed Brown   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
74706ebdd98SJed Brown   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
74806ebdd98SJed Brown   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
7498865f1eaSKarl Rupp 
75047c6ae99SBarry Smith   mine->n      = n;
75106ebdd98SJed Brown   mine->nlocal = nlocal;
75247c6ae99SBarry Smith   mine->dm     = dm;
7530298fd71SBarry Smith   mine->next   = NULL;
75447c6ae99SBarry Smith   com->n      += n;
7557ac2b803SAlex Fikl   com->nghost += nlocal;
75647c6ae99SBarry Smith 
75747c6ae99SBarry Smith   /* add to end of list */
7588865f1eaSKarl Rupp   if (!next) com->next = mine;
7598865f1eaSKarl Rupp   else {
76047c6ae99SBarry Smith     while (next->next) next = next->next;
76147c6ae99SBarry Smith     next->next = mine;
76247c6ae99SBarry Smith   }
76347c6ae99SBarry Smith   com->nDM++;
76447c6ae99SBarry Smith   com->nmine++;
76547c6ae99SBarry Smith   PetscFunctionReturn(0);
76647c6ae99SBarry Smith }
76747c6ae99SBarry Smith 
7689804daf3SBarry Smith #include <petscdraw.h>
76926887b52SJed Brown PETSC_EXTERN PetscErrorCode  VecView_MPI(Vec,PetscViewer);
7707087cfbeSBarry Smith PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
77147c6ae99SBarry Smith {
77247c6ae99SBarry Smith   DM                     dm;
77347c6ae99SBarry Smith   PetscErrorCode         ierr;
77447c6ae99SBarry Smith   struct DMCompositeLink *next;
77547c6ae99SBarry Smith   PetscBool              isdraw;
776cef07954SSatish Balay   DM_Composite           *com;
77747c6ae99SBarry Smith 
77847c6ae99SBarry Smith   PetscFunctionBegin;
779c688c046SMatthew G Knepley   ierr = VecGetDM(gvec, &dm);CHKERRQ(ierr);
780ce94432eSBarry Smith   if (!dm) SETERRQ(PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
78147c6ae99SBarry Smith   com  = (DM_Composite*)dm->data;
78247c6ae99SBarry Smith   next = com->next;
78347c6ae99SBarry Smith 
784251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
78547c6ae99SBarry Smith   if (!isdraw) {
78647c6ae99SBarry Smith     /* do I really want to call this? */
78747c6ae99SBarry Smith     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
78847c6ae99SBarry Smith   } else {
78947c6ae99SBarry Smith     PetscInt cnt = 0;
79047c6ae99SBarry Smith 
79147c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
79247c6ae99SBarry Smith     while (next) {
79347c6ae99SBarry Smith       Vec               vec;
794ca4278abSLisandro Dalcin       const PetscScalar *array;
79547c6ae99SBarry Smith       PetscInt          bs;
79647c6ae99SBarry Smith 
7979ae5db72SJed Brown       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
7989ae5db72SJed Brown       ierr = DMGetGlobalVector(next->dm,&vec);CHKERRQ(ierr);
799ca4278abSLisandro Dalcin       ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr);
800ca4278abSLisandro Dalcin       ierr = VecPlaceArray(vec,(PetscScalar*)array+next->rstart);CHKERRQ(ierr);
801ca4278abSLisandro Dalcin       ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr);
80247c6ae99SBarry Smith       ierr = VecView(vec,viewer);CHKERRQ(ierr);
8039ae5db72SJed Brown       ierr = VecResetArray(vec);CHKERRQ(ierr);
804ca4278abSLisandro Dalcin       ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
8059ae5db72SJed Brown       ierr = DMRestoreGlobalVector(next->dm,&vec);CHKERRQ(ierr);
80647c6ae99SBarry Smith       ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
80747c6ae99SBarry Smith       cnt += bs;
80847c6ae99SBarry Smith       next = next->next;
80947c6ae99SBarry Smith     }
81047c6ae99SBarry Smith     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
81147c6ae99SBarry Smith   }
81247c6ae99SBarry Smith   PetscFunctionReturn(0);
81347c6ae99SBarry Smith }
81447c6ae99SBarry Smith 
8157087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
81647c6ae99SBarry Smith {
81747c6ae99SBarry Smith   PetscErrorCode ierr;
81847c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
81947c6ae99SBarry Smith 
82047c6ae99SBarry Smith   PetscFunctionBegin;
82147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
822d7bf68aeSBarry Smith   ierr = DMSetUp(dm);CHKERRQ(ierr);
823ce94432eSBarry Smith   ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),com->n,com->N,gvec);CHKERRQ(ierr);
824c688c046SMatthew G Knepley   ierr = VecSetDM(*gvec, dm);CHKERRQ(ierr);
82547c6ae99SBarry Smith   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);CHKERRQ(ierr);
82647c6ae99SBarry Smith   PetscFunctionReturn(0);
82747c6ae99SBarry Smith }
82847c6ae99SBarry Smith 
8297087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
83047c6ae99SBarry Smith {
83147c6ae99SBarry Smith   PetscErrorCode ierr;
83247c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
83347c6ae99SBarry Smith 
83447c6ae99SBarry Smith   PetscFunctionBegin;
83547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
83647c6ae99SBarry Smith   if (!com->setup) {
837d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
83847c6ae99SBarry Smith   }
839f0e01b1fSVincent Le Chenadec   ierr = VecCreateSeq(PETSC_COMM_SELF,com->nghost,lvec);CHKERRQ(ierr);
840c688c046SMatthew G Knepley   ierr = VecSetDM(*lvec, dm);CHKERRQ(ierr);
84147c6ae99SBarry Smith   PetscFunctionReturn(0);
84247c6ae99SBarry Smith }
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith /*@C
8459ae5db72SJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
84647c6ae99SBarry Smith 
84706ebdd98SJed Brown     Collective on DM
84847c6ae99SBarry Smith 
84947c6ae99SBarry Smith     Input Parameter:
85047c6ae99SBarry Smith .    dm - the packer object
85147c6ae99SBarry Smith 
85247c6ae99SBarry Smith     Output Parameters:
8539ae5db72SJed Brown .    ltogs - the individual mappings for each packed vector. Note that this includes
8549ae5db72SJed Brown            all the ghost points that individual ghosted DMDA's may have.
85547c6ae99SBarry Smith 
85647c6ae99SBarry Smith     Level: advanced
85747c6ae99SBarry Smith 
85847c6ae99SBarry Smith     Notes:
8596eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
86047c6ae99SBarry Smith 
861*f5f57ec0SBarry Smith     Not available from Fortran
862*f5f57ec0SBarry Smith 
8639ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
86447c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
86547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
86647c6ae99SBarry Smith 
86747c6ae99SBarry Smith @*/
8687087cfbeSBarry Smith PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
86947c6ae99SBarry Smith {
87047c6ae99SBarry Smith   PetscErrorCode         ierr;
87147c6ae99SBarry Smith   PetscInt               i,*idx,n,cnt;
87247c6ae99SBarry Smith   struct DMCompositeLink *next;
87347c6ae99SBarry Smith   PetscMPIInt            rank;
87447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
87547c6ae99SBarry Smith 
87647c6ae99SBarry Smith   PetscFunctionBegin;
87747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
878728e99d6SJed Brown   ierr = DMSetUp(dm);CHKERRQ(ierr);
879854ce69bSBarry Smith   ierr = PetscMalloc1(com->nDM,ltogs);CHKERRQ(ierr);
88047c6ae99SBarry Smith   next = com->next;
881ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
88247c6ae99SBarry Smith 
88347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
88447c6ae99SBarry Smith   cnt = 0;
88547c6ae99SBarry Smith   while (next) {
8866eb61c8cSJed Brown     ISLocalToGlobalMapping ltog;
8876eb61c8cSJed Brown     PetscMPIInt            size;
88886994e45SJed Brown     const PetscInt         *suboff,*indices;
8896eb61c8cSJed Brown     Vec                    global;
89047c6ae99SBarry Smith 
8916eb61c8cSJed Brown     /* Get sub-DM global indices for each local dof */
8921411c6eeSJed Brown     ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
8936eb61c8cSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
89486994e45SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
895785e854fSJed Brown     ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr);
89647c6ae99SBarry Smith 
8976eb61c8cSJed Brown     /* Get the offsets for the sub-DM global vector */
8986eb61c8cSJed Brown     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
8996eb61c8cSJed Brown     ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
900ce94432eSBarry Smith     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)global),&size);CHKERRQ(ierr);
9016eb61c8cSJed Brown 
9026eb61c8cSJed Brown     /* Shift the sub-DM definition of the global space to the composite global space */
9036eb61c8cSJed Brown     for (i=0; i<n; i++) {
90486994e45SJed Brown       PetscInt subi = indices[i],lo = 0,hi = size,t;
9056eb61c8cSJed Brown       /* Binary search to find which rank owns subi */
9066eb61c8cSJed Brown       while (hi-lo > 1) {
9076eb61c8cSJed Brown         t = lo + (hi-lo)/2;
9086eb61c8cSJed Brown         if (suboff[t] > subi) hi = t;
9096eb61c8cSJed Brown         else                  lo = t;
9106eb61c8cSJed Brown       }
9116eb61c8cSJed Brown       idx[i] = subi - suboff[lo] + next->grstarts[lo];
9126eb61c8cSJed Brown     }
91386994e45SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
914f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
9156eb61c8cSJed Brown     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
91647c6ae99SBarry Smith     next = next->next;
91747c6ae99SBarry Smith     cnt++;
91847c6ae99SBarry Smith   }
91947c6ae99SBarry Smith   PetscFunctionReturn(0);
92047c6ae99SBarry Smith }
92147c6ae99SBarry Smith 
92287c85e80SJed Brown /*@C
9239ae5db72SJed Brown    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
92487c85e80SJed Brown 
92587c85e80SJed Brown    Not Collective
92687c85e80SJed Brown 
92787c85e80SJed Brown    Input Arguments:
92887c85e80SJed Brown . dm - composite DM
92987c85e80SJed Brown 
93087c85e80SJed Brown    Output Arguments:
93187c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
93287c85e80SJed Brown 
93387c85e80SJed Brown    Level: intermediate
93487c85e80SJed Brown 
93587c85e80SJed Brown    Notes:
93687c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
93787c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
9389ae5db72SJed Brown    scatter to a composite local vector.  The user should not typically need to know which is being done.
93987c85e80SJed Brown 
94087c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
94187c85e80SJed Brown 
94287c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
94387c85e80SJed Brown 
94487c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
94587c85e80SJed Brown 
946*f5f57ec0SBarry Smith    Not available from Fortran
947*f5f57ec0SBarry Smith 
94887c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
94987c85e80SJed Brown @*/
9507087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
95187c85e80SJed Brown {
95287c85e80SJed Brown   PetscErrorCode         ierr;
95387c85e80SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
95487c85e80SJed Brown   struct DMCompositeLink *link;
95587c85e80SJed Brown   PetscInt               cnt,start;
95687c85e80SJed Brown 
95787c85e80SJed Brown   PetscFunctionBegin;
95887c85e80SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
95987c85e80SJed Brown   PetscValidPointer(is,2);
960785e854fSJed Brown   ierr = PetscMalloc1(com->nmine,is);CHKERRQ(ierr);
96106ebdd98SJed Brown   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
962520db06cSJed Brown     PetscInt bs;
9639ae5db72SJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
9641411c6eeSJed Brown     ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
965520db06cSJed Brown     ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
966520db06cSJed Brown   }
96787c85e80SJed Brown   PetscFunctionReturn(0);
96887c85e80SJed Brown }
96987c85e80SJed Brown 
97047c6ae99SBarry Smith /*@C
97147c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
97247c6ae99SBarry Smith 
97347c6ae99SBarry Smith     Collective on DMComposite
97447c6ae99SBarry Smith 
97547c6ae99SBarry Smith     Input Parameter:
97647c6ae99SBarry Smith .    dm - the packer object
97747c6ae99SBarry Smith 
97847c6ae99SBarry Smith     Output Parameters:
97947c6ae99SBarry Smith .    is - the array of index sets
98047c6ae99SBarry Smith 
98147c6ae99SBarry Smith     Level: advanced
98247c6ae99SBarry Smith 
98347c6ae99SBarry Smith     Notes:
98447c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
98547c6ae99SBarry Smith 
98647c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
98747c6ae99SBarry Smith 
9886eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
9896eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
9906eb61c8cSJed Brown        indices.
99147c6ae99SBarry Smith 
992f3cb0f7eSJed Brown     Fortran Notes:
993f3cb0f7eSJed Brown 
994f3cb0f7eSJed Brown        The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM().
995f3cb0f7eSJed Brown 
9969ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
99747c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
99847c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
99947c6ae99SBarry Smith 
100047c6ae99SBarry Smith @*/
10017087cfbeSBarry Smith PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
100247c6ae99SBarry Smith {
100347c6ae99SBarry Smith   PetscErrorCode         ierr;
100466bb578eSMark Adams   PetscInt               cnt = 0;
100547c6ae99SBarry Smith   struct DMCompositeLink *next;
100647c6ae99SBarry Smith   PetscMPIInt            rank;
100747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
100847c6ae99SBarry Smith 
100947c6ae99SBarry Smith   PetscFunctionBegin;
101047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1011854ce69bSBarry Smith   ierr = PetscMalloc1(com->nDM,is);CHKERRQ(ierr);
101247c6ae99SBarry Smith   next = com->next;
1013ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
101447c6ae99SBarry Smith 
101547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
101647c6ae99SBarry Smith   while (next) {
101766bb578eSMark Adams     ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt]);CHKERRQ(ierr);
10180f21e855SMatthew G. Knepley     if (dm->prob) {
101965c226d8SMatthew G. Knepley       MatNullSpace space;
102065c226d8SMatthew G. Knepley       Mat          pmat;
10210f21e855SMatthew G. Knepley       PetscObject  disc;
10220f21e855SMatthew G. Knepley       PetscInt     Nf;
102365c226d8SMatthew G. Knepley 
10242764a2aaSMatthew G. Knepley       ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
1025f24dd8d2SMatthew G. Knepley       if (cnt < Nf) {
10262764a2aaSMatthew G. Knepley         ierr = PetscDSGetDiscretization(dm->prob, cnt, &disc);CHKERRQ(ierr);
10270f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nullspace", (PetscObject*) &space);CHKERRQ(ierr);
1028aac2dd2dSMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space);CHKERRQ(ierr);}
10290f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space);CHKERRQ(ierr);
1030aac2dd2dSMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space);CHKERRQ(ierr);}
10310f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr);
1032aac2dd2dSMatthew G. Knepley         if (pmat) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat);CHKERRQ(ierr);}
103365c226d8SMatthew G. Knepley       }
1034f24dd8d2SMatthew G. Knepley     }
103547c6ae99SBarry Smith     cnt++;
103647c6ae99SBarry Smith     next = next->next;
103747c6ae99SBarry Smith   }
103847c6ae99SBarry Smith   PetscFunctionReturn(0);
103947c6ae99SBarry Smith }
104047c6ae99SBarry Smith 
104121c9b008SJed Brown PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields)
10424d343eeaSMatthew G Knepley {
10434d343eeaSMatthew G Knepley   PetscInt       nDM;
10444d343eeaSMatthew G Knepley   DM             *dms;
10454d343eeaSMatthew G Knepley   PetscInt       i;
10464d343eeaSMatthew G Knepley   PetscErrorCode ierr;
10474d343eeaSMatthew G Knepley 
10484d343eeaSMatthew G Knepley   PetscFunctionBegin;
10494d343eeaSMatthew G Knepley   ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
10508865f1eaSKarl Rupp   if (numFields) *numFields = nDM;
10514d343eeaSMatthew G Knepley   ierr = DMCompositeGetGlobalISs(dm, fields);CHKERRQ(ierr);
10524d343eeaSMatthew G Knepley   if (fieldNames) {
1053785e854fSJed Brown     ierr = PetscMalloc1(nDM, &dms);CHKERRQ(ierr);
1054785e854fSJed Brown     ierr = PetscMalloc1(nDM, fieldNames);CHKERRQ(ierr);
10554d343eeaSMatthew G Knepley     ierr = DMCompositeGetEntriesArray(dm, dms);CHKERRQ(ierr);
10564d343eeaSMatthew G Knepley     for (i=0; i<nDM; i++) {
10574d343eeaSMatthew G Knepley       char       buf[256];
10584d343eeaSMatthew G Knepley       const char *splitname;
10594d343eeaSMatthew G Knepley 
10604d343eeaSMatthew G Knepley       /* Split naming precedence: object name, prefix, number */
10614d343eeaSMatthew G Knepley       splitname = ((PetscObject) dm)->name;
10624d343eeaSMatthew G Knepley       if (!splitname) {
10634d343eeaSMatthew G Knepley         ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr);
10644d343eeaSMatthew G Knepley         if (splitname) {
10654d343eeaSMatthew G Knepley           size_t len;
10668caf3d72SBarry Smith           ierr                 = PetscStrncpy(buf,splitname,sizeof(buf));CHKERRQ(ierr);
10678caf3d72SBarry Smith           buf[sizeof(buf) - 1] = 0;
10684d343eeaSMatthew G Knepley           ierr                 = PetscStrlen(buf,&len);CHKERRQ(ierr);
10694d343eeaSMatthew G Knepley           if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
10704d343eeaSMatthew G Knepley           splitname = buf;
10714d343eeaSMatthew G Knepley         }
10724d343eeaSMatthew G Knepley       }
10734d343eeaSMatthew G Knepley       if (!splitname) {
10748caf3d72SBarry Smith         ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",i);CHKERRQ(ierr);
10754d343eeaSMatthew G Knepley         splitname = buf;
10764d343eeaSMatthew G Knepley       }
107721c9b008SJed Brown       ierr = PetscStrallocpy(splitname,&(*fieldNames)[i]);CHKERRQ(ierr);
10784d343eeaSMatthew G Knepley     }
10794d343eeaSMatthew G Knepley     ierr = PetscFree(dms);CHKERRQ(ierr);
10804d343eeaSMatthew G Knepley   }
10814d343eeaSMatthew G Knepley   PetscFunctionReturn(0);
10824d343eeaSMatthew G Knepley }
10834d343eeaSMatthew G Knepley 
1084e7c4fc90SDmitry Karpeev /*
1085e7c4fc90SDmitry Karpeev  This could take over from DMCreateFieldIS(), as it is more general,
10860298fd71SBarry Smith  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1087e7c4fc90SDmitry Karpeev  At this point it's probably best to be less intrusive, however.
1088e7c4fc90SDmitry Karpeev  */
108916621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
1090e7c4fc90SDmitry Karpeev {
1091e7c4fc90SDmitry Karpeev   PetscInt       nDM;
1092e7c4fc90SDmitry Karpeev   PetscInt       i;
1093e7c4fc90SDmitry Karpeev   PetscErrorCode ierr;
1094e7c4fc90SDmitry Karpeev 
1095e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
1096e7c4fc90SDmitry Karpeev   ierr = DMCreateFieldIS_Composite(dm, len, namelist, islist);CHKERRQ(ierr);
1097e7c4fc90SDmitry Karpeev   if (dmlist) {
1098e7c4fc90SDmitry Karpeev     ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr);
1099785e854fSJed Brown     ierr = PetscMalloc1(nDM, dmlist);CHKERRQ(ierr);
1100e7c4fc90SDmitry Karpeev     ierr = DMCompositeGetEntriesArray(dm, *dmlist);CHKERRQ(ierr);
1101e7c4fc90SDmitry Karpeev     for (i=0; i<nDM; i++) {
1102e7c4fc90SDmitry Karpeev       ierr = PetscObjectReference((PetscObject)((*dmlist)[i]));CHKERRQ(ierr);
1103e7c4fc90SDmitry Karpeev     }
1104e7c4fc90SDmitry Karpeev   }
1105e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
1106e7c4fc90SDmitry Karpeev }
1107e7c4fc90SDmitry Karpeev 
1108e7c4fc90SDmitry Karpeev 
1109e7c4fc90SDmitry Karpeev 
111047c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
111147c6ae99SBarry Smith /*@C
11129ae5db72SJed Brown     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
111347c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
111447c6ae99SBarry Smith 
111547c6ae99SBarry Smith     Not Collective
111647c6ae99SBarry Smith 
111747c6ae99SBarry Smith     Input Parameter:
111847c6ae99SBarry Smith .    dm - the packer object
111947c6ae99SBarry Smith 
112047c6ae99SBarry Smith     Output Parameter:
11219ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
112247c6ae99SBarry Smith 
112347c6ae99SBarry Smith     Level: advanced
112447c6ae99SBarry Smith 
1125*f5f57ec0SBarry Smith     Not available from Fortran
1126*f5f57ec0SBarry Smith 
11279ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
11286eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
112947c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
113047c6ae99SBarry Smith 
113147c6ae99SBarry Smith @*/
11327087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
113347c6ae99SBarry Smith {
113447c6ae99SBarry Smith   va_list                Argp;
113547c6ae99SBarry Smith   PetscErrorCode         ierr;
113647c6ae99SBarry Smith   struct DMCompositeLink *next;
113747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
113847c6ae99SBarry Smith 
113947c6ae99SBarry Smith   PetscFunctionBegin;
114047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
114147c6ae99SBarry Smith   next = com->next;
114247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
114347c6ae99SBarry Smith   va_start(Argp,dm);
114447c6ae99SBarry Smith   while (next) {
114547c6ae99SBarry Smith     Vec *vec;
114647c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
114706930112SJed Brown     if (vec) {ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);}
114847c6ae99SBarry Smith     next = next->next;
114947c6ae99SBarry Smith   }
115047c6ae99SBarry Smith   va_end(Argp);
115147c6ae99SBarry Smith   PetscFunctionReturn(0);
115247c6ae99SBarry Smith }
115347c6ae99SBarry Smith 
115447c6ae99SBarry Smith /*@C
11559ae5db72SJed Brown     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
115647c6ae99SBarry Smith 
115747c6ae99SBarry Smith     Not Collective
115847c6ae99SBarry Smith 
115947c6ae99SBarry Smith     Input Parameter:
116047c6ae99SBarry Smith .    dm - the packer object
116147c6ae99SBarry Smith 
116247c6ae99SBarry Smith     Output Parameter:
11639ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
116447c6ae99SBarry Smith 
116547c6ae99SBarry Smith     Level: advanced
116647c6ae99SBarry Smith 
1167*f5f57ec0SBarry Smith     Not available from Fortran
1168*f5f57ec0SBarry Smith 
11699ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(),
11706eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
117147c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
117247c6ae99SBarry Smith 
117347c6ae99SBarry Smith @*/
11747087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
117547c6ae99SBarry Smith {
117647c6ae99SBarry Smith   va_list                Argp;
117747c6ae99SBarry Smith   PetscErrorCode         ierr;
117847c6ae99SBarry Smith   struct DMCompositeLink *next;
117947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
118047c6ae99SBarry Smith 
118147c6ae99SBarry Smith   PetscFunctionBegin;
118247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
118347c6ae99SBarry Smith   next = com->next;
118447c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
118547c6ae99SBarry Smith   va_start(Argp,dm);
118647c6ae99SBarry Smith   while (next) {
118747c6ae99SBarry Smith     Vec *vec;
118847c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
118906930112SJed Brown     if (vec) {ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);}
119047c6ae99SBarry Smith     next = next->next;
119147c6ae99SBarry Smith   }
119247c6ae99SBarry Smith   va_end(Argp);
119347c6ae99SBarry Smith   PetscFunctionReturn(0);
119447c6ae99SBarry Smith }
119547c6ae99SBarry Smith 
119647c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
119747c6ae99SBarry Smith /*@C
11989ae5db72SJed Brown     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
119947c6ae99SBarry Smith 
120047c6ae99SBarry Smith     Not Collective
120147c6ae99SBarry Smith 
120247c6ae99SBarry Smith     Input Parameter:
120347c6ae99SBarry Smith .    dm - the packer object
120447c6ae99SBarry Smith 
120547c6ae99SBarry Smith     Output Parameter:
12069ae5db72SJed Brown .   DM ... - the individual entries (DMs)
120747c6ae99SBarry Smith 
120847c6ae99SBarry Smith     Level: advanced
120947c6ae99SBarry Smith 
1210*f5f57ec0SBarry Smith     Fortran Notes: Available as DMCompositeGetEntries() for one output DM, DMCompositeGetEntries2() for 2, etc
1211*f5f57ec0SBarry Smith 
12122fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray()
12136eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
121447c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
121547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
121647c6ae99SBarry Smith 
121747c6ae99SBarry Smith @*/
12187087cfbeSBarry Smith PetscErrorCode  DMCompositeGetEntries(DM dm,...)
121947c6ae99SBarry Smith {
122047c6ae99SBarry Smith   va_list                Argp;
122147c6ae99SBarry Smith   struct DMCompositeLink *next;
122247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
122347c6ae99SBarry Smith 
122447c6ae99SBarry Smith   PetscFunctionBegin;
122547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
122647c6ae99SBarry Smith   next = com->next;
122747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
122847c6ae99SBarry Smith   va_start(Argp,dm);
122947c6ae99SBarry Smith   while (next) {
123047c6ae99SBarry Smith     DM *dmn;
123147c6ae99SBarry Smith     dmn = va_arg(Argp, DM*);
12329ae5db72SJed Brown     if (dmn) *dmn = next->dm;
123347c6ae99SBarry Smith     next = next->next;
123447c6ae99SBarry Smith   }
123547c6ae99SBarry Smith   va_end(Argp);
123647c6ae99SBarry Smith   PetscFunctionReturn(0);
123747c6ae99SBarry Smith }
123847c6ae99SBarry Smith 
1239dbab29e1SMark F. Adams /*@C
12402fa5ba8aSJed Brown     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
12412fa5ba8aSJed Brown 
12422fa5ba8aSJed Brown     Not Collective
12432fa5ba8aSJed Brown 
12442fa5ba8aSJed Brown     Input Parameter:
1245907376e6SBarry Smith .    dm - the packer object
1246907376e6SBarry Smith 
1247907376e6SBarry Smith     Output Parameter:
1248907376e6SBarry Smith .    dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs
12492fa5ba8aSJed Brown 
12502fa5ba8aSJed Brown     Level: advanced
12512fa5ba8aSJed Brown 
12522fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries()
12532fa5ba8aSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
12542fa5ba8aSJed Brown          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
12552fa5ba8aSJed Brown          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
12562fa5ba8aSJed Brown 
12572fa5ba8aSJed Brown @*/
12582fa5ba8aSJed Brown PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
12592fa5ba8aSJed Brown {
12602fa5ba8aSJed Brown   struct DMCompositeLink *next;
12612fa5ba8aSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
12622fa5ba8aSJed Brown   PetscInt               i;
12632fa5ba8aSJed Brown 
12642fa5ba8aSJed Brown   PetscFunctionBegin;
12652fa5ba8aSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12662fa5ba8aSJed Brown   /* loop over packed objects, handling one at at time */
12672fa5ba8aSJed Brown   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
12682fa5ba8aSJed Brown   PetscFunctionReturn(0);
12692fa5ba8aSJed Brown }
12702fa5ba8aSJed Brown 
1271e10fd815SStefano Zampini typedef struct {
1272e10fd815SStefano Zampini   DM          dm;
1273e10fd815SStefano Zampini   PetscViewer *subv;
1274e10fd815SStefano Zampini   Vec         *vecs;
1275e10fd815SStefano Zampini } GLVisViewerCtx;
1276e10fd815SStefano Zampini 
1277e10fd815SStefano Zampini static PetscErrorCode  DestroyGLVisViewerCtx_Private(void *vctx)
1278e10fd815SStefano Zampini {
1279e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
1280e10fd815SStefano Zampini   PetscInt       i,n;
1281e10fd815SStefano Zampini   PetscErrorCode ierr;
1282e10fd815SStefano Zampini 
1283e10fd815SStefano Zampini   PetscFunctionBegin;
1284e10fd815SStefano Zampini   ierr = DMCompositeGetNumberDM(ctx->dm,&n);CHKERRQ(ierr);
1285e10fd815SStefano Zampini   for (i = 0; i < n; i++) {
1286e10fd815SStefano Zampini     ierr = PetscViewerDestroy(&ctx->subv[i]);CHKERRQ(ierr);
1287e10fd815SStefano Zampini   }
1288e10fd815SStefano Zampini   ierr = PetscFree2(ctx->subv,ctx->vecs);CHKERRQ(ierr);
1289e10fd815SStefano Zampini   ierr = DMDestroy(&ctx->dm);CHKERRQ(ierr);
1290e10fd815SStefano Zampini   ierr = PetscFree(ctx);CHKERRQ(ierr);
1291e10fd815SStefano Zampini   PetscFunctionReturn(0);
1292e10fd815SStefano Zampini }
1293e10fd815SStefano Zampini 
1294e10fd815SStefano Zampini static PetscErrorCode  DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
1295e10fd815SStefano Zampini {
1296e10fd815SStefano Zampini   Vec            X = (Vec)oX;
1297e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
1298e10fd815SStefano Zampini   PetscInt       i,n,cumf;
1299e10fd815SStefano Zampini   PetscErrorCode ierr;
1300e10fd815SStefano Zampini 
1301e10fd815SStefano Zampini   PetscFunctionBegin;
1302e10fd815SStefano Zampini   ierr = DMCompositeGetNumberDM(ctx->dm,&n);CHKERRQ(ierr);
1303e10fd815SStefano Zampini   ierr = DMCompositeGetAccessArray(ctx->dm,X,n,NULL,ctx->vecs);CHKERRQ(ierr);
1304e10fd815SStefano Zampini   for (i = 0, cumf = 0; i < n; i++) {
1305e10fd815SStefano Zampini     PetscErrorCode (*g2l)(PetscObject,PetscInt,PetscObject[],void*);
1306e10fd815SStefano Zampini     void           *fctx;
1307e10fd815SStefano Zampini     PetscInt       nfi;
1308e10fd815SStefano Zampini 
1309e10fd815SStefano Zampini     ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nfi,NULL,NULL,&g2l,NULL,&fctx);CHKERRQ(ierr);
1310e10fd815SStefano Zampini     if (!nfi) continue;
1311e10fd815SStefano Zampini     if (g2l) {
1312e10fd815SStefano Zampini       ierr = (*g2l)((PetscObject)ctx->vecs[i],nfi,oXfield+cumf,fctx);CHKERRQ(ierr);
1313e10fd815SStefano Zampini     } else {
1314e10fd815SStefano Zampini       ierr = VecCopy(ctx->vecs[i],(Vec)(oXfield[cumf]));CHKERRQ(ierr);
1315e10fd815SStefano Zampini     }
1316e10fd815SStefano Zampini     cumf += nfi;
1317e10fd815SStefano Zampini   }
1318e10fd815SStefano Zampini   ierr = DMCompositeRestoreAccessArray(ctx->dm,X,n,NULL,ctx->vecs);CHKERRQ(ierr);
1319e10fd815SStefano Zampini   PetscFunctionReturn(0);
1320e10fd815SStefano Zampini }
1321e10fd815SStefano Zampini 
1322e10fd815SStefano Zampini static PetscErrorCode  DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer)
1323e10fd815SStefano Zampini {
1324e10fd815SStefano Zampini   DM             dm = (DM)odm, *dms;
1325e10fd815SStefano Zampini   Vec            *Ufds;
1326e10fd815SStefano Zampini   GLVisViewerCtx *ctx;
1327e10fd815SStefano Zampini   PetscInt       i,n,tnf,*sdim;
1328e10fd815SStefano Zampini   char           **fecs;
1329e10fd815SStefano Zampini   PetscErrorCode ierr;
1330e10fd815SStefano Zampini 
1331e10fd815SStefano Zampini   PetscFunctionBegin;
1332e10fd815SStefano Zampini   ierr = PetscNew(&ctx);CHKERRQ(ierr);
1333e10fd815SStefano Zampini   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
1334e10fd815SStefano Zampini   ctx->dm = dm;
1335e10fd815SStefano Zampini   ierr = DMCompositeGetNumberDM(dm,&n);CHKERRQ(ierr);
1336e10fd815SStefano Zampini   ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr);
1337e10fd815SStefano Zampini   ierr = DMCompositeGetEntriesArray(dm,dms);CHKERRQ(ierr);
1338e10fd815SStefano Zampini   ierr = PetscMalloc2(n,&ctx->subv,n,&ctx->vecs);CHKERRQ(ierr);
1339e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1340e10fd815SStefano Zampini     PetscInt nf;
1341e10fd815SStefano Zampini 
1342e10fd815SStefano Zampini     ierr = PetscViewerCreate(PetscObjectComm(odm),&ctx->subv[i]);CHKERRQ(ierr);
1343e10fd815SStefano Zampini     ierr = PetscViewerSetType(ctx->subv[i],PETSCVIEWERGLVIS);CHKERRQ(ierr);
1344e10fd815SStefano Zampini     ierr = PetscViewerGLVisSetDM_Private(ctx->subv[i],(PetscObject)dms[i]);CHKERRQ(ierr);
1345e10fd815SStefano Zampini     ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1346e10fd815SStefano Zampini     tnf += nf;
1347e10fd815SStefano Zampini   }
1348e10fd815SStefano Zampini   ierr = PetscFree(dms);CHKERRQ(ierr);
1349e10fd815SStefano Zampini   ierr = PetscMalloc3(tnf,&fecs,tnf,&sdim,tnf,&Ufds);CHKERRQ(ierr);
1350e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1351e10fd815SStefano Zampini     PetscInt   *sd,nf,f;
1352e10fd815SStefano Zampini     const char **fec;
1353e10fd815SStefano Zampini     Vec        *Uf;
1354e10fd815SStefano Zampini 
1355e10fd815SStefano Zampini     ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,&fec,&sd,NULL,(PetscObject**)&Uf,NULL);CHKERRQ(ierr);
1356e10fd815SStefano Zampini     for (f = 0; f < nf; f++) {
1357e10fd815SStefano Zampini       ierr = PetscStrallocpy(fec[f],&fecs[tnf+f]);CHKERRQ(ierr);
1358e10fd815SStefano Zampini       Ufds[tnf+f] = Uf[f];
1359e10fd815SStefano Zampini       sdim[tnf+f] = sd[f];
1360e10fd815SStefano Zampini     }
1361e10fd815SStefano Zampini     tnf += nf;
1362e10fd815SStefano Zampini   }
1363e10fd815SStefano Zampini   ierr = PetscViewerGLVisSetFields(viewer,tnf,(const char**)fecs,sdim,DMCompositeSampleGLVisFields_Private,(PetscObject*)Ufds,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
1364e10fd815SStefano Zampini   for (i = 0; i < tnf; i++) {
1365e10fd815SStefano Zampini     ierr = PetscFree(fecs[i]);CHKERRQ(ierr);
1366e10fd815SStefano Zampini   }
1367e10fd815SStefano Zampini   ierr = PetscFree3(fecs,sdim,Ufds);CHKERRQ(ierr);
1368e10fd815SStefano Zampini   PetscFunctionReturn(0);
1369e10fd815SStefano Zampini }
1370e10fd815SStefano Zampini 
13717087cfbeSBarry Smith PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
137247c6ae99SBarry Smith {
137347c6ae99SBarry Smith   PetscErrorCode         ierr;
137447c6ae99SBarry Smith   struct DMCompositeLink *next;
137547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmi->data;
137647c6ae99SBarry Smith   DM                     dm;
137747c6ae99SBarry Smith 
137847c6ae99SBarry Smith   PetscFunctionBegin;
137947c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1380ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
1381ce94432eSBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
1382ce94432eSBarry Smith   }
13832ce3a92bSJed Brown   ierr = DMSetUp(dmi);CHKERRQ(ierr);
138447c6ae99SBarry Smith   next = com->next;
138547c6ae99SBarry Smith   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
138647c6ae99SBarry Smith 
138747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
138847c6ae99SBarry Smith   while (next) {
138947c6ae99SBarry Smith     ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
139047c6ae99SBarry Smith     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
139147c6ae99SBarry Smith     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
139247c6ae99SBarry Smith     next = next->next;
139347c6ae99SBarry Smith   }
139447c6ae99SBarry Smith   PetscFunctionReturn(0);
139547c6ae99SBarry Smith }
139647c6ae99SBarry Smith 
139714354c39SJed Brown PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
139814354c39SJed Brown {
139914354c39SJed Brown   PetscErrorCode         ierr;
140014354c39SJed Brown   struct DMCompositeLink *next;
140114354c39SJed Brown   DM_Composite           *com = (DM_Composite*)dmi->data;
140214354c39SJed Brown   DM                     dm;
140314354c39SJed Brown 
140414354c39SJed Brown   PetscFunctionBegin;
140514354c39SJed Brown   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
14062ce3a92bSJed Brown   ierr = DMSetUp(dmi);CHKERRQ(ierr);
14072ee06e3bSJed Brown   if (comm == MPI_COMM_NULL) {
140825296bd5SBarry Smith     ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr);
140925296bd5SBarry Smith   }
141014354c39SJed Brown   next = com->next;
141114354c39SJed Brown   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
141214354c39SJed Brown 
141314354c39SJed Brown   /* loop over packed objects, handling one at at time */
141414354c39SJed Brown   while (next) {
141514354c39SJed Brown     ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr);
141614354c39SJed Brown     ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
141714354c39SJed Brown     ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
141814354c39SJed Brown     next = next->next;
141914354c39SJed Brown   }
142014354c39SJed Brown   PetscFunctionReturn(0);
142114354c39SJed Brown }
142247c6ae99SBarry Smith 
1423e727c939SJed Brown PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
142447c6ae99SBarry Smith {
142547c6ae99SBarry Smith   PetscErrorCode         ierr;
14269ae5db72SJed Brown   PetscInt               m,n,M,N,nDM,i;
142747c6ae99SBarry Smith   struct DMCompositeLink *nextc;
142847c6ae99SBarry Smith   struct DMCompositeLink *nextf;
142925296bd5SBarry Smith   Vec                    gcoarse,gfine,*vecs;
143047c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
143147c6ae99SBarry Smith   DM_Composite           *comfine   = (DM_Composite*)fine->data;
14329ae5db72SJed Brown   Mat                    *mats;
143347c6ae99SBarry Smith 
143447c6ae99SBarry Smith   PetscFunctionBegin;
143547c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
143647c6ae99SBarry Smith   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
1437f692024eSJed Brown   ierr = DMSetUp(coarse);CHKERRQ(ierr);
1438f692024eSJed Brown   ierr = DMSetUp(fine);CHKERRQ(ierr);
143947c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14409ae5db72SJed Brown   ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
14419ae5db72SJed Brown   ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr);
144247c6ae99SBarry Smith   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
144347c6ae99SBarry Smith   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
144447c6ae99SBarry Smith   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
144547c6ae99SBarry Smith   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
14469ae5db72SJed Brown   ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
14479ae5db72SJed Brown   ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr);
144847c6ae99SBarry Smith 
14499ae5db72SJed Brown   nDM = comfine->nDM;
1450ce94432eSBarry 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);
14511795a4d1SJed Brown   ierr = PetscCalloc1(nDM*nDM,&mats);CHKERRQ(ierr);
145225296bd5SBarry Smith   if (v) {
14531795a4d1SJed Brown     ierr = PetscCalloc1(nDM,&vecs);CHKERRQ(ierr);
145425296bd5SBarry Smith   }
145547c6ae99SBarry Smith 
145647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
14579ae5db72SJed Brown   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
145825296bd5SBarry Smith     if (!v) {
14590298fd71SBarry Smith       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL);CHKERRQ(ierr);
146025296bd5SBarry Smith     } else {
146125296bd5SBarry Smith       ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr);
146225296bd5SBarry Smith     }
146347c6ae99SBarry Smith   }
1464ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A);CHKERRQ(ierr);
146525296bd5SBarry Smith   if (v) {
1466ce94432eSBarry Smith     ierr = VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v);CHKERRQ(ierr);
146725296bd5SBarry Smith   }
14689ae5db72SJed Brown   for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);}
14699ae5db72SJed Brown   ierr = PetscFree(mats);CHKERRQ(ierr);
147025296bd5SBarry Smith   if (v) {
147125296bd5SBarry Smith     for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);}
147225296bd5SBarry Smith     ierr = PetscFree(vecs);CHKERRQ(ierr);
147325296bd5SBarry Smith   }
147447c6ae99SBarry Smith   PetscFunctionReturn(0);
147547c6ae99SBarry Smith }
147647c6ae99SBarry Smith 
1477184d77edSJed Brown static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
14781411c6eeSJed Brown {
14791411c6eeSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
14801411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1481f7efa3c7SJed Brown   PetscInt               i;
14821411c6eeSJed Brown   PetscErrorCode         ierr;
14831411c6eeSJed Brown 
14841411c6eeSJed Brown   PetscFunctionBegin;
14851411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
14861411c6eeSJed Brown   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
1487ce94432eSBarry Smith   ierr = ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr);
14889ae5db72SJed Brown   for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
14891411c6eeSJed Brown   ierr = PetscFree(ltogs);CHKERRQ(ierr);
14901411c6eeSJed Brown   PetscFunctionReturn(0);
14911411c6eeSJed Brown }
14921411c6eeSJed Brown 
14931411c6eeSJed Brown 
1494b412c318SBarry Smith PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring)
149547c6ae99SBarry Smith {
149647c6ae99SBarry Smith   PetscErrorCode  ierr;
149747c6ae99SBarry Smith   PetscInt        n,i,cnt;
149847c6ae99SBarry Smith   ISColoringValue *colors;
149947c6ae99SBarry Smith   PetscBool       dense  = PETSC_FALSE;
150047c6ae99SBarry Smith   ISColoringValue maxcol = 0;
150147c6ae99SBarry Smith   DM_Composite    *com   = (DM_Composite*)dm->data;
150247c6ae99SBarry Smith 
150347c6ae99SBarry Smith   PetscFunctionBegin;
150447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
15055bdb020cSBarry Smith   if (ctype == IS_COLORING_LOCAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported");
1506e3247f34SBarry Smith   else if (ctype == IS_COLORING_GLOBAL) {
150747c6ae99SBarry Smith     n = com->n;
1508ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
1509785e854fSJed Brown   ierr = PetscMalloc1(n,&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
151047c6ae99SBarry Smith 
1511c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL);CHKERRQ(ierr);
151247c6ae99SBarry Smith   if (dense) {
151347c6ae99SBarry Smith     for (i=0; i<n; i++) {
151447c6ae99SBarry Smith       colors[i] = (ISColoringValue)(com->rstart + i);
151547c6ae99SBarry Smith     }
151647c6ae99SBarry Smith     maxcol = com->N;
151747c6ae99SBarry Smith   } else {
151847c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
151947c6ae99SBarry Smith     PetscMPIInt            rank;
152047c6ae99SBarry Smith 
1521ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
152247c6ae99SBarry Smith     cnt  = 0;
152347c6ae99SBarry Smith     while (next) {
152447c6ae99SBarry Smith       ISColoring lcoloring;
152547c6ae99SBarry Smith 
1526b412c318SBarry Smith       ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring);CHKERRQ(ierr);
152747c6ae99SBarry Smith       for (i=0; i<lcoloring->N; i++) {
152847c6ae99SBarry Smith         colors[cnt++] = maxcol + lcoloring->colors[i];
152947c6ae99SBarry Smith       }
153047c6ae99SBarry Smith       maxcol += lcoloring->n;
1531fcfd50ebSBarry Smith       ierr    = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
153247c6ae99SBarry Smith       next    = next->next;
153347c6ae99SBarry Smith     }
153447c6ae99SBarry Smith   }
1535aaf3ff59SMatthew G. Knepley   ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,PETSC_OWN_POINTER,coloring);CHKERRQ(ierr);
153647c6ae99SBarry Smith   PetscFunctionReturn(0);
153747c6ae99SBarry Smith }
153847c6ae99SBarry Smith 
15397087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
154047c6ae99SBarry Smith {
154147c6ae99SBarry Smith   PetscErrorCode         ierr;
154247c6ae99SBarry Smith   struct DMCompositeLink *next;
154347c6ae99SBarry Smith   PetscScalar            *garray,*larray;
154447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
154547c6ae99SBarry Smith 
154647c6ae99SBarry Smith   PetscFunctionBegin;
154747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
154847c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
154939d35262SVincent Le Chenadec 
155047c6ae99SBarry Smith   if (!com->setup) {
1551d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
155247c6ae99SBarry Smith   }
155339d35262SVincent Le Chenadec 
155447c6ae99SBarry Smith   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
155547c6ae99SBarry Smith   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
155647c6ae99SBarry Smith 
155747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
155839d35262SVincent Le Chenadec   next = com->next;
155947c6ae99SBarry Smith   while (next) {
156047c6ae99SBarry Smith     Vec      local,global;
156147c6ae99SBarry Smith     PetscInt N;
156247c6ae99SBarry Smith 
156347c6ae99SBarry Smith     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
156447c6ae99SBarry Smith     ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
156547c6ae99SBarry Smith     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
156647c6ae99SBarry Smith     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
156747c6ae99SBarry Smith     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
156847c6ae99SBarry Smith     ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
156947c6ae99SBarry Smith     ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
157047c6ae99SBarry Smith     ierr = VecResetArray(global);CHKERRQ(ierr);
157147c6ae99SBarry Smith     ierr = VecResetArray(local);CHKERRQ(ierr);
157247c6ae99SBarry Smith     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
157339d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr);
157439d35262SVincent Le Chenadec 
157506ebdd98SJed Brown     larray += next->nlocal;
157639d35262SVincent Le Chenadec     garray += next->n;
157747c6ae99SBarry Smith     next    = next->next;
157847c6ae99SBarry Smith   }
157947c6ae99SBarry Smith 
15800298fd71SBarry Smith   ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr);
15810298fd71SBarry Smith   ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr);
158247c6ae99SBarry Smith   PetscFunctionReturn(0);
158347c6ae99SBarry Smith }
158447c6ae99SBarry Smith 
15857087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
15860c010503SBarry Smith {
15870c010503SBarry Smith   PetscFunctionBegin;
158839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
158939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
159039d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,4);
159139d35262SVincent Le Chenadec   PetscFunctionReturn(0);
159239d35262SVincent Le Chenadec }
159339d35262SVincent Le Chenadec 
159439d35262SVincent Le Chenadec PetscErrorCode  DMLocalToGlobalBegin_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
159539d35262SVincent Le Chenadec {
159639d35262SVincent Le Chenadec   PetscErrorCode         ierr;
159739d35262SVincent Le Chenadec   struct DMCompositeLink *next;
159839d35262SVincent Le Chenadec   PetscScalar            *larray,*garray;
159939d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite*)dm->data;
160039d35262SVincent Le Chenadec 
160139d35262SVincent Le Chenadec   PetscFunctionBegin;
160239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
160339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
160439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
160539d35262SVincent Le Chenadec 
160639d35262SVincent Le Chenadec   if (!com->setup) {
160739d35262SVincent Le Chenadec     ierr = DMSetUp(dm);CHKERRQ(ierr);
160839d35262SVincent Le Chenadec   }
160939d35262SVincent Le Chenadec 
161039d35262SVincent Le Chenadec   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
161139d35262SVincent Le Chenadec   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
161239d35262SVincent Le Chenadec 
161339d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
161439d35262SVincent Le Chenadec   next = com->next;
161539d35262SVincent Le Chenadec   while (next) {
161639d35262SVincent Le Chenadec     Vec      global,local;
161739d35262SVincent Le Chenadec 
161839d35262SVincent Le Chenadec     ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
161939d35262SVincent Le Chenadec     ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
162039d35262SVincent Le Chenadec     ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
162139d35262SVincent Le Chenadec     ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
162239d35262SVincent Le Chenadec     ierr = DMLocalToGlobalBegin(next->dm,local,mode,global);CHKERRQ(ierr);
162339d35262SVincent Le Chenadec     ierr = DMLocalToGlobalEnd(next->dm,local,mode,global);CHKERRQ(ierr);
162439d35262SVincent Le Chenadec     ierr = VecResetArray(local);CHKERRQ(ierr);
162539d35262SVincent Le Chenadec     ierr = VecResetArray(global);CHKERRQ(ierr);
162639d35262SVincent Le Chenadec     ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
162739d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr);
162839d35262SVincent Le Chenadec 
162939d35262SVincent Le Chenadec     garray += next->n;
163039d35262SVincent Le Chenadec     larray += next->nlocal;
163139d35262SVincent Le Chenadec     next    = next->next;
163239d35262SVincent Le Chenadec   }
163339d35262SVincent Le Chenadec 
163439d35262SVincent Le Chenadec   ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr);
163539d35262SVincent Le Chenadec   ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr);
163639d35262SVincent Le Chenadec 
163739d35262SVincent Le Chenadec   PetscFunctionReturn(0);
163839d35262SVincent Le Chenadec }
163939d35262SVincent Le Chenadec 
164039d35262SVincent Le Chenadec PetscErrorCode  DMLocalToGlobalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
164139d35262SVincent Le Chenadec {
164239d35262SVincent Le Chenadec   PetscFunctionBegin;
164339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
164439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
164539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
164639d35262SVincent Le Chenadec   PetscFunctionReturn(0);
164739d35262SVincent Le Chenadec }
164839d35262SVincent Le Chenadec 
164939d35262SVincent Le Chenadec PetscErrorCode  DMLocalToLocalBegin_Composite(DM dm,Vec vec1,InsertMode mode,Vec vec2)
165039d35262SVincent Le Chenadec {
165139d35262SVincent Le Chenadec   PetscErrorCode         ierr;
165239d35262SVincent Le Chenadec   struct DMCompositeLink *next;
165339d35262SVincent Le Chenadec   PetscScalar            *array1,*array2;
165439d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite*)dm->data;
165539d35262SVincent Le Chenadec 
165639d35262SVincent Le Chenadec   PetscFunctionBegin;
165739d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
165839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec1,VEC_CLASSID,2);
165939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec2,VEC_CLASSID,4);
166039d35262SVincent Le Chenadec 
166139d35262SVincent Le Chenadec   if (!com->setup) {
166239d35262SVincent Le Chenadec     ierr = DMSetUp(dm);CHKERRQ(ierr);
166339d35262SVincent Le Chenadec   }
166439d35262SVincent Le Chenadec 
166539d35262SVincent Le Chenadec   ierr = VecGetArray(vec1,&array1);CHKERRQ(ierr);
166639d35262SVincent Le Chenadec   ierr = VecGetArray(vec2,&array2);CHKERRQ(ierr);
166739d35262SVincent Le Chenadec 
166839d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
166939d35262SVincent Le Chenadec   next = com->next;
167039d35262SVincent Le Chenadec   while (next) {
167139d35262SVincent Le Chenadec     Vec      local1,local2;
167239d35262SVincent Le Chenadec 
167339d35262SVincent Le Chenadec     ierr = DMGetLocalVector(next->dm,&local1);CHKERRQ(ierr);
167439d35262SVincent Le Chenadec     ierr = VecPlaceArray(local1,array1);CHKERRQ(ierr);
167539d35262SVincent Le Chenadec     ierr = DMGetLocalVector(next->dm,&local2);CHKERRQ(ierr);
167639d35262SVincent Le Chenadec     ierr = VecPlaceArray(local2,array2);CHKERRQ(ierr);
167739d35262SVincent Le Chenadec     ierr = DMLocalToLocalBegin(next->dm,local1,mode,local2);CHKERRQ(ierr);
167839d35262SVincent Le Chenadec     ierr = DMLocalToLocalEnd(next->dm,local1,mode,local2);CHKERRQ(ierr);
167939d35262SVincent Le Chenadec     ierr = VecResetArray(local2);CHKERRQ(ierr);
168039d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local2);CHKERRQ(ierr);
168139d35262SVincent Le Chenadec     ierr = VecResetArray(local1);CHKERRQ(ierr);
168239d35262SVincent Le Chenadec     ierr = DMRestoreLocalVector(next->dm,&local1);CHKERRQ(ierr);
168339d35262SVincent Le Chenadec 
168439d35262SVincent Le Chenadec     array1 += next->nlocal;
168539d35262SVincent Le Chenadec     array2 += next->nlocal;
168639d35262SVincent Le Chenadec     next    = next->next;
168739d35262SVincent Le Chenadec   }
168839d35262SVincent Le Chenadec 
168939d35262SVincent Le Chenadec   ierr = VecRestoreArray(vec1,NULL);CHKERRQ(ierr);
169039d35262SVincent Le Chenadec   ierr = VecRestoreArray(vec2,NULL);CHKERRQ(ierr);
169139d35262SVincent Le Chenadec 
169239d35262SVincent Le Chenadec   PetscFunctionReturn(0);
169339d35262SVincent Le Chenadec }
169439d35262SVincent Le Chenadec 
169539d35262SVincent Le Chenadec PetscErrorCode  DMLocalToLocalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
169639d35262SVincent Le Chenadec {
169739d35262SVincent Le Chenadec   PetscFunctionBegin;
169839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
169939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
170039d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
17010c010503SBarry Smith   PetscFunctionReturn(0);
17020c010503SBarry Smith }
170347c6ae99SBarry Smith 
17046ae3a549SBarry Smith /*MC
17056ae3a549SBarry Smith    DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs
17066ae3a549SBarry Smith 
17076ae3a549SBarry Smith   Level: intermediate
17086ae3a549SBarry Smith 
17091abcec8cSBarry Smith .seealso: DMType, DM, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate()
17106ae3a549SBarry Smith M*/
17116ae3a549SBarry Smith 
17126ae3a549SBarry Smith 
17138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1714a4121054SBarry Smith {
1715a4121054SBarry Smith   PetscErrorCode ierr;
1716a4121054SBarry Smith   DM_Composite   *com;
1717a4121054SBarry Smith 
1718a4121054SBarry Smith   PetscFunctionBegin;
1719b00a9115SJed Brown   ierr          = PetscNewLog(p,&com);CHKERRQ(ierr);
1720a4121054SBarry Smith   p->data       = com;
1721a4121054SBarry Smith   ierr          = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1722a4121054SBarry Smith   com->n        = 0;
17237ac2b803SAlex Fikl   com->nghost   = 0;
17240298fd71SBarry Smith   com->next     = NULL;
1725a4121054SBarry Smith   com->nDM      = 0;
1726a4121054SBarry Smith 
1727a4121054SBarry Smith   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1728a4121054SBarry Smith   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1729184d77edSJed Brown   p->ops->getlocaltoglobalmapping         = DMGetLocalToGlobalMapping_Composite;
17304d343eeaSMatthew G Knepley   p->ops->createfieldis                   = DMCreateFieldIS_Composite;
173116621825SDmitry Karpeev   p->ops->createfielddecomposition        = DMCreateFieldDecomposition_Composite;
1732a4121054SBarry Smith   p->ops->refine                          = DMRefine_Composite;
173314354c39SJed Brown   p->ops->coarsen                         = DMCoarsen_Composite;
173425296bd5SBarry Smith   p->ops->createinterpolation             = DMCreateInterpolation_Composite;
173525296bd5SBarry Smith   p->ops->creatematrix                    = DMCreateMatrix_Composite;
1736e727c939SJed Brown   p->ops->getcoloring                     = DMCreateColoring_Composite;
1737a4121054SBarry Smith   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1738a4121054SBarry Smith   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
173939d35262SVincent Le Chenadec   p->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Composite;
174039d35262SVincent Le Chenadec   p->ops->localtoglobalend                = DMLocalToGlobalEnd_Composite;
174139d35262SVincent Le Chenadec   p->ops->localtolocalbegin               = DMLocalToLocalBegin_Composite;
174239d35262SVincent Le Chenadec   p->ops->localtolocalend                 = DMLocalToLocalEnd_Composite;
1743a4121054SBarry Smith   p->ops->destroy                         = DMDestroy_Composite;
1744a4121054SBarry Smith   p->ops->view                            = DMView_Composite;
1745a4121054SBarry Smith   p->ops->setup                           = DMSetUp_Composite;
1746e10fd815SStefano Zampini 
1747e10fd815SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)p,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Composite);CHKERRQ(ierr);
1748a4121054SBarry Smith   PetscFunctionReturn(0);
1749a4121054SBarry Smith }
1750a4121054SBarry Smith 
17511fd49c25SBarry Smith /*@
17520c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
17530c010503SBarry Smith       vectors made up of several subvectors.
17540c010503SBarry Smith 
17550c010503SBarry Smith     Collective on MPI_Comm
175647c6ae99SBarry Smith 
175747c6ae99SBarry Smith     Input Parameter:
17580c010503SBarry Smith .   comm - the processors that will share the global vector
17590c010503SBarry Smith 
17600c010503SBarry Smith     Output Parameters:
17610c010503SBarry Smith .   packer - the packer object
176247c6ae99SBarry Smith 
176347c6ae99SBarry Smith     Level: advanced
176447c6ae99SBarry Smith 
17651abcec8cSBarry Smith .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(), DMCOMPOSITE,DMCreate()
17666eb61c8cSJed Brown          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
176747c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
176847c6ae99SBarry Smith 
176947c6ae99SBarry Smith @*/
17707087cfbeSBarry Smith PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
177147c6ae99SBarry Smith {
17720c010503SBarry Smith   PetscErrorCode ierr;
17730c010503SBarry Smith 
177447c6ae99SBarry Smith   PetscFunctionBegin;
17750c010503SBarry Smith   PetscValidPointer(packer,2);
1776a4121054SBarry Smith   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1777a4121054SBarry Smith   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
177847c6ae99SBarry Smith   PetscFunctionReturn(0);
177947c6ae99SBarry Smith }
1780