xref: /petsc/src/dm/impls/composite/pack.c (revision 966bd95a39c2334d2e2ce17ad22128f3c1861eeb)
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 {
589f196a02SMartin Diehl   PetscBool     isascii;
5947c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
6047c6ae99SBarry Smith 
6147c6ae99SBarry Smith   PetscFunctionBegin;
629f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii));
639f196a02SMartin Diehl   if (isascii) {
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 
163dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetEntries()`, `DMCompositeScatter()`
16447c6ae99SBarry Smith @*/
165d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccess(DM dm, Vec gvec, ...)
166d71ae5a4SJacob Faibussowitsch {
16747c6ae99SBarry Smith   va_list                 Argp;
16847c6ae99SBarry Smith   struct DMCompositeLink *next;
16947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
1705edff71fSBarry Smith   PetscInt                readonly;
17171b14b3eSStefano Zampini   PetscBool               flg;
17247c6ae99SBarry Smith 
17347c6ae99SBarry Smith   PetscFunctionBegin;
17447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
17547c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
1769566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1777a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
17847c6ae99SBarry Smith   next = com->next;
17948a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
18047c6ae99SBarry Smith 
1819566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec, &readonly));
18215229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
18347c6ae99SBarry Smith   va_start(Argp, gvec);
18447c6ae99SBarry Smith   while (next) {
18547c6ae99SBarry Smith     Vec *vec;
18647c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
1879ae5db72SJed Brown     if (vec) {
1889566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, vec));
1895edff71fSBarry Smith       if (readonly) {
1905edff71fSBarry Smith         const PetscScalar *array;
1919566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(gvec, &array));
1929566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec, array + next->rstart));
1939566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(*vec));
1949566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(gvec, &array));
1955edff71fSBarry Smith       } else {
1965edff71fSBarry Smith         PetscScalar *array;
1979566063dSJacob Faibussowitsch         PetscCall(VecGetArray(gvec, &array));
1989566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec, array + next->rstart));
1999566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(gvec, &array));
20047c6ae99SBarry Smith       }
2015edff71fSBarry Smith     }
20247c6ae99SBarry Smith     next = next->next;
20347c6ae99SBarry Smith   }
20447c6ae99SBarry Smith   va_end(Argp);
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20647c6ae99SBarry Smith }
20747c6ae99SBarry Smith 
2085d83a8b1SBarry Smith /*@
209f73e5cebSJed Brown   DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
210f73e5cebSJed Brown   representation.
211f73e5cebSJed Brown 
21220f4b53cSBarry Smith   Collective
213f73e5cebSJed Brown 
214f73e5cebSJed Brown   Input Parameters:
215dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE`
216f73e5cebSJed Brown . pvec    - packed vector
217f73e5cebSJed Brown . nwanted - number of vectors wanted
218ce78bad3SBarry Smith - wanted  - sorted array of integers indicating thde vectors wanted, or `NULL` to get all vectors, length `nwanted`
219f73e5cebSJed Brown 
2202fe279fdSBarry Smith   Output Parameter:
221f13dfd9eSBarry Smith . vecs - array of requested global vectors (must be previously allocated and of length `nwanted`)
222f73e5cebSJed Brown 
223f73e5cebSJed Brown   Level: advanced
224f73e5cebSJed Brown 
225dce8aebaSBarry Smith   Note:
226dce8aebaSBarry Smith   Use `DMCompositeRestoreAccessArray()` to return the vectors when you no longer need them
227dce8aebaSBarry Smith 
228dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
229f73e5cebSJed Brown @*/
2305d83a8b1SBarry Smith PetscErrorCode DMCompositeGetAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[])
231d71ae5a4SJacob Faibussowitsch {
232f73e5cebSJed Brown   struct DMCompositeLink *link;
233f73e5cebSJed Brown   PetscInt                i, wnum;
234f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
235bee642f7SBarry Smith   PetscInt                readonly;
23671b14b3eSStefano Zampini   PetscBool               flg;
237f73e5cebSJed Brown 
238f73e5cebSJed Brown   PetscFunctionBegin;
239f73e5cebSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
240f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
2419566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
2427a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
24348a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
244f73e5cebSJed Brown 
2459566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
246f73e5cebSJed Brown   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
247f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
248f73e5cebSJed Brown       Vec v;
2499566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(link->dm, &v));
250bee642f7SBarry Smith       if (readonly) {
251bee642f7SBarry Smith         const PetscScalar *array;
2529566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec, &array));
2539566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + link->rstart));
2549566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
2559566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec, &array));
256bee642f7SBarry Smith       } else {
257bee642f7SBarry Smith         PetscScalar *array;
2589566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec, &array));
2599566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + link->rstart));
2609566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec, &array));
261bee642f7SBarry Smith       }
262f73e5cebSJed Brown       vecs[wnum++] = v;
263f73e5cebSJed Brown     }
264f73e5cebSJed Brown   }
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
266f73e5cebSJed Brown }
267f73e5cebSJed Brown 
2685d83a8b1SBarry Smith /*@
2697ac2b803SAlex Fikl   DMCompositeGetLocalAccessArray - Allows one to access the individual
2707ac2b803SAlex Fikl   packed vectors in their local representation.
2717ac2b803SAlex Fikl 
27220f4b53cSBarry Smith   Collective
2737ac2b803SAlex Fikl 
2747ac2b803SAlex Fikl   Input Parameters:
275dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE`
2767ac2b803SAlex Fikl . pvec    - packed vector
2777ac2b803SAlex Fikl . nwanted - number of vectors wanted
278f13dfd9eSBarry Smith - wanted  - sorted array of vectors wanted, or `NULL` to get all vectors, length `nwanted`
2797ac2b803SAlex Fikl 
2802fe279fdSBarry Smith   Output Parameter:
281f13dfd9eSBarry Smith . vecs - array of requested local vectors (must be allocated and of length `nwanted`)
2827ac2b803SAlex Fikl 
2837ac2b803SAlex Fikl   Level: advanced
2847ac2b803SAlex Fikl 
285dce8aebaSBarry Smith   Note:
286dce8aebaSBarry Smith   Use `DMCompositeRestoreLocalAccessArray()` to return the vectors
287dce8aebaSBarry Smith   when you no longer need them.
288dce8aebaSBarry Smith 
289dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`,
290db781477SPatrick Sanan           `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
2917ac2b803SAlex Fikl @*/
2925d83a8b1SBarry Smith PetscErrorCode DMCompositeGetLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[])
293d71ae5a4SJacob Faibussowitsch {
2947ac2b803SAlex Fikl   struct DMCompositeLink *link;
2957ac2b803SAlex Fikl   PetscInt                i, wnum;
2967ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite *)dm->data;
2977ac2b803SAlex Fikl   PetscInt                readonly;
2987ac2b803SAlex Fikl   PetscInt                nlocal = 0;
29971b14b3eSStefano Zampini   PetscBool               flg;
3007ac2b803SAlex Fikl 
3017ac2b803SAlex Fikl   PetscFunctionBegin;
3027ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3037ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
3049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3057a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
30648a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
3077ac2b803SAlex Fikl 
3089566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
3097ac2b803SAlex Fikl   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
3107ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
3117ac2b803SAlex Fikl       Vec v;
3129566063dSJacob Faibussowitsch       PetscCall(DMGetLocalVector(link->dm, &v));
3137ac2b803SAlex Fikl       if (readonly) {
3147ac2b803SAlex Fikl         const PetscScalar *array;
3159566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec, &array));
3169566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + nlocal));
317b1c3483dSMark 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
3189566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
3199566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec, &array));
3207ac2b803SAlex Fikl       } else {
3217ac2b803SAlex Fikl         PetscScalar *array;
3229566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec, &array));
3239566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + nlocal));
3249566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec, &array));
3257ac2b803SAlex Fikl       }
3267ac2b803SAlex Fikl       vecs[wnum++] = v;
3277ac2b803SAlex Fikl     }
3287ac2b803SAlex Fikl 
3297ac2b803SAlex Fikl     nlocal += link->nlocal;
3307ac2b803SAlex Fikl   }
3313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3327ac2b803SAlex Fikl }
3337ac2b803SAlex Fikl 
33447c6ae99SBarry Smith /*@C
335dce8aebaSBarry Smith   DMCompositeRestoreAccess - Returns the vectors obtained with `DMCompositeGetAccess()`
33647c6ae99SBarry Smith   representation.
33747c6ae99SBarry Smith 
33820f4b53cSBarry Smith   Collective
33947c6ae99SBarry Smith 
3409ae5db72SJed Brown   Input Parameters:
341dce8aebaSBarry Smith + dm   - the `DMCOMPOSITE` object
34247c6ae99SBarry Smith . gvec - the global vector
343a4e35b19SJacob Faibussowitsch - ...  - the individual parallel vectors, `NULL` for those that are not needed
34447c6ae99SBarry Smith 
34547c6ae99SBarry Smith   Level: advanced
34647c6ae99SBarry Smith 
347dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
348db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`,
34942747ad1SJacob Faibussowitsch          `DMCompositeGetAccess()`
35047c6ae99SBarry Smith @*/
351d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...)
352d71ae5a4SJacob Faibussowitsch {
35347c6ae99SBarry Smith   va_list                 Argp;
35447c6ae99SBarry Smith   struct DMCompositeLink *next;
35547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
3565edff71fSBarry Smith   PetscInt                readonly;
35771b14b3eSStefano Zampini   PetscBool               flg;
35847c6ae99SBarry Smith 
35947c6ae99SBarry Smith   PetscFunctionBegin;
36047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36147c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
3629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3637a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
36447c6ae99SBarry Smith   next = com->next;
36548a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
36647c6ae99SBarry Smith 
3679566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec, &readonly));
36815229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
36947c6ae99SBarry Smith   va_start(Argp, gvec);
37047c6ae99SBarry Smith   while (next) {
37147c6ae99SBarry Smith     Vec *vec;
37247c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
3739ae5db72SJed Brown     if (vec) {
3749566063dSJacob Faibussowitsch       PetscCall(VecResetArray(*vec));
3751baa6e33SBarry Smith       if (readonly) PetscCall(VecLockReadPop(*vec));
3769566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, vec));
37747c6ae99SBarry Smith     }
37847c6ae99SBarry Smith     next = next->next;
37947c6ae99SBarry Smith   }
38047c6ae99SBarry Smith   va_end(Argp);
3813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38247c6ae99SBarry Smith }
38347c6ae99SBarry Smith 
3845d83a8b1SBarry Smith /*@
385dce8aebaSBarry Smith   DMCompositeRestoreAccessArray - Returns the vectors obtained with `DMCompositeGetAccessArray()`
386f73e5cebSJed Brown 
38720f4b53cSBarry Smith   Collective
388f73e5cebSJed Brown 
389f73e5cebSJed Brown   Input Parameters:
390dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE` object
391f73e5cebSJed Brown . pvec    - packed vector
392f73e5cebSJed Brown . nwanted - number of vectors wanted
393f13dfd9eSBarry Smith . wanted  - sorted array of vectors wanted, or `NULL` to restore all vectors
394f13dfd9eSBarry Smith - vecs    - array of global vectors
395f73e5cebSJed Brown 
396f73e5cebSJed Brown   Level: advanced
397f73e5cebSJed Brown 
398dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
399f73e5cebSJed Brown @*/
4005d83a8b1SBarry Smith PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[])
401d71ae5a4SJacob Faibussowitsch {
402f73e5cebSJed Brown   struct DMCompositeLink *link;
403f73e5cebSJed Brown   PetscInt                i, wnum;
404f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
405bee642f7SBarry Smith   PetscInt                readonly;
40671b14b3eSStefano Zampini   PetscBool               flg;
407f73e5cebSJed Brown 
408f73e5cebSJed Brown   PetscFunctionBegin;
409f73e5cebSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
410f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4127a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
41348a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
414f73e5cebSJed Brown 
4159566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
416f73e5cebSJed Brown   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
417f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
4189566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
41948a46eb9SPierre Jolivet       if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4209566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum]));
421f73e5cebSJed Brown       wnum++;
422f73e5cebSJed Brown     }
423f73e5cebSJed Brown   }
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
425f73e5cebSJed Brown }
426f73e5cebSJed Brown 
4275d83a8b1SBarry Smith /*@
428dce8aebaSBarry Smith   DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with `DMCompositeGetLocalAccessArray()`.
4297ac2b803SAlex Fikl 
43020f4b53cSBarry Smith   Collective
4317ac2b803SAlex Fikl 
4327ac2b803SAlex Fikl   Input Parameters:
433dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE` object
4347ac2b803SAlex Fikl . pvec    - packed vector
4357ac2b803SAlex Fikl . nwanted - number of vectors wanted
436f13dfd9eSBarry Smith . wanted  - sorted array of vectors wanted, or `NULL` to restore all vectors
437f13dfd9eSBarry Smith - vecs    - array of local vectors
4387ac2b803SAlex Fikl 
4397ac2b803SAlex Fikl   Level: advanced
4407ac2b803SAlex Fikl 
441dce8aebaSBarry Smith   Note:
442f13dfd9eSBarry Smith   `nwanted` and `wanted` must match the values given to `DMCompositeGetLocalAccessArray()`
4437ac2b803SAlex Fikl   otherwise the call will fail.
4447ac2b803SAlex Fikl 
445dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`,
446db781477SPatrick Sanan           `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`,
447db781477SPatrick Sanan           `DMCompositeScatter()`, `DMCompositeGather()`
4487ac2b803SAlex Fikl @*/
4495d83a8b1SBarry Smith PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec *vecs)
450d71ae5a4SJacob Faibussowitsch {
4517ac2b803SAlex Fikl   struct DMCompositeLink *link;
4527ac2b803SAlex Fikl   PetscInt                i, wnum;
4537ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite *)dm->data;
4547ac2b803SAlex Fikl   PetscInt                readonly;
45571b14b3eSStefano Zampini   PetscBool               flg;
4567ac2b803SAlex Fikl 
4577ac2b803SAlex Fikl   PetscFunctionBegin;
4587ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4597ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4609566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4617a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
46248a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
4637ac2b803SAlex Fikl 
4649566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
4657ac2b803SAlex Fikl   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
4667ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
4679566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
46848a46eb9SPierre Jolivet       if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4699566063dSJacob Faibussowitsch       PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum]));
4707ac2b803SAlex Fikl       wnum++;
4717ac2b803SAlex Fikl     }
4727ac2b803SAlex Fikl   }
4733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4747ac2b803SAlex Fikl }
4757ac2b803SAlex Fikl 
47647c6ae99SBarry Smith /*@C
47747c6ae99SBarry Smith   DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
47847c6ae99SBarry Smith 
47920f4b53cSBarry Smith   Collective
48047c6ae99SBarry Smith 
4819ae5db72SJed Brown   Input Parameters:
482dce8aebaSBarry Smith + dm   - the `DMCOMPOSITE` object
48347c6ae99SBarry Smith . gvec - the global vector
484a4e35b19SJacob Faibussowitsch - ...  - the individual sequential vectors, `NULL` for those that are not needed
48547c6ae99SBarry Smith 
48647c6ae99SBarry Smith   Level: advanced
48747c6ae99SBarry Smith 
488dce8aebaSBarry Smith   Note:
489dce8aebaSBarry Smith   `DMCompositeScatterArray()` is a non-variadic alternative that is often more convenient for library callers and is
4906f3c3dcfSJed Brown   accessible from Fortran.
4916f3c3dcfSJed Brown 
492dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
493db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
494db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
495db781477SPatrick Sanan          `DMCompositeScatterArray()`
49647c6ae99SBarry Smith @*/
497d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...)
498d71ae5a4SJacob Faibussowitsch {
49947c6ae99SBarry Smith   va_list                 Argp;
50047c6ae99SBarry Smith   struct DMCompositeLink *next;
501c599c493SJunchao Zhang   PETSC_UNUSED PetscInt   cnt;
50247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
50371b14b3eSStefano Zampini   PetscBool               flg;
50447c6ae99SBarry Smith 
50547c6ae99SBarry Smith   PetscFunctionBegin;
50647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
50747c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5097a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
51048a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
51147c6ae99SBarry Smith 
51215229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
51347c6ae99SBarry Smith   va_start(Argp, gvec);
5148fd8f222SJed Brown   for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
5159ae5db72SJed Brown     Vec local;
5169ae5db72SJed Brown     local = va_arg(Argp, Vec);
5179ae5db72SJed Brown     if (local) {
5189ae5db72SJed Brown       Vec                global;
5195edff71fSBarry Smith       const PetscScalar *array;
52063a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
5219566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
5229566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
5239566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
5249566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local));
5259566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local));
5269566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
5279566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5289566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
52947c6ae99SBarry Smith     }
53047c6ae99SBarry Smith   }
53147c6ae99SBarry Smith   va_end(Argp);
5323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53347c6ae99SBarry Smith }
53447c6ae99SBarry Smith 
5356f3c3dcfSJed Brown /*@
5366f3c3dcfSJed Brown   DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5376f3c3dcfSJed Brown 
53820f4b53cSBarry Smith   Collective
5396f3c3dcfSJed Brown 
5406f3c3dcfSJed Brown   Input Parameters:
541dce8aebaSBarry Smith + dm    - the `DMCOMPOSITE` object
5426f3c3dcfSJed Brown . gvec  - the global vector
543a2b725a8SWilliam Gropp - lvecs - array of local vectors, NULL for any that are not needed
5446f3c3dcfSJed Brown 
5456f3c3dcfSJed Brown   Level: advanced
5466f3c3dcfSJed Brown 
5476f3c3dcfSJed Brown   Note:
548dce8aebaSBarry Smith   This is a non-variadic alternative to `DMCompositeScatter()`
5496f3c3dcfSJed Brown 
550dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`
551db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
552db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
5536f3c3dcfSJed Brown @*/
554d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs)
555d71ae5a4SJacob Faibussowitsch {
5566f3c3dcfSJed Brown   struct DMCompositeLink *next;
5576f3c3dcfSJed Brown   PetscInt                i;
5586f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
55971b14b3eSStefano Zampini   PetscBool               flg;
5606f3c3dcfSJed Brown 
5616f3c3dcfSJed Brown   PetscFunctionBegin;
5626f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5636f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5657a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
56648a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
5676f3c3dcfSJed Brown 
56815229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
5696f3c3dcfSJed Brown   for (i = 0, next = com->next; next; next = next->next, i++) {
5706f3c3dcfSJed Brown     if (lvecs[i]) {
5716f3c3dcfSJed Brown       Vec                global;
572c5d31e75SLisandro Dalcin       const PetscScalar *array;
5736f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3);
5749566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
5759566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
5769566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart));
5779566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i]));
5789566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i]));
5799566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
5809566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5819566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
5826f3c3dcfSJed Brown     }
5836f3c3dcfSJed Brown   }
5843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5856f3c3dcfSJed Brown }
5866f3c3dcfSJed Brown 
58747c6ae99SBarry Smith /*@C
58847c6ae99SBarry Smith   DMCompositeGather - Gathers into a global packed vector from its individual local vectors
58947c6ae99SBarry Smith 
59020f4b53cSBarry Smith   Collective
59147c6ae99SBarry Smith 
592d8d19677SJose E. Roman   Input Parameters:
593dce8aebaSBarry Smith + dm    - the `DMCOMPOSITE` object
594dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES`
595a4e35b19SJacob Faibussowitsch . gvec  - the global vector
596a4e35b19SJacob Faibussowitsch - ...   - the individual sequential vectors, `NULL` for any that are not needed
59747c6ae99SBarry Smith 
59847c6ae99SBarry Smith   Level: advanced
59947c6ae99SBarry Smith 
60060225df5SJacob Faibussowitsch   Fortran Notes:
601dce8aebaSBarry Smith   Fortran users should use `DMCompositeGatherArray()`
602f5f57ec0SBarry Smith 
603dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
604db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
605db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
60647c6ae99SBarry Smith @*/
607d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...)
608d71ae5a4SJacob Faibussowitsch {
60947c6ae99SBarry Smith   va_list                 Argp;
61047c6ae99SBarry Smith   struct DMCompositeLink *next;
61147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
612c599c493SJunchao Zhang   PETSC_UNUSED PetscInt   cnt;
61371b14b3eSStefano Zampini   PetscBool               flg;
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith   PetscFunctionBegin;
61647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
617064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6189566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6197a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
62048a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
62147c6ae99SBarry Smith 
62215229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
6231dac896bSSatish Balay   va_start(Argp, gvec);
6248fd8f222SJed Brown   for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
6259ae5db72SJed Brown     Vec local;
6269ae5db72SJed Brown     local = va_arg(Argp, Vec);
6279ae5db72SJed Brown     if (local) {
62847c6ae99SBarry Smith       PetscScalar *array;
6299ae5db72SJed Brown       Vec          global;
63063a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
6319566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
6329566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec, &array));
6339566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
6349566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global));
6359566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global));
6369566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec, &array));
6379566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6389566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
63947c6ae99SBarry Smith     }
64047c6ae99SBarry Smith   }
64147c6ae99SBarry Smith   va_end(Argp);
6423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64347c6ae99SBarry Smith }
64447c6ae99SBarry Smith 
6456f3c3dcfSJed Brown /*@
6466f3c3dcfSJed Brown   DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6476f3c3dcfSJed Brown 
64820f4b53cSBarry Smith   Collective
6496f3c3dcfSJed Brown 
650d8d19677SJose E. Roman   Input Parameters:
651dce8aebaSBarry Smith + dm    - the `DMCOMPOSITE` object
6526f3c3dcfSJed Brown . gvec  - the global vector
653dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES`
6546f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed
6556f3c3dcfSJed Brown 
6566f3c3dcfSJed Brown   Level: advanced
6576f3c3dcfSJed Brown 
658dce8aebaSBarry Smith   Note:
659dce8aebaSBarry Smith   This is a non-variadic alternative to `DMCompositeGather()`.
6606f3c3dcfSJed Brown 
661dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
662db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
663db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`,
6646f3c3dcfSJed Brown @*/
665d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs)
666d71ae5a4SJacob Faibussowitsch {
6676f3c3dcfSJed Brown   struct DMCompositeLink *next;
6686f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
6696f3c3dcfSJed Brown   PetscInt                i;
67071b14b3eSStefano Zampini   PetscBool               flg;
6716f3c3dcfSJed Brown 
6726f3c3dcfSJed Brown   PetscFunctionBegin;
6736f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
674064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6759566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6767a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
67748a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
6786f3c3dcfSJed Brown 
67915229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
6806f3c3dcfSJed Brown   for (next = com->next, i = 0; next; next = next->next, i++) {
6816f3c3dcfSJed Brown     if (lvecs[i]) {
6826f3c3dcfSJed Brown       PetscScalar *array;
6836f3c3dcfSJed Brown       Vec          global;
684064a246eSJacob Faibussowitsch       PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4);
6859566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
6869566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec, &array));
6879566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
6889566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global));
6899566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global));
6909566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec, &array));
6919566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6929566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
6936f3c3dcfSJed Brown     }
6946f3c3dcfSJed Brown   }
6953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6966f3c3dcfSJed Brown }
6976f3c3dcfSJed Brown 
698f5f57ec0SBarry Smith /*@
699dce8aebaSBarry Smith   DMCompositeAddDM - adds a `DM` vector to a `DMCOMPOSITE`
70047c6ae99SBarry Smith 
70120f4b53cSBarry Smith   Collective
70247c6ae99SBarry Smith 
703d8d19677SJose E. Roman   Input Parameters:
704dce8aebaSBarry Smith + dmc - the  `DMCOMPOSITE` object
705dce8aebaSBarry Smith - dm  - the `DM` object
70647c6ae99SBarry Smith 
70747c6ae99SBarry Smith   Level: advanced
70847c6ae99SBarry Smith 
70942747ad1SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeGather()`, `DMCreateGlobalVector()`,
710db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
711db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
71247c6ae99SBarry Smith @*/
713d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeAddDM(DM dmc, DM dm)
714d71ae5a4SJacob Faibussowitsch {
71506ebdd98SJed Brown   PetscInt                n, nlocal;
71647c6ae99SBarry Smith   struct DMCompositeLink *mine, *next;
71706ebdd98SJed Brown   Vec                     global, local;
71847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dmc->data;
71971b14b3eSStefano Zampini   PetscBool               flg;
72047c6ae99SBarry Smith 
72147c6ae99SBarry Smith   PetscFunctionBegin;
72247c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc, DM_CLASSID, 1);
72347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
7249566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg));
7257a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
72647c6ae99SBarry Smith   next = com->next;
72728b400f6SJacob Faibussowitsch   PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite");
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith   /* create new link */
7309566063dSJacob Faibussowitsch   PetscCall(PetscNew(&mine));
7319566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
732cac35740SJunchao Zhang   PetscCall(DMCreateGlobalVector(dm, &global)); // Not using DMGetGlobalVector() since it will stash the vector with a type decided by dm,
733cac35740SJunchao Zhang   PetscCall(VecGetLocalSize(global, &n));       // while we may want to set dmc's vectype later, say via DMSetFromOptions(dmc), and we
734cac35740SJunchao Zhang   PetscCall(VecDestroy(&global));               // want to propagate the type to dm.
735cac35740SJunchao Zhang   PetscCall(DMCreateLocalVector(dm, &local));   // Not using DMGetLocalVector(), same reason as above.
7369566063dSJacob Faibussowitsch   PetscCall(VecGetSize(local, &nlocal));
737cac35740SJunchao Zhang   PetscCall(VecDestroy(&local));
7388865f1eaSKarl Rupp 
73947c6ae99SBarry Smith   mine->n      = n;
74006ebdd98SJed Brown   mine->nlocal = nlocal;
74147c6ae99SBarry Smith   mine->dm     = dm;
7420298fd71SBarry Smith   mine->next   = NULL;
74347c6ae99SBarry Smith   com->n += n;
7447ac2b803SAlex Fikl   com->nghost += nlocal;
74547c6ae99SBarry Smith 
74647c6ae99SBarry Smith   /* add to end of list */
7478865f1eaSKarl Rupp   if (!next) com->next = mine;
7488865f1eaSKarl Rupp   else {
74947c6ae99SBarry Smith     while (next->next) next = next->next;
75047c6ae99SBarry Smith     next->next = mine;
75147c6ae99SBarry Smith   }
75247c6ae99SBarry Smith   com->nDM++;
75347c6ae99SBarry Smith   com->nmine++;
7543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75547c6ae99SBarry Smith }
75647c6ae99SBarry Smith 
7579804daf3SBarry Smith #include <petscdraw.h>
758d6acfc2dSPierre Jolivet PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
759d6acfc2dSPierre Jolivet 
76066976f2fSJacob Faibussowitsch static PetscErrorCode VecView_DMComposite(Vec gvec, PetscViewer viewer)
761d71ae5a4SJacob Faibussowitsch {
76247c6ae99SBarry Smith   DM                      dm;
76347c6ae99SBarry Smith   struct DMCompositeLink *next;
76447c6ae99SBarry Smith   PetscBool               isdraw;
765cef07954SSatish Balay   DM_Composite           *com;
76647c6ae99SBarry Smith 
76747c6ae99SBarry Smith   PetscFunctionBegin;
7689566063dSJacob Faibussowitsch   PetscCall(VecGetDM(gvec, &dm));
7697a8be351SBarry Smith   PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite");
77047c6ae99SBarry Smith   com  = (DM_Composite *)dm->data;
77147c6ae99SBarry Smith   next = com->next;
77247c6ae99SBarry Smith 
7739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
77447c6ae99SBarry Smith   if (!isdraw) {
77547c6ae99SBarry Smith     /* do I really want to call this? */
7769566063dSJacob Faibussowitsch     PetscCall(VecView_MPI(gvec, viewer));
77747c6ae99SBarry Smith   } else {
77847c6ae99SBarry Smith     PetscInt cnt = 0;
77947c6ae99SBarry Smith 
78015229ffcSPierre Jolivet     /* loop over packed objects, handling one at a time */
78147c6ae99SBarry Smith     while (next) {
78247c6ae99SBarry Smith       Vec                vec;
783ca4278abSLisandro Dalcin       const PetscScalar *array;
78447c6ae99SBarry Smith       PetscInt           bs;
78547c6ae99SBarry Smith 
7869ae5db72SJed Brown       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
7879566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &vec));
7889566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
7899566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart));
7909566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
7919566063dSJacob Faibussowitsch       PetscCall(VecView(vec, viewer));
7929566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vec));
7939566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(vec, &bs));
7949566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &vec));
7959566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawBaseAdd(viewer, bs));
79647c6ae99SBarry Smith       cnt += bs;
79747c6ae99SBarry Smith       next = next->next;
79847c6ae99SBarry Smith     }
7999566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt));
80047c6ae99SBarry Smith   }
8013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80247c6ae99SBarry Smith }
80347c6ae99SBarry Smith 
80466976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec)
805d71ae5a4SJacob Faibussowitsch {
80647c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith   PetscFunctionBegin;
80947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8109566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
8119566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8129566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec));
8139566063dSJacob Faibussowitsch   PetscCall(VecSetType(*gvec, dm->vectype));
8149566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*gvec, com->n, com->N));
8159566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*gvec, dm));
8169566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (void (*)(void))VecView_DMComposite));
8173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
81847c6ae99SBarry Smith }
81947c6ae99SBarry Smith 
82066976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec)
821d71ae5a4SJacob Faibussowitsch {
82247c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
82347c6ae99SBarry Smith 
82447c6ae99SBarry Smith   PetscFunctionBegin;
82547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
82647c6ae99SBarry Smith   if (!com->setup) {
8279566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(dm));
8289566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
82947c6ae99SBarry Smith   }
8309566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, lvec));
8319566063dSJacob Faibussowitsch   PetscCall(VecSetType(*lvec, dm->vectype));
8329566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE));
8339566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*lvec, dm));
8343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83547c6ae99SBarry Smith }
83647c6ae99SBarry Smith 
83747c6ae99SBarry Smith /*@C
838dce8aebaSBarry Smith   DMCompositeGetISLocalToGlobalMappings - gets an `ISLocalToGlobalMapping` for each `DM` in the `DMCOMPOSITE`, maps to the composite global space
83947c6ae99SBarry Smith 
84020f4b53cSBarry Smith   Collective; No Fortran Support
84147c6ae99SBarry Smith 
84247c6ae99SBarry Smith   Input Parameter:
843dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
84447c6ae99SBarry Smith 
8452fe279fdSBarry Smith   Output Parameter:
8469ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes
847dce8aebaSBarry Smith            all the ghost points that individual ghosted `DMDA` may have.
84847c6ae99SBarry Smith 
84947c6ae99SBarry Smith   Level: advanced
85047c6ae99SBarry Smith 
851dce8aebaSBarry Smith   Note:
852f13dfd9eSBarry Smith   Each entry of `ltogs` should be destroyed with `ISLocalToGlobalMappingDestroy()`, `ltogs` should be freed with `PetscFree()`.
85347c6ae99SBarry Smith 
854dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
855db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
856c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
85747c6ae99SBarry Smith @*/
858cc4c1da9SBarry Smith PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping *ltogs[])
859d71ae5a4SJacob Faibussowitsch {
86047c6ae99SBarry Smith   PetscInt                i, *idx, n, cnt;
86147c6ae99SBarry Smith   struct DMCompositeLink *next;
86247c6ae99SBarry Smith   PetscMPIInt             rank;
86347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
86471b14b3eSStefano Zampini   PetscBool               flg;
86547c6ae99SBarry Smith 
86647c6ae99SBarry Smith   PetscFunctionBegin;
86747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
8697a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
8709566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM, ltogs));
87247c6ae99SBarry Smith   next = com->next;
8739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
87447c6ae99SBarry Smith 
87515229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
87647c6ae99SBarry Smith   cnt = 0;
87747c6ae99SBarry Smith   while (next) {
8786eb61c8cSJed Brown     ISLocalToGlobalMapping ltog;
8796eb61c8cSJed Brown     PetscMPIInt            size;
88086994e45SJed Brown     const PetscInt        *suboff, *indices;
8816eb61c8cSJed Brown     Vec                    global;
88247c6ae99SBarry Smith 
8836eb61c8cSJed Brown     /* Get sub-DM global indices for each local dof */
8849566063dSJacob Faibussowitsch     PetscCall(DMGetLocalToGlobalMapping(next->dm, &ltog));
8859566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n));
8869566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices));
8879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &idx));
88847c6ae99SBarry Smith 
8896eb61c8cSJed Brown     /* Get the offsets for the sub-DM global vector */
8909566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
8919566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRanges(global, &suboff));
8929566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size));
8936eb61c8cSJed Brown 
8946eb61c8cSJed Brown     /* Shift the sub-DM definition of the global space to the composite global space */
8956eb61c8cSJed Brown     for (i = 0; i < n; i++) {
89686994e45SJed Brown       PetscInt subi = indices[i], lo = 0, hi = size, t;
897e333b035SStefano Zampini       /* There's no consensus on what a negative index means,
898e333b035SStefano Zampini          except for skipping when setting the values in vectors and matrices */
8999371c9d4SSatish Balay       if (subi < 0) {
9009371c9d4SSatish Balay         idx[i] = subi - next->grstarts[rank];
9019371c9d4SSatish Balay         continue;
9029371c9d4SSatish Balay       }
9036eb61c8cSJed Brown       /* Binary search to find which rank owns subi */
9046eb61c8cSJed Brown       while (hi - lo > 1) {
9056eb61c8cSJed Brown         t = lo + (hi - lo) / 2;
9066eb61c8cSJed Brown         if (suboff[t] > subi) hi = t;
9076eb61c8cSJed Brown         else lo = t;
9086eb61c8cSJed Brown       }
9096eb61c8cSJed Brown       idx[i] = subi - suboff[lo] + next->grstarts[lo];
9106eb61c8cSJed Brown     }
9119566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices));
9129566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt]));
9139566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
91447c6ae99SBarry Smith     next = next->next;
91547c6ae99SBarry Smith     cnt++;
91647c6ae99SBarry Smith   }
9173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91847c6ae99SBarry Smith }
91947c6ae99SBarry Smith 
92087c85e80SJed Brown /*@C
9219ae5db72SJed Brown   DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
92287c85e80SJed Brown 
92320f4b53cSBarry Smith   Not Collective; No Fortran Support
92487c85e80SJed Brown 
9254165533cSJose E. Roman   Input Parameter:
926dce8aebaSBarry Smith . dm - the `DMCOMPOSITE`
92787c85e80SJed Brown 
9284165533cSJose E. Roman   Output Parameter:
92915229ffcSPierre Jolivet . is - array of serial index sets for each component of the `DMCOMPOSITE`
93087c85e80SJed Brown 
93187c85e80SJed Brown   Level: intermediate
93287c85e80SJed Brown 
93387c85e80SJed Brown   Notes:
93487c85e80SJed Brown   At present, a composite local vector does not normally exist.  This function is used to provide index sets for
93515229ffcSPierre Jolivet   `MatGetLocalSubMatrix()`.  In the future, the scatters for each entry in the `DMCOMPOSITE` may be merged into a single
9369ae5db72SJed Brown   scatter to a composite local vector.  The user should not typically need to know which is being done.
93787c85e80SJed Brown 
938dce8aebaSBarry Smith   To get the composite global indices at all local points (including ghosts), use `DMCompositeGetISLocalToGlobalMappings()`.
93987c85e80SJed Brown 
940dce8aebaSBarry Smith   To get index sets for pieces of the composite global vector, use `DMCompositeGetGlobalISs()`.
94187c85e80SJed Brown 
942dce8aebaSBarry Smith   Each returned `IS` should be destroyed with `ISDestroy()`, the array should be freed with `PetscFree()`.
94387c85e80SJed Brown 
944f13dfd9eSBarry Smith   Fortran Note:
945ce78bad3SBarry Smith   Use `DMCompositeRestoreLocalISs()` to release the `is`.
946f13dfd9eSBarry Smith 
947f13dfd9eSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`,
948f13dfd9eSBarry Smith           `MatCreateLocalRef()`, `DMCompositeGetNumberDM()`
94987c85e80SJed Brown @*/
950cc4c1da9SBarry Smith PetscErrorCode DMCompositeGetLocalISs(DM dm, IS *is[])
951d71ae5a4SJacob Faibussowitsch {
95287c85e80SJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
95387c85e80SJed Brown   struct DMCompositeLink *link;
95487c85e80SJed Brown   PetscInt                cnt, start;
95571b14b3eSStefano Zampini   PetscBool               flg;
95687c85e80SJed Brown 
95787c85e80SJed Brown   PetscFunctionBegin;
95887c85e80SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9594f572ea9SToby Isaac   PetscAssertPointer(is, 2);
9609566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
9617a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
9629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nmine, is));
96306ebdd98SJed Brown   for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) {
964520db06cSJed Brown     PetscInt bs;
9659566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt]));
9669566063dSJacob Faibussowitsch     PetscCall(DMGetBlockSize(link->dm, &bs));
9679566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize((*is)[cnt], bs));
968520db06cSJed Brown   }
9693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97087c85e80SJed Brown }
97187c85e80SJed Brown 
97247c6ae99SBarry Smith /*@C
973dce8aebaSBarry Smith   DMCompositeGetGlobalISs - Gets the index sets for each composed object in a `DMCOMPOSITE`
97447c6ae99SBarry Smith 
97520f4b53cSBarry Smith   Collective
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith   Input Parameter:
978dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
97947c6ae99SBarry Smith 
9802fe279fdSBarry Smith   Output Parameter:
98147c6ae99SBarry Smith . is - the array of index sets
98247c6ae99SBarry Smith 
98347c6ae99SBarry Smith   Level: advanced
98447c6ae99SBarry Smith 
98547c6ae99SBarry Smith   Notes:
986f13dfd9eSBarry Smith   The `is` entries should be destroyed with `ISDestroy()`, `is` should be freed with `PetscFree()`
98747c6ae99SBarry Smith 
98847c6ae99SBarry Smith   These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
98947c6ae99SBarry Smith 
990dce8aebaSBarry Smith   Use `DMCompositeGetLocalISs()` for index sets in the packed local numbering, and
991dce8aebaSBarry Smith   `DMCompositeGetISLocalToGlobalMappings()` for to map local sub-`DM` (including ghost) indices to packed global
9926eb61c8cSJed Brown   indices.
99347c6ae99SBarry Smith 
994ce78bad3SBarry Smith   Fortran Note:
995ce78bad3SBarry Smith   Use `DMCompositeRestoreGlobalISs()` to release the `is`.
996f3cb0f7eSJed Brown 
997dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
998db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
999c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
100047c6ae99SBarry Smith @*/
1001d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[])
1002d71ae5a4SJacob Faibussowitsch {
100366bb578eSMark Adams   PetscInt                cnt = 0;
100447c6ae99SBarry Smith   struct DMCompositeLink *next;
100547c6ae99SBarry Smith   PetscMPIInt             rank;
100647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
100771b14b3eSStefano Zampini   PetscBool               flg;
100847c6ae99SBarry Smith 
100947c6ae99SBarry Smith   PetscFunctionBegin;
101047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
10127a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
10137a8be351SBarry Smith   PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before");
10149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM, is));
101547c6ae99SBarry Smith   next = com->next;
10169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
101747c6ae99SBarry Smith 
101815229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
101947c6ae99SBarry Smith   while (next) {
1020e5e52638SMatthew G. Knepley     PetscDS prob;
1021e5e52638SMatthew G. Knepley 
10229566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt]));
10239566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &prob));
1024e5e52638SMatthew G. Knepley     if (prob) {
102565c226d8SMatthew G. Knepley       MatNullSpace space;
102665c226d8SMatthew G. Knepley       Mat          pmat;
10270f21e855SMatthew G. Knepley       PetscObject  disc;
10280f21e855SMatthew G. Knepley       PetscInt     Nf;
102965c226d8SMatthew G. Knepley 
10309566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(prob, &Nf));
1031f24dd8d2SMatthew G. Knepley       if (cnt < Nf) {
10329566063dSJacob Faibussowitsch         PetscCall(PetscDSGetDiscretization(prob, cnt, &disc));
10339566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space));
10349566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space));
10359566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space));
10369566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space));
10379566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat));
10389566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat));
103965c226d8SMatthew G. Knepley       }
1040f24dd8d2SMatthew G. Knepley     }
104147c6ae99SBarry Smith     cnt++;
104247c6ae99SBarry Smith     next = next->next;
104347c6ae99SBarry Smith   }
10443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104547c6ae99SBarry Smith }
104647c6ae99SBarry Smith 
104766976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
1048d71ae5a4SJacob Faibussowitsch {
10494d343eeaSMatthew G Knepley   PetscInt nDM;
10504d343eeaSMatthew G Knepley   DM      *dms;
10514d343eeaSMatthew G Knepley   PetscInt i;
10524d343eeaSMatthew G Knepley 
10534d343eeaSMatthew G Knepley   PetscFunctionBegin;
10549566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm, &nDM));
10558865f1eaSKarl Rupp   if (numFields) *numFields = nDM;
10569566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetGlobalISs(dm, fields));
10574d343eeaSMatthew G Knepley   if (fieldNames) {
10589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, &dms));
10599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, fieldNames));
10609566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, dms));
10614d343eeaSMatthew G Knepley     for (i = 0; i < nDM; i++) {
10624d343eeaSMatthew G Knepley       char        buf[256];
10634d343eeaSMatthew G Knepley       const char *splitname;
10644d343eeaSMatthew G Knepley 
10654d343eeaSMatthew G Knepley       /* Split naming precedence: object name, prefix, number */
10664d343eeaSMatthew G Knepley       splitname = ((PetscObject)dm)->name;
10674d343eeaSMatthew G Knepley       if (!splitname) {
10689566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname));
10694d343eeaSMatthew G Knepley         if (splitname) {
10704d343eeaSMatthew G Knepley           size_t len;
10719566063dSJacob Faibussowitsch           PetscCall(PetscStrncpy(buf, splitname, sizeof(buf)));
10728caf3d72SBarry Smith           buf[sizeof(buf) - 1] = 0;
10739566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(buf, &len));
10744d343eeaSMatthew G Knepley           if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */
10754d343eeaSMatthew G Knepley           splitname = buf;
10764d343eeaSMatthew G Knepley         }
10774d343eeaSMatthew G Knepley       }
10784d343eeaSMatthew G Knepley       if (!splitname) {
107963a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i));
10804d343eeaSMatthew G Knepley         splitname = buf;
10814d343eeaSMatthew G Knepley       }
10829566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i]));
10834d343eeaSMatthew G Knepley     }
10849566063dSJacob Faibussowitsch     PetscCall(PetscFree(dms));
10854d343eeaSMatthew G Knepley   }
10863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10874d343eeaSMatthew G Knepley }
10884d343eeaSMatthew G Knepley 
1089e7c4fc90SDmitry Karpeev /*
1090e7c4fc90SDmitry Karpeev  This could take over from DMCreateFieldIS(), as it is more general,
10910298fd71SBarry Smith  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1092e7c4fc90SDmitry Karpeev  At this point it's probably best to be less intrusive, however.
1093e7c4fc90SDmitry Karpeev  */
109466976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
1095d71ae5a4SJacob Faibussowitsch {
1096e7c4fc90SDmitry Karpeev   PetscInt nDM;
1097e7c4fc90SDmitry Karpeev   PetscInt i;
1098e7c4fc90SDmitry Karpeev 
1099e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
11009566063dSJacob Faibussowitsch   PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist));
1101e7c4fc90SDmitry Karpeev   if (dmlist) {
11029566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetNumberDM(dm, &nDM));
11039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, dmlist));
11049566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, *dmlist));
110548a46eb9SPierre Jolivet     for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i])));
1106e7c4fc90SDmitry Karpeev   }
11073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1108e7c4fc90SDmitry Karpeev }
1109e7c4fc90SDmitry Karpeev 
111047c6ae99SBarry Smith /*@C
1111dce8aebaSBarry Smith   DMCompositeGetLocalVectors - Gets local vectors for each part of a `DMCOMPOSITE`
1112dce8aebaSBarry Smith   Use `DMCompositeRestoreLocalVectors()` to return them.
111347c6ae99SBarry Smith 
111420f4b53cSBarry Smith   Not Collective; No Fortran Support
111547c6ae99SBarry Smith 
111647c6ae99SBarry Smith   Input Parameter:
1117dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
111847c6ae99SBarry Smith 
111947c6ae99SBarry Smith   Output Parameter:
1120a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`s
112147c6ae99SBarry Smith 
112247c6ae99SBarry Smith   Level: advanced
112347c6ae99SBarry Smith 
1124dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1125db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1126db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
112747c6ae99SBarry Smith @*/
1128d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...)
1129d71ae5a4SJacob Faibussowitsch {
113047c6ae99SBarry Smith   va_list                 Argp;
113147c6ae99SBarry Smith   struct DMCompositeLink *next;
113247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
113371b14b3eSStefano Zampini   PetscBool               flg;
113447c6ae99SBarry Smith 
113547c6ae99SBarry Smith   PetscFunctionBegin;
113647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11379566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11387a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
113947c6ae99SBarry Smith   next = com->next;
114015229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
114147c6ae99SBarry Smith   va_start(Argp, dm);
114247c6ae99SBarry Smith   while (next) {
114347c6ae99SBarry Smith     Vec *vec;
114447c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
11459566063dSJacob Faibussowitsch     if (vec) PetscCall(DMGetLocalVector(next->dm, vec));
114647c6ae99SBarry Smith     next = next->next;
114747c6ae99SBarry Smith   }
114847c6ae99SBarry Smith   va_end(Argp);
11493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115047c6ae99SBarry Smith }
115147c6ae99SBarry Smith 
115247c6ae99SBarry Smith /*@C
1153dce8aebaSBarry Smith   DMCompositeRestoreLocalVectors - Restores local vectors for each part of a `DMCOMPOSITE`
115447c6ae99SBarry Smith 
115520f4b53cSBarry Smith   Not Collective; No Fortran Support
115647c6ae99SBarry Smith 
115747c6ae99SBarry Smith   Input Parameter:
1158dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
115947c6ae99SBarry Smith 
116047c6ae99SBarry Smith   Output Parameter:
1161a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`
116247c6ae99SBarry Smith 
116347c6ae99SBarry Smith   Level: advanced
116447c6ae99SBarry Smith 
1165dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1166db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1167db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
116847c6ae99SBarry Smith @*/
1169d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...)
1170d71ae5a4SJacob Faibussowitsch {
117147c6ae99SBarry Smith   va_list                 Argp;
117247c6ae99SBarry Smith   struct DMCompositeLink *next;
117347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
117471b14b3eSStefano Zampini   PetscBool               flg;
117547c6ae99SBarry Smith 
117647c6ae99SBarry Smith   PetscFunctionBegin;
117747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11789566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11797a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
118047c6ae99SBarry Smith   next = com->next;
118115229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
118247c6ae99SBarry Smith   va_start(Argp, dm);
118347c6ae99SBarry Smith   while (next) {
118447c6ae99SBarry Smith     Vec *vec;
118547c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
11869566063dSJacob Faibussowitsch     if (vec) PetscCall(DMRestoreLocalVector(next->dm, vec));
118747c6ae99SBarry Smith     next = next->next;
118847c6ae99SBarry Smith   }
118947c6ae99SBarry Smith   va_end(Argp);
11903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119147c6ae99SBarry Smith }
119247c6ae99SBarry Smith 
119347c6ae99SBarry Smith /*@C
1194dce8aebaSBarry Smith   DMCompositeGetEntries - Gets the `DM` for each entry in a `DMCOMPOSITE`.
119547c6ae99SBarry Smith 
119647c6ae99SBarry Smith   Not Collective
119747c6ae99SBarry Smith 
119847c6ae99SBarry Smith   Input Parameter:
1199dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
120047c6ae99SBarry Smith 
120147c6ae99SBarry Smith   Output Parameter:
1202a4e35b19SJacob Faibussowitsch . ... - the individual entries `DM`
120347c6ae99SBarry Smith 
120447c6ae99SBarry Smith   Level: advanced
120547c6ae99SBarry Smith 
120660225df5SJacob Faibussowitsch   Fortran Notes:
12075d83a8b1SBarry Smith   Use `DMCompositeGetEntriesArray()`
1208f5f57ec0SBarry Smith 
1209dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()`
1210db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
121160225df5SJacob Faibussowitsch          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`
121247c6ae99SBarry Smith @*/
1213d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntries(DM dm, ...)
1214d71ae5a4SJacob Faibussowitsch {
121547c6ae99SBarry Smith   va_list                 Argp;
121647c6ae99SBarry Smith   struct DMCompositeLink *next;
121747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
121871b14b3eSStefano Zampini   PetscBool               flg;
121947c6ae99SBarry Smith 
122047c6ae99SBarry Smith   PetscFunctionBegin;
122147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12229566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
12237a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
122447c6ae99SBarry Smith   next = com->next;
122515229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
122647c6ae99SBarry Smith   va_start(Argp, dm);
122747c6ae99SBarry Smith   while (next) {
122847c6ae99SBarry Smith     DM *dmn;
122947c6ae99SBarry Smith     dmn = va_arg(Argp, DM *);
12309ae5db72SJed Brown     if (dmn) *dmn = next->dm;
123147c6ae99SBarry Smith     next = next->next;
123247c6ae99SBarry Smith   }
123347c6ae99SBarry Smith   va_end(Argp);
12343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
123547c6ae99SBarry Smith }
123647c6ae99SBarry Smith 
1237cc4c1da9SBarry Smith /*@
1238dce8aebaSBarry Smith   DMCompositeGetEntriesArray - Gets the DM for each entry in a `DMCOMPOSITE`
12392fa5ba8aSJed Brown 
12402fa5ba8aSJed Brown   Not Collective
12412fa5ba8aSJed Brown 
12422fa5ba8aSJed Brown   Input Parameter:
1243dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
1244907376e6SBarry Smith 
1245907376e6SBarry Smith   Output Parameter:
1246dce8aebaSBarry Smith . dms - array of sufficient length (see `DMCompositeGetNumberDM()`) to hold the individual `DM`
12472fa5ba8aSJed Brown 
12482fa5ba8aSJed Brown   Level: advanced
12492fa5ba8aSJed Brown 
1250dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()`
1251db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
125260225df5SJacob Faibussowitsch          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`
12532fa5ba8aSJed Brown @*/
1254d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntriesArray(DM dm, DM dms[])
1255d71ae5a4SJacob Faibussowitsch {
12562fa5ba8aSJed Brown   struct DMCompositeLink *next;
12572fa5ba8aSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
12582fa5ba8aSJed Brown   PetscInt                i;
125971b14b3eSStefano Zampini   PetscBool               flg;
12602fa5ba8aSJed Brown 
12612fa5ba8aSJed Brown   PetscFunctionBegin;
12622fa5ba8aSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
12647a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
126515229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
12662fa5ba8aSJed Brown   for (next = com->next, i = 0; next; next = next->next, i++) dms[i] = next->dm;
12673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12682fa5ba8aSJed Brown }
12692fa5ba8aSJed Brown 
1270e10fd815SStefano Zampini typedef struct {
1271e10fd815SStefano Zampini   DM           dm;
1272e10fd815SStefano Zampini   PetscViewer *subv;
1273e10fd815SStefano Zampini   Vec         *vecs;
1274e10fd815SStefano Zampini } GLVisViewerCtx;
1275e10fd815SStefano Zampini 
1276d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx)
1277d71ae5a4SJacob Faibussowitsch {
1278e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
1279e10fd815SStefano Zampini   PetscInt        i, n;
1280e10fd815SStefano Zampini 
1281e10fd815SStefano Zampini   PetscFunctionBegin;
12829566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(ctx->dm, &n));
128348a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscViewerDestroy(&ctx->subv[i]));
12849566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ctx->subv, ctx->vecs));
12859566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dm));
12869566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
12873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1288e10fd815SStefano Zampini }
1289e10fd815SStefano Zampini 
1290d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
1291d71ae5a4SJacob Faibussowitsch {
1292e10fd815SStefano Zampini   Vec             X   = (Vec)oX;
1293e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
1294e10fd815SStefano Zampini   PetscInt        i, n, cumf;
1295e10fd815SStefano Zampini 
1296e10fd815SStefano Zampini   PetscFunctionBegin;
12979566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(ctx->dm, &n));
12989566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccessArray(ctx->dm, X, n, NULL, ctx->vecs));
1299e10fd815SStefano Zampini   for (i = 0, cumf = 0; i < n; i++) {
1300e10fd815SStefano Zampini     PetscErrorCode (*g2l)(PetscObject, PetscInt, PetscObject[], void *);
1301e10fd815SStefano Zampini     void    *fctx;
1302e10fd815SStefano Zampini     PetscInt nfi;
1303e10fd815SStefano Zampini 
130434e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nfi, NULL, NULL, &g2l, NULL, &fctx));
1305e10fd815SStefano Zampini     if (!nfi) continue;
13061baa6e33SBarry Smith     if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i], nfi, oXfield + cumf, fctx));
1307f4f49eeaSPierre Jolivet     else PetscCall(VecCopy(ctx->vecs[i], (Vec)oXfield[cumf]));
1308e10fd815SStefano Zampini     cumf += nfi;
1309e10fd815SStefano Zampini   }
13109566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccessArray(ctx->dm, X, n, NULL, ctx->vecs));
13113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1312e10fd815SStefano Zampini }
1313e10fd815SStefano Zampini 
1314d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer)
1315d71ae5a4SJacob Faibussowitsch {
1316e10fd815SStefano Zampini   DM              dm = (DM)odm, *dms;
1317e10fd815SStefano Zampini   Vec            *Ufds;
1318e10fd815SStefano Zampini   GLVisViewerCtx *ctx;
1319e10fd815SStefano Zampini   PetscInt        i, n, tnf, *sdim;
1320e10fd815SStefano Zampini   char          **fecs;
1321e10fd815SStefano Zampini 
1322e10fd815SStefano Zampini   PetscFunctionBegin;
13239566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
13249566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
1325e10fd815SStefano Zampini   ctx->dm = dm;
13269566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm, &n));
13279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &dms));
13289566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetEntriesArray(dm, dms));
13299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n, &ctx->subv, n, &ctx->vecs));
1330e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1331e10fd815SStefano Zampini     PetscInt nf;
1332e10fd815SStefano Zampini 
13339566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PetscObjectComm(odm), &ctx->subv[i]));
13349566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(ctx->subv[i], PETSCVIEWERGLVIS));
133534e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisSetDM_Internal(ctx->subv[i], (PetscObject)dms[i]));
133634e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, NULL, NULL, NULL, NULL, NULL));
1337e10fd815SStefano Zampini     tnf += nf;
1338e10fd815SStefano Zampini   }
13399566063dSJacob Faibussowitsch   PetscCall(PetscFree(dms));
13409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(tnf, &fecs, tnf, &sdim, tnf, &Ufds));
1341e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1342e10fd815SStefano Zampini     PetscInt    *sd, nf, f;
1343e10fd815SStefano Zampini     const char **fec;
1344e10fd815SStefano Zampini     Vec         *Uf;
1345e10fd815SStefano Zampini 
134634e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, &fec, &sd, NULL, (PetscObject **)&Uf, NULL));
1347e10fd815SStefano Zampini     for (f = 0; f < nf; f++) {
13489566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(fec[f], &fecs[tnf + f]));
1349e10fd815SStefano Zampini       Ufds[tnf + f] = Uf[f];
1350e10fd815SStefano Zampini       sdim[tnf + f] = sd[f];
1351e10fd815SStefano Zampini     }
1352e10fd815SStefano Zampini     tnf += nf;
1353e10fd815SStefano Zampini   }
13549566063dSJacob Faibussowitsch   PetscCall(PetscViewerGLVisSetFields(viewer, tnf, (const char **)fecs, sdim, DMCompositeSampleGLVisFields_Private, (PetscObject *)Ufds, ctx, DestroyGLVisViewerCtx_Private));
135548a46eb9SPierre Jolivet   for (i = 0; i < tnf; i++) PetscCall(PetscFree(fecs[i]));
13569566063dSJacob Faibussowitsch   PetscCall(PetscFree3(fecs, sdim, Ufds));
13573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1358e10fd815SStefano Zampini }
1359e10fd815SStefano Zampini 
136066976f2fSJacob Faibussowitsch static PetscErrorCode DMRefine_Composite(DM dmi, MPI_Comm comm, DM *fine)
1361d71ae5a4SJacob Faibussowitsch {
136247c6ae99SBarry Smith   struct DMCompositeLink *next;
136347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dmi->data;
136447c6ae99SBarry Smith   DM                      dm;
136547c6ae99SBarry Smith 
136647c6ae99SBarry Smith   PetscFunctionBegin;
136747c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi, DM_CLASSID, 1);
136848a46eb9SPierre Jolivet   if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm));
13699566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmi));
137047c6ae99SBarry Smith   next = com->next;
13719566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm, fine));
137247c6ae99SBarry Smith 
137315229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
137447c6ae99SBarry Smith   while (next) {
13759566063dSJacob Faibussowitsch     PetscCall(DMRefine(next->dm, comm, &dm));
13769566063dSJacob Faibussowitsch     PetscCall(DMCompositeAddDM(*fine, dm));
13779566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dm));
137847c6ae99SBarry Smith     next = next->next;
137947c6ae99SBarry Smith   }
13803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138147c6ae99SBarry Smith }
138247c6ae99SBarry Smith 
138366976f2fSJacob Faibussowitsch static PetscErrorCode DMCoarsen_Composite(DM dmi, MPI_Comm comm, DM *fine)
1384d71ae5a4SJacob Faibussowitsch {
138514354c39SJed Brown   struct DMCompositeLink *next;
138614354c39SJed Brown   DM_Composite           *com = (DM_Composite *)dmi->data;
138714354c39SJed Brown   DM                      dm;
138814354c39SJed Brown 
138914354c39SJed Brown   PetscFunctionBegin;
139014354c39SJed Brown   PetscValidHeaderSpecific(dmi, DM_CLASSID, 1);
13919566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmi));
139248a46eb9SPierre Jolivet   if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm));
139314354c39SJed Brown   next = com->next;
13949566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm, fine));
139514354c39SJed Brown 
139615229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
139714354c39SJed Brown   while (next) {
13989566063dSJacob Faibussowitsch     PetscCall(DMCoarsen(next->dm, comm, &dm));
13999566063dSJacob Faibussowitsch     PetscCall(DMCompositeAddDM(*fine, dm));
14009566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dm));
140114354c39SJed Brown     next = next->next;
140214354c39SJed Brown   }
14033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
140414354c39SJed Brown }
140547c6ae99SBarry Smith 
140666976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_Composite(DM coarse, DM fine, Mat *A, Vec *v)
1407d71ae5a4SJacob Faibussowitsch {
14089ae5db72SJed Brown   PetscInt                m, n, M, N, nDM, i;
140947c6ae99SBarry Smith   struct DMCompositeLink *nextc;
141047c6ae99SBarry Smith   struct DMCompositeLink *nextf;
141125296bd5SBarry Smith   Vec                     gcoarse, gfine, *vecs;
141247c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite *)coarse->data;
141347c6ae99SBarry Smith   DM_Composite           *comfine   = (DM_Composite *)fine->data;
14149ae5db72SJed Brown   Mat                    *mats;
141547c6ae99SBarry Smith 
141647c6ae99SBarry Smith   PetscFunctionBegin;
141747c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse, DM_CLASSID, 1);
141847c6ae99SBarry Smith   PetscValidHeaderSpecific(fine, DM_CLASSID, 2);
14199566063dSJacob Faibussowitsch   PetscCall(DMSetUp(coarse));
14209566063dSJacob Faibussowitsch   PetscCall(DMSetUp(fine));
142147c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14229566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(coarse, &gcoarse));
14239566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(fine, &gfine));
14249566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gcoarse, &n));
14259566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gfine, &m));
14269566063dSJacob Faibussowitsch   PetscCall(VecGetSize(gcoarse, &N));
14279566063dSJacob Faibussowitsch   PetscCall(VecGetSize(gfine, &M));
14289566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(coarse, &gcoarse));
14299566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(fine, &gfine));
143047c6ae99SBarry Smith 
14319ae5db72SJed Brown   nDM = comfine->nDM;
143263a3b9bcSJacob 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);
14339566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nDM * nDM, &mats));
143448a46eb9SPierre Jolivet   if (v) PetscCall(PetscCalloc1(nDM, &vecs));
143547c6ae99SBarry Smith 
143615229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
14379ae5db72SJed Brown   for (nextc = comcoarse->next, nextf = comfine->next, i = 0; nextc; nextc = nextc->next, nextf = nextf->next, i++) {
14381baa6e33SBarry Smith     if (!v) PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], NULL));
14391baa6e33SBarry Smith     else PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], &vecs[i]));
144047c6ae99SBarry Smith   }
14419566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, nDM, NULL, mats, A));
14421baa6e33SBarry Smith   if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, vecs, v));
14439566063dSJacob Faibussowitsch   for (i = 0; i < nDM * nDM; i++) PetscCall(MatDestroy(&mats[i]));
14449566063dSJacob Faibussowitsch   PetscCall(PetscFree(mats));
144525296bd5SBarry Smith   if (v) {
14469566063dSJacob Faibussowitsch     for (i = 0; i < nDM; i++) PetscCall(VecDestroy(&vecs[i]));
14479566063dSJacob Faibussowitsch     PetscCall(PetscFree(vecs));
144825296bd5SBarry Smith   }
14493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
145047c6ae99SBarry Smith }
145147c6ae99SBarry Smith 
1452d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
1453d71ae5a4SJacob Faibussowitsch {
14541411c6eeSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
14551411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1456f7efa3c7SJed Brown   PetscInt                i;
14571411c6eeSJed Brown 
14581411c6eeSJed Brown   PetscFunctionBegin;
14591411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
14609566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetISLocalToGlobalMappings(dm, &ltogs));
14619566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm), com->nDM, ltogs, &dm->ltogmap));
14629566063dSJacob Faibussowitsch   for (i = 0; i < com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(&ltogs[i]));
14639566063dSJacob Faibussowitsch   PetscCall(PetscFree(ltogs));
14643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14651411c6eeSJed Brown }
14661411c6eeSJed Brown 
146766976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateColoring_Composite(DM dm, ISColoringType ctype, ISColoring *coloring)
1468d71ae5a4SJacob Faibussowitsch {
146947c6ae99SBarry Smith   PetscInt         n, i, cnt;
147047c6ae99SBarry Smith   ISColoringValue *colors;
147147c6ae99SBarry Smith   PetscBool        dense  = PETSC_FALSE;
147247c6ae99SBarry Smith   ISColoringValue  maxcol = 0;
147347c6ae99SBarry Smith   DM_Composite    *com    = (DM_Composite *)dm->data;
147447c6ae99SBarry Smith 
147547c6ae99SBarry Smith   PetscFunctionBegin;
147647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
147708401ef6SPierre Jolivet   PetscCheck(ctype != IS_COLORING_LOCAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only global coloring supported");
1478*966bd95aSPierre Jolivet   PetscCheck(ctype == IS_COLORING_GLOBAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unknown ISColoringType");
147947c6ae99SBarry Smith   n = com->n;
14809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &colors)); /* freed in ISColoringDestroy() */
148147c6ae99SBarry Smith 
14829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL));
148347c6ae99SBarry Smith   if (dense) {
14846497c311SBarry Smith     PetscCall(ISColoringValueCast(com->N, &maxcol));
14856497c311SBarry Smith     for (i = 0; i < n; i++) PetscCall(ISColoringValueCast(com->rstart + i, colors + i));
148647c6ae99SBarry Smith   } else {
148747c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
148847c6ae99SBarry Smith     PetscMPIInt             rank;
148947c6ae99SBarry Smith 
14909566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
149147c6ae99SBarry Smith     cnt = 0;
149247c6ae99SBarry Smith     while (next) {
149347c6ae99SBarry Smith       ISColoring lcoloring;
149447c6ae99SBarry Smith 
14959566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(next->dm, IS_COLORING_GLOBAL, &lcoloring));
14966497c311SBarry Smith       for (i = 0; i < lcoloring->N; i++) PetscCall(ISColoringValueCast(maxcol + lcoloring->colors[i], colors + cnt++));
149747c6ae99SBarry Smith       maxcol += lcoloring->n;
14989566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&lcoloring));
149947c6ae99SBarry Smith       next = next->next;
150047c6ae99SBarry Smith     }
150147c6ae99SBarry Smith   }
15029566063dSJacob Faibussowitsch   PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), maxcol, n, colors, PETSC_OWN_POINTER, coloring));
15033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150447c6ae99SBarry Smith }
150547c6ae99SBarry Smith 
150666976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec)
1507d71ae5a4SJacob Faibussowitsch {
150847c6ae99SBarry Smith   struct DMCompositeLink *next;
150947c6ae99SBarry Smith   PetscScalar            *garray, *larray;
151047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
151147c6ae99SBarry Smith 
151247c6ae99SBarry Smith   PetscFunctionBegin;
151347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
151447c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
151539d35262SVincent Le Chenadec 
151648a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
151739d35262SVincent Le Chenadec 
15189566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gvec, &garray));
15199566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lvec, &larray));
152047c6ae99SBarry Smith 
152115229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
152239d35262SVincent Le Chenadec   next = com->next;
152347c6ae99SBarry Smith   while (next) {
152447c6ae99SBarry Smith     Vec      local, global;
152547c6ae99SBarry Smith     PetscInt N;
152647c6ae99SBarry Smith 
15279566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
15289566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(global, &N));
15299566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(global, garray));
15309566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local));
15319566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local, larray));
15329566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(next->dm, global, mode, local));
15339566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(next->dm, global, mode, local));
15349566063dSJacob Faibussowitsch     PetscCall(VecResetArray(global));
15359566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local));
15369566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
15379566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local));
153839d35262SVincent Le Chenadec 
153906ebdd98SJed Brown     larray += next->nlocal;
154039d35262SVincent Le Chenadec     garray += next->n;
154147c6ae99SBarry Smith     next = next->next;
154247c6ae99SBarry Smith   }
154347c6ae99SBarry Smith 
15449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gvec, NULL));
15459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lvec, NULL));
15463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
154747c6ae99SBarry Smith }
154847c6ae99SBarry Smith 
154966976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec)
1550d71ae5a4SJacob Faibussowitsch {
15510c010503SBarry Smith   PetscFunctionBegin;
155239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
155339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
155439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 4);
15553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
155639d35262SVincent Le Chenadec }
155739d35262SVincent Le Chenadec 
155866976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1559d71ae5a4SJacob Faibussowitsch {
156039d35262SVincent Le Chenadec   struct DMCompositeLink *next;
156139d35262SVincent Le Chenadec   PetscScalar            *larray, *garray;
156239d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite *)dm->data;
156339d35262SVincent Le Chenadec 
156439d35262SVincent Le Chenadec   PetscFunctionBegin;
156539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
156639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
156739d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
156839d35262SVincent Le Chenadec 
156948a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
157039d35262SVincent Le Chenadec 
15719566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lvec, &larray));
15729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gvec, &garray));
157339d35262SVincent Le Chenadec 
157415229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
157539d35262SVincent Le Chenadec   next = com->next;
157639d35262SVincent Le Chenadec   while (next) {
157739d35262SVincent Le Chenadec     Vec global, local;
157839d35262SVincent Le Chenadec 
15799566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local));
15809566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local, larray));
15819566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
15829566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(global, garray));
15839566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(next->dm, local, mode, global));
15849566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(next->dm, local, mode, global));
15859566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local));
15869566063dSJacob Faibussowitsch     PetscCall(VecResetArray(global));
15879566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
15889566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local));
158939d35262SVincent Le Chenadec 
159039d35262SVincent Le Chenadec     garray += next->n;
159139d35262SVincent Le Chenadec     larray += next->nlocal;
159239d35262SVincent Le Chenadec     next = next->next;
159339d35262SVincent Le Chenadec   }
159439d35262SVincent Le Chenadec 
15959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gvec, NULL));
15969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lvec, NULL));
15973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159839d35262SVincent Le Chenadec }
159939d35262SVincent Le Chenadec 
160066976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1601d71ae5a4SJacob Faibussowitsch {
160239d35262SVincent Le Chenadec   PetscFunctionBegin;
160339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
160439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
160539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
160739d35262SVincent Le Chenadec }
160839d35262SVincent Le Chenadec 
160966976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2)
1610d71ae5a4SJacob Faibussowitsch {
161139d35262SVincent Le Chenadec   struct DMCompositeLink *next;
161239d35262SVincent Le Chenadec   PetscScalar            *array1, *array2;
161339d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite *)dm->data;
161439d35262SVincent Le Chenadec 
161539d35262SVincent Le Chenadec   PetscFunctionBegin;
161639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
161739d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2);
161839d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4);
161939d35262SVincent Le Chenadec 
162048a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
162139d35262SVincent Le Chenadec 
16229566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec1, &array1));
16239566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec2, &array2));
162439d35262SVincent Le Chenadec 
162515229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
162639d35262SVincent Le Chenadec   next = com->next;
162739d35262SVincent Le Chenadec   while (next) {
162839d35262SVincent Le Chenadec     Vec local1, local2;
162939d35262SVincent Le Chenadec 
16309566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local1));
16319566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local1, array1));
16329566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local2));
16339566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local2, array2));
16349566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2));
16359566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2));
16369566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local2));
16379566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local2));
16389566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local1));
16399566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local1));
164039d35262SVincent Le Chenadec 
164139d35262SVincent Le Chenadec     array1 += next->nlocal;
164239d35262SVincent Le Chenadec     array2 += next->nlocal;
164339d35262SVincent Le Chenadec     next = next->next;
164439d35262SVincent Le Chenadec   }
164539d35262SVincent Le Chenadec 
16469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec1, NULL));
16479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec2, NULL));
16483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
164939d35262SVincent Le Chenadec }
165039d35262SVincent Le Chenadec 
165166976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1652d71ae5a4SJacob Faibussowitsch {
165339d35262SVincent Le Chenadec   PetscFunctionBegin;
165439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
165539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
165639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16580c010503SBarry Smith }
165947c6ae99SBarry Smith 
16606ae3a549SBarry Smith /*MC
1661dce8aebaSBarry Smith    DMCOMPOSITE = "composite" - A `DM` object that is used to manage data for a collection of `DM`
16626ae3a549SBarry Smith 
16636ae3a549SBarry Smith   Level: intermediate
16646ae3a549SBarry Smith 
1665db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()`
16666ae3a549SBarry Smith M*/
16676ae3a549SBarry Smith 
1668d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1669d71ae5a4SJacob Faibussowitsch {
1670a4121054SBarry Smith   DM_Composite *com;
1671a4121054SBarry Smith 
1672a4121054SBarry Smith   PetscFunctionBegin;
16734dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&com));
1674a4121054SBarry Smith   p->data     = com;
1675a4121054SBarry Smith   com->n      = 0;
16767ac2b803SAlex Fikl   com->nghost = 0;
16770298fd71SBarry Smith   com->next   = NULL;
1678a4121054SBarry Smith   com->nDM    = 0;
1679a4121054SBarry Smith 
1680a4121054SBarry Smith   p->ops->createglobalvector       = DMCreateGlobalVector_Composite;
1681a4121054SBarry Smith   p->ops->createlocalvector        = DMCreateLocalVector_Composite;
1682184d77edSJed Brown   p->ops->getlocaltoglobalmapping  = DMGetLocalToGlobalMapping_Composite;
16834d343eeaSMatthew G Knepley   p->ops->createfieldis            = DMCreateFieldIS_Composite;
168416621825SDmitry Karpeev   p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite;
1685a4121054SBarry Smith   p->ops->refine                   = DMRefine_Composite;
168614354c39SJed Brown   p->ops->coarsen                  = DMCoarsen_Composite;
168725296bd5SBarry Smith   p->ops->createinterpolation      = DMCreateInterpolation_Composite;
168825296bd5SBarry Smith   p->ops->creatematrix             = DMCreateMatrix_Composite;
1689e727c939SJed Brown   p->ops->getcoloring              = DMCreateColoring_Composite;
1690a4121054SBarry Smith   p->ops->globaltolocalbegin       = DMGlobalToLocalBegin_Composite;
1691a4121054SBarry Smith   p->ops->globaltolocalend         = DMGlobalToLocalEnd_Composite;
169239d35262SVincent Le Chenadec   p->ops->localtoglobalbegin       = DMLocalToGlobalBegin_Composite;
169339d35262SVincent Le Chenadec   p->ops->localtoglobalend         = DMLocalToGlobalEnd_Composite;
169439d35262SVincent Le Chenadec   p->ops->localtolocalbegin        = DMLocalToLocalBegin_Composite;
169539d35262SVincent Le Chenadec   p->ops->localtolocalend          = DMLocalToLocalEnd_Composite;
1696a4121054SBarry Smith   p->ops->destroy                  = DMDestroy_Composite;
1697a4121054SBarry Smith   p->ops->view                     = DMView_Composite;
1698a4121054SBarry Smith   p->ops->setup                    = DMSetUp_Composite;
1699e10fd815SStefano Zampini 
17009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite));
17013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1702a4121054SBarry Smith }
1703a4121054SBarry Smith 
17041fd49c25SBarry Smith /*@
1705dce8aebaSBarry Smith   DMCompositeCreate - Creates a `DMCOMPOSITE`, used to generate "composite"
17060c010503SBarry Smith   vectors made up of several subvectors.
17070c010503SBarry Smith 
1708d083f849SBarry Smith   Collective
170947c6ae99SBarry Smith 
171047c6ae99SBarry Smith   Input Parameter:
17110c010503SBarry Smith . comm - the processors that will share the global vector
17120c010503SBarry Smith 
17132fe279fdSBarry Smith   Output Parameter:
1714dce8aebaSBarry Smith . packer - the `DMCOMPOSITE` object
171547c6ae99SBarry Smith 
171647c6ae99SBarry Smith   Level: advanced
171747c6ae99SBarry Smith 
171860225df5SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCreate()`
1719db781477SPatrick Sanan           `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`
1720db781477SPatrick Sanan           `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
172147c6ae99SBarry Smith @*/
1722d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer)
1723d71ae5a4SJacob Faibussowitsch {
172447c6ae99SBarry Smith   PetscFunctionBegin;
17254f572ea9SToby Isaac   PetscAssertPointer(packer, 2);
17269566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, packer));
17279566063dSJacob Faibussowitsch   PetscCall(DMSetType(*packer, DMCOMPOSITE));
17283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172947c6ae99SBarry Smith }
1730