1ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h> /*I "petscdmcomposite.h" I*/
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>
3e10fd815SStefano Zampini #include <petsc/private/glvisviewerimpl.h>
42764a2aaSMatthew G. Knepley #include <petscds.h>
547c6ae99SBarry Smith
647c6ae99SBarry Smith /*@C
747c6ae99SBarry Smith DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
8dce8aebaSBarry Smith separate components `DM` in a `DMCOMPOSITE` to build the correct matrix nonzero structure.
947c6ae99SBarry Smith
1020f4b53cSBarry Smith Logically Collective; No Fortran Support
1147c6ae99SBarry Smith
12d8d19677SJose E. Roman Input Parameters:
1347c6ae99SBarry Smith + dm - the composite object
1460225df5SJacob Faibussowitsch - FormCoupleLocations - routine to set the nonzero locations in the matrix
1547c6ae99SBarry Smith
1647c6ae99SBarry Smith Level: advanced
1747c6ae99SBarry Smith
18761ec508SZach Atkins Notes:
19dce8aebaSBarry Smith See `DMSetApplicationContext()` and `DMGetApplicationContext()` for how to get user information into
2047c6ae99SBarry Smith this routine
2147c6ae99SBarry Smith
22761ec508SZach Atkins The provided function should have a signature matching
23761ec508SZach Atkins .vb
24761ec508SZach Atkins PetscErrorCode your_form_couple_method(DM dm, Mat J, PetscInt *dnz, PetscInt *onz, PetscInt rstart, PetscInt nrows, PetscInt start, PetscInt end);
25761ec508SZach Atkins .ve
26761ec508SZach Atkins where
27761ec508SZach Atkins + dm - the composite object
28761ec508SZach Atkins . J - the constructed matrix, or `NULL`. If provided, the function should fill the existing nonzero pattern with zeros (only `dm` and `rstart` are valid in this case).
29761ec508SZach Atkins . dnz - array counting the number of on-diagonal non-zero entries per row, where on-diagonal means that this process owns both the row and column
30761ec508SZach Atkins . onz - array counting the number of off-diagonal non-zero entries per row, where off-diagonal means that this process owns the row
31761ec508SZach Atkins . rstart - offset into `*nz` arrays, for local row index `r`, update `onz[r - rstart]` or `dnz[r - rstart]`
32761ec508SZach Atkins . nrows - number of owned global rows
33761ec508SZach Atkins . start - the first owned global index
34761ec508SZach Atkins - end - the last owned global index + 1
35761ec508SZach Atkins
36761ec508SZach Atkins If `J` is not `NULL`, then the only other valid parameter is `rstart`
37761ec508SZach Atkins
38761ec508SZach Atkins The user coupling function has a weird and poorly documented interface and is not tested, it should be removed
39761ec508SZach Atkins
40dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`
4147c6ae99SBarry Smith @*/
DMCompositeSetCoupling(DM dm,PetscErrorCode (* FormCoupleLocations)(DM,Mat,PetscInt *,PetscInt *,PetscInt,PetscInt,PetscInt,PetscInt))42d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeSetCoupling(DM dm, PetscErrorCode (*FormCoupleLocations)(DM, Mat, PetscInt *, PetscInt *, PetscInt, PetscInt, PetscInt, PetscInt))
43d71ae5a4SJacob Faibussowitsch {
4447c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
4571b14b3eSStefano Zampini PetscBool flg;
4647c6ae99SBarry Smith
4747c6ae99SBarry Smith PetscFunctionBegin;
489566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
497a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
5047c6ae99SBarry Smith com->FormCoupleLocations = FormCoupleLocations;
513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5247c6ae99SBarry Smith }
5347c6ae99SBarry Smith
DMDestroy_Composite(DM dm)5466976f2fSJacob Faibussowitsch static PetscErrorCode DMDestroy_Composite(DM dm)
55d71ae5a4SJacob Faibussowitsch {
5647c6ae99SBarry Smith struct DMCompositeLink *next, *prev;
5747c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
5847c6ae99SBarry Smith
5947c6ae99SBarry Smith PetscFunctionBegin;
6047c6ae99SBarry Smith next = com->next;
6147c6ae99SBarry Smith while (next) {
6247c6ae99SBarry Smith prev = next;
6347c6ae99SBarry Smith next = next->next;
649566063dSJacob Faibussowitsch PetscCall(DMDestroy(&prev->dm));
659566063dSJacob Faibussowitsch PetscCall(PetscFree(prev->grstarts));
669566063dSJacob Faibussowitsch PetscCall(PetscFree(prev));
6747c6ae99SBarry Smith }
689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL));
69435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
709566063dSJacob Faibussowitsch PetscCall(PetscFree(com));
713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7247c6ae99SBarry Smith }
7347c6ae99SBarry Smith
DMView_Composite(DM dm,PetscViewer v)7466976f2fSJacob Faibussowitsch static PetscErrorCode DMView_Composite(DM dm, PetscViewer v)
75d71ae5a4SJacob Faibussowitsch {
769f196a02SMartin Diehl PetscBool isascii;
7747c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
7847c6ae99SBarry Smith
7947c6ae99SBarry Smith PetscFunctionBegin;
809f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii));
819f196a02SMartin Diehl if (isascii) {
8247c6ae99SBarry Smith struct DMCompositeLink *lnk = com->next;
8347c6ae99SBarry Smith PetscInt i;
8447c6ae99SBarry Smith
859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "DM (%s)\n", ((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix"));
8663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, " contains %" PetscInt_FMT " DMs\n", com->nDM));
879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v));
8847c6ae99SBarry Smith for (i = 0; lnk; lnk = lnk->next, i++) {
8963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "Link %" PetscInt_FMT ": DM of type %s\n", i, ((PetscObject)lnk->dm)->type_name));
909566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v));
919566063dSJacob Faibussowitsch PetscCall(DMView(lnk->dm, v));
929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v));
9347c6ae99SBarry Smith }
949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v));
9547c6ae99SBarry Smith }
963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9747c6ae99SBarry Smith }
9847c6ae99SBarry Smith
DMSetUp_Composite(DM dm)9966976f2fSJacob Faibussowitsch static PetscErrorCode DMSetUp_Composite(DM dm)
100d71ae5a4SJacob Faibussowitsch {
10147c6ae99SBarry Smith PetscInt nprev = 0;
10247c6ae99SBarry Smith PetscMPIInt rank, size;
10347c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
10447c6ae99SBarry Smith struct DMCompositeLink *next = com->next;
10547c6ae99SBarry Smith PetscLayout map;
10647c6ae99SBarry Smith
10747c6ae99SBarry Smith PetscFunctionBegin;
10828b400f6SJacob Faibussowitsch PetscCheck(!com->setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Packer has already been setup");
1099566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &map));
1109566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(map, com->n));
1119566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetSize(map, PETSC_DETERMINE));
1129566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(map, 1));
1139566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(map));
1149566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &com->N));
1159566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(map, &com->rstart, NULL));
1169566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&map));
11747c6ae99SBarry Smith
1189ae5db72SJed Brown /* now set the rstart for each linked vector */
1199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
1209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
12147c6ae99SBarry Smith while (next) {
12247c6ae99SBarry Smith next->rstart = nprev;
12306ebdd98SJed Brown nprev += next->n;
12447c6ae99SBarry Smith next->grstart = com->rstart + next->rstart;
1259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &next->grstarts));
1269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&next->grstart, 1, MPIU_INT, next->grstarts, 1, MPIU_INT, PetscObjectComm((PetscObject)dm)));
12747c6ae99SBarry Smith next = next->next;
12847c6ae99SBarry Smith }
12947c6ae99SBarry Smith com->setup = PETSC_TRUE;
1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13147c6ae99SBarry Smith }
13247c6ae99SBarry Smith
13373e31fe2SJed Brown /*@
134aaa8cc7dSPierre Jolivet DMCompositeGetNumberDM - Gets the number of `DM` objects in the `DMCOMPOSITE`
13547c6ae99SBarry Smith representation.
13647c6ae99SBarry Smith
13747c6ae99SBarry Smith Not Collective
13847c6ae99SBarry Smith
13947c6ae99SBarry Smith Input Parameter:
140dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
14147c6ae99SBarry Smith
14247c6ae99SBarry Smith Output Parameter:
143dce8aebaSBarry Smith . nDM - the number of `DM`
14447c6ae99SBarry Smith
14547c6ae99SBarry Smith Level: beginner
14647c6ae99SBarry Smith
147dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`
14847c6ae99SBarry Smith @*/
DMCompositeGetNumberDM(DM dm,PetscInt * nDM)149d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetNumberDM(DM dm, PetscInt *nDM)
150d71ae5a4SJacob Faibussowitsch {
15147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
15271b14b3eSStefano Zampini PetscBool flg;
1535fd66863SKarl Rupp
15447c6ae99SBarry Smith PetscFunctionBegin;
15547c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1569566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1577a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
15847c6ae99SBarry Smith *nDM = com->nDM;
1593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
16047c6ae99SBarry Smith }
16147c6ae99SBarry Smith
16247c6ae99SBarry Smith /*@C
16347c6ae99SBarry Smith DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
16447c6ae99SBarry Smith representation.
16547c6ae99SBarry Smith
16620f4b53cSBarry Smith Collective
16747c6ae99SBarry Smith
1689ae5db72SJed Brown Input Parameters:
169dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object
1709ae5db72SJed Brown - gvec - the global vector
1719ae5db72SJed Brown
1722fe279fdSBarry Smith Output Parameter:
173a4e35b19SJacob Faibussowitsch . ... - the packed parallel vectors, `NULL` for those that are not needed
17447c6ae99SBarry Smith
17547c6ae99SBarry Smith Level: advanced
17647c6ae99SBarry Smith
177dce8aebaSBarry Smith Note:
178dce8aebaSBarry Smith Use `DMCompositeRestoreAccess()` to return the vectors when you no longer need them
179dce8aebaSBarry Smith
180dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetEntries()`, `DMCompositeScatter()`
18147c6ae99SBarry Smith @*/
DMCompositeGetAccess(DM dm,Vec gvec,...)182d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccess(DM dm, Vec gvec, ...)
183d71ae5a4SJacob Faibussowitsch {
18447c6ae99SBarry Smith va_list Argp;
18547c6ae99SBarry Smith struct DMCompositeLink *next;
18647c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
1875edff71fSBarry Smith PetscInt readonly;
18871b14b3eSStefano Zampini PetscBool flg;
18947c6ae99SBarry Smith
19047c6ae99SBarry Smith PetscFunctionBegin;
19147c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
19247c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
1939566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
1947a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
19547c6ae99SBarry Smith next = com->next;
19648a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
19747c6ae99SBarry Smith
1989566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec, &readonly));
19915229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
20047c6ae99SBarry Smith va_start(Argp, gvec);
20147c6ae99SBarry Smith while (next) {
20247c6ae99SBarry Smith Vec *vec;
20347c6ae99SBarry Smith vec = va_arg(Argp, Vec *);
2049ae5db72SJed Brown if (vec) {
2059566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, vec));
2065edff71fSBarry Smith if (readonly) {
2075edff71fSBarry Smith const PetscScalar *array;
2089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array));
2099566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(*vec, array + next->rstart));
2109566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(*vec));
2119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array));
2125edff71fSBarry Smith } else {
2135edff71fSBarry Smith PetscScalar *array;
2149566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array));
2159566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(*vec, array + next->rstart));
2169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array));
21747c6ae99SBarry Smith }
2185edff71fSBarry Smith }
21947c6ae99SBarry Smith next = next->next;
22047c6ae99SBarry Smith }
22147c6ae99SBarry Smith va_end(Argp);
2223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
22347c6ae99SBarry Smith }
22447c6ae99SBarry Smith
2255d83a8b1SBarry Smith /*@
226f73e5cebSJed Brown DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global
227f73e5cebSJed Brown representation.
228f73e5cebSJed Brown
22920f4b53cSBarry Smith Collective
230f73e5cebSJed Brown
231f73e5cebSJed Brown Input Parameters:
232dce8aebaSBarry Smith + dm - the `DMCOMPOSITE`
233f73e5cebSJed Brown . pvec - packed vector
234f73e5cebSJed Brown . nwanted - number of vectors wanted
235ce78bad3SBarry Smith - wanted - sorted array of integers indicating thde vectors wanted, or `NULL` to get all vectors, length `nwanted`
236f73e5cebSJed Brown
2372fe279fdSBarry Smith Output Parameter:
238f13dfd9eSBarry Smith . vecs - array of requested global vectors (must be previously allocated and of length `nwanted`)
239f73e5cebSJed Brown
240f73e5cebSJed Brown Level: advanced
241f73e5cebSJed Brown
242dce8aebaSBarry Smith Note:
243dce8aebaSBarry Smith Use `DMCompositeRestoreAccessArray()` to return the vectors when you no longer need them
244dce8aebaSBarry Smith
245dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
246f73e5cebSJed Brown @*/
DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt wanted[],Vec vecs[])2475d83a8b1SBarry Smith PetscErrorCode DMCompositeGetAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[])
248d71ae5a4SJacob Faibussowitsch {
249f73e5cebSJed Brown struct DMCompositeLink *link;
250f73e5cebSJed Brown PetscInt i, wnum;
251f73e5cebSJed Brown DM_Composite *com = (DM_Composite *)dm->data;
252bee642f7SBarry Smith PetscInt readonly;
25371b14b3eSStefano Zampini PetscBool flg;
254f73e5cebSJed Brown
255f73e5cebSJed Brown PetscFunctionBegin;
256f73e5cebSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
257f73e5cebSJed Brown PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
2589566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
2597a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
26048a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
261f73e5cebSJed Brown
2629566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly));
263f73e5cebSJed Brown for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
264f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) {
265f73e5cebSJed Brown Vec v;
2669566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(link->dm, &v));
267bee642f7SBarry Smith if (readonly) {
268bee642f7SBarry Smith const PetscScalar *array;
2699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvec, &array));
2709566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + link->rstart));
2719566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(v));
2729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvec, &array));
273bee642f7SBarry Smith } else {
274bee642f7SBarry Smith PetscScalar *array;
2759566063dSJacob Faibussowitsch PetscCall(VecGetArray(pvec, &array));
2769566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + link->rstart));
2779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pvec, &array));
278bee642f7SBarry Smith }
279f73e5cebSJed Brown vecs[wnum++] = v;
280f73e5cebSJed Brown }
281f73e5cebSJed Brown }
2823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
283f73e5cebSJed Brown }
284f73e5cebSJed Brown
2855d83a8b1SBarry Smith /*@
2867ac2b803SAlex Fikl DMCompositeGetLocalAccessArray - Allows one to access the individual
2877ac2b803SAlex Fikl packed vectors in their local representation.
2887ac2b803SAlex Fikl
28920f4b53cSBarry Smith Collective
2907ac2b803SAlex Fikl
2917ac2b803SAlex Fikl Input Parameters:
292dce8aebaSBarry Smith + dm - the `DMCOMPOSITE`
2937ac2b803SAlex Fikl . pvec - packed vector
2947ac2b803SAlex Fikl . nwanted - number of vectors wanted
295f13dfd9eSBarry Smith - wanted - sorted array of vectors wanted, or `NULL` to get all vectors, length `nwanted`
2967ac2b803SAlex Fikl
2972fe279fdSBarry Smith Output Parameter:
298f13dfd9eSBarry Smith . vecs - array of requested local vectors (must be allocated and of length `nwanted`)
2997ac2b803SAlex Fikl
3007ac2b803SAlex Fikl Level: advanced
3017ac2b803SAlex Fikl
302dce8aebaSBarry Smith Note:
303dce8aebaSBarry Smith Use `DMCompositeRestoreLocalAccessArray()` to return the vectors
304dce8aebaSBarry Smith when you no longer need them.
305dce8aebaSBarry Smith
306dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`,
307db781477SPatrick Sanan `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
3087ac2b803SAlex Fikl @*/
DMCompositeGetLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt wanted[],Vec vecs[])3095d83a8b1SBarry Smith PetscErrorCode DMCompositeGetLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[])
310d71ae5a4SJacob Faibussowitsch {
3117ac2b803SAlex Fikl struct DMCompositeLink *link;
3127ac2b803SAlex Fikl PetscInt i, wnum;
3137ac2b803SAlex Fikl DM_Composite *com = (DM_Composite *)dm->data;
3147ac2b803SAlex Fikl PetscInt readonly;
3157ac2b803SAlex Fikl PetscInt nlocal = 0;
31671b14b3eSStefano Zampini PetscBool flg;
3177ac2b803SAlex Fikl
3187ac2b803SAlex Fikl PetscFunctionBegin;
3197ac2b803SAlex Fikl PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3207ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
3219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3227a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
32348a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
3247ac2b803SAlex Fikl
3259566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly));
3267ac2b803SAlex Fikl for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
3277ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) {
3287ac2b803SAlex Fikl Vec v;
3299566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(link->dm, &v));
3307ac2b803SAlex Fikl if (readonly) {
3317ac2b803SAlex Fikl const PetscScalar *array;
3329566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvec, &array));
3339566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + nlocal));
334b1c3483dSMark 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
3359566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(v));
3369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvec, &array));
3377ac2b803SAlex Fikl } else {
3387ac2b803SAlex Fikl PetscScalar *array;
3399566063dSJacob Faibussowitsch PetscCall(VecGetArray(pvec, &array));
3409566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + nlocal));
3419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pvec, &array));
3427ac2b803SAlex Fikl }
3437ac2b803SAlex Fikl vecs[wnum++] = v;
3447ac2b803SAlex Fikl }
3457ac2b803SAlex Fikl
3467ac2b803SAlex Fikl nlocal += link->nlocal;
3477ac2b803SAlex Fikl }
3483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3497ac2b803SAlex Fikl }
3507ac2b803SAlex Fikl
35147c6ae99SBarry Smith /*@C
352dce8aebaSBarry Smith DMCompositeRestoreAccess - Returns the vectors obtained with `DMCompositeGetAccess()`
35347c6ae99SBarry Smith representation.
35447c6ae99SBarry Smith
35520f4b53cSBarry Smith Collective
35647c6ae99SBarry Smith
3579ae5db72SJed Brown Input Parameters:
358dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object
35947c6ae99SBarry Smith . gvec - the global vector
360a4e35b19SJacob Faibussowitsch - ... - the individual parallel vectors, `NULL` for those that are not needed
36147c6ae99SBarry Smith
36247c6ae99SBarry Smith Level: advanced
36347c6ae99SBarry Smith
364dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
365db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`,
36642747ad1SJacob Faibussowitsch `DMCompositeGetAccess()`
36747c6ae99SBarry Smith @*/
DMCompositeRestoreAccess(DM dm,Vec gvec,...)368d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...)
369d71ae5a4SJacob Faibussowitsch {
37047c6ae99SBarry Smith va_list Argp;
37147c6ae99SBarry Smith struct DMCompositeLink *next;
37247c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
3735edff71fSBarry Smith PetscInt readonly;
37471b14b3eSStefano Zampini PetscBool flg;
37547c6ae99SBarry Smith
37647c6ae99SBarry Smith PetscFunctionBegin;
37747c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
37847c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
3799566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
3807a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
38147c6ae99SBarry Smith next = com->next;
38248a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
38347c6ae99SBarry Smith
3849566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec, &readonly));
38515229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
38647c6ae99SBarry Smith va_start(Argp, gvec);
38747c6ae99SBarry Smith while (next) {
38847c6ae99SBarry Smith Vec *vec;
38947c6ae99SBarry Smith vec = va_arg(Argp, Vec *);
3909ae5db72SJed Brown if (vec) {
3919566063dSJacob Faibussowitsch PetscCall(VecResetArray(*vec));
3921baa6e33SBarry Smith if (readonly) PetscCall(VecLockReadPop(*vec));
3939566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, vec));
39447c6ae99SBarry Smith }
39547c6ae99SBarry Smith next = next->next;
39647c6ae99SBarry Smith }
39747c6ae99SBarry Smith va_end(Argp);
3983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
39947c6ae99SBarry Smith }
40047c6ae99SBarry Smith
4015d83a8b1SBarry Smith /*@
402dce8aebaSBarry Smith DMCompositeRestoreAccessArray - Returns the vectors obtained with `DMCompositeGetAccessArray()`
403f73e5cebSJed Brown
40420f4b53cSBarry Smith Collective
405f73e5cebSJed Brown
406f73e5cebSJed Brown Input Parameters:
407dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object
408f73e5cebSJed Brown . pvec - packed vector
409f73e5cebSJed Brown . nwanted - number of vectors wanted
410f13dfd9eSBarry Smith . wanted - sorted array of vectors wanted, or `NULL` to restore all vectors
411f13dfd9eSBarry Smith - vecs - array of global vectors
412f73e5cebSJed Brown
413f73e5cebSJed Brown Level: advanced
414f73e5cebSJed Brown
415dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()`
416f73e5cebSJed Brown @*/
DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt wanted[],Vec vecs[])4175d83a8b1SBarry Smith PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[])
418d71ae5a4SJacob Faibussowitsch {
419f73e5cebSJed Brown struct DMCompositeLink *link;
420f73e5cebSJed Brown PetscInt i, wnum;
421f73e5cebSJed Brown DM_Composite *com = (DM_Composite *)dm->data;
422bee642f7SBarry Smith PetscInt readonly;
42371b14b3eSStefano Zampini PetscBool flg;
424f73e5cebSJed Brown
425f73e5cebSJed Brown PetscFunctionBegin;
426f73e5cebSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
427f73e5cebSJed Brown PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4297a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
43048a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
431f73e5cebSJed Brown
4329566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly));
433f73e5cebSJed Brown for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
434f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) {
4359566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum]));
43648a46eb9SPierre Jolivet if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4379566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum]));
438f73e5cebSJed Brown wnum++;
439f73e5cebSJed Brown }
440f73e5cebSJed Brown }
4413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
442f73e5cebSJed Brown }
443f73e5cebSJed Brown
4445d83a8b1SBarry Smith /*@
445dce8aebaSBarry Smith DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with `DMCompositeGetLocalAccessArray()`.
4467ac2b803SAlex Fikl
44720f4b53cSBarry Smith Collective
4487ac2b803SAlex Fikl
4497ac2b803SAlex Fikl Input Parameters:
450dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object
4517ac2b803SAlex Fikl . pvec - packed vector
4527ac2b803SAlex Fikl . nwanted - number of vectors wanted
453f13dfd9eSBarry Smith . wanted - sorted array of vectors wanted, or `NULL` to restore all vectors
454f13dfd9eSBarry Smith - vecs - array of local vectors
4557ac2b803SAlex Fikl
4567ac2b803SAlex Fikl Level: advanced
4577ac2b803SAlex Fikl
458dce8aebaSBarry Smith Note:
459f13dfd9eSBarry Smith `nwanted` and `wanted` must match the values given to `DMCompositeGetLocalAccessArray()`
4607ac2b803SAlex Fikl otherwise the call will fail.
4617ac2b803SAlex Fikl
462dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`,
463db781477SPatrick Sanan `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`,
464db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeGather()`
4657ac2b803SAlex Fikl @*/
DMCompositeRestoreLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt wanted[],Vec * vecs)4665d83a8b1SBarry Smith PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec *vecs)
467d71ae5a4SJacob Faibussowitsch {
4687ac2b803SAlex Fikl struct DMCompositeLink *link;
4697ac2b803SAlex Fikl PetscInt i, wnum;
4707ac2b803SAlex Fikl DM_Composite *com = (DM_Composite *)dm->data;
4717ac2b803SAlex Fikl PetscInt readonly;
47271b14b3eSStefano Zampini PetscBool flg;
4737ac2b803SAlex Fikl
4747ac2b803SAlex Fikl PetscFunctionBegin;
4757ac2b803SAlex Fikl PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4767ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2);
4779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
4787a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
47948a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
4807ac2b803SAlex Fikl
4819566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly));
4827ac2b803SAlex Fikl for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) {
4837ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) {
4849566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum]));
48548a46eb9SPierre Jolivet if (readonly) PetscCall(VecLockReadPop(vecs[wnum]));
4869566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum]));
4877ac2b803SAlex Fikl wnum++;
4887ac2b803SAlex Fikl }
4897ac2b803SAlex Fikl }
4903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4917ac2b803SAlex Fikl }
4927ac2b803SAlex Fikl
49347c6ae99SBarry Smith /*@C
49447c6ae99SBarry Smith DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
49547c6ae99SBarry Smith
49620f4b53cSBarry Smith Collective
49747c6ae99SBarry Smith
4989ae5db72SJed Brown Input Parameters:
499dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object
50047c6ae99SBarry Smith . gvec - the global vector
501a4e35b19SJacob Faibussowitsch - ... - the individual sequential vectors, `NULL` for those that are not needed
50247c6ae99SBarry Smith
50347c6ae99SBarry Smith Level: advanced
50447c6ae99SBarry Smith
505dce8aebaSBarry Smith Note:
506dce8aebaSBarry Smith `DMCompositeScatterArray()` is a non-variadic alternative that is often more convenient for library callers and is
5076f3c3dcfSJed Brown accessible from Fortran.
5086f3c3dcfSJed Brown
509dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
510db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
511db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
512db781477SPatrick Sanan `DMCompositeScatterArray()`
51347c6ae99SBarry Smith @*/
DMCompositeScatter(DM dm,Vec gvec,...)514d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...)
515d71ae5a4SJacob Faibussowitsch {
51647c6ae99SBarry Smith va_list Argp;
51747c6ae99SBarry Smith struct DMCompositeLink *next;
518c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt;
51947c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
52071b14b3eSStefano Zampini PetscBool flg;
52147c6ae99SBarry Smith
52247c6ae99SBarry Smith PetscFunctionBegin;
52347c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
52447c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5267a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
52748a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
52847c6ae99SBarry Smith
52915229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
53047c6ae99SBarry Smith va_start(Argp, gvec);
5318fd8f222SJed Brown for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
5329ae5db72SJed Brown Vec local;
5339ae5db72SJed Brown local = va_arg(Argp, Vec);
5349ae5db72SJed Brown if (local) {
5359ae5db72SJed Brown Vec global;
5365edff71fSBarry Smith const PetscScalar *array;
53763a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
5389566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global));
5399566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array));
5409566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart));
5419566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local));
5429566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local));
5439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array));
5449566063dSJacob Faibussowitsch PetscCall(VecResetArray(global));
5459566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global));
54647c6ae99SBarry Smith }
54747c6ae99SBarry Smith }
54847c6ae99SBarry Smith va_end(Argp);
5493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
55047c6ae99SBarry Smith }
55147c6ae99SBarry Smith
5526f3c3dcfSJed Brown /*@
5536f3c3dcfSJed Brown DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors
5546f3c3dcfSJed Brown
55520f4b53cSBarry Smith Collective
5566f3c3dcfSJed Brown
5576f3c3dcfSJed Brown Input Parameters:
558dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object
5596f3c3dcfSJed Brown . gvec - the global vector
560a2b725a8SWilliam Gropp - lvecs - array of local vectors, NULL for any that are not needed
5616f3c3dcfSJed Brown
5626f3c3dcfSJed Brown Level: advanced
5636f3c3dcfSJed Brown
5646f3c3dcfSJed Brown Note:
565dce8aebaSBarry Smith This is a non-variadic alternative to `DMCompositeScatter()`
5666f3c3dcfSJed Brown
567dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`
568db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
569db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
5706f3c3dcfSJed Brown @*/
DMCompositeScatterArray(DM dm,Vec gvec,Vec * lvecs)571d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs)
572d71ae5a4SJacob Faibussowitsch {
5736f3c3dcfSJed Brown struct DMCompositeLink *next;
5746f3c3dcfSJed Brown PetscInt i;
5756f3c3dcfSJed Brown DM_Composite *com = (DM_Composite *)dm->data;
57671b14b3eSStefano Zampini PetscBool flg;
5776f3c3dcfSJed Brown
5786f3c3dcfSJed Brown PetscFunctionBegin;
5796f3c3dcfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5806f3c3dcfSJed Brown PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
5819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
5827a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
58348a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
5846f3c3dcfSJed Brown
58515229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
5866f3c3dcfSJed Brown for (i = 0, next = com->next; next; next = next->next, i++) {
5876f3c3dcfSJed Brown if (lvecs[i]) {
5886f3c3dcfSJed Brown Vec global;
589c5d31e75SLisandro Dalcin const PetscScalar *array;
5906f3c3dcfSJed Brown PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3);
5919566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global));
5929566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array));
5939566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart));
5949566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i]));
5959566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i]));
5969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array));
5979566063dSJacob Faibussowitsch PetscCall(VecResetArray(global));
5989566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global));
5996f3c3dcfSJed Brown }
6006f3c3dcfSJed Brown }
6013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6026f3c3dcfSJed Brown }
6036f3c3dcfSJed Brown
60447c6ae99SBarry Smith /*@C
60547c6ae99SBarry Smith DMCompositeGather - Gathers into a global packed vector from its individual local vectors
60647c6ae99SBarry Smith
60720f4b53cSBarry Smith Collective
60847c6ae99SBarry Smith
609d8d19677SJose E. Roman Input Parameters:
610dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object
611dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES`
612a4e35b19SJacob Faibussowitsch . gvec - the global vector
613a4e35b19SJacob Faibussowitsch - ... - the individual sequential vectors, `NULL` for any that are not needed
61447c6ae99SBarry Smith
61547c6ae99SBarry Smith Level: advanced
61647c6ae99SBarry Smith
61760225df5SJacob Faibussowitsch Fortran Notes:
618dce8aebaSBarry Smith Fortran users should use `DMCompositeGatherArray()`
619f5f57ec0SBarry Smith
620dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
621db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
622db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
62347c6ae99SBarry Smith @*/
DMCompositeGather(DM dm,InsertMode imode,Vec gvec,...)624d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...)
625d71ae5a4SJacob Faibussowitsch {
62647c6ae99SBarry Smith va_list Argp;
62747c6ae99SBarry Smith struct DMCompositeLink *next;
62847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
629c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt;
63071b14b3eSStefano Zampini PetscBool flg;
63147c6ae99SBarry Smith
63247c6ae99SBarry Smith PetscFunctionBegin;
63347c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
634064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6359566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6367a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
63748a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
63847c6ae99SBarry Smith
63915229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
6401dac896bSSatish Balay va_start(Argp, gvec);
6418fd8f222SJed Brown for (cnt = 3, next = com->next; next; cnt++, next = next->next) {
6429ae5db72SJed Brown Vec local;
6439ae5db72SJed Brown local = va_arg(Argp, Vec);
6449ae5db72SJed Brown if (local) {
64547c6ae99SBarry Smith PetscScalar *array;
6469ae5db72SJed Brown Vec global;
64763a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt));
6489566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global));
6499566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array));
6509566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart));
6519566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global));
6529566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global));
6539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array));
6549566063dSJacob Faibussowitsch PetscCall(VecResetArray(global));
6559566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global));
65647c6ae99SBarry Smith }
65747c6ae99SBarry Smith }
65847c6ae99SBarry Smith va_end(Argp);
6593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
66047c6ae99SBarry Smith }
66147c6ae99SBarry Smith
6626f3c3dcfSJed Brown /*@
6636f3c3dcfSJed Brown DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors
6646f3c3dcfSJed Brown
66520f4b53cSBarry Smith Collective
6666f3c3dcfSJed Brown
667d8d19677SJose E. Roman Input Parameters:
668dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object
6696f3c3dcfSJed Brown . gvec - the global vector
670dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES`
6716f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed
6726f3c3dcfSJed Brown
6736f3c3dcfSJed Brown Level: advanced
6746f3c3dcfSJed Brown
675dce8aebaSBarry Smith Note:
676dce8aebaSBarry Smith This is a non-variadic alternative to `DMCompositeGather()`.
6776f3c3dcfSJed Brown
678dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
679db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
680db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`,
6816f3c3dcfSJed Brown @*/
DMCompositeGatherArray(DM dm,InsertMode imode,Vec gvec,Vec * lvecs)682d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs)
683d71ae5a4SJacob Faibussowitsch {
6846f3c3dcfSJed Brown struct DMCompositeLink *next;
6856f3c3dcfSJed Brown DM_Composite *com = (DM_Composite *)dm->data;
6866f3c3dcfSJed Brown PetscInt i;
68771b14b3eSStefano Zampini PetscBool flg;
6886f3c3dcfSJed Brown
6896f3c3dcfSJed Brown PetscFunctionBegin;
6906f3c3dcfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
691064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3);
6929566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
6937a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
69448a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
6956f3c3dcfSJed Brown
69615229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
6976f3c3dcfSJed Brown for (next = com->next, i = 0; next; next = next->next, i++) {
6986f3c3dcfSJed Brown if (lvecs[i]) {
6996f3c3dcfSJed Brown PetscScalar *array;
7006f3c3dcfSJed Brown Vec global;
701064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4);
7029566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global));
7039566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array));
7049566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart));
7059566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global));
7069566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global));
7079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array));
7089566063dSJacob Faibussowitsch PetscCall(VecResetArray(global));
7099566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global));
7106f3c3dcfSJed Brown }
7116f3c3dcfSJed Brown }
7123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7136f3c3dcfSJed Brown }
7146f3c3dcfSJed Brown
715f5f57ec0SBarry Smith /*@
716dce8aebaSBarry Smith DMCompositeAddDM - adds a `DM` vector to a `DMCOMPOSITE`
71747c6ae99SBarry Smith
71820f4b53cSBarry Smith Collective
71947c6ae99SBarry Smith
720d8d19677SJose E. Roman Input Parameters:
721dce8aebaSBarry Smith + dmc - the `DMCOMPOSITE` object
722dce8aebaSBarry Smith - dm - the `DM` object
72347c6ae99SBarry Smith
72447c6ae99SBarry Smith Level: advanced
72547c6ae99SBarry Smith
72642747ad1SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeGather()`, `DMCreateGlobalVector()`,
727db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
728db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
72947c6ae99SBarry Smith @*/
DMCompositeAddDM(DM dmc,DM dm)730d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeAddDM(DM dmc, DM dm)
731d71ae5a4SJacob Faibussowitsch {
73206ebdd98SJed Brown PetscInt n, nlocal;
73347c6ae99SBarry Smith struct DMCompositeLink *mine, *next;
73406ebdd98SJed Brown Vec global, local;
73547c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dmc->data;
73671b14b3eSStefano Zampini PetscBool flg;
73747c6ae99SBarry Smith
73847c6ae99SBarry Smith PetscFunctionBegin;
73947c6ae99SBarry Smith PetscValidHeaderSpecific(dmc, DM_CLASSID, 1);
74047c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
7419566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg));
7427a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
74347c6ae99SBarry Smith next = com->next;
74428b400f6SJacob Faibussowitsch PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite");
74547c6ae99SBarry Smith
74647c6ae99SBarry Smith /* create new link */
7479566063dSJacob Faibussowitsch PetscCall(PetscNew(&mine));
7489566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm));
749cac35740SJunchao Zhang PetscCall(DMCreateGlobalVector(dm, &global)); // Not using DMGetGlobalVector() since it will stash the vector with a type decided by dm,
750cac35740SJunchao Zhang PetscCall(VecGetLocalSize(global, &n)); // while we may want to set dmc's vectype later, say via DMSetFromOptions(dmc), and we
751cac35740SJunchao Zhang PetscCall(VecDestroy(&global)); // want to propagate the type to dm.
752cac35740SJunchao Zhang PetscCall(DMCreateLocalVector(dm, &local)); // Not using DMGetLocalVector(), same reason as above.
7539566063dSJacob Faibussowitsch PetscCall(VecGetSize(local, &nlocal));
754cac35740SJunchao Zhang PetscCall(VecDestroy(&local));
7558865f1eaSKarl Rupp
75647c6ae99SBarry Smith mine->n = n;
75706ebdd98SJed Brown mine->nlocal = nlocal;
75847c6ae99SBarry Smith mine->dm = dm;
7590298fd71SBarry Smith mine->next = NULL;
76047c6ae99SBarry Smith com->n += n;
7617ac2b803SAlex Fikl com->nghost += nlocal;
76247c6ae99SBarry Smith
76347c6ae99SBarry Smith /* add to end of list */
7648865f1eaSKarl Rupp if (!next) com->next = mine;
7658865f1eaSKarl Rupp else {
76647c6ae99SBarry Smith while (next->next) next = next->next;
76747c6ae99SBarry Smith next->next = mine;
76847c6ae99SBarry Smith }
76947c6ae99SBarry Smith com->nDM++;
77047c6ae99SBarry Smith com->nmine++;
7713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
77247c6ae99SBarry Smith }
77347c6ae99SBarry Smith
7749804daf3SBarry Smith #include <petscdraw.h>
775d6acfc2dSPierre Jolivet PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
776d6acfc2dSPierre Jolivet
VecView_DMComposite(Vec gvec,PetscViewer viewer)77766976f2fSJacob Faibussowitsch static PetscErrorCode VecView_DMComposite(Vec gvec, PetscViewer viewer)
778d71ae5a4SJacob Faibussowitsch {
77947c6ae99SBarry Smith DM dm;
78047c6ae99SBarry Smith struct DMCompositeLink *next;
78147c6ae99SBarry Smith PetscBool isdraw;
782cef07954SSatish Balay DM_Composite *com;
78347c6ae99SBarry Smith
78447c6ae99SBarry Smith PetscFunctionBegin;
7859566063dSJacob Faibussowitsch PetscCall(VecGetDM(gvec, &dm));
7867a8be351SBarry Smith PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite");
78747c6ae99SBarry Smith com = (DM_Composite *)dm->data;
78847c6ae99SBarry Smith next = com->next;
78947c6ae99SBarry Smith
7909566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
79147c6ae99SBarry Smith if (!isdraw) {
79247c6ae99SBarry Smith /* do I really want to call this? */
7939566063dSJacob Faibussowitsch PetscCall(VecView_MPI(gvec, viewer));
79447c6ae99SBarry Smith } else {
79547c6ae99SBarry Smith PetscInt cnt = 0;
79647c6ae99SBarry Smith
79715229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
79847c6ae99SBarry Smith while (next) {
79947c6ae99SBarry Smith Vec vec;
800ca4278abSLisandro Dalcin const PetscScalar *array;
80147c6ae99SBarry Smith PetscInt bs;
80247c6ae99SBarry Smith
8039ae5db72SJed Brown /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */
8049566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &vec));
8059566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array));
8069566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart));
8079566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array));
8089566063dSJacob Faibussowitsch PetscCall(VecView(vec, viewer));
8099566063dSJacob Faibussowitsch PetscCall(VecResetArray(vec));
8109566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vec, &bs));
8119566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &vec));
8129566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer, bs));
81347c6ae99SBarry Smith cnt += bs;
81447c6ae99SBarry Smith next = next->next;
81547c6ae99SBarry Smith }
8169566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt));
81747c6ae99SBarry Smith }
8183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
81947c6ae99SBarry Smith }
82047c6ae99SBarry Smith
DMCreateGlobalVector_Composite(DM dm,Vec * gvec)82166976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec)
822d71ae5a4SJacob Faibussowitsch {
82347c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
82447c6ae99SBarry Smith
82547c6ae99SBarry Smith PetscFunctionBegin;
82647c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8279566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm));
8289566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm));
8299566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec));
8309566063dSJacob Faibussowitsch PetscCall(VecSetType(*gvec, dm->vectype));
8319566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*gvec, com->n, com->N));
8329566063dSJacob Faibussowitsch PetscCall(VecSetDM(*gvec, dm));
83357d50842SBarry Smith PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (PetscErrorCodeFn *)VecView_DMComposite));
8343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
83547c6ae99SBarry Smith }
83647c6ae99SBarry Smith
DMCreateLocalVector_Composite(DM dm,Vec * lvec)83766976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec)
838d71ae5a4SJacob Faibussowitsch {
83947c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
84047c6ae99SBarry Smith
84147c6ae99SBarry Smith PetscFunctionBegin;
84247c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
84347c6ae99SBarry Smith if (!com->setup) {
8449566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm));
8459566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm));
84647c6ae99SBarry Smith }
8479566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, lvec));
8489566063dSJacob Faibussowitsch PetscCall(VecSetType(*lvec, dm->vectype));
8499566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE));
8509566063dSJacob Faibussowitsch PetscCall(VecSetDM(*lvec, dm));
8513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
85247c6ae99SBarry Smith }
85347c6ae99SBarry Smith
85447c6ae99SBarry Smith /*@C
855dce8aebaSBarry Smith DMCompositeGetISLocalToGlobalMappings - gets an `ISLocalToGlobalMapping` for each `DM` in the `DMCOMPOSITE`, maps to the composite global space
85647c6ae99SBarry Smith
85720f4b53cSBarry Smith Collective; No Fortran Support
85847c6ae99SBarry Smith
85947c6ae99SBarry Smith Input Parameter:
860dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
86147c6ae99SBarry Smith
8622fe279fdSBarry Smith Output Parameter:
8639ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes
864dce8aebaSBarry Smith all the ghost points that individual ghosted `DMDA` may have.
86547c6ae99SBarry Smith
86647c6ae99SBarry Smith Level: advanced
86747c6ae99SBarry Smith
868dce8aebaSBarry Smith Note:
869f13dfd9eSBarry Smith Each entry of `ltogs` should be destroyed with `ISLocalToGlobalMappingDestroy()`, `ltogs` should be freed with `PetscFree()`.
87047c6ae99SBarry Smith
871dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
872db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
873c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
87447c6ae99SBarry Smith @*/
DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping * ltogs[])875cc4c1da9SBarry Smith PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping *ltogs[])
876d71ae5a4SJacob Faibussowitsch {
87747c6ae99SBarry Smith PetscInt i, *idx, n, cnt;
87847c6ae99SBarry Smith struct DMCompositeLink *next;
87947c6ae99SBarry Smith PetscMPIInt rank;
88047c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
88171b14b3eSStefano Zampini PetscBool flg;
88247c6ae99SBarry Smith
88347c6ae99SBarry Smith PetscFunctionBegin;
88447c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8859566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
8867a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
8879566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm));
8889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM, ltogs));
88947c6ae99SBarry Smith next = com->next;
8909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
89147c6ae99SBarry Smith
89215229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
89347c6ae99SBarry Smith cnt = 0;
89447c6ae99SBarry Smith while (next) {
8956eb61c8cSJed Brown ISLocalToGlobalMapping ltog;
8966eb61c8cSJed Brown PetscMPIInt size;
89786994e45SJed Brown const PetscInt *suboff, *indices;
8986eb61c8cSJed Brown Vec global;
89947c6ae99SBarry Smith
9006eb61c8cSJed Brown /* Get sub-DM global indices for each local dof */
9019566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(next->dm, <og));
9029566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n));
9039566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices));
9049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idx));
90547c6ae99SBarry Smith
9066eb61c8cSJed Brown /* Get the offsets for the sub-DM global vector */
9079566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global));
9089566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(global, &suboff));
9099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size));
9106eb61c8cSJed Brown
9116eb61c8cSJed Brown /* Shift the sub-DM definition of the global space to the composite global space */
9126eb61c8cSJed Brown for (i = 0; i < n; i++) {
91386994e45SJed Brown PetscInt subi = indices[i], lo = 0, hi = size, t;
914e333b035SStefano Zampini /* There's no consensus on what a negative index means,
915e333b035SStefano Zampini except for skipping when setting the values in vectors and matrices */
9169371c9d4SSatish Balay if (subi < 0) {
9179371c9d4SSatish Balay idx[i] = subi - next->grstarts[rank];
9189371c9d4SSatish Balay continue;
9199371c9d4SSatish Balay }
9206eb61c8cSJed Brown /* Binary search to find which rank owns subi */
9216eb61c8cSJed Brown while (hi - lo > 1) {
9226eb61c8cSJed Brown t = lo + (hi - lo) / 2;
9236eb61c8cSJed Brown if (suboff[t] > subi) hi = t;
9246eb61c8cSJed Brown else lo = t;
9256eb61c8cSJed Brown }
9266eb61c8cSJed Brown idx[i] = subi - suboff[lo] + next->grstarts[lo];
9276eb61c8cSJed Brown }
9289566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices));
9299566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt]));
9309566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global));
93147c6ae99SBarry Smith next = next->next;
93247c6ae99SBarry Smith cnt++;
93347c6ae99SBarry Smith }
9343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
93547c6ae99SBarry Smith }
93647c6ae99SBarry Smith
93787c85e80SJed Brown /*@C
9389ae5db72SJed Brown DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector
93987c85e80SJed Brown
94020f4b53cSBarry Smith Not Collective; No Fortran Support
94187c85e80SJed Brown
9424165533cSJose E. Roman Input Parameter:
943dce8aebaSBarry Smith . dm - the `DMCOMPOSITE`
94487c85e80SJed Brown
9454165533cSJose E. Roman Output Parameter:
94615229ffcSPierre Jolivet . is - array of serial index sets for each component of the `DMCOMPOSITE`
94787c85e80SJed Brown
94887c85e80SJed Brown Level: intermediate
94987c85e80SJed Brown
95087c85e80SJed Brown Notes:
95187c85e80SJed Brown At present, a composite local vector does not normally exist. This function is used to provide index sets for
95215229ffcSPierre Jolivet `MatGetLocalSubMatrix()`. In the future, the scatters for each entry in the `DMCOMPOSITE` may be merged into a single
9539ae5db72SJed Brown scatter to a composite local vector. The user should not typically need to know which is being done.
95487c85e80SJed Brown
955dce8aebaSBarry Smith To get the composite global indices at all local points (including ghosts), use `DMCompositeGetISLocalToGlobalMappings()`.
95687c85e80SJed Brown
957dce8aebaSBarry Smith To get index sets for pieces of the composite global vector, use `DMCompositeGetGlobalISs()`.
95887c85e80SJed Brown
959dce8aebaSBarry Smith Each returned `IS` should be destroyed with `ISDestroy()`, the array should be freed with `PetscFree()`.
96087c85e80SJed Brown
961f13dfd9eSBarry Smith Fortran Note:
962ce78bad3SBarry Smith Use `DMCompositeRestoreLocalISs()` to release the `is`.
963f13dfd9eSBarry Smith
964f13dfd9eSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`,
965f13dfd9eSBarry Smith `MatCreateLocalRef()`, `DMCompositeGetNumberDM()`
96687c85e80SJed Brown @*/
DMCompositeGetLocalISs(DM dm,IS * is[])967cc4c1da9SBarry Smith PetscErrorCode DMCompositeGetLocalISs(DM dm, IS *is[])
968d71ae5a4SJacob Faibussowitsch {
96987c85e80SJed Brown DM_Composite *com = (DM_Composite *)dm->data;
97087c85e80SJed Brown struct DMCompositeLink *link;
97187c85e80SJed Brown PetscInt cnt, start;
97271b14b3eSStefano Zampini PetscBool flg;
97387c85e80SJed Brown
97487c85e80SJed Brown PetscFunctionBegin;
97587c85e80SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9764f572ea9SToby Isaac PetscAssertPointer(is, 2);
9779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
9787a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
9799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nmine, is));
98006ebdd98SJed Brown for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) {
981520db06cSJed Brown PetscInt bs;
9829566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt]));
9839566063dSJacob Faibussowitsch PetscCall(DMGetBlockSize(link->dm, &bs));
9849566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize((*is)[cnt], bs));
985520db06cSJed Brown }
9863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
98787c85e80SJed Brown }
98887c85e80SJed Brown
98947c6ae99SBarry Smith /*@C
990dce8aebaSBarry Smith DMCompositeGetGlobalISs - Gets the index sets for each composed object in a `DMCOMPOSITE`
99147c6ae99SBarry Smith
99220f4b53cSBarry Smith Collective
99347c6ae99SBarry Smith
99447c6ae99SBarry Smith Input Parameter:
995dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
99647c6ae99SBarry Smith
9972fe279fdSBarry Smith Output Parameter:
99847c6ae99SBarry Smith . is - the array of index sets
99947c6ae99SBarry Smith
100047c6ae99SBarry Smith Level: advanced
100147c6ae99SBarry Smith
100247c6ae99SBarry Smith Notes:
1003f13dfd9eSBarry Smith The `is` entries should be destroyed with `ISDestroy()`, `is` should be freed with `PetscFree()`
100447c6ae99SBarry Smith
100547c6ae99SBarry Smith These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
100647c6ae99SBarry Smith
1007dce8aebaSBarry Smith Use `DMCompositeGetLocalISs()` for index sets in the packed local numbering, and
1008dce8aebaSBarry Smith `DMCompositeGetISLocalToGlobalMappings()` for to map local sub-`DM` (including ghost) indices to packed global
10096eb61c8cSJed Brown indices.
101047c6ae99SBarry Smith
1011ce78bad3SBarry Smith Fortran Note:
1012ce78bad3SBarry Smith Use `DMCompositeRestoreGlobalISs()` to release the `is`.
1013f3cb0f7eSJed Brown
1014dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1015db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`,
1016c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
101747c6ae99SBarry Smith @*/
DMCompositeGetGlobalISs(DM dm,IS * is[])1018d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[])
1019d71ae5a4SJacob Faibussowitsch {
102066bb578eSMark Adams PetscInt cnt = 0;
102147c6ae99SBarry Smith struct DMCompositeLink *next;
102247c6ae99SBarry Smith PetscMPIInt rank;
102347c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
102471b14b3eSStefano Zampini PetscBool flg;
102547c6ae99SBarry Smith
102647c6ae99SBarry Smith PetscFunctionBegin;
102747c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
10297a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
10307a8be351SBarry Smith PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before");
10319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM, is));
103247c6ae99SBarry Smith next = com->next;
10339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
103447c6ae99SBarry Smith
103515229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
103647c6ae99SBarry Smith while (next) {
1037e5e52638SMatthew G. Knepley PetscDS prob;
1038e5e52638SMatthew G. Knepley
10399566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt]));
10409566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob));
1041e5e52638SMatthew G. Knepley if (prob) {
104265c226d8SMatthew G. Knepley MatNullSpace space;
104365c226d8SMatthew G. Knepley Mat pmat;
10440f21e855SMatthew G. Knepley PetscObject disc;
10450f21e855SMatthew G. Knepley PetscInt Nf;
104665c226d8SMatthew G. Knepley
10479566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf));
1048f24dd8d2SMatthew G. Knepley if (cnt < Nf) {
10499566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(prob, cnt, &disc));
10509566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space));
10519566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space));
10529566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space));
10539566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space));
10549566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat));
10559566063dSJacob Faibussowitsch if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat));
105665c226d8SMatthew G. Knepley }
1057f24dd8d2SMatthew G. Knepley }
105847c6ae99SBarry Smith cnt++;
105947c6ae99SBarry Smith next = next->next;
106047c6ae99SBarry Smith }
10613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
106247c6ae99SBarry Smith }
106347c6ae99SBarry Smith
DMCreateFieldIS_Composite(DM dm,PetscInt * numFields,char *** fieldNames,IS ** fields)106466976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
1065d71ae5a4SJacob Faibussowitsch {
10664d343eeaSMatthew G Knepley PetscInt nDM;
10674d343eeaSMatthew G Knepley DM *dms;
10684d343eeaSMatthew G Knepley PetscInt i;
10694d343eeaSMatthew G Knepley
10704d343eeaSMatthew G Knepley PetscFunctionBegin;
10719566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM));
10728865f1eaSKarl Rupp if (numFields) *numFields = nDM;
10739566063dSJacob Faibussowitsch PetscCall(DMCompositeGetGlobalISs(dm, fields));
10744d343eeaSMatthew G Knepley if (fieldNames) {
10759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, &dms));
10769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, fieldNames));
10779566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms));
10784d343eeaSMatthew G Knepley for (i = 0; i < nDM; i++) {
10794d343eeaSMatthew G Knepley char buf[256];
10804d343eeaSMatthew G Knepley const char *splitname;
10814d343eeaSMatthew G Knepley
10824d343eeaSMatthew G Knepley /* Split naming precedence: object name, prefix, number */
10834d343eeaSMatthew G Knepley splitname = ((PetscObject)dm)->name;
10844d343eeaSMatthew G Knepley if (!splitname) {
10859566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname));
10864d343eeaSMatthew G Knepley if (splitname) {
10874d343eeaSMatthew G Knepley size_t len;
10889566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(buf, splitname, sizeof(buf)));
10898caf3d72SBarry Smith buf[sizeof(buf) - 1] = 0;
10909566063dSJacob Faibussowitsch PetscCall(PetscStrlen(buf, &len));
10914d343eeaSMatthew G Knepley if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */
10924d343eeaSMatthew G Knepley splitname = buf;
10934d343eeaSMatthew G Knepley }
10944d343eeaSMatthew G Knepley }
10954d343eeaSMatthew G Knepley if (!splitname) {
109663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i));
10974d343eeaSMatthew G Knepley splitname = buf;
10984d343eeaSMatthew G Knepley }
10999566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i]));
11004d343eeaSMatthew G Knepley }
11019566063dSJacob Faibussowitsch PetscCall(PetscFree(dms));
11024d343eeaSMatthew G Knepley }
11033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11044d343eeaSMatthew G Knepley }
11054d343eeaSMatthew G Knepley
1106e7c4fc90SDmitry Karpeev /*
1107e7c4fc90SDmitry Karpeev This could take over from DMCreateFieldIS(), as it is more general,
11080298fd71SBarry Smith making DMCreateFieldIS() a special case -- calling with dmlist == NULL;
1109e7c4fc90SDmitry Karpeev At this point it's probably best to be less intrusive, however.
1110e7c4fc90SDmitry Karpeev */
DMCreateFieldDecomposition_Composite(DM dm,PetscInt * len,char *** namelist,IS ** islist,DM ** dmlist)111166976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist)
1112d71ae5a4SJacob Faibussowitsch {
1113e7c4fc90SDmitry Karpeev PetscInt nDM;
1114e7c4fc90SDmitry Karpeev PetscInt i;
1115e7c4fc90SDmitry Karpeev
1116e7c4fc90SDmitry Karpeev PetscFunctionBegin;
11179566063dSJacob Faibussowitsch PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist));
1118e7c4fc90SDmitry Karpeev if (dmlist) {
11199566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM));
11209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, dmlist));
11219566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, *dmlist));
112248a46eb9SPierre Jolivet for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i])));
1123e7c4fc90SDmitry Karpeev }
11243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1125e7c4fc90SDmitry Karpeev }
1126e7c4fc90SDmitry Karpeev
112747c6ae99SBarry Smith /*@C
1128dce8aebaSBarry Smith DMCompositeGetLocalVectors - Gets local vectors for each part of a `DMCOMPOSITE`
1129dce8aebaSBarry Smith Use `DMCompositeRestoreLocalVectors()` to return them.
113047c6ae99SBarry Smith
113120f4b53cSBarry Smith Not Collective; No Fortran Support
113247c6ae99SBarry Smith
113347c6ae99SBarry Smith Input Parameter:
1134dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
113547c6ae99SBarry Smith
113647c6ae99SBarry Smith Output Parameter:
1137a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`s
113847c6ae99SBarry Smith
113947c6ae99SBarry Smith Level: advanced
114047c6ae99SBarry Smith
1141dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1142db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1143db781477SPatrick Sanan `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
114447c6ae99SBarry Smith @*/
DMCompositeGetLocalVectors(DM dm,...)1145d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...)
1146d71ae5a4SJacob Faibussowitsch {
114747c6ae99SBarry Smith va_list Argp;
114847c6ae99SBarry Smith struct DMCompositeLink *next;
114947c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
115071b14b3eSStefano Zampini PetscBool flg;
115147c6ae99SBarry Smith
115247c6ae99SBarry Smith PetscFunctionBegin;
115347c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11557a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
115647c6ae99SBarry Smith next = com->next;
115715229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
115847c6ae99SBarry Smith va_start(Argp, dm);
115947c6ae99SBarry Smith while (next) {
116047c6ae99SBarry Smith Vec *vec;
116147c6ae99SBarry Smith vec = va_arg(Argp, Vec *);
11629566063dSJacob Faibussowitsch if (vec) PetscCall(DMGetLocalVector(next->dm, vec));
116347c6ae99SBarry Smith next = next->next;
116447c6ae99SBarry Smith }
116547c6ae99SBarry Smith va_end(Argp);
11663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
116747c6ae99SBarry Smith }
116847c6ae99SBarry Smith
116947c6ae99SBarry Smith /*@C
1170dce8aebaSBarry Smith DMCompositeRestoreLocalVectors - Restores local vectors for each part of a `DMCOMPOSITE`
117147c6ae99SBarry Smith
117220f4b53cSBarry Smith Not Collective; No Fortran Support
117347c6ae99SBarry Smith
117447c6ae99SBarry Smith Input Parameter:
1175dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
117647c6ae99SBarry Smith
117747c6ae99SBarry Smith Output Parameter:
1178a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`
117947c6ae99SBarry Smith
118047c6ae99SBarry Smith Level: advanced
118147c6ae99SBarry Smith
1182dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`,
1183db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
1184db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()`
118547c6ae99SBarry Smith @*/
DMCompositeRestoreLocalVectors(DM dm,...)1186d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...)
1187d71ae5a4SJacob Faibussowitsch {
118847c6ae99SBarry Smith va_list Argp;
118947c6ae99SBarry Smith struct DMCompositeLink *next;
119047c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
119171b14b3eSStefano Zampini PetscBool flg;
119247c6ae99SBarry Smith
119347c6ae99SBarry Smith PetscFunctionBegin;
119447c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11959566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
11967a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
119747c6ae99SBarry Smith next = com->next;
119815229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
119947c6ae99SBarry Smith va_start(Argp, dm);
120047c6ae99SBarry Smith while (next) {
120147c6ae99SBarry Smith Vec *vec;
120247c6ae99SBarry Smith vec = va_arg(Argp, Vec *);
12039566063dSJacob Faibussowitsch if (vec) PetscCall(DMRestoreLocalVector(next->dm, vec));
120447c6ae99SBarry Smith next = next->next;
120547c6ae99SBarry Smith }
120647c6ae99SBarry Smith va_end(Argp);
12073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
120847c6ae99SBarry Smith }
120947c6ae99SBarry Smith
121047c6ae99SBarry Smith /*@C
1211dce8aebaSBarry Smith DMCompositeGetEntries - Gets the `DM` for each entry in a `DMCOMPOSITE`.
121247c6ae99SBarry Smith
121347c6ae99SBarry Smith Not Collective
121447c6ae99SBarry Smith
121547c6ae99SBarry Smith Input Parameter:
1216dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
121747c6ae99SBarry Smith
121847c6ae99SBarry Smith Output Parameter:
1219a4e35b19SJacob Faibussowitsch . ... - the individual entries `DM`
122047c6ae99SBarry Smith
122147c6ae99SBarry Smith Level: advanced
122247c6ae99SBarry Smith
122360225df5SJacob Faibussowitsch Fortran Notes:
12245d83a8b1SBarry Smith Use `DMCompositeGetEntriesArray()`
1225f5f57ec0SBarry Smith
1226dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()`
1227db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
122860225df5SJacob Faibussowitsch `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`
122947c6ae99SBarry Smith @*/
DMCompositeGetEntries(DM dm,...)1230d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntries(DM dm, ...)
1231d71ae5a4SJacob Faibussowitsch {
123247c6ae99SBarry Smith va_list Argp;
123347c6ae99SBarry Smith struct DMCompositeLink *next;
123447c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
123571b14b3eSStefano Zampini PetscBool flg;
123647c6ae99SBarry Smith
123747c6ae99SBarry Smith PetscFunctionBegin;
123847c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12399566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
12407a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
124147c6ae99SBarry Smith next = com->next;
124215229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
124347c6ae99SBarry Smith va_start(Argp, dm);
124447c6ae99SBarry Smith while (next) {
124547c6ae99SBarry Smith DM *dmn;
124647c6ae99SBarry Smith dmn = va_arg(Argp, DM *);
12479ae5db72SJed Brown if (dmn) *dmn = next->dm;
124847c6ae99SBarry Smith next = next->next;
124947c6ae99SBarry Smith }
125047c6ae99SBarry Smith va_end(Argp);
12513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
125247c6ae99SBarry Smith }
125347c6ae99SBarry Smith
1254cc4c1da9SBarry Smith /*@
1255dce8aebaSBarry Smith DMCompositeGetEntriesArray - Gets the DM for each entry in a `DMCOMPOSITE`
12562fa5ba8aSJed Brown
12572fa5ba8aSJed Brown Not Collective
12582fa5ba8aSJed Brown
12592fa5ba8aSJed Brown Input Parameter:
1260dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object
1261907376e6SBarry Smith
1262907376e6SBarry Smith Output Parameter:
1263dce8aebaSBarry Smith . dms - array of sufficient length (see `DMCompositeGetNumberDM()`) to hold the individual `DM`
12642fa5ba8aSJed Brown
12652fa5ba8aSJed Brown Level: advanced
12662fa5ba8aSJed Brown
1267dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()`
1268db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`,
126960225df5SJacob Faibussowitsch `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`
12702fa5ba8aSJed Brown @*/
DMCompositeGetEntriesArray(DM dm,DM dms[])1271d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntriesArray(DM dm, DM dms[])
1272d71ae5a4SJacob Faibussowitsch {
12732fa5ba8aSJed Brown struct DMCompositeLink *next;
12742fa5ba8aSJed Brown DM_Composite *com = (DM_Composite *)dm->data;
12752fa5ba8aSJed Brown PetscInt i;
127671b14b3eSStefano Zampini PetscBool flg;
12772fa5ba8aSJed Brown
12782fa5ba8aSJed Brown PetscFunctionBegin;
12792fa5ba8aSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12809566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg));
12817a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name);
128215229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
12832fa5ba8aSJed Brown for (next = com->next, i = 0; next; next = next->next, i++) dms[i] = next->dm;
12843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12852fa5ba8aSJed Brown }
12862fa5ba8aSJed Brown
1287e10fd815SStefano Zampini typedef struct {
1288e10fd815SStefano Zampini DM dm;
1289e10fd815SStefano Zampini PetscViewer *subv;
1290e10fd815SStefano Zampini Vec *vecs;
1291e10fd815SStefano Zampini } GLVisViewerCtx;
1292e10fd815SStefano Zampini
DestroyGLVisViewerCtx_Private(PetscCtxRt vctx)1293*2a8381b2SBarry Smith static PetscErrorCode DestroyGLVisViewerCtx_Private(PetscCtxRt vctx)
1294d71ae5a4SJacob Faibussowitsch {
1295e6aa7a3bSBarry Smith GLVisViewerCtx *ctx = *(GLVisViewerCtx **)vctx;
1296e10fd815SStefano Zampini PetscInt i, n;
1297e10fd815SStefano Zampini
1298e10fd815SStefano Zampini PetscFunctionBegin;
12999566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm, &n));
130048a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(PetscViewerDestroy(&ctx->subv[i]));
13019566063dSJacob Faibussowitsch PetscCall(PetscFree2(ctx->subv, ctx->vecs));
13029566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx->dm));
13039566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx));
13043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1305e10fd815SStefano Zampini }
1306e10fd815SStefano Zampini
DMCompositeSampleGLVisFields_Private(PetscObject oX,PetscInt nf,PetscObject oXfield[],void * vctx)1307d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
1308d71ae5a4SJacob Faibussowitsch {
1309e10fd815SStefano Zampini Vec X = (Vec)oX;
1310e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx;
1311e10fd815SStefano Zampini PetscInt i, n, cumf;
1312e10fd815SStefano Zampini
1313e10fd815SStefano Zampini PetscFunctionBegin;
13149566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm, &n));
13159566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(ctx->dm, X, n, NULL, ctx->vecs));
1316e10fd815SStefano Zampini for (i = 0, cumf = 0; i < n; i++) {
1317e10fd815SStefano Zampini PetscErrorCode (*g2l)(PetscObject, PetscInt, PetscObject[], void *);
1318e10fd815SStefano Zampini void *fctx;
1319e10fd815SStefano Zampini PetscInt nfi;
1320e10fd815SStefano Zampini
132134e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nfi, NULL, NULL, &g2l, NULL, &fctx));
1322e10fd815SStefano Zampini if (!nfi) continue;
13231baa6e33SBarry Smith if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i], nfi, oXfield + cumf, fctx));
1324f4f49eeaSPierre Jolivet else PetscCall(VecCopy(ctx->vecs[i], (Vec)oXfield[cumf]));
1325e10fd815SStefano Zampini cumf += nfi;
1326e10fd815SStefano Zampini }
13279566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(ctx->dm, X, n, NULL, ctx->vecs));
13283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1329e10fd815SStefano Zampini }
1330e10fd815SStefano Zampini
DMSetUpGLVisViewer_Composite(PetscObject odm,PetscViewer viewer)1331d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer)
1332d71ae5a4SJacob Faibussowitsch {
1333e10fd815SStefano Zampini DM dm = (DM)odm, *dms;
1334e10fd815SStefano Zampini Vec *Ufds;
1335e10fd815SStefano Zampini GLVisViewerCtx *ctx;
1336e10fd815SStefano Zampini PetscInt i, n, tnf, *sdim;
1337e10fd815SStefano Zampini char **fecs;
1338e10fd815SStefano Zampini
1339e10fd815SStefano Zampini PetscFunctionBegin;
13409566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx));
13419566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm));
1342e10fd815SStefano Zampini ctx->dm = dm;
13439566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &n));
13449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dms));
13459566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms));
13469566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &ctx->subv, n, &ctx->vecs));
1347e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) {
1348e10fd815SStefano Zampini PetscInt nf;
1349e10fd815SStefano Zampini
13509566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PetscObjectComm(odm), &ctx->subv[i]));
13519566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(ctx->subv[i], PETSCVIEWERGLVIS));
135234e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisSetDM_Internal(ctx->subv[i], (PetscObject)dms[i]));
135334e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, NULL, NULL, NULL, NULL, NULL));
1354e10fd815SStefano Zampini tnf += nf;
1355e10fd815SStefano Zampini }
13569566063dSJacob Faibussowitsch PetscCall(PetscFree(dms));
13579566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(tnf, &fecs, tnf, &sdim, tnf, &Ufds));
1358e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) {
1359e10fd815SStefano Zampini PetscInt *sd, nf, f;
1360e10fd815SStefano Zampini const char **fec;
1361e10fd815SStefano Zampini Vec *Uf;
1362e10fd815SStefano Zampini
136334e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, &fec, &sd, NULL, (PetscObject **)&Uf, NULL));
1364e10fd815SStefano Zampini for (f = 0; f < nf; f++) {
13659566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec[f], &fecs[tnf + f]));
1366e10fd815SStefano Zampini Ufds[tnf + f] = Uf[f];
1367e10fd815SStefano Zampini sdim[tnf + f] = sd[f];
1368e10fd815SStefano Zampini }
1369e10fd815SStefano Zampini tnf += nf;
1370e10fd815SStefano Zampini }
13719566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetFields(viewer, tnf, (const char **)fecs, sdim, DMCompositeSampleGLVisFields_Private, (PetscObject *)Ufds, ctx, DestroyGLVisViewerCtx_Private));
137248a46eb9SPierre Jolivet for (i = 0; i < tnf; i++) PetscCall(PetscFree(fecs[i]));
13739566063dSJacob Faibussowitsch PetscCall(PetscFree3(fecs, sdim, Ufds));
13743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1375e10fd815SStefano Zampini }
1376e10fd815SStefano Zampini
DMRefine_Composite(DM dmi,MPI_Comm comm,DM * fine)137766976f2fSJacob Faibussowitsch static PetscErrorCode DMRefine_Composite(DM dmi, MPI_Comm comm, DM *fine)
1378d71ae5a4SJacob Faibussowitsch {
137947c6ae99SBarry Smith struct DMCompositeLink *next;
138047c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dmi->data;
138147c6ae99SBarry Smith DM dm;
138247c6ae99SBarry Smith
138347c6ae99SBarry Smith PetscFunctionBegin;
138447c6ae99SBarry Smith PetscValidHeaderSpecific(dmi, DM_CLASSID, 1);
138548a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm));
13869566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi));
138747c6ae99SBarry Smith next = com->next;
13889566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm, fine));
138947c6ae99SBarry Smith
139015229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
139147c6ae99SBarry Smith while (next) {
13929566063dSJacob Faibussowitsch PetscCall(DMRefine(next->dm, comm, &dm));
13939566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine, dm));
13949566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm));
139547c6ae99SBarry Smith next = next->next;
139647c6ae99SBarry Smith }
13973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
139847c6ae99SBarry Smith }
139947c6ae99SBarry Smith
DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM * fine)140066976f2fSJacob Faibussowitsch static PetscErrorCode DMCoarsen_Composite(DM dmi, MPI_Comm comm, DM *fine)
1401d71ae5a4SJacob Faibussowitsch {
140214354c39SJed Brown struct DMCompositeLink *next;
140314354c39SJed Brown DM_Composite *com = (DM_Composite *)dmi->data;
140414354c39SJed Brown DM dm;
140514354c39SJed Brown
140614354c39SJed Brown PetscFunctionBegin;
140714354c39SJed Brown PetscValidHeaderSpecific(dmi, DM_CLASSID, 1);
14089566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi));
140948a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm));
141014354c39SJed Brown next = com->next;
14119566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm, fine));
141214354c39SJed Brown
141315229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
141414354c39SJed Brown while (next) {
14159566063dSJacob Faibussowitsch PetscCall(DMCoarsen(next->dm, comm, &dm));
14169566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine, dm));
14179566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm));
141814354c39SJed Brown next = next->next;
141914354c39SJed Brown }
14203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
142114354c39SJed Brown }
142247c6ae99SBarry Smith
DMCreateInterpolation_Composite(DM coarse,DM fine,Mat * A,Vec * v)142366976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_Composite(DM coarse, DM fine, Mat *A, Vec *v)
1424d71ae5a4SJacob Faibussowitsch {
14259ae5db72SJed Brown PetscInt m, n, M, N, nDM, i;
142647c6ae99SBarry Smith struct DMCompositeLink *nextc;
142747c6ae99SBarry Smith struct DMCompositeLink *nextf;
142825296bd5SBarry Smith Vec gcoarse, gfine, *vecs;
142947c6ae99SBarry Smith DM_Composite *comcoarse = (DM_Composite *)coarse->data;
143047c6ae99SBarry Smith DM_Composite *comfine = (DM_Composite *)fine->data;
14319ae5db72SJed Brown Mat *mats;
143247c6ae99SBarry Smith
143347c6ae99SBarry Smith PetscFunctionBegin;
143447c6ae99SBarry Smith PetscValidHeaderSpecific(coarse, DM_CLASSID, 1);
143547c6ae99SBarry Smith PetscValidHeaderSpecific(fine, DM_CLASSID, 2);
14369566063dSJacob Faibussowitsch PetscCall(DMSetUp(coarse));
14379566063dSJacob Faibussowitsch PetscCall(DMSetUp(fine));
143847c6ae99SBarry Smith /* use global vectors only for determining matrix layout */
14399566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(coarse, &gcoarse));
14409566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fine, &gfine));
14419566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gcoarse, &n));
14429566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gfine, &m));
14439566063dSJacob Faibussowitsch PetscCall(VecGetSize(gcoarse, &N));
14449566063dSJacob Faibussowitsch PetscCall(VecGetSize(gfine, &M));
14459566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(coarse, &gcoarse));
14469566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fine, &gfine));
144747c6ae99SBarry Smith
14489ae5db72SJed Brown nDM = comfine->nDM;
144963a3b9bcSJacob 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);
14509566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nDM * nDM, &mats));
145148a46eb9SPierre Jolivet if (v) PetscCall(PetscCalloc1(nDM, &vecs));
145247c6ae99SBarry Smith
145315229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
14549ae5db72SJed Brown for (nextc = comcoarse->next, nextf = comfine->next, i = 0; nextc; nextc = nextc->next, nextf = nextf->next, i++) {
14551baa6e33SBarry Smith if (!v) PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], NULL));
14561baa6e33SBarry Smith else PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], &vecs[i]));
145747c6ae99SBarry Smith }
14589566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, nDM, NULL, mats, A));
14591baa6e33SBarry Smith if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, vecs, v));
14609566063dSJacob Faibussowitsch for (i = 0; i < nDM * nDM; i++) PetscCall(MatDestroy(&mats[i]));
14619566063dSJacob Faibussowitsch PetscCall(PetscFree(mats));
146225296bd5SBarry Smith if (v) {
14639566063dSJacob Faibussowitsch for (i = 0; i < nDM; i++) PetscCall(VecDestroy(&vecs[i]));
14649566063dSJacob Faibussowitsch PetscCall(PetscFree(vecs));
146525296bd5SBarry Smith }
14663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
146747c6ae99SBarry Smith }
146847c6ae99SBarry Smith
DMGetLocalToGlobalMapping_Composite(DM dm)1469d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm)
1470d71ae5a4SJacob Faibussowitsch {
14711411c6eeSJed Brown DM_Composite *com = (DM_Composite *)dm->data;
14721411c6eeSJed Brown ISLocalToGlobalMapping *ltogs;
1473f7efa3c7SJed Brown PetscInt i;
14741411c6eeSJed Brown
14751411c6eeSJed Brown PetscFunctionBegin;
14761411c6eeSJed Brown /* Set the ISLocalToGlobalMapping on the new matrix */
14779566063dSJacob Faibussowitsch PetscCall(DMCompositeGetISLocalToGlobalMappings(dm, <ogs));
14789566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm), com->nDM, ltogs, &dm->ltogmap));
14799566063dSJacob Faibussowitsch for (i = 0; i < com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(<ogs[i]));
14809566063dSJacob Faibussowitsch PetscCall(PetscFree(ltogs));
14813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14821411c6eeSJed Brown }
14831411c6eeSJed Brown
DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring * coloring)148466976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateColoring_Composite(DM dm, ISColoringType ctype, ISColoring *coloring)
1485d71ae5a4SJacob Faibussowitsch {
148647c6ae99SBarry Smith PetscInt n, i, cnt;
148747c6ae99SBarry Smith ISColoringValue *colors;
148847c6ae99SBarry Smith PetscBool dense = PETSC_FALSE;
148947c6ae99SBarry Smith ISColoringValue maxcol = 0;
149047c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
149147c6ae99SBarry Smith
149247c6ae99SBarry Smith PetscFunctionBegin;
149347c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
149408401ef6SPierre Jolivet PetscCheck(ctype != IS_COLORING_LOCAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only global coloring supported");
1495966bd95aSPierre Jolivet PetscCheck(ctype == IS_COLORING_GLOBAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unknown ISColoringType");
149647c6ae99SBarry Smith n = com->n;
14979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &colors)); /* freed in ISColoringDestroy() */
149847c6ae99SBarry Smith
14999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL));
150047c6ae99SBarry Smith if (dense) {
15016497c311SBarry Smith PetscCall(ISColoringValueCast(com->N, &maxcol));
15026497c311SBarry Smith for (i = 0; i < n; i++) PetscCall(ISColoringValueCast(com->rstart + i, colors + i));
150347c6ae99SBarry Smith } else {
150447c6ae99SBarry Smith struct DMCompositeLink *next = com->next;
150547c6ae99SBarry Smith PetscMPIInt rank;
150647c6ae99SBarry Smith
15079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
150847c6ae99SBarry Smith cnt = 0;
150947c6ae99SBarry Smith while (next) {
151047c6ae99SBarry Smith ISColoring lcoloring;
151147c6ae99SBarry Smith
15129566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(next->dm, IS_COLORING_GLOBAL, &lcoloring));
15136497c311SBarry Smith for (i = 0; i < lcoloring->N; i++) PetscCall(ISColoringValueCast(maxcol + lcoloring->colors[i], colors + cnt++));
151447c6ae99SBarry Smith maxcol += lcoloring->n;
15159566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&lcoloring));
151647c6ae99SBarry Smith next = next->next;
151747c6ae99SBarry Smith }
151847c6ae99SBarry Smith }
15199566063dSJacob Faibussowitsch PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), maxcol, n, colors, PETSC_OWN_POINTER, coloring));
15203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
152147c6ae99SBarry Smith }
152247c6ae99SBarry Smith
DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)152366976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec)
1524d71ae5a4SJacob Faibussowitsch {
152547c6ae99SBarry Smith struct DMCompositeLink *next;
152651a5d4e4SZach Atkins const PetscScalar *garray, *garrayhead;
152751a5d4e4SZach Atkins PetscScalar *larray, *larrayhead;
152847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data;
152947c6ae99SBarry Smith
153047c6ae99SBarry Smith PetscFunctionBegin;
153147c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
153247c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
153339d35262SVincent Le Chenadec
153448a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
153539d35262SVincent Le Chenadec
153651a5d4e4SZach Atkins PetscCall(VecGetArrayRead(gvec, &garray));
153751a5d4e4SZach Atkins garrayhead = garray;
15389566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec, &larray));
153951a5d4e4SZach Atkins larrayhead = larray;
154047c6ae99SBarry Smith
154115229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
154239d35262SVincent Le Chenadec next = com->next;
154347c6ae99SBarry Smith while (next) {
154447c6ae99SBarry Smith Vec local, global;
154547c6ae99SBarry Smith PetscInt N;
154647c6ae99SBarry Smith
15479566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global));
15489566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global, &N));
154951a5d4e4SZach Atkins PetscCall(VecPlaceArray(global, garrayhead));
15509566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local));
155151a5d4e4SZach Atkins PetscCall(VecPlaceArray(local, larrayhead));
15529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, mode, local));
15539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, mode, local));
15549566063dSJacob Faibussowitsch PetscCall(VecResetArray(global));
15559566063dSJacob Faibussowitsch PetscCall(VecResetArray(local));
15569566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global));
15579566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local));
155839d35262SVincent Le Chenadec
155951a5d4e4SZach Atkins larrayhead += next->nlocal;
156051a5d4e4SZach Atkins garrayhead += next->n;
156147c6ae99SBarry Smith next = next->next;
156247c6ae99SBarry Smith }
156351a5d4e4SZach Atkins PetscCall(VecRestoreArrayRead(gvec, &garray));
156451a5d4e4SZach Atkins PetscCall(VecRestoreArray(lvec, &larray));
15653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
156647c6ae99SBarry Smith }
156747c6ae99SBarry Smith
DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)156866976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec)
1569d71ae5a4SJacob Faibussowitsch {
15700c010503SBarry Smith PetscFunctionBegin;
157139d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
157239d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2);
157339d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 4);
15743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
157539d35262SVincent Le Chenadec }
157639d35262SVincent Le Chenadec
DMLocalToGlobalBegin_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)157766976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1578d71ae5a4SJacob Faibussowitsch {
157939d35262SVincent Le Chenadec struct DMCompositeLink *next;
158051a5d4e4SZach Atkins const PetscScalar *larray, *larrayhead;
158151a5d4e4SZach Atkins PetscScalar *garray, *garrayhead;
158239d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite *)dm->data;
158339d35262SVincent Le Chenadec
158439d35262SVincent Le Chenadec PetscFunctionBegin;
158539d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
158639d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
158739d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
158839d35262SVincent Le Chenadec
158948a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
159039d35262SVincent Le Chenadec
159151a5d4e4SZach Atkins PetscCall(VecGetArrayRead(lvec, &larray));
159251a5d4e4SZach Atkins larrayhead = larray;
15939566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &garray));
159451a5d4e4SZach Atkins garrayhead = garray;
159539d35262SVincent Le Chenadec
159615229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
159739d35262SVincent Le Chenadec next = com->next;
159839d35262SVincent Le Chenadec while (next) {
159939d35262SVincent Le Chenadec Vec global, local;
160039d35262SVincent Le Chenadec
16019566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local));
160251a5d4e4SZach Atkins PetscCall(VecPlaceArray(local, larrayhead));
16039566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global));
160451a5d4e4SZach Atkins PetscCall(VecPlaceArray(global, garrayhead));
16059566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, local, mode, global));
16069566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, local, mode, global));
16079566063dSJacob Faibussowitsch PetscCall(VecResetArray(local));
16089566063dSJacob Faibussowitsch PetscCall(VecResetArray(global));
16099566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global));
16109566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local));
161139d35262SVincent Le Chenadec
161251a5d4e4SZach Atkins garrayhead += next->n;
161351a5d4e4SZach Atkins larrayhead += next->nlocal;
161439d35262SVincent Le Chenadec next = next->next;
161539d35262SVincent Le Chenadec }
161639d35262SVincent Le Chenadec
161751a5d4e4SZach Atkins PetscCall(VecRestoreArray(gvec, &garray));
161851a5d4e4SZach Atkins PetscCall(VecRestoreArrayRead(lvec, &larray));
16193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
162039d35262SVincent Le Chenadec }
162139d35262SVincent Le Chenadec
DMLocalToGlobalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)162266976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1623d71ae5a4SJacob Faibussowitsch {
162439d35262SVincent Le Chenadec PetscFunctionBegin;
162539d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
162639d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
162739d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
162939d35262SVincent Le Chenadec }
163039d35262SVincent Le Chenadec
DMLocalToLocalBegin_Composite(DM dm,Vec vec1,InsertMode mode,Vec vec2)163166976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2)
1632d71ae5a4SJacob Faibussowitsch {
163339d35262SVincent Le Chenadec struct DMCompositeLink *next;
163451a5d4e4SZach Atkins const PetscScalar *array1, *array1head;
163551a5d4e4SZach Atkins PetscScalar *array2, *array2head;
163639d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite *)dm->data;
163739d35262SVincent Le Chenadec
163839d35262SVincent Le Chenadec PetscFunctionBegin;
163939d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
164039d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2);
164139d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4);
164239d35262SVincent Le Chenadec
164348a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm));
164439d35262SVincent Le Chenadec
164551a5d4e4SZach Atkins PetscCall(VecGetArrayRead(vec1, &array1));
164651a5d4e4SZach Atkins array1head = array1;
16479566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec2, &array2));
164851a5d4e4SZach Atkins array2head = array2;
164939d35262SVincent Le Chenadec
165015229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */
165139d35262SVincent Le Chenadec next = com->next;
165239d35262SVincent Le Chenadec while (next) {
165339d35262SVincent Le Chenadec Vec local1, local2;
165439d35262SVincent Le Chenadec
16559566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local1));
165651a5d4e4SZach Atkins PetscCall(VecPlaceArray(local1, array1head));
16579566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local2));
165851a5d4e4SZach Atkins PetscCall(VecPlaceArray(local2, array2head));
16599566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2));
16609566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2));
16619566063dSJacob Faibussowitsch PetscCall(VecResetArray(local2));
16629566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local2));
16639566063dSJacob Faibussowitsch PetscCall(VecResetArray(local1));
16649566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local1));
166539d35262SVincent Le Chenadec
166651a5d4e4SZach Atkins array1head += next->nlocal;
166751a5d4e4SZach Atkins array2head += next->nlocal;
166839d35262SVincent Le Chenadec next = next->next;
166939d35262SVincent Le Chenadec }
167039d35262SVincent Le Chenadec
167151a5d4e4SZach Atkins PetscCall(VecRestoreArrayRead(vec1, &array1));
167251a5d4e4SZach Atkins PetscCall(VecRestoreArray(vec2, &array2));
16733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
167439d35262SVincent Le Chenadec }
167539d35262SVincent Le Chenadec
DMLocalToLocalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec)167666976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec)
1677d71ae5a4SJacob Faibussowitsch {
167839d35262SVincent Le Chenadec PetscFunctionBegin;
167939d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
168039d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2);
168139d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4);
16823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
16830c010503SBarry Smith }
168447c6ae99SBarry Smith
16856ae3a549SBarry Smith /*MC
1686dce8aebaSBarry Smith DMCOMPOSITE = "composite" - A `DM` object that is used to manage data for a collection of `DM`
16876ae3a549SBarry Smith
16886ae3a549SBarry Smith Level: intermediate
16896ae3a549SBarry Smith
1690db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()`
16916ae3a549SBarry Smith M*/
16926ae3a549SBarry Smith
DMCreate_Composite(DM p)1693d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p)
1694d71ae5a4SJacob Faibussowitsch {
1695a4121054SBarry Smith DM_Composite *com;
1696a4121054SBarry Smith
1697a4121054SBarry Smith PetscFunctionBegin;
16984dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&com));
1699a4121054SBarry Smith p->data = com;
1700a4121054SBarry Smith com->n = 0;
17017ac2b803SAlex Fikl com->nghost = 0;
17020298fd71SBarry Smith com->next = NULL;
1703a4121054SBarry Smith com->nDM = 0;
1704a4121054SBarry Smith
1705a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Composite;
1706a4121054SBarry Smith p->ops->createlocalvector = DMCreateLocalVector_Composite;
1707184d77edSJed Brown p->ops->getlocaltoglobalmapping = DMGetLocalToGlobalMapping_Composite;
17084d343eeaSMatthew G Knepley p->ops->createfieldis = DMCreateFieldIS_Composite;
170916621825SDmitry Karpeev p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite;
1710a4121054SBarry Smith p->ops->refine = DMRefine_Composite;
171114354c39SJed Brown p->ops->coarsen = DMCoarsen_Composite;
171225296bd5SBarry Smith p->ops->createinterpolation = DMCreateInterpolation_Composite;
171325296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Composite;
1714e727c939SJed Brown p->ops->getcoloring = DMCreateColoring_Composite;
1715a4121054SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite;
1716a4121054SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite;
171739d35262SVincent Le Chenadec p->ops->localtoglobalbegin = DMLocalToGlobalBegin_Composite;
171839d35262SVincent Le Chenadec p->ops->localtoglobalend = DMLocalToGlobalEnd_Composite;
171939d35262SVincent Le Chenadec p->ops->localtolocalbegin = DMLocalToLocalBegin_Composite;
172039d35262SVincent Le Chenadec p->ops->localtolocalend = DMLocalToLocalEnd_Composite;
1721a4121054SBarry Smith p->ops->destroy = DMDestroy_Composite;
1722a4121054SBarry Smith p->ops->view = DMView_Composite;
1723a4121054SBarry Smith p->ops->setup = DMSetUp_Composite;
1724e10fd815SStefano Zampini
17259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite));
17263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1727a4121054SBarry Smith }
1728a4121054SBarry Smith
17291fd49c25SBarry Smith /*@
1730dce8aebaSBarry Smith DMCompositeCreate - Creates a `DMCOMPOSITE`, used to generate "composite"
17310c010503SBarry Smith vectors made up of several subvectors.
17320c010503SBarry Smith
1733d083f849SBarry Smith Collective
173447c6ae99SBarry Smith
173547c6ae99SBarry Smith Input Parameter:
17360c010503SBarry Smith . comm - the processors that will share the global vector
17370c010503SBarry Smith
17382fe279fdSBarry Smith Output Parameter:
1739dce8aebaSBarry Smith . packer - the `DMCOMPOSITE` object
174047c6ae99SBarry Smith
174147c6ae99SBarry Smith Level: advanced
174247c6ae99SBarry Smith
174360225df5SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCreate()`
1744db781477SPatrick Sanan `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`
1745db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`
174647c6ae99SBarry Smith @*/
DMCompositeCreate(MPI_Comm comm,DM * packer)1747d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer)
1748d71ae5a4SJacob Faibussowitsch {
174947c6ae99SBarry Smith PetscFunctionBegin;
17504f572ea9SToby Isaac PetscAssertPointer(packer, 2);
17519566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, packer));
17529566063dSJacob Faibussowitsch PetscCall(DMSetType(*packer, DMCOMPOSITE));
17533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
175447c6ae99SBarry Smith }
1755