xref: /petsc/src/dm/impls/composite/pack.c (revision c2e3fba1fe1cda7e6350bbca19c4ed35ce95940a)
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 
11d083f849SBarry Smith     Logically Collective
1247c6ae99SBarry Smith 
13d8d19677SJose E. Roman     Input Parameters:
1447c6ae99SBarry Smith +   dm - the composite object
1547c6ae99SBarry Smith -   formcouplelocations - routine to set the nonzero locations in the matrix
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith     Level: advanced
1847c6ae99SBarry Smith 
19f5f57ec0SBarry Smith     Not available from Fortran
20f5f57ec0SBarry Smith 
2195452b02SPatrick Sanan     Notes:
2295452b02SPatrick Sanan     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;
2971b14b3eSStefano Zampini   PetscBool      flg;
3047c6ae99SBarry Smith 
3147c6ae99SBarry Smith   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
337a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
3447c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
3547c6ae99SBarry Smith   PetscFunctionReturn(0);
3647c6ae99SBarry Smith }
3747c6ae99SBarry Smith 
386bf464f9SBarry Smith PetscErrorCode  DMDestroy_Composite(DM dm)
3947c6ae99SBarry Smith {
4047c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
4147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
4247c6ae99SBarry Smith 
4347c6ae99SBarry Smith   PetscFunctionBegin;
4447c6ae99SBarry Smith   next = com->next;
4547c6ae99SBarry Smith   while (next) {
4647c6ae99SBarry Smith     prev = next;
4747c6ae99SBarry Smith     next = next->next;
489566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&prev->dm));
499566063dSJacob Faibussowitsch     PetscCall(PetscFree(prev->grstarts));
509566063dSJacob Faibussowitsch     PetscCall(PetscFree(prev));
5147c6ae99SBarry Smith   }
529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL));
53435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
549566063dSJacob Faibussowitsch   PetscCall(PetscFree(com));
5547c6ae99SBarry Smith   PetscFunctionReturn(0);
5647c6ae99SBarry Smith }
5747c6ae99SBarry Smith 
587087cfbeSBarry Smith PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
5947c6ae99SBarry Smith {
6047c6ae99SBarry Smith   PetscBool      iascii;
6147c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
6247c6ae99SBarry Smith 
6347c6ae99SBarry Smith   PetscFunctionBegin;
649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii));
6547c6ae99SBarry Smith   if (iascii) {
6647c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
6747c6ae99SBarry Smith     PetscInt               i;
6847c6ae99SBarry Smith 
699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix"));
7063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v,"  contains %" PetscInt_FMT " DMs\n",com->nDM));
719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(v));
7247c6ae99SBarry Smith     for (i=0; lnk; lnk=lnk->next,i++) {
7363a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(v,"Link %" PetscInt_FMT ": DM of type %s\n",i,((PetscObject)lnk->dm)->type_name));
749566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(v));
759566063dSJacob Faibussowitsch       PetscCall(DMView(lnk->dm,v));
769566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(v));
7747c6ae99SBarry Smith     }
789566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(v));
7947c6ae99SBarry Smith   }
8047c6ae99SBarry Smith   PetscFunctionReturn(0);
8147c6ae99SBarry Smith }
8247c6ae99SBarry Smith 
8347c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
847087cfbeSBarry Smith PetscErrorCode  DMSetUp_Composite(DM dm)
8547c6ae99SBarry Smith {
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;
9328b400f6SJacob Faibussowitsch   PetscCheck(!com->setup,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
949566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm),&map));
959566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(map,com->n));
969566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(map,PETSC_DETERMINE));
979566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(map,1));
989566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(map));
999566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(map,&com->N));
1009566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(map,&com->rstart,NULL));
1019566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&map));
10247c6ae99SBarry Smith 
1039ae5db72SJed Brown   /* now set the rstart for each linked vector */
1049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
1059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size));
10647c6ae99SBarry Smith   while (next) {
10747c6ae99SBarry Smith     next->rstart  = nprev;
10806ebdd98SJed Brown     nprev        += next->n;
10947c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
1109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size,&next->grstarts));
1119566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,PetscObjectComm((PetscObject)dm)));
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;
13871b14b3eSStefano Zampini   PetscBool      flg;
1395fd66863SKarl Rupp 
14047c6ae99SBarry Smith   PetscFunctionBegin;
14147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1429566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
1437a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
14447c6ae99SBarry Smith   *nDM = com->nDM;
14547c6ae99SBarry Smith   PetscFunctionReturn(0);
14647c6ae99SBarry Smith }
14747c6ae99SBarry Smith 
14847c6ae99SBarry Smith /*@C
14947c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
15047c6ae99SBarry Smith        representation.
15147c6ae99SBarry Smith 
152d083f849SBarry Smith     Collective on dm
15347c6ae99SBarry Smith 
1549ae5db72SJed Brown     Input Parameters:
15547c6ae99SBarry Smith +    dm - the packer object
1569ae5db72SJed Brown -    gvec - the global vector
1579ae5db72SJed Brown 
1589ae5db72SJed Brown     Output Parameters:
1590298fd71SBarry Smith .    Vec* ... - the packed parallel vectors, NULL for those that are not needed
16047c6ae99SBarry Smith 
16195452b02SPatrick Sanan     Notes:
16295452b02SPatrick Sanan     Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
16347c6ae99SBarry Smith 
164f73e5cebSJed Brown     Fortran Notes:
165f73e5cebSJed Brown 
166f73e5cebSJed Brown     Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
167f73e5cebSJed Brown     or use the alternative interface DMCompositeGetAccessArray().
168f73e5cebSJed Brown 
16947c6ae99SBarry Smith     Level: advanced
17047c6ae99SBarry Smith 
171db781477SPatrick Sanan .seealso: `DMCompositeGetEntries()`, `DMCompositeScatter()`
17247c6ae99SBarry Smith @*/
1737087cfbeSBarry Smith PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
17447c6ae99SBarry Smith {
17547c6ae99SBarry Smith   va_list                Argp;
17647c6ae99SBarry Smith   struct DMCompositeLink *next;
17747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
1785edff71fSBarry Smith   PetscInt               readonly;
17971b14b3eSStefano Zampini   PetscBool              flg;
18047c6ae99SBarry Smith 
18147c6ae99SBarry Smith   PetscFunctionBegin;
18247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
18347c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
1849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
1857a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
18647c6ae99SBarry Smith   next = com->next;
18747c6ae99SBarry Smith   if (!com->setup) {
1889566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
18947c6ae99SBarry Smith   }
19047c6ae99SBarry Smith 
1919566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec,&readonly));
19247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
19347c6ae99SBarry Smith   va_start(Argp,gvec);
19447c6ae99SBarry Smith   while (next) {
19547c6ae99SBarry Smith     Vec *vec;
19647c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
1979ae5db72SJed Brown     if (vec) {
1989566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm,vec));
1995edff71fSBarry Smith       if (readonly) {
2005edff71fSBarry Smith         const PetscScalar *array;
2019566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(gvec,&array));
2029566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec,array+next->rstart));
2039566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(*vec));
2049566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(gvec,&array));
2055edff71fSBarry Smith       } else {
2065edff71fSBarry Smith         PetscScalar *array;
2079566063dSJacob Faibussowitsch         PetscCall(VecGetArray(gvec,&array));
2089566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec,array+next->rstart));
2099566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(gvec,&array));
21047c6ae99SBarry Smith       }
2115edff71fSBarry Smith     }
21247c6ae99SBarry Smith     next = next->next;
21347c6ae99SBarry Smith   }
21447c6ae99SBarry Smith   va_end(Argp);
21547c6ae99SBarry Smith   PetscFunctionReturn(0);
21647c6ae99SBarry Smith }
21747c6ae99SBarry Smith 
218f73e5cebSJed Brown /*@C
219f73e5cebSJed Brown     DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
220f73e5cebSJed Brown        representation.
221f73e5cebSJed Brown 
222d083f849SBarry Smith     Collective on dm
223f73e5cebSJed Brown 
224f73e5cebSJed Brown     Input Parameters:
225f73e5cebSJed Brown +    dm - the packer object
226f73e5cebSJed Brown .    pvec - packed vector
227f73e5cebSJed Brown .    nwanted - number of vectors wanted
2280298fd71SBarry Smith -    wanted - sorted array of vectors wanted, or NULL to get all vectors
229f73e5cebSJed Brown 
230f73e5cebSJed Brown     Output Parameters:
231f73e5cebSJed Brown .    vecs - array of requested global vectors (must be allocated)
232f73e5cebSJed Brown 
23395452b02SPatrick Sanan     Notes:
23495452b02SPatrick Sanan     Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them
235f73e5cebSJed Brown 
236f73e5cebSJed Brown     Level: advanced
237f73e5cebSJed Brown 
238db781477SPatrick Sanan .seealso: `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
239f73e5cebSJed Brown @*/
240f73e5cebSJed Brown PetscErrorCode  DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
241f73e5cebSJed Brown {
242f73e5cebSJed Brown   struct DMCompositeLink *link;
243f73e5cebSJed Brown   PetscInt               i,wnum;
244f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
245bee642f7SBarry Smith   PetscInt               readonly;
24671b14b3eSStefano Zampini   PetscBool              flg;
247f73e5cebSJed Brown 
248f73e5cebSJed Brown   PetscFunctionBegin;
249f73e5cebSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
250f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
2519566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
2527a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
253f73e5cebSJed Brown   if (!com->setup) {
2549566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
255f73e5cebSJed Brown   }
256f73e5cebSJed Brown 
2579566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec,&readonly));
258f73e5cebSJed Brown   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
259f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
260f73e5cebSJed Brown       Vec v;
2619566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(link->dm,&v));
262bee642f7SBarry Smith       if (readonly) {
263bee642f7SBarry Smith         const PetscScalar *array;
2649566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec,&array));
2659566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v,array+link->rstart));
2669566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
2679566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec,&array));
268bee642f7SBarry Smith       } else {
269bee642f7SBarry Smith         PetscScalar *array;
2709566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec,&array));
2719566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v,array+link->rstart));
2729566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec,&array));
273bee642f7SBarry Smith       }
274f73e5cebSJed Brown       vecs[wnum++] = v;
275f73e5cebSJed Brown     }
276f73e5cebSJed Brown   }
277f73e5cebSJed Brown   PetscFunctionReturn(0);
278f73e5cebSJed Brown }
279f73e5cebSJed Brown 
2807ac2b803SAlex Fikl /*@C
2817ac2b803SAlex Fikl     DMCompositeGetLocalAccessArray - Allows one to access the individual
2827ac2b803SAlex Fikl     packed vectors in their local representation.
2837ac2b803SAlex Fikl 
284d083f849SBarry Smith     Collective on dm.
2857ac2b803SAlex Fikl 
2867ac2b803SAlex Fikl     Input Parameters:
2877ac2b803SAlex Fikl +    dm - the packer object
2887ac2b803SAlex Fikl .    pvec - packed vector
2897ac2b803SAlex Fikl .    nwanted - number of vectors wanted
2907ac2b803SAlex Fikl -    wanted - sorted array of vectors wanted, or NULL to get all vectors
2917ac2b803SAlex Fikl 
2927ac2b803SAlex Fikl     Output Parameters:
2937ac2b803SAlex Fikl .    vecs - array of requested local vectors (must be allocated)
2947ac2b803SAlex Fikl 
29595452b02SPatrick Sanan     Notes:
29695452b02SPatrick Sanan     Use DMCompositeRestoreLocalAccessArray() to return the vectors
2977ac2b803SAlex Fikl     when you no longer need them.
2987ac2b803SAlex Fikl 
2997ac2b803SAlex Fikl     Level: advanced
3007ac2b803SAlex Fikl 
301db781477SPatrick Sanan .seealso: `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`,
302db781477SPatrick Sanan           `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
3037ac2b803SAlex Fikl @*/
3047ac2b803SAlex Fikl PetscErrorCode  DMCompositeGetLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
3057ac2b803SAlex Fikl {
3067ac2b803SAlex Fikl   struct DMCompositeLink *link;
3077ac2b803SAlex Fikl   PetscInt               i,wnum;
3087ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite*)dm->data;
3097ac2b803SAlex Fikl   PetscInt               readonly;
3107ac2b803SAlex Fikl   PetscInt               nlocal = 0;
31171b14b3eSStefano Zampini   PetscBool              flg;
3127ac2b803SAlex Fikl 
3137ac2b803SAlex Fikl   PetscFunctionBegin;
3147ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3157ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
3169566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
3177a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
3187ac2b803SAlex Fikl   if (!com->setup) {
3199566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
3207ac2b803SAlex Fikl   }
3217ac2b803SAlex Fikl 
3229566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec,&readonly));
3237ac2b803SAlex Fikl   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
3247ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
3257ac2b803SAlex Fikl       Vec v;
3269566063dSJacob Faibussowitsch       PetscCall(DMGetLocalVector(link->dm,&v));
3277ac2b803SAlex Fikl       if (readonly) {
3287ac2b803SAlex Fikl         const PetscScalar *array;
3299566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec,&array));
3309566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v,array+nlocal));
331b1c3483dSMark Adams         // this method does not make sense. The local vectors are not updated with a global-to-local and the user can not do it because it is locked
3329566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
3339566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec,&array));
3347ac2b803SAlex Fikl       } else {
3357ac2b803SAlex Fikl         PetscScalar *array;
3369566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec,&array));
3379566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v,array+nlocal));
3389566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec,&array));
3397ac2b803SAlex Fikl       }
3407ac2b803SAlex Fikl       vecs[wnum++] = v;
3417ac2b803SAlex Fikl     }
3427ac2b803SAlex Fikl 
3437ac2b803SAlex Fikl     nlocal += link->nlocal;
3447ac2b803SAlex Fikl   }
3457ac2b803SAlex Fikl 
3467ac2b803SAlex Fikl   PetscFunctionReturn(0);
3477ac2b803SAlex Fikl }
3487ac2b803SAlex Fikl 
34947c6ae99SBarry Smith /*@C
350aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
35147c6ae99SBarry Smith        representation.
35247c6ae99SBarry Smith 
353d083f849SBarry Smith     Collective on dm
35447c6ae99SBarry Smith 
3559ae5db72SJed Brown     Input Parameters:
35647c6ae99SBarry Smith +    dm - the packer object
35747c6ae99SBarry Smith .    gvec - the global vector
3580298fd71SBarry Smith -    Vec* ... - the individual parallel vectors, NULL for those that are not needed
35947c6ae99SBarry Smith 
36047c6ae99SBarry Smith     Level: advanced
36147c6ae99SBarry Smith 
362db781477SPatrick Sanan .seealso `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
363db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`,
364db781477SPatrick Sanan          `DMCompositeRestoreAccess()`, `DMCompositeGetAccess()`
36547c6ae99SBarry Smith 
36647c6ae99SBarry Smith @*/
3677087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
36847c6ae99SBarry Smith {
36947c6ae99SBarry Smith   va_list                Argp;
37047c6ae99SBarry Smith   struct DMCompositeLink *next;
37147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
3725edff71fSBarry Smith   PetscInt               readonly;
37371b14b3eSStefano Zampini   PetscBool              flg;
37447c6ae99SBarry Smith 
37547c6ae99SBarry Smith   PetscFunctionBegin;
37647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
37747c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
3789566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
3797a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
38047c6ae99SBarry Smith   next = com->next;
38147c6ae99SBarry Smith   if (!com->setup) {
3829566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
38347c6ae99SBarry Smith   }
38447c6ae99SBarry Smith 
3859566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec,&readonly));
38647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
38747c6ae99SBarry Smith   va_start(Argp,gvec);
38847c6ae99SBarry Smith   while (next) {
38947c6ae99SBarry Smith     Vec *vec;
39047c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
3919ae5db72SJed Brown     if (vec) {
3929566063dSJacob Faibussowitsch       PetscCall(VecResetArray(*vec));
3935edff71fSBarry Smith       if (readonly) {
3949566063dSJacob Faibussowitsch         PetscCall(VecLockReadPop(*vec));
3955edff71fSBarry Smith       }
3969566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm,vec));
39747c6ae99SBarry Smith     }
39847c6ae99SBarry Smith     next = next->next;
39947c6ae99SBarry Smith   }
40047c6ae99SBarry Smith   va_end(Argp);
40147c6ae99SBarry Smith   PetscFunctionReturn(0);
40247c6ae99SBarry Smith }
40347c6ae99SBarry Smith 
404f73e5cebSJed Brown /*@C
405f73e5cebSJed Brown     DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray()
406f73e5cebSJed Brown 
407d083f849SBarry Smith     Collective on dm
408f73e5cebSJed Brown 
409f73e5cebSJed Brown     Input Parameters:
410f73e5cebSJed Brown +    dm - the packer object
411f73e5cebSJed Brown .    pvec - packed vector
412f73e5cebSJed Brown .    nwanted - number of vectors wanted
4130298fd71SBarry Smith .    wanted - sorted array of vectors wanted, or NULL to get all vectors
414f73e5cebSJed Brown -    vecs - array of global vectors to return
415f73e5cebSJed Brown 
416f73e5cebSJed Brown     Level: advanced
417f73e5cebSJed Brown 
418db781477SPatrick Sanan .seealso: `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
419f73e5cebSJed Brown @*/
420f73e5cebSJed Brown PetscErrorCode  DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
421f73e5cebSJed Brown {
422f73e5cebSJed Brown   struct DMCompositeLink *link;
423f73e5cebSJed Brown   PetscInt               i,wnum;
424f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
425bee642f7SBarry Smith   PetscInt               readonly;
42671b14b3eSStefano Zampini   PetscBool              flg;
427f73e5cebSJed Brown 
428f73e5cebSJed Brown   PetscFunctionBegin;
429f73e5cebSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
430f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
4319566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
4327a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
433f73e5cebSJed Brown   if (!com->setup) {
4349566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
435f73e5cebSJed Brown   }
436f73e5cebSJed Brown 
4379566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec,&readonly));
438f73e5cebSJed Brown   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
439f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
4409566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
441bee642f7SBarry Smith       if (readonly) {
4429566063dSJacob Faibussowitsch         PetscCall(VecLockReadPop(vecs[wnum]));
443bee642f7SBarry Smith       }
4449566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(link->dm,&vecs[wnum]));
445f73e5cebSJed Brown       wnum++;
446f73e5cebSJed Brown     }
447f73e5cebSJed Brown   }
448f73e5cebSJed Brown   PetscFunctionReturn(0);
449f73e5cebSJed Brown }
450f73e5cebSJed Brown 
4517ac2b803SAlex Fikl /*@C
4527ac2b803SAlex Fikl     DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray().
4537ac2b803SAlex Fikl 
454d083f849SBarry Smith     Collective on dm.
4557ac2b803SAlex Fikl 
4567ac2b803SAlex Fikl     Input Parameters:
4577ac2b803SAlex Fikl +    dm - the packer object
4587ac2b803SAlex Fikl .    pvec - packed vector
4597ac2b803SAlex Fikl .    nwanted - number of vectors wanted
4607ac2b803SAlex Fikl .    wanted - sorted array of vectors wanted, or NULL to restore all vectors
4617ac2b803SAlex Fikl -    vecs - array of local vectors to return
4627ac2b803SAlex Fikl 
4637ac2b803SAlex Fikl     Level: advanced
4647ac2b803SAlex Fikl 
4657ac2b803SAlex Fikl     Notes:
4667ac2b803SAlex Fikl     nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray()
4677ac2b803SAlex Fikl     otherwise the call will fail.
4687ac2b803SAlex Fikl 
469db781477SPatrick Sanan .seealso: `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`,
470db781477SPatrick Sanan           `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`,
471db781477SPatrick Sanan           `DMCompositeScatter()`, `DMCompositeGather()`
4727ac2b803SAlex Fikl @*/
4737ac2b803SAlex Fikl PetscErrorCode  DMCompositeRestoreLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs)
4747ac2b803SAlex Fikl {
4757ac2b803SAlex Fikl   struct DMCompositeLink *link;
4767ac2b803SAlex Fikl   PetscInt               i,wnum;
4777ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite*)dm->data;
4787ac2b803SAlex Fikl   PetscInt               readonly;
47971b14b3eSStefano Zampini   PetscBool              flg;
4807ac2b803SAlex Fikl 
4817ac2b803SAlex Fikl   PetscFunctionBegin;
4827ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4837ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec,VEC_CLASSID,2);
4849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
4857a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
4867ac2b803SAlex Fikl   if (!com->setup) {
4879566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
4887ac2b803SAlex Fikl   }
4897ac2b803SAlex Fikl 
4909566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec,&readonly));
4917ac2b803SAlex Fikl   for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) {
4927ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
4939566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
4947ac2b803SAlex Fikl       if (readonly) {
4959566063dSJacob Faibussowitsch         PetscCall(VecLockReadPop(vecs[wnum]));
4967ac2b803SAlex Fikl       }
4979566063dSJacob Faibussowitsch       PetscCall(DMRestoreLocalVector(link->dm,&vecs[wnum]));
4987ac2b803SAlex Fikl       wnum++;
4997ac2b803SAlex Fikl     }
5007ac2b803SAlex Fikl   }
5017ac2b803SAlex Fikl   PetscFunctionReturn(0);
5027ac2b803SAlex Fikl }
5037ac2b803SAlex Fikl 
50447c6ae99SBarry Smith /*@C
50547c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
50647c6ae99SBarry Smith 
507d083f849SBarry Smith     Collective on dm
50847c6ae99SBarry Smith 
5099ae5db72SJed Brown     Input Parameters:
51047c6ae99SBarry Smith +    dm - the packer object
51147c6ae99SBarry Smith .    gvec - the global vector
5120298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for those that are not needed
51347c6ae99SBarry Smith 
51447c6ae99SBarry Smith     Level: advanced
51547c6ae99SBarry Smith 
5166f3c3dcfSJed Brown     Notes:
5176f3c3dcfSJed Brown     DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is
5186f3c3dcfSJed Brown     accessible from Fortran.
5196f3c3dcfSJed Brown 
520db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
521db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
522db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
523db781477SPatrick Sanan          `DMCompositeScatterArray()`
52447c6ae99SBarry Smith 
52547c6ae99SBarry Smith @*/
5267087cfbeSBarry Smith PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
52747c6ae99SBarry Smith {
52847c6ae99SBarry Smith   va_list                Argp;
52947c6ae99SBarry Smith   struct DMCompositeLink *next;
5308fd8f222SJed Brown   PetscInt               cnt;
53147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
53271b14b3eSStefano Zampini   PetscBool              flg;
53347c6ae99SBarry Smith 
53447c6ae99SBarry Smith   PetscFunctionBegin;
53547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
53647c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
5379566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
5387a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
53947c6ae99SBarry Smith   if (!com->setup) {
5409566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
54147c6ae99SBarry Smith   }
54247c6ae99SBarry Smith 
54347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
54447c6ae99SBarry Smith   va_start(Argp,gvec);
5458fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
5469ae5db72SJed Brown     Vec local;
5479ae5db72SJed Brown     local = va_arg(Argp, Vec);
5489ae5db72SJed Brown     if (local) {
5499ae5db72SJed Brown       Vec               global;
5505edff71fSBarry Smith       const PetscScalar *array;
55163a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local,VEC_CLASSID,(int)cnt));
5529566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm,&global));
5539566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec,&array));
5549566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global,array+next->rstart));
5559566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local));
5569566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local));
5579566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec,&array));
5589566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5599566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm,&global));
56047c6ae99SBarry Smith     }
56147c6ae99SBarry Smith   }
56247c6ae99SBarry Smith   va_end(Argp);
56347c6ae99SBarry Smith   PetscFunctionReturn(0);
56447c6ae99SBarry Smith }
56547c6ae99SBarry Smith 
5666f3c3dcfSJed Brown /*@
5676f3c3dcfSJed Brown     DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5686f3c3dcfSJed Brown 
569d083f849SBarry Smith     Collective on dm
5706f3c3dcfSJed Brown 
5716f3c3dcfSJed Brown     Input Parameters:
5726f3c3dcfSJed Brown +    dm - the packer object
5736f3c3dcfSJed Brown .    gvec - the global vector
574a2b725a8SWilliam Gropp -    lvecs - array of local vectors, NULL for any that are not needed
5756f3c3dcfSJed Brown 
5766f3c3dcfSJed Brown     Level: advanced
5776f3c3dcfSJed Brown 
5786f3c3dcfSJed Brown     Note:
579907376e6SBarry Smith     This is a non-variadic alternative to DMCompositeScatter()
5806f3c3dcfSJed Brown 
581db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`
582db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
583db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
5846f3c3dcfSJed Brown 
5856f3c3dcfSJed Brown @*/
5866f3c3dcfSJed Brown PetscErrorCode  DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs)
5876f3c3dcfSJed Brown {
5886f3c3dcfSJed Brown   struct DMCompositeLink *next;
5896f3c3dcfSJed Brown   PetscInt               i;
5906f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
59171b14b3eSStefano Zampini   PetscBool              flg;
5926f3c3dcfSJed Brown 
5936f3c3dcfSJed Brown   PetscFunctionBegin;
5946f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5956f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
5969566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
5977a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
5986f3c3dcfSJed Brown   if (!com->setup) {
5999566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
6006f3c3dcfSJed Brown   }
6016f3c3dcfSJed Brown 
6026f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
6036f3c3dcfSJed Brown   for (i=0,next=com->next; next; next=next->next,i++) {
6046f3c3dcfSJed Brown     if (lvecs[i]) {
6056f3c3dcfSJed Brown       Vec         global;
606c5d31e75SLisandro Dalcin       const PetscScalar *array;
6076f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3);
6089566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm,&global));
6099566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec,&array));
6109566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global,(PetscScalar*)array+next->rstart));
6119566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]));
6129566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]));
6139566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec,&array));
6149566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6159566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm,&global));
6166f3c3dcfSJed Brown     }
6176f3c3dcfSJed Brown   }
6186f3c3dcfSJed Brown   PetscFunctionReturn(0);
6196f3c3dcfSJed Brown }
6206f3c3dcfSJed Brown 
62147c6ae99SBarry Smith /*@C
62247c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
62347c6ae99SBarry Smith 
624d083f849SBarry Smith     Collective on dm
62547c6ae99SBarry Smith 
626d8d19677SJose E. Roman     Input Parameters:
62747c6ae99SBarry Smith +    dm - the packer object
62847c6ae99SBarry Smith .    gvec - the global vector
629907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6300298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for any that are not needed
63147c6ae99SBarry Smith 
63247c6ae99SBarry Smith     Level: advanced
63347c6ae99SBarry Smith 
634f5f57ec0SBarry Smith     Not available from Fortran, Fortran users can use DMCompositeGatherArray()
635f5f57ec0SBarry Smith 
636db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
637db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
638db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
63947c6ae99SBarry Smith 
64047c6ae99SBarry Smith @*/
6411dac896bSSatish Balay PetscErrorCode  DMCompositeGather(DM dm,InsertMode imode,Vec gvec,...)
64247c6ae99SBarry Smith {
64347c6ae99SBarry Smith   va_list                Argp;
64447c6ae99SBarry Smith   struct DMCompositeLink *next;
64547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
6468fd8f222SJed Brown   PetscInt               cnt;
64771b14b3eSStefano Zampini   PetscBool              flg;
64847c6ae99SBarry Smith 
64947c6ae99SBarry Smith   PetscFunctionBegin;
65047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
651064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec,VEC_CLASSID,3);
6529566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
6537a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
65447c6ae99SBarry Smith   if (!com->setup) {
6559566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
65647c6ae99SBarry Smith   }
65747c6ae99SBarry Smith 
65847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
6591dac896bSSatish Balay   va_start(Argp,gvec);
6608fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
6619ae5db72SJed Brown     Vec local;
6629ae5db72SJed Brown     local = va_arg(Argp, Vec);
6639ae5db72SJed Brown     if (local) {
66447c6ae99SBarry Smith       PetscScalar *array;
6659ae5db72SJed Brown       Vec         global;
66663a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local,VEC_CLASSID,(int)cnt));
6679566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm,&global));
6689566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec,&array));
6699566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global,array+next->rstart));
6709566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm,local,imode,global));
6719566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm,local,imode,global));
6729566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec,&array));
6739566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6749566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm,&global));
67547c6ae99SBarry Smith     }
67647c6ae99SBarry Smith   }
67747c6ae99SBarry Smith   va_end(Argp);
67847c6ae99SBarry Smith   PetscFunctionReturn(0);
67947c6ae99SBarry Smith }
68047c6ae99SBarry Smith 
6816f3c3dcfSJed Brown /*@
6826f3c3dcfSJed Brown     DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6836f3c3dcfSJed Brown 
684d083f849SBarry Smith     Collective on dm
6856f3c3dcfSJed Brown 
686d8d19677SJose E. Roman     Input Parameters:
6876f3c3dcfSJed Brown +    dm - the packer object
6886f3c3dcfSJed Brown .    gvec - the global vector
689907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6906f3c3dcfSJed Brown -    lvecs - the individual sequential vectors, NULL for any that are not needed
6916f3c3dcfSJed Brown 
6926f3c3dcfSJed Brown     Level: advanced
6936f3c3dcfSJed Brown 
6946f3c3dcfSJed Brown     Notes:
6956f3c3dcfSJed Brown     This is a non-variadic alternative to DMCompositeGather().
6966f3c3dcfSJed Brown 
697db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
698db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
699db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`,
7006f3c3dcfSJed Brown @*/
7011dac896bSSatish Balay PetscErrorCode  DMCompositeGatherArray(DM dm,InsertMode imode,Vec gvec,Vec *lvecs)
7026f3c3dcfSJed Brown {
7036f3c3dcfSJed Brown   struct DMCompositeLink *next;
7046f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
7056f3c3dcfSJed Brown   PetscInt               i;
70671b14b3eSStefano Zampini   PetscBool              flg;
7076f3c3dcfSJed Brown 
7086f3c3dcfSJed Brown   PetscFunctionBegin;
7096f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
710064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec,VEC_CLASSID,3);
7119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
7127a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
7136f3c3dcfSJed Brown   if (!com->setup) {
7149566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
7156f3c3dcfSJed Brown   }
7166f3c3dcfSJed Brown 
7176f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
7186f3c3dcfSJed Brown   for (next=com->next,i=0; next; next=next->next,i++) {
7196f3c3dcfSJed Brown     if (lvecs[i]) {
7206f3c3dcfSJed Brown       PetscScalar *array;
7216f3c3dcfSJed Brown       Vec         global;
722064a246eSJacob Faibussowitsch       PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,4);
7239566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm,&global));
7249566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec,&array));
7259566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global,array+next->rstart));
7269566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global));
7279566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global));
7289566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec,&array));
7299566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
7309566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm,&global));
7316f3c3dcfSJed Brown     }
7326f3c3dcfSJed Brown   }
7336f3c3dcfSJed Brown   PetscFunctionReturn(0);
7346f3c3dcfSJed Brown }
7356f3c3dcfSJed Brown 
736f5f57ec0SBarry Smith /*@
737aa219208SBarry Smith     DMCompositeAddDM - adds a DM vector to a DMComposite
73847c6ae99SBarry Smith 
739d083f849SBarry Smith     Collective on dm
74047c6ae99SBarry Smith 
741d8d19677SJose E. Roman     Input Parameters:
7429b52a9b5SPatrick Sanan +    dmc - the DMComposite (packer) object
743f5f57ec0SBarry Smith -    dm - the DM object
74447c6ae99SBarry Smith 
74547c6ae99SBarry Smith     Level: advanced
74647c6ae99SBarry Smith 
747db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeGather()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
748db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
749db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
75047c6ae99SBarry Smith 
75147c6ae99SBarry Smith @*/
7527087cfbeSBarry Smith PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
75347c6ae99SBarry Smith {
75406ebdd98SJed Brown   PetscInt               n,nlocal;
75547c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
75606ebdd98SJed Brown   Vec                    global,local;
75747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmc->data;
75871b14b3eSStefano Zampini   PetscBool              flg;
75947c6ae99SBarry Smith 
76047c6ae99SBarry Smith   PetscFunctionBegin;
76147c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
76247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
7639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dmc,DMCOMPOSITE,&flg));
7647a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
76547c6ae99SBarry Smith   next = com->next;
76628b400f6SJacob Faibussowitsch   PetscCheck(!com->setup,PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
76747c6ae99SBarry Smith 
76847c6ae99SBarry Smith   /* create new link */
7699566063dSJacob Faibussowitsch   PetscCall(PetscNew(&mine));
7709566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
7719566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm,&global));
7729566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global,&n));
7739566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm,&global));
7749566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&local));
7759566063dSJacob Faibussowitsch   PetscCall(VecGetSize(local,&nlocal));
7769566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&local));
7778865f1eaSKarl Rupp 
77847c6ae99SBarry Smith   mine->n      = n;
77906ebdd98SJed Brown   mine->nlocal = nlocal;
78047c6ae99SBarry Smith   mine->dm     = dm;
7810298fd71SBarry Smith   mine->next   = NULL;
78247c6ae99SBarry Smith   com->n      += n;
7837ac2b803SAlex Fikl   com->nghost += nlocal;
78447c6ae99SBarry Smith 
78547c6ae99SBarry Smith   /* add to end of list */
7868865f1eaSKarl Rupp   if (!next) com->next = mine;
7878865f1eaSKarl Rupp   else {
78847c6ae99SBarry Smith     while (next->next) next = next->next;
78947c6ae99SBarry Smith     next->next = mine;
79047c6ae99SBarry Smith   }
79147c6ae99SBarry Smith   com->nDM++;
79247c6ae99SBarry Smith   com->nmine++;
79347c6ae99SBarry Smith   PetscFunctionReturn(0);
79447c6ae99SBarry Smith }
79547c6ae99SBarry Smith 
7969804daf3SBarry Smith #include <petscdraw.h>
79726887b52SJed Brown PETSC_EXTERN PetscErrorCode  VecView_MPI(Vec,PetscViewer);
7987087cfbeSBarry Smith PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
79947c6ae99SBarry Smith {
80047c6ae99SBarry Smith   DM                     dm;
80147c6ae99SBarry Smith   struct DMCompositeLink *next;
80247c6ae99SBarry Smith   PetscBool              isdraw;
803cef07954SSatish Balay   DM_Composite           *com;
80447c6ae99SBarry Smith 
80547c6ae99SBarry Smith   PetscFunctionBegin;
8069566063dSJacob Faibussowitsch   PetscCall(VecGetDM(gvec, &dm));
8077a8be351SBarry Smith   PetscCheck(dm,PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
80847c6ae99SBarry Smith   com  = (DM_Composite*)dm->data;
80947c6ae99SBarry Smith   next = com->next;
81047c6ae99SBarry Smith 
8119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
81247c6ae99SBarry Smith   if (!isdraw) {
81347c6ae99SBarry Smith     /* do I really want to call this? */
8149566063dSJacob Faibussowitsch     PetscCall(VecView_MPI(gvec,viewer));
81547c6ae99SBarry Smith   } else {
81647c6ae99SBarry Smith     PetscInt cnt = 0;
81747c6ae99SBarry Smith 
81847c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
81947c6ae99SBarry Smith     while (next) {
82047c6ae99SBarry Smith       Vec               vec;
821ca4278abSLisandro Dalcin       const PetscScalar *array;
82247c6ae99SBarry Smith       PetscInt          bs;
82347c6ae99SBarry Smith 
8249ae5db72SJed Brown       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
8259566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm,&vec));
8269566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec,&array));
8279566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(vec,(PetscScalar*)array+next->rstart));
8289566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec,&array));
8299566063dSJacob Faibussowitsch       PetscCall(VecView(vec,viewer));
8309566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vec));
8319566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(vec,&bs));
8329566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm,&vec));
8339566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawBaseAdd(viewer,bs));
83447c6ae99SBarry Smith       cnt += bs;
83547c6ae99SBarry Smith       next = next->next;
83647c6ae99SBarry Smith     }
8379566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawBaseAdd(viewer,-cnt));
83847c6ae99SBarry Smith   }
83947c6ae99SBarry Smith   PetscFunctionReturn(0);
84047c6ae99SBarry Smith }
84147c6ae99SBarry Smith 
8427087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
84347c6ae99SBarry Smith {
84447c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
84547c6ae99SBarry Smith 
84647c6ae99SBarry Smith   PetscFunctionBegin;
84747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
8489566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
8499566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8509566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm),gvec));
8519566063dSJacob Faibussowitsch   PetscCall(VecSetType(*gvec,dm->vectype));
8529566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*gvec,com->n,com->N));
8539566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*gvec, dm));
8549566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite));
85547c6ae99SBarry Smith   PetscFunctionReturn(0);
85647c6ae99SBarry Smith }
85747c6ae99SBarry Smith 
8587087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
85947c6ae99SBarry Smith {
86047c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite*)dm->data;
86147c6ae99SBarry Smith 
86247c6ae99SBarry Smith   PetscFunctionBegin;
86347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
86447c6ae99SBarry Smith   if (!com->setup) {
8659566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(dm));
8669566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
86747c6ae99SBarry Smith   }
8689566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF,lvec));
8699566063dSJacob Faibussowitsch   PetscCall(VecSetType(*lvec,dm->vectype));
8709566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*lvec,com->nghost,PETSC_DECIDE));
8719566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*lvec, dm));
87247c6ae99SBarry Smith   PetscFunctionReturn(0);
87347c6ae99SBarry Smith }
87447c6ae99SBarry Smith 
87547c6ae99SBarry Smith /*@C
8769ae5db72SJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
87747c6ae99SBarry Smith 
87806ebdd98SJed Brown     Collective on DM
87947c6ae99SBarry Smith 
88047c6ae99SBarry Smith     Input Parameter:
88147c6ae99SBarry Smith .    dm - the packer object
88247c6ae99SBarry Smith 
88347c6ae99SBarry Smith     Output Parameters:
8849ae5db72SJed Brown .    ltogs - the individual mappings for each packed vector. Note that this includes
8859ae5db72SJed Brown            all the ghost points that individual ghosted DMDA's may have.
88647c6ae99SBarry Smith 
88747c6ae99SBarry Smith     Level: advanced
88847c6ae99SBarry Smith 
88947c6ae99SBarry Smith     Notes:
8906eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
89147c6ae99SBarry Smith 
892f5f57ec0SBarry Smith     Not available from Fortran
893f5f57ec0SBarry Smith 
894db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
895db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
896*c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
89747c6ae99SBarry Smith 
89847c6ae99SBarry Smith @*/
8997087cfbeSBarry Smith PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
90047c6ae99SBarry Smith {
90147c6ae99SBarry Smith   PetscInt               i,*idx,n,cnt;
90247c6ae99SBarry Smith   struct DMCompositeLink *next;
90347c6ae99SBarry Smith   PetscMPIInt            rank;
90447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
90571b14b3eSStefano Zampini   PetscBool              flg;
90647c6ae99SBarry Smith 
90747c6ae99SBarry Smith   PetscFunctionBegin;
90847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
9099566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
9107a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
9119566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
9129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM,ltogs));
91347c6ae99SBarry Smith   next = com->next;
9149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
91547c6ae99SBarry Smith 
91647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
91747c6ae99SBarry Smith   cnt = 0;
91847c6ae99SBarry Smith   while (next) {
9196eb61c8cSJed Brown     ISLocalToGlobalMapping ltog;
9206eb61c8cSJed Brown     PetscMPIInt            size;
92186994e45SJed Brown     const PetscInt         *suboff,*indices;
9226eb61c8cSJed Brown     Vec                    global;
92347c6ae99SBarry Smith 
9246eb61c8cSJed Brown     /* Get sub-DM global indices for each local dof */
9259566063dSJacob Faibussowitsch     PetscCall(DMGetLocalToGlobalMapping(next->dm,&ltog));
9269566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltog,&n));
9279566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltog,&indices));
9289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&idx));
92947c6ae99SBarry Smith 
9306eb61c8cSJed Brown     /* Get the offsets for the sub-DM global vector */
9319566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm,&global));
9329566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRanges(global,&suboff));
9339566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global),&size));
9346eb61c8cSJed Brown 
9356eb61c8cSJed Brown     /* Shift the sub-DM definition of the global space to the composite global space */
9366eb61c8cSJed Brown     for (i=0; i<n; i++) {
93786994e45SJed Brown       PetscInt subi = indices[i],lo = 0,hi = size,t;
938e333b035SStefano Zampini       /* There's no consensus on what a negative index means,
939e333b035SStefano Zampini          except for skipping when setting the values in vectors and matrices */
940e333b035SStefano Zampini       if (subi < 0) { idx[i] = subi - next->grstarts[rank]; continue; }
9416eb61c8cSJed Brown       /* Binary search to find which rank owns subi */
9426eb61c8cSJed Brown       while (hi-lo > 1) {
9436eb61c8cSJed Brown         t = lo + (hi-lo)/2;
9446eb61c8cSJed Brown         if (suboff[t] > subi) hi = t;
9456eb61c8cSJed Brown         else                  lo = t;
9466eb61c8cSJed Brown       }
9476eb61c8cSJed Brown       idx[i] = subi - suboff[lo] + next->grstarts[lo];
9486eb61c8cSJed Brown     }
9499566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog,&indices));
9509566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]));
9519566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm,&global));
95247c6ae99SBarry Smith     next = next->next;
95347c6ae99SBarry Smith     cnt++;
95447c6ae99SBarry Smith   }
95547c6ae99SBarry Smith   PetscFunctionReturn(0);
95647c6ae99SBarry Smith }
95747c6ae99SBarry Smith 
95887c85e80SJed Brown /*@C
9599ae5db72SJed Brown    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
96087c85e80SJed Brown 
96187c85e80SJed Brown    Not Collective
96287c85e80SJed Brown 
9634165533cSJose E. Roman    Input Parameter:
96487c85e80SJed Brown . dm - composite DM
96587c85e80SJed Brown 
9664165533cSJose E. Roman    Output Parameter:
96787c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
96887c85e80SJed Brown 
96987c85e80SJed Brown    Level: intermediate
97087c85e80SJed Brown 
97187c85e80SJed Brown    Notes:
97287c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
97387c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
9749ae5db72SJed Brown    scatter to a composite local vector.  The user should not typically need to know which is being done.
97587c85e80SJed Brown 
97687c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
97787c85e80SJed Brown 
97887c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
97987c85e80SJed Brown 
98087c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
98187c85e80SJed Brown 
982f5f57ec0SBarry Smith    Not available from Fortran
983f5f57ec0SBarry Smith 
984db781477SPatrick Sanan .seealso: `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`, `MatCreateLocalRef()`
98587c85e80SJed Brown @*/
9867087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
98787c85e80SJed Brown {
98887c85e80SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
98987c85e80SJed Brown   struct DMCompositeLink *link;
99087c85e80SJed Brown   PetscInt               cnt,start;
99171b14b3eSStefano Zampini   PetscBool              flg;
99287c85e80SJed Brown 
99387c85e80SJed Brown   PetscFunctionBegin;
99487c85e80SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
99587c85e80SJed Brown   PetscValidPointer(is,2);
9969566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
9977a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
9989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nmine,is));
99906ebdd98SJed Brown   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
1000520db06cSJed Brown     PetscInt bs;
10019566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]));
10029566063dSJacob Faibussowitsch     PetscCall(DMGetBlockSize(link->dm,&bs));
10039566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize((*is)[cnt],bs));
1004520db06cSJed Brown   }
100587c85e80SJed Brown   PetscFunctionReturn(0);
100687c85e80SJed Brown }
100787c85e80SJed Brown 
100847c6ae99SBarry Smith /*@C
100947c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
101047c6ae99SBarry Smith 
1011d083f849SBarry Smith     Collective on dm
101247c6ae99SBarry Smith 
101347c6ae99SBarry Smith     Input Parameter:
101447c6ae99SBarry Smith .    dm - the packer object
101547c6ae99SBarry Smith 
101647c6ae99SBarry Smith     Output Parameters:
101747c6ae99SBarry Smith .    is - the array of index sets
101847c6ae99SBarry Smith 
101947c6ae99SBarry Smith     Level: advanced
102047c6ae99SBarry Smith 
102147c6ae99SBarry Smith     Notes:
102247c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
102347c6ae99SBarry Smith 
102447c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
102547c6ae99SBarry Smith 
10266eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
10276eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
10286eb61c8cSJed Brown        indices.
102947c6ae99SBarry Smith 
1030f3cb0f7eSJed Brown     Fortran Notes:
1031f3cb0f7eSJed Brown 
1032f3cb0f7eSJed Brown        The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM().
1033f3cb0f7eSJed Brown 
1034db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1035db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
1036*c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
103747c6ae99SBarry Smith 
103847c6ae99SBarry Smith @*/
10397087cfbeSBarry Smith PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
104047c6ae99SBarry Smith {
104166bb578eSMark Adams   PetscInt               cnt = 0;
104247c6ae99SBarry Smith   struct DMCompositeLink *next;
104347c6ae99SBarry Smith   PetscMPIInt            rank;
104447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
104571b14b3eSStefano Zampini   PetscBool              flg;
104647c6ae99SBarry Smith 
104747c6ae99SBarry Smith   PetscFunctionBegin;
104847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
10499566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
10507a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
10517a8be351SBarry Smith   PetscCheck(dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() before");
10529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM,is));
105347c6ae99SBarry Smith   next = com->next;
10549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
105547c6ae99SBarry Smith 
105647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
105747c6ae99SBarry Smith   while (next) {
1058e5e52638SMatthew G. Knepley     PetscDS prob;
1059e5e52638SMatthew G. Knepley 
10609566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt]));
10619566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &prob));
1062e5e52638SMatthew G. Knepley     if (prob) {
106365c226d8SMatthew G. Knepley       MatNullSpace space;
106465c226d8SMatthew G. Knepley       Mat          pmat;
10650f21e855SMatthew G. Knepley       PetscObject  disc;
10660f21e855SMatthew G. Knepley       PetscInt     Nf;
106765c226d8SMatthew G. Knepley 
10689566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(prob, &Nf));
1069f24dd8d2SMatthew G. Knepley       if (cnt < Nf) {
10709566063dSJacob Faibussowitsch         PetscCall(PetscDSGetDiscretization(prob, cnt, &disc));
10719566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject*) &space));
10729566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space));
10739566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space));
10749566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space));
10759566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat));
10769566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat));
107765c226d8SMatthew G. Knepley       }
1078f24dd8d2SMatthew G. Knepley     }
107947c6ae99SBarry Smith     cnt++;
108047c6ae99SBarry Smith     next = next->next;
108147c6ae99SBarry Smith   }
108247c6ae99SBarry Smith   PetscFunctionReturn(0);
108347c6ae99SBarry Smith }
108447c6ae99SBarry Smith 
108521c9b008SJed Brown PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields)
10864d343eeaSMatthew G Knepley {
10874d343eeaSMatthew G Knepley   PetscInt       nDM;
10884d343eeaSMatthew G Knepley   DM             *dms;
10894d343eeaSMatthew G Knepley   PetscInt       i;
10904d343eeaSMatthew G Knepley 
10914d343eeaSMatthew G Knepley   PetscFunctionBegin;
10929566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm, &nDM));
10938865f1eaSKarl Rupp   if (numFields) *numFields = nDM;
10949566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetGlobalISs(dm, fields));
10954d343eeaSMatthew G Knepley   if (fieldNames) {
10969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, &dms));
10979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, fieldNames));
10989566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, dms));
10994d343eeaSMatthew G Knepley     for (i=0; i<nDM; i++) {
11004d343eeaSMatthew G Knepley       char       buf[256];
11014d343eeaSMatthew G Knepley       const char *splitname;
11024d343eeaSMatthew G Knepley 
11034d343eeaSMatthew G Knepley       /* Split naming precedence: object name, prefix, number */
11044d343eeaSMatthew G Knepley       splitname = ((PetscObject) dm)->name;
11054d343eeaSMatthew G Knepley       if (!splitname) {
11069566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname));
11074d343eeaSMatthew G Knepley         if (splitname) {
11084d343eeaSMatthew G Knepley           size_t len;
11099566063dSJacob Faibussowitsch           PetscCall(PetscStrncpy(buf,splitname,sizeof(buf)));
11108caf3d72SBarry Smith           buf[sizeof(buf) - 1] = 0;
11119566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(buf,&len));
11124d343eeaSMatthew G Knepley           if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */
11134d343eeaSMatthew G Knepley           splitname = buf;
11144d343eeaSMatthew G Knepley         }
11154d343eeaSMatthew G Knepley       }
11164d343eeaSMatthew G Knepley       if (!splitname) {
111763a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(buf,sizeof(buf),"%" PetscInt_FMT,i));
11184d343eeaSMatthew G Knepley         splitname = buf;
11194d343eeaSMatthew G Knepley       }
11209566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(splitname,&(*fieldNames)[i]));
11214d343eeaSMatthew G Knepley     }
11229566063dSJacob Faibussowitsch     PetscCall(PetscFree(dms));
11234d343eeaSMatthew G Knepley   }
11244d343eeaSMatthew G Knepley   PetscFunctionReturn(0);
11254d343eeaSMatthew G Knepley }
11264d343eeaSMatthew G Knepley 
1127e7c4fc90SDmitry Karpeev /*
1128e7c4fc90SDmitry Karpeev  This could take over from DMCreateFieldIS(), as it is more general,
11290298fd71SBarry Smith  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1130e7c4fc90SDmitry Karpeev  At this point it's probably best to be less intrusive, however.
1131e7c4fc90SDmitry Karpeev  */
113216621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
1133e7c4fc90SDmitry Karpeev {
1134e7c4fc90SDmitry Karpeev   PetscInt       nDM;
1135e7c4fc90SDmitry Karpeev   PetscInt       i;
1136e7c4fc90SDmitry Karpeev 
1137e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
11389566063dSJacob Faibussowitsch   PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist));
1139e7c4fc90SDmitry Karpeev   if (dmlist) {
11409566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetNumberDM(dm, &nDM));
11419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, dmlist));
11429566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, *dmlist));
1143e7c4fc90SDmitry Karpeev     for (i=0; i<nDM; i++) {
11449566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i])));
1145e7c4fc90SDmitry Karpeev     }
1146e7c4fc90SDmitry Karpeev   }
1147e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
1148e7c4fc90SDmitry Karpeev }
1149e7c4fc90SDmitry Karpeev 
115047c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
115147c6ae99SBarry Smith /*@C
11529ae5db72SJed Brown     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
115347c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
115447c6ae99SBarry Smith 
115547c6ae99SBarry Smith     Not Collective
115647c6ae99SBarry Smith 
115747c6ae99SBarry Smith     Input Parameter:
115847c6ae99SBarry Smith .    dm - the packer object
115947c6ae99SBarry Smith 
116047c6ae99SBarry Smith     Output Parameter:
11619ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
116247c6ae99SBarry Smith 
116347c6ae99SBarry Smith     Level: advanced
116447c6ae99SBarry Smith 
1165f5f57ec0SBarry Smith     Not available from Fortran
1166f5f57ec0SBarry Smith 
1167db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1168db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1169db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
117047c6ae99SBarry Smith 
117147c6ae99SBarry Smith @*/
11727087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
117347c6ae99SBarry Smith {
117447c6ae99SBarry Smith   va_list                Argp;
117547c6ae99SBarry Smith   struct DMCompositeLink *next;
117647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
117771b14b3eSStefano Zampini   PetscBool              flg;
117847c6ae99SBarry Smith 
117947c6ae99SBarry Smith   PetscFunctionBegin;
118047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
11819566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
11827a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
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*);
11899566063dSJacob Faibussowitsch     if (vec) PetscCall(DMGetLocalVector(next->dm,vec));
119047c6ae99SBarry Smith     next = next->next;
119147c6ae99SBarry Smith   }
119247c6ae99SBarry Smith   va_end(Argp);
119347c6ae99SBarry Smith   PetscFunctionReturn(0);
119447c6ae99SBarry Smith }
119547c6ae99SBarry Smith 
119647c6ae99SBarry Smith /*@C
11979ae5db72SJed Brown     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
119847c6ae99SBarry Smith 
119947c6ae99SBarry Smith     Not Collective
120047c6ae99SBarry Smith 
120147c6ae99SBarry Smith     Input Parameter:
120247c6ae99SBarry Smith .    dm - the packer object
120347c6ae99SBarry Smith 
120447c6ae99SBarry Smith     Output Parameter:
12059ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
120647c6ae99SBarry Smith 
120747c6ae99SBarry Smith     Level: advanced
120847c6ae99SBarry Smith 
1209f5f57ec0SBarry Smith     Not available from Fortran
1210f5f57ec0SBarry Smith 
1211db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1212db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1213db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
121447c6ae99SBarry Smith 
121547c6ae99SBarry Smith @*/
12167087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
121747c6ae99SBarry Smith {
121847c6ae99SBarry Smith   va_list                Argp;
121947c6ae99SBarry Smith   struct DMCompositeLink *next;
122047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
122171b14b3eSStefano Zampini   PetscBool              flg;
122247c6ae99SBarry Smith 
122347c6ae99SBarry Smith   PetscFunctionBegin;
122447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
12267a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
122747c6ae99SBarry Smith   next = com->next;
122847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
122947c6ae99SBarry Smith   va_start(Argp,dm);
123047c6ae99SBarry Smith   while (next) {
123147c6ae99SBarry Smith     Vec *vec;
123247c6ae99SBarry Smith     vec = va_arg(Argp, Vec*);
12339566063dSJacob Faibussowitsch     if (vec) PetscCall(DMRestoreLocalVector(next->dm,vec));
123447c6ae99SBarry Smith     next = next->next;
123547c6ae99SBarry Smith   }
123647c6ae99SBarry Smith   va_end(Argp);
123747c6ae99SBarry Smith   PetscFunctionReturn(0);
123847c6ae99SBarry Smith }
123947c6ae99SBarry Smith 
124047c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
124147c6ae99SBarry Smith /*@C
12429ae5db72SJed Brown     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
124347c6ae99SBarry Smith 
124447c6ae99SBarry Smith     Not Collective
124547c6ae99SBarry Smith 
124647c6ae99SBarry Smith     Input Parameter:
124747c6ae99SBarry Smith .    dm - the packer object
124847c6ae99SBarry Smith 
124947c6ae99SBarry Smith     Output Parameter:
12509ae5db72SJed Brown .   DM ... - the individual entries (DMs)
125147c6ae99SBarry Smith 
125247c6ae99SBarry Smith     Level: advanced
125347c6ae99SBarry Smith 
125495452b02SPatrick Sanan     Fortran Notes:
125595452b02SPatrick Sanan     Available as DMCompositeGetEntries() for one output DM, DMCompositeGetEntries2() for 2, etc
1256f5f57ec0SBarry Smith 
1257db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()`
1258db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1259db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`,
1260db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`
126147c6ae99SBarry Smith 
126247c6ae99SBarry Smith @*/
12637087cfbeSBarry Smith PetscErrorCode  DMCompositeGetEntries(DM dm,...)
126447c6ae99SBarry Smith {
126547c6ae99SBarry Smith   va_list                Argp;
126647c6ae99SBarry Smith   struct DMCompositeLink *next;
126747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
126871b14b3eSStefano Zampini   PetscBool              flg;
126947c6ae99SBarry Smith 
127047c6ae99SBarry Smith   PetscFunctionBegin;
127147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
12737a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
127447c6ae99SBarry Smith   next = com->next;
127547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
127647c6ae99SBarry Smith   va_start(Argp,dm);
127747c6ae99SBarry Smith   while (next) {
127847c6ae99SBarry Smith     DM *dmn;
127947c6ae99SBarry Smith     dmn = va_arg(Argp, DM*);
12809ae5db72SJed Brown     if (dmn) *dmn = next->dm;
128147c6ae99SBarry Smith     next = next->next;
128247c6ae99SBarry Smith   }
128347c6ae99SBarry Smith   va_end(Argp);
128447c6ae99SBarry Smith   PetscFunctionReturn(0);
128547c6ae99SBarry Smith }
128647c6ae99SBarry Smith 
1287dbab29e1SMark F. Adams /*@C
12882fa5ba8aSJed Brown     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
12892fa5ba8aSJed Brown 
12902fa5ba8aSJed Brown     Not Collective
12912fa5ba8aSJed Brown 
12922fa5ba8aSJed Brown     Input Parameter:
1293907376e6SBarry Smith .    dm - the packer object
1294907376e6SBarry Smith 
1295907376e6SBarry Smith     Output Parameter:
1296907376e6SBarry Smith .    dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs
12972fa5ba8aSJed Brown 
12982fa5ba8aSJed Brown     Level: advanced
12992fa5ba8aSJed Brown 
1300db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()`
1301db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1302db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`,
1303db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`
13042fa5ba8aSJed Brown 
13052fa5ba8aSJed Brown @*/
13062fa5ba8aSJed Brown PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[])
13072fa5ba8aSJed Brown {
13082fa5ba8aSJed Brown   struct DMCompositeLink *next;
13092fa5ba8aSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
13102fa5ba8aSJed Brown   PetscInt               i;
131171b14b3eSStefano Zampini   PetscBool              flg;
13122fa5ba8aSJed Brown 
13132fa5ba8aSJed Brown   PetscFunctionBegin;
13142fa5ba8aSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
13159566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg));
13167a8be351SBarry Smith   PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name);
13172fa5ba8aSJed Brown   /* loop over packed objects, handling one at at time */
13182fa5ba8aSJed Brown   for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm;
13192fa5ba8aSJed Brown   PetscFunctionReturn(0);
13202fa5ba8aSJed Brown }
13212fa5ba8aSJed Brown 
1322e10fd815SStefano Zampini typedef struct {
1323e10fd815SStefano Zampini   DM          dm;
1324e10fd815SStefano Zampini   PetscViewer *subv;
1325e10fd815SStefano Zampini   Vec         *vecs;
1326e10fd815SStefano Zampini } GLVisViewerCtx;
1327e10fd815SStefano Zampini 
1328e10fd815SStefano Zampini static PetscErrorCode  DestroyGLVisViewerCtx_Private(void *vctx)
1329e10fd815SStefano Zampini {
1330e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
1331e10fd815SStefano Zampini   PetscInt       i,n;
1332e10fd815SStefano Zampini 
1333e10fd815SStefano Zampini   PetscFunctionBegin;
13349566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(ctx->dm,&n));
1335e10fd815SStefano Zampini   for (i = 0; i < n; i++) {
13369566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&ctx->subv[i]));
1337e10fd815SStefano Zampini   }
13389566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ctx->subv,ctx->vecs));
13399566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dm));
13409566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
1341e10fd815SStefano Zampini   PetscFunctionReturn(0);
1342e10fd815SStefano Zampini }
1343e10fd815SStefano Zampini 
1344e10fd815SStefano Zampini static PetscErrorCode  DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
1345e10fd815SStefano Zampini {
1346e10fd815SStefano Zampini   Vec            X = (Vec)oX;
1347e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
1348e10fd815SStefano Zampini   PetscInt       i,n,cumf;
1349e10fd815SStefano Zampini 
1350e10fd815SStefano Zampini   PetscFunctionBegin;
13519566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(ctx->dm,&n));
13529566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccessArray(ctx->dm,X,n,NULL,ctx->vecs));
1353e10fd815SStefano Zampini   for (i = 0, cumf = 0; i < n; i++) {
1354e10fd815SStefano Zampini     PetscErrorCode (*g2l)(PetscObject,PetscInt,PetscObject[],void*);
1355e10fd815SStefano Zampini     void           *fctx;
1356e10fd815SStefano Zampini     PetscInt       nfi;
1357e10fd815SStefano Zampini 
13589566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i],&nfi,NULL,NULL,&g2l,NULL,&fctx));
1359e10fd815SStefano Zampini     if (!nfi) continue;
1360e10fd815SStefano Zampini     if (g2l) {
13619566063dSJacob Faibussowitsch       PetscCall((*g2l)((PetscObject)ctx->vecs[i],nfi,oXfield+cumf,fctx));
1362e10fd815SStefano Zampini     } else {
13639566063dSJacob Faibussowitsch       PetscCall(VecCopy(ctx->vecs[i],(Vec)(oXfield[cumf])));
1364e10fd815SStefano Zampini     }
1365e10fd815SStefano Zampini     cumf += nfi;
1366e10fd815SStefano Zampini   }
13679566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccessArray(ctx->dm,X,n,NULL,ctx->vecs));
1368e10fd815SStefano Zampini   PetscFunctionReturn(0);
1369e10fd815SStefano Zampini }
1370e10fd815SStefano Zampini 
1371e10fd815SStefano Zampini static PetscErrorCode  DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer)
1372e10fd815SStefano Zampini {
1373e10fd815SStefano Zampini   DM             dm = (DM)odm, *dms;
1374e10fd815SStefano Zampini   Vec            *Ufds;
1375e10fd815SStefano Zampini   GLVisViewerCtx *ctx;
1376e10fd815SStefano Zampini   PetscInt       i,n,tnf,*sdim;
1377e10fd815SStefano Zampini   char           **fecs;
1378e10fd815SStefano Zampini 
1379e10fd815SStefano Zampini   PetscFunctionBegin;
13809566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
13819566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
1382e10fd815SStefano Zampini   ctx->dm = dm;
13839566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm,&n));
13849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&dms));
13859566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetEntriesArray(dm,dms));
13869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n,&ctx->subv,n,&ctx->vecs));
1387e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1388e10fd815SStefano Zampini     PetscInt nf;
1389e10fd815SStefano Zampini 
13909566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PetscObjectComm(odm),&ctx->subv[i]));
13919566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(ctx->subv[i],PETSCVIEWERGLVIS));
13929566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisSetDM_Private(ctx->subv[i],(PetscObject)dms[i]));
13939566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,NULL,NULL,NULL,NULL,NULL));
1394e10fd815SStefano Zampini     tnf += nf;
1395e10fd815SStefano Zampini   }
13969566063dSJacob Faibussowitsch   PetscCall(PetscFree(dms));
13979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(tnf,&fecs,tnf,&sdim,tnf,&Ufds));
1398e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1399e10fd815SStefano Zampini     PetscInt   *sd,nf,f;
1400e10fd815SStefano Zampini     const char **fec;
1401e10fd815SStefano Zampini     Vec        *Uf;
1402e10fd815SStefano Zampini 
14039566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,&fec,&sd,NULL,(PetscObject**)&Uf,NULL));
1404e10fd815SStefano Zampini     for (f = 0; f < nf; f++) {
14059566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(fec[f],&fecs[tnf+f]));
1406e10fd815SStefano Zampini       Ufds[tnf+f] = Uf[f];
1407e10fd815SStefano Zampini       sdim[tnf+f] = sd[f];
1408e10fd815SStefano Zampini     }
1409e10fd815SStefano Zampini     tnf += nf;
1410e10fd815SStefano Zampini   }
14119566063dSJacob Faibussowitsch   PetscCall(PetscViewerGLVisSetFields(viewer,tnf,(const char**)fecs,sdim,DMCompositeSampleGLVisFields_Private,(PetscObject*)Ufds,ctx,DestroyGLVisViewerCtx_Private));
1412e10fd815SStefano Zampini   for (i = 0; i < tnf; i++) {
14139566063dSJacob Faibussowitsch     PetscCall(PetscFree(fecs[i]));
1414e10fd815SStefano Zampini   }
14159566063dSJacob Faibussowitsch   PetscCall(PetscFree3(fecs,sdim,Ufds));
1416e10fd815SStefano Zampini   PetscFunctionReturn(0);
1417e10fd815SStefano Zampini }
1418e10fd815SStefano Zampini 
14197087cfbeSBarry Smith PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
142047c6ae99SBarry Smith {
142147c6ae99SBarry Smith   struct DMCompositeLink *next;
142247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmi->data;
142347c6ae99SBarry Smith   DM                     dm;
142447c6ae99SBarry Smith 
142547c6ae99SBarry Smith   PetscFunctionBegin;
142647c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
1427ce94432eSBarry Smith   if (comm == MPI_COMM_NULL) {
14289566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)dmi,&comm));
1429ce94432eSBarry Smith   }
14309566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmi));
143147c6ae99SBarry Smith   next = com->next;
14329566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm,fine));
143347c6ae99SBarry Smith 
143447c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
143547c6ae99SBarry Smith   while (next) {
14369566063dSJacob Faibussowitsch     PetscCall(DMRefine(next->dm,comm,&dm));
14379566063dSJacob Faibussowitsch     PetscCall(DMCompositeAddDM(*fine,dm));
14389566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dm));
143947c6ae99SBarry Smith     next = next->next;
144047c6ae99SBarry Smith   }
144147c6ae99SBarry Smith   PetscFunctionReturn(0);
144247c6ae99SBarry Smith }
144347c6ae99SBarry Smith 
144414354c39SJed Brown PetscErrorCode  DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine)
144514354c39SJed Brown {
144614354c39SJed Brown   struct DMCompositeLink *next;
144714354c39SJed Brown   DM_Composite           *com = (DM_Composite*)dmi->data;
144814354c39SJed Brown   DM                     dm;
144914354c39SJed Brown 
145014354c39SJed Brown   PetscFunctionBegin;
145114354c39SJed Brown   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
14529566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmi));
14532ee06e3bSJed Brown   if (comm == MPI_COMM_NULL) {
14549566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)dmi,&comm));
145525296bd5SBarry Smith   }
145614354c39SJed Brown   next = com->next;
14579566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm,fine));
145814354c39SJed Brown 
145914354c39SJed Brown   /* loop over packed objects, handling one at at time */
146014354c39SJed Brown   while (next) {
14619566063dSJacob Faibussowitsch     PetscCall(DMCoarsen(next->dm,comm,&dm));
14629566063dSJacob Faibussowitsch     PetscCall(DMCompositeAddDM(*fine,dm));
14639566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dm));
146414354c39SJed Brown     next = next->next;
146514354c39SJed Brown   }
146614354c39SJed Brown   PetscFunctionReturn(0);
146714354c39SJed Brown }
146847c6ae99SBarry Smith 
1469e727c939SJed Brown PetscErrorCode  DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
147047c6ae99SBarry Smith {
14719ae5db72SJed Brown   PetscInt               m,n,M,N,nDM,i;
147247c6ae99SBarry Smith   struct DMCompositeLink *nextc;
147347c6ae99SBarry Smith   struct DMCompositeLink *nextf;
147425296bd5SBarry Smith   Vec                    gcoarse,gfine,*vecs;
147547c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
147647c6ae99SBarry Smith   DM_Composite           *comfine   = (DM_Composite*)fine->data;
14779ae5db72SJed Brown   Mat                    *mats;
147847c6ae99SBarry Smith 
147947c6ae99SBarry Smith   PetscFunctionBegin;
148047c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
148147c6ae99SBarry Smith   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
14829566063dSJacob Faibussowitsch   PetscCall(DMSetUp(coarse));
14839566063dSJacob Faibussowitsch   PetscCall(DMSetUp(fine));
148447c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14859566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(coarse,&gcoarse));
14869566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(fine,&gfine));
14879566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gcoarse,&n));
14889566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gfine,&m));
14899566063dSJacob Faibussowitsch   PetscCall(VecGetSize(gcoarse,&N));
14909566063dSJacob Faibussowitsch   PetscCall(VecGetSize(gfine,&M));
14919566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(coarse,&gcoarse));
14929566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(fine,&gfine));
149347c6ae99SBarry Smith 
14949ae5db72SJed Brown   nDM = comfine->nDM;
149563a3b9bcSJacob Faibussowitsch   PetscCheck(nDM == comcoarse->nDM,PetscObjectComm((PetscObject)fine),PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %" PetscInt_FMT " entries, but coarse has %" PetscInt_FMT,nDM,comcoarse->nDM);
14969566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nDM*nDM,&mats));
149725296bd5SBarry Smith   if (v) {
14989566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nDM,&vecs));
149925296bd5SBarry Smith   }
150047c6ae99SBarry Smith 
150147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
15029ae5db72SJed Brown   for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) {
150325296bd5SBarry Smith     if (!v) {
15049566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL));
150525296bd5SBarry Smith     } else {
15069566063dSJacob Faibussowitsch       PetscCall(DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]));
150725296bd5SBarry Smith     }
150847c6ae99SBarry Smith   }
15099566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A));
151025296bd5SBarry Smith   if (v) {
15119566063dSJacob Faibussowitsch     PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v));
151225296bd5SBarry Smith   }
15139566063dSJacob Faibussowitsch   for (i=0; i<nDM*nDM; i++) PetscCall(MatDestroy(&mats[i]));
15149566063dSJacob Faibussowitsch   PetscCall(PetscFree(mats));
151525296bd5SBarry Smith   if (v) {
15169566063dSJacob Faibussowitsch     for (i=0; i<nDM; i++) PetscCall(VecDestroy(&vecs[i]));
15179566063dSJacob Faibussowitsch     PetscCall(PetscFree(vecs));
151825296bd5SBarry Smith   }
151947c6ae99SBarry Smith   PetscFunctionReturn(0);
152047c6ae99SBarry Smith }
152147c6ae99SBarry Smith 
1522184d77edSJed Brown static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
15231411c6eeSJed Brown {
15241411c6eeSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
15251411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1526f7efa3c7SJed Brown   PetscInt               i;
15271411c6eeSJed Brown 
15281411c6eeSJed Brown   PetscFunctionBegin;
15291411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
15309566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetISLocalToGlobalMappings(dm,&ltogs));
15319566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap));
15329566063dSJacob Faibussowitsch   for (i=0; i<com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(&ltogs[i]));
15339566063dSJacob Faibussowitsch   PetscCall(PetscFree(ltogs));
15341411c6eeSJed Brown   PetscFunctionReturn(0);
15351411c6eeSJed Brown }
15361411c6eeSJed Brown 
1537b412c318SBarry Smith PetscErrorCode  DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring)
153847c6ae99SBarry Smith {
153947c6ae99SBarry Smith   PetscInt        n,i,cnt;
154047c6ae99SBarry Smith   ISColoringValue *colors;
154147c6ae99SBarry Smith   PetscBool       dense  = PETSC_FALSE;
154247c6ae99SBarry Smith   ISColoringValue maxcol = 0;
154347c6ae99SBarry Smith   DM_Composite    *com   = (DM_Composite*)dm->data;
154447c6ae99SBarry Smith 
154547c6ae99SBarry Smith   PetscFunctionBegin;
154647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
154708401ef6SPierre Jolivet   PetscCheck(ctype != IS_COLORING_LOCAL,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported");
1548e3247f34SBarry Smith   else if (ctype == IS_COLORING_GLOBAL) {
154947c6ae99SBarry Smith     n = com->n;
1550ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
15519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&colors)); /* freed in ISColoringDestroy() */
155247c6ae99SBarry Smith 
15539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL));
155447c6ae99SBarry Smith   if (dense) {
155547c6ae99SBarry Smith     for (i=0; i<n; i++) {
155647c6ae99SBarry Smith       colors[i] = (ISColoringValue)(com->rstart + i);
155747c6ae99SBarry Smith     }
155847c6ae99SBarry Smith     maxcol = com->N;
155947c6ae99SBarry Smith   } else {
156047c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
156147c6ae99SBarry Smith     PetscMPIInt            rank;
156247c6ae99SBarry Smith 
15639566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
156447c6ae99SBarry Smith     cnt  = 0;
156547c6ae99SBarry Smith     while (next) {
156647c6ae99SBarry Smith       ISColoring lcoloring;
156747c6ae99SBarry Smith 
15689566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring));
156947c6ae99SBarry Smith       for (i=0; i<lcoloring->N; i++) {
157047c6ae99SBarry Smith         colors[cnt++] = maxcol + lcoloring->colors[i];
157147c6ae99SBarry Smith       }
157247c6ae99SBarry Smith       maxcol += lcoloring->n;
15739566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&lcoloring));
157447c6ae99SBarry Smith       next    = next->next;
157547c6ae99SBarry Smith     }
157647c6ae99SBarry Smith   }
15779566063dSJacob Faibussowitsch   PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,PETSC_OWN_POINTER,coloring));
157847c6ae99SBarry Smith   PetscFunctionReturn(0);
157947c6ae99SBarry Smith }
158047c6ae99SBarry Smith 
15817087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
158247c6ae99SBarry Smith {
158347c6ae99SBarry Smith   struct DMCompositeLink *next;
158447c6ae99SBarry Smith   PetscScalar            *garray,*larray;
158547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
158647c6ae99SBarry Smith 
158747c6ae99SBarry Smith   PetscFunctionBegin;
158847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
158947c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
159039d35262SVincent Le Chenadec 
159147c6ae99SBarry Smith   if (!com->setup) {
15929566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
159347c6ae99SBarry Smith   }
159439d35262SVincent Le Chenadec 
15959566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gvec,&garray));
15969566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lvec,&larray));
159747c6ae99SBarry Smith 
159847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
159939d35262SVincent Le Chenadec   next = com->next;
160047c6ae99SBarry Smith   while (next) {
160147c6ae99SBarry Smith     Vec      local,global;
160247c6ae99SBarry Smith     PetscInt N;
160347c6ae99SBarry Smith 
16049566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm,&global));
16059566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(global,&N));
16069566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(global,garray));
16079566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm,&local));
16089566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local,larray));
16099566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(next->dm,global,mode,local));
16109566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(next->dm,global,mode,local));
16119566063dSJacob Faibussowitsch     PetscCall(VecResetArray(global));
16129566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local));
16139566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm,&global));
16149566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm,&local));
161539d35262SVincent Le Chenadec 
161606ebdd98SJed Brown     larray += next->nlocal;
161739d35262SVincent Le Chenadec     garray += next->n;
161847c6ae99SBarry Smith     next    = next->next;
161947c6ae99SBarry Smith   }
162047c6ae99SBarry Smith 
16219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gvec,NULL));
16229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lvec,NULL));
162347c6ae99SBarry Smith   PetscFunctionReturn(0);
162447c6ae99SBarry Smith }
162547c6ae99SBarry Smith 
16267087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
16270c010503SBarry Smith {
16280c010503SBarry Smith   PetscFunctionBegin;
162939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
163039d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
163139d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,4);
163239d35262SVincent Le Chenadec   PetscFunctionReturn(0);
163339d35262SVincent Le Chenadec }
163439d35262SVincent Le Chenadec 
163539d35262SVincent Le Chenadec PetscErrorCode  DMLocalToGlobalBegin_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
163639d35262SVincent Le Chenadec {
163739d35262SVincent Le Chenadec   struct DMCompositeLink *next;
163839d35262SVincent Le Chenadec   PetscScalar            *larray,*garray;
163939d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite*)dm->data;
164039d35262SVincent Le Chenadec 
164139d35262SVincent Le Chenadec   PetscFunctionBegin;
164239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
164339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
164439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
164539d35262SVincent Le Chenadec 
164639d35262SVincent Le Chenadec   if (!com->setup) {
16479566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
164839d35262SVincent Le Chenadec   }
164939d35262SVincent Le Chenadec 
16509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lvec,&larray));
16519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gvec,&garray));
165239d35262SVincent Le Chenadec 
165339d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
165439d35262SVincent Le Chenadec   next = com->next;
165539d35262SVincent Le Chenadec   while (next) {
165639d35262SVincent Le Chenadec     Vec      global,local;
165739d35262SVincent Le Chenadec 
16589566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm,&local));
16599566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local,larray));
16609566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm,&global));
16619566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(global,garray));
16629566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(next->dm,local,mode,global));
16639566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(next->dm,local,mode,global));
16649566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local));
16659566063dSJacob Faibussowitsch     PetscCall(VecResetArray(global));
16669566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm,&global));
16679566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm,&local));
166839d35262SVincent Le Chenadec 
166939d35262SVincent Le Chenadec     garray += next->n;
167039d35262SVincent Le Chenadec     larray += next->nlocal;
167139d35262SVincent Le Chenadec     next    = next->next;
167239d35262SVincent Le Chenadec   }
167339d35262SVincent Le Chenadec 
16749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gvec,NULL));
16759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lvec,NULL));
167639d35262SVincent Le Chenadec 
167739d35262SVincent Le Chenadec   PetscFunctionReturn(0);
167839d35262SVincent Le Chenadec }
167939d35262SVincent Le Chenadec 
168039d35262SVincent Le Chenadec PetscErrorCode  DMLocalToGlobalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
168139d35262SVincent Le Chenadec {
168239d35262SVincent Le Chenadec   PetscFunctionBegin;
168339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
168439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
168539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
168639d35262SVincent Le Chenadec   PetscFunctionReturn(0);
168739d35262SVincent Le Chenadec }
168839d35262SVincent Le Chenadec 
168939d35262SVincent Le Chenadec PetscErrorCode  DMLocalToLocalBegin_Composite(DM dm,Vec vec1,InsertMode mode,Vec vec2)
169039d35262SVincent Le Chenadec {
169139d35262SVincent Le Chenadec   struct DMCompositeLink *next;
169239d35262SVincent Le Chenadec   PetscScalar            *array1,*array2;
169339d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite*)dm->data;
169439d35262SVincent Le Chenadec 
169539d35262SVincent Le Chenadec   PetscFunctionBegin;
169639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
169739d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec1,VEC_CLASSID,2);
169839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec2,VEC_CLASSID,4);
169939d35262SVincent Le Chenadec 
170039d35262SVincent Le Chenadec   if (!com->setup) {
17019566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
170239d35262SVincent Le Chenadec   }
170339d35262SVincent Le Chenadec 
17049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec1,&array1));
17059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec2,&array2));
170639d35262SVincent Le Chenadec 
170739d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
170839d35262SVincent Le Chenadec   next = com->next;
170939d35262SVincent Le Chenadec   while (next) {
171039d35262SVincent Le Chenadec     Vec      local1,local2;
171139d35262SVincent Le Chenadec 
17129566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm,&local1));
17139566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local1,array1));
17149566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm,&local2));
17159566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local2,array2));
17169566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalBegin(next->dm,local1,mode,local2));
17179566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalEnd(next->dm,local1,mode,local2));
17189566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local2));
17199566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm,&local2));
17209566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local1));
17219566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm,&local1));
172239d35262SVincent Le Chenadec 
172339d35262SVincent Le Chenadec     array1 += next->nlocal;
172439d35262SVincent Le Chenadec     array2 += next->nlocal;
172539d35262SVincent Le Chenadec     next    = next->next;
172639d35262SVincent Le Chenadec   }
172739d35262SVincent Le Chenadec 
17289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec1,NULL));
17299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec2,NULL));
173039d35262SVincent Le Chenadec 
173139d35262SVincent Le Chenadec   PetscFunctionReturn(0);
173239d35262SVincent Le Chenadec }
173339d35262SVincent Le Chenadec 
173439d35262SVincent Le Chenadec PetscErrorCode  DMLocalToLocalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)
173539d35262SVincent Le Chenadec {
173639d35262SVincent Le Chenadec   PetscFunctionBegin;
173739d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
173839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec,VEC_CLASSID,2);
173939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec,VEC_CLASSID,4);
17400c010503SBarry Smith   PetscFunctionReturn(0);
17410c010503SBarry Smith }
174247c6ae99SBarry Smith 
17436ae3a549SBarry Smith /*MC
17446ae3a549SBarry Smith    DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs
17456ae3a549SBarry Smith 
17466ae3a549SBarry Smith   Level: intermediate
17476ae3a549SBarry Smith 
1748db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()`
17496ae3a549SBarry Smith M*/
17506ae3a549SBarry Smith 
17518cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1752a4121054SBarry Smith {
1753a4121054SBarry Smith   DM_Composite   *com;
1754a4121054SBarry Smith 
1755a4121054SBarry Smith   PetscFunctionBegin;
17569566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(p,&com));
1757a4121054SBarry Smith   p->data       = com;
1758a4121054SBarry Smith   com->n        = 0;
17597ac2b803SAlex Fikl   com->nghost   = 0;
17600298fd71SBarry Smith   com->next     = NULL;
1761a4121054SBarry Smith   com->nDM      = 0;
1762a4121054SBarry Smith 
1763a4121054SBarry Smith   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1764a4121054SBarry Smith   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
1765184d77edSJed Brown   p->ops->getlocaltoglobalmapping         = DMGetLocalToGlobalMapping_Composite;
17664d343eeaSMatthew G Knepley   p->ops->createfieldis                   = DMCreateFieldIS_Composite;
176716621825SDmitry Karpeev   p->ops->createfielddecomposition        = DMCreateFieldDecomposition_Composite;
1768a4121054SBarry Smith   p->ops->refine                          = DMRefine_Composite;
176914354c39SJed Brown   p->ops->coarsen                         = DMCoarsen_Composite;
177025296bd5SBarry Smith   p->ops->createinterpolation             = DMCreateInterpolation_Composite;
177125296bd5SBarry Smith   p->ops->creatematrix                    = DMCreateMatrix_Composite;
1772e727c939SJed Brown   p->ops->getcoloring                     = DMCreateColoring_Composite;
1773a4121054SBarry Smith   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1774a4121054SBarry Smith   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
177539d35262SVincent Le Chenadec   p->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Composite;
177639d35262SVincent Le Chenadec   p->ops->localtoglobalend                = DMLocalToGlobalEnd_Composite;
177739d35262SVincent Le Chenadec   p->ops->localtolocalbegin               = DMLocalToLocalBegin_Composite;
177839d35262SVincent Le Chenadec   p->ops->localtolocalend                 = DMLocalToLocalEnd_Composite;
1779a4121054SBarry Smith   p->ops->destroy                         = DMDestroy_Composite;
1780a4121054SBarry Smith   p->ops->view                            = DMView_Composite;
1781a4121054SBarry Smith   p->ops->setup                           = DMSetUp_Composite;
1782e10fd815SStefano Zampini 
17839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)p,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Composite));
1784a4121054SBarry Smith   PetscFunctionReturn(0);
1785a4121054SBarry Smith }
1786a4121054SBarry Smith 
17871fd49c25SBarry Smith /*@
17880c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
17890c010503SBarry Smith       vectors made up of several subvectors.
17900c010503SBarry Smith 
1791d083f849SBarry Smith     Collective
179247c6ae99SBarry Smith 
179347c6ae99SBarry Smith     Input Parameter:
17940c010503SBarry Smith .   comm - the processors that will share the global vector
17950c010503SBarry Smith 
17960c010503SBarry Smith     Output Parameters:
17970c010503SBarry Smith .   packer - the packer object
179847c6ae99SBarry Smith 
179947c6ae99SBarry Smith     Level: advanced
180047c6ae99SBarry Smith 
1801*c2e3fba1SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCOMPOSITE`, `DMCreate()`
1802db781477SPatrick Sanan          `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`
1803db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
180447c6ae99SBarry Smith 
180547c6ae99SBarry Smith @*/
18067087cfbeSBarry Smith PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
180747c6ae99SBarry Smith {
180847c6ae99SBarry Smith   PetscFunctionBegin;
18090c010503SBarry Smith   PetscValidPointer(packer,2);
18109566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm,packer));
18119566063dSJacob Faibussowitsch   PetscCall(DMSetType(*packer,DMCOMPOSITE));
181247c6ae99SBarry Smith   PetscFunctionReturn(0);
181347c6ae99SBarry Smith }
1814