xref: /petsc/src/dm/impls/composite/pack.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
147c6ae99SBarry Smith 
2ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h> /*I  "petscdmcomposite.h"  I*/
3af0996ceSBarry Smith #include <petsc/private/isimpl.h>
4e10fd815SStefano Zampini #include <petsc/private/glvisviewerimpl.h>
52764a2aaSMatthew G. Knepley #include <petscds.h>
647c6ae99SBarry Smith 
747c6ae99SBarry Smith /*@C
847c6ae99SBarry Smith     DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
9bebe2cf6SSatish Balay       separate components (DMs) in a DMto build the correct matrix nonzero structure.
1047c6ae99SBarry Smith 
11d083f849SBarry Smith     Logically Collective
1247c6ae99SBarry Smith 
13d8d19677SJose E. Roman     Input Parameters:
1447c6ae99SBarry Smith +   dm - the composite object
1547c6ae99SBarry Smith -   formcouplelocations - routine to set the nonzero locations in the matrix
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith     Level: advanced
1847c6ae99SBarry Smith 
19f5f57ec0SBarry Smith     Not available from Fortran
20f5f57ec0SBarry Smith 
2195452b02SPatrick Sanan     Notes:
2295452b02SPatrick Sanan     See DMSetApplicationContext() and DMGetApplicationContext() for how to get user information into
2347c6ae99SBarry Smith         this routine
2447c6ae99SBarry Smith 
2547c6ae99SBarry Smith @*/
269371c9d4SSatish Balay PetscErrorCode DMCompositeSetCoupling(DM dm, PetscErrorCode (*FormCoupleLocations)(DM, Mat, PetscInt *, PetscInt *, PetscInt, PetscInt, PetscInt, PetscInt)) {
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;
3447c6ae99SBarry Smith   PetscFunctionReturn(0);
3547c6ae99SBarry Smith }
3647c6ae99SBarry Smith 
379371c9d4SSatish Balay PetscErrorCode DMDestroy_Composite(DM dm) {
3847c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
3947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
4047c6ae99SBarry Smith 
4147c6ae99SBarry Smith   PetscFunctionBegin;
4247c6ae99SBarry Smith   next = com->next;
4347c6ae99SBarry Smith   while (next) {
4447c6ae99SBarry Smith     prev = next;
4547c6ae99SBarry Smith     next = next->next;
469566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&prev->dm));
479566063dSJacob Faibussowitsch     PetscCall(PetscFree(prev->grstarts));
489566063dSJacob Faibussowitsch     PetscCall(PetscFree(prev));
4947c6ae99SBarry Smith   }
509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL));
51435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
529566063dSJacob Faibussowitsch   PetscCall(PetscFree(com));
5347c6ae99SBarry Smith   PetscFunctionReturn(0);
5447c6ae99SBarry Smith }
5547c6ae99SBarry Smith 
569371c9d4SSatish Balay PetscErrorCode DMView_Composite(DM dm, PetscViewer v) {
5747c6ae99SBarry Smith   PetscBool     iascii;
5847c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
5947c6ae99SBarry Smith 
6047c6ae99SBarry Smith   PetscFunctionBegin;
619566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii));
6247c6ae99SBarry Smith   if (iascii) {
6347c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
6447c6ae99SBarry Smith     PetscInt                i;
6547c6ae99SBarry Smith 
669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "DM (%s)\n", ((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix"));
6763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "  contains %" PetscInt_FMT " DMs\n", com->nDM));
689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(v));
6947c6ae99SBarry Smith     for (i = 0; lnk; lnk = lnk->next, i++) {
7063a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(v, "Link %" PetscInt_FMT ": DM of type %s\n", i, ((PetscObject)lnk->dm)->type_name));
719566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(v));
729566063dSJacob Faibussowitsch       PetscCall(DMView(lnk->dm, v));
739566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(v));
7447c6ae99SBarry Smith     }
759566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(v));
7647c6ae99SBarry Smith   }
7747c6ae99SBarry Smith   PetscFunctionReturn(0);
7847c6ae99SBarry Smith }
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
819371c9d4SSatish Balay PetscErrorCode DMSetUp_Composite(DM dm) {
8247c6ae99SBarry Smith   PetscInt                nprev = 0;
8347c6ae99SBarry Smith   PetscMPIInt             rank, size;
8447c6ae99SBarry Smith   DM_Composite           *com  = (DM_Composite *)dm->data;
8547c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
8647c6ae99SBarry Smith   PetscLayout             map;
8747c6ae99SBarry Smith 
8847c6ae99SBarry Smith   PetscFunctionBegin;
8928b400f6SJacob Faibussowitsch   PetscCheck(!com->setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Packer has already been setup");
909566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &map));
919566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(map, com->n));
929566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(map, PETSC_DETERMINE));
939566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(map, 1));
949566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(map));
959566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetSize(map, &com->N));
969566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(map, &com->rstart, NULL));
979566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&map));
9847c6ae99SBarry Smith 
999ae5db72SJed Brown   /* now set the rstart for each linked vector */
1009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
10247c6ae99SBarry Smith   while (next) {
10347c6ae99SBarry Smith     next->rstart = nprev;
10406ebdd98SJed Brown     nprev += next->n;
10547c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
1069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size, &next->grstarts));
1079566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgather(&next->grstart, 1, MPIU_INT, next->grstarts, 1, MPIU_INT, PetscObjectComm((PetscObject)dm)));
10847c6ae99SBarry Smith     next = next->next;
10947c6ae99SBarry Smith   }
11047c6ae99SBarry Smith   com->setup = PETSC_TRUE;
11147c6ae99SBarry Smith   PetscFunctionReturn(0);
11247c6ae99SBarry Smith }
11347c6ae99SBarry Smith 
11447c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
11547c6ae99SBarry Smith 
11673e31fe2SJed Brown /*@
11747c6ae99SBarry Smith     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
11847c6ae99SBarry Smith        representation.
11947c6ae99SBarry Smith 
12047c6ae99SBarry Smith     Not Collective
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith     Input Parameter:
12347c6ae99SBarry Smith .    dm - the packer object
12447c6ae99SBarry Smith 
12547c6ae99SBarry Smith     Output Parameter:
12647c6ae99SBarry Smith .     nDM - the number of DMs
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith     Level: beginner
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith @*/
1319371c9d4SSatish Balay PetscErrorCode DMCompositeGetNumberDM(DM dm, PetscInt *nDM) {
13247c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
13371b14b3eSStefano Zampini   PetscBool     flg;
1345fd66863SKarl Rupp 
13547c6ae99SBarry Smith   PetscFunctionBegin;
13647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1379566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1387a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
13947c6ae99SBarry Smith   *nDM = com->nDM;
14047c6ae99SBarry Smith   PetscFunctionReturn(0);
14147c6ae99SBarry Smith }
14247c6ae99SBarry Smith 
14347c6ae99SBarry Smith /*@C
14447c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
14547c6ae99SBarry Smith        representation.
14647c6ae99SBarry Smith 
147d083f849SBarry Smith     Collective on dm
14847c6ae99SBarry Smith 
1499ae5db72SJed Brown     Input Parameters:
15047c6ae99SBarry Smith +    dm - the packer object
1519ae5db72SJed Brown -    gvec - the global vector
1529ae5db72SJed Brown 
1539ae5db72SJed Brown     Output Parameters:
1540298fd71SBarry Smith .    Vec* ... - the packed parallel vectors, NULL for those that are not needed
15547c6ae99SBarry Smith 
15695452b02SPatrick Sanan     Notes:
15795452b02SPatrick Sanan     Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
15847c6ae99SBarry Smith 
159f73e5cebSJed Brown     Fortran Notes:
160f73e5cebSJed Brown 
161f73e5cebSJed Brown     Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4)
162f73e5cebSJed Brown     or use the alternative interface DMCompositeGetAccessArray().
163f73e5cebSJed Brown 
16447c6ae99SBarry Smith     Level: advanced
16547c6ae99SBarry Smith 
166db781477SPatrick Sanan .seealso: `DMCompositeGetEntries()`, `DMCompositeScatter()`
16747c6ae99SBarry Smith @*/
1689371c9d4SSatish Balay PetscErrorCode DMCompositeGetAccess(DM dm, Vec gvec, ...) {
16947c6ae99SBarry Smith   va_list                 Argp;
17047c6ae99SBarry Smith   struct DMCompositeLink *next;
17147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
1725edff71fSBarry Smith   PetscInt                readonly;
17371b14b3eSStefano Zampini   PetscBool               flg;
17447c6ae99SBarry Smith 
17547c6ae99SBarry Smith   PetscFunctionBegin;
17647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
17747c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
1789566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1797a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
18047c6ae99SBarry Smith   next = com->next;
18148a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
18247c6ae99SBarry Smith 
1839566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec, &readonly));
18447c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
18547c6ae99SBarry Smith   va_start(Argp, gvec);
18647c6ae99SBarry Smith   while (next) {
18747c6ae99SBarry Smith     Vec *vec;
18847c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
1899ae5db72SJed Brown     if (vec) {
1909566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, vec));
1915edff71fSBarry Smith       if (readonly) {
1925edff71fSBarry Smith         const PetscScalar *array;
1939566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(gvec, &array));
1949566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec, array + next->rstart));
1959566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(*vec));
1969566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(gvec, &array));
1975edff71fSBarry Smith       } else {
1985edff71fSBarry Smith         PetscScalar *array;
1999566063dSJacob Faibussowitsch         PetscCall(VecGetArray(gvec, &array));
2009566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(*vec, array + next->rstart));
2019566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(gvec, &array));
20247c6ae99SBarry Smith       }
2035edff71fSBarry Smith     }
20447c6ae99SBarry Smith     next = next->next;
20547c6ae99SBarry Smith   }
20647c6ae99SBarry Smith   va_end(Argp);
20747c6ae99SBarry Smith   PetscFunctionReturn(0);
20847c6ae99SBarry Smith }
20947c6ae99SBarry Smith 
210f73e5cebSJed Brown /*@C
211f73e5cebSJed Brown     DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
212f73e5cebSJed Brown        representation.
213f73e5cebSJed Brown 
214d083f849SBarry Smith     Collective on dm
215f73e5cebSJed Brown 
216f73e5cebSJed Brown     Input Parameters:
217f73e5cebSJed Brown +    dm - the packer object
218f73e5cebSJed Brown .    pvec - packed vector
219f73e5cebSJed Brown .    nwanted - number of vectors wanted
2200298fd71SBarry Smith -    wanted - sorted array of vectors wanted, or NULL to get all vectors
221f73e5cebSJed Brown 
222f73e5cebSJed Brown     Output Parameters:
223f73e5cebSJed Brown .    vecs - array of requested global vectors (must be allocated)
224f73e5cebSJed Brown 
22595452b02SPatrick Sanan     Notes:
22695452b02SPatrick Sanan     Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them
227f73e5cebSJed Brown 
228f73e5cebSJed Brown     Level: advanced
229f73e5cebSJed Brown 
230db781477SPatrick Sanan .seealso: `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
231f73e5cebSJed Brown @*/
2329371c9d4SSatish Balay PetscErrorCode DMCompositeGetAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) {
233f73e5cebSJed Brown   struct DMCompositeLink *link;
234f73e5cebSJed Brown   PetscInt                i, wnum;
235f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
236bee642f7SBarry Smith   PetscInt                readonly;
23771b14b3eSStefano Zampini   PetscBool               flg;
238f73e5cebSJed Brown 
239f73e5cebSJed Brown   PetscFunctionBegin;
240f73e5cebSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
241f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
2429566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
2437a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
24448a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
245f73e5cebSJed Brown 
2469566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
247f73e5cebSJed Brown   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
248f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
249f73e5cebSJed Brown       Vec v;
2509566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(link->dm, &v));
251bee642f7SBarry Smith       if (readonly) {
252bee642f7SBarry Smith         const PetscScalar *array;
2539566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec, &array));
2549566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + link->rstart));
2559566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
2569566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec, &array));
257bee642f7SBarry Smith       } else {
258bee642f7SBarry Smith         PetscScalar *array;
2599566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec, &array));
2609566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + link->rstart));
2619566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec, &array));
262bee642f7SBarry Smith       }
263f73e5cebSJed Brown       vecs[wnum++] = v;
264f73e5cebSJed Brown     }
265f73e5cebSJed Brown   }
266f73e5cebSJed Brown   PetscFunctionReturn(0);
267f73e5cebSJed Brown }
268f73e5cebSJed Brown 
2697ac2b803SAlex Fikl /*@C
2707ac2b803SAlex Fikl     DMCompositeGetLocalAccessArray - Allows one to access the individual
2717ac2b803SAlex Fikl     packed vectors in their local representation.
2727ac2b803SAlex Fikl 
273d083f849SBarry Smith     Collective on dm.
2747ac2b803SAlex Fikl 
2757ac2b803SAlex Fikl     Input Parameters:
2767ac2b803SAlex Fikl +    dm - the packer object
2777ac2b803SAlex Fikl .    pvec - packed vector
2787ac2b803SAlex Fikl .    nwanted - number of vectors wanted
2797ac2b803SAlex Fikl -    wanted - sorted array of vectors wanted, or NULL to get all vectors
2807ac2b803SAlex Fikl 
2817ac2b803SAlex Fikl     Output Parameters:
2827ac2b803SAlex Fikl .    vecs - array of requested local vectors (must be allocated)
2837ac2b803SAlex Fikl 
28495452b02SPatrick Sanan     Notes:
28595452b02SPatrick Sanan     Use DMCompositeRestoreLocalAccessArray() to return the vectors
2867ac2b803SAlex Fikl     when you no longer need them.
2877ac2b803SAlex Fikl 
2887ac2b803SAlex Fikl     Level: advanced
2897ac2b803SAlex Fikl 
290db781477SPatrick Sanan .seealso: `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`,
291db781477SPatrick Sanan           `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
2927ac2b803SAlex Fikl @*/
2939371c9d4SSatish Balay PetscErrorCode DMCompositeGetLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) {
2947ac2b803SAlex Fikl   struct DMCompositeLink *link;
2957ac2b803SAlex Fikl   PetscInt                i, wnum;
2967ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite *)dm->data;
2977ac2b803SAlex Fikl   PetscInt                readonly;
2987ac2b803SAlex Fikl   PetscInt                nlocal = 0;
29971b14b3eSStefano Zampini   PetscBool               flg;
3007ac2b803SAlex Fikl 
3017ac2b803SAlex Fikl   PetscFunctionBegin;
3027ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3037ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
3049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3057a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
30648a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
3077ac2b803SAlex Fikl 
3089566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
3097ac2b803SAlex Fikl   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
3107ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
3117ac2b803SAlex Fikl       Vec v;
3129566063dSJacob Faibussowitsch       PetscCall(DMGetLocalVector(link->dm, &v));
3137ac2b803SAlex Fikl       if (readonly) {
3147ac2b803SAlex Fikl         const PetscScalar *array;
3159566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(pvec, &array));
3169566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + nlocal));
317b1c3483dSMark Adams         // this method does not make sense. The local vectors are not updated with a global-to-local and the user can not do it because it is locked
3189566063dSJacob Faibussowitsch         PetscCall(VecLockReadPush(v));
3199566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(pvec, &array));
3207ac2b803SAlex Fikl       } else {
3217ac2b803SAlex Fikl         PetscScalar *array;
3229566063dSJacob Faibussowitsch         PetscCall(VecGetArray(pvec, &array));
3239566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(v, array + nlocal));
3249566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(pvec, &array));
3257ac2b803SAlex Fikl       }
3267ac2b803SAlex Fikl       vecs[wnum++] = v;
3277ac2b803SAlex Fikl     }
3287ac2b803SAlex Fikl 
3297ac2b803SAlex Fikl     nlocal += link->nlocal;
3307ac2b803SAlex Fikl   }
3317ac2b803SAlex Fikl 
3327ac2b803SAlex Fikl   PetscFunctionReturn(0);
3337ac2b803SAlex Fikl }
3347ac2b803SAlex Fikl 
33547c6ae99SBarry Smith /*@C
336aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
33747c6ae99SBarry Smith        representation.
33847c6ae99SBarry Smith 
339d083f849SBarry Smith     Collective on dm
34047c6ae99SBarry Smith 
3419ae5db72SJed Brown     Input Parameters:
34247c6ae99SBarry Smith +    dm - the packer object
34347c6ae99SBarry Smith .    gvec - the global vector
3440298fd71SBarry Smith -    Vec* ... - the individual parallel vectors, NULL for those that are not needed
34547c6ae99SBarry Smith 
34647c6ae99SBarry Smith     Level: advanced
34747c6ae99SBarry Smith 
348db781477SPatrick Sanan .seealso `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
349db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`,
350db781477SPatrick Sanan          `DMCompositeRestoreAccess()`, `DMCompositeGetAccess()`
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith @*/
3539371c9d4SSatish Balay PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...) {
35447c6ae99SBarry Smith   va_list                 Argp;
35547c6ae99SBarry Smith   struct DMCompositeLink *next;
35647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
3575edff71fSBarry Smith   PetscInt                readonly;
35871b14b3eSStefano Zampini   PetscBool               flg;
35947c6ae99SBarry Smith 
36047c6ae99SBarry Smith   PetscFunctionBegin;
36147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36247c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
3639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3647a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
36547c6ae99SBarry Smith   next = com->next;
36648a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
36747c6ae99SBarry Smith 
3689566063dSJacob Faibussowitsch   PetscCall(VecLockGet(gvec, &readonly));
36947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
37047c6ae99SBarry Smith   va_start(Argp, gvec);
37147c6ae99SBarry Smith   while (next) {
37247c6ae99SBarry Smith     Vec *vec;
37347c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
3749ae5db72SJed Brown     if (vec) {
3759566063dSJacob Faibussowitsch       PetscCall(VecResetArray(*vec));
3761baa6e33SBarry Smith       if (readonly) PetscCall(VecLockReadPop(*vec));
3779566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, vec));
37847c6ae99SBarry Smith     }
37947c6ae99SBarry Smith     next = next->next;
38047c6ae99SBarry Smith   }
38147c6ae99SBarry Smith   va_end(Argp);
38247c6ae99SBarry Smith   PetscFunctionReturn(0);
38347c6ae99SBarry Smith }
38447c6ae99SBarry Smith 
385f73e5cebSJed Brown /*@C
386f73e5cebSJed Brown     DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray()
387f73e5cebSJed Brown 
388d083f849SBarry Smith     Collective on dm
389f73e5cebSJed Brown 
390f73e5cebSJed Brown     Input Parameters:
391f73e5cebSJed Brown +    dm - the packer object
392f73e5cebSJed Brown .    pvec - packed vector
393f73e5cebSJed Brown .    nwanted - number of vectors wanted
3940298fd71SBarry Smith .    wanted - sorted array of vectors wanted, or NULL to get all vectors
395f73e5cebSJed Brown -    vecs - array of global vectors to return
396f73e5cebSJed Brown 
397f73e5cebSJed Brown     Level: advanced
398f73e5cebSJed Brown 
399db781477SPatrick Sanan .seealso: `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
400f73e5cebSJed Brown @*/
4019371c9d4SSatish Balay PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) {
402f73e5cebSJed Brown   struct DMCompositeLink *link;
403f73e5cebSJed Brown   PetscInt                i, wnum;
404f73e5cebSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
405bee642f7SBarry Smith   PetscInt                readonly;
40671b14b3eSStefano Zampini   PetscBool               flg;
407f73e5cebSJed Brown 
408f73e5cebSJed Brown   PetscFunctionBegin;
409f73e5cebSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
410f73e5cebSJed Brown   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4127a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
41348a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
414f73e5cebSJed Brown 
4159566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
416f73e5cebSJed Brown   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
417f73e5cebSJed Brown     if (!wanted || i == wanted[wnum]) {
4189566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
41948a46eb9SPierre Jolivet       if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4209566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum]));
421f73e5cebSJed Brown       wnum++;
422f73e5cebSJed Brown     }
423f73e5cebSJed Brown   }
424f73e5cebSJed Brown   PetscFunctionReturn(0);
425f73e5cebSJed Brown }
426f73e5cebSJed Brown 
4277ac2b803SAlex Fikl /*@C
4287ac2b803SAlex Fikl     DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray().
4297ac2b803SAlex Fikl 
430d083f849SBarry Smith     Collective on dm.
4317ac2b803SAlex Fikl 
4327ac2b803SAlex Fikl     Input Parameters:
4337ac2b803SAlex Fikl +    dm - the packer object
4347ac2b803SAlex Fikl .    pvec - packed vector
4357ac2b803SAlex Fikl .    nwanted - number of vectors wanted
4367ac2b803SAlex Fikl .    wanted - sorted array of vectors wanted, or NULL to restore all vectors
4377ac2b803SAlex Fikl -    vecs - array of local vectors to return
4387ac2b803SAlex Fikl 
4397ac2b803SAlex Fikl     Level: advanced
4407ac2b803SAlex Fikl 
4417ac2b803SAlex Fikl     Notes:
4427ac2b803SAlex Fikl     nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray()
4437ac2b803SAlex Fikl     otherwise the call will fail.
4447ac2b803SAlex Fikl 
445db781477SPatrick Sanan .seealso: `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`,
446db781477SPatrick Sanan           `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`,
447db781477SPatrick Sanan           `DMCompositeScatter()`, `DMCompositeGather()`
4487ac2b803SAlex Fikl @*/
4499371c9d4SSatish Balay PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) {
4507ac2b803SAlex Fikl   struct DMCompositeLink *link;
4517ac2b803SAlex Fikl   PetscInt                i, wnum;
4527ac2b803SAlex Fikl   DM_Composite           *com = (DM_Composite *)dm->data;
4537ac2b803SAlex Fikl   PetscInt                readonly;
45471b14b3eSStefano Zampini   PetscBool               flg;
4557ac2b803SAlex Fikl 
4567ac2b803SAlex Fikl   PetscFunctionBegin;
4577ac2b803SAlex Fikl   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4587ac2b803SAlex Fikl   PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4599566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4607a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
46148a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
4627ac2b803SAlex Fikl 
4639566063dSJacob Faibussowitsch   PetscCall(VecLockGet(pvec, &readonly));
4647ac2b803SAlex Fikl   for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
4657ac2b803SAlex Fikl     if (!wanted || i == wanted[wnum]) {
4669566063dSJacob Faibussowitsch       PetscCall(VecResetArray(vecs[wnum]));
46748a46eb9SPierre Jolivet       if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4689566063dSJacob Faibussowitsch       PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum]));
4697ac2b803SAlex Fikl       wnum++;
4707ac2b803SAlex Fikl     }
4717ac2b803SAlex Fikl   }
4727ac2b803SAlex Fikl   PetscFunctionReturn(0);
4737ac2b803SAlex Fikl }
4747ac2b803SAlex Fikl 
47547c6ae99SBarry Smith /*@C
47647c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
47747c6ae99SBarry Smith 
478d083f849SBarry Smith     Collective on dm
47947c6ae99SBarry Smith 
4809ae5db72SJed Brown     Input Parameters:
48147c6ae99SBarry Smith +    dm - the packer object
48247c6ae99SBarry Smith .    gvec - the global vector
4830298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for those that are not needed
48447c6ae99SBarry Smith 
48547c6ae99SBarry Smith     Level: advanced
48647c6ae99SBarry Smith 
4876f3c3dcfSJed Brown     Notes:
4886f3c3dcfSJed Brown     DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is
4896f3c3dcfSJed Brown     accessible from Fortran.
4906f3c3dcfSJed Brown 
491db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
492db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
493db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
494db781477SPatrick Sanan          `DMCompositeScatterArray()`
49547c6ae99SBarry Smith 
49647c6ae99SBarry Smith @*/
4979371c9d4SSatish Balay PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...) {
49847c6ae99SBarry Smith   va_list                 Argp;
49947c6ae99SBarry Smith   struct DMCompositeLink *next;
500c599c493SJunchao Zhang   PETSC_UNUSED PetscInt   cnt;
50147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
50271b14b3eSStefano Zampini   PetscBool               flg;
50347c6ae99SBarry Smith 
50447c6ae99SBarry Smith   PetscFunctionBegin;
50547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
50647c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5079566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5087a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
50948a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
51247c6ae99SBarry Smith   va_start(Argp, gvec);
5138fd8f222SJed Brown   for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
5149ae5db72SJed Brown     Vec local;
5159ae5db72SJed Brown     local = va_arg(Argp, Vec);
5169ae5db72SJed Brown     if (local) {
5179ae5db72SJed Brown       Vec                global;
5185edff71fSBarry Smith       const PetscScalar *array;
51963a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
5209566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
5219566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
5229566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
5239566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local));
5249566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local));
5259566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
5269566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5279566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
52847c6ae99SBarry Smith     }
52947c6ae99SBarry Smith   }
53047c6ae99SBarry Smith   va_end(Argp);
53147c6ae99SBarry Smith   PetscFunctionReturn(0);
53247c6ae99SBarry Smith }
53347c6ae99SBarry Smith 
5346f3c3dcfSJed Brown /*@
5356f3c3dcfSJed Brown     DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5366f3c3dcfSJed Brown 
537d083f849SBarry Smith     Collective on dm
5386f3c3dcfSJed Brown 
5396f3c3dcfSJed Brown     Input Parameters:
5406f3c3dcfSJed Brown +    dm - the packer object
5416f3c3dcfSJed Brown .    gvec - the global vector
542a2b725a8SWilliam Gropp -    lvecs - array of local vectors, NULL for any that are not needed
5436f3c3dcfSJed Brown 
5446f3c3dcfSJed Brown     Level: advanced
5456f3c3dcfSJed Brown 
5466f3c3dcfSJed Brown     Note:
547907376e6SBarry Smith     This is a non-variadic alternative to DMCompositeScatter()
5486f3c3dcfSJed Brown 
549db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`
550db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
551db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
5526f3c3dcfSJed Brown 
5536f3c3dcfSJed Brown @*/
5549371c9d4SSatish Balay PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs) {
5556f3c3dcfSJed Brown   struct DMCompositeLink *next;
5566f3c3dcfSJed Brown   PetscInt                i;
5576f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
55871b14b3eSStefano Zampini   PetscBool               flg;
5596f3c3dcfSJed Brown 
5606f3c3dcfSJed Brown   PetscFunctionBegin;
5616f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5626f3c3dcfSJed Brown   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5647a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
56548a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
5666f3c3dcfSJed Brown 
5676f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
5686f3c3dcfSJed Brown   for (i = 0, next = com->next; next; next = next->next, i++) {
5696f3c3dcfSJed Brown     if (lvecs[i]) {
5706f3c3dcfSJed Brown       Vec                global;
571c5d31e75SLisandro Dalcin       const PetscScalar *array;
5726f3c3dcfSJed Brown       PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3);
5739566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
5749566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(gvec, &array));
5759566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart));
5769566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i]));
5779566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i]));
5789566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(gvec, &array));
5799566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
5809566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
5816f3c3dcfSJed Brown     }
5826f3c3dcfSJed Brown   }
5836f3c3dcfSJed Brown   PetscFunctionReturn(0);
5846f3c3dcfSJed Brown }
5856f3c3dcfSJed Brown 
58647c6ae99SBarry Smith /*@C
58747c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
58847c6ae99SBarry Smith 
589d083f849SBarry Smith     Collective on dm
59047c6ae99SBarry Smith 
591d8d19677SJose E. Roman     Input Parameters:
59247c6ae99SBarry Smith +    dm - the packer object
59347c6ae99SBarry Smith .    gvec - the global vector
594907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
5950298fd71SBarry Smith -    Vec ... - the individual sequential vectors, NULL for any that are not needed
59647c6ae99SBarry Smith 
59747c6ae99SBarry Smith     Level: advanced
59847c6ae99SBarry Smith 
599f5f57ec0SBarry Smith     Not available from Fortran, Fortran users can use DMCompositeGatherArray()
600f5f57ec0SBarry Smith 
601db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
602db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
603db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
60447c6ae99SBarry Smith 
60547c6ae99SBarry Smith @*/
6069371c9d4SSatish Balay PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...) {
60747c6ae99SBarry Smith   va_list                 Argp;
60847c6ae99SBarry Smith   struct DMCompositeLink *next;
60947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
610c599c493SJunchao Zhang   PETSC_UNUSED PetscInt   cnt;
61171b14b3eSStefano Zampini   PetscBool               flg;
61247c6ae99SBarry Smith 
61347c6ae99SBarry Smith   PetscFunctionBegin;
61447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
615064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6169566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6177a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
61848a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
61947c6ae99SBarry Smith 
62047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
6211dac896bSSatish Balay   va_start(Argp, gvec);
6228fd8f222SJed Brown   for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
6239ae5db72SJed Brown     Vec local;
6249ae5db72SJed Brown     local = va_arg(Argp, Vec);
6259ae5db72SJed Brown     if (local) {
62647c6ae99SBarry Smith       PetscScalar *array;
6279ae5db72SJed Brown       Vec          global;
62863a3b9bcSJacob Faibussowitsch       PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
6299566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
6309566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec, &array));
6319566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
6329566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global));
6339566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global));
6349566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec, &array));
6359566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6369566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
63747c6ae99SBarry Smith     }
63847c6ae99SBarry Smith   }
63947c6ae99SBarry Smith   va_end(Argp);
64047c6ae99SBarry Smith   PetscFunctionReturn(0);
64147c6ae99SBarry Smith }
64247c6ae99SBarry Smith 
6436f3c3dcfSJed Brown /*@
6446f3c3dcfSJed Brown     DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6456f3c3dcfSJed Brown 
646d083f849SBarry Smith     Collective on dm
6476f3c3dcfSJed Brown 
648d8d19677SJose E. Roman     Input Parameters:
6496f3c3dcfSJed Brown +    dm - the packer object
6506f3c3dcfSJed Brown .    gvec - the global vector
651907376e6SBarry Smith .    imode - INSERT_VALUES or ADD_VALUES
6526f3c3dcfSJed Brown -    lvecs - the individual sequential vectors, NULL for any that are not needed
6536f3c3dcfSJed Brown 
6546f3c3dcfSJed Brown     Level: advanced
6556f3c3dcfSJed Brown 
6566f3c3dcfSJed Brown     Notes:
6576f3c3dcfSJed Brown     This is a non-variadic alternative to DMCompositeGather().
6586f3c3dcfSJed Brown 
659db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
660db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
661db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`,
6626f3c3dcfSJed Brown @*/
6639371c9d4SSatish Balay PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs) {
6646f3c3dcfSJed Brown   struct DMCompositeLink *next;
6656f3c3dcfSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
6666f3c3dcfSJed Brown   PetscInt                i;
66771b14b3eSStefano Zampini   PetscBool               flg;
6686f3c3dcfSJed Brown 
6696f3c3dcfSJed Brown   PetscFunctionBegin;
6706f3c3dcfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
671064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6737a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
67448a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
6756f3c3dcfSJed Brown 
6766f3c3dcfSJed Brown   /* loop over packed objects, handling one at at time */
6776f3c3dcfSJed Brown   for (next = com->next, i = 0; next; next = next->next, i++) {
6786f3c3dcfSJed Brown     if (lvecs[i]) {
6796f3c3dcfSJed Brown       PetscScalar *array;
6806f3c3dcfSJed Brown       Vec          global;
681064a246eSJacob Faibussowitsch       PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4);
6829566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(next->dm, &global));
6839566063dSJacob Faibussowitsch       PetscCall(VecGetArray(gvec, &array));
6849566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(global, array + next->rstart));
6859566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global));
6869566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global));
6879566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(gvec, &array));
6889566063dSJacob Faibussowitsch       PetscCall(VecResetArray(global));
6899566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(next->dm, &global));
6906f3c3dcfSJed Brown     }
6916f3c3dcfSJed Brown   }
6926f3c3dcfSJed Brown   PetscFunctionReturn(0);
6936f3c3dcfSJed Brown }
6946f3c3dcfSJed Brown 
695f5f57ec0SBarry Smith /*@
696aa219208SBarry Smith     DMCompositeAddDM - adds a DM vector to a DMComposite
69747c6ae99SBarry Smith 
698d083f849SBarry Smith     Collective on dm
69947c6ae99SBarry Smith 
700d8d19677SJose E. Roman     Input Parameters:
7019b52a9b5SPatrick Sanan +    dmc - the DMComposite (packer) object
702f5f57ec0SBarry Smith -    dm - the DM object
70347c6ae99SBarry Smith 
70447c6ae99SBarry Smith     Level: advanced
70547c6ae99SBarry Smith 
706db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeGather()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
707db781477SPatrick Sanan          `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
708db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
70947c6ae99SBarry Smith 
71047c6ae99SBarry Smith @*/
7119371c9d4SSatish Balay PetscErrorCode DMCompositeAddDM(DM dmc, DM dm) {
71206ebdd98SJed Brown   PetscInt                n, nlocal;
71347c6ae99SBarry Smith   struct DMCompositeLink *mine, *next;
71406ebdd98SJed Brown   Vec                     global, local;
71547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dmc->data;
71671b14b3eSStefano Zampini   PetscBool               flg;
71747c6ae99SBarry Smith 
71847c6ae99SBarry Smith   PetscFunctionBegin;
71947c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc, DM_CLASSID, 1);
72047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
7219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg));
7227a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
72347c6ae99SBarry Smith   next = com->next;
72428b400f6SJacob Faibussowitsch   PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite");
72547c6ae99SBarry Smith 
72647c6ae99SBarry Smith   /* create new link */
7279566063dSJacob Faibussowitsch   PetscCall(PetscNew(&mine));
7289566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
7299566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &global));
7309566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global, &n));
7319566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &global));
7329566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &local));
7339566063dSJacob Faibussowitsch   PetscCall(VecGetSize(local, &nlocal));
7349566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &local));
7358865f1eaSKarl Rupp 
73647c6ae99SBarry Smith   mine->n      = n;
73706ebdd98SJed Brown   mine->nlocal = nlocal;
73847c6ae99SBarry Smith   mine->dm     = dm;
7390298fd71SBarry Smith   mine->next   = NULL;
74047c6ae99SBarry Smith   com->n += n;
7417ac2b803SAlex Fikl   com->nghost += nlocal;
74247c6ae99SBarry Smith 
74347c6ae99SBarry Smith   /* add to end of list */
7448865f1eaSKarl Rupp   if (!next) com->next = mine;
7458865f1eaSKarl Rupp   else {
74647c6ae99SBarry Smith     while (next->next) next = next->next;
74747c6ae99SBarry Smith     next->next = mine;
74847c6ae99SBarry Smith   }
74947c6ae99SBarry Smith   com->nDM++;
75047c6ae99SBarry Smith   com->nmine++;
75147c6ae99SBarry Smith   PetscFunctionReturn(0);
75247c6ae99SBarry Smith }
75347c6ae99SBarry Smith 
7549804daf3SBarry Smith #include <petscdraw.h>
75526887b52SJed Brown PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
7569371c9d4SSatish Balay PetscErrorCode              VecView_DMComposite(Vec gvec, PetscViewer viewer) {
75747c6ae99SBarry Smith                DM                      dm;
75847c6ae99SBarry Smith                struct DMCompositeLink *next;
75947c6ae99SBarry Smith                PetscBool               isdraw;
760cef07954SSatish Balay                DM_Composite           *com;
76147c6ae99SBarry Smith 
76247c6ae99SBarry Smith                PetscFunctionBegin;
7639566063dSJacob Faibussowitsch                PetscCall(VecGetDM(gvec, &dm));
7647a8be351SBarry Smith                PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite");
76547c6ae99SBarry Smith                com  = (DM_Composite *)dm->data;
76647c6ae99SBarry Smith                next = com->next;
76747c6ae99SBarry Smith 
7689566063dSJacob Faibussowitsch                PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
76947c6ae99SBarry Smith                if (!isdraw) {
77047c6ae99SBarry Smith                  /* do I really want to call this? */
7719566063dSJacob Faibussowitsch     PetscCall(VecView_MPI(gvec, viewer));
77247c6ae99SBarry Smith   } else {
77347c6ae99SBarry Smith                  PetscInt cnt = 0;
77447c6ae99SBarry Smith 
77547c6ae99SBarry Smith                  /* loop over packed objects, handling one at at time */
77647c6ae99SBarry Smith                  while (next) {
77747c6ae99SBarry Smith                    Vec                vec;
778ca4278abSLisandro Dalcin                    const PetscScalar *array;
77947c6ae99SBarry Smith                    PetscInt           bs;
78047c6ae99SBarry Smith 
7819ae5db72SJed Brown                    /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
7829566063dSJacob Faibussowitsch                    PetscCall(DMGetGlobalVector(next->dm, &vec));
7839566063dSJacob Faibussowitsch                    PetscCall(VecGetArrayRead(gvec, &array));
7849566063dSJacob Faibussowitsch                    PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart));
7859566063dSJacob Faibussowitsch                    PetscCall(VecRestoreArrayRead(gvec, &array));
7869566063dSJacob Faibussowitsch                    PetscCall(VecView(vec, viewer));
7879566063dSJacob Faibussowitsch                    PetscCall(VecResetArray(vec));
7889566063dSJacob Faibussowitsch                    PetscCall(VecGetBlockSize(vec, &bs));
7899566063dSJacob Faibussowitsch                    PetscCall(DMRestoreGlobalVector(next->dm, &vec));
7909566063dSJacob Faibussowitsch                    PetscCall(PetscViewerDrawBaseAdd(viewer, bs));
79147c6ae99SBarry Smith                    cnt += bs;
79247c6ae99SBarry Smith                    next = next->next;
79347c6ae99SBarry Smith     }
7949566063dSJacob Faibussowitsch                  PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt));
79547c6ae99SBarry Smith   }
79647c6ae99SBarry Smith                PetscFunctionReturn(0);
79747c6ae99SBarry Smith }
79847c6ae99SBarry Smith 
7999371c9d4SSatish Balay PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec) {
80047c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
80147c6ae99SBarry Smith 
80247c6ae99SBarry Smith   PetscFunctionBegin;
80347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8049566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
8059566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8069566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec));
8079566063dSJacob Faibussowitsch   PetscCall(VecSetType(*gvec, dm->vectype));
8089566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*gvec, com->n, com->N));
8099566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*gvec, dm));
8109566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (void (*)(void))VecView_DMComposite));
81147c6ae99SBarry Smith   PetscFunctionReturn(0);
81247c6ae99SBarry Smith }
81347c6ae99SBarry Smith 
8149371c9d4SSatish Balay PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec) {
81547c6ae99SBarry Smith   DM_Composite *com = (DM_Composite *)dm->data;
81647c6ae99SBarry Smith 
81747c6ae99SBarry Smith   PetscFunctionBegin;
81847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
81947c6ae99SBarry Smith   if (!com->setup) {
8209566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(dm));
8219566063dSJacob Faibussowitsch     PetscCall(DMSetUp(dm));
82247c6ae99SBarry Smith   }
8239566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, lvec));
8249566063dSJacob Faibussowitsch   PetscCall(VecSetType(*lvec, dm->vectype));
8259566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE));
8269566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*lvec, dm));
82747c6ae99SBarry Smith   PetscFunctionReturn(0);
82847c6ae99SBarry Smith }
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith /*@C
8319ae5db72SJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space
83247c6ae99SBarry Smith 
83306ebdd98SJed Brown     Collective on DM
83447c6ae99SBarry Smith 
83547c6ae99SBarry Smith     Input Parameter:
83647c6ae99SBarry Smith .    dm - the packer object
83747c6ae99SBarry Smith 
83847c6ae99SBarry Smith     Output Parameters:
8399ae5db72SJed Brown .    ltogs - the individual mappings for each packed vector. Note that this includes
8409ae5db72SJed Brown            all the ghost points that individual ghosted DMDA's may have.
84147c6ae99SBarry Smith 
84247c6ae99SBarry Smith     Level: advanced
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith     Notes:
8456eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
84647c6ae99SBarry Smith 
847f5f57ec0SBarry Smith     Not available from Fortran
848f5f57ec0SBarry Smith 
849db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
850db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
851c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
85247c6ae99SBarry Smith 
85347c6ae99SBarry Smith @*/
8549371c9d4SSatish Balay PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping **ltogs) {
85547c6ae99SBarry Smith   PetscInt                i, *idx, n, cnt;
85647c6ae99SBarry Smith   struct DMCompositeLink *next;
85747c6ae99SBarry Smith   PetscMPIInt             rank;
85847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
85971b14b3eSStefano Zampini   PetscBool               flg;
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith   PetscFunctionBegin;
86247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
8647a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
8659566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
8669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM, ltogs));
86747c6ae99SBarry Smith   next = com->next;
8689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
86947c6ae99SBarry Smith 
87047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
87147c6ae99SBarry Smith   cnt = 0;
87247c6ae99SBarry Smith   while (next) {
8736eb61c8cSJed Brown     ISLocalToGlobalMapping ltog;
8746eb61c8cSJed Brown     PetscMPIInt            size;
87586994e45SJed Brown     const PetscInt        *suboff, *indices;
8766eb61c8cSJed Brown     Vec                    global;
87747c6ae99SBarry Smith 
8786eb61c8cSJed Brown     /* Get sub-DM global indices for each local dof */
8799566063dSJacob Faibussowitsch     PetscCall(DMGetLocalToGlobalMapping(next->dm, &ltog));
8809566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n));
8819566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices));
8829566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &idx));
88347c6ae99SBarry Smith 
8846eb61c8cSJed Brown     /* Get the offsets for the sub-DM global vector */
8859566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
8869566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRanges(global, &suboff));
8879566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size));
8886eb61c8cSJed Brown 
8896eb61c8cSJed Brown     /* Shift the sub-DM definition of the global space to the composite global space */
8906eb61c8cSJed Brown     for (i = 0; i < n; i++) {
89186994e45SJed Brown       PetscInt subi = indices[i], lo = 0, hi = size, t;
892e333b035SStefano Zampini       /* There's no consensus on what a negative index means,
893e333b035SStefano Zampini          except for skipping when setting the values in vectors and matrices */
8949371c9d4SSatish Balay       if (subi < 0) {
8959371c9d4SSatish Balay         idx[i] = subi - next->grstarts[rank];
8969371c9d4SSatish Balay         continue;
8979371c9d4SSatish Balay       }
8986eb61c8cSJed Brown       /* Binary search to find which rank owns subi */
8996eb61c8cSJed Brown       while (hi - lo > 1) {
9006eb61c8cSJed Brown         t = lo + (hi - lo) / 2;
9016eb61c8cSJed Brown         if (suboff[t] > subi) hi = t;
9026eb61c8cSJed Brown         else lo = t;
9036eb61c8cSJed Brown       }
9046eb61c8cSJed Brown       idx[i] = subi - suboff[lo] + next->grstarts[lo];
9056eb61c8cSJed Brown     }
9069566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices));
9079566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt]));
9089566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
90947c6ae99SBarry Smith     next = next->next;
91047c6ae99SBarry Smith     cnt++;
91147c6ae99SBarry Smith   }
91247c6ae99SBarry Smith   PetscFunctionReturn(0);
91347c6ae99SBarry Smith }
91447c6ae99SBarry Smith 
91587c85e80SJed Brown /*@C
9169ae5db72SJed Brown    DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
91787c85e80SJed Brown 
91887c85e80SJed Brown    Not Collective
91987c85e80SJed Brown 
9204165533cSJose E. Roman    Input Parameter:
92187c85e80SJed Brown . dm - composite DM
92287c85e80SJed Brown 
9234165533cSJose E. Roman    Output Parameter:
92487c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
92587c85e80SJed Brown 
92687c85e80SJed Brown    Level: intermediate
92787c85e80SJed Brown 
92887c85e80SJed Brown    Notes:
92987c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
93087c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
9319ae5db72SJed Brown    scatter to a composite local vector.  The user should not typically need to know which is being done.
93287c85e80SJed Brown 
93387c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
93487c85e80SJed Brown 
93587c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
93687c85e80SJed Brown 
93787c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
93887c85e80SJed Brown 
939f5f57ec0SBarry Smith    Not available from Fortran
940f5f57ec0SBarry Smith 
941db781477SPatrick Sanan .seealso: `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`, `MatCreateLocalRef()`
94287c85e80SJed Brown @*/
9439371c9d4SSatish Balay PetscErrorCode DMCompositeGetLocalISs(DM dm, IS **is) {
94487c85e80SJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
94587c85e80SJed Brown   struct DMCompositeLink *link;
94687c85e80SJed Brown   PetscInt                cnt, start;
94771b14b3eSStefano Zampini   PetscBool               flg;
94887c85e80SJed Brown 
94987c85e80SJed Brown   PetscFunctionBegin;
95087c85e80SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
95187c85e80SJed Brown   PetscValidPointer(is, 2);
9529566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
9537a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
9549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nmine, is));
95506ebdd98SJed Brown   for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) {
956520db06cSJed Brown     PetscInt bs;
9579566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt]));
9589566063dSJacob Faibussowitsch     PetscCall(DMGetBlockSize(link->dm, &bs));
9599566063dSJacob Faibussowitsch     PetscCall(ISSetBlockSize((*is)[cnt], bs));
960520db06cSJed Brown   }
96187c85e80SJed Brown   PetscFunctionReturn(0);
96287c85e80SJed Brown }
96387c85e80SJed Brown 
96447c6ae99SBarry Smith /*@C
96547c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
96647c6ae99SBarry Smith 
967d083f849SBarry Smith     Collective on dm
96847c6ae99SBarry Smith 
96947c6ae99SBarry Smith     Input Parameter:
97047c6ae99SBarry Smith .    dm - the packer object
97147c6ae99SBarry Smith 
97247c6ae99SBarry Smith     Output Parameters:
97347c6ae99SBarry Smith .    is - the array of index sets
97447c6ae99SBarry Smith 
97547c6ae99SBarry Smith     Level: advanced
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith     Notes:
97847c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
97947c6ae99SBarry Smith 
98047c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
98147c6ae99SBarry Smith 
9826eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
9836eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
9846eb61c8cSJed Brown        indices.
98547c6ae99SBarry Smith 
986f3cb0f7eSJed Brown     Fortran Notes:
987f3cb0f7eSJed Brown 
988f3cb0f7eSJed Brown        The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM().
989f3cb0f7eSJed Brown 
990db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
991db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
992c2e3fba1SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
99347c6ae99SBarry Smith 
99447c6ae99SBarry Smith @*/
9959371c9d4SSatish Balay PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[]) {
99666bb578eSMark Adams   PetscInt                cnt = 0;
99747c6ae99SBarry Smith   struct DMCompositeLink *next;
99847c6ae99SBarry Smith   PetscMPIInt             rank;
99947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
100071b14b3eSStefano Zampini   PetscBool               flg;
100147c6ae99SBarry Smith 
100247c6ae99SBarry Smith   PetscFunctionBegin;
100347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
10057a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
10067a8be351SBarry Smith   PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before");
10079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(com->nDM, is));
100847c6ae99SBarry Smith   next = com->next;
10099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
101047c6ae99SBarry Smith 
101147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
101247c6ae99SBarry Smith   while (next) {
1013e5e52638SMatthew G. Knepley     PetscDS prob;
1014e5e52638SMatthew G. Knepley 
10159566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt]));
10169566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &prob));
1017e5e52638SMatthew G. Knepley     if (prob) {
101865c226d8SMatthew G. Knepley       MatNullSpace space;
101965c226d8SMatthew G. Knepley       Mat          pmat;
10200f21e855SMatthew G. Knepley       PetscObject  disc;
10210f21e855SMatthew G. Knepley       PetscInt     Nf;
102265c226d8SMatthew G. Knepley 
10239566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(prob, &Nf));
1024f24dd8d2SMatthew G. Knepley       if (cnt < Nf) {
10259566063dSJacob Faibussowitsch         PetscCall(PetscDSGetDiscretization(prob, cnt, &disc));
10269566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space));
10279566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space));
10289566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space));
10299566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space));
10309566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat));
10319566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat));
103265c226d8SMatthew G. Knepley       }
1033f24dd8d2SMatthew G. Knepley     }
103447c6ae99SBarry Smith     cnt++;
103547c6ae99SBarry Smith     next = next->next;
103647c6ae99SBarry Smith   }
103747c6ae99SBarry Smith   PetscFunctionReturn(0);
103847c6ae99SBarry Smith }
103947c6ae99SBarry Smith 
10409371c9d4SSatish Balay PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) {
10414d343eeaSMatthew G Knepley   PetscInt nDM;
10424d343eeaSMatthew G Knepley   DM      *dms;
10434d343eeaSMatthew G Knepley   PetscInt i;
10444d343eeaSMatthew G Knepley 
10454d343eeaSMatthew G Knepley   PetscFunctionBegin;
10469566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm, &nDM));
10478865f1eaSKarl Rupp   if (numFields) *numFields = nDM;
10489566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetGlobalISs(dm, fields));
10494d343eeaSMatthew G Knepley   if (fieldNames) {
10509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, &dms));
10519566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, fieldNames));
10529566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, dms));
10534d343eeaSMatthew G Knepley     for (i = 0; i < nDM; i++) {
10544d343eeaSMatthew G Knepley       char        buf[256];
10554d343eeaSMatthew G Knepley       const char *splitname;
10564d343eeaSMatthew G Knepley 
10574d343eeaSMatthew G Knepley       /* Split naming precedence: object name, prefix, number */
10584d343eeaSMatthew G Knepley       splitname = ((PetscObject)dm)->name;
10594d343eeaSMatthew G Knepley       if (!splitname) {
10609566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname));
10614d343eeaSMatthew G Knepley         if (splitname) {
10624d343eeaSMatthew G Knepley           size_t len;
10639566063dSJacob Faibussowitsch           PetscCall(PetscStrncpy(buf, splitname, sizeof(buf)));
10648caf3d72SBarry Smith           buf[sizeof(buf) - 1] = 0;
10659566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(buf, &len));
10664d343eeaSMatthew G Knepley           if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */
10674d343eeaSMatthew G Knepley           splitname = buf;
10684d343eeaSMatthew G Knepley         }
10694d343eeaSMatthew G Knepley       }
10704d343eeaSMatthew G Knepley       if (!splitname) {
107163a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i));
10724d343eeaSMatthew G Knepley         splitname = buf;
10734d343eeaSMatthew G Knepley       }
10749566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i]));
10754d343eeaSMatthew G Knepley     }
10769566063dSJacob Faibussowitsch     PetscCall(PetscFree(dms));
10774d343eeaSMatthew G Knepley   }
10784d343eeaSMatthew G Knepley   PetscFunctionReturn(0);
10794d343eeaSMatthew G Knepley }
10804d343eeaSMatthew G Knepley 
1081e7c4fc90SDmitry Karpeev /*
1082e7c4fc90SDmitry Karpeev  This could take over from DMCreateFieldIS(), as it is more general,
10830298fd71SBarry Smith  making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1084e7c4fc90SDmitry Karpeev  At this point it's probably best to be less intrusive, however.
1085e7c4fc90SDmitry Karpeev  */
10869371c9d4SSatish Balay PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) {
1087e7c4fc90SDmitry Karpeev   PetscInt nDM;
1088e7c4fc90SDmitry Karpeev   PetscInt i;
1089e7c4fc90SDmitry Karpeev 
1090e7c4fc90SDmitry Karpeev   PetscFunctionBegin;
10919566063dSJacob Faibussowitsch   PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist));
1092e7c4fc90SDmitry Karpeev   if (dmlist) {
10939566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetNumberDM(dm, &nDM));
10949566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nDM, dmlist));
10959566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetEntriesArray(dm, *dmlist));
109648a46eb9SPierre Jolivet     for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i])));
1097e7c4fc90SDmitry Karpeev   }
1098e7c4fc90SDmitry Karpeev   PetscFunctionReturn(0);
1099e7c4fc90SDmitry Karpeev }
1100e7c4fc90SDmitry Karpeev 
110147c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
110247c6ae99SBarry Smith /*@C
11039ae5db72SJed Brown     DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite.
110447c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
110547c6ae99SBarry Smith 
110647c6ae99SBarry Smith     Not Collective
110747c6ae99SBarry Smith 
110847c6ae99SBarry Smith     Input Parameter:
110947c6ae99SBarry Smith .    dm - the packer object
111047c6ae99SBarry Smith 
111147c6ae99SBarry Smith     Output Parameter:
11129ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
111347c6ae99SBarry Smith 
111447c6ae99SBarry Smith     Level: advanced
111547c6ae99SBarry Smith 
1116f5f57ec0SBarry Smith     Not available from Fortran
1117f5f57ec0SBarry Smith 
1118db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1119db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1120db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
112147c6ae99SBarry Smith 
112247c6ae99SBarry Smith @*/
11239371c9d4SSatish Balay PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...) {
112447c6ae99SBarry Smith   va_list                 Argp;
112547c6ae99SBarry Smith   struct DMCompositeLink *next;
112647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
112771b14b3eSStefano Zampini   PetscBool               flg;
112847c6ae99SBarry Smith 
112947c6ae99SBarry Smith   PetscFunctionBegin;
113047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11319566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11327a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
113347c6ae99SBarry Smith   next = com->next;
113447c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
113547c6ae99SBarry Smith   va_start(Argp, dm);
113647c6ae99SBarry Smith   while (next) {
113747c6ae99SBarry Smith     Vec *vec;
113847c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
11399566063dSJacob Faibussowitsch     if (vec) PetscCall(DMGetLocalVector(next->dm, vec));
114047c6ae99SBarry Smith     next = next->next;
114147c6ae99SBarry Smith   }
114247c6ae99SBarry Smith   va_end(Argp);
114347c6ae99SBarry Smith   PetscFunctionReturn(0);
114447c6ae99SBarry Smith }
114547c6ae99SBarry Smith 
114647c6ae99SBarry Smith /*@C
11479ae5db72SJed Brown     DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite.
114847c6ae99SBarry Smith 
114947c6ae99SBarry Smith     Not Collective
115047c6ae99SBarry Smith 
115147c6ae99SBarry Smith     Input Parameter:
115247c6ae99SBarry Smith .    dm - the packer object
115347c6ae99SBarry Smith 
115447c6ae99SBarry Smith     Output Parameter:
11559ae5db72SJed Brown .   Vec ... - the individual sequential Vecs
115647c6ae99SBarry Smith 
115747c6ae99SBarry Smith     Level: advanced
115847c6ae99SBarry Smith 
1159f5f57ec0SBarry Smith     Not available from Fortran
1160f5f57ec0SBarry Smith 
1161db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1162db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1163db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
116447c6ae99SBarry Smith 
116547c6ae99SBarry Smith @*/
11669371c9d4SSatish Balay PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...) {
116747c6ae99SBarry Smith   va_list                 Argp;
116847c6ae99SBarry Smith   struct DMCompositeLink *next;
116947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
117071b14b3eSStefano Zampini   PetscBool               flg;
117147c6ae99SBarry Smith 
117247c6ae99SBarry Smith   PetscFunctionBegin;
117347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11749566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11757a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
117647c6ae99SBarry Smith   next = com->next;
117747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
117847c6ae99SBarry Smith   va_start(Argp, dm);
117947c6ae99SBarry Smith   while (next) {
118047c6ae99SBarry Smith     Vec *vec;
118147c6ae99SBarry Smith     vec = va_arg(Argp, Vec *);
11829566063dSJacob Faibussowitsch     if (vec) PetscCall(DMRestoreLocalVector(next->dm, vec));
118347c6ae99SBarry Smith     next = next->next;
118447c6ae99SBarry Smith   }
118547c6ae99SBarry Smith   va_end(Argp);
118647c6ae99SBarry Smith   PetscFunctionReturn(0);
118747c6ae99SBarry Smith }
118847c6ae99SBarry Smith 
118947c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
119047c6ae99SBarry Smith /*@C
11919ae5db72SJed Brown     DMCompositeGetEntries - Gets the DM for each entry in a DMComposite.
119247c6ae99SBarry Smith 
119347c6ae99SBarry Smith     Not Collective
119447c6ae99SBarry Smith 
119547c6ae99SBarry Smith     Input Parameter:
119647c6ae99SBarry Smith .    dm - the packer object
119747c6ae99SBarry Smith 
119847c6ae99SBarry Smith     Output Parameter:
11999ae5db72SJed Brown .   DM ... - the individual entries (DMs)
120047c6ae99SBarry Smith 
120147c6ae99SBarry Smith     Level: advanced
120247c6ae99SBarry Smith 
120395452b02SPatrick Sanan     Fortran Notes:
120495452b02SPatrick Sanan     Available as DMCompositeGetEntries() for one output DM, DMCompositeGetEntries2() for 2, etc
1205f5f57ec0SBarry Smith 
1206db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()`
1207db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1208db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`,
1209db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`
121047c6ae99SBarry Smith 
121147c6ae99SBarry Smith @*/
12129371c9d4SSatish Balay PetscErrorCode DMCompositeGetEntries(DM dm, ...) {
121347c6ae99SBarry Smith   va_list                 Argp;
121447c6ae99SBarry Smith   struct DMCompositeLink *next;
121547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
121671b14b3eSStefano Zampini   PetscBool               flg;
121747c6ae99SBarry Smith 
121847c6ae99SBarry Smith   PetscFunctionBegin;
121947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
12217a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
122247c6ae99SBarry Smith   next = com->next;
122347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
122447c6ae99SBarry Smith   va_start(Argp, dm);
122547c6ae99SBarry Smith   while (next) {
122647c6ae99SBarry Smith     DM *dmn;
122747c6ae99SBarry Smith     dmn = va_arg(Argp, DM *);
12289ae5db72SJed Brown     if (dmn) *dmn = next->dm;
122947c6ae99SBarry Smith     next = next->next;
123047c6ae99SBarry Smith   }
123147c6ae99SBarry Smith   va_end(Argp);
123247c6ae99SBarry Smith   PetscFunctionReturn(0);
123347c6ae99SBarry Smith }
123447c6ae99SBarry Smith 
1235dbab29e1SMark F. Adams /*@C
12362fa5ba8aSJed Brown     DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite.
12372fa5ba8aSJed Brown 
12382fa5ba8aSJed Brown     Not Collective
12392fa5ba8aSJed Brown 
12402fa5ba8aSJed Brown     Input Parameter:
1241907376e6SBarry Smith .    dm - the packer object
1242907376e6SBarry Smith 
1243907376e6SBarry Smith     Output Parameter:
1244907376e6SBarry Smith .    dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs
12452fa5ba8aSJed Brown 
12462fa5ba8aSJed Brown     Level: advanced
12472fa5ba8aSJed Brown 
1248db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()`
1249db781477SPatrick Sanan          `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1250db781477SPatrick Sanan          `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`,
1251db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`
12522fa5ba8aSJed Brown 
12532fa5ba8aSJed Brown @*/
12549371c9d4SSatish Balay PetscErrorCode DMCompositeGetEntriesArray(DM dm, DM dms[]) {
12552fa5ba8aSJed Brown   struct DMCompositeLink *next;
12562fa5ba8aSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
12572fa5ba8aSJed Brown   PetscInt                i;
125871b14b3eSStefano Zampini   PetscBool               flg;
12592fa5ba8aSJed Brown 
12602fa5ba8aSJed Brown   PetscFunctionBegin;
12612fa5ba8aSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
12637a8be351SBarry Smith   PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
12642fa5ba8aSJed Brown   /* loop over packed objects, handling one at at time */
12652fa5ba8aSJed Brown   for (next = com->next, i = 0; next; next = next->next, i++) dms[i] = next->dm;
12662fa5ba8aSJed Brown   PetscFunctionReturn(0);
12672fa5ba8aSJed Brown }
12682fa5ba8aSJed Brown 
1269e10fd815SStefano Zampini typedef struct {
1270e10fd815SStefano Zampini   DM           dm;
1271e10fd815SStefano Zampini   PetscViewer *subv;
1272e10fd815SStefano Zampini   Vec         *vecs;
1273e10fd815SStefano Zampini } GLVisViewerCtx;
1274e10fd815SStefano Zampini 
12759371c9d4SSatish Balay static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) {
1276e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
1277e10fd815SStefano Zampini   PetscInt        i, n;
1278e10fd815SStefano Zampini 
1279e10fd815SStefano Zampini   PetscFunctionBegin;
12809566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(ctx->dm, &n));
128148a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscViewerDestroy(&ctx->subv[i]));
12829566063dSJacob Faibussowitsch   PetscCall(PetscFree2(ctx->subv, ctx->vecs));
12839566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dm));
12849566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
1285e10fd815SStefano Zampini   PetscFunctionReturn(0);
1286e10fd815SStefano Zampini }
1287e10fd815SStefano Zampini 
12889371c9d4SSatish Balay static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) {
1289e10fd815SStefano Zampini   Vec             X   = (Vec)oX;
1290e10fd815SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
1291e10fd815SStefano Zampini   PetscInt        i, n, cumf;
1292e10fd815SStefano Zampini 
1293e10fd815SStefano Zampini   PetscFunctionBegin;
12949566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(ctx->dm, &n));
12959566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccessArray(ctx->dm, X, n, NULL, ctx->vecs));
1296e10fd815SStefano Zampini   for (i = 0, cumf = 0; i < n; i++) {
1297e10fd815SStefano Zampini     PetscErrorCode (*g2l)(PetscObject, PetscInt, PetscObject[], void *);
1298e10fd815SStefano Zampini     void    *fctx;
1299e10fd815SStefano Zampini     PetscInt nfi;
1300e10fd815SStefano Zampini 
13019566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nfi, NULL, NULL, &g2l, NULL, &fctx));
1302e10fd815SStefano Zampini     if (!nfi) continue;
13031baa6e33SBarry Smith     if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i], nfi, oXfield + cumf, fctx));
13041baa6e33SBarry Smith     else PetscCall(VecCopy(ctx->vecs[i], (Vec)(oXfield[cumf])));
1305e10fd815SStefano Zampini     cumf += nfi;
1306e10fd815SStefano Zampini   }
13079566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccessArray(ctx->dm, X, n, NULL, ctx->vecs));
1308e10fd815SStefano Zampini   PetscFunctionReturn(0);
1309e10fd815SStefano Zampini }
1310e10fd815SStefano Zampini 
13119371c9d4SSatish Balay static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer) {
1312e10fd815SStefano Zampini   DM              dm = (DM)odm, *dms;
1313e10fd815SStefano Zampini   Vec            *Ufds;
1314e10fd815SStefano Zampini   GLVisViewerCtx *ctx;
1315e10fd815SStefano Zampini   PetscInt        i, n, tnf, *sdim;
1316e10fd815SStefano Zampini   char          **fecs;
1317e10fd815SStefano Zampini 
1318e10fd815SStefano Zampini   PetscFunctionBegin;
13199566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
13209566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
1321e10fd815SStefano Zampini   ctx->dm = dm;
13229566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(dm, &n));
13239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &dms));
13249566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetEntriesArray(dm, dms));
13259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n, &ctx->subv, n, &ctx->vecs));
1326e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1327e10fd815SStefano Zampini     PetscInt nf;
1328e10fd815SStefano Zampini 
13299566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PetscObjectComm(odm), &ctx->subv[i]));
13309566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(ctx->subv[i], PETSCVIEWERGLVIS));
13319566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisSetDM_Private(ctx->subv[i], (PetscObject)dms[i]));
13329566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nf, NULL, NULL, NULL, NULL, NULL));
1333e10fd815SStefano Zampini     tnf += nf;
1334e10fd815SStefano Zampini   }
13359566063dSJacob Faibussowitsch   PetscCall(PetscFree(dms));
13369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(tnf, &fecs, tnf, &sdim, tnf, &Ufds));
1337e10fd815SStefano Zampini   for (i = 0, tnf = 0; i < n; i++) {
1338e10fd815SStefano Zampini     PetscInt    *sd, nf, f;
1339e10fd815SStefano Zampini     const char **fec;
1340e10fd815SStefano Zampini     Vec         *Uf;
1341e10fd815SStefano Zampini 
13429566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nf, &fec, &sd, NULL, (PetscObject **)&Uf, NULL));
1343e10fd815SStefano Zampini     for (f = 0; f < nf; f++) {
13449566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(fec[f], &fecs[tnf + f]));
1345e10fd815SStefano Zampini       Ufds[tnf + f] = Uf[f];
1346e10fd815SStefano Zampini       sdim[tnf + f] = sd[f];
1347e10fd815SStefano Zampini     }
1348e10fd815SStefano Zampini     tnf += nf;
1349e10fd815SStefano Zampini   }
13509566063dSJacob Faibussowitsch   PetscCall(PetscViewerGLVisSetFields(viewer, tnf, (const char **)fecs, sdim, DMCompositeSampleGLVisFields_Private, (PetscObject *)Ufds, ctx, DestroyGLVisViewerCtx_Private));
135148a46eb9SPierre Jolivet   for (i = 0; i < tnf; i++) PetscCall(PetscFree(fecs[i]));
13529566063dSJacob Faibussowitsch   PetscCall(PetscFree3(fecs, sdim, Ufds));
1353e10fd815SStefano Zampini   PetscFunctionReturn(0);
1354e10fd815SStefano Zampini }
1355e10fd815SStefano Zampini 
13569371c9d4SSatish Balay PetscErrorCode DMRefine_Composite(DM dmi, MPI_Comm comm, DM *fine) {
135747c6ae99SBarry Smith   struct DMCompositeLink *next;
135847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dmi->data;
135947c6ae99SBarry Smith   DM                      dm;
136047c6ae99SBarry Smith 
136147c6ae99SBarry Smith   PetscFunctionBegin;
136247c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi, DM_CLASSID, 1);
136348a46eb9SPierre Jolivet   if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm));
13649566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmi));
136547c6ae99SBarry Smith   next = com->next;
13669566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm, fine));
136747c6ae99SBarry Smith 
136847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
136947c6ae99SBarry Smith   while (next) {
13709566063dSJacob Faibussowitsch     PetscCall(DMRefine(next->dm, comm, &dm));
13719566063dSJacob Faibussowitsch     PetscCall(DMCompositeAddDM(*fine, dm));
13729566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dm));
137347c6ae99SBarry Smith     next = next->next;
137447c6ae99SBarry Smith   }
137547c6ae99SBarry Smith   PetscFunctionReturn(0);
137647c6ae99SBarry Smith }
137747c6ae99SBarry Smith 
13789371c9d4SSatish Balay PetscErrorCode DMCoarsen_Composite(DM dmi, MPI_Comm comm, DM *fine) {
137914354c39SJed Brown   struct DMCompositeLink *next;
138014354c39SJed Brown   DM_Composite           *com = (DM_Composite *)dmi->data;
138114354c39SJed Brown   DM                      dm;
138214354c39SJed Brown 
138314354c39SJed Brown   PetscFunctionBegin;
138414354c39SJed Brown   PetscValidHeaderSpecific(dmi, DM_CLASSID, 1);
13859566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmi));
138648a46eb9SPierre Jolivet   if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm));
138714354c39SJed Brown   next = com->next;
13889566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(comm, fine));
138914354c39SJed Brown 
139014354c39SJed Brown   /* loop over packed objects, handling one at at time */
139114354c39SJed Brown   while (next) {
13929566063dSJacob Faibussowitsch     PetscCall(DMCoarsen(next->dm, comm, &dm));
13939566063dSJacob Faibussowitsch     PetscCall(DMCompositeAddDM(*fine, dm));
13949566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)dm));
139514354c39SJed Brown     next = next->next;
139614354c39SJed Brown   }
139714354c39SJed Brown   PetscFunctionReturn(0);
139814354c39SJed Brown }
139947c6ae99SBarry Smith 
14009371c9d4SSatish Balay PetscErrorCode DMCreateInterpolation_Composite(DM coarse, DM fine, Mat *A, Vec *v) {
14019ae5db72SJed Brown   PetscInt                m, n, M, N, nDM, i;
140247c6ae99SBarry Smith   struct DMCompositeLink *nextc;
140347c6ae99SBarry Smith   struct DMCompositeLink *nextf;
140425296bd5SBarry Smith   Vec                     gcoarse, gfine, *vecs;
140547c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite *)coarse->data;
140647c6ae99SBarry Smith   DM_Composite           *comfine   = (DM_Composite *)fine->data;
14079ae5db72SJed Brown   Mat                    *mats;
140847c6ae99SBarry Smith 
140947c6ae99SBarry Smith   PetscFunctionBegin;
141047c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse, DM_CLASSID, 1);
141147c6ae99SBarry Smith   PetscValidHeaderSpecific(fine, DM_CLASSID, 2);
14129566063dSJacob Faibussowitsch   PetscCall(DMSetUp(coarse));
14139566063dSJacob Faibussowitsch   PetscCall(DMSetUp(fine));
141447c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14159566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(coarse, &gcoarse));
14169566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(fine, &gfine));
14179566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gcoarse, &n));
14189566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(gfine, &m));
14199566063dSJacob Faibussowitsch   PetscCall(VecGetSize(gcoarse, &N));
14209566063dSJacob Faibussowitsch   PetscCall(VecGetSize(gfine, &M));
14219566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(coarse, &gcoarse));
14229566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(fine, &gfine));
142347c6ae99SBarry Smith 
14249ae5db72SJed Brown   nDM = comfine->nDM;
142563a3b9bcSJacob 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);
14269566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nDM * nDM, &mats));
142748a46eb9SPierre Jolivet   if (v) PetscCall(PetscCalloc1(nDM, &vecs));
142847c6ae99SBarry Smith 
142947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
14309ae5db72SJed Brown   for (nextc = comcoarse->next, nextf = comfine->next, i = 0; nextc; nextc = nextc->next, nextf = nextf->next, i++) {
14311baa6e33SBarry Smith     if (!v) PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], NULL));
14321baa6e33SBarry Smith     else PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], &vecs[i]));
143347c6ae99SBarry Smith   }
14349566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, nDM, NULL, mats, A));
14351baa6e33SBarry Smith   if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, vecs, v));
14369566063dSJacob Faibussowitsch   for (i = 0; i < nDM * nDM; i++) PetscCall(MatDestroy(&mats[i]));
14379566063dSJacob Faibussowitsch   PetscCall(PetscFree(mats));
143825296bd5SBarry Smith   if (v) {
14399566063dSJacob Faibussowitsch     for (i = 0; i < nDM; i++) PetscCall(VecDestroy(&vecs[i]));
14409566063dSJacob Faibussowitsch     PetscCall(PetscFree(vecs));
144125296bd5SBarry Smith   }
144247c6ae99SBarry Smith   PetscFunctionReturn(0);
144347c6ae99SBarry Smith }
144447c6ae99SBarry Smith 
14459371c9d4SSatish Balay static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm) {
14461411c6eeSJed Brown   DM_Composite           *com = (DM_Composite *)dm->data;
14471411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1448f7efa3c7SJed Brown   PetscInt                i;
14491411c6eeSJed Brown 
14501411c6eeSJed Brown   PetscFunctionBegin;
14511411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
14529566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetISLocalToGlobalMappings(dm, &ltogs));
14539566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm), com->nDM, ltogs, &dm->ltogmap));
14549566063dSJacob Faibussowitsch   for (i = 0; i < com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(&ltogs[i]));
14559566063dSJacob Faibussowitsch   PetscCall(PetscFree(ltogs));
14561411c6eeSJed Brown   PetscFunctionReturn(0);
14571411c6eeSJed Brown }
14581411c6eeSJed Brown 
14599371c9d4SSatish Balay PetscErrorCode DMCreateColoring_Composite(DM dm, ISColoringType ctype, ISColoring *coloring) {
146047c6ae99SBarry Smith   PetscInt         n, i, cnt;
146147c6ae99SBarry Smith   ISColoringValue *colors;
146247c6ae99SBarry Smith   PetscBool        dense  = PETSC_FALSE;
146347c6ae99SBarry Smith   ISColoringValue  maxcol = 0;
146447c6ae99SBarry Smith   DM_Composite    *com    = (DM_Composite *)dm->data;
146547c6ae99SBarry Smith 
146647c6ae99SBarry Smith   PetscFunctionBegin;
146747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
146808401ef6SPierre Jolivet   PetscCheck(ctype != IS_COLORING_LOCAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only global coloring supported");
1469f7d195e4SLawrence Mitchell   if (ctype == IS_COLORING_GLOBAL) {
147047c6ae99SBarry Smith     n = com->n;
1471ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unknown ISColoringType");
14729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &colors)); /* freed in ISColoringDestroy() */
147347c6ae99SBarry Smith 
14749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL));
147547c6ae99SBarry Smith   if (dense) {
1476ad540459SPierre Jolivet     for (i = 0; i < n; i++) colors[i] = (ISColoringValue)(com->rstart + i);
147747c6ae99SBarry Smith     maxcol = com->N;
147847c6ae99SBarry Smith   } else {
147947c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
148047c6ae99SBarry Smith     PetscMPIInt             rank;
148147c6ae99SBarry Smith 
14829566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
148347c6ae99SBarry Smith     cnt = 0;
148447c6ae99SBarry Smith     while (next) {
148547c6ae99SBarry Smith       ISColoring lcoloring;
148647c6ae99SBarry Smith 
14879566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(next->dm, IS_COLORING_GLOBAL, &lcoloring));
1488ad540459SPierre Jolivet       for (i = 0; i < lcoloring->N; i++) colors[cnt++] = maxcol + lcoloring->colors[i];
148947c6ae99SBarry Smith       maxcol += lcoloring->n;
14909566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&lcoloring));
149147c6ae99SBarry Smith       next = next->next;
149247c6ae99SBarry Smith     }
149347c6ae99SBarry Smith   }
14949566063dSJacob Faibussowitsch   PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), maxcol, n, colors, PETSC_OWN_POINTER, coloring));
149547c6ae99SBarry Smith   PetscFunctionReturn(0);
149647c6ae99SBarry Smith }
149747c6ae99SBarry Smith 
14989371c9d4SSatish Balay PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) {
149947c6ae99SBarry Smith   struct DMCompositeLink *next;
150047c6ae99SBarry Smith   PetscScalar            *garray, *larray;
150147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
150247c6ae99SBarry Smith 
150347c6ae99SBarry Smith   PetscFunctionBegin;
150447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
150547c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
150639d35262SVincent Le Chenadec 
150748a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
150839d35262SVincent Le Chenadec 
15099566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gvec, &garray));
15109566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lvec, &larray));
151147c6ae99SBarry Smith 
151247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
151339d35262SVincent Le Chenadec   next = com->next;
151447c6ae99SBarry Smith   while (next) {
151547c6ae99SBarry Smith     Vec      local, global;
151647c6ae99SBarry Smith     PetscInt N;
151747c6ae99SBarry Smith 
15189566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
15199566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(global, &N));
15209566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(global, garray));
15219566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local));
15229566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local, larray));
15239566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(next->dm, global, mode, local));
15249566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(next->dm, global, mode, local));
15259566063dSJacob Faibussowitsch     PetscCall(VecResetArray(global));
15269566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local));
15279566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
15289566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local));
152939d35262SVincent Le Chenadec 
153006ebdd98SJed Brown     larray += next->nlocal;
153139d35262SVincent Le Chenadec     garray += next->n;
153247c6ae99SBarry Smith     next = next->next;
153347c6ae99SBarry Smith   }
153447c6ae99SBarry Smith 
15359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gvec, NULL));
15369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lvec, NULL));
153747c6ae99SBarry Smith   PetscFunctionReturn(0);
153847c6ae99SBarry Smith }
153947c6ae99SBarry Smith 
15409371c9d4SSatish Balay PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) {
15410c010503SBarry Smith   PetscFunctionBegin;
154239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
154339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
154439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 4);
154539d35262SVincent Le Chenadec   PetscFunctionReturn(0);
154639d35262SVincent Le Chenadec }
154739d35262SVincent Le Chenadec 
15489371c9d4SSatish Balay PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) {
154939d35262SVincent Le Chenadec   struct DMCompositeLink *next;
155039d35262SVincent Le Chenadec   PetscScalar            *larray, *garray;
155139d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite *)dm->data;
155239d35262SVincent Le Chenadec 
155339d35262SVincent Le Chenadec   PetscFunctionBegin;
155439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
155539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
155639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
155739d35262SVincent Le Chenadec 
155848a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
155939d35262SVincent Le Chenadec 
15609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lvec, &larray));
15619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gvec, &garray));
156239d35262SVincent Le Chenadec 
156339d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
156439d35262SVincent Le Chenadec   next = com->next;
156539d35262SVincent Le Chenadec   while (next) {
156639d35262SVincent Le Chenadec     Vec global, local;
156739d35262SVincent Le Chenadec 
15689566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local));
15699566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local, larray));
15709566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(next->dm, &global));
15719566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(global, garray));
15729566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(next->dm, local, mode, global));
15739566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(next->dm, local, mode, global));
15749566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local));
15759566063dSJacob Faibussowitsch     PetscCall(VecResetArray(global));
15769566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(next->dm, &global));
15779566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local));
157839d35262SVincent Le Chenadec 
157939d35262SVincent Le Chenadec     garray += next->n;
158039d35262SVincent Le Chenadec     larray += next->nlocal;
158139d35262SVincent Le Chenadec     next = next->next;
158239d35262SVincent Le Chenadec   }
158339d35262SVincent Le Chenadec 
15849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gvec, NULL));
15859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lvec, NULL));
158639d35262SVincent Le Chenadec 
158739d35262SVincent Le Chenadec   PetscFunctionReturn(0);
158839d35262SVincent Le Chenadec }
158939d35262SVincent Le Chenadec 
15909371c9d4SSatish Balay PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) {
159139d35262SVincent Le Chenadec   PetscFunctionBegin;
159239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
159339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
159439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
159539d35262SVincent Le Chenadec   PetscFunctionReturn(0);
159639d35262SVincent Le Chenadec }
159739d35262SVincent Le Chenadec 
15989371c9d4SSatish Balay PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2) {
159939d35262SVincent Le Chenadec   struct DMCompositeLink *next;
160039d35262SVincent Le Chenadec   PetscScalar            *array1, *array2;
160139d35262SVincent Le Chenadec   DM_Composite           *com = (DM_Composite *)dm->data;
160239d35262SVincent Le Chenadec 
160339d35262SVincent Le Chenadec   PetscFunctionBegin;
160439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
160539d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2);
160639d35262SVincent Le Chenadec   PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4);
160739d35262SVincent Le Chenadec 
160848a46eb9SPierre Jolivet   if (!com->setup) PetscCall(DMSetUp(dm));
160939d35262SVincent Le Chenadec 
16109566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec1, &array1));
16119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec2, &array2));
161239d35262SVincent Le Chenadec 
161339d35262SVincent Le Chenadec   /* loop over packed objects, handling one at at time */
161439d35262SVincent Le Chenadec   next = com->next;
161539d35262SVincent Le Chenadec   while (next) {
161639d35262SVincent Le Chenadec     Vec local1, local2;
161739d35262SVincent Le Chenadec 
16189566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local1));
16199566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local1, array1));
16209566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(next->dm, &local2));
16219566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(local2, array2));
16229566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2));
16239566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2));
16249566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local2));
16259566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local2));
16269566063dSJacob Faibussowitsch     PetscCall(VecResetArray(local1));
16279566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(next->dm, &local1));
162839d35262SVincent Le Chenadec 
162939d35262SVincent Le Chenadec     array1 += next->nlocal;
163039d35262SVincent Le Chenadec     array2 += next->nlocal;
163139d35262SVincent Le Chenadec     next = next->next;
163239d35262SVincent Le Chenadec   }
163339d35262SVincent Le Chenadec 
16349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec1, NULL));
16359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec2, NULL));
163639d35262SVincent Le Chenadec 
163739d35262SVincent Le Chenadec   PetscFunctionReturn(0);
163839d35262SVincent Le Chenadec }
163939d35262SVincent Le Chenadec 
16409371c9d4SSatish Balay PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) {
164139d35262SVincent Le Chenadec   PetscFunctionBegin;
164239d35262SVincent Le Chenadec   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
164339d35262SVincent Le Chenadec   PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
164439d35262SVincent Le Chenadec   PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16450c010503SBarry Smith   PetscFunctionReturn(0);
16460c010503SBarry Smith }
164747c6ae99SBarry Smith 
16486ae3a549SBarry Smith /*MC
16496ae3a549SBarry Smith    DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs
16506ae3a549SBarry Smith 
16516ae3a549SBarry Smith   Level: intermediate
16526ae3a549SBarry Smith 
1653db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()`
16546ae3a549SBarry Smith M*/
16556ae3a549SBarry Smith 
16569371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p) {
1657a4121054SBarry Smith   DM_Composite *com;
1658a4121054SBarry Smith 
1659a4121054SBarry Smith   PetscFunctionBegin;
1660*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&com));
1661a4121054SBarry Smith   p->data     = com;
1662a4121054SBarry Smith   com->n      = 0;
16637ac2b803SAlex Fikl   com->nghost = 0;
16640298fd71SBarry Smith   com->next   = NULL;
1665a4121054SBarry Smith   com->nDM    = 0;
1666a4121054SBarry Smith 
1667a4121054SBarry Smith   p->ops->createglobalvector       = DMCreateGlobalVector_Composite;
1668a4121054SBarry Smith   p->ops->createlocalvector        = DMCreateLocalVector_Composite;
1669184d77edSJed Brown   p->ops->getlocaltoglobalmapping  = DMGetLocalToGlobalMapping_Composite;
16704d343eeaSMatthew G Knepley   p->ops->createfieldis            = DMCreateFieldIS_Composite;
167116621825SDmitry Karpeev   p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite;
1672a4121054SBarry Smith   p->ops->refine                   = DMRefine_Composite;
167314354c39SJed Brown   p->ops->coarsen                  = DMCoarsen_Composite;
167425296bd5SBarry Smith   p->ops->createinterpolation      = DMCreateInterpolation_Composite;
167525296bd5SBarry Smith   p->ops->creatematrix             = DMCreateMatrix_Composite;
1676e727c939SJed Brown   p->ops->getcoloring              = DMCreateColoring_Composite;
1677a4121054SBarry Smith   p->ops->globaltolocalbegin       = DMGlobalToLocalBegin_Composite;
1678a4121054SBarry Smith   p->ops->globaltolocalend         = DMGlobalToLocalEnd_Composite;
167939d35262SVincent Le Chenadec   p->ops->localtoglobalbegin       = DMLocalToGlobalBegin_Composite;
168039d35262SVincent Le Chenadec   p->ops->localtoglobalend         = DMLocalToGlobalEnd_Composite;
168139d35262SVincent Le Chenadec   p->ops->localtolocalbegin        = DMLocalToLocalBegin_Composite;
168239d35262SVincent Le Chenadec   p->ops->localtolocalend          = DMLocalToLocalEnd_Composite;
1683a4121054SBarry Smith   p->ops->destroy                  = DMDestroy_Composite;
1684a4121054SBarry Smith   p->ops->view                     = DMView_Composite;
1685a4121054SBarry Smith   p->ops->setup                    = DMSetUp_Composite;
1686e10fd815SStefano Zampini 
16879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite));
1688a4121054SBarry Smith   PetscFunctionReturn(0);
1689a4121054SBarry Smith }
1690a4121054SBarry Smith 
16911fd49c25SBarry Smith /*@
16920c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
16930c010503SBarry Smith       vectors made up of several subvectors.
16940c010503SBarry Smith 
1695d083f849SBarry Smith     Collective
169647c6ae99SBarry Smith 
169747c6ae99SBarry Smith     Input Parameter:
16980c010503SBarry Smith .   comm - the processors that will share the global vector
16990c010503SBarry Smith 
17000c010503SBarry Smith     Output Parameters:
17010c010503SBarry Smith .   packer - the packer object
170247c6ae99SBarry Smith 
170347c6ae99SBarry Smith     Level: advanced
170447c6ae99SBarry Smith 
1705c2e3fba1SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCOMPOSITE`, `DMCreate()`
1706db781477SPatrick Sanan          `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`
1707db781477SPatrick Sanan          `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
170847c6ae99SBarry Smith 
170947c6ae99SBarry Smith @*/
17109371c9d4SSatish Balay PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer) {
171147c6ae99SBarry Smith   PetscFunctionBegin;
17120c010503SBarry Smith   PetscValidPointer(packer, 2);
17139566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, packer));
17149566063dSJacob Faibussowitsch   PetscCall(DMSetType(*packer, DMCOMPOSITE));
171547c6ae99SBarry Smith   PetscFunctionReturn(0);
171647c6ae99SBarry Smith }
1717