xref: /petsc/src/dm/impls/composite/pack.c (revision aaa8cc7d2a5c3913edcbb923e20f154fe9c4aa65)
147c6ae99SBarry Smith 
2ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h> /*I  "petscdmcomposite.h"  I*/
3af0996ceSBarry Smith #include <petsc/private/isimpl.h>
4e10fd815SStefano Zampini #include <petsc/private/glvisviewerimpl.h>
52764a2aaSMatthew G. Knepley #include <petscds.h>
647c6ae99SBarry Smith 
747c6ae99SBarry Smith /*@C
847c6ae99SBarry Smith     DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
9dce8aebaSBarry Smith       separate components `DM` in a `DMCOMPOSITE` to build the correct matrix nonzero structure.
1047c6ae99SBarry Smith 
1120f4b53cSBarry Smith     Logically Collective; No Fortran Support
1247c6ae99SBarry Smith 
13d8d19677SJose E. Roman     Input Parameters:
1447c6ae99SBarry Smith +   dm - the composite object
1547c6ae99SBarry Smith -   formcouplelocations - routine to set the nonzero locations in the matrix
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith     Level: advanced
1847c6ae99SBarry Smith 
19dce8aebaSBarry Smith     Note:
20dce8aebaSBarry Smith     See `DMSetApplicationContext()` and `DMGetApplicationContext()` for how to get user information into
2147c6ae99SBarry Smith     this routine
2247c6ae99SBarry Smith 
23dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`
2447c6ae99SBarry Smith @*/
25d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeSetCoupling(DM dm, PetscErrorCode (*FormCoupleLocations)(DM, Mat, PetscInt *, PetscInt *, PetscInt, PetscInt, PetscInt, PetscInt))
26d71ae5a4SJacob Faibussowitsch {
2747c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
2871b14b3eSStefano Zampini   PetscBool     flg;
2947c6ae99SBarry Smith 
3047c6ae99SBarry Smith   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
327a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
3347c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3547c6ae99SBarry Smith }
3647c6ae99SBarry Smith 
37d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDestroy_Composite(DM dm)
38d71ae5a4SJacob Faibussowitsch {
3947c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
4047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
4147c6ae99SBarry Smith 
4247c6ae99SBarry Smith   PetscFunctionBegin;
4347c6ae99SBarry Smith   next = com->next;
4447c6ae99SBarry Smith   while (next) {
4547c6ae99SBarry Smith     prev = next;
4647c6ae99SBarry Smith     next = next->next;
479566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&prev->dm));
489566063dSJacob Faibussowitsch     PetscCall(PetscFree(prev->grstarts));
499566063dSJacob Faibussowitsch     PetscCall(PetscFree(prev));
5047c6ae99SBarry Smith   }
519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL));
52435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
539566063dSJacob Faibussowitsch   PetscCall(PetscFree(com));
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5547c6ae99SBarry Smith }
5647c6ae99SBarry Smith 
57d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_Composite(DM dm, PetscViewer v)
58d71ae5a4SJacob Faibussowitsch {
5947c6ae99SBarry Smith   PetscBool     iascii;
6047c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
6147c6ae99SBarry Smith 
6247c6ae99SBarry Smith   PetscFunctionBegin;
639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii));
6447c6ae99SBarry Smith   if (iascii) {
6547c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
6647c6ae99SBarry Smith     PetscInt                i;
6747c6ae99SBarry Smith 
689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "DM (%s)\n", ((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix"));
6963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "  contains %" PetscInt_FMT " DMs\n", com->nDM));
709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(v));
7147c6ae99SBarry Smith     for (i = 0; lnk; lnk = lnk->next, i++) {
7263a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(v, "Link %" PetscInt_FMT ": DM of type %s\n", i, ((PetscObject)lnk->dm)->type_name));
739566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(v));
749566063dSJacob Faibussowitsch       PetscCall(DMView(lnk->dm, v));
759566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(v));
7647c6ae99SBarry Smith     }
779566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(v));
7847c6ae99SBarry Smith   }
793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8047c6ae99SBarry Smith }
8147c6ae99SBarry Smith 
8247c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
83d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUp_Composite(DM dm)
84d71ae5a4SJacob Faibussowitsch {
8547c6ae99SBarry Smith   PetscInt                nprev = 0;
8647c6ae99SBarry Smith   PetscMPIInt             rank, size;
8747c6ae99SBarry Smith   DM_Composite           *com  = (DM_Composite *)dm->data;
8847c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
8947c6ae99SBarry Smith   PetscLayout             map;
9047c6ae99SBarry Smith 
9147c6ae99SBarry Smith   PetscFunctionBegin;
9228b400f6SJacob Faibussowitsch   PetscCheck(!com->setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Packer has already been setup");
939566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &map));
949566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(map, com->n));
959566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(map, PETSC_DETERMINE));
969566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(map, 1));
979566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(map));
989566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(map, &com->N));
999566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(map, &com->rstart, NULL));
1009566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&map));
10147c6ae99SBarry Smith 
1029ae5db72SJed Brown   /* now set the rstart for each linked vector */
1039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
10547c6ae99SBarry Smith   while (next) {
10647c6ae99SBarry Smith     next->rstart = nprev;
10706ebdd98SJed Brown     nprev += next->n;
10847c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
1099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &next->grstarts));
1109566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgather(&next->grstart, 1, MPIU_INT, next->grstarts, 1, MPIU_INT, PetscObjectComm((PetscObject)dm)));
11147c6ae99SBarry Smith     next = next->next;
11247c6ae99SBarry Smith   }
11347c6ae99SBarry Smith   com->setup = PETSC_TRUE;
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11547c6ae99SBarry Smith }
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
11847c6ae99SBarry Smith 
11973e31fe2SJed Brown /*@
120*aaa8cc7dSPierre Jolivet     DMCompositeGetNumberDM - Gets the number of `DM` objects in the `DMCOMPOSITE`
12147c6ae99SBarry Smith        representation.
12247c6ae99SBarry Smith 
12347c6ae99SBarry Smith     Not Collective
12447c6ae99SBarry Smith 
12547c6ae99SBarry Smith     Input Parameter:
126dce8aebaSBarry Smith .    dm - the `DMCOMPOSITE` object
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith     Output Parameter:
129dce8aebaSBarry Smith .     nDM - the number of `DM`
13047c6ae99SBarry Smith 
13147c6ae99SBarry Smith     Level: beginner
13247c6ae99SBarry Smith 
133dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`
13447c6ae99SBarry Smith @*/
135d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetNumberDM(DM dm, PetscInt *nDM)
136d71ae5a4SJacob Faibussowitsch {
13747c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
13871b14b3eSStefano Zampini   PetscBool     flg;
1395fd66863SKarl Rupp 
14047c6ae99SBarry Smith   PetscFunctionBegin;
14147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1429566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1437a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
14447c6ae99SBarry Smith   *nDM = com->nDM;
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14647c6ae99SBarry Smith }
14747c6ae99SBarry Smith 
14847c6ae99SBarry Smith /*@C
14947c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
15047c6ae99SBarry Smith        representation.
15147c6ae99SBarry Smith 
15220f4b53cSBarry Smith     Collective
15347c6ae99SBarry Smith 
1549ae5db72SJed Brown     Input Parameters:
155dce8aebaSBarry Smith +    dm - the `DMCOMPOSITE` object
1569ae5db72SJed Brown -    gvec - the global vector
1579ae5db72SJed Brown 
1582fe279fdSBarry Smith     Output Parameter:
1590298fd71SBarry Smith .    Vec* ... - the packed parallel vectors, NULL for those that are not needed
16047c6ae99SBarry Smith 
16147c6ae99SBarry Smith     Level: advanced
16247c6ae99SBarry Smith 
163dce8aebaSBarry Smith     Note:
164dce8aebaSBarry Smith     Use `DMCompositeRestoreAccess()` to return the vectors when you no longer need them
165dce8aebaSBarry Smith 
166dce8aebaSBarry Smith     Fortran Note:
167dce8aebaSBarry Smith     Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
168dce8aebaSBarry Smith     or use the alternative interface `DMCompositeGetAccessArray()`.
169dce8aebaSBarry Smith 
170dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetEntries()`, `DMCompositeScatter()`
17147c6ae99SBarry Smith @*/
172d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccess(DM dm, Vec gvec, ...)
173d71ae5a4SJacob Faibussowitsch {
17447c6ae99SBarry Smith   va_list                 Argp;
17547c6ae99SBarry Smith   struct DMCompositeLink *next;
17647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
1775edff71fSBarry Smith   PetscInt                readonly;
17871b14b3eSStefano Zampini   PetscBool               flg;
17947c6ae99SBarry Smith 
18047c6ae99SBarry Smith   PetscFunctionBegin;
18147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
18247c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
1839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1847a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
18547c6ae99SBarry Smith   next = com->next;
18648a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
18747c6ae99SBarry Smith 
1889566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec, &readonly));
18947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
19047c6ae99SBarry Smith   va_start(Argp, gvec);
19147c6ae99SBarry Smith   while (next) {
19247c6ae99SBarry Smith     Vec *vec;
19347c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
1949ae5db72SJed Brown     if (vec) {
1959566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, vec));
1965edff71fSBarry Smith       if (readonly) {
1975edff71fSBarry Smith         const PetscScalar *array;
1989566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(gvec, &array));
1999566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec, array + next->rstart));
2009566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(*vec));
2019566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(gvec, &array));
2025edff71fSBarry Smith       } else {
2035edff71fSBarry Smith         PetscScalar *array;
2049566063dSJacob Faibussowitsch         PetscCall(VecGetArray(gvec, &array));
2059566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec, array + next->rstart));
2069566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(gvec, &array));
20747c6ae99SBarry Smith       }
2085edff71fSBarry Smith     }
20947c6ae99SBarry Smith     next = next->next;
21047c6ae99SBarry Smith   }
21147c6ae99SBarry Smith   va_end(Argp);
2123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21347c6ae99SBarry Smith }
21447c6ae99SBarry Smith 
215f73e5cebSJed Brown /*@C
216f73e5cebSJed Brown     DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
217f73e5cebSJed Brown        representation.
218f73e5cebSJed Brown 
21920f4b53cSBarry Smith     Collective
220f73e5cebSJed Brown 
221f73e5cebSJed Brown     Input Parameters:
222dce8aebaSBarry Smith +    dm - the `DMCOMPOSITE`
223f73e5cebSJed Brown .    pvec - packed vector
224f73e5cebSJed Brown .    nwanted - number of vectors wanted
2250298fd71SBarry Smith -    wanted - sorted array of vectors wanted, or NULL to get all vectors
226f73e5cebSJed Brown 
2272fe279fdSBarry Smith     Output Parameter:
228f73e5cebSJed Brown .    vecs - array of requested global vectors (must be allocated)
229f73e5cebSJed Brown 
230f73e5cebSJed Brown     Level: advanced
231f73e5cebSJed Brown 
232dce8aebaSBarry Smith     Note:
233dce8aebaSBarry Smith     Use `DMCompositeRestoreAccessArray()` to return the vectors when you no longer need them
234dce8aebaSBarry Smith 
235dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
236f73e5cebSJed Brown @*/
237d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs)
238d71ae5a4SJacob Faibussowitsch {
239f73e5cebSJed Brown   struct DMCompositeLink *link;
240f73e5cebSJed Brown   PetscInt                i, wnum;
241f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
242bee642f7SBarry Smith   PetscInt                readonly;
24371b14b3eSStefano Zampini   PetscBool               flg;
244f73e5cebSJed Brown 
245f73e5cebSJed Brown   PetscFunctionBegin;
246f73e5cebSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
247f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
2489566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
2497a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
25048a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
251f73e5cebSJed Brown 
2529566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
253f73e5cebSJed Brown   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
254f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
255f73e5cebSJed Brown       Vec v;
2569566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(link->dm, &v));
257bee642f7SBarry Smith       if (readonly) {
258bee642f7SBarry Smith         const PetscScalar *array;
2599566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec, &array));
2609566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + link->rstart));
2619566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
2629566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec, &array));
263bee642f7SBarry Smith       } else {
264bee642f7SBarry Smith         PetscScalar *array;
2659566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec, &array));
2669566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + link->rstart));
2679566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec, &array));
268bee642f7SBarry Smith       }
269f73e5cebSJed Brown       vecs[wnum++] = v;
270f73e5cebSJed Brown     }
271f73e5cebSJed Brown   }
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
273f73e5cebSJed Brown }
274f73e5cebSJed Brown 
2757ac2b803SAlex Fikl /*@C
2767ac2b803SAlex Fikl     DMCompositeGetLocalAccessArray - Allows one to access the individual
2777ac2b803SAlex Fikl     packed vectors in their local representation.
2787ac2b803SAlex Fikl 
27920f4b53cSBarry Smith     Collective
2807ac2b803SAlex Fikl 
2817ac2b803SAlex Fikl     Input Parameters:
282dce8aebaSBarry Smith +    dm - the `DMCOMPOSITE`
2837ac2b803SAlex Fikl .    pvec - packed vector
2847ac2b803SAlex Fikl .    nwanted - number of vectors wanted
2857ac2b803SAlex Fikl -    wanted - sorted array of vectors wanted, or NULL to get all vectors
2867ac2b803SAlex Fikl 
2872fe279fdSBarry Smith     Output Parameter:
2887ac2b803SAlex Fikl .    vecs - array of requested local vectors (must be allocated)
2897ac2b803SAlex Fikl 
2907ac2b803SAlex Fikl     Level: advanced
2917ac2b803SAlex Fikl 
292dce8aebaSBarry Smith     Note:
293dce8aebaSBarry Smith     Use `DMCompositeRestoreLocalAccessArray()` to return the vectors
294dce8aebaSBarry Smith     when you no longer need them.
295dce8aebaSBarry Smith 
296dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`,
297db781477SPatrick Sanan           `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
2987ac2b803SAlex Fikl @*/
299d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs)
300d71ae5a4SJacob Faibussowitsch {
3017ac2b803SAlex Fikl   struct DMCompositeLink *link;
3027ac2b803SAlex Fikl   PetscInt                i, wnum;
3037ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite *)dm->data;
3047ac2b803SAlex Fikl   PetscInt                readonly;
3057ac2b803SAlex Fikl   PetscInt                nlocal = 0;
30671b14b3eSStefano Zampini   PetscBool               flg;
3077ac2b803SAlex Fikl 
3087ac2b803SAlex Fikl   PetscFunctionBegin;
3097ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3107ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
3119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3127a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
31348a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
3147ac2b803SAlex Fikl 
3159566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
3167ac2b803SAlex Fikl   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
3177ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
3187ac2b803SAlex Fikl       Vec v;
3199566063dSJacob Faibussowitsch       PetscCall(DMGetLocalVector(link->dm, &v));
3207ac2b803SAlex Fikl       if (readonly) {
3217ac2b803SAlex Fikl         const PetscScalar *array;
3229566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec, &array));
3239566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + nlocal));
324b1c3483dSMark 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
3259566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
3269566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec, &array));
3277ac2b803SAlex Fikl       } else {
3287ac2b803SAlex Fikl         PetscScalar *array;
3299566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec, &array));
3309566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + nlocal));
3319566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec, &array));
3327ac2b803SAlex Fikl       }
3337ac2b803SAlex Fikl       vecs[wnum++] = v;
3347ac2b803SAlex Fikl     }
3357ac2b803SAlex Fikl 
3367ac2b803SAlex Fikl     nlocal += link->nlocal;
3377ac2b803SAlex Fikl   }
3387ac2b803SAlex Fikl 
3393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3407ac2b803SAlex Fikl }
3417ac2b803SAlex Fikl 
34247c6ae99SBarry Smith /*@C
343dce8aebaSBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with `DMCompositeGetAccess()`
34447c6ae99SBarry Smith        representation.
34547c6ae99SBarry Smith 
34620f4b53cSBarry Smith     Collective
34747c6ae99SBarry Smith 
3489ae5db72SJed Brown     Input Parameters:
349dce8aebaSBarry Smith +    dm - the `DMCOMPOSITE` object
35047c6ae99SBarry Smith .    gvec - the global vector
3510298fd71SBarry Smith -    Vec* ... - the individual parallel vectors, NULL for those that are not needed
35247c6ae99SBarry Smith 
35347c6ae99SBarry Smith     Level: advanced
35447c6ae99SBarry Smith 
355dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
356db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`,
357db781477SPatrick Sanan          `DMCompositeRestoreAccess()`, `DMCompositeGetAccess()`
35847c6ae99SBarry Smith @*/
359d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...)
360d71ae5a4SJacob Faibussowitsch {
36147c6ae99SBarry Smith   va_list                 Argp;
36247c6ae99SBarry Smith   struct DMCompositeLink *next;
36347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
3645edff71fSBarry Smith   PetscInt                readonly;
36571b14b3eSStefano Zampini   PetscBool               flg;
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith   PetscFunctionBegin;
36847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36947c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
3709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3717a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
37247c6ae99SBarry Smith   next = com->next;
37348a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
37447c6ae99SBarry Smith 
3759566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec, &readonly));
37647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
37747c6ae99SBarry Smith   va_start(Argp, gvec);
37847c6ae99SBarry Smith   while (next) {
37947c6ae99SBarry Smith     Vec *vec;
38047c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
3819ae5db72SJed Brown     if (vec) {
3829566063dSJacob Faibussowitsch       PetscCall(VecResetArray(*vec));
3831baa6e33SBarry Smith       if (readonly) PetscCall(VecLockReadPop(*vec));
3849566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, vec));
38547c6ae99SBarry Smith     }
38647c6ae99SBarry Smith     next = next->next;
38747c6ae99SBarry Smith   }
38847c6ae99SBarry Smith   va_end(Argp);
3893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39047c6ae99SBarry Smith }
39147c6ae99SBarry Smith 
392f73e5cebSJed Brown /*@C
393dce8aebaSBarry Smith     DMCompositeRestoreAccessArray - Returns the vectors obtained with `DMCompositeGetAccessArray()`
394f73e5cebSJed Brown 
39520f4b53cSBarry Smith     Collective
396f73e5cebSJed Brown 
397f73e5cebSJed Brown     Input Parameters:
398dce8aebaSBarry Smith +    dm - the `DMCOMPOSITE` object
399f73e5cebSJed Brown .    pvec - packed vector
400f73e5cebSJed Brown .    nwanted - number of vectors wanted
4010298fd71SBarry Smith .    wanted - sorted array of vectors wanted, or NULL to get all vectors
402f73e5cebSJed Brown -    vecs - array of global vectors to return
403f73e5cebSJed Brown 
404f73e5cebSJed Brown     Level: advanced
405f73e5cebSJed Brown 
406dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
407f73e5cebSJed Brown @*/
408d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs)
409d71ae5a4SJacob Faibussowitsch {
410f73e5cebSJed Brown   struct DMCompositeLink *link;
411f73e5cebSJed Brown   PetscInt                i, wnum;
412f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
413bee642f7SBarry Smith   PetscInt                readonly;
41471b14b3eSStefano Zampini   PetscBool               flg;
415f73e5cebSJed Brown 
416f73e5cebSJed Brown   PetscFunctionBegin;
417f73e5cebSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
418f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4199566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4207a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
42148a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
422f73e5cebSJed Brown 
4239566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
424f73e5cebSJed Brown   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
425f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
4269566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
42748a46eb9SPierre Jolivet       if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4289566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum]));
429f73e5cebSJed Brown       wnum++;
430f73e5cebSJed Brown     }
431f73e5cebSJed Brown   }
4323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
433f73e5cebSJed Brown }
434f73e5cebSJed Brown 
4357ac2b803SAlex Fikl /*@C
436dce8aebaSBarry Smith     DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with `DMCompositeGetLocalAccessArray()`.
4377ac2b803SAlex Fikl 
43820f4b53cSBarry Smith     Collective
4397ac2b803SAlex Fikl 
4407ac2b803SAlex Fikl     Input Parameters:
441dce8aebaSBarry Smith +    dm - the `DMCOMPOSITE` object
4427ac2b803SAlex Fikl .    pvec - packed vector
4437ac2b803SAlex Fikl .    nwanted - number of vectors wanted
4447ac2b803SAlex Fikl .    wanted - sorted array of vectors wanted, or NULL to restore all vectors
4457ac2b803SAlex Fikl -    vecs - array of local vectors to return
4467ac2b803SAlex Fikl 
4477ac2b803SAlex Fikl     Level: advanced
4487ac2b803SAlex Fikl 
449dce8aebaSBarry Smith     Note:
450dce8aebaSBarry Smith     nwanted and wanted must match the values given to `DMCompositeGetLocalAccessArray()`
4517ac2b803SAlex Fikl     otherwise the call will fail.
4527ac2b803SAlex Fikl 
453dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`,
454db781477SPatrick Sanan           `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`,
455db781477SPatrick Sanan           `DMCompositeScatter()`, `DMCompositeGather()`
4567ac2b803SAlex Fikl @*/
457d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs)
458d71ae5a4SJacob Faibussowitsch {
4597ac2b803SAlex Fikl   struct DMCompositeLink *link;
4607ac2b803SAlex Fikl   PetscInt                i, wnum;
4617ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite *)dm->data;
4627ac2b803SAlex Fikl   PetscInt                readonly;
46371b14b3eSStefano Zampini   PetscBool               flg;
4647ac2b803SAlex Fikl 
4657ac2b803SAlex Fikl   PetscFunctionBegin;
4667ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4677ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4697a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
47048a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
4717ac2b803SAlex Fikl 
4729566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
4737ac2b803SAlex Fikl   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
4747ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
4759566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
47648a46eb9SPierre Jolivet       if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4779566063dSJacob Faibussowitsch       PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum]));
4787ac2b803SAlex Fikl       wnum++;
4797ac2b803SAlex Fikl     }
4807ac2b803SAlex Fikl   }
4813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4827ac2b803SAlex Fikl }
4837ac2b803SAlex Fikl 
48447c6ae99SBarry Smith /*@C
48547c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
48647c6ae99SBarry Smith 
48720f4b53cSBarry Smith     Collective
48847c6ae99SBarry Smith 
4899ae5db72SJed Brown     Input Parameters:
490dce8aebaSBarry Smith +    dm - the `DMCOMPOSITE` object
49147c6ae99SBarry Smith .    gvec - the global vector
4920298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for those that are not needed
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith     Level: advanced
49547c6ae99SBarry Smith 
496dce8aebaSBarry Smith     Note:
497dce8aebaSBarry Smith     `DMCompositeScatterArray()` is a non-variadic alternative that is often more convenient for library callers and is
4986f3c3dcfSJed Brown     accessible from Fortran.
4996f3c3dcfSJed Brown 
500dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
501db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
502db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
503db781477SPatrick Sanan          `DMCompositeScatterArray()`
50447c6ae99SBarry Smith @*/
505d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...)
506d71ae5a4SJacob Faibussowitsch {
50747c6ae99SBarry Smith   va_list                 Argp;
50847c6ae99SBarry Smith   struct DMCompositeLink *next;
509c599c493SJunchao Zhang   PETSC_UNUSED PetscInt   cnt;
51047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
51171b14b3eSStefano Zampini   PetscBool               flg;
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith   PetscFunctionBegin;
51447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
51547c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5169566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5177a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
51848a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
52147c6ae99SBarry Smith   va_start(Argp, gvec);
5228fd8f222SJed Brown   for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
5239ae5db72SJed Brown     Vec local;
5249ae5db72SJed Brown     local = va_arg(Argp, Vec);
5259ae5db72SJed Brown     if (local) {
5269ae5db72SJed Brown       Vec                global;
5275edff71fSBarry Smith       const PetscScalar *array;
52863a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
5299566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
5309566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
5319566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
5329566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local));
5339566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local));
5349566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
5359566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5369566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
53747c6ae99SBarry Smith     }
53847c6ae99SBarry Smith   }
53947c6ae99SBarry Smith   va_end(Argp);
5403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54147c6ae99SBarry Smith }
54247c6ae99SBarry Smith 
5436f3c3dcfSJed Brown /*@
5446f3c3dcfSJed Brown     DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5456f3c3dcfSJed Brown 
54620f4b53cSBarry Smith     Collective
5476f3c3dcfSJed Brown 
5486f3c3dcfSJed Brown     Input Parameters:
549dce8aebaSBarry Smith +    dm - the `DMCOMPOSITE` object
5506f3c3dcfSJed Brown .    gvec - the global vector
551a2b725a8SWilliam Gropp -    lvecs - array of local vectors, NULL for any that are not needed
5526f3c3dcfSJed Brown 
5536f3c3dcfSJed Brown     Level: advanced
5546f3c3dcfSJed Brown 
5556f3c3dcfSJed Brown     Note:
556dce8aebaSBarry Smith     This is a non-variadic alternative to `DMCompositeScatter()`
5576f3c3dcfSJed Brown 
558dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`
559db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
560db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
5616f3c3dcfSJed Brown @*/
562d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs)
563d71ae5a4SJacob Faibussowitsch {
5646f3c3dcfSJed Brown   struct DMCompositeLink *next;
5656f3c3dcfSJed Brown   PetscInt                i;
5666f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
56771b14b3eSStefano Zampini   PetscBool               flg;
5686f3c3dcfSJed Brown 
5696f3c3dcfSJed Brown   PetscFunctionBegin;
5706f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5716f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5737a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
57448a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
5756f3c3dcfSJed Brown 
5766f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
5776f3c3dcfSJed Brown   for (i = 0, next = com->next; next; next = next->next, i++) {
5786f3c3dcfSJed Brown     if (lvecs[i]) {
5796f3c3dcfSJed Brown       Vec                global;
580c5d31e75SLisandro Dalcin       const PetscScalar *array;
5816f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3);
5829566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
5839566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
5849566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart));
5859566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i]));
5869566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i]));
5879566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
5889566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5899566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
5906f3c3dcfSJed Brown     }
5916f3c3dcfSJed Brown   }
5923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5936f3c3dcfSJed Brown }
5946f3c3dcfSJed Brown 
59547c6ae99SBarry Smith /*@C
59647c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
59747c6ae99SBarry Smith 
59820f4b53cSBarry Smith     Collective
59947c6ae99SBarry Smith 
600d8d19677SJose E. Roman     Input Parameters:
601dce8aebaSBarry Smith +    dm - the `DMCOMPOSITE` object
60247c6ae99SBarry Smith .    gvec - the global vector
603dce8aebaSBarry Smith .    imode - `INSERT_VALUES` or `ADD_VALUES`
6040298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for any that are not needed
60547c6ae99SBarry Smith 
60647c6ae99SBarry Smith     Level: advanced
60747c6ae99SBarry Smith 
608dce8aebaSBarry Smith     Fortran Note:
609dce8aebaSBarry Smith     Fortran users should use `DMCompositeGatherArray()`
610f5f57ec0SBarry Smith 
611dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
612db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
613db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
61447c6ae99SBarry Smith @*/
615d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...)
616d71ae5a4SJacob Faibussowitsch {
61747c6ae99SBarry Smith   va_list                 Argp;
61847c6ae99SBarry Smith   struct DMCompositeLink *next;
61947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
620c599c493SJunchao Zhang   PETSC_UNUSED PetscInt   cnt;
62171b14b3eSStefano Zampini   PetscBool               flg;
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith   PetscFunctionBegin;
62447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
625064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6269566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6277a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
62848a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
62947c6ae99SBarry Smith 
63047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
6311dac896bSSatish Balay   va_start(Argp, gvec);
6328fd8f222SJed Brown   for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
6339ae5db72SJed Brown     Vec local;
6349ae5db72SJed Brown     local = va_arg(Argp, Vec);
6359ae5db72SJed Brown     if (local) {
63647c6ae99SBarry Smith       PetscScalar *array;
6379ae5db72SJed Brown       Vec          global;
63863a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
6399566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
6409566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec, &array));
6419566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
6429566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global));
6439566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global));
6449566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec, &array));
6459566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6469566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
64747c6ae99SBarry Smith     }
64847c6ae99SBarry Smith   }
64947c6ae99SBarry Smith   va_end(Argp);
6503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
65147c6ae99SBarry Smith }
65247c6ae99SBarry Smith 
6536f3c3dcfSJed Brown /*@
6546f3c3dcfSJed Brown     DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6556f3c3dcfSJed Brown 
65620f4b53cSBarry Smith     Collective
6576f3c3dcfSJed Brown 
658d8d19677SJose E. Roman     Input Parameters:
659dce8aebaSBarry Smith +    dm - the `DMCOMPOSITE` object
6606f3c3dcfSJed Brown .    gvec - the global vector
661dce8aebaSBarry Smith .    imode - `INSERT_VALUES` or `ADD_VALUES`
6626f3c3dcfSJed Brown -    lvecs - the individual sequential vectors, NULL for any that are not needed
6636f3c3dcfSJed Brown 
6646f3c3dcfSJed Brown     Level: advanced
6656f3c3dcfSJed Brown 
666dce8aebaSBarry Smith     Note:
667dce8aebaSBarry Smith     This is a non-variadic alternative to `DMCompositeGather()`.
6686f3c3dcfSJed Brown 
669dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
670db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
671db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`,
6726f3c3dcfSJed Brown @*/
673d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs)
674d71ae5a4SJacob Faibussowitsch {
6756f3c3dcfSJed Brown   struct DMCompositeLink *next;
6766f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
6776f3c3dcfSJed Brown   PetscInt                i;
67871b14b3eSStefano Zampini   PetscBool               flg;
6796f3c3dcfSJed Brown 
6806f3c3dcfSJed Brown   PetscFunctionBegin;
6816f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
682064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6847a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
68548a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
6866f3c3dcfSJed Brown 
6876f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
6886f3c3dcfSJed Brown   for (next = com->next, i = 0; next; next = next->next, i++) {
6896f3c3dcfSJed Brown     if (lvecs[i]) {
6906f3c3dcfSJed Brown       PetscScalar *array;
6916f3c3dcfSJed Brown       Vec          global;
692064a246eSJacob Faibussowitsch       PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4);
6939566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
6949566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec, &array));
6959566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
6969566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global));
6979566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global));
6989566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec, &array));
6999566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
7009566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
7016f3c3dcfSJed Brown     }
7026f3c3dcfSJed Brown   }
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7046f3c3dcfSJed Brown }
7056f3c3dcfSJed Brown 
706f5f57ec0SBarry Smith /*@
707dce8aebaSBarry Smith     DMCompositeAddDM - adds a `DM` vector to a `DMCOMPOSITE`
70847c6ae99SBarry Smith 
70920f4b53cSBarry Smith     Collective
71047c6ae99SBarry Smith 
711d8d19677SJose E. Roman     Input Parameters:
712dce8aebaSBarry Smith +    dmc - the  `DMCOMPOSITE` object
713dce8aebaSBarry Smith -    dm - the `DM` object
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith     Level: advanced
71647c6ae99SBarry Smith 
717dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeGather()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
718db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
719db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
72047c6ae99SBarry Smith @*/
721d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeAddDM(DM dmc, DM dm)
722d71ae5a4SJacob Faibussowitsch {
72306ebdd98SJed Brown   PetscInt                n, nlocal;
72447c6ae99SBarry Smith   struct DMCompositeLink *mine, *next;
72506ebdd98SJed Brown   Vec                     global, local;
72647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dmc->data;
72771b14b3eSStefano Zampini   PetscBool               flg;
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith   PetscFunctionBegin;
73047c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc, DM_CLASSID, 1);
73147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
7329566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg));
7337a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
73447c6ae99SBarry Smith   next = com->next;
73528b400f6SJacob Faibussowitsch   PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite");
73647c6ae99SBarry Smith 
73747c6ae99SBarry Smith   /* create new link */
7389566063dSJacob Faibussowitsch   PetscCall(PetscNew(&mine));
7399566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
7409566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &global));
7419566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global, &n));
7429566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &global));
7439566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &local));
7449566063dSJacob Faibussowitsch   PetscCall(VecGetSize(local, &nlocal));
7459566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &local));
7468865f1eaSKarl Rupp 
74747c6ae99SBarry Smith   mine->n      = n;
74806ebdd98SJed Brown   mine->nlocal = nlocal;
74947c6ae99SBarry Smith   mine->dm     = dm;
7500298fd71SBarry Smith   mine->next   = NULL;
75147c6ae99SBarry Smith   com->n += n;
7527ac2b803SAlex Fikl   com->nghost += nlocal;
75347c6ae99SBarry Smith 
75447c6ae99SBarry Smith   /* add to end of list */
7558865f1eaSKarl Rupp   if (!next) com->next = mine;
7568865f1eaSKarl Rupp   else {
75747c6ae99SBarry Smith     while (next->next) next = next->next;
75847c6ae99SBarry Smith     next->next = mine;
75947c6ae99SBarry Smith   }
76047c6ae99SBarry Smith   com->nDM++;
76147c6ae99SBarry Smith   com->nmine++;
7623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76347c6ae99SBarry Smith }
76447c6ae99SBarry Smith 
7659804daf3SBarry Smith #include <petscdraw.h>
76626887b52SJed Brown PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
767d71ae5a4SJacob Faibussowitsch PetscErrorCode              VecView_DMComposite(Vec gvec, PetscViewer viewer)
768d71ae5a4SJacob Faibussowitsch {
76947c6ae99SBarry Smith   DM                      dm;
77047c6ae99SBarry Smith   struct DMCompositeLink *next;
77147c6ae99SBarry Smith   PetscBool               isdraw;
772cef07954SSatish Balay   DM_Composite           *com;
77347c6ae99SBarry Smith 
77447c6ae99SBarry Smith   PetscFunctionBegin;
7759566063dSJacob Faibussowitsch   PetscCall(VecGetDM(gvec, &dm));
7767a8be351SBarry Smith   PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite");
77747c6ae99SBarry Smith   com  = (DM_Composite *)dm->data;
77847c6ae99SBarry Smith   next = com->next;
77947c6ae99SBarry Smith 
7809566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
78147c6ae99SBarry Smith   if (!isdraw) {
78247c6ae99SBarry Smith     /* do I really want to call this? */
7839566063dSJacob Faibussowitsch     PetscCall(VecView_MPI(gvec, viewer));
78447c6ae99SBarry Smith   } else {
78547c6ae99SBarry Smith     PetscInt cnt = 0;
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
78847c6ae99SBarry Smith     while (next) {
78947c6ae99SBarry Smith       Vec                vec;
790ca4278abSLisandro Dalcin       const PetscScalar *array;
79147c6ae99SBarry Smith       PetscInt           bs;
79247c6ae99SBarry Smith 
7939ae5db72SJed Brown       /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
7949566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &vec));
7959566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
7969566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart));
7979566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
7989566063dSJacob Faibussowitsch       PetscCall(VecView(vec, viewer));
7999566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vec));
8009566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(vec, &bs));
8019566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &vec));
8029566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawBaseAdd(viewer, bs));
80347c6ae99SBarry Smith       cnt += bs;
80447c6ae99SBarry Smith       next = next->next;
80547c6ae99SBarry Smith     }
8069566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt));
80747c6ae99SBarry Smith   }
8083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80947c6ae99SBarry Smith }
81047c6ae99SBarry Smith 
811d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec)
812d71ae5a4SJacob Faibussowitsch {
81347c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
81447c6ae99SBarry Smith 
81547c6ae99SBarry Smith   PetscFunctionBegin;
81647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8179566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
8189566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8199566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec));
8209566063dSJacob Faibussowitsch   PetscCall(VecSetType(*gvec, dm->vectype));
8219566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*gvec, com->n, com->N));
8229566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*gvec, dm));
8239566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (void (*)(void))VecView_DMComposite));
8243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82547c6ae99SBarry Smith }
82647c6ae99SBarry Smith 
827d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec)
828d71ae5a4SJacob Faibussowitsch {
82947c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
83047c6ae99SBarry Smith 
83147c6ae99SBarry Smith   PetscFunctionBegin;
83247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
83347c6ae99SBarry Smith   if (!com->setup) {
8349566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(dm));
8359566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
83647c6ae99SBarry Smith   }
8379566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, lvec));
8389566063dSJacob Faibussowitsch   PetscCall(VecSetType(*lvec, dm->vectype));
8399566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE));
8409566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*lvec, dm));
8413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
84247c6ae99SBarry Smith }
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith /*@C
845dce8aebaSBarry Smith     DMCompositeGetISLocalToGlobalMappings - gets an `ISLocalToGlobalMapping` for each `DM` in the `DMCOMPOSITE`, maps to the composite global space
84647c6ae99SBarry Smith 
84720f4b53cSBarry Smith     Collective; No Fortran Support
84847c6ae99SBarry Smith 
84947c6ae99SBarry Smith     Input Parameter:
850dce8aebaSBarry Smith .    dm - the `DMCOMPOSITE` object
85147c6ae99SBarry Smith 
8522fe279fdSBarry Smith     Output Parameter:
8539ae5db72SJed Brown .    ltogs - the individual mappings for each packed vector. Note that this includes
854dce8aebaSBarry Smith            all the ghost points that individual ghosted `DMDA` may have.
85547c6ae99SBarry Smith 
85647c6ae99SBarry Smith     Level: advanced
85747c6ae99SBarry Smith 
858dce8aebaSBarry Smith     Note:
859dce8aebaSBarry Smith     Each entry of ltogs should be destroyed with `ISLocalToGlobalMappingDestroy()`, the ltogs array should be freed with `PetscFree()`.
86047c6ae99SBarry Smith 
861dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
862db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
863c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
86447c6ae99SBarry Smith @*/
865d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping **ltogs)
866d71ae5a4SJacob Faibussowitsch {
86747c6ae99SBarry Smith   PetscInt                i, *idx, n, cnt;
86847c6ae99SBarry Smith   struct DMCompositeLink *next;
86947c6ae99SBarry Smith   PetscMPIInt             rank;
87047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
87171b14b3eSStefano Zampini   PetscBool               flg;
87247c6ae99SBarry Smith 
87347c6ae99SBarry Smith   PetscFunctionBegin;
87447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8759566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
8767a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
8779566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM, ltogs));
87947c6ae99SBarry Smith   next = com->next;
8809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
88147c6ae99SBarry Smith 
88247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
88347c6ae99SBarry Smith   cnt = 0;
88447c6ae99SBarry Smith   while (next) {
8856eb61c8cSJed Brown     ISLocalToGlobalMapping ltog;
8866eb61c8cSJed Brown     PetscMPIInt            size;
88786994e45SJed Brown     const PetscInt        *suboff, *indices;
8886eb61c8cSJed Brown     Vec                    global;
88947c6ae99SBarry Smith 
8906eb61c8cSJed Brown     /* Get sub-DM global indices for each local dof */
8919566063dSJacob Faibussowitsch     PetscCall(DMGetLocalToGlobalMapping(next->dm, &ltog));
8929566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n));
8939566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices));
8949566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &idx));
89547c6ae99SBarry Smith 
8966eb61c8cSJed Brown     /* Get the offsets for the sub-DM global vector */
8979566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
8989566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRanges(global, &suboff));
8999566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size));
9006eb61c8cSJed Brown 
9016eb61c8cSJed Brown     /* Shift the sub-DM definition of the global space to the composite global space */
9026eb61c8cSJed Brown     for (i = 0; i < n; i++) {
90386994e45SJed Brown       PetscInt subi = indices[i], lo = 0, hi = size, t;
904e333b035SStefano Zampini       /* There's no consensus on what a negative index means,
905e333b035SStefano Zampini          except for skipping when setting the values in vectors and matrices */
9069371c9d4SSatish Balay       if (subi < 0) {
9079371c9d4SSatish Balay         idx[i] = subi - next->grstarts[rank];
9089371c9d4SSatish Balay         continue;
9099371c9d4SSatish Balay       }
9106eb61c8cSJed Brown       /* Binary search to find which rank owns subi */
9116eb61c8cSJed Brown       while (hi - lo > 1) {
9126eb61c8cSJed Brown         t = lo + (hi - lo) / 2;
9136eb61c8cSJed Brown         if (suboff[t] > subi) hi = t;
9146eb61c8cSJed Brown         else lo = t;
9156eb61c8cSJed Brown       }
9166eb61c8cSJed Brown       idx[i] = subi - suboff[lo] + next->grstarts[lo];
9176eb61c8cSJed Brown     }
9189566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices));
9199566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt]));
9209566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
92147c6ae99SBarry Smith     next = next->next;
92247c6ae99SBarry Smith     cnt++;
92347c6ae99SBarry Smith   }
9243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92547c6ae99SBarry Smith }
92647c6ae99SBarry Smith 
92787c85e80SJed Brown /*@C
9289ae5db72SJed Brown    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
92987c85e80SJed Brown 
93020f4b53cSBarry Smith    Not Collective; No Fortran Support
93187c85e80SJed Brown 
9324165533cSJose E. Roman    Input Parameter:
933dce8aebaSBarry Smith . dm - the `DMCOMPOSITE`
93487c85e80SJed Brown 
9354165533cSJose E. Roman    Output Parameter:
936dce8aebaSBarry Smith . is - array of serial index sets for each each component of the `DMCOMPOSITE`
93787c85e80SJed Brown 
93887c85e80SJed Brown    Level: intermediate
93987c85e80SJed Brown 
94087c85e80SJed Brown    Notes:
94187c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
942dce8aebaSBarry Smith    `MatGetLocalSubMatrix()`.  In the future, the scatters for each entry in the `DMCOMPOSITE` may be be merged into a single
9439ae5db72SJed Brown    scatter to a composite local vector.  The user should not typically need to know which is being done.
94487c85e80SJed Brown 
945dce8aebaSBarry Smith    To get the composite global indices at all local points (including ghosts), use `DMCompositeGetISLocalToGlobalMappings()`.
94687c85e80SJed Brown 
947dce8aebaSBarry Smith    To get index sets for pieces of the composite global vector, use `DMCompositeGetGlobalISs()`.
94887c85e80SJed Brown 
949dce8aebaSBarry Smith    Each returned `IS` should be destroyed with `ISDestroy()`, the array should be freed with `PetscFree()`.
95087c85e80SJed Brown 
951dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`, `MatCreateLocalRef()`
95287c85e80SJed Brown @*/
953d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalISs(DM dm, IS **is)
954d71ae5a4SJacob Faibussowitsch {
95587c85e80SJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
95687c85e80SJed Brown   struct DMCompositeLink *link;
95787c85e80SJed Brown   PetscInt                cnt, start;
95871b14b3eSStefano Zampini   PetscBool               flg;
95987c85e80SJed Brown 
96087c85e80SJed Brown   PetscFunctionBegin;
96187c85e80SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
96287c85e80SJed Brown   PetscValidPointer(is, 2);
9639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
9647a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
9659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nmine, is));
96606ebdd98SJed Brown   for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) {
967520db06cSJed Brown     PetscInt bs;
9689566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt]));
9699566063dSJacob Faibussowitsch     PetscCall(DMGetBlockSize(link->dm, &bs));
9709566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize((*is)[cnt], bs));
971520db06cSJed Brown   }
9723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97387c85e80SJed Brown }
97487c85e80SJed Brown 
97547c6ae99SBarry Smith /*@C
976dce8aebaSBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object in a `DMCOMPOSITE`
97747c6ae99SBarry Smith 
97820f4b53cSBarry Smith     Collective
97947c6ae99SBarry Smith 
98047c6ae99SBarry Smith     Input Parameter:
981dce8aebaSBarry Smith .    dm - the `DMCOMPOSITE` object
98247c6ae99SBarry Smith 
9832fe279fdSBarry Smith     Output Parameter:
98447c6ae99SBarry Smith .    is - the array of index sets
98547c6ae99SBarry Smith 
98647c6ae99SBarry Smith     Level: advanced
98747c6ae99SBarry Smith 
98847c6ae99SBarry Smith     Notes:
989dce8aebaSBarry Smith        The is entries should be destroyed with `ISDestroy()`, the is array should be freed with `PetscFree()`
99047c6ae99SBarry Smith 
99147c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
99247c6ae99SBarry Smith 
993dce8aebaSBarry Smith        Use `DMCompositeGetLocalISs()` for index sets in the packed local numbering, and
994dce8aebaSBarry Smith        `DMCompositeGetISLocalToGlobalMappings()` for to map local sub-`DM` (including ghost) indices to packed global
9956eb61c8cSJed Brown        indices.
99647c6ae99SBarry Smith 
997dce8aebaSBarry Smith     Fortran Note:
998dce8aebaSBarry Smith     The output argument 'is' must be an allocated array of sufficient length, which can be learned using `DMCompositeGetNumberDM()`.
999f3cb0f7eSJed Brown 
1000dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1001db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
1002c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
100347c6ae99SBarry Smith @*/
1004d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[])
1005d71ae5a4SJacob Faibussowitsch {
100666bb578eSMark Adams   PetscInt                cnt = 0;
100747c6ae99SBarry Smith   struct DMCompositeLink *next;
100847c6ae99SBarry Smith   PetscMPIInt             rank;
100947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
101071b14b3eSStefano Zampini   PetscBool               flg;
101147c6ae99SBarry Smith 
101247c6ae99SBarry Smith   PetscFunctionBegin;
101347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
10157a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
10167a8be351SBarry Smith   PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before");
10179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM, is));
101847c6ae99SBarry Smith   next = com->next;
10199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
102047c6ae99SBarry Smith 
102147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
102247c6ae99SBarry Smith   while (next) {
1023e5e52638SMatthew G. Knepley     PetscDS prob;
1024e5e52638SMatthew G. Knepley 
10259566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt]));
10269566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &prob));
1027e5e52638SMatthew G. Knepley     if (prob) {
102865c226d8SMatthew G. Knepley       MatNullSpace space;
102965c226d8SMatthew G. Knepley       Mat          pmat;
10300f21e855SMatthew G. Knepley       PetscObject  disc;
10310f21e855SMatthew G. Knepley       PetscInt     Nf;
103265c226d8SMatthew G. Knepley 
10339566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(prob, &Nf));
1034f24dd8d2SMatthew G. Knepley       if (cnt < Nf) {
10359566063dSJacob Faibussowitsch         PetscCall(PetscDSGetDiscretization(prob, cnt, &disc));
10369566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space));
10379566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space));
10389566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space));
10399566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space));
10409566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat));
10419566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat));
104265c226d8SMatthew G. Knepley       }
1043f24dd8d2SMatthew G. Knepley     }
104447c6ae99SBarry Smith     cnt++;
104547c6ae99SBarry Smith     next = next->next;
104647c6ae99SBarry Smith   }
10473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104847c6ae99SBarry Smith }
104947c6ae99SBarry Smith 
1050d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
1051d71ae5a4SJacob Faibussowitsch {
10524d343eeaSMatthew G Knepley   PetscInt nDM;
10534d343eeaSMatthew G Knepley   DM      *dms;
10544d343eeaSMatthew G Knepley   PetscInt i;
10554d343eeaSMatthew G Knepley 
10564d343eeaSMatthew G Knepley   PetscFunctionBegin;
10579566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm, &nDM));
10588865f1eaSKarl Rupp   if (numFields) *numFields = nDM;
10599566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetGlobalISs(dm, fields));
10604d343eeaSMatthew G Knepley   if (fieldNames) {
10619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, &dms));
10629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, fieldNames));
10639566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, dms));
10644d343eeaSMatthew G Knepley     for (i = 0; i < nDM; i++) {
10654d343eeaSMatthew G Knepley       char        buf[256];
10664d343eeaSMatthew G Knepley       const char *splitname;
10674d343eeaSMatthew G Knepley 
10684d343eeaSMatthew G Knepley       /* Split naming precedence: object name, prefix, number */
10694d343eeaSMatthew G Knepley       splitname = ((PetscObject)dm)->name;
10704d343eeaSMatthew G Knepley       if (!splitname) {
10719566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname));
10724d343eeaSMatthew G Knepley         if (splitname) {
10734d343eeaSMatthew G Knepley           size_t len;
10749566063dSJacob Faibussowitsch           PetscCall(PetscStrncpy(buf, splitname, sizeof(buf)));
10758caf3d72SBarry Smith           buf[sizeof(buf) - 1] = 0;
10769566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(buf, &len));
10774d343eeaSMatthew G Knepley           if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */
10784d343eeaSMatthew G Knepley           splitname = buf;
10794d343eeaSMatthew G Knepley         }
10804d343eeaSMatthew G Knepley       }
10814d343eeaSMatthew G Knepley       if (!splitname) {
108263a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i));
10834d343eeaSMatthew G Knepley         splitname = buf;
10844d343eeaSMatthew G Knepley       }
10859566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i]));
10864d343eeaSMatthew G Knepley     }
10879566063dSJacob Faibussowitsch     PetscCall(PetscFree(dms));
10884d343eeaSMatthew G Knepley   }
10893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10904d343eeaSMatthew G Knepley }
10914d343eeaSMatthew G Knepley 
1092e7c4fc90SDmitry Karpeev /*
1093e7c4fc90SDmitry Karpeev  This could take over from DMCreateFieldIS(), as it is more general,
10940298fd71SBarry Smith  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1095e7c4fc90SDmitry Karpeev  At this point it's probably best to be less intrusive, however.
1096e7c4fc90SDmitry Karpeev  */
1097d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
1098d71ae5a4SJacob Faibussowitsch {
1099e7c4fc90SDmitry Karpeev   PetscInt nDM;
1100e7c4fc90SDmitry Karpeev   PetscInt i;
1101e7c4fc90SDmitry Karpeev 
1102e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
11039566063dSJacob Faibussowitsch   PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist));
1104e7c4fc90SDmitry Karpeev   if (dmlist) {
11059566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetNumberDM(dm, &nDM));
11069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, dmlist));
11079566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, *dmlist));
110848a46eb9SPierre Jolivet     for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i])));
1109e7c4fc90SDmitry Karpeev   }
11103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1111e7c4fc90SDmitry Karpeev }
1112e7c4fc90SDmitry Karpeev 
111347c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
111447c6ae99SBarry Smith /*@C
1115dce8aebaSBarry Smith     DMCompositeGetLocalVectors - Gets local vectors for each part of a `DMCOMPOSITE`
1116dce8aebaSBarry Smith        Use `DMCompositeRestoreLocalVectors()` to return them.
111747c6ae99SBarry Smith 
111820f4b53cSBarry Smith     Not Collective; No Fortran Support
111947c6ae99SBarry Smith 
112047c6ae99SBarry Smith     Input Parameter:
1121dce8aebaSBarry Smith .    dm - the `DMCOMPOSITE` object
112247c6ae99SBarry Smith 
112347c6ae99SBarry Smith     Output Parameter:
112420f4b53cSBarry Smith .   Vec ... - the individual sequential `Vec`s
112547c6ae99SBarry Smith 
112647c6ae99SBarry Smith     Level: advanced
112747c6ae99SBarry Smith 
1128dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1129db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1130db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
113147c6ae99SBarry Smith @*/
1132d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...)
1133d71ae5a4SJacob Faibussowitsch {
113447c6ae99SBarry Smith   va_list                 Argp;
113547c6ae99SBarry Smith   struct DMCompositeLink *next;
113647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
113771b14b3eSStefano Zampini   PetscBool               flg;
113847c6ae99SBarry Smith 
113947c6ae99SBarry Smith   PetscFunctionBegin;
114047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11419566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11427a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
114347c6ae99SBarry Smith   next = com->next;
114447c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
114547c6ae99SBarry Smith   va_start(Argp, dm);
114647c6ae99SBarry Smith   while (next) {
114747c6ae99SBarry Smith     Vec *vec;
114847c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
11499566063dSJacob Faibussowitsch     if (vec) PetscCall(DMGetLocalVector(next->dm, vec));
115047c6ae99SBarry Smith     next = next->next;
115147c6ae99SBarry Smith   }
115247c6ae99SBarry Smith   va_end(Argp);
11533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115447c6ae99SBarry Smith }
115547c6ae99SBarry Smith 
115647c6ae99SBarry Smith /*@C
1157dce8aebaSBarry Smith     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a `DMCOMPOSITE`
115847c6ae99SBarry Smith 
115920f4b53cSBarry Smith     Not Collective; No Fortran Support
116047c6ae99SBarry Smith 
116147c6ae99SBarry Smith     Input Parameter:
1162dce8aebaSBarry Smith .    dm - the `DMCOMPOSITE` object
116347c6ae99SBarry Smith 
116447c6ae99SBarry Smith     Output Parameter:
1165dce8aebaSBarry Smith .   Vec ... - the individual sequential `Vec`
116647c6ae99SBarry Smith 
116747c6ae99SBarry Smith     Level: advanced
116847c6ae99SBarry Smith 
1169dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1170db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1171db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
117247c6ae99SBarry Smith @*/
1173d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...)
1174d71ae5a4SJacob Faibussowitsch {
117547c6ae99SBarry Smith   va_list                 Argp;
117647c6ae99SBarry Smith   struct DMCompositeLink *next;
117747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
117871b14b3eSStefano Zampini   PetscBool               flg;
117947c6ae99SBarry Smith 
118047c6ae99SBarry Smith   PetscFunctionBegin;
118147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11829566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11837a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
118447c6ae99SBarry Smith   next = com->next;
118547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
118647c6ae99SBarry Smith   va_start(Argp, dm);
118747c6ae99SBarry Smith   while (next) {
118847c6ae99SBarry Smith     Vec *vec;
118947c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
11909566063dSJacob Faibussowitsch     if (vec) PetscCall(DMRestoreLocalVector(next->dm, vec));
119147c6ae99SBarry Smith     next = next->next;
119247c6ae99SBarry Smith   }
119347c6ae99SBarry Smith   va_end(Argp);
11943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119547c6ae99SBarry Smith }
119647c6ae99SBarry Smith 
119747c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
119847c6ae99SBarry Smith /*@C
1199dce8aebaSBarry Smith     DMCompositeGetEntries - Gets the `DM` for each entry in a `DMCOMPOSITE`.
120047c6ae99SBarry Smith 
120147c6ae99SBarry Smith     Not Collective
120247c6ae99SBarry Smith 
120347c6ae99SBarry Smith     Input Parameter:
1204dce8aebaSBarry Smith .    dm - the `DMCOMPOSITE` object
120547c6ae99SBarry Smith 
120647c6ae99SBarry Smith     Output Parameter:
1207dce8aebaSBarry Smith .   DM ... - the individual entries `DM`
120847c6ae99SBarry Smith 
120947c6ae99SBarry Smith     Level: advanced
121047c6ae99SBarry Smith 
1211dce8aebaSBarry Smith     Fortran Note:
1212dce8aebaSBarry Smith     Available as `DMCompositeGetEntries()` for one output `DM`, DMCompositeGetEntries2() for 2, etc
1213f5f57ec0SBarry Smith 
1214dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()`
1215db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1216db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`,
1217db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`
121847c6ae99SBarry Smith @*/
1219d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntries(DM dm, ...)
1220d71ae5a4SJacob Faibussowitsch {
122147c6ae99SBarry Smith   va_list                 Argp;
122247c6ae99SBarry Smith   struct DMCompositeLink *next;
122347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
122471b14b3eSStefano Zampini   PetscBool               flg;
122547c6ae99SBarry Smith 
122647c6ae99SBarry Smith   PetscFunctionBegin;
122747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12289566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
12297a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
123047c6ae99SBarry Smith   next = com->next;
123147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
123247c6ae99SBarry Smith   va_start(Argp, dm);
123347c6ae99SBarry Smith   while (next) {
123447c6ae99SBarry Smith     DM *dmn;
123547c6ae99SBarry Smith     dmn = va_arg(Argp, DM *);
12369ae5db72SJed Brown     if (dmn) *dmn = next->dm;
123747c6ae99SBarry Smith     next = next->next;
123847c6ae99SBarry Smith   }
123947c6ae99SBarry Smith   va_end(Argp);
12403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124147c6ae99SBarry Smith }
124247c6ae99SBarry Smith 
1243dbab29e1SMark F. Adams /*@C
1244dce8aebaSBarry Smith     DMCompositeGetEntriesArray - Gets the DM for each entry in a `DMCOMPOSITE`
12452fa5ba8aSJed Brown 
12462fa5ba8aSJed Brown     Not Collective
12472fa5ba8aSJed Brown 
12482fa5ba8aSJed Brown     Input Parameter:
1249dce8aebaSBarry Smith .    dm - the `DMCOMPOSITE` object
1250907376e6SBarry Smith 
1251907376e6SBarry Smith     Output Parameter:
1252dce8aebaSBarry Smith .    dms - array of sufficient length (see `DMCompositeGetNumberDM()`) to hold the individual `DM`
12532fa5ba8aSJed Brown 
12542fa5ba8aSJed Brown     Level: advanced
12552fa5ba8aSJed Brown 
1256dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()`
1257db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1258db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`,
1259db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`
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);
12722fa5ba8aSJed Brown   /* loop over packed objects, handling one at at 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 
13119566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Private(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));
13141baa6e33SBarry Smith     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));
13429566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisSetDM_Private(ctx->subv[i], (PetscObject)dms[i]));
13439566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Private(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 
13539566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Private(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 
1367d71ae5a4SJacob Faibussowitsch 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 
138047c6ae99SBarry Smith   /* loop over packed objects, handling one at at 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 
1390d71ae5a4SJacob Faibussowitsch 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 
140314354c39SJed Brown   /* loop over packed objects, handling one at at 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 
1413d71ae5a4SJacob Faibussowitsch 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 
144347c6ae99SBarry Smith   /* loop over packed objects, handling one at at 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 
1474d71ae5a4SJacob Faibussowitsch 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 
1514d71ae5a4SJacob Faibussowitsch 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 
152947c6ae99SBarry Smith   /* loop over packed objects, handling one at at 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 
1557d71ae5a4SJacob Faibussowitsch 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 
1566d71ae5a4SJacob Faibussowitsch 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 
158239d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at 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));
160539d35262SVincent Le Chenadec 
16063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
160739d35262SVincent Le Chenadec }
160839d35262SVincent Le Chenadec 
1609d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1610d71ae5a4SJacob Faibussowitsch {
161139d35262SVincent Le Chenadec   PetscFunctionBegin;
161239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
161339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
161439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
161639d35262SVincent Le Chenadec }
161739d35262SVincent Le Chenadec 
1618d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2)
1619d71ae5a4SJacob Faibussowitsch {
162039d35262SVincent Le Chenadec   struct DMCompositeLink *next;
162139d35262SVincent Le Chenadec   PetscScalar            *array1, *array2;
162239d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite *)dm->data;
162339d35262SVincent Le Chenadec 
162439d35262SVincent Le Chenadec   PetscFunctionBegin;
162539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
162639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2);
162739d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4);
162839d35262SVincent Le Chenadec 
162948a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
163039d35262SVincent Le Chenadec 
16319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec1, &array1));
16329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec2, &array2));
163339d35262SVincent Le Chenadec 
163439d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
163539d35262SVincent Le Chenadec   next = com->next;
163639d35262SVincent Le Chenadec   while (next) {
163739d35262SVincent Le Chenadec     Vec local1, local2;
163839d35262SVincent Le Chenadec 
16399566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local1));
16409566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local1, array1));
16419566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local2));
16429566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local2, array2));
16439566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2));
16449566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2));
16459566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local2));
16469566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local2));
16479566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local1));
16489566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local1));
164939d35262SVincent Le Chenadec 
165039d35262SVincent Le Chenadec     array1 += next->nlocal;
165139d35262SVincent Le Chenadec     array2 += next->nlocal;
165239d35262SVincent Le Chenadec     next = next->next;
165339d35262SVincent Le Chenadec   }
165439d35262SVincent Le Chenadec 
16559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec1, NULL));
16569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec2, NULL));
165739d35262SVincent Le Chenadec 
16583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165939d35262SVincent Le Chenadec }
166039d35262SVincent Le Chenadec 
1661d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1662d71ae5a4SJacob Faibussowitsch {
166339d35262SVincent Le Chenadec   PetscFunctionBegin;
166439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
166539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
166639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16680c010503SBarry Smith }
166947c6ae99SBarry Smith 
16706ae3a549SBarry Smith /*MC
1671dce8aebaSBarry Smith    DMCOMPOSITE = "composite" - A `DM` object that is used to manage data for a collection of `DM`
16726ae3a549SBarry Smith 
16736ae3a549SBarry Smith   Level: intermediate
16746ae3a549SBarry Smith 
1675db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()`
16766ae3a549SBarry Smith M*/
16776ae3a549SBarry Smith 
1678d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1679d71ae5a4SJacob Faibussowitsch {
1680a4121054SBarry Smith   DM_Composite *com;
1681a4121054SBarry Smith 
1682a4121054SBarry Smith   PetscFunctionBegin;
16834dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&com));
1684a4121054SBarry Smith   p->data     = com;
1685a4121054SBarry Smith   com->n      = 0;
16867ac2b803SAlex Fikl   com->nghost = 0;
16870298fd71SBarry Smith   com->next   = NULL;
1688a4121054SBarry Smith   com->nDM    = 0;
1689a4121054SBarry Smith 
1690a4121054SBarry Smith   p->ops->createglobalvector       = DMCreateGlobalVector_Composite;
1691a4121054SBarry Smith   p->ops->createlocalvector        = DMCreateLocalVector_Composite;
1692184d77edSJed Brown   p->ops->getlocaltoglobalmapping  = DMGetLocalToGlobalMapping_Composite;
16934d343eeaSMatthew G Knepley   p->ops->createfieldis            = DMCreateFieldIS_Composite;
169416621825SDmitry Karpeev   p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite;
1695a4121054SBarry Smith   p->ops->refine                   = DMRefine_Composite;
169614354c39SJed Brown   p->ops->coarsen                  = DMCoarsen_Composite;
169725296bd5SBarry Smith   p->ops->createinterpolation      = DMCreateInterpolation_Composite;
169825296bd5SBarry Smith   p->ops->creatematrix             = DMCreateMatrix_Composite;
1699e727c939SJed Brown   p->ops->getcoloring              = DMCreateColoring_Composite;
1700a4121054SBarry Smith   p->ops->globaltolocalbegin       = DMGlobalToLocalBegin_Composite;
1701a4121054SBarry Smith   p->ops->globaltolocalend         = DMGlobalToLocalEnd_Composite;
170239d35262SVincent Le Chenadec   p->ops->localtoglobalbegin       = DMLocalToGlobalBegin_Composite;
170339d35262SVincent Le Chenadec   p->ops->localtoglobalend         = DMLocalToGlobalEnd_Composite;
170439d35262SVincent Le Chenadec   p->ops->localtolocalbegin        = DMLocalToLocalBegin_Composite;
170539d35262SVincent Le Chenadec   p->ops->localtolocalend          = DMLocalToLocalEnd_Composite;
1706a4121054SBarry Smith   p->ops->destroy                  = DMDestroy_Composite;
1707a4121054SBarry Smith   p->ops->view                     = DMView_Composite;
1708a4121054SBarry Smith   p->ops->setup                    = DMSetUp_Composite;
1709e10fd815SStefano Zampini 
17109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite));
17113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1712a4121054SBarry Smith }
1713a4121054SBarry Smith 
17141fd49c25SBarry Smith /*@
1715dce8aebaSBarry Smith     DMCompositeCreate - Creates a `DMCOMPOSITE`, used to generate "composite"
17160c010503SBarry Smith       vectors made up of several subvectors.
17170c010503SBarry Smith 
1718d083f849SBarry Smith     Collective
171947c6ae99SBarry Smith 
172047c6ae99SBarry Smith     Input Parameter:
17210c010503SBarry Smith .   comm - the processors that will share the global vector
17220c010503SBarry Smith 
17232fe279fdSBarry Smith     Output Parameter:
1724dce8aebaSBarry Smith .   packer - the `DMCOMPOSITE` object
172547c6ae99SBarry Smith 
172647c6ae99SBarry Smith     Level: advanced
172747c6ae99SBarry Smith 
1728dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCOMPOSITE`, `DMCreate()`
1729db781477SPatrick Sanan           `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`
1730db781477SPatrick Sanan           `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
173147c6ae99SBarry Smith @*/
1732d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer)
1733d71ae5a4SJacob Faibussowitsch {
173447c6ae99SBarry Smith   PetscFunctionBegin;
17350c010503SBarry Smith   PetscValidPointer(packer, 2);
17369566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, packer));
17379566063dSJacob Faibussowitsch   PetscCall(DMSetType(*packer, DMCOMPOSITE));
17383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173947c6ae99SBarry Smith }
1740