xref: /petsc/src/dm/impls/composite/pack.c (revision 6497c311e7b976d467be1503c1effce92a60525c)
1ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h> /*I  "petscdmcomposite.h"  I*/
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>
3e10fd815SStefano Zampini #include <petsc/private/glvisviewerimpl.h>
42764a2aaSMatthew G. Knepley #include <petscds.h>
547c6ae99SBarry Smith 
647c6ae99SBarry Smith /*@C
747c6ae99SBarry Smith   DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
8dce8aebaSBarry Smith   separate components `DM` in a `DMCOMPOSITE` to build the correct matrix nonzero structure.
947c6ae99SBarry Smith 
1020f4b53cSBarry Smith   Logically Collective; No Fortran Support
1147c6ae99SBarry Smith 
12d8d19677SJose E. Roman   Input Parameters:
1347c6ae99SBarry Smith + dm                  - the composite object
1460225df5SJacob Faibussowitsch - FormCoupleLocations - routine to set the nonzero locations in the matrix
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   Level: advanced
1747c6ae99SBarry Smith 
18dce8aebaSBarry Smith   Note:
19dce8aebaSBarry Smith   See `DMSetApplicationContext()` and `DMGetApplicationContext()` for how to get user information into
2047c6ae99SBarry Smith   this routine
2147c6ae99SBarry Smith 
22dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`
2347c6ae99SBarry Smith @*/
24d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeSetCoupling(DM dm, PetscErrorCode (*FormCoupleLocations)(DM, Mat, PetscInt *, PetscInt *, PetscInt, PetscInt, PetscInt, PetscInt))
25d71ae5a4SJacob Faibussowitsch {
2647c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
2771b14b3eSStefano Zampini   PetscBool     flg;
2847c6ae99SBarry Smith 
2947c6ae99SBarry Smith   PetscFunctionBegin;
309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
317a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
3247c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3447c6ae99SBarry Smith }
3547c6ae99SBarry Smith 
3666976f2fSJacob Faibussowitsch static PetscErrorCode DMDestroy_Composite(DM dm)
37d71ae5a4SJacob Faibussowitsch {
3847c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
3947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
4047c6ae99SBarry Smith 
4147c6ae99SBarry Smith   PetscFunctionBegin;
4247c6ae99SBarry Smith   next = com->next;
4347c6ae99SBarry Smith   while (next) {
4447c6ae99SBarry Smith     prev = next;
4547c6ae99SBarry Smith     next = next->next;
469566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&prev->dm));
479566063dSJacob Faibussowitsch     PetscCall(PetscFree(prev->grstarts));
489566063dSJacob Faibussowitsch     PetscCall(PetscFree(prev));
4947c6ae99SBarry Smith   }
509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL));
51435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
529566063dSJacob Faibussowitsch   PetscCall(PetscFree(com));
533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5447c6ae99SBarry Smith }
5547c6ae99SBarry Smith 
5666976f2fSJacob Faibussowitsch static PetscErrorCode DMView_Composite(DM dm, PetscViewer v)
57d71ae5a4SJacob Faibussowitsch {
5847c6ae99SBarry Smith   PetscBool     iascii;
5947c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
6047c6ae99SBarry Smith 
6147c6ae99SBarry Smith   PetscFunctionBegin;
629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii));
6347c6ae99SBarry Smith   if (iascii) {
6447c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
6547c6ae99SBarry Smith     PetscInt                i;
6647c6ae99SBarry Smith 
679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "DM (%s)\n", ((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix"));
6863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "  contains %" PetscInt_FMT " DMs\n", com->nDM));
699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(v));
7047c6ae99SBarry Smith     for (i = 0; lnk; lnk = lnk->next, i++) {
7163a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(v, "Link %" PetscInt_FMT ": DM of type %s\n", i, ((PetscObject)lnk->dm)->type_name));
729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(v));
739566063dSJacob Faibussowitsch       PetscCall(DMView(lnk->dm, v));
749566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(v));
7547c6ae99SBarry Smith     }
769566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(v));
7747c6ae99SBarry Smith   }
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7947c6ae99SBarry Smith }
8047c6ae99SBarry Smith 
8147c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
8266976f2fSJacob Faibussowitsch static PetscErrorCode DMSetUp_Composite(DM dm)
83d71ae5a4SJacob Faibussowitsch {
8447c6ae99SBarry Smith   PetscInt                nprev = 0;
8547c6ae99SBarry Smith   PetscMPIInt             rank, size;
8647c6ae99SBarry Smith   DM_Composite           *com  = (DM_Composite *)dm->data;
8747c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
8847c6ae99SBarry Smith   PetscLayout             map;
8947c6ae99SBarry Smith 
9047c6ae99SBarry Smith   PetscFunctionBegin;
9128b400f6SJacob Faibussowitsch   PetscCheck(!com->setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Packer has already been setup");
929566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &map));
939566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(map, com->n));
949566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(map, PETSC_DETERMINE));
959566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(map, 1));
969566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(map));
979566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(map, &com->N));
989566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(map, &com->rstart, NULL));
999566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&map));
10047c6ae99SBarry Smith 
1019ae5db72SJed Brown   /* now set the rstart for each linked vector */
1029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
10447c6ae99SBarry Smith   while (next) {
10547c6ae99SBarry Smith     next->rstart = nprev;
10606ebdd98SJed Brown     nprev += next->n;
10747c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
1089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &next->grstarts));
1099566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgather(&next->grstart, 1, MPIU_INT, next->grstarts, 1, MPIU_INT, PetscObjectComm((PetscObject)dm)));
11047c6ae99SBarry Smith     next = next->next;
11147c6ae99SBarry Smith   }
11247c6ae99SBarry Smith   com->setup = PETSC_TRUE;
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11447c6ae99SBarry Smith }
11547c6ae99SBarry Smith 
11673e31fe2SJed Brown /*@
117aaa8cc7dSPierre Jolivet   DMCompositeGetNumberDM - Gets the number of `DM` objects in the `DMCOMPOSITE`
11847c6ae99SBarry Smith   representation.
11947c6ae99SBarry Smith 
12047c6ae99SBarry Smith   Not Collective
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith   Input Parameter:
123dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
12447c6ae99SBarry Smith 
12547c6ae99SBarry Smith   Output Parameter:
126dce8aebaSBarry Smith . nDM - the number of `DM`
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith   Level: beginner
12947c6ae99SBarry Smith 
130dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`
13147c6ae99SBarry Smith @*/
132d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetNumberDM(DM dm, PetscInt *nDM)
133d71ae5a4SJacob Faibussowitsch {
13447c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
13571b14b3eSStefano Zampini   PetscBool     flg;
1365fd66863SKarl Rupp 
13747c6ae99SBarry Smith   PetscFunctionBegin;
13847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1399566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1407a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
14147c6ae99SBarry Smith   *nDM = com->nDM;
1423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14347c6ae99SBarry Smith }
14447c6ae99SBarry Smith 
14547c6ae99SBarry Smith /*@C
14647c6ae99SBarry Smith   DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
14747c6ae99SBarry Smith   representation.
14847c6ae99SBarry Smith 
14920f4b53cSBarry Smith   Collective
15047c6ae99SBarry Smith 
1519ae5db72SJed Brown   Input Parameters:
152dce8aebaSBarry Smith + dm   - the `DMCOMPOSITE` object
1539ae5db72SJed Brown - gvec - the global vector
1549ae5db72SJed Brown 
1552fe279fdSBarry Smith   Output Parameter:
156a4e35b19SJacob Faibussowitsch . ... - the packed parallel vectors, `NULL` for those that are not needed
15747c6ae99SBarry Smith 
15847c6ae99SBarry Smith   Level: advanced
15947c6ae99SBarry Smith 
160dce8aebaSBarry Smith   Note:
161dce8aebaSBarry Smith   Use `DMCompositeRestoreAccess()` to return the vectors when you no longer need them
162dce8aebaSBarry Smith 
16360225df5SJacob Faibussowitsch   Fortran Notes:
164dce8aebaSBarry Smith   Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
165dce8aebaSBarry Smith   or use the alternative interface `DMCompositeGetAccessArray()`.
166dce8aebaSBarry Smith 
167dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetEntries()`, `DMCompositeScatter()`
16847c6ae99SBarry Smith @*/
169d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccess(DM dm, Vec gvec, ...)
170d71ae5a4SJacob Faibussowitsch {
17147c6ae99SBarry Smith   va_list                 Argp;
17247c6ae99SBarry Smith   struct DMCompositeLink *next;
17347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
1745edff71fSBarry Smith   PetscInt                readonly;
17571b14b3eSStefano Zampini   PetscBool               flg;
17647c6ae99SBarry Smith 
17747c6ae99SBarry Smith   PetscFunctionBegin;
17847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
17947c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
1809566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1817a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
18247c6ae99SBarry Smith   next = com->next;
18348a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
18447c6ae99SBarry Smith 
1859566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec, &readonly));
18615229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
18747c6ae99SBarry Smith   va_start(Argp, gvec);
18847c6ae99SBarry Smith   while (next) {
18947c6ae99SBarry Smith     Vec *vec;
19047c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
1919ae5db72SJed Brown     if (vec) {
1929566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, vec));
1935edff71fSBarry Smith       if (readonly) {
1945edff71fSBarry Smith         const PetscScalar *array;
1959566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(gvec, &array));
1969566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec, array + next->rstart));
1979566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(*vec));
1989566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(gvec, &array));
1995edff71fSBarry Smith       } else {
2005edff71fSBarry Smith         PetscScalar *array;
2019566063dSJacob Faibussowitsch         PetscCall(VecGetArray(gvec, &array));
2029566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec, array + next->rstart));
2039566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(gvec, &array));
20447c6ae99SBarry Smith       }
2055edff71fSBarry Smith     }
20647c6ae99SBarry Smith     next = next->next;
20747c6ae99SBarry Smith   }
20847c6ae99SBarry Smith   va_end(Argp);
2093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21047c6ae99SBarry Smith }
21147c6ae99SBarry Smith 
2125d83a8b1SBarry Smith /*@
213f73e5cebSJed Brown   DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
214f73e5cebSJed Brown   representation.
215f73e5cebSJed Brown 
21620f4b53cSBarry Smith   Collective
217f73e5cebSJed Brown 
218f73e5cebSJed Brown   Input Parameters:
219dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE`
220f73e5cebSJed Brown . pvec    - packed vector
221f73e5cebSJed Brown . nwanted - number of vectors wanted
222f13dfd9eSBarry Smith - wanted  - sorted array of vectors wanted, or `NULL` to get all vectors, length `nwanted`
223f73e5cebSJed Brown 
2242fe279fdSBarry Smith   Output Parameter:
225f13dfd9eSBarry Smith . vecs - array of requested global vectors (must be previously allocated and of length `nwanted`)
226f73e5cebSJed Brown 
227f73e5cebSJed Brown   Level: advanced
228f73e5cebSJed Brown 
229dce8aebaSBarry Smith   Note:
230dce8aebaSBarry Smith   Use `DMCompositeRestoreAccessArray()` to return the vectors when you no longer need them
231dce8aebaSBarry Smith 
232dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
233f73e5cebSJed Brown @*/
2345d83a8b1SBarry Smith PetscErrorCode DMCompositeGetAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[])
235d71ae5a4SJacob Faibussowitsch {
236f73e5cebSJed Brown   struct DMCompositeLink *link;
237f73e5cebSJed Brown   PetscInt                i, wnum;
238f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
239bee642f7SBarry Smith   PetscInt                readonly;
24071b14b3eSStefano Zampini   PetscBool               flg;
241f73e5cebSJed Brown 
242f73e5cebSJed Brown   PetscFunctionBegin;
243f73e5cebSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
244f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
2459566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
2467a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
24748a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
248f73e5cebSJed Brown 
2499566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
250f73e5cebSJed Brown   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
251f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
252f73e5cebSJed Brown       Vec v;
2539566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(link->dm, &v));
254bee642f7SBarry Smith       if (readonly) {
255bee642f7SBarry Smith         const PetscScalar *array;
2569566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec, &array));
2579566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + link->rstart));
2589566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
2599566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec, &array));
260bee642f7SBarry Smith       } else {
261bee642f7SBarry Smith         PetscScalar *array;
2629566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec, &array));
2639566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + link->rstart));
2649566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec, &array));
265bee642f7SBarry Smith       }
266f73e5cebSJed Brown       vecs[wnum++] = v;
267f73e5cebSJed Brown     }
268f73e5cebSJed Brown   }
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
270f73e5cebSJed Brown }
271f73e5cebSJed Brown 
2725d83a8b1SBarry Smith /*@
2737ac2b803SAlex Fikl   DMCompositeGetLocalAccessArray - Allows one to access the individual
2747ac2b803SAlex Fikl   packed vectors in their local representation.
2757ac2b803SAlex Fikl 
27620f4b53cSBarry Smith   Collective
2777ac2b803SAlex Fikl 
2787ac2b803SAlex Fikl   Input Parameters:
279dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE`
2807ac2b803SAlex Fikl . pvec    - packed vector
2817ac2b803SAlex Fikl . nwanted - number of vectors wanted
282f13dfd9eSBarry Smith - wanted  - sorted array of vectors wanted, or `NULL` to get all vectors, length `nwanted`
2837ac2b803SAlex Fikl 
2842fe279fdSBarry Smith   Output Parameter:
285f13dfd9eSBarry Smith . vecs - array of requested local vectors (must be allocated and of length `nwanted`)
2867ac2b803SAlex Fikl 
2877ac2b803SAlex Fikl   Level: advanced
2887ac2b803SAlex Fikl 
289dce8aebaSBarry Smith   Note:
290dce8aebaSBarry Smith   Use `DMCompositeRestoreLocalAccessArray()` to return the vectors
291dce8aebaSBarry Smith   when you no longer need them.
292dce8aebaSBarry Smith 
293dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`,
294db781477SPatrick Sanan           `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
2957ac2b803SAlex Fikl @*/
2965d83a8b1SBarry Smith PetscErrorCode DMCompositeGetLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[])
297d71ae5a4SJacob Faibussowitsch {
2987ac2b803SAlex Fikl   struct DMCompositeLink *link;
2997ac2b803SAlex Fikl   PetscInt                i, wnum;
3007ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite *)dm->data;
3017ac2b803SAlex Fikl   PetscInt                readonly;
3027ac2b803SAlex Fikl   PetscInt                nlocal = 0;
30371b14b3eSStefano Zampini   PetscBool               flg;
3047ac2b803SAlex Fikl 
3057ac2b803SAlex Fikl   PetscFunctionBegin;
3067ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3077ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
3089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3097a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
31048a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
3117ac2b803SAlex Fikl 
3129566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
3137ac2b803SAlex Fikl   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
3147ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
3157ac2b803SAlex Fikl       Vec v;
3169566063dSJacob Faibussowitsch       PetscCall(DMGetLocalVector(link->dm, &v));
3177ac2b803SAlex Fikl       if (readonly) {
3187ac2b803SAlex Fikl         const PetscScalar *array;
3199566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec, &array));
3209566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + nlocal));
321b1c3483dSMark 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
3229566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
3239566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec, &array));
3247ac2b803SAlex Fikl       } else {
3257ac2b803SAlex Fikl         PetscScalar *array;
3269566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec, &array));
3279566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + nlocal));
3289566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec, &array));
3297ac2b803SAlex Fikl       }
3307ac2b803SAlex Fikl       vecs[wnum++] = v;
3317ac2b803SAlex Fikl     }
3327ac2b803SAlex Fikl 
3337ac2b803SAlex Fikl     nlocal += link->nlocal;
3347ac2b803SAlex Fikl   }
3353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3367ac2b803SAlex Fikl }
3377ac2b803SAlex Fikl 
33847c6ae99SBarry Smith /*@C
339dce8aebaSBarry Smith   DMCompositeRestoreAccess - Returns the vectors obtained with `DMCompositeGetAccess()`
34047c6ae99SBarry Smith   representation.
34147c6ae99SBarry Smith 
34220f4b53cSBarry Smith   Collective
34347c6ae99SBarry Smith 
3449ae5db72SJed Brown   Input Parameters:
345dce8aebaSBarry Smith + dm   - the `DMCOMPOSITE` object
34647c6ae99SBarry Smith . gvec - the global vector
347a4e35b19SJacob Faibussowitsch - ...  - the individual parallel vectors, `NULL` for those that are not needed
34847c6ae99SBarry Smith 
34947c6ae99SBarry Smith   Level: advanced
35047c6ae99SBarry Smith 
351dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
352db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`,
35342747ad1SJacob Faibussowitsch          `DMCompositeGetAccess()`
35447c6ae99SBarry Smith @*/
355d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...)
356d71ae5a4SJacob Faibussowitsch {
35747c6ae99SBarry Smith   va_list                 Argp;
35847c6ae99SBarry Smith   struct DMCompositeLink *next;
35947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
3605edff71fSBarry Smith   PetscInt                readonly;
36171b14b3eSStefano Zampini   PetscBool               flg;
36247c6ae99SBarry Smith 
36347c6ae99SBarry Smith   PetscFunctionBegin;
36447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36547c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
3669566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3677a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
36847c6ae99SBarry Smith   next = com->next;
36948a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
37047c6ae99SBarry Smith 
3719566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec, &readonly));
37215229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
37347c6ae99SBarry Smith   va_start(Argp, gvec);
37447c6ae99SBarry Smith   while (next) {
37547c6ae99SBarry Smith     Vec *vec;
37647c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
3779ae5db72SJed Brown     if (vec) {
3789566063dSJacob Faibussowitsch       PetscCall(VecResetArray(*vec));
3791baa6e33SBarry Smith       if (readonly) PetscCall(VecLockReadPop(*vec));
3809566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, vec));
38147c6ae99SBarry Smith     }
38247c6ae99SBarry Smith     next = next->next;
38347c6ae99SBarry Smith   }
38447c6ae99SBarry Smith   va_end(Argp);
3853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38647c6ae99SBarry Smith }
38747c6ae99SBarry Smith 
3885d83a8b1SBarry Smith /*@
389dce8aebaSBarry Smith   DMCompositeRestoreAccessArray - Returns the vectors obtained with `DMCompositeGetAccessArray()`
390f73e5cebSJed Brown 
39120f4b53cSBarry Smith   Collective
392f73e5cebSJed Brown 
393f73e5cebSJed Brown   Input Parameters:
394dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE` object
395f73e5cebSJed Brown . pvec    - packed vector
396f73e5cebSJed Brown . nwanted - number of vectors wanted
397f13dfd9eSBarry Smith . wanted  - sorted array of vectors wanted, or `NULL` to restore all vectors
398f13dfd9eSBarry Smith - vecs    - array of global vectors
399f73e5cebSJed Brown 
400f73e5cebSJed Brown   Level: advanced
401f73e5cebSJed Brown 
402dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
403f73e5cebSJed Brown @*/
4045d83a8b1SBarry Smith PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[])
405d71ae5a4SJacob Faibussowitsch {
406f73e5cebSJed Brown   struct DMCompositeLink *link;
407f73e5cebSJed Brown   PetscInt                i, wnum;
408f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
409bee642f7SBarry Smith   PetscInt                readonly;
41071b14b3eSStefano Zampini   PetscBool               flg;
411f73e5cebSJed Brown 
412f73e5cebSJed Brown   PetscFunctionBegin;
413f73e5cebSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
414f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4159566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4167a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
41748a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
418f73e5cebSJed Brown 
4199566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
420f73e5cebSJed Brown   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
421f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
4229566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
42348a46eb9SPierre Jolivet       if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4249566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum]));
425f73e5cebSJed Brown       wnum++;
426f73e5cebSJed Brown     }
427f73e5cebSJed Brown   }
4283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
429f73e5cebSJed Brown }
430f73e5cebSJed Brown 
4315d83a8b1SBarry Smith /*@
432dce8aebaSBarry Smith   DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with `DMCompositeGetLocalAccessArray()`.
4337ac2b803SAlex Fikl 
43420f4b53cSBarry Smith   Collective
4357ac2b803SAlex Fikl 
4367ac2b803SAlex Fikl   Input Parameters:
437dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE` object
4387ac2b803SAlex Fikl . pvec    - packed vector
4397ac2b803SAlex Fikl . nwanted - number of vectors wanted
440f13dfd9eSBarry Smith . wanted  - sorted array of vectors wanted, or `NULL` to restore all vectors
441f13dfd9eSBarry Smith - vecs    - array of local vectors
4427ac2b803SAlex Fikl 
4437ac2b803SAlex Fikl   Level: advanced
4447ac2b803SAlex Fikl 
445dce8aebaSBarry Smith   Note:
446f13dfd9eSBarry Smith   `nwanted` and `wanted` must match the values given to `DMCompositeGetLocalAccessArray()`
4477ac2b803SAlex Fikl   otherwise the call will fail.
4487ac2b803SAlex Fikl 
449dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`,
450db781477SPatrick Sanan           `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`,
451db781477SPatrick Sanan           `DMCompositeScatter()`, `DMCompositeGather()`
4527ac2b803SAlex Fikl @*/
4535d83a8b1SBarry Smith PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec *vecs)
454d71ae5a4SJacob Faibussowitsch {
4557ac2b803SAlex Fikl   struct DMCompositeLink *link;
4567ac2b803SAlex Fikl   PetscInt                i, wnum;
4577ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite *)dm->data;
4587ac2b803SAlex Fikl   PetscInt                readonly;
45971b14b3eSStefano Zampini   PetscBool               flg;
4607ac2b803SAlex Fikl 
4617ac2b803SAlex Fikl   PetscFunctionBegin;
4627ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4637ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4657a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
46648a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
4677ac2b803SAlex Fikl 
4689566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
4697ac2b803SAlex Fikl   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
4707ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
4719566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
47248a46eb9SPierre Jolivet       if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4739566063dSJacob Faibussowitsch       PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum]));
4747ac2b803SAlex Fikl       wnum++;
4757ac2b803SAlex Fikl     }
4767ac2b803SAlex Fikl   }
4773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4787ac2b803SAlex Fikl }
4797ac2b803SAlex Fikl 
48047c6ae99SBarry Smith /*@C
48147c6ae99SBarry Smith   DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
48247c6ae99SBarry Smith 
48320f4b53cSBarry Smith   Collective
48447c6ae99SBarry Smith 
4859ae5db72SJed Brown   Input Parameters:
486dce8aebaSBarry Smith + dm   - the `DMCOMPOSITE` object
48747c6ae99SBarry Smith . gvec - the global vector
488a4e35b19SJacob Faibussowitsch - ...  - the individual sequential vectors, `NULL` for those that are not needed
48947c6ae99SBarry Smith 
49047c6ae99SBarry Smith   Level: advanced
49147c6ae99SBarry Smith 
492dce8aebaSBarry Smith   Note:
493dce8aebaSBarry Smith   `DMCompositeScatterArray()` is a non-variadic alternative that is often more convenient for library callers and is
4946f3c3dcfSJed Brown   accessible from Fortran.
4956f3c3dcfSJed Brown 
496dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
497db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
498db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
499db781477SPatrick Sanan          `DMCompositeScatterArray()`
50047c6ae99SBarry Smith @*/
501d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...)
502d71ae5a4SJacob Faibussowitsch {
50347c6ae99SBarry Smith   va_list                 Argp;
50447c6ae99SBarry Smith   struct DMCompositeLink *next;
505c599c493SJunchao Zhang   PETSC_UNUSED PetscInt   cnt;
50647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
50771b14b3eSStefano Zampini   PetscBool               flg;
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith   PetscFunctionBegin;
51047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
51147c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5129566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5137a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
51448a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
51547c6ae99SBarry Smith 
51615229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
51747c6ae99SBarry Smith   va_start(Argp, gvec);
5188fd8f222SJed Brown   for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
5199ae5db72SJed Brown     Vec local;
5209ae5db72SJed Brown     local = va_arg(Argp, Vec);
5219ae5db72SJed Brown     if (local) {
5229ae5db72SJed Brown       Vec                global;
5235edff71fSBarry Smith       const PetscScalar *array;
52463a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
5259566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
5269566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
5279566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
5289566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local));
5299566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local));
5309566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
5319566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5329566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
53347c6ae99SBarry Smith     }
53447c6ae99SBarry Smith   }
53547c6ae99SBarry Smith   va_end(Argp);
5363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53747c6ae99SBarry Smith }
53847c6ae99SBarry Smith 
5396f3c3dcfSJed Brown /*@
5406f3c3dcfSJed Brown   DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5416f3c3dcfSJed Brown 
54220f4b53cSBarry Smith   Collective
5436f3c3dcfSJed Brown 
5446f3c3dcfSJed Brown   Input Parameters:
545dce8aebaSBarry Smith + dm    - the `DMCOMPOSITE` object
5466f3c3dcfSJed Brown . gvec  - the global vector
547a2b725a8SWilliam Gropp - lvecs - array of local vectors, NULL for any that are not needed
5486f3c3dcfSJed Brown 
5496f3c3dcfSJed Brown   Level: advanced
5506f3c3dcfSJed Brown 
5516f3c3dcfSJed Brown   Note:
552dce8aebaSBarry Smith   This is a non-variadic alternative to `DMCompositeScatter()`
5536f3c3dcfSJed Brown 
554dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`
555db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
556db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
5576f3c3dcfSJed Brown @*/
558d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs)
559d71ae5a4SJacob Faibussowitsch {
5606f3c3dcfSJed Brown   struct DMCompositeLink *next;
5616f3c3dcfSJed Brown   PetscInt                i;
5626f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
56371b14b3eSStefano Zampini   PetscBool               flg;
5646f3c3dcfSJed Brown 
5656f3c3dcfSJed Brown   PetscFunctionBegin;
5666f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5676f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5697a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
57048a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
5716f3c3dcfSJed Brown 
57215229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
5736f3c3dcfSJed Brown   for (i = 0, next = com->next; next; next = next->next, i++) {
5746f3c3dcfSJed Brown     if (lvecs[i]) {
5756f3c3dcfSJed Brown       Vec                global;
576c5d31e75SLisandro Dalcin       const PetscScalar *array;
5776f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3);
5789566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
5799566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
5809566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart));
5819566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i]));
5829566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i]));
5839566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
5849566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5859566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
5866f3c3dcfSJed Brown     }
5876f3c3dcfSJed Brown   }
5883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5896f3c3dcfSJed Brown }
5906f3c3dcfSJed Brown 
59147c6ae99SBarry Smith /*@C
59247c6ae99SBarry Smith   DMCompositeGather - Gathers into a global packed vector from its individual local vectors
59347c6ae99SBarry Smith 
59420f4b53cSBarry Smith   Collective
59547c6ae99SBarry Smith 
596d8d19677SJose E. Roman   Input Parameters:
597dce8aebaSBarry Smith + dm    - the `DMCOMPOSITE` object
598dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES`
599a4e35b19SJacob Faibussowitsch . gvec  - the global vector
600a4e35b19SJacob Faibussowitsch - ...   - the individual sequential vectors, `NULL` for any that are not needed
60147c6ae99SBarry Smith 
60247c6ae99SBarry Smith   Level: advanced
60347c6ae99SBarry Smith 
60460225df5SJacob Faibussowitsch   Fortran Notes:
605dce8aebaSBarry Smith   Fortran users should use `DMCompositeGatherArray()`
606f5f57ec0SBarry Smith 
607dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
608db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
609db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
61047c6ae99SBarry Smith @*/
611d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...)
612d71ae5a4SJacob Faibussowitsch {
61347c6ae99SBarry Smith   va_list                 Argp;
61447c6ae99SBarry Smith   struct DMCompositeLink *next;
61547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
616c599c493SJunchao Zhang   PETSC_UNUSED PetscInt   cnt;
61771b14b3eSStefano Zampini   PetscBool               flg;
61847c6ae99SBarry Smith 
61947c6ae99SBarry Smith   PetscFunctionBegin;
62047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
621064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6229566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6237a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
62448a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
62547c6ae99SBarry Smith 
62615229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
6271dac896bSSatish Balay   va_start(Argp, gvec);
6288fd8f222SJed Brown   for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
6299ae5db72SJed Brown     Vec local;
6309ae5db72SJed Brown     local = va_arg(Argp, Vec);
6319ae5db72SJed Brown     if (local) {
63247c6ae99SBarry Smith       PetscScalar *array;
6339ae5db72SJed Brown       Vec          global;
63463a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
6359566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
6369566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec, &array));
6379566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
6389566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global));
6399566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global));
6409566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec, &array));
6419566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6429566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
64347c6ae99SBarry Smith     }
64447c6ae99SBarry Smith   }
64547c6ae99SBarry Smith   va_end(Argp);
6463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64747c6ae99SBarry Smith }
64847c6ae99SBarry Smith 
6496f3c3dcfSJed Brown /*@
6506f3c3dcfSJed Brown   DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6516f3c3dcfSJed Brown 
65220f4b53cSBarry Smith   Collective
6536f3c3dcfSJed Brown 
654d8d19677SJose E. Roman   Input Parameters:
655dce8aebaSBarry Smith + dm    - the `DMCOMPOSITE` object
6566f3c3dcfSJed Brown . gvec  - the global vector
657dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES`
6586f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed
6596f3c3dcfSJed Brown 
6606f3c3dcfSJed Brown   Level: advanced
6616f3c3dcfSJed Brown 
662dce8aebaSBarry Smith   Note:
663dce8aebaSBarry Smith   This is a non-variadic alternative to `DMCompositeGather()`.
6646f3c3dcfSJed Brown 
665dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
666db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
667db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`,
6686f3c3dcfSJed Brown @*/
669d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs)
670d71ae5a4SJacob Faibussowitsch {
6716f3c3dcfSJed Brown   struct DMCompositeLink *next;
6726f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
6736f3c3dcfSJed Brown   PetscInt                i;
67471b14b3eSStefano Zampini   PetscBool               flg;
6756f3c3dcfSJed Brown 
6766f3c3dcfSJed Brown   PetscFunctionBegin;
6776f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
678064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6799566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6807a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
68148a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
6826f3c3dcfSJed Brown 
68315229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
6846f3c3dcfSJed Brown   for (next = com->next, i = 0; next; next = next->next, i++) {
6856f3c3dcfSJed Brown     if (lvecs[i]) {
6866f3c3dcfSJed Brown       PetscScalar *array;
6876f3c3dcfSJed Brown       Vec          global;
688064a246eSJacob Faibussowitsch       PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4);
6899566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
6909566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec, &array));
6919566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
6929566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global));
6939566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global));
6949566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec, &array));
6959566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6969566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
6976f3c3dcfSJed Brown     }
6986f3c3dcfSJed Brown   }
6993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7006f3c3dcfSJed Brown }
7016f3c3dcfSJed Brown 
702f5f57ec0SBarry Smith /*@
703dce8aebaSBarry Smith   DMCompositeAddDM - adds a `DM` vector to a `DMCOMPOSITE`
70447c6ae99SBarry Smith 
70520f4b53cSBarry Smith   Collective
70647c6ae99SBarry Smith 
707d8d19677SJose E. Roman   Input Parameters:
708dce8aebaSBarry Smith + dmc - the  `DMCOMPOSITE` object
709dce8aebaSBarry Smith - dm  - the `DM` object
71047c6ae99SBarry Smith 
71147c6ae99SBarry Smith   Level: advanced
71247c6ae99SBarry Smith 
71342747ad1SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeGather()`, `DMCreateGlobalVector()`,
714db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
715db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
71647c6ae99SBarry Smith @*/
717d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeAddDM(DM dmc, DM dm)
718d71ae5a4SJacob Faibussowitsch {
71906ebdd98SJed Brown   PetscInt                n, nlocal;
72047c6ae99SBarry Smith   struct DMCompositeLink *mine, *next;
72106ebdd98SJed Brown   Vec                     global, local;
72247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dmc->data;
72371b14b3eSStefano Zampini   PetscBool               flg;
72447c6ae99SBarry Smith 
72547c6ae99SBarry Smith   PetscFunctionBegin;
72647c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc, DM_CLASSID, 1);
72747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
7289566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg));
7297a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
73047c6ae99SBarry Smith   next = com->next;
73128b400f6SJacob Faibussowitsch   PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite");
73247c6ae99SBarry Smith 
73347c6ae99SBarry Smith   /* create new link */
7349566063dSJacob Faibussowitsch   PetscCall(PetscNew(&mine));
7359566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
736cac35740SJunchao Zhang   PetscCall(DMCreateGlobalVector(dm, &global)); // Not using DMGetGlobalVector() since it will stash the vector with a type decided by dm,
737cac35740SJunchao Zhang   PetscCall(VecGetLocalSize(global, &n));       // while we may want to set dmc's vectype later, say via DMSetFromOptions(dmc), and we
738cac35740SJunchao Zhang   PetscCall(VecDestroy(&global));               // want to propagate the type to dm.
739cac35740SJunchao Zhang   PetscCall(DMCreateLocalVector(dm, &local));   // Not using DMGetLocalVector(), same reason as above.
7409566063dSJacob Faibussowitsch   PetscCall(VecGetSize(local, &nlocal));
741cac35740SJunchao Zhang   PetscCall(VecDestroy(&local));
7428865f1eaSKarl Rupp 
74347c6ae99SBarry Smith   mine->n      = n;
74406ebdd98SJed Brown   mine->nlocal = nlocal;
74547c6ae99SBarry Smith   mine->dm     = dm;
7460298fd71SBarry Smith   mine->next   = NULL;
74747c6ae99SBarry Smith   com->n += n;
7487ac2b803SAlex Fikl   com->nghost += nlocal;
74947c6ae99SBarry Smith 
75047c6ae99SBarry Smith   /* add to end of list */
7518865f1eaSKarl Rupp   if (!next) com->next = mine;
7528865f1eaSKarl Rupp   else {
75347c6ae99SBarry Smith     while (next->next) next = next->next;
75447c6ae99SBarry Smith     next->next = mine;
75547c6ae99SBarry Smith   }
75647c6ae99SBarry Smith   com->nDM++;
75747c6ae99SBarry Smith   com->nmine++;
7583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75947c6ae99SBarry Smith }
76047c6ae99SBarry Smith 
7619804daf3SBarry Smith #include <petscdraw.h>
76226887b52SJed Brown PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
76366976f2fSJacob Faibussowitsch static PetscErrorCode       VecView_DMComposite(Vec gvec, PetscViewer viewer)
764d71ae5a4SJacob Faibussowitsch {
76547c6ae99SBarry Smith   DM                      dm;
76647c6ae99SBarry Smith   struct DMCompositeLink *next;
76747c6ae99SBarry Smith   PetscBool               isdraw;
768cef07954SSatish Balay   DM_Composite           *com;
76947c6ae99SBarry Smith 
77047c6ae99SBarry Smith   PetscFunctionBegin;
7719566063dSJacob Faibussowitsch   PetscCall(VecGetDM(gvec, &dm));
7727a8be351SBarry Smith   PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite");
77347c6ae99SBarry Smith   com  = (DM_Composite *)dm->data;
77447c6ae99SBarry Smith   next = com->next;
77547c6ae99SBarry Smith 
7769566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
77747c6ae99SBarry Smith   if (!isdraw) {
77847c6ae99SBarry Smith     /* do I really want to call this? */
7799566063dSJacob Faibussowitsch     PetscCall(VecView_MPI(gvec, viewer));
78047c6ae99SBarry Smith   } else {
78147c6ae99SBarry Smith     PetscInt cnt = 0;
78247c6ae99SBarry Smith 
78315229ffcSPierre Jolivet     /* loop over packed objects, handling one at a time */
78447c6ae99SBarry Smith     while (next) {
78547c6ae99SBarry Smith       Vec                vec;
786ca4278abSLisandro Dalcin       const PetscScalar *array;
78747c6ae99SBarry Smith       PetscInt           bs;
78847c6ae99SBarry Smith 
7899ae5db72SJed Brown       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
7909566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &vec));
7919566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
7929566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart));
7939566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
7949566063dSJacob Faibussowitsch       PetscCall(VecView(vec, viewer));
7959566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vec));
7969566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(vec, &bs));
7979566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &vec));
7989566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawBaseAdd(viewer, bs));
79947c6ae99SBarry Smith       cnt += bs;
80047c6ae99SBarry Smith       next = next->next;
80147c6ae99SBarry Smith     }
8029566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt));
80347c6ae99SBarry Smith   }
8043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80547c6ae99SBarry Smith }
80647c6ae99SBarry Smith 
80766976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec)
808d71ae5a4SJacob Faibussowitsch {
80947c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
81047c6ae99SBarry Smith 
81147c6ae99SBarry Smith   PetscFunctionBegin;
81247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8139566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
8149566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8159566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec));
8169566063dSJacob Faibussowitsch   PetscCall(VecSetType(*gvec, dm->vectype));
8179566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*gvec, com->n, com->N));
8189566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*gvec, dm));
8199566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (void (*)(void))VecView_DMComposite));
8203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82147c6ae99SBarry Smith }
82247c6ae99SBarry Smith 
82366976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec)
824d71ae5a4SJacob Faibussowitsch {
82547c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
82647c6ae99SBarry Smith 
82747c6ae99SBarry Smith   PetscFunctionBegin;
82847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
82947c6ae99SBarry Smith   if (!com->setup) {
8309566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(dm));
8319566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
83247c6ae99SBarry Smith   }
8339566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, lvec));
8349566063dSJacob Faibussowitsch   PetscCall(VecSetType(*lvec, dm->vectype));
8359566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE));
8369566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*lvec, dm));
8373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83847c6ae99SBarry Smith }
83947c6ae99SBarry Smith 
84047c6ae99SBarry Smith /*@C
841dce8aebaSBarry Smith   DMCompositeGetISLocalToGlobalMappings - gets an `ISLocalToGlobalMapping` for each `DM` in the `DMCOMPOSITE`, maps to the composite global space
84247c6ae99SBarry Smith 
84320f4b53cSBarry Smith   Collective; No Fortran Support
84447c6ae99SBarry Smith 
84547c6ae99SBarry Smith   Input Parameter:
846dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
84747c6ae99SBarry Smith 
8482fe279fdSBarry Smith   Output Parameter:
8499ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes
850dce8aebaSBarry Smith            all the ghost points that individual ghosted `DMDA` may have.
85147c6ae99SBarry Smith 
85247c6ae99SBarry Smith   Level: advanced
85347c6ae99SBarry Smith 
854dce8aebaSBarry Smith   Note:
855f13dfd9eSBarry Smith   Each entry of `ltogs` should be destroyed with `ISLocalToGlobalMappingDestroy()`, `ltogs` should be freed with `PetscFree()`.
85647c6ae99SBarry Smith 
857dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
858db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
859c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
86047c6ae99SBarry Smith @*/
861cc4c1da9SBarry Smith PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping *ltogs[])
862d71ae5a4SJacob Faibussowitsch {
86347c6ae99SBarry Smith   PetscInt                i, *idx, n, cnt;
86447c6ae99SBarry Smith   struct DMCompositeLink *next;
86547c6ae99SBarry Smith   PetscMPIInt             rank;
86647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
86771b14b3eSStefano Zampini   PetscBool               flg;
86847c6ae99SBarry Smith 
86947c6ae99SBarry Smith   PetscFunctionBegin;
87047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8719566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
8727a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
8739566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM, ltogs));
87547c6ae99SBarry Smith   next = com->next;
8769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
87747c6ae99SBarry Smith 
87815229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
87947c6ae99SBarry Smith   cnt = 0;
88047c6ae99SBarry Smith   while (next) {
8816eb61c8cSJed Brown     ISLocalToGlobalMapping ltog;
8826eb61c8cSJed Brown     PetscMPIInt            size;
88386994e45SJed Brown     const PetscInt        *suboff, *indices;
8846eb61c8cSJed Brown     Vec                    global;
88547c6ae99SBarry Smith 
8866eb61c8cSJed Brown     /* Get sub-DM global indices for each local dof */
8879566063dSJacob Faibussowitsch     PetscCall(DMGetLocalToGlobalMapping(next->dm, &ltog));
8889566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n));
8899566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices));
8909566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &idx));
89147c6ae99SBarry Smith 
8926eb61c8cSJed Brown     /* Get the offsets for the sub-DM global vector */
8939566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
8949566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRanges(global, &suboff));
8959566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size));
8966eb61c8cSJed Brown 
8976eb61c8cSJed Brown     /* Shift the sub-DM definition of the global space to the composite global space */
8986eb61c8cSJed Brown     for (i = 0; i < n; i++) {
89986994e45SJed Brown       PetscInt subi = indices[i], lo = 0, hi = size, t;
900e333b035SStefano Zampini       /* There's no consensus on what a negative index means,
901e333b035SStefano Zampini          except for skipping when setting the values in vectors and matrices */
9029371c9d4SSatish Balay       if (subi < 0) {
9039371c9d4SSatish Balay         idx[i] = subi - next->grstarts[rank];
9049371c9d4SSatish Balay         continue;
9059371c9d4SSatish Balay       }
9066eb61c8cSJed Brown       /* Binary search to find which rank owns subi */
9076eb61c8cSJed Brown       while (hi - lo > 1) {
9086eb61c8cSJed Brown         t = lo + (hi - lo) / 2;
9096eb61c8cSJed Brown         if (suboff[t] > subi) hi = t;
9106eb61c8cSJed Brown         else lo = t;
9116eb61c8cSJed Brown       }
9126eb61c8cSJed Brown       idx[i] = subi - suboff[lo] + next->grstarts[lo];
9136eb61c8cSJed Brown     }
9149566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices));
9159566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt]));
9169566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
91747c6ae99SBarry Smith     next = next->next;
91847c6ae99SBarry Smith     cnt++;
91947c6ae99SBarry Smith   }
9203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92147c6ae99SBarry Smith }
92247c6ae99SBarry Smith 
92387c85e80SJed Brown /*@C
9249ae5db72SJed Brown   DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
92587c85e80SJed Brown 
92620f4b53cSBarry Smith   Not Collective; No Fortran Support
92787c85e80SJed Brown 
9284165533cSJose E. Roman   Input Parameter:
929dce8aebaSBarry Smith . dm - the `DMCOMPOSITE`
93087c85e80SJed Brown 
9314165533cSJose E. Roman   Output Parameter:
93215229ffcSPierre Jolivet . is - array of serial index sets for each component of the `DMCOMPOSITE`
93387c85e80SJed Brown 
93487c85e80SJed Brown   Level: intermediate
93587c85e80SJed Brown 
93687c85e80SJed Brown   Notes:
93787c85e80SJed Brown   At present, a composite local vector does not normally exist.  This function is used to provide index sets for
93815229ffcSPierre Jolivet   `MatGetLocalSubMatrix()`.  In the future, the scatters for each entry in the `DMCOMPOSITE` may be merged into a single
9399ae5db72SJed Brown   scatter to a composite local vector.  The user should not typically need to know which is being done.
94087c85e80SJed Brown 
941dce8aebaSBarry Smith   To get the composite global indices at all local points (including ghosts), use `DMCompositeGetISLocalToGlobalMappings()`.
94287c85e80SJed Brown 
943dce8aebaSBarry Smith   To get index sets for pieces of the composite global vector, use `DMCompositeGetGlobalISs()`.
94487c85e80SJed Brown 
945dce8aebaSBarry Smith   Each returned `IS` should be destroyed with `ISDestroy()`, the array should be freed with `PetscFree()`.
94687c85e80SJed Brown 
947f13dfd9eSBarry Smith   Fortran Note:
948f13dfd9eSBarry Smith   Pass in an array long enough to hold all the `IS`, see `DMCompositeGetNumberDM()`
949f13dfd9eSBarry Smith 
950f13dfd9eSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`,
951f13dfd9eSBarry Smith           `MatCreateLocalRef()`, `DMCompositeGetNumberDM()`
95287c85e80SJed Brown @*/
953cc4c1da9SBarry Smith PetscErrorCode DMCompositeGetLocalISs(DM dm, IS *is[])
954d71ae5a4SJacob Faibussowitsch {
95587c85e80SJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
95687c85e80SJed Brown   struct DMCompositeLink *link;
95787c85e80SJed Brown   PetscInt                cnt, start;
95871b14b3eSStefano Zampini   PetscBool               flg;
95987c85e80SJed Brown 
96087c85e80SJed Brown   PetscFunctionBegin;
96187c85e80SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9624f572ea9SToby Isaac   PetscAssertPointer(is, 2);
9639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
9647a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
9659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nmine, is));
96606ebdd98SJed Brown   for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) {
967520db06cSJed Brown     PetscInt bs;
9689566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt]));
9699566063dSJacob Faibussowitsch     PetscCall(DMGetBlockSize(link->dm, &bs));
9709566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize((*is)[cnt], bs));
971520db06cSJed Brown   }
9723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97387c85e80SJed Brown }
97487c85e80SJed Brown 
97547c6ae99SBarry Smith /*@C
976dce8aebaSBarry Smith   DMCompositeGetGlobalISs - Gets the index sets for each composed object in a `DMCOMPOSITE`
97747c6ae99SBarry Smith 
97820f4b53cSBarry Smith   Collective
97947c6ae99SBarry Smith 
98047c6ae99SBarry Smith   Input Parameter:
981dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
98247c6ae99SBarry Smith 
9832fe279fdSBarry Smith   Output Parameter:
98447c6ae99SBarry Smith . is - the array of index sets
98547c6ae99SBarry Smith 
98647c6ae99SBarry Smith   Level: advanced
98747c6ae99SBarry Smith 
98847c6ae99SBarry Smith   Notes:
989f13dfd9eSBarry Smith   The `is` entries should be destroyed with `ISDestroy()`, `is` should be freed with `PetscFree()`
99047c6ae99SBarry Smith 
99147c6ae99SBarry Smith   These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
99247c6ae99SBarry Smith 
993dce8aebaSBarry Smith   Use `DMCompositeGetLocalISs()` for index sets in the packed local numbering, and
994dce8aebaSBarry Smith   `DMCompositeGetISLocalToGlobalMappings()` for to map local sub-`DM` (including ghost) indices to packed global
9956eb61c8cSJed Brown   indices.
99647c6ae99SBarry Smith 
99760225df5SJacob Faibussowitsch   Fortran Notes:
998dce8aebaSBarry Smith   The output argument 'is' must be an allocated array of sufficient length, which can be learned using `DMCompositeGetNumberDM()`.
999f3cb0f7eSJed Brown 
1000dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1001db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
1002c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
100347c6ae99SBarry Smith @*/
1004d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[])
1005d71ae5a4SJacob Faibussowitsch {
100666bb578eSMark Adams   PetscInt                cnt = 0;
100747c6ae99SBarry Smith   struct DMCompositeLink *next;
100847c6ae99SBarry Smith   PetscMPIInt             rank;
100947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
101071b14b3eSStefano Zampini   PetscBool               flg;
101147c6ae99SBarry Smith 
101247c6ae99SBarry Smith   PetscFunctionBegin;
101347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
10157a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
10167a8be351SBarry Smith   PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before");
10179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM, is));
101847c6ae99SBarry Smith   next = com->next;
10199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
102047c6ae99SBarry Smith 
102115229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
102247c6ae99SBarry Smith   while (next) {
1023e5e52638SMatthew G. Knepley     PetscDS prob;
1024e5e52638SMatthew G. Knepley 
10259566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt]));
10269566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &prob));
1027e5e52638SMatthew G. Knepley     if (prob) {
102865c226d8SMatthew G. Knepley       MatNullSpace space;
102965c226d8SMatthew G. Knepley       Mat          pmat;
10300f21e855SMatthew G. Knepley       PetscObject  disc;
10310f21e855SMatthew G. Knepley       PetscInt     Nf;
103265c226d8SMatthew G. Knepley 
10339566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(prob, &Nf));
1034f24dd8d2SMatthew G. Knepley       if (cnt < Nf) {
10359566063dSJacob Faibussowitsch         PetscCall(PetscDSGetDiscretization(prob, cnt, &disc));
10369566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space));
10379566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space));
10389566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space));
10399566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space));
10409566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat));
10419566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat));
104265c226d8SMatthew G. Knepley       }
1043f24dd8d2SMatthew G. Knepley     }
104447c6ae99SBarry Smith     cnt++;
104547c6ae99SBarry Smith     next = next->next;
104647c6ae99SBarry Smith   }
10473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104847c6ae99SBarry Smith }
104947c6ae99SBarry Smith 
105066976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
1051d71ae5a4SJacob Faibussowitsch {
10524d343eeaSMatthew G Knepley   PetscInt nDM;
10534d343eeaSMatthew G Knepley   DM      *dms;
10544d343eeaSMatthew G Knepley   PetscInt i;
10554d343eeaSMatthew G Knepley 
10564d343eeaSMatthew G Knepley   PetscFunctionBegin;
10579566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm, &nDM));
10588865f1eaSKarl Rupp   if (numFields) *numFields = nDM;
10599566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetGlobalISs(dm, fields));
10604d343eeaSMatthew G Knepley   if (fieldNames) {
10619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, &dms));
10629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, fieldNames));
10639566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, dms));
10644d343eeaSMatthew G Knepley     for (i = 0; i < nDM; i++) {
10654d343eeaSMatthew G Knepley       char        buf[256];
10664d343eeaSMatthew G Knepley       const char *splitname;
10674d343eeaSMatthew G Knepley 
10684d343eeaSMatthew G Knepley       /* Split naming precedence: object name, prefix, number */
10694d343eeaSMatthew G Knepley       splitname = ((PetscObject)dm)->name;
10704d343eeaSMatthew G Knepley       if (!splitname) {
10719566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname));
10724d343eeaSMatthew G Knepley         if (splitname) {
10734d343eeaSMatthew G Knepley           size_t len;
10749566063dSJacob Faibussowitsch           PetscCall(PetscStrncpy(buf, splitname, sizeof(buf)));
10758caf3d72SBarry Smith           buf[sizeof(buf) - 1] = 0;
10769566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(buf, &len));
10774d343eeaSMatthew G Knepley           if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */
10784d343eeaSMatthew G Knepley           splitname = buf;
10794d343eeaSMatthew G Knepley         }
10804d343eeaSMatthew G Knepley       }
10814d343eeaSMatthew G Knepley       if (!splitname) {
108263a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i));
10834d343eeaSMatthew G Knepley         splitname = buf;
10844d343eeaSMatthew G Knepley       }
10859566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i]));
10864d343eeaSMatthew G Knepley     }
10879566063dSJacob Faibussowitsch     PetscCall(PetscFree(dms));
10884d343eeaSMatthew G Knepley   }
10893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10904d343eeaSMatthew G Knepley }
10914d343eeaSMatthew G Knepley 
1092e7c4fc90SDmitry Karpeev /*
1093e7c4fc90SDmitry Karpeev  This could take over from DMCreateFieldIS(), as it is more general,
10940298fd71SBarry Smith  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1095e7c4fc90SDmitry Karpeev  At this point it's probably best to be less intrusive, however.
1096e7c4fc90SDmitry Karpeev  */
109766976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
1098d71ae5a4SJacob Faibussowitsch {
1099e7c4fc90SDmitry Karpeev   PetscInt nDM;
1100e7c4fc90SDmitry Karpeev   PetscInt i;
1101e7c4fc90SDmitry Karpeev 
1102e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
11039566063dSJacob Faibussowitsch   PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist));
1104e7c4fc90SDmitry Karpeev   if (dmlist) {
11059566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetNumberDM(dm, &nDM));
11069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, dmlist));
11079566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, *dmlist));
110848a46eb9SPierre Jolivet     for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i])));
1109e7c4fc90SDmitry Karpeev   }
11103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1111e7c4fc90SDmitry Karpeev }
1112e7c4fc90SDmitry Karpeev 
111347c6ae99SBarry Smith /*@C
1114dce8aebaSBarry Smith   DMCompositeGetLocalVectors - Gets local vectors for each part of a `DMCOMPOSITE`
1115dce8aebaSBarry Smith   Use `DMCompositeRestoreLocalVectors()` to return them.
111647c6ae99SBarry Smith 
111720f4b53cSBarry Smith   Not Collective; No Fortran Support
111847c6ae99SBarry Smith 
111947c6ae99SBarry Smith   Input Parameter:
1120dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
112147c6ae99SBarry Smith 
112247c6ae99SBarry Smith   Output Parameter:
1123a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`s
112447c6ae99SBarry Smith 
112547c6ae99SBarry Smith   Level: advanced
112647c6ae99SBarry Smith 
1127dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1128db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1129db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
113047c6ae99SBarry Smith @*/
1131d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...)
1132d71ae5a4SJacob Faibussowitsch {
113347c6ae99SBarry Smith   va_list                 Argp;
113447c6ae99SBarry Smith   struct DMCompositeLink *next;
113547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
113671b14b3eSStefano Zampini   PetscBool               flg;
113747c6ae99SBarry Smith 
113847c6ae99SBarry Smith   PetscFunctionBegin;
113947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11409566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11417a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
114247c6ae99SBarry Smith   next = com->next;
114315229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
114447c6ae99SBarry Smith   va_start(Argp, dm);
114547c6ae99SBarry Smith   while (next) {
114647c6ae99SBarry Smith     Vec *vec;
114747c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
11489566063dSJacob Faibussowitsch     if (vec) PetscCall(DMGetLocalVector(next->dm, vec));
114947c6ae99SBarry Smith     next = next->next;
115047c6ae99SBarry Smith   }
115147c6ae99SBarry Smith   va_end(Argp);
11523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115347c6ae99SBarry Smith }
115447c6ae99SBarry Smith 
115547c6ae99SBarry Smith /*@C
1156dce8aebaSBarry Smith   DMCompositeRestoreLocalVectors - Restores local vectors for each part of a `DMCOMPOSITE`
115747c6ae99SBarry Smith 
115820f4b53cSBarry Smith   Not Collective; No Fortran Support
115947c6ae99SBarry Smith 
116047c6ae99SBarry Smith   Input Parameter:
1161dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
116247c6ae99SBarry Smith 
116347c6ae99SBarry Smith   Output Parameter:
1164a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`
116547c6ae99SBarry Smith 
116647c6ae99SBarry Smith   Level: advanced
116747c6ae99SBarry Smith 
1168dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1169db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1170db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
117147c6ae99SBarry Smith @*/
1172d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...)
1173d71ae5a4SJacob Faibussowitsch {
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;
118415229ffcSPierre Jolivet   /* loop over packed objects, handling one at a 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(DMRestoreLocalVector(next->dm, vec));
119047c6ae99SBarry Smith     next = next->next;
119147c6ae99SBarry Smith   }
119247c6ae99SBarry Smith   va_end(Argp);
11933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119447c6ae99SBarry Smith }
119547c6ae99SBarry Smith 
119647c6ae99SBarry Smith /*@C
1197dce8aebaSBarry Smith   DMCompositeGetEntries - Gets the `DM` for each entry in a `DMCOMPOSITE`.
119847c6ae99SBarry Smith 
119947c6ae99SBarry Smith   Not Collective
120047c6ae99SBarry Smith 
120147c6ae99SBarry Smith   Input Parameter:
1202dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
120347c6ae99SBarry Smith 
120447c6ae99SBarry Smith   Output Parameter:
1205a4e35b19SJacob Faibussowitsch . ... - the individual entries `DM`
120647c6ae99SBarry Smith 
120747c6ae99SBarry Smith   Level: advanced
120847c6ae99SBarry Smith 
120960225df5SJacob Faibussowitsch   Fortran Notes:
12105d83a8b1SBarry Smith   Use `DMCompositeGetEntriesArray()`
1211f5f57ec0SBarry Smith 
1212dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()`
1213db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
121460225df5SJacob Faibussowitsch          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`
121547c6ae99SBarry Smith @*/
1216d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntries(DM dm, ...)
1217d71ae5a4SJacob Faibussowitsch {
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;
122815229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
122947c6ae99SBarry Smith   va_start(Argp, dm);
123047c6ae99SBarry Smith   while (next) {
123147c6ae99SBarry Smith     DM *dmn;
123247c6ae99SBarry Smith     dmn = va_arg(Argp, DM *);
12339ae5db72SJed Brown     if (dmn) *dmn = next->dm;
123447c6ae99SBarry Smith     next = next->next;
123547c6ae99SBarry Smith   }
123647c6ae99SBarry Smith   va_end(Argp);
12373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
123847c6ae99SBarry Smith }
123947c6ae99SBarry Smith 
1240cc4c1da9SBarry Smith /*@
1241dce8aebaSBarry Smith   DMCompositeGetEntriesArray - Gets the DM for each entry in a `DMCOMPOSITE`
12422fa5ba8aSJed Brown 
12432fa5ba8aSJed Brown   Not Collective
12442fa5ba8aSJed Brown 
12452fa5ba8aSJed Brown   Input Parameter:
1246dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
1247907376e6SBarry Smith 
1248907376e6SBarry Smith   Output Parameter:
1249dce8aebaSBarry Smith . dms - array of sufficient length (see `DMCompositeGetNumberDM()`) to hold the individual `DM`
12502fa5ba8aSJed Brown 
12512fa5ba8aSJed Brown   Level: advanced
12522fa5ba8aSJed Brown 
1253dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()`
1254db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
125560225df5SJacob Faibussowitsch          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`
12562fa5ba8aSJed Brown @*/
1257d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntriesArray(DM dm, DM dms[])
1258d71ae5a4SJacob Faibussowitsch {
12592fa5ba8aSJed Brown   struct DMCompositeLink *next;
12602fa5ba8aSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
12612fa5ba8aSJed Brown   PetscInt                i;
126271b14b3eSStefano Zampini   PetscBool               flg;
12632fa5ba8aSJed Brown 
12642fa5ba8aSJed Brown   PetscFunctionBegin;
12652fa5ba8aSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12669566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
12677a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
126815229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
12692fa5ba8aSJed Brown   for (next = com->next, i = 0; next; next = next->next, i++) dms[i] = next->dm;
12703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12712fa5ba8aSJed Brown }
12722fa5ba8aSJed Brown 
1273e10fd815SStefano Zampini typedef struct {
1274e10fd815SStefano Zampini   DM           dm;
1275e10fd815SStefano Zampini   PetscViewer *subv;
1276e10fd815SStefano Zampini   Vec         *vecs;
1277e10fd815SStefano Zampini } GLVisViewerCtx;
1278e10fd815SStefano Zampini 
1279d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx)
1280d71ae5a4SJacob Faibussowitsch {
1281e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
1282e10fd815SStefano Zampini   PetscInt        i, n;
1283e10fd815SStefano Zampini 
1284e10fd815SStefano Zampini   PetscFunctionBegin;
12859566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(ctx->dm, &n));
128648a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscViewerDestroy(&ctx->subv[i]));
12879566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ctx->subv, ctx->vecs));
12889566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dm));
12899566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
12903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1291e10fd815SStefano Zampini }
1292e10fd815SStefano Zampini 
1293d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
1294d71ae5a4SJacob Faibussowitsch {
1295e10fd815SStefano Zampini   Vec             X   = (Vec)oX;
1296e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
1297e10fd815SStefano Zampini   PetscInt        i, n, cumf;
1298e10fd815SStefano Zampini 
1299e10fd815SStefano Zampini   PetscFunctionBegin;
13009566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(ctx->dm, &n));
13019566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccessArray(ctx->dm, X, n, NULL, ctx->vecs));
1302e10fd815SStefano Zampini   for (i = 0, cumf = 0; i < n; i++) {
1303e10fd815SStefano Zampini     PetscErrorCode (*g2l)(PetscObject, PetscInt, PetscObject[], void *);
1304e10fd815SStefano Zampini     void    *fctx;
1305e10fd815SStefano Zampini     PetscInt nfi;
1306e10fd815SStefano Zampini 
130734e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nfi, NULL, NULL, &g2l, NULL, &fctx));
1308e10fd815SStefano Zampini     if (!nfi) continue;
13091baa6e33SBarry Smith     if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i], nfi, oXfield + cumf, fctx));
1310f4f49eeaSPierre Jolivet     else PetscCall(VecCopy(ctx->vecs[i], (Vec)oXfield[cumf]));
1311e10fd815SStefano Zampini     cumf += nfi;
1312e10fd815SStefano Zampini   }
13139566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccessArray(ctx->dm, X, n, NULL, ctx->vecs));
13143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1315e10fd815SStefano Zampini }
1316e10fd815SStefano Zampini 
1317d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer)
1318d71ae5a4SJacob Faibussowitsch {
1319e10fd815SStefano Zampini   DM              dm = (DM)odm, *dms;
1320e10fd815SStefano Zampini   Vec            *Ufds;
1321e10fd815SStefano Zampini   GLVisViewerCtx *ctx;
1322e10fd815SStefano Zampini   PetscInt        i, n, tnf, *sdim;
1323e10fd815SStefano Zampini   char          **fecs;
1324e10fd815SStefano Zampini 
1325e10fd815SStefano Zampini   PetscFunctionBegin;
13269566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
13279566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
1328e10fd815SStefano Zampini   ctx->dm = dm;
13299566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm, &n));
13309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &dms));
13319566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetEntriesArray(dm, dms));
13329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n, &ctx->subv, n, &ctx->vecs));
1333e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1334e10fd815SStefano Zampini     PetscInt nf;
1335e10fd815SStefano Zampini 
13369566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PetscObjectComm(odm), &ctx->subv[i]));
13379566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(ctx->subv[i], PETSCVIEWERGLVIS));
133834e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisSetDM_Internal(ctx->subv[i], (PetscObject)dms[i]));
133934e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, NULL, NULL, NULL, NULL, NULL));
1340e10fd815SStefano Zampini     tnf += nf;
1341e10fd815SStefano Zampini   }
13429566063dSJacob Faibussowitsch   PetscCall(PetscFree(dms));
13439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(tnf, &fecs, tnf, &sdim, tnf, &Ufds));
1344e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1345e10fd815SStefano Zampini     PetscInt    *sd, nf, f;
1346e10fd815SStefano Zampini     const char **fec;
1347e10fd815SStefano Zampini     Vec         *Uf;
1348e10fd815SStefano Zampini 
134934e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, &fec, &sd, NULL, (PetscObject **)&Uf, NULL));
1350e10fd815SStefano Zampini     for (f = 0; f < nf; f++) {
13519566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(fec[f], &fecs[tnf + f]));
1352e10fd815SStefano Zampini       Ufds[tnf + f] = Uf[f];
1353e10fd815SStefano Zampini       sdim[tnf + f] = sd[f];
1354e10fd815SStefano Zampini     }
1355e10fd815SStefano Zampini     tnf += nf;
1356e10fd815SStefano Zampini   }
13579566063dSJacob Faibussowitsch   PetscCall(PetscViewerGLVisSetFields(viewer, tnf, (const char **)fecs, sdim, DMCompositeSampleGLVisFields_Private, (PetscObject *)Ufds, ctx, DestroyGLVisViewerCtx_Private));
135848a46eb9SPierre Jolivet   for (i = 0; i < tnf; i++) PetscCall(PetscFree(fecs[i]));
13599566063dSJacob Faibussowitsch   PetscCall(PetscFree3(fecs, sdim, Ufds));
13603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1361e10fd815SStefano Zampini }
1362e10fd815SStefano Zampini 
136366976f2fSJacob Faibussowitsch static PetscErrorCode DMRefine_Composite(DM dmi, MPI_Comm comm, DM *fine)
1364d71ae5a4SJacob Faibussowitsch {
136547c6ae99SBarry Smith   struct DMCompositeLink *next;
136647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dmi->data;
136747c6ae99SBarry Smith   DM                      dm;
136847c6ae99SBarry Smith 
136947c6ae99SBarry Smith   PetscFunctionBegin;
137047c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi, DM_CLASSID, 1);
137148a46eb9SPierre Jolivet   if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm));
13729566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmi));
137347c6ae99SBarry Smith   next = com->next;
13749566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm, fine));
137547c6ae99SBarry Smith 
137615229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
137747c6ae99SBarry Smith   while (next) {
13789566063dSJacob Faibussowitsch     PetscCall(DMRefine(next->dm, comm, &dm));
13799566063dSJacob Faibussowitsch     PetscCall(DMCompositeAddDM(*fine, dm));
13809566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dm));
138147c6ae99SBarry Smith     next = next->next;
138247c6ae99SBarry Smith   }
13833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138447c6ae99SBarry Smith }
138547c6ae99SBarry Smith 
138666976f2fSJacob Faibussowitsch static PetscErrorCode DMCoarsen_Composite(DM dmi, MPI_Comm comm, DM *fine)
1387d71ae5a4SJacob Faibussowitsch {
138814354c39SJed Brown   struct DMCompositeLink *next;
138914354c39SJed Brown   DM_Composite           *com = (DM_Composite *)dmi->data;
139014354c39SJed Brown   DM                      dm;
139114354c39SJed Brown 
139214354c39SJed Brown   PetscFunctionBegin;
139314354c39SJed Brown   PetscValidHeaderSpecific(dmi, DM_CLASSID, 1);
13949566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmi));
139548a46eb9SPierre Jolivet   if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm));
139614354c39SJed Brown   next = com->next;
13979566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm, fine));
139814354c39SJed Brown 
139915229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
140014354c39SJed Brown   while (next) {
14019566063dSJacob Faibussowitsch     PetscCall(DMCoarsen(next->dm, comm, &dm));
14029566063dSJacob Faibussowitsch     PetscCall(DMCompositeAddDM(*fine, dm));
14039566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dm));
140414354c39SJed Brown     next = next->next;
140514354c39SJed Brown   }
14063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
140714354c39SJed Brown }
140847c6ae99SBarry Smith 
140966976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_Composite(DM coarse, DM fine, Mat *A, Vec *v)
1410d71ae5a4SJacob Faibussowitsch {
14119ae5db72SJed Brown   PetscInt                m, n, M, N, nDM, i;
141247c6ae99SBarry Smith   struct DMCompositeLink *nextc;
141347c6ae99SBarry Smith   struct DMCompositeLink *nextf;
141425296bd5SBarry Smith   Vec                     gcoarse, gfine, *vecs;
141547c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite *)coarse->data;
141647c6ae99SBarry Smith   DM_Composite           *comfine   = (DM_Composite *)fine->data;
14179ae5db72SJed Brown   Mat                    *mats;
141847c6ae99SBarry Smith 
141947c6ae99SBarry Smith   PetscFunctionBegin;
142047c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse, DM_CLASSID, 1);
142147c6ae99SBarry Smith   PetscValidHeaderSpecific(fine, DM_CLASSID, 2);
14229566063dSJacob Faibussowitsch   PetscCall(DMSetUp(coarse));
14239566063dSJacob Faibussowitsch   PetscCall(DMSetUp(fine));
142447c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14259566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(coarse, &gcoarse));
14269566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(fine, &gfine));
14279566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gcoarse, &n));
14289566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gfine, &m));
14299566063dSJacob Faibussowitsch   PetscCall(VecGetSize(gcoarse, &N));
14309566063dSJacob Faibussowitsch   PetscCall(VecGetSize(gfine, &M));
14319566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(coarse, &gcoarse));
14329566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(fine, &gfine));
143347c6ae99SBarry Smith 
14349ae5db72SJed Brown   nDM = comfine->nDM;
143563a3b9bcSJacob 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);
14369566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nDM * nDM, &mats));
143748a46eb9SPierre Jolivet   if (v) PetscCall(PetscCalloc1(nDM, &vecs));
143847c6ae99SBarry Smith 
143915229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
14409ae5db72SJed Brown   for (nextc = comcoarse->next, nextf = comfine->next, i = 0; nextc; nextc = nextc->next, nextf = nextf->next, i++) {
14411baa6e33SBarry Smith     if (!v) PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], NULL));
14421baa6e33SBarry Smith     else PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], &vecs[i]));
144347c6ae99SBarry Smith   }
14449566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, nDM, NULL, mats, A));
14451baa6e33SBarry Smith   if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, vecs, v));
14469566063dSJacob Faibussowitsch   for (i = 0; i < nDM * nDM; i++) PetscCall(MatDestroy(&mats[i]));
14479566063dSJacob Faibussowitsch   PetscCall(PetscFree(mats));
144825296bd5SBarry Smith   if (v) {
14499566063dSJacob Faibussowitsch     for (i = 0; i < nDM; i++) PetscCall(VecDestroy(&vecs[i]));
14509566063dSJacob Faibussowitsch     PetscCall(PetscFree(vecs));
145125296bd5SBarry Smith   }
14523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
145347c6ae99SBarry Smith }
145447c6ae99SBarry Smith 
1455d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
1456d71ae5a4SJacob Faibussowitsch {
14571411c6eeSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
14581411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1459f7efa3c7SJed Brown   PetscInt                i;
14601411c6eeSJed Brown 
14611411c6eeSJed Brown   PetscFunctionBegin;
14621411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
14639566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetISLocalToGlobalMappings(dm, &ltogs));
14649566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm), com->nDM, ltogs, &dm->ltogmap));
14659566063dSJacob Faibussowitsch   for (i = 0; i < com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(&ltogs[i]));
14669566063dSJacob Faibussowitsch   PetscCall(PetscFree(ltogs));
14673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14681411c6eeSJed Brown }
14691411c6eeSJed Brown 
147066976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateColoring_Composite(DM dm, ISColoringType ctype, ISColoring *coloring)
1471d71ae5a4SJacob Faibussowitsch {
147247c6ae99SBarry Smith   PetscInt         n, i, cnt;
147347c6ae99SBarry Smith   ISColoringValue *colors;
147447c6ae99SBarry Smith   PetscBool        dense  = PETSC_FALSE;
147547c6ae99SBarry Smith   ISColoringValue  maxcol = 0;
147647c6ae99SBarry Smith   DM_Composite    *com    = (DM_Composite *)dm->data;
147747c6ae99SBarry Smith 
147847c6ae99SBarry Smith   PetscFunctionBegin;
147947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
148008401ef6SPierre Jolivet   PetscCheck(ctype != IS_COLORING_LOCAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only global coloring supported");
1481f7d195e4SLawrence Mitchell   if (ctype == IS_COLORING_GLOBAL) {
148247c6ae99SBarry Smith     n = com->n;
1483ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unknown ISColoringType");
14849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &colors)); /* freed in ISColoringDestroy() */
148547c6ae99SBarry Smith 
14869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL));
148747c6ae99SBarry Smith   if (dense) {
1488*6497c311SBarry Smith     PetscCall(ISColoringValueCast(com->N, &maxcol));
1489*6497c311SBarry Smith     for (i = 0; i < n; i++) PetscCall(ISColoringValueCast(com->rstart + i, colors + i));
149047c6ae99SBarry Smith   } else {
149147c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
149247c6ae99SBarry Smith     PetscMPIInt             rank;
149347c6ae99SBarry Smith 
14949566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
149547c6ae99SBarry Smith     cnt = 0;
149647c6ae99SBarry Smith     while (next) {
149747c6ae99SBarry Smith       ISColoring lcoloring;
149847c6ae99SBarry Smith 
14999566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(next->dm, IS_COLORING_GLOBAL, &lcoloring));
1500*6497c311SBarry Smith       for (i = 0; i < lcoloring->N; i++) PetscCall(ISColoringValueCast(maxcol + lcoloring->colors[i], colors + cnt++));
150147c6ae99SBarry Smith       maxcol += lcoloring->n;
15029566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&lcoloring));
150347c6ae99SBarry Smith       next = next->next;
150447c6ae99SBarry Smith     }
150547c6ae99SBarry Smith   }
15069566063dSJacob Faibussowitsch   PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), maxcol, n, colors, PETSC_OWN_POINTER, coloring));
15073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150847c6ae99SBarry Smith }
150947c6ae99SBarry Smith 
151066976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec)
1511d71ae5a4SJacob Faibussowitsch {
151247c6ae99SBarry Smith   struct DMCompositeLink *next;
151347c6ae99SBarry Smith   PetscScalar            *garray, *larray;
151447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
151547c6ae99SBarry Smith 
151647c6ae99SBarry Smith   PetscFunctionBegin;
151747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
151847c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
151939d35262SVincent Le Chenadec 
152048a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
152139d35262SVincent Le Chenadec 
15229566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gvec, &garray));
15239566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lvec, &larray));
152447c6ae99SBarry Smith 
152515229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
152639d35262SVincent Le Chenadec   next = com->next;
152747c6ae99SBarry Smith   while (next) {
152847c6ae99SBarry Smith     Vec      local, global;
152947c6ae99SBarry Smith     PetscInt N;
153047c6ae99SBarry Smith 
15319566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
15329566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(global, &N));
15339566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(global, garray));
15349566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local));
15359566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local, larray));
15369566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(next->dm, global, mode, local));
15379566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(next->dm, global, mode, local));
15389566063dSJacob Faibussowitsch     PetscCall(VecResetArray(global));
15399566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local));
15409566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
15419566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local));
154239d35262SVincent Le Chenadec 
154306ebdd98SJed Brown     larray += next->nlocal;
154439d35262SVincent Le Chenadec     garray += next->n;
154547c6ae99SBarry Smith     next = next->next;
154647c6ae99SBarry Smith   }
154747c6ae99SBarry Smith 
15489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gvec, NULL));
15499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lvec, NULL));
15503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
155147c6ae99SBarry Smith }
155247c6ae99SBarry Smith 
155366976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec)
1554d71ae5a4SJacob Faibussowitsch {
15550c010503SBarry Smith   PetscFunctionBegin;
155639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
155739d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
155839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 4);
15593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
156039d35262SVincent Le Chenadec }
156139d35262SVincent Le Chenadec 
156266976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1563d71ae5a4SJacob Faibussowitsch {
156439d35262SVincent Le Chenadec   struct DMCompositeLink *next;
156539d35262SVincent Le Chenadec   PetscScalar            *larray, *garray;
156639d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite *)dm->data;
156739d35262SVincent Le Chenadec 
156839d35262SVincent Le Chenadec   PetscFunctionBegin;
156939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
157039d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
157139d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
157239d35262SVincent Le Chenadec 
157348a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
157439d35262SVincent Le Chenadec 
15759566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lvec, &larray));
15769566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gvec, &garray));
157739d35262SVincent Le Chenadec 
157815229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
157939d35262SVincent Le Chenadec   next = com->next;
158039d35262SVincent Le Chenadec   while (next) {
158139d35262SVincent Le Chenadec     Vec global, local;
158239d35262SVincent Le Chenadec 
15839566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local));
15849566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local, larray));
15859566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
15869566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(global, garray));
15879566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(next->dm, local, mode, global));
15889566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(next->dm, local, mode, global));
15899566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local));
15909566063dSJacob Faibussowitsch     PetscCall(VecResetArray(global));
15919566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
15929566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local));
159339d35262SVincent Le Chenadec 
159439d35262SVincent Le Chenadec     garray += next->n;
159539d35262SVincent Le Chenadec     larray += next->nlocal;
159639d35262SVincent Le Chenadec     next = next->next;
159739d35262SVincent Le Chenadec   }
159839d35262SVincent Le Chenadec 
15999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gvec, NULL));
16009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lvec, NULL));
16013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
160239d35262SVincent Le Chenadec }
160339d35262SVincent Le Chenadec 
160466976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1605d71ae5a4SJacob Faibussowitsch {
160639d35262SVincent Le Chenadec   PetscFunctionBegin;
160739d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
160839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
160939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
161139d35262SVincent Le Chenadec }
161239d35262SVincent Le Chenadec 
161366976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2)
1614d71ae5a4SJacob Faibussowitsch {
161539d35262SVincent Le Chenadec   struct DMCompositeLink *next;
161639d35262SVincent Le Chenadec   PetscScalar            *array1, *array2;
161739d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite *)dm->data;
161839d35262SVincent Le Chenadec 
161939d35262SVincent Le Chenadec   PetscFunctionBegin;
162039d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
162139d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2);
162239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4);
162339d35262SVincent Le Chenadec 
162448a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
162539d35262SVincent Le Chenadec 
16269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec1, &array1));
16279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec2, &array2));
162839d35262SVincent Le Chenadec 
162915229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
163039d35262SVincent Le Chenadec   next = com->next;
163139d35262SVincent Le Chenadec   while (next) {
163239d35262SVincent Le Chenadec     Vec local1, local2;
163339d35262SVincent Le Chenadec 
16349566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local1));
16359566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local1, array1));
16369566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local2));
16379566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local2, array2));
16389566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2));
16399566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2));
16409566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local2));
16419566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local2));
16429566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local1));
16439566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local1));
164439d35262SVincent Le Chenadec 
164539d35262SVincent Le Chenadec     array1 += next->nlocal;
164639d35262SVincent Le Chenadec     array2 += next->nlocal;
164739d35262SVincent Le Chenadec     next = next->next;
164839d35262SVincent Le Chenadec   }
164939d35262SVincent Le Chenadec 
16509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec1, NULL));
16519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec2, NULL));
16523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165339d35262SVincent Le Chenadec }
165439d35262SVincent Le Chenadec 
165566976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1656d71ae5a4SJacob Faibussowitsch {
165739d35262SVincent Le Chenadec   PetscFunctionBegin;
165839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
165939d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
166039d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16620c010503SBarry Smith }
166347c6ae99SBarry Smith 
16646ae3a549SBarry Smith /*MC
1665dce8aebaSBarry Smith    DMCOMPOSITE = "composite" - A `DM` object that is used to manage data for a collection of `DM`
16666ae3a549SBarry Smith 
16676ae3a549SBarry Smith   Level: intermediate
16686ae3a549SBarry Smith 
1669db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()`
16706ae3a549SBarry Smith M*/
16716ae3a549SBarry Smith 
1672d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1673d71ae5a4SJacob Faibussowitsch {
1674a4121054SBarry Smith   DM_Composite *com;
1675a4121054SBarry Smith 
1676a4121054SBarry Smith   PetscFunctionBegin;
16774dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&com));
1678a4121054SBarry Smith   p->data     = com;
1679a4121054SBarry Smith   com->n      = 0;
16807ac2b803SAlex Fikl   com->nghost = 0;
16810298fd71SBarry Smith   com->next   = NULL;
1682a4121054SBarry Smith   com->nDM    = 0;
1683a4121054SBarry Smith 
1684a4121054SBarry Smith   p->ops->createglobalvector       = DMCreateGlobalVector_Composite;
1685a4121054SBarry Smith   p->ops->createlocalvector        = DMCreateLocalVector_Composite;
1686184d77edSJed Brown   p->ops->getlocaltoglobalmapping  = DMGetLocalToGlobalMapping_Composite;
16874d343eeaSMatthew G Knepley   p->ops->createfieldis            = DMCreateFieldIS_Composite;
168816621825SDmitry Karpeev   p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite;
1689a4121054SBarry Smith   p->ops->refine                   = DMRefine_Composite;
169014354c39SJed Brown   p->ops->coarsen                  = DMCoarsen_Composite;
169125296bd5SBarry Smith   p->ops->createinterpolation      = DMCreateInterpolation_Composite;
169225296bd5SBarry Smith   p->ops->creatematrix             = DMCreateMatrix_Composite;
1693e727c939SJed Brown   p->ops->getcoloring              = DMCreateColoring_Composite;
1694a4121054SBarry Smith   p->ops->globaltolocalbegin       = DMGlobalToLocalBegin_Composite;
1695a4121054SBarry Smith   p->ops->globaltolocalend         = DMGlobalToLocalEnd_Composite;
169639d35262SVincent Le Chenadec   p->ops->localtoglobalbegin       = DMLocalToGlobalBegin_Composite;
169739d35262SVincent Le Chenadec   p->ops->localtoglobalend         = DMLocalToGlobalEnd_Composite;
169839d35262SVincent Le Chenadec   p->ops->localtolocalbegin        = DMLocalToLocalBegin_Composite;
169939d35262SVincent Le Chenadec   p->ops->localtolocalend          = DMLocalToLocalEnd_Composite;
1700a4121054SBarry Smith   p->ops->destroy                  = DMDestroy_Composite;
1701a4121054SBarry Smith   p->ops->view                     = DMView_Composite;
1702a4121054SBarry Smith   p->ops->setup                    = DMSetUp_Composite;
1703e10fd815SStefano Zampini 
17049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite));
17053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1706a4121054SBarry Smith }
1707a4121054SBarry Smith 
17081fd49c25SBarry Smith /*@
1709dce8aebaSBarry Smith   DMCompositeCreate - Creates a `DMCOMPOSITE`, used to generate "composite"
17100c010503SBarry Smith   vectors made up of several subvectors.
17110c010503SBarry Smith 
1712d083f849SBarry Smith   Collective
171347c6ae99SBarry Smith 
171447c6ae99SBarry Smith   Input Parameter:
17150c010503SBarry Smith . comm - the processors that will share the global vector
17160c010503SBarry Smith 
17172fe279fdSBarry Smith   Output Parameter:
1718dce8aebaSBarry Smith . packer - the `DMCOMPOSITE` object
171947c6ae99SBarry Smith 
172047c6ae99SBarry Smith   Level: advanced
172147c6ae99SBarry Smith 
172260225df5SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCreate()`
1723db781477SPatrick Sanan           `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`
1724db781477SPatrick Sanan           `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
172547c6ae99SBarry Smith @*/
1726d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer)
1727d71ae5a4SJacob Faibussowitsch {
172847c6ae99SBarry Smith   PetscFunctionBegin;
17294f572ea9SToby Isaac   PetscAssertPointer(packer, 2);
17309566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, packer));
17319566063dSJacob Faibussowitsch   PetscCall(DMSetType(*packer, DMCOMPOSITE));
17323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173347c6ae99SBarry Smith }
1734