xref: /petsc/src/dm/impls/composite/pack.c (revision f13dfd9ea68e0ddeee984e65c377a1819eab8a8a)
1ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h> /*I  "petscdmcomposite.h"  I*/
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>
3e10fd815SStefano Zampini #include <petsc/private/glvisviewerimpl.h>
42764a2aaSMatthew G. Knepley #include <petscds.h>
547c6ae99SBarry Smith 
647c6ae99SBarry Smith /*@C
747c6ae99SBarry Smith   DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
8dce8aebaSBarry Smith   separate components `DM` in a `DMCOMPOSITE` to build the correct matrix nonzero structure.
947c6ae99SBarry Smith 
1020f4b53cSBarry Smith   Logically Collective; No Fortran Support
1147c6ae99SBarry Smith 
12d8d19677SJose E. Roman   Input Parameters:
1347c6ae99SBarry Smith + dm                  - the composite object
1460225df5SJacob Faibussowitsch - FormCoupleLocations - routine to set the nonzero locations in the matrix
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   Level: advanced
1747c6ae99SBarry Smith 
18dce8aebaSBarry Smith   Note:
19dce8aebaSBarry Smith   See `DMSetApplicationContext()` and `DMGetApplicationContext()` for how to get user information into
2047c6ae99SBarry Smith   this routine
2147c6ae99SBarry Smith 
22dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`
2347c6ae99SBarry Smith @*/
24d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeSetCoupling(DM dm, PetscErrorCode (*FormCoupleLocations)(DM, Mat, PetscInt *, PetscInt *, PetscInt, PetscInt, PetscInt, PetscInt))
25d71ae5a4SJacob Faibussowitsch {
2647c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
2771b14b3eSStefano Zampini   PetscBool     flg;
2847c6ae99SBarry Smith 
2947c6ae99SBarry Smith   PetscFunctionBegin;
309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
317a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
3247c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3447c6ae99SBarry Smith }
3547c6ae99SBarry Smith 
3666976f2fSJacob Faibussowitsch static PetscErrorCode DMDestroy_Composite(DM dm)
37d71ae5a4SJacob Faibussowitsch {
3847c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
3947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
4047c6ae99SBarry Smith 
4147c6ae99SBarry Smith   PetscFunctionBegin;
4247c6ae99SBarry Smith   next = com->next;
4347c6ae99SBarry Smith   while (next) {
4447c6ae99SBarry Smith     prev = next;
4547c6ae99SBarry Smith     next = next->next;
469566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&prev->dm));
479566063dSJacob Faibussowitsch     PetscCall(PetscFree(prev->grstarts));
489566063dSJacob Faibussowitsch     PetscCall(PetscFree(prev));
4947c6ae99SBarry Smith   }
509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL));
51435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
529566063dSJacob Faibussowitsch   PetscCall(PetscFree(com));
533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5447c6ae99SBarry Smith }
5547c6ae99SBarry Smith 
5666976f2fSJacob Faibussowitsch static PetscErrorCode DMView_Composite(DM dm, PetscViewer v)
57d71ae5a4SJacob Faibussowitsch {
5847c6ae99SBarry Smith   PetscBool     iascii;
5947c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
6047c6ae99SBarry Smith 
6147c6ae99SBarry Smith   PetscFunctionBegin;
629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii));
6347c6ae99SBarry Smith   if (iascii) {
6447c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
6547c6ae99SBarry Smith     PetscInt                i;
6647c6ae99SBarry Smith 
679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "DM (%s)\n", ((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix"));
6863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "  contains %" PetscInt_FMT " DMs\n", com->nDM));
699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(v));
7047c6ae99SBarry Smith     for (i = 0; lnk; lnk = lnk->next, i++) {
7163a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(v, "Link %" PetscInt_FMT ": DM of type %s\n", i, ((PetscObject)lnk->dm)->type_name));
729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(v));
739566063dSJacob Faibussowitsch       PetscCall(DMView(lnk->dm, v));
749566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(v));
7547c6ae99SBarry Smith     }
769566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(v));
7747c6ae99SBarry Smith   }
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7947c6ae99SBarry Smith }
8047c6ae99SBarry Smith 
8147c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
8266976f2fSJacob Faibussowitsch static PetscErrorCode DMSetUp_Composite(DM dm)
83d71ae5a4SJacob Faibussowitsch {
8447c6ae99SBarry Smith   PetscInt                nprev = 0;
8547c6ae99SBarry Smith   PetscMPIInt             rank, size;
8647c6ae99SBarry Smith   DM_Composite           *com  = (DM_Composite *)dm->data;
8747c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
8847c6ae99SBarry Smith   PetscLayout             map;
8947c6ae99SBarry Smith 
9047c6ae99SBarry Smith   PetscFunctionBegin;
9128b400f6SJacob Faibussowitsch   PetscCheck(!com->setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Packer has already been setup");
929566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &map));
939566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(map, com->n));
949566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(map, PETSC_DETERMINE));
959566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(map, 1));
969566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(map));
979566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(map, &com->N));
989566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(map, &com->rstart, NULL));
999566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&map));
10047c6ae99SBarry Smith 
1019ae5db72SJed Brown   /* now set the rstart for each linked vector */
1029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
10447c6ae99SBarry Smith   while (next) {
10547c6ae99SBarry Smith     next->rstart = nprev;
10606ebdd98SJed Brown     nprev += next->n;
10747c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
1089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &next->grstarts));
1099566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgather(&next->grstart, 1, MPIU_INT, next->grstarts, 1, MPIU_INT, PetscObjectComm((PetscObject)dm)));
11047c6ae99SBarry Smith     next = next->next;
11147c6ae99SBarry Smith   }
11247c6ae99SBarry Smith   com->setup = PETSC_TRUE;
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11447c6ae99SBarry Smith }
11547c6ae99SBarry Smith 
11647c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
11747c6ae99SBarry Smith 
11873e31fe2SJed Brown /*@
119aaa8cc7dSPierre Jolivet   DMCompositeGetNumberDM - Gets the number of `DM` objects in the `DMCOMPOSITE`
12047c6ae99SBarry Smith   representation.
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith   Not Collective
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith   Input Parameter:
125dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith   Output Parameter:
128dce8aebaSBarry Smith . nDM - the number of `DM`
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith   Level: beginner
13147c6ae99SBarry Smith 
132dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`
13347c6ae99SBarry Smith @*/
134d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetNumberDM(DM dm, PetscInt *nDM)
135d71ae5a4SJacob Faibussowitsch {
13647c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
13771b14b3eSStefano Zampini   PetscBool     flg;
1385fd66863SKarl Rupp 
13947c6ae99SBarry Smith   PetscFunctionBegin;
14047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1419566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1427a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
14347c6ae99SBarry Smith   *nDM = com->nDM;
1443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14547c6ae99SBarry Smith }
14647c6ae99SBarry Smith 
14747c6ae99SBarry Smith /*@C
14847c6ae99SBarry Smith   DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
14947c6ae99SBarry Smith   representation.
15047c6ae99SBarry Smith 
15120f4b53cSBarry Smith   Collective
15247c6ae99SBarry Smith 
1539ae5db72SJed Brown   Input Parameters:
154dce8aebaSBarry Smith + dm   - the `DMCOMPOSITE` object
1559ae5db72SJed Brown - gvec - the global vector
1569ae5db72SJed Brown 
1572fe279fdSBarry Smith   Output Parameter:
158a4e35b19SJacob Faibussowitsch . ... - the packed parallel vectors, `NULL` for those that are not needed
15947c6ae99SBarry Smith 
16047c6ae99SBarry Smith   Level: advanced
16147c6ae99SBarry Smith 
162dce8aebaSBarry Smith   Note:
163dce8aebaSBarry Smith   Use `DMCompositeRestoreAccess()` to return the vectors when you no longer need them
164dce8aebaSBarry Smith 
16560225df5SJacob Faibussowitsch   Fortran Notes:
166dce8aebaSBarry Smith   Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
167dce8aebaSBarry Smith   or use the alternative interface `DMCompositeGetAccessArray()`.
168dce8aebaSBarry Smith 
169dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetEntries()`, `DMCompositeScatter()`
17047c6ae99SBarry Smith @*/
171d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccess(DM dm, Vec gvec, ...)
172d71ae5a4SJacob Faibussowitsch {
17347c6ae99SBarry Smith   va_list                 Argp;
17447c6ae99SBarry Smith   struct DMCompositeLink *next;
17547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
1765edff71fSBarry Smith   PetscInt                readonly;
17771b14b3eSStefano Zampini   PetscBool               flg;
17847c6ae99SBarry Smith 
17947c6ae99SBarry Smith   PetscFunctionBegin;
18047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
18147c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
1829566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1837a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
18447c6ae99SBarry Smith   next = com->next;
18548a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
18647c6ae99SBarry Smith 
1879566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec, &readonly));
18815229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
18947c6ae99SBarry Smith   va_start(Argp, gvec);
19047c6ae99SBarry Smith   while (next) {
19147c6ae99SBarry Smith     Vec *vec;
19247c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
1939ae5db72SJed Brown     if (vec) {
1949566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, vec));
1955edff71fSBarry Smith       if (readonly) {
1965edff71fSBarry Smith         const PetscScalar *array;
1979566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(gvec, &array));
1989566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec, array + next->rstart));
1999566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(*vec));
2009566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(gvec, &array));
2015edff71fSBarry Smith       } else {
2025edff71fSBarry Smith         PetscScalar *array;
2039566063dSJacob Faibussowitsch         PetscCall(VecGetArray(gvec, &array));
2049566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec, array + next->rstart));
2059566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(gvec, &array));
20647c6ae99SBarry Smith       }
2075edff71fSBarry Smith     }
20847c6ae99SBarry Smith     next = next->next;
20947c6ae99SBarry Smith   }
21047c6ae99SBarry Smith   va_end(Argp);
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21247c6ae99SBarry Smith }
21347c6ae99SBarry Smith 
214f73e5cebSJed Brown /*@C
215f73e5cebSJed Brown   DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
216f73e5cebSJed Brown   representation.
217f73e5cebSJed Brown 
21820f4b53cSBarry Smith   Collective
219f73e5cebSJed Brown 
220f73e5cebSJed Brown   Input Parameters:
221dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE`
222f73e5cebSJed Brown . pvec    - packed vector
223f73e5cebSJed Brown . nwanted - number of vectors wanted
224*f13dfd9eSBarry Smith - wanted  - sorted array of vectors wanted, or `NULL` to get all vectors, length `nwanted`
225f73e5cebSJed Brown 
2262fe279fdSBarry Smith   Output Parameter:
227*f13dfd9eSBarry Smith . vecs - array of requested global vectors (must be previously allocated and of length `nwanted`)
228f73e5cebSJed Brown 
229f73e5cebSJed Brown   Level: advanced
230f73e5cebSJed Brown 
231dce8aebaSBarry Smith   Note:
232dce8aebaSBarry Smith   Use `DMCompositeRestoreAccessArray()` to return the vectors when you no longer need them
233dce8aebaSBarry Smith 
234dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
235f73e5cebSJed Brown @*/
236d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs)
237d71ae5a4SJacob Faibussowitsch {
238f73e5cebSJed Brown   struct DMCompositeLink *link;
239f73e5cebSJed Brown   PetscInt                i, wnum;
240f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
241bee642f7SBarry Smith   PetscInt                readonly;
24271b14b3eSStefano Zampini   PetscBool               flg;
243f73e5cebSJed Brown 
244f73e5cebSJed Brown   PetscFunctionBegin;
245f73e5cebSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
246f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
2479566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
2487a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
24948a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
250f73e5cebSJed Brown 
2519566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
252f73e5cebSJed Brown   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
253f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
254f73e5cebSJed Brown       Vec v;
2559566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(link->dm, &v));
256bee642f7SBarry Smith       if (readonly) {
257bee642f7SBarry Smith         const PetscScalar *array;
2589566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec, &array));
2599566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + link->rstart));
2609566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
2619566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec, &array));
262bee642f7SBarry Smith       } else {
263bee642f7SBarry Smith         PetscScalar *array;
2649566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec, &array));
2659566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + link->rstart));
2669566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec, &array));
267bee642f7SBarry Smith       }
268f73e5cebSJed Brown       vecs[wnum++] = v;
269f73e5cebSJed Brown     }
270f73e5cebSJed Brown   }
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272f73e5cebSJed Brown }
273f73e5cebSJed Brown 
2747ac2b803SAlex Fikl /*@C
2757ac2b803SAlex Fikl   DMCompositeGetLocalAccessArray - Allows one to access the individual
2767ac2b803SAlex Fikl   packed vectors in their local representation.
2777ac2b803SAlex Fikl 
27820f4b53cSBarry Smith   Collective
2797ac2b803SAlex Fikl 
2807ac2b803SAlex Fikl   Input Parameters:
281dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE`
2827ac2b803SAlex Fikl . pvec    - packed vector
2837ac2b803SAlex Fikl . nwanted - number of vectors wanted
284*f13dfd9eSBarry Smith - wanted  - sorted array of vectors wanted, or `NULL` to get all vectors, length `nwanted`
2857ac2b803SAlex Fikl 
2862fe279fdSBarry Smith   Output Parameter:
287*f13dfd9eSBarry Smith . vecs - array of requested local vectors (must be allocated and of length `nwanted`)
2887ac2b803SAlex Fikl 
2897ac2b803SAlex Fikl   Level: advanced
2907ac2b803SAlex Fikl 
291dce8aebaSBarry Smith   Note:
292dce8aebaSBarry Smith   Use `DMCompositeRestoreLocalAccessArray()` to return the vectors
293dce8aebaSBarry Smith   when you no longer need them.
294dce8aebaSBarry Smith 
295dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`,
296db781477SPatrick Sanan           `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
2977ac2b803SAlex Fikl @*/
298d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs)
299d71ae5a4SJacob Faibussowitsch {
3007ac2b803SAlex Fikl   struct DMCompositeLink *link;
3017ac2b803SAlex Fikl   PetscInt                i, wnum;
3027ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite *)dm->data;
3037ac2b803SAlex Fikl   PetscInt                readonly;
3047ac2b803SAlex Fikl   PetscInt                nlocal = 0;
30571b14b3eSStefano Zampini   PetscBool               flg;
3067ac2b803SAlex Fikl 
3077ac2b803SAlex Fikl   PetscFunctionBegin;
3087ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3097ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
3109566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3117a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
31248a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
3137ac2b803SAlex Fikl 
3149566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
3157ac2b803SAlex Fikl   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
3167ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
3177ac2b803SAlex Fikl       Vec v;
3189566063dSJacob Faibussowitsch       PetscCall(DMGetLocalVector(link->dm, &v));
3197ac2b803SAlex Fikl       if (readonly) {
3207ac2b803SAlex Fikl         const PetscScalar *array;
3219566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec, &array));
3229566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + nlocal));
323b1c3483dSMark 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
3249566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
3259566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec, &array));
3267ac2b803SAlex Fikl       } else {
3277ac2b803SAlex Fikl         PetscScalar *array;
3289566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec, &array));
3299566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + nlocal));
3309566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec, &array));
3317ac2b803SAlex Fikl       }
3327ac2b803SAlex Fikl       vecs[wnum++] = v;
3337ac2b803SAlex Fikl     }
3347ac2b803SAlex Fikl 
3357ac2b803SAlex Fikl     nlocal += link->nlocal;
3367ac2b803SAlex Fikl   }
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3387ac2b803SAlex Fikl }
3397ac2b803SAlex Fikl 
34047c6ae99SBarry Smith /*@C
341dce8aebaSBarry Smith   DMCompositeRestoreAccess - Returns the vectors obtained with `DMCompositeGetAccess()`
34247c6ae99SBarry Smith   representation.
34347c6ae99SBarry Smith 
34420f4b53cSBarry Smith   Collective
34547c6ae99SBarry Smith 
3469ae5db72SJed Brown   Input Parameters:
347dce8aebaSBarry Smith + dm   - the `DMCOMPOSITE` object
34847c6ae99SBarry Smith . gvec - the global vector
349a4e35b19SJacob Faibussowitsch - ...  - the individual parallel vectors, `NULL` for those that are not needed
35047c6ae99SBarry Smith 
35147c6ae99SBarry Smith   Level: advanced
35247c6ae99SBarry Smith 
353dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
354db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`,
35542747ad1SJacob Faibussowitsch          `DMCompositeGetAccess()`
35647c6ae99SBarry Smith @*/
357d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...)
358d71ae5a4SJacob Faibussowitsch {
35947c6ae99SBarry Smith   va_list                 Argp;
36047c6ae99SBarry Smith   struct DMCompositeLink *next;
36147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
3625edff71fSBarry Smith   PetscInt                readonly;
36371b14b3eSStefano Zampini   PetscBool               flg;
36447c6ae99SBarry Smith 
36547c6ae99SBarry Smith   PetscFunctionBegin;
36647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36747c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
3689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3697a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
37047c6ae99SBarry Smith   next = com->next;
37148a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
37247c6ae99SBarry Smith 
3739566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec, &readonly));
37415229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
37547c6ae99SBarry Smith   va_start(Argp, gvec);
37647c6ae99SBarry Smith   while (next) {
37747c6ae99SBarry Smith     Vec *vec;
37847c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
3799ae5db72SJed Brown     if (vec) {
3809566063dSJacob Faibussowitsch       PetscCall(VecResetArray(*vec));
3811baa6e33SBarry Smith       if (readonly) PetscCall(VecLockReadPop(*vec));
3829566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, vec));
38347c6ae99SBarry Smith     }
38447c6ae99SBarry Smith     next = next->next;
38547c6ae99SBarry Smith   }
38647c6ae99SBarry Smith   va_end(Argp);
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38847c6ae99SBarry Smith }
38947c6ae99SBarry Smith 
390f73e5cebSJed Brown /*@C
391dce8aebaSBarry Smith   DMCompositeRestoreAccessArray - Returns the vectors obtained with `DMCompositeGetAccessArray()`
392f73e5cebSJed Brown 
39320f4b53cSBarry Smith   Collective
394f73e5cebSJed Brown 
395f73e5cebSJed Brown   Input Parameters:
396dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE` object
397f73e5cebSJed Brown . pvec    - packed vector
398f73e5cebSJed Brown . nwanted - number of vectors wanted
399*f13dfd9eSBarry Smith . wanted  - sorted array of vectors wanted, or `NULL` to restore all vectors
400*f13dfd9eSBarry Smith - vecs    - array of global vectors
401f73e5cebSJed Brown 
402f73e5cebSJed Brown   Level: advanced
403f73e5cebSJed Brown 
404dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
405f73e5cebSJed Brown @*/
406d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs)
407d71ae5a4SJacob Faibussowitsch {
408f73e5cebSJed Brown   struct DMCompositeLink *link;
409f73e5cebSJed Brown   PetscInt                i, wnum;
410f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
411bee642f7SBarry Smith   PetscInt                readonly;
41271b14b3eSStefano Zampini   PetscBool               flg;
413f73e5cebSJed Brown 
414f73e5cebSJed Brown   PetscFunctionBegin;
415f73e5cebSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
416f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4179566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4187a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
41948a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
420f73e5cebSJed Brown 
4219566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
422f73e5cebSJed Brown   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
423f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
4249566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
42548a46eb9SPierre Jolivet       if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4269566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum]));
427f73e5cebSJed Brown       wnum++;
428f73e5cebSJed Brown     }
429f73e5cebSJed Brown   }
4303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
431f73e5cebSJed Brown }
432f73e5cebSJed Brown 
4337ac2b803SAlex Fikl /*@C
434dce8aebaSBarry Smith   DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with `DMCompositeGetLocalAccessArray()`.
4357ac2b803SAlex Fikl 
43620f4b53cSBarry Smith   Collective
4377ac2b803SAlex Fikl 
4387ac2b803SAlex Fikl   Input Parameters:
439dce8aebaSBarry Smith + dm      - the `DMCOMPOSITE` object
4407ac2b803SAlex Fikl . pvec    - packed vector
4417ac2b803SAlex Fikl . nwanted - number of vectors wanted
442*f13dfd9eSBarry Smith . wanted  - sorted array of vectors wanted, or `NULL` to restore all vectors
443*f13dfd9eSBarry Smith - vecs    - array of local vectors
4447ac2b803SAlex Fikl 
4457ac2b803SAlex Fikl   Level: advanced
4467ac2b803SAlex Fikl 
447dce8aebaSBarry Smith   Note:
448*f13dfd9eSBarry Smith   `nwanted` and `wanted` must match the values given to `DMCompositeGetLocalAccessArray()`
4497ac2b803SAlex Fikl   otherwise the call will fail.
4507ac2b803SAlex Fikl 
451dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`,
452db781477SPatrick Sanan           `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`,
453db781477SPatrick Sanan           `DMCompositeScatter()`, `DMCompositeGather()`
4547ac2b803SAlex Fikl @*/
455d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs)
456d71ae5a4SJacob Faibussowitsch {
4577ac2b803SAlex Fikl   struct DMCompositeLink *link;
4587ac2b803SAlex Fikl   PetscInt                i, wnum;
4597ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite *)dm->data;
4607ac2b803SAlex Fikl   PetscInt                readonly;
46171b14b3eSStefano Zampini   PetscBool               flg;
4627ac2b803SAlex Fikl 
4637ac2b803SAlex Fikl   PetscFunctionBegin;
4647ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4657ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4669566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4677a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
46848a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
4697ac2b803SAlex Fikl 
4709566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
4717ac2b803SAlex Fikl   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
4727ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
4739566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
47448a46eb9SPierre Jolivet       if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4759566063dSJacob Faibussowitsch       PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum]));
4767ac2b803SAlex Fikl       wnum++;
4777ac2b803SAlex Fikl     }
4787ac2b803SAlex Fikl   }
4793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4807ac2b803SAlex Fikl }
4817ac2b803SAlex Fikl 
48247c6ae99SBarry Smith /*@C
48347c6ae99SBarry Smith   DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
48447c6ae99SBarry Smith 
48520f4b53cSBarry Smith   Collective
48647c6ae99SBarry Smith 
4879ae5db72SJed Brown   Input Parameters:
488dce8aebaSBarry Smith + dm   - the `DMCOMPOSITE` object
48947c6ae99SBarry Smith . gvec - the global vector
490a4e35b19SJacob Faibussowitsch - ...  - the individual sequential vectors, `NULL` for those that are not needed
49147c6ae99SBarry Smith 
49247c6ae99SBarry Smith   Level: advanced
49347c6ae99SBarry Smith 
494dce8aebaSBarry Smith   Note:
495dce8aebaSBarry Smith   `DMCompositeScatterArray()` is a non-variadic alternative that is often more convenient for library callers and is
4966f3c3dcfSJed Brown   accessible from Fortran.
4976f3c3dcfSJed Brown 
498dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
499db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
500db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
501db781477SPatrick Sanan          `DMCompositeScatterArray()`
50247c6ae99SBarry Smith @*/
503d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...)
504d71ae5a4SJacob Faibussowitsch {
50547c6ae99SBarry Smith   va_list                 Argp;
50647c6ae99SBarry Smith   struct DMCompositeLink *next;
507c599c493SJunchao Zhang   PETSC_UNUSED PetscInt   cnt;
50847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
50971b14b3eSStefano Zampini   PetscBool               flg;
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith   PetscFunctionBegin;
51247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
51347c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5157a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
51648a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
51747c6ae99SBarry Smith 
51815229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
51947c6ae99SBarry Smith   va_start(Argp, gvec);
5208fd8f222SJed Brown   for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
5219ae5db72SJed Brown     Vec local;
5229ae5db72SJed Brown     local = va_arg(Argp, Vec);
5239ae5db72SJed Brown     if (local) {
5249ae5db72SJed Brown       Vec                global;
5255edff71fSBarry Smith       const PetscScalar *array;
52663a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
5279566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
5289566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
5299566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
5309566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local));
5319566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local));
5329566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
5339566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5349566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
53547c6ae99SBarry Smith     }
53647c6ae99SBarry Smith   }
53747c6ae99SBarry Smith   va_end(Argp);
5383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53947c6ae99SBarry Smith }
54047c6ae99SBarry Smith 
5416f3c3dcfSJed Brown /*@
5426f3c3dcfSJed Brown   DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5436f3c3dcfSJed Brown 
54420f4b53cSBarry Smith   Collective
5456f3c3dcfSJed Brown 
5466f3c3dcfSJed Brown   Input Parameters:
547dce8aebaSBarry Smith + dm    - the `DMCOMPOSITE` object
5486f3c3dcfSJed Brown . gvec  - the global vector
549a2b725a8SWilliam Gropp - lvecs - array of local vectors, NULL for any that are not needed
5506f3c3dcfSJed Brown 
5516f3c3dcfSJed Brown   Level: advanced
5526f3c3dcfSJed Brown 
5536f3c3dcfSJed Brown   Note:
554dce8aebaSBarry Smith   This is a non-variadic alternative to `DMCompositeScatter()`
5556f3c3dcfSJed Brown 
556dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`
557db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
558db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
5596f3c3dcfSJed Brown @*/
560d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs)
561d71ae5a4SJacob Faibussowitsch {
5626f3c3dcfSJed Brown   struct DMCompositeLink *next;
5636f3c3dcfSJed Brown   PetscInt                i;
5646f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
56571b14b3eSStefano Zampini   PetscBool               flg;
5666f3c3dcfSJed Brown 
5676f3c3dcfSJed Brown   PetscFunctionBegin;
5686f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5696f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5717a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
57248a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
5736f3c3dcfSJed Brown 
57415229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
5756f3c3dcfSJed Brown   for (i = 0, next = com->next; next; next = next->next, i++) {
5766f3c3dcfSJed Brown     if (lvecs[i]) {
5776f3c3dcfSJed Brown       Vec                global;
578c5d31e75SLisandro Dalcin       const PetscScalar *array;
5796f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3);
5809566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
5819566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
5829566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart));
5839566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i]));
5849566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i]));
5859566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
5869566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5879566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
5886f3c3dcfSJed Brown     }
5896f3c3dcfSJed Brown   }
5903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5916f3c3dcfSJed Brown }
5926f3c3dcfSJed Brown 
59347c6ae99SBarry Smith /*@C
59447c6ae99SBarry Smith   DMCompositeGather - Gathers into a global packed vector from its individual local vectors
59547c6ae99SBarry Smith 
59620f4b53cSBarry Smith   Collective
59747c6ae99SBarry Smith 
598d8d19677SJose E. Roman   Input Parameters:
599dce8aebaSBarry Smith + dm    - the `DMCOMPOSITE` object
600dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES`
601a4e35b19SJacob Faibussowitsch . gvec  - the global vector
602a4e35b19SJacob Faibussowitsch - ...   - the individual sequential vectors, `NULL` for any that are not needed
60347c6ae99SBarry Smith 
60447c6ae99SBarry Smith   Level: advanced
60547c6ae99SBarry Smith 
60660225df5SJacob Faibussowitsch   Fortran Notes:
607dce8aebaSBarry Smith   Fortran users should use `DMCompositeGatherArray()`
608f5f57ec0SBarry Smith 
609dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
610db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
611db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
61247c6ae99SBarry Smith @*/
613d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...)
614d71ae5a4SJacob Faibussowitsch {
61547c6ae99SBarry Smith   va_list                 Argp;
61647c6ae99SBarry Smith   struct DMCompositeLink *next;
61747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
618c599c493SJunchao Zhang   PETSC_UNUSED PetscInt   cnt;
61971b14b3eSStefano Zampini   PetscBool               flg;
62047c6ae99SBarry Smith 
62147c6ae99SBarry Smith   PetscFunctionBegin;
62247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
623064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6249566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6257a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
62648a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
62747c6ae99SBarry Smith 
62815229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
6291dac896bSSatish Balay   va_start(Argp, gvec);
6308fd8f222SJed Brown   for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
6319ae5db72SJed Brown     Vec local;
6329ae5db72SJed Brown     local = va_arg(Argp, Vec);
6339ae5db72SJed Brown     if (local) {
63447c6ae99SBarry Smith       PetscScalar *array;
6359ae5db72SJed Brown       Vec          global;
63663a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
6379566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
6389566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec, &array));
6399566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
6409566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global));
6419566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global));
6429566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec, &array));
6439566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6449566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
64547c6ae99SBarry Smith     }
64647c6ae99SBarry Smith   }
64747c6ae99SBarry Smith   va_end(Argp);
6483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64947c6ae99SBarry Smith }
65047c6ae99SBarry Smith 
6516f3c3dcfSJed Brown /*@
6526f3c3dcfSJed Brown   DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6536f3c3dcfSJed Brown 
65420f4b53cSBarry Smith   Collective
6556f3c3dcfSJed Brown 
656d8d19677SJose E. Roman   Input Parameters:
657dce8aebaSBarry Smith + dm    - the `DMCOMPOSITE` object
6586f3c3dcfSJed Brown . gvec  - the global vector
659dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES`
6606f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed
6616f3c3dcfSJed Brown 
6626f3c3dcfSJed Brown   Level: advanced
6636f3c3dcfSJed Brown 
664dce8aebaSBarry Smith   Note:
665dce8aebaSBarry Smith   This is a non-variadic alternative to `DMCompositeGather()`.
6666f3c3dcfSJed Brown 
667dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
668db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
669db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`,
6706f3c3dcfSJed Brown @*/
671d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs)
672d71ae5a4SJacob Faibussowitsch {
6736f3c3dcfSJed Brown   struct DMCompositeLink *next;
6746f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
6756f3c3dcfSJed Brown   PetscInt                i;
67671b14b3eSStefano Zampini   PetscBool               flg;
6776f3c3dcfSJed Brown 
6786f3c3dcfSJed Brown   PetscFunctionBegin;
6796f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
680064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6819566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6827a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
68348a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
6846f3c3dcfSJed Brown 
68515229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
6866f3c3dcfSJed Brown   for (next = com->next, i = 0; next; next = next->next, i++) {
6876f3c3dcfSJed Brown     if (lvecs[i]) {
6886f3c3dcfSJed Brown       PetscScalar *array;
6896f3c3dcfSJed Brown       Vec          global;
690064a246eSJacob Faibussowitsch       PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4);
6919566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
6929566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec, &array));
6939566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
6949566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global));
6959566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global));
6969566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec, &array));
6979566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6989566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
6996f3c3dcfSJed Brown     }
7006f3c3dcfSJed Brown   }
7013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7026f3c3dcfSJed Brown }
7036f3c3dcfSJed Brown 
704f5f57ec0SBarry Smith /*@
705dce8aebaSBarry Smith   DMCompositeAddDM - adds a `DM` vector to a `DMCOMPOSITE`
70647c6ae99SBarry Smith 
70720f4b53cSBarry Smith   Collective
70847c6ae99SBarry Smith 
709d8d19677SJose E. Roman   Input Parameters:
710dce8aebaSBarry Smith + dmc - the  `DMCOMPOSITE` object
711dce8aebaSBarry Smith - dm  - the `DM` object
71247c6ae99SBarry Smith 
71347c6ae99SBarry Smith   Level: advanced
71447c6ae99SBarry Smith 
71542747ad1SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeGather()`, `DMCreateGlobalVector()`,
716db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
717db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
71847c6ae99SBarry Smith @*/
719d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeAddDM(DM dmc, DM dm)
720d71ae5a4SJacob Faibussowitsch {
72106ebdd98SJed Brown   PetscInt                n, nlocal;
72247c6ae99SBarry Smith   struct DMCompositeLink *mine, *next;
72306ebdd98SJed Brown   Vec                     global, local;
72447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dmc->data;
72571b14b3eSStefano Zampini   PetscBool               flg;
72647c6ae99SBarry Smith 
72747c6ae99SBarry Smith   PetscFunctionBegin;
72847c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc, DM_CLASSID, 1);
72947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
7309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg));
7317a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
73247c6ae99SBarry Smith   next = com->next;
73328b400f6SJacob Faibussowitsch   PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite");
73447c6ae99SBarry Smith 
73547c6ae99SBarry Smith   /* create new link */
7369566063dSJacob Faibussowitsch   PetscCall(PetscNew(&mine));
7379566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
7389566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &global));
7399566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global, &n));
7409566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &global));
7419566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &local));
7429566063dSJacob Faibussowitsch   PetscCall(VecGetSize(local, &nlocal));
7439566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &local));
7448865f1eaSKarl Rupp 
74547c6ae99SBarry Smith   mine->n      = n;
74606ebdd98SJed Brown   mine->nlocal = nlocal;
74747c6ae99SBarry Smith   mine->dm     = dm;
7480298fd71SBarry Smith   mine->next   = NULL;
74947c6ae99SBarry Smith   com->n += n;
7507ac2b803SAlex Fikl   com->nghost += nlocal;
75147c6ae99SBarry Smith 
75247c6ae99SBarry Smith   /* add to end of list */
7538865f1eaSKarl Rupp   if (!next) com->next = mine;
7548865f1eaSKarl Rupp   else {
75547c6ae99SBarry Smith     while (next->next) next = next->next;
75647c6ae99SBarry Smith     next->next = mine;
75747c6ae99SBarry Smith   }
75847c6ae99SBarry Smith   com->nDM++;
75947c6ae99SBarry Smith   com->nmine++;
7603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76147c6ae99SBarry Smith }
76247c6ae99SBarry Smith 
7639804daf3SBarry Smith #include <petscdraw.h>
76426887b52SJed Brown PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
76566976f2fSJacob Faibussowitsch static PetscErrorCode       VecView_DMComposite(Vec gvec, PetscViewer viewer)
766d71ae5a4SJacob Faibussowitsch {
76747c6ae99SBarry Smith   DM                      dm;
76847c6ae99SBarry Smith   struct DMCompositeLink *next;
76947c6ae99SBarry Smith   PetscBool               isdraw;
770cef07954SSatish Balay   DM_Composite           *com;
77147c6ae99SBarry Smith 
77247c6ae99SBarry Smith   PetscFunctionBegin;
7739566063dSJacob Faibussowitsch   PetscCall(VecGetDM(gvec, &dm));
7747a8be351SBarry Smith   PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite");
77547c6ae99SBarry Smith   com  = (DM_Composite *)dm->data;
77647c6ae99SBarry Smith   next = com->next;
77747c6ae99SBarry Smith 
7789566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
77947c6ae99SBarry Smith   if (!isdraw) {
78047c6ae99SBarry Smith     /* do I really want to call this? */
7819566063dSJacob Faibussowitsch     PetscCall(VecView_MPI(gvec, viewer));
78247c6ae99SBarry Smith   } else {
78347c6ae99SBarry Smith     PetscInt cnt = 0;
78447c6ae99SBarry Smith 
78515229ffcSPierre Jolivet     /* loop over packed objects, handling one at a time */
78647c6ae99SBarry Smith     while (next) {
78747c6ae99SBarry Smith       Vec                vec;
788ca4278abSLisandro Dalcin       const PetscScalar *array;
78947c6ae99SBarry Smith       PetscInt           bs;
79047c6ae99SBarry Smith 
7919ae5db72SJed Brown       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
7929566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &vec));
7939566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
7949566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart));
7959566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
7969566063dSJacob Faibussowitsch       PetscCall(VecView(vec, viewer));
7979566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vec));
7989566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(vec, &bs));
7999566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &vec));
8009566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawBaseAdd(viewer, bs));
80147c6ae99SBarry Smith       cnt += bs;
80247c6ae99SBarry Smith       next = next->next;
80347c6ae99SBarry Smith     }
8049566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt));
80547c6ae99SBarry Smith   }
8063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80747c6ae99SBarry Smith }
80847c6ae99SBarry Smith 
80966976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec)
810d71ae5a4SJacob Faibussowitsch {
81147c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
81247c6ae99SBarry Smith 
81347c6ae99SBarry Smith   PetscFunctionBegin;
81447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8159566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
8169566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8179566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec));
8189566063dSJacob Faibussowitsch   PetscCall(VecSetType(*gvec, dm->vectype));
8199566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*gvec, com->n, com->N));
8209566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*gvec, dm));
8219566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (void (*)(void))VecView_DMComposite));
8223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82347c6ae99SBarry Smith }
82447c6ae99SBarry Smith 
82566976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec)
826d71ae5a4SJacob Faibussowitsch {
82747c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
82847c6ae99SBarry Smith 
82947c6ae99SBarry Smith   PetscFunctionBegin;
83047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
83147c6ae99SBarry Smith   if (!com->setup) {
8329566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(dm));
8339566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
83447c6ae99SBarry Smith   }
8359566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, lvec));
8369566063dSJacob Faibussowitsch   PetscCall(VecSetType(*lvec, dm->vectype));
8379566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE));
8389566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*lvec, dm));
8393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
84047c6ae99SBarry Smith }
84147c6ae99SBarry Smith 
84247c6ae99SBarry Smith /*@C
843dce8aebaSBarry Smith   DMCompositeGetISLocalToGlobalMappings - gets an `ISLocalToGlobalMapping` for each `DM` in the `DMCOMPOSITE`, maps to the composite global space
84447c6ae99SBarry Smith 
84520f4b53cSBarry Smith   Collective; No Fortran Support
84647c6ae99SBarry Smith 
84747c6ae99SBarry Smith   Input Parameter:
848dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
84947c6ae99SBarry Smith 
8502fe279fdSBarry Smith   Output Parameter:
8519ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes
852dce8aebaSBarry Smith            all the ghost points that individual ghosted `DMDA` may have.
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith   Level: advanced
85547c6ae99SBarry Smith 
856dce8aebaSBarry Smith   Note:
857*f13dfd9eSBarry Smith   Each entry of `ltogs` should be destroyed with `ISLocalToGlobalMappingDestroy()`, `ltogs` should be freed with `PetscFree()`.
85847c6ae99SBarry Smith 
859dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
860db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
861c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
86247c6ae99SBarry Smith @*/
863d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping **ltogs)
864d71ae5a4SJacob Faibussowitsch {
86547c6ae99SBarry Smith   PetscInt                i, *idx, n, cnt;
86647c6ae99SBarry Smith   struct DMCompositeLink *next;
86747c6ae99SBarry Smith   PetscMPIInt             rank;
86847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
86971b14b3eSStefano Zampini   PetscBool               flg;
87047c6ae99SBarry Smith 
87147c6ae99SBarry Smith   PetscFunctionBegin;
87247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
8747a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
8759566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM, ltogs));
87747c6ae99SBarry Smith   next = com->next;
8789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
87947c6ae99SBarry Smith 
88015229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
88147c6ae99SBarry Smith   cnt = 0;
88247c6ae99SBarry Smith   while (next) {
8836eb61c8cSJed Brown     ISLocalToGlobalMapping ltog;
8846eb61c8cSJed Brown     PetscMPIInt            size;
88586994e45SJed Brown     const PetscInt        *suboff, *indices;
8866eb61c8cSJed Brown     Vec                    global;
88747c6ae99SBarry Smith 
8886eb61c8cSJed Brown     /* Get sub-DM global indices for each local dof */
8899566063dSJacob Faibussowitsch     PetscCall(DMGetLocalToGlobalMapping(next->dm, &ltog));
8909566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n));
8919566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices));
8929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &idx));
89347c6ae99SBarry Smith 
8946eb61c8cSJed Brown     /* Get the offsets for the sub-DM global vector */
8959566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
8969566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRanges(global, &suboff));
8979566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size));
8986eb61c8cSJed Brown 
8996eb61c8cSJed Brown     /* Shift the sub-DM definition of the global space to the composite global space */
9006eb61c8cSJed Brown     for (i = 0; i < n; i++) {
90186994e45SJed Brown       PetscInt subi = indices[i], lo = 0, hi = size, t;
902e333b035SStefano Zampini       /* There's no consensus on what a negative index means,
903e333b035SStefano Zampini          except for skipping when setting the values in vectors and matrices */
9049371c9d4SSatish Balay       if (subi < 0) {
9059371c9d4SSatish Balay         idx[i] = subi - next->grstarts[rank];
9069371c9d4SSatish Balay         continue;
9079371c9d4SSatish Balay       }
9086eb61c8cSJed Brown       /* Binary search to find which rank owns subi */
9096eb61c8cSJed Brown       while (hi - lo > 1) {
9106eb61c8cSJed Brown         t = lo + (hi - lo) / 2;
9116eb61c8cSJed Brown         if (suboff[t] > subi) hi = t;
9126eb61c8cSJed Brown         else lo = t;
9136eb61c8cSJed Brown       }
9146eb61c8cSJed Brown       idx[i] = subi - suboff[lo] + next->grstarts[lo];
9156eb61c8cSJed Brown     }
9169566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices));
9179566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt]));
9189566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
91947c6ae99SBarry Smith     next = next->next;
92047c6ae99SBarry Smith     cnt++;
92147c6ae99SBarry Smith   }
9223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92347c6ae99SBarry Smith }
92447c6ae99SBarry Smith 
92587c85e80SJed Brown /*@C
9269ae5db72SJed Brown   DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
92787c85e80SJed Brown 
92820f4b53cSBarry Smith   Not Collective; No Fortran Support
92987c85e80SJed Brown 
9304165533cSJose E. Roman   Input Parameter:
931dce8aebaSBarry Smith . dm - the `DMCOMPOSITE`
93287c85e80SJed Brown 
9334165533cSJose E. Roman   Output Parameter:
93415229ffcSPierre Jolivet . is - array of serial index sets for each component of the `DMCOMPOSITE`
93587c85e80SJed Brown 
93687c85e80SJed Brown   Level: intermediate
93787c85e80SJed Brown 
93887c85e80SJed Brown   Notes:
93987c85e80SJed Brown   At present, a composite local vector does not normally exist.  This function is used to provide index sets for
94015229ffcSPierre Jolivet   `MatGetLocalSubMatrix()`.  In the future, the scatters for each entry in the `DMCOMPOSITE` may be merged into a single
9419ae5db72SJed Brown   scatter to a composite local vector.  The user should not typically need to know which is being done.
94287c85e80SJed Brown 
943dce8aebaSBarry Smith   To get the composite global indices at all local points (including ghosts), use `DMCompositeGetISLocalToGlobalMappings()`.
94487c85e80SJed Brown 
945dce8aebaSBarry Smith   To get index sets for pieces of the composite global vector, use `DMCompositeGetGlobalISs()`.
94687c85e80SJed Brown 
947dce8aebaSBarry Smith   Each returned `IS` should be destroyed with `ISDestroy()`, the array should be freed with `PetscFree()`.
94887c85e80SJed Brown 
949*f13dfd9eSBarry Smith   Fortran Note:
950*f13dfd9eSBarry Smith   Pass in an array long enough to hold all the `IS`, see `DMCompositeGetNumberDM()`
951*f13dfd9eSBarry Smith 
952*f13dfd9eSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`,
953*f13dfd9eSBarry Smith           `MatCreateLocalRef()`, `DMCompositeGetNumberDM()`
95487c85e80SJed Brown @*/
955d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalISs(DM dm, IS **is)
956d71ae5a4SJacob Faibussowitsch {
95787c85e80SJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
95887c85e80SJed Brown   struct DMCompositeLink *link;
95987c85e80SJed Brown   PetscInt                cnt, start;
96071b14b3eSStefano Zampini   PetscBool               flg;
96187c85e80SJed Brown 
96287c85e80SJed Brown   PetscFunctionBegin;
96387c85e80SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9644f572ea9SToby Isaac   PetscAssertPointer(is, 2);
9659566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
9667a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
9679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nmine, is));
96806ebdd98SJed Brown   for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) {
969520db06cSJed Brown     PetscInt bs;
9709566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt]));
9719566063dSJacob Faibussowitsch     PetscCall(DMGetBlockSize(link->dm, &bs));
9729566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize((*is)[cnt], bs));
973520db06cSJed Brown   }
9743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97587c85e80SJed Brown }
97687c85e80SJed Brown 
97747c6ae99SBarry Smith /*@C
978dce8aebaSBarry Smith   DMCompositeGetGlobalISs - Gets the index sets for each composed object in a `DMCOMPOSITE`
97947c6ae99SBarry Smith 
98020f4b53cSBarry Smith   Collective
98147c6ae99SBarry Smith 
98247c6ae99SBarry Smith   Input Parameter:
983dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
98447c6ae99SBarry Smith 
9852fe279fdSBarry Smith   Output Parameter:
98647c6ae99SBarry Smith . is - the array of index sets
98747c6ae99SBarry Smith 
98847c6ae99SBarry Smith   Level: advanced
98947c6ae99SBarry Smith 
99047c6ae99SBarry Smith   Notes:
991*f13dfd9eSBarry Smith   The `is` entries should be destroyed with `ISDestroy()`, `is` should be freed with `PetscFree()`
99247c6ae99SBarry Smith 
99347c6ae99SBarry Smith   These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
99447c6ae99SBarry Smith 
995dce8aebaSBarry Smith   Use `DMCompositeGetLocalISs()` for index sets in the packed local numbering, and
996dce8aebaSBarry Smith   `DMCompositeGetISLocalToGlobalMappings()` for to map local sub-`DM` (including ghost) indices to packed global
9976eb61c8cSJed Brown   indices.
99847c6ae99SBarry Smith 
99960225df5SJacob Faibussowitsch   Fortran Notes:
1000dce8aebaSBarry Smith   The output argument 'is' must be an allocated array of sufficient length, which can be learned using `DMCompositeGetNumberDM()`.
1001f3cb0f7eSJed Brown 
1002dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1003db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
1004c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
100547c6ae99SBarry Smith @*/
1006d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[])
1007d71ae5a4SJacob Faibussowitsch {
100866bb578eSMark Adams   PetscInt                cnt = 0;
100947c6ae99SBarry Smith   struct DMCompositeLink *next;
101047c6ae99SBarry Smith   PetscMPIInt             rank;
101147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
101271b14b3eSStefano Zampini   PetscBool               flg;
101347c6ae99SBarry Smith 
101447c6ae99SBarry Smith   PetscFunctionBegin;
101547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10169566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
10177a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
10187a8be351SBarry Smith   PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before");
10199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM, is));
102047c6ae99SBarry Smith   next = com->next;
10219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
102247c6ae99SBarry Smith 
102315229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
102447c6ae99SBarry Smith   while (next) {
1025e5e52638SMatthew G. Knepley     PetscDS prob;
1026e5e52638SMatthew G. Knepley 
10279566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt]));
10289566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &prob));
1029e5e52638SMatthew G. Knepley     if (prob) {
103065c226d8SMatthew G. Knepley       MatNullSpace space;
103165c226d8SMatthew G. Knepley       Mat          pmat;
10320f21e855SMatthew G. Knepley       PetscObject  disc;
10330f21e855SMatthew G. Knepley       PetscInt     Nf;
103465c226d8SMatthew G. Knepley 
10359566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(prob, &Nf));
1036f24dd8d2SMatthew G. Knepley       if (cnt < Nf) {
10379566063dSJacob Faibussowitsch         PetscCall(PetscDSGetDiscretization(prob, cnt, &disc));
10389566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space));
10399566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space));
10409566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space));
10419566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space));
10429566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat));
10439566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat));
104465c226d8SMatthew G. Knepley       }
1045f24dd8d2SMatthew G. Knepley     }
104647c6ae99SBarry Smith     cnt++;
104747c6ae99SBarry Smith     next = next->next;
104847c6ae99SBarry Smith   }
10493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105047c6ae99SBarry Smith }
105147c6ae99SBarry Smith 
105266976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
1053d71ae5a4SJacob Faibussowitsch {
10544d343eeaSMatthew G Knepley   PetscInt nDM;
10554d343eeaSMatthew G Knepley   DM      *dms;
10564d343eeaSMatthew G Knepley   PetscInt i;
10574d343eeaSMatthew G Knepley 
10584d343eeaSMatthew G Knepley   PetscFunctionBegin;
10599566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm, &nDM));
10608865f1eaSKarl Rupp   if (numFields) *numFields = nDM;
10619566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetGlobalISs(dm, fields));
10624d343eeaSMatthew G Knepley   if (fieldNames) {
10639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, &dms));
10649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, fieldNames));
10659566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, dms));
10664d343eeaSMatthew G Knepley     for (i = 0; i < nDM; i++) {
10674d343eeaSMatthew G Knepley       char        buf[256];
10684d343eeaSMatthew G Knepley       const char *splitname;
10694d343eeaSMatthew G Knepley 
10704d343eeaSMatthew G Knepley       /* Split naming precedence: object name, prefix, number */
10714d343eeaSMatthew G Knepley       splitname = ((PetscObject)dm)->name;
10724d343eeaSMatthew G Knepley       if (!splitname) {
10739566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname));
10744d343eeaSMatthew G Knepley         if (splitname) {
10754d343eeaSMatthew G Knepley           size_t len;
10769566063dSJacob Faibussowitsch           PetscCall(PetscStrncpy(buf, splitname, sizeof(buf)));
10778caf3d72SBarry Smith           buf[sizeof(buf) - 1] = 0;
10789566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(buf, &len));
10794d343eeaSMatthew G Knepley           if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */
10804d343eeaSMatthew G Knepley           splitname = buf;
10814d343eeaSMatthew G Knepley         }
10824d343eeaSMatthew G Knepley       }
10834d343eeaSMatthew G Knepley       if (!splitname) {
108463a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i));
10854d343eeaSMatthew G Knepley         splitname = buf;
10864d343eeaSMatthew G Knepley       }
10879566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i]));
10884d343eeaSMatthew G Knepley     }
10899566063dSJacob Faibussowitsch     PetscCall(PetscFree(dms));
10904d343eeaSMatthew G Knepley   }
10913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10924d343eeaSMatthew G Knepley }
10934d343eeaSMatthew G Knepley 
1094e7c4fc90SDmitry Karpeev /*
1095e7c4fc90SDmitry Karpeev  This could take over from DMCreateFieldIS(), as it is more general,
10960298fd71SBarry Smith  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1097e7c4fc90SDmitry Karpeev  At this point it's probably best to be less intrusive, however.
1098e7c4fc90SDmitry Karpeev  */
109966976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
1100d71ae5a4SJacob Faibussowitsch {
1101e7c4fc90SDmitry Karpeev   PetscInt nDM;
1102e7c4fc90SDmitry Karpeev   PetscInt i;
1103e7c4fc90SDmitry Karpeev 
1104e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
11059566063dSJacob Faibussowitsch   PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist));
1106e7c4fc90SDmitry Karpeev   if (dmlist) {
11079566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetNumberDM(dm, &nDM));
11089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, dmlist));
11099566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, *dmlist));
111048a46eb9SPierre Jolivet     for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i])));
1111e7c4fc90SDmitry Karpeev   }
11123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1113e7c4fc90SDmitry Karpeev }
1114e7c4fc90SDmitry Karpeev 
111547c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
111647c6ae99SBarry Smith /*@C
1117dce8aebaSBarry Smith   DMCompositeGetLocalVectors - Gets local vectors for each part of a `DMCOMPOSITE`
1118dce8aebaSBarry Smith   Use `DMCompositeRestoreLocalVectors()` to return them.
111947c6ae99SBarry Smith 
112020f4b53cSBarry Smith   Not Collective; No Fortran Support
112147c6ae99SBarry Smith 
112247c6ae99SBarry Smith   Input Parameter:
1123dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
112447c6ae99SBarry Smith 
112547c6ae99SBarry Smith   Output Parameter:
1126a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`s
112747c6ae99SBarry Smith 
112847c6ae99SBarry Smith   Level: advanced
112947c6ae99SBarry Smith 
1130dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1131db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1132db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
113347c6ae99SBarry Smith @*/
1134d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...)
1135d71ae5a4SJacob Faibussowitsch {
113647c6ae99SBarry Smith   va_list                 Argp;
113747c6ae99SBarry Smith   struct DMCompositeLink *next;
113847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
113971b14b3eSStefano Zampini   PetscBool               flg;
114047c6ae99SBarry Smith 
114147c6ae99SBarry Smith   PetscFunctionBegin;
114247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11447a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
114547c6ae99SBarry Smith   next = com->next;
114615229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
114747c6ae99SBarry Smith   va_start(Argp, dm);
114847c6ae99SBarry Smith   while (next) {
114947c6ae99SBarry Smith     Vec *vec;
115047c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
11519566063dSJacob Faibussowitsch     if (vec) PetscCall(DMGetLocalVector(next->dm, vec));
115247c6ae99SBarry Smith     next = next->next;
115347c6ae99SBarry Smith   }
115447c6ae99SBarry Smith   va_end(Argp);
11553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115647c6ae99SBarry Smith }
115747c6ae99SBarry Smith 
115847c6ae99SBarry Smith /*@C
1159dce8aebaSBarry Smith   DMCompositeRestoreLocalVectors - Restores local vectors for each part of a `DMCOMPOSITE`
116047c6ae99SBarry Smith 
116120f4b53cSBarry Smith   Not Collective; No Fortran Support
116247c6ae99SBarry Smith 
116347c6ae99SBarry Smith   Input Parameter:
1164dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
116547c6ae99SBarry Smith 
116647c6ae99SBarry Smith   Output Parameter:
1167a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`
116847c6ae99SBarry Smith 
116947c6ae99SBarry Smith   Level: advanced
117047c6ae99SBarry Smith 
1171dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1172db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1173db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
117447c6ae99SBarry Smith @*/
1175d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...)
1176d71ae5a4SJacob Faibussowitsch {
117747c6ae99SBarry Smith   va_list                 Argp;
117847c6ae99SBarry Smith   struct DMCompositeLink *next;
117947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
118071b14b3eSStefano Zampini   PetscBool               flg;
118147c6ae99SBarry Smith 
118247c6ae99SBarry Smith   PetscFunctionBegin;
118347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11857a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
118647c6ae99SBarry Smith   next = com->next;
118715229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
118847c6ae99SBarry Smith   va_start(Argp, dm);
118947c6ae99SBarry Smith   while (next) {
119047c6ae99SBarry Smith     Vec *vec;
119147c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
11929566063dSJacob Faibussowitsch     if (vec) PetscCall(DMRestoreLocalVector(next->dm, vec));
119347c6ae99SBarry Smith     next = next->next;
119447c6ae99SBarry Smith   }
119547c6ae99SBarry Smith   va_end(Argp);
11963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119747c6ae99SBarry Smith }
119847c6ae99SBarry Smith 
119947c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
120047c6ae99SBarry Smith /*@C
1201dce8aebaSBarry Smith   DMCompositeGetEntries - Gets the `DM` for each entry in a `DMCOMPOSITE`.
120247c6ae99SBarry Smith 
120347c6ae99SBarry Smith   Not Collective
120447c6ae99SBarry Smith 
120547c6ae99SBarry Smith   Input Parameter:
1206dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
120747c6ae99SBarry Smith 
120847c6ae99SBarry Smith   Output Parameter:
1209a4e35b19SJacob Faibussowitsch . ... - the individual entries `DM`
121047c6ae99SBarry Smith 
121147c6ae99SBarry Smith   Level: advanced
121247c6ae99SBarry Smith 
121360225df5SJacob Faibussowitsch   Fortran Notes:
1214dce8aebaSBarry Smith   Available as `DMCompositeGetEntries()` for one output `DM`, DMCompositeGetEntries2() for 2, etc
1215f5f57ec0SBarry Smith 
1216dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()`
1217db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
121860225df5SJacob Faibussowitsch          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`
121947c6ae99SBarry Smith @*/
1220d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntries(DM dm, ...)
1221d71ae5a4SJacob Faibussowitsch {
122247c6ae99SBarry Smith   va_list                 Argp;
122347c6ae99SBarry Smith   struct DMCompositeLink *next;
122447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
122571b14b3eSStefano Zampini   PetscBool               flg;
122647c6ae99SBarry Smith 
122747c6ae99SBarry Smith   PetscFunctionBegin;
122847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12299566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
12307a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
123147c6ae99SBarry Smith   next = com->next;
123215229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
123347c6ae99SBarry Smith   va_start(Argp, dm);
123447c6ae99SBarry Smith   while (next) {
123547c6ae99SBarry Smith     DM *dmn;
123647c6ae99SBarry Smith     dmn = va_arg(Argp, DM *);
12379ae5db72SJed Brown     if (dmn) *dmn = next->dm;
123847c6ae99SBarry Smith     next = next->next;
123947c6ae99SBarry Smith   }
124047c6ae99SBarry Smith   va_end(Argp);
12413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124247c6ae99SBarry Smith }
124347c6ae99SBarry Smith 
1244dbab29e1SMark F. Adams /*@C
1245dce8aebaSBarry Smith   DMCompositeGetEntriesArray - Gets the DM for each entry in a `DMCOMPOSITE`
12462fa5ba8aSJed Brown 
12472fa5ba8aSJed Brown   Not Collective
12482fa5ba8aSJed Brown 
12492fa5ba8aSJed Brown   Input Parameter:
1250dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
1251907376e6SBarry Smith 
1252907376e6SBarry Smith   Output Parameter:
1253dce8aebaSBarry Smith . dms - array of sufficient length (see `DMCompositeGetNumberDM()`) to hold the individual `DM`
12542fa5ba8aSJed Brown 
12552fa5ba8aSJed Brown   Level: advanced
12562fa5ba8aSJed Brown 
1257dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()`
1258db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
125960225df5SJacob Faibussowitsch          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`
12602fa5ba8aSJed Brown @*/
1261d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntriesArray(DM dm, DM dms[])
1262d71ae5a4SJacob Faibussowitsch {
12632fa5ba8aSJed Brown   struct DMCompositeLink *next;
12642fa5ba8aSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
12652fa5ba8aSJed Brown   PetscInt                i;
126671b14b3eSStefano Zampini   PetscBool               flg;
12672fa5ba8aSJed Brown 
12682fa5ba8aSJed Brown   PetscFunctionBegin;
12692fa5ba8aSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
12717a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
127215229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
12732fa5ba8aSJed Brown   for (next = com->next, i = 0; next; next = next->next, i++) dms[i] = next->dm;
12743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12752fa5ba8aSJed Brown }
12762fa5ba8aSJed Brown 
1277e10fd815SStefano Zampini typedef struct {
1278e10fd815SStefano Zampini   DM           dm;
1279e10fd815SStefano Zampini   PetscViewer *subv;
1280e10fd815SStefano Zampini   Vec         *vecs;
1281e10fd815SStefano Zampini } GLVisViewerCtx;
1282e10fd815SStefano Zampini 
1283d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx)
1284d71ae5a4SJacob Faibussowitsch {
1285e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
1286e10fd815SStefano Zampini   PetscInt        i, n;
1287e10fd815SStefano Zampini 
1288e10fd815SStefano Zampini   PetscFunctionBegin;
12899566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(ctx->dm, &n));
129048a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscViewerDestroy(&ctx->subv[i]));
12919566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ctx->subv, ctx->vecs));
12929566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dm));
12939566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
12943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1295e10fd815SStefano Zampini }
1296e10fd815SStefano Zampini 
1297d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
1298d71ae5a4SJacob Faibussowitsch {
1299e10fd815SStefano Zampini   Vec             X   = (Vec)oX;
1300e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
1301e10fd815SStefano Zampini   PetscInt        i, n, cumf;
1302e10fd815SStefano Zampini 
1303e10fd815SStefano Zampini   PetscFunctionBegin;
13049566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(ctx->dm, &n));
13059566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccessArray(ctx->dm, X, n, NULL, ctx->vecs));
1306e10fd815SStefano Zampini   for (i = 0, cumf = 0; i < n; i++) {
1307e10fd815SStefano Zampini     PetscErrorCode (*g2l)(PetscObject, PetscInt, PetscObject[], void *);
1308e10fd815SStefano Zampini     void    *fctx;
1309e10fd815SStefano Zampini     PetscInt nfi;
1310e10fd815SStefano Zampini 
131134e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nfi, NULL, NULL, &g2l, NULL, &fctx));
1312e10fd815SStefano Zampini     if (!nfi) continue;
13131baa6e33SBarry Smith     if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i], nfi, oXfield + cumf, fctx));
1314f4f49eeaSPierre Jolivet     else PetscCall(VecCopy(ctx->vecs[i], (Vec)oXfield[cumf]));
1315e10fd815SStefano Zampini     cumf += nfi;
1316e10fd815SStefano Zampini   }
13179566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccessArray(ctx->dm, X, n, NULL, ctx->vecs));
13183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1319e10fd815SStefano Zampini }
1320e10fd815SStefano Zampini 
1321d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer)
1322d71ae5a4SJacob Faibussowitsch {
1323e10fd815SStefano Zampini   DM              dm = (DM)odm, *dms;
1324e10fd815SStefano Zampini   Vec            *Ufds;
1325e10fd815SStefano Zampini   GLVisViewerCtx *ctx;
1326e10fd815SStefano Zampini   PetscInt        i, n, tnf, *sdim;
1327e10fd815SStefano Zampini   char          **fecs;
1328e10fd815SStefano Zampini 
1329e10fd815SStefano Zampini   PetscFunctionBegin;
13309566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
13319566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
1332e10fd815SStefano Zampini   ctx->dm = dm;
13339566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm, &n));
13349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &dms));
13359566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetEntriesArray(dm, dms));
13369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n, &ctx->subv, n, &ctx->vecs));
1337e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1338e10fd815SStefano Zampini     PetscInt nf;
1339e10fd815SStefano Zampini 
13409566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PetscObjectComm(odm), &ctx->subv[i]));
13419566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(ctx->subv[i], PETSCVIEWERGLVIS));
134234e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisSetDM_Internal(ctx->subv[i], (PetscObject)dms[i]));
134334e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, NULL, NULL, NULL, NULL, NULL));
1344e10fd815SStefano Zampini     tnf += nf;
1345e10fd815SStefano Zampini   }
13469566063dSJacob Faibussowitsch   PetscCall(PetscFree(dms));
13479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(tnf, &fecs, tnf, &sdim, tnf, &Ufds));
1348e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1349e10fd815SStefano Zampini     PetscInt    *sd, nf, f;
1350e10fd815SStefano Zampini     const char **fec;
1351e10fd815SStefano Zampini     Vec         *Uf;
1352e10fd815SStefano Zampini 
135334e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, &fec, &sd, NULL, (PetscObject **)&Uf, NULL));
1354e10fd815SStefano Zampini     for (f = 0; f < nf; f++) {
13559566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(fec[f], &fecs[tnf + f]));
1356e10fd815SStefano Zampini       Ufds[tnf + f] = Uf[f];
1357e10fd815SStefano Zampini       sdim[tnf + f] = sd[f];
1358e10fd815SStefano Zampini     }
1359e10fd815SStefano Zampini     tnf += nf;
1360e10fd815SStefano Zampini   }
13619566063dSJacob Faibussowitsch   PetscCall(PetscViewerGLVisSetFields(viewer, tnf, (const char **)fecs, sdim, DMCompositeSampleGLVisFields_Private, (PetscObject *)Ufds, ctx, DestroyGLVisViewerCtx_Private));
136248a46eb9SPierre Jolivet   for (i = 0; i < tnf; i++) PetscCall(PetscFree(fecs[i]));
13639566063dSJacob Faibussowitsch   PetscCall(PetscFree3(fecs, sdim, Ufds));
13643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1365e10fd815SStefano Zampini }
1366e10fd815SStefano Zampini 
136766976f2fSJacob Faibussowitsch static PetscErrorCode DMRefine_Composite(DM dmi, MPI_Comm comm, DM *fine)
1368d71ae5a4SJacob Faibussowitsch {
136947c6ae99SBarry Smith   struct DMCompositeLink *next;
137047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dmi->data;
137147c6ae99SBarry Smith   DM                      dm;
137247c6ae99SBarry Smith 
137347c6ae99SBarry Smith   PetscFunctionBegin;
137447c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi, DM_CLASSID, 1);
137548a46eb9SPierre Jolivet   if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm));
13769566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmi));
137747c6ae99SBarry Smith   next = com->next;
13789566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm, fine));
137947c6ae99SBarry Smith 
138015229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
138147c6ae99SBarry Smith   while (next) {
13829566063dSJacob Faibussowitsch     PetscCall(DMRefine(next->dm, comm, &dm));
13839566063dSJacob Faibussowitsch     PetscCall(DMCompositeAddDM(*fine, dm));
13849566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dm));
138547c6ae99SBarry Smith     next = next->next;
138647c6ae99SBarry Smith   }
13873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138847c6ae99SBarry Smith }
138947c6ae99SBarry Smith 
139066976f2fSJacob Faibussowitsch static PetscErrorCode DMCoarsen_Composite(DM dmi, MPI_Comm comm, DM *fine)
1391d71ae5a4SJacob Faibussowitsch {
139214354c39SJed Brown   struct DMCompositeLink *next;
139314354c39SJed Brown   DM_Composite           *com = (DM_Composite *)dmi->data;
139414354c39SJed Brown   DM                      dm;
139514354c39SJed Brown 
139614354c39SJed Brown   PetscFunctionBegin;
139714354c39SJed Brown   PetscValidHeaderSpecific(dmi, DM_CLASSID, 1);
13989566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmi));
139948a46eb9SPierre Jolivet   if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm));
140014354c39SJed Brown   next = com->next;
14019566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm, fine));
140214354c39SJed Brown 
140315229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
140414354c39SJed Brown   while (next) {
14059566063dSJacob Faibussowitsch     PetscCall(DMCoarsen(next->dm, comm, &dm));
14069566063dSJacob Faibussowitsch     PetscCall(DMCompositeAddDM(*fine, dm));
14079566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dm));
140814354c39SJed Brown     next = next->next;
140914354c39SJed Brown   }
14103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
141114354c39SJed Brown }
141247c6ae99SBarry Smith 
141366976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_Composite(DM coarse, DM fine, Mat *A, Vec *v)
1414d71ae5a4SJacob Faibussowitsch {
14159ae5db72SJed Brown   PetscInt                m, n, M, N, nDM, i;
141647c6ae99SBarry Smith   struct DMCompositeLink *nextc;
141747c6ae99SBarry Smith   struct DMCompositeLink *nextf;
141825296bd5SBarry Smith   Vec                     gcoarse, gfine, *vecs;
141947c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite *)coarse->data;
142047c6ae99SBarry Smith   DM_Composite           *comfine   = (DM_Composite *)fine->data;
14219ae5db72SJed Brown   Mat                    *mats;
142247c6ae99SBarry Smith 
142347c6ae99SBarry Smith   PetscFunctionBegin;
142447c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse, DM_CLASSID, 1);
142547c6ae99SBarry Smith   PetscValidHeaderSpecific(fine, DM_CLASSID, 2);
14269566063dSJacob Faibussowitsch   PetscCall(DMSetUp(coarse));
14279566063dSJacob Faibussowitsch   PetscCall(DMSetUp(fine));
142847c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14299566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(coarse, &gcoarse));
14309566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(fine, &gfine));
14319566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gcoarse, &n));
14329566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gfine, &m));
14339566063dSJacob Faibussowitsch   PetscCall(VecGetSize(gcoarse, &N));
14349566063dSJacob Faibussowitsch   PetscCall(VecGetSize(gfine, &M));
14359566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(coarse, &gcoarse));
14369566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(fine, &gfine));
143747c6ae99SBarry Smith 
14389ae5db72SJed Brown   nDM = comfine->nDM;
143963a3b9bcSJacob 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);
14409566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nDM * nDM, &mats));
144148a46eb9SPierre Jolivet   if (v) PetscCall(PetscCalloc1(nDM, &vecs));
144247c6ae99SBarry Smith 
144315229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
14449ae5db72SJed Brown   for (nextc = comcoarse->next, nextf = comfine->next, i = 0; nextc; nextc = nextc->next, nextf = nextf->next, i++) {
14451baa6e33SBarry Smith     if (!v) PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], NULL));
14461baa6e33SBarry Smith     else PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], &vecs[i]));
144747c6ae99SBarry Smith   }
14489566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, nDM, NULL, mats, A));
14491baa6e33SBarry Smith   if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, vecs, v));
14509566063dSJacob Faibussowitsch   for (i = 0; i < nDM * nDM; i++) PetscCall(MatDestroy(&mats[i]));
14519566063dSJacob Faibussowitsch   PetscCall(PetscFree(mats));
145225296bd5SBarry Smith   if (v) {
14539566063dSJacob Faibussowitsch     for (i = 0; i < nDM; i++) PetscCall(VecDestroy(&vecs[i]));
14549566063dSJacob Faibussowitsch     PetscCall(PetscFree(vecs));
145525296bd5SBarry Smith   }
14563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
145747c6ae99SBarry Smith }
145847c6ae99SBarry Smith 
1459d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
1460d71ae5a4SJacob Faibussowitsch {
14611411c6eeSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
14621411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1463f7efa3c7SJed Brown   PetscInt                i;
14641411c6eeSJed Brown 
14651411c6eeSJed Brown   PetscFunctionBegin;
14661411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
14679566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetISLocalToGlobalMappings(dm, &ltogs));
14689566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm), com->nDM, ltogs, &dm->ltogmap));
14699566063dSJacob Faibussowitsch   for (i = 0; i < com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(&ltogs[i]));
14709566063dSJacob Faibussowitsch   PetscCall(PetscFree(ltogs));
14713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14721411c6eeSJed Brown }
14731411c6eeSJed Brown 
147466976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateColoring_Composite(DM dm, ISColoringType ctype, ISColoring *coloring)
1475d71ae5a4SJacob Faibussowitsch {
147647c6ae99SBarry Smith   PetscInt         n, i, cnt;
147747c6ae99SBarry Smith   ISColoringValue *colors;
147847c6ae99SBarry Smith   PetscBool        dense  = PETSC_FALSE;
147947c6ae99SBarry Smith   ISColoringValue  maxcol = 0;
148047c6ae99SBarry Smith   DM_Composite    *com    = (DM_Composite *)dm->data;
148147c6ae99SBarry Smith 
148247c6ae99SBarry Smith   PetscFunctionBegin;
148347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
148408401ef6SPierre Jolivet   PetscCheck(ctype != IS_COLORING_LOCAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only global coloring supported");
1485f7d195e4SLawrence Mitchell   if (ctype == IS_COLORING_GLOBAL) {
148647c6ae99SBarry Smith     n = com->n;
1487ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unknown ISColoringType");
14889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &colors)); /* freed in ISColoringDestroy() */
148947c6ae99SBarry Smith 
14909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL));
149147c6ae99SBarry Smith   if (dense) {
1492ad540459SPierre Jolivet     for (i = 0; i < n; i++) colors[i] = (ISColoringValue)(com->rstart + i);
149347c6ae99SBarry Smith     maxcol = com->N;
149447c6ae99SBarry Smith   } else {
149547c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
149647c6ae99SBarry Smith     PetscMPIInt             rank;
149747c6ae99SBarry Smith 
14989566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
149947c6ae99SBarry Smith     cnt = 0;
150047c6ae99SBarry Smith     while (next) {
150147c6ae99SBarry Smith       ISColoring lcoloring;
150247c6ae99SBarry Smith 
15039566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(next->dm, IS_COLORING_GLOBAL, &lcoloring));
1504ad540459SPierre Jolivet       for (i = 0; i < lcoloring->N; i++) colors[cnt++] = maxcol + lcoloring->colors[i];
150547c6ae99SBarry Smith       maxcol += lcoloring->n;
15069566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&lcoloring));
150747c6ae99SBarry Smith       next = next->next;
150847c6ae99SBarry Smith     }
150947c6ae99SBarry Smith   }
15109566063dSJacob Faibussowitsch   PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), maxcol, n, colors, PETSC_OWN_POINTER, coloring));
15113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
151247c6ae99SBarry Smith }
151347c6ae99SBarry Smith 
151466976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec)
1515d71ae5a4SJacob Faibussowitsch {
151647c6ae99SBarry Smith   struct DMCompositeLink *next;
151747c6ae99SBarry Smith   PetscScalar            *garray, *larray;
151847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
151947c6ae99SBarry Smith 
152047c6ae99SBarry Smith   PetscFunctionBegin;
152147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
152247c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
152339d35262SVincent Le Chenadec 
152448a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
152539d35262SVincent Le Chenadec 
15269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gvec, &garray));
15279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lvec, &larray));
152847c6ae99SBarry Smith 
152915229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
153039d35262SVincent Le Chenadec   next = com->next;
153147c6ae99SBarry Smith   while (next) {
153247c6ae99SBarry Smith     Vec      local, global;
153347c6ae99SBarry Smith     PetscInt N;
153447c6ae99SBarry Smith 
15359566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
15369566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(global, &N));
15379566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(global, garray));
15389566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local));
15399566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local, larray));
15409566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(next->dm, global, mode, local));
15419566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(next->dm, global, mode, local));
15429566063dSJacob Faibussowitsch     PetscCall(VecResetArray(global));
15439566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local));
15449566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
15459566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local));
154639d35262SVincent Le Chenadec 
154706ebdd98SJed Brown     larray += next->nlocal;
154839d35262SVincent Le Chenadec     garray += next->n;
154947c6ae99SBarry Smith     next = next->next;
155047c6ae99SBarry Smith   }
155147c6ae99SBarry Smith 
15529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gvec, NULL));
15539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lvec, NULL));
15543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
155547c6ae99SBarry Smith }
155647c6ae99SBarry Smith 
155766976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec)
1558d71ae5a4SJacob Faibussowitsch {
15590c010503SBarry Smith   PetscFunctionBegin;
156039d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
156139d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
156239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 4);
15633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
156439d35262SVincent Le Chenadec }
156539d35262SVincent Le Chenadec 
156666976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1567d71ae5a4SJacob Faibussowitsch {
156839d35262SVincent Le Chenadec   struct DMCompositeLink *next;
156939d35262SVincent Le Chenadec   PetscScalar            *larray, *garray;
157039d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite *)dm->data;
157139d35262SVincent Le Chenadec 
157239d35262SVincent Le Chenadec   PetscFunctionBegin;
157339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
157439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
157539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
157639d35262SVincent Le Chenadec 
157748a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
157839d35262SVincent Le Chenadec 
15799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lvec, &larray));
15809566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gvec, &garray));
158139d35262SVincent Le Chenadec 
158215229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
158339d35262SVincent Le Chenadec   next = com->next;
158439d35262SVincent Le Chenadec   while (next) {
158539d35262SVincent Le Chenadec     Vec global, local;
158639d35262SVincent Le Chenadec 
15879566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local));
15889566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local, larray));
15899566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
15909566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(global, garray));
15919566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(next->dm, local, mode, global));
15929566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(next->dm, local, mode, global));
15939566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local));
15949566063dSJacob Faibussowitsch     PetscCall(VecResetArray(global));
15959566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
15969566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local));
159739d35262SVincent Le Chenadec 
159839d35262SVincent Le Chenadec     garray += next->n;
159939d35262SVincent Le Chenadec     larray += next->nlocal;
160039d35262SVincent Le Chenadec     next = next->next;
160139d35262SVincent Le Chenadec   }
160239d35262SVincent Le Chenadec 
16039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gvec, NULL));
16049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lvec, NULL));
16053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
160639d35262SVincent Le Chenadec }
160739d35262SVincent Le Chenadec 
160866976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1609d71ae5a4SJacob Faibussowitsch {
161039d35262SVincent Le Chenadec   PetscFunctionBegin;
161139d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
161239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
161339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
161539d35262SVincent Le Chenadec }
161639d35262SVincent Le Chenadec 
161766976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2)
1618d71ae5a4SJacob Faibussowitsch {
161939d35262SVincent Le Chenadec   struct DMCompositeLink *next;
162039d35262SVincent Le Chenadec   PetscScalar            *array1, *array2;
162139d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite *)dm->data;
162239d35262SVincent Le Chenadec 
162339d35262SVincent Le Chenadec   PetscFunctionBegin;
162439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
162539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2);
162639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4);
162739d35262SVincent Le Chenadec 
162848a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
162939d35262SVincent Le Chenadec 
16309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec1, &array1));
16319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec2, &array2));
163239d35262SVincent Le Chenadec 
163315229ffcSPierre Jolivet   /* loop over packed objects, handling one at a time */
163439d35262SVincent Le Chenadec   next = com->next;
163539d35262SVincent Le Chenadec   while (next) {
163639d35262SVincent Le Chenadec     Vec local1, local2;
163739d35262SVincent Le Chenadec 
16389566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local1));
16399566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local1, array1));
16409566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local2));
16419566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local2, array2));
16429566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2));
16439566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2));
16449566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local2));
16459566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local2));
16469566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local1));
16479566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local1));
164839d35262SVincent Le Chenadec 
164939d35262SVincent Le Chenadec     array1 += next->nlocal;
165039d35262SVincent Le Chenadec     array2 += next->nlocal;
165139d35262SVincent Le Chenadec     next = next->next;
165239d35262SVincent Le Chenadec   }
165339d35262SVincent Le Chenadec 
16549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec1, NULL));
16559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec2, NULL));
16563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165739d35262SVincent Le Chenadec }
165839d35262SVincent Le Chenadec 
165966976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1660d71ae5a4SJacob Faibussowitsch {
166139d35262SVincent Le Chenadec   PetscFunctionBegin;
166239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
166339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
166439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16660c010503SBarry Smith }
166747c6ae99SBarry Smith 
16686ae3a549SBarry Smith /*MC
1669dce8aebaSBarry Smith    DMCOMPOSITE = "composite" - A `DM` object that is used to manage data for a collection of `DM`
16706ae3a549SBarry Smith 
16716ae3a549SBarry Smith   Level: intermediate
16726ae3a549SBarry Smith 
1673db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()`
16746ae3a549SBarry Smith M*/
16756ae3a549SBarry Smith 
1676d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1677d71ae5a4SJacob Faibussowitsch {
1678a4121054SBarry Smith   DM_Composite *com;
1679a4121054SBarry Smith 
1680a4121054SBarry Smith   PetscFunctionBegin;
16814dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&com));
1682a4121054SBarry Smith   p->data     = com;
1683a4121054SBarry Smith   com->n      = 0;
16847ac2b803SAlex Fikl   com->nghost = 0;
16850298fd71SBarry Smith   com->next   = NULL;
1686a4121054SBarry Smith   com->nDM    = 0;
1687a4121054SBarry Smith 
1688a4121054SBarry Smith   p->ops->createglobalvector       = DMCreateGlobalVector_Composite;
1689a4121054SBarry Smith   p->ops->createlocalvector        = DMCreateLocalVector_Composite;
1690184d77edSJed Brown   p->ops->getlocaltoglobalmapping  = DMGetLocalToGlobalMapping_Composite;
16914d343eeaSMatthew G Knepley   p->ops->createfieldis            = DMCreateFieldIS_Composite;
169216621825SDmitry Karpeev   p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite;
1693a4121054SBarry Smith   p->ops->refine                   = DMRefine_Composite;
169414354c39SJed Brown   p->ops->coarsen                  = DMCoarsen_Composite;
169525296bd5SBarry Smith   p->ops->createinterpolation      = DMCreateInterpolation_Composite;
169625296bd5SBarry Smith   p->ops->creatematrix             = DMCreateMatrix_Composite;
1697e727c939SJed Brown   p->ops->getcoloring              = DMCreateColoring_Composite;
1698a4121054SBarry Smith   p->ops->globaltolocalbegin       = DMGlobalToLocalBegin_Composite;
1699a4121054SBarry Smith   p->ops->globaltolocalend         = DMGlobalToLocalEnd_Composite;
170039d35262SVincent Le Chenadec   p->ops->localtoglobalbegin       = DMLocalToGlobalBegin_Composite;
170139d35262SVincent Le Chenadec   p->ops->localtoglobalend         = DMLocalToGlobalEnd_Composite;
170239d35262SVincent Le Chenadec   p->ops->localtolocalbegin        = DMLocalToLocalBegin_Composite;
170339d35262SVincent Le Chenadec   p->ops->localtolocalend          = DMLocalToLocalEnd_Composite;
1704a4121054SBarry Smith   p->ops->destroy                  = DMDestroy_Composite;
1705a4121054SBarry Smith   p->ops->view                     = DMView_Composite;
1706a4121054SBarry Smith   p->ops->setup                    = DMSetUp_Composite;
1707e10fd815SStefano Zampini 
17089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite));
17093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1710a4121054SBarry Smith }
1711a4121054SBarry Smith 
17121fd49c25SBarry Smith /*@
1713dce8aebaSBarry Smith   DMCompositeCreate - Creates a `DMCOMPOSITE`, used to generate "composite"
17140c010503SBarry Smith   vectors made up of several subvectors.
17150c010503SBarry Smith 
1716d083f849SBarry Smith   Collective
171747c6ae99SBarry Smith 
171847c6ae99SBarry Smith   Input Parameter:
17190c010503SBarry Smith . comm - the processors that will share the global vector
17200c010503SBarry Smith 
17212fe279fdSBarry Smith   Output Parameter:
1722dce8aebaSBarry Smith . packer - the `DMCOMPOSITE` object
172347c6ae99SBarry Smith 
172447c6ae99SBarry Smith   Level: advanced
172547c6ae99SBarry Smith 
172660225df5SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCreate()`
1727db781477SPatrick Sanan           `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`
1728db781477SPatrick Sanan           `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
172947c6ae99SBarry Smith @*/
1730d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer)
1731d71ae5a4SJacob Faibussowitsch {
173247c6ae99SBarry Smith   PetscFunctionBegin;
17334f572ea9SToby Isaac   PetscAssertPointer(packer, 2);
17349566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, packer));
17359566063dSJacob Faibussowitsch   PetscCall(DMSetType(*packer, DMCOMPOSITE));
17363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173747c6ae99SBarry Smith }
1738