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 @*/ 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 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 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 9947c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/ 10066976f2fSJacob Faibussowitsch static PetscErrorCode DMSetUp_Composite(DM dm) 101d71ae5a4SJacob Faibussowitsch { 10247c6ae99SBarry Smith PetscInt nprev = 0; 10347c6ae99SBarry Smith PetscMPIInt rank, size; 10447c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 10547c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 10647c6ae99SBarry Smith PetscLayout map; 10747c6ae99SBarry Smith 10847c6ae99SBarry Smith PetscFunctionBegin; 10928b400f6SJacob Faibussowitsch PetscCheck(!com->setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Packer has already been setup"); 1109566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &map)); 1119566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(map, com->n)); 1129566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetSize(map, PETSC_DETERMINE)); 1139566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(map, 1)); 1149566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(map)); 1159566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &com->N)); 1169566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(map, &com->rstart, NULL)); 1179566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&map)); 11847c6ae99SBarry Smith 1199ae5db72SJed Brown /* now set the rstart for each linked vector */ 1209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 12247c6ae99SBarry Smith while (next) { 12347c6ae99SBarry Smith next->rstart = nprev; 12406ebdd98SJed Brown nprev += next->n; 12547c6ae99SBarry Smith next->grstart = com->rstart + next->rstart; 1269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &next->grstarts)); 1279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&next->grstart, 1, MPIU_INT, next->grstarts, 1, MPIU_INT, PetscObjectComm((PetscObject)dm))); 12847c6ae99SBarry Smith next = next->next; 12947c6ae99SBarry Smith } 13047c6ae99SBarry Smith com->setup = PETSC_TRUE; 1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13247c6ae99SBarry Smith } 13347c6ae99SBarry Smith 13473e31fe2SJed Brown /*@ 135aaa8cc7dSPierre Jolivet DMCompositeGetNumberDM - Gets the number of `DM` objects in the `DMCOMPOSITE` 13647c6ae99SBarry Smith representation. 13747c6ae99SBarry Smith 13847c6ae99SBarry Smith Not Collective 13947c6ae99SBarry Smith 14047c6ae99SBarry Smith Input Parameter: 141dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 14247c6ae99SBarry Smith 14347c6ae99SBarry Smith Output Parameter: 144dce8aebaSBarry Smith . nDM - the number of `DM` 14547c6ae99SBarry Smith 14647c6ae99SBarry Smith Level: beginner 14747c6ae99SBarry Smith 148dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM` 14947c6ae99SBarry Smith @*/ 150d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetNumberDM(DM dm, PetscInt *nDM) 151d71ae5a4SJacob Faibussowitsch { 15247c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 15371b14b3eSStefano Zampini PetscBool flg; 1545fd66863SKarl Rupp 15547c6ae99SBarry Smith PetscFunctionBegin; 15647c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1579566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1587a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 15947c6ae99SBarry Smith *nDM = com->nDM; 1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16147c6ae99SBarry Smith } 16247c6ae99SBarry Smith 16347c6ae99SBarry Smith /*@C 16447c6ae99SBarry Smith DMCompositeGetAccess - Allows one to access the individual packed vectors in their global 16547c6ae99SBarry Smith representation. 16647c6ae99SBarry Smith 16720f4b53cSBarry Smith Collective 16847c6ae99SBarry Smith 1699ae5db72SJed Brown Input Parameters: 170dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 1719ae5db72SJed Brown - gvec - the global vector 1729ae5db72SJed Brown 1732fe279fdSBarry Smith Output Parameter: 174a4e35b19SJacob Faibussowitsch . ... - the packed parallel vectors, `NULL` for those that are not needed 17547c6ae99SBarry Smith 17647c6ae99SBarry Smith Level: advanced 17747c6ae99SBarry Smith 178dce8aebaSBarry Smith Note: 179dce8aebaSBarry Smith Use `DMCompositeRestoreAccess()` to return the vectors when you no longer need them 180dce8aebaSBarry Smith 181dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetEntries()`, `DMCompositeScatter()` 18247c6ae99SBarry Smith @*/ 183d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccess(DM dm, Vec gvec, ...) 184d71ae5a4SJacob Faibussowitsch { 18547c6ae99SBarry Smith va_list Argp; 18647c6ae99SBarry Smith struct DMCompositeLink *next; 18747c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 1885edff71fSBarry Smith PetscInt readonly; 18971b14b3eSStefano Zampini PetscBool flg; 19047c6ae99SBarry Smith 19147c6ae99SBarry Smith PetscFunctionBegin; 19247c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 19347c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 1949566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1957a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 19647c6ae99SBarry Smith next = com->next; 19748a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 19847c6ae99SBarry Smith 1999566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec, &readonly)); 20015229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 20147c6ae99SBarry Smith va_start(Argp, gvec); 20247c6ae99SBarry Smith while (next) { 20347c6ae99SBarry Smith Vec *vec; 20447c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 2059ae5db72SJed Brown if (vec) { 2069566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, vec)); 2075edff71fSBarry Smith if (readonly) { 2085edff71fSBarry Smith const PetscScalar *array; 2099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 2109566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(*vec, array + next->rstart)); 2119566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(*vec)); 2129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 2135edff71fSBarry Smith } else { 2145edff71fSBarry Smith PetscScalar *array; 2159566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 2169566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(*vec, array + next->rstart)); 2179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 21847c6ae99SBarry Smith } 2195edff71fSBarry Smith } 22047c6ae99SBarry Smith next = next->next; 22147c6ae99SBarry Smith } 22247c6ae99SBarry Smith va_end(Argp); 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22447c6ae99SBarry Smith } 22547c6ae99SBarry Smith 2265d83a8b1SBarry Smith /*@ 227f73e5cebSJed Brown DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global 228f73e5cebSJed Brown representation. 229f73e5cebSJed Brown 23020f4b53cSBarry Smith Collective 231f73e5cebSJed Brown 232f73e5cebSJed Brown Input Parameters: 233dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` 234f73e5cebSJed Brown . pvec - packed vector 235f73e5cebSJed Brown . nwanted - number of vectors wanted 236ce78bad3SBarry Smith - wanted - sorted array of integers indicating thde vectors wanted, or `NULL` to get all vectors, length `nwanted` 237f73e5cebSJed Brown 2382fe279fdSBarry Smith Output Parameter: 239f13dfd9eSBarry Smith . vecs - array of requested global vectors (must be previously allocated and of length `nwanted`) 240f73e5cebSJed Brown 241f73e5cebSJed Brown Level: advanced 242f73e5cebSJed Brown 243dce8aebaSBarry Smith Note: 244dce8aebaSBarry Smith Use `DMCompositeRestoreAccessArray()` to return the vectors when you no longer need them 245dce8aebaSBarry Smith 246dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 247f73e5cebSJed Brown @*/ 2485d83a8b1SBarry Smith PetscErrorCode DMCompositeGetAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[]) 249d71ae5a4SJacob Faibussowitsch { 250f73e5cebSJed Brown struct DMCompositeLink *link; 251f73e5cebSJed Brown PetscInt i, wnum; 252f73e5cebSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 253bee642f7SBarry Smith PetscInt readonly; 25471b14b3eSStefano Zampini PetscBool flg; 255f73e5cebSJed Brown 256f73e5cebSJed Brown PetscFunctionBegin; 257f73e5cebSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 258f73e5cebSJed Brown PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 2599566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 2607a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 26148a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 262f73e5cebSJed Brown 2639566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 264f73e5cebSJed Brown for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 265f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 266f73e5cebSJed Brown Vec v; 2679566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(link->dm, &v)); 268bee642f7SBarry Smith if (readonly) { 269bee642f7SBarry Smith const PetscScalar *array; 2709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvec, &array)); 2719566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + link->rstart)); 2729566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(v)); 2739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvec, &array)); 274bee642f7SBarry Smith } else { 275bee642f7SBarry Smith PetscScalar *array; 2769566063dSJacob Faibussowitsch PetscCall(VecGetArray(pvec, &array)); 2779566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + link->rstart)); 2789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pvec, &array)); 279bee642f7SBarry Smith } 280f73e5cebSJed Brown vecs[wnum++] = v; 281f73e5cebSJed Brown } 282f73e5cebSJed Brown } 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 284f73e5cebSJed Brown } 285f73e5cebSJed Brown 2865d83a8b1SBarry Smith /*@ 2877ac2b803SAlex Fikl DMCompositeGetLocalAccessArray - Allows one to access the individual 2887ac2b803SAlex Fikl packed vectors in their local representation. 2897ac2b803SAlex Fikl 29020f4b53cSBarry Smith Collective 2917ac2b803SAlex Fikl 2927ac2b803SAlex Fikl Input Parameters: 293dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` 2947ac2b803SAlex Fikl . pvec - packed vector 2957ac2b803SAlex Fikl . nwanted - number of vectors wanted 296f13dfd9eSBarry Smith - wanted - sorted array of vectors wanted, or `NULL` to get all vectors, length `nwanted` 2977ac2b803SAlex Fikl 2982fe279fdSBarry Smith Output Parameter: 299f13dfd9eSBarry Smith . vecs - array of requested local vectors (must be allocated and of length `nwanted`) 3007ac2b803SAlex Fikl 3017ac2b803SAlex Fikl Level: advanced 3027ac2b803SAlex Fikl 303dce8aebaSBarry Smith Note: 304dce8aebaSBarry Smith Use `DMCompositeRestoreLocalAccessArray()` to return the vectors 305dce8aebaSBarry Smith when you no longer need them. 306dce8aebaSBarry Smith 307dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`, 308db781477SPatrick Sanan `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 3097ac2b803SAlex Fikl @*/ 3105d83a8b1SBarry Smith PetscErrorCode DMCompositeGetLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[]) 311d71ae5a4SJacob Faibussowitsch { 3127ac2b803SAlex Fikl struct DMCompositeLink *link; 3137ac2b803SAlex Fikl PetscInt i, wnum; 3147ac2b803SAlex Fikl DM_Composite *com = (DM_Composite *)dm->data; 3157ac2b803SAlex Fikl PetscInt readonly; 3167ac2b803SAlex Fikl PetscInt nlocal = 0; 31771b14b3eSStefano Zampini PetscBool flg; 3187ac2b803SAlex Fikl 3197ac2b803SAlex Fikl PetscFunctionBegin; 3207ac2b803SAlex Fikl PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3217ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 3229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 3237a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 32448a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 3257ac2b803SAlex Fikl 3269566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 3277ac2b803SAlex Fikl for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 3287ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 3297ac2b803SAlex Fikl Vec v; 3309566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(link->dm, &v)); 3317ac2b803SAlex Fikl if (readonly) { 3327ac2b803SAlex Fikl const PetscScalar *array; 3339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvec, &array)); 3349566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + nlocal)); 335b1c3483dSMark 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 3369566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(v)); 3379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvec, &array)); 3387ac2b803SAlex Fikl } else { 3397ac2b803SAlex Fikl PetscScalar *array; 3409566063dSJacob Faibussowitsch PetscCall(VecGetArray(pvec, &array)); 3419566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + nlocal)); 3429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pvec, &array)); 3437ac2b803SAlex Fikl } 3447ac2b803SAlex Fikl vecs[wnum++] = v; 3457ac2b803SAlex Fikl } 3467ac2b803SAlex Fikl 3477ac2b803SAlex Fikl nlocal += link->nlocal; 3487ac2b803SAlex Fikl } 3493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3507ac2b803SAlex Fikl } 3517ac2b803SAlex Fikl 35247c6ae99SBarry Smith /*@C 353dce8aebaSBarry Smith DMCompositeRestoreAccess - Returns the vectors obtained with `DMCompositeGetAccess()` 35447c6ae99SBarry Smith representation. 35547c6ae99SBarry Smith 35620f4b53cSBarry Smith Collective 35747c6ae99SBarry Smith 3589ae5db72SJed Brown Input Parameters: 359dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 36047c6ae99SBarry Smith . gvec - the global vector 361a4e35b19SJacob Faibussowitsch - ... - the individual parallel vectors, `NULL` for those that are not needed 36247c6ae99SBarry Smith 36347c6ae99SBarry Smith Level: advanced 36447c6ae99SBarry Smith 365dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 366db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`, 36742747ad1SJacob Faibussowitsch `DMCompositeGetAccess()` 36847c6ae99SBarry Smith @*/ 369d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...) 370d71ae5a4SJacob Faibussowitsch { 37147c6ae99SBarry Smith va_list Argp; 37247c6ae99SBarry Smith struct DMCompositeLink *next; 37347c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 3745edff71fSBarry Smith PetscInt readonly; 37571b14b3eSStefano Zampini PetscBool flg; 37647c6ae99SBarry Smith 37747c6ae99SBarry Smith PetscFunctionBegin; 37847c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 37947c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 3809566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 3817a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 38247c6ae99SBarry Smith next = com->next; 38348a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 38447c6ae99SBarry Smith 3859566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec, &readonly)); 38615229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 38747c6ae99SBarry Smith va_start(Argp, gvec); 38847c6ae99SBarry Smith while (next) { 38947c6ae99SBarry Smith Vec *vec; 39047c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 3919ae5db72SJed Brown if (vec) { 3929566063dSJacob Faibussowitsch PetscCall(VecResetArray(*vec)); 3931baa6e33SBarry Smith if (readonly) PetscCall(VecLockReadPop(*vec)); 3949566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, vec)); 39547c6ae99SBarry Smith } 39647c6ae99SBarry Smith next = next->next; 39747c6ae99SBarry Smith } 39847c6ae99SBarry Smith va_end(Argp); 3993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40047c6ae99SBarry Smith } 40147c6ae99SBarry Smith 4025d83a8b1SBarry Smith /*@ 403dce8aebaSBarry Smith DMCompositeRestoreAccessArray - Returns the vectors obtained with `DMCompositeGetAccessArray()` 404f73e5cebSJed Brown 40520f4b53cSBarry Smith Collective 406f73e5cebSJed Brown 407f73e5cebSJed Brown Input Parameters: 408dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 409f73e5cebSJed Brown . pvec - packed vector 410f73e5cebSJed Brown . nwanted - number of vectors wanted 411f13dfd9eSBarry Smith . wanted - sorted array of vectors wanted, or `NULL` to restore all vectors 412f13dfd9eSBarry Smith - vecs - array of global vectors 413f73e5cebSJed Brown 414f73e5cebSJed Brown Level: advanced 415f73e5cebSJed Brown 416dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 417f73e5cebSJed Brown @*/ 4185d83a8b1SBarry Smith PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec vecs[]) 419d71ae5a4SJacob Faibussowitsch { 420f73e5cebSJed Brown struct DMCompositeLink *link; 421f73e5cebSJed Brown PetscInt i, wnum; 422f73e5cebSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 423bee642f7SBarry Smith PetscInt readonly; 42471b14b3eSStefano Zampini PetscBool flg; 425f73e5cebSJed Brown 426f73e5cebSJed Brown PetscFunctionBegin; 427f73e5cebSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 428f73e5cebSJed Brown PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 4299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 4307a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 43148a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 432f73e5cebSJed Brown 4339566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 434f73e5cebSJed Brown for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 435f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 4369566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum])); 43748a46eb9SPierre Jolivet if (readonly) PetscCall(VecLockReadPop(vecs[wnum])); 4389566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum])); 439f73e5cebSJed Brown wnum++; 440f73e5cebSJed Brown } 441f73e5cebSJed Brown } 4423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 443f73e5cebSJed Brown } 444f73e5cebSJed Brown 4455d83a8b1SBarry Smith /*@ 446dce8aebaSBarry Smith DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with `DMCompositeGetLocalAccessArray()`. 4477ac2b803SAlex Fikl 44820f4b53cSBarry Smith Collective 4497ac2b803SAlex Fikl 4507ac2b803SAlex Fikl Input Parameters: 451dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 4527ac2b803SAlex Fikl . pvec - packed vector 4537ac2b803SAlex Fikl . nwanted - number of vectors wanted 454f13dfd9eSBarry Smith . wanted - sorted array of vectors wanted, or `NULL` to restore all vectors 455f13dfd9eSBarry Smith - vecs - array of local vectors 4567ac2b803SAlex Fikl 4577ac2b803SAlex Fikl Level: advanced 4587ac2b803SAlex Fikl 459dce8aebaSBarry Smith Note: 460f13dfd9eSBarry Smith `nwanted` and `wanted` must match the values given to `DMCompositeGetLocalAccessArray()` 4617ac2b803SAlex Fikl otherwise the call will fail. 4627ac2b803SAlex Fikl 463dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`, 464db781477SPatrick Sanan `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, 465db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeGather()` 4667ac2b803SAlex Fikl @*/ 4675d83a8b1SBarry Smith PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt wanted[], Vec *vecs) 468d71ae5a4SJacob Faibussowitsch { 4697ac2b803SAlex Fikl struct DMCompositeLink *link; 4707ac2b803SAlex Fikl PetscInt i, wnum; 4717ac2b803SAlex Fikl DM_Composite *com = (DM_Composite *)dm->data; 4727ac2b803SAlex Fikl PetscInt readonly; 47371b14b3eSStefano Zampini PetscBool flg; 4747ac2b803SAlex Fikl 4757ac2b803SAlex Fikl PetscFunctionBegin; 4767ac2b803SAlex Fikl PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4777ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 4789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 4797a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 48048a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 4817ac2b803SAlex Fikl 4829566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 4837ac2b803SAlex Fikl for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 4847ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 4859566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum])); 48648a46eb9SPierre Jolivet if (readonly) PetscCall(VecLockReadPop(vecs[wnum])); 4879566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum])); 4887ac2b803SAlex Fikl wnum++; 4897ac2b803SAlex Fikl } 4907ac2b803SAlex Fikl } 4913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4927ac2b803SAlex Fikl } 4937ac2b803SAlex Fikl 49447c6ae99SBarry Smith /*@C 49547c6ae99SBarry Smith DMCompositeScatter - Scatters from a global packed vector into its individual local vectors 49647c6ae99SBarry Smith 49720f4b53cSBarry Smith Collective 49847c6ae99SBarry Smith 4999ae5db72SJed Brown Input Parameters: 500dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 50147c6ae99SBarry Smith . gvec - the global vector 502a4e35b19SJacob Faibussowitsch - ... - the individual sequential vectors, `NULL` for those that are not needed 50347c6ae99SBarry Smith 50447c6ae99SBarry Smith Level: advanced 50547c6ae99SBarry Smith 506dce8aebaSBarry Smith Note: 507dce8aebaSBarry Smith `DMCompositeScatterArray()` is a non-variadic alternative that is often more convenient for library callers and is 5086f3c3dcfSJed Brown accessible from Fortran. 5096f3c3dcfSJed Brown 510dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 511db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 512db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 513db781477SPatrick Sanan `DMCompositeScatterArray()` 51447c6ae99SBarry Smith @*/ 515d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...) 516d71ae5a4SJacob Faibussowitsch { 51747c6ae99SBarry Smith va_list Argp; 51847c6ae99SBarry Smith struct DMCompositeLink *next; 519c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt; 52047c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 52171b14b3eSStefano Zampini PetscBool flg; 52247c6ae99SBarry Smith 52347c6ae99SBarry Smith PetscFunctionBegin; 52447c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 52547c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 5269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 5277a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 52848a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 52947c6ae99SBarry Smith 53015229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 53147c6ae99SBarry Smith va_start(Argp, gvec); 5328fd8f222SJed Brown for (cnt = 3, next = com->next; next; cnt++, next = next->next) { 5339ae5db72SJed Brown Vec local; 5349ae5db72SJed Brown local = va_arg(Argp, Vec); 5359ae5db72SJed Brown if (local) { 5369ae5db72SJed Brown Vec global; 5375edff71fSBarry Smith const PetscScalar *array; 53863a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt)); 5399566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 5409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 5419566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 5429566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local)); 5439566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local)); 5449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 5459566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 5469566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 54747c6ae99SBarry Smith } 54847c6ae99SBarry Smith } 54947c6ae99SBarry Smith va_end(Argp); 5503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55147c6ae99SBarry Smith } 55247c6ae99SBarry Smith 5536f3c3dcfSJed Brown /*@ 5546f3c3dcfSJed Brown DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors 5556f3c3dcfSJed Brown 55620f4b53cSBarry Smith Collective 5576f3c3dcfSJed Brown 5586f3c3dcfSJed Brown Input Parameters: 559dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 5606f3c3dcfSJed Brown . gvec - the global vector 561a2b725a8SWilliam Gropp - lvecs - array of local vectors, NULL for any that are not needed 5626f3c3dcfSJed Brown 5636f3c3dcfSJed Brown Level: advanced 5646f3c3dcfSJed Brown 5656f3c3dcfSJed Brown Note: 566dce8aebaSBarry Smith This is a non-variadic alternative to `DMCompositeScatter()` 5676f3c3dcfSJed Brown 568dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()` 569db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 570db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 5716f3c3dcfSJed Brown @*/ 572d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs) 573d71ae5a4SJacob Faibussowitsch { 5746f3c3dcfSJed Brown struct DMCompositeLink *next; 5756f3c3dcfSJed Brown PetscInt i; 5766f3c3dcfSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 57771b14b3eSStefano Zampini PetscBool flg; 5786f3c3dcfSJed Brown 5796f3c3dcfSJed Brown PetscFunctionBegin; 5806f3c3dcfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5816f3c3dcfSJed Brown PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 5829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 5837a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 58448a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 5856f3c3dcfSJed Brown 58615229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 5876f3c3dcfSJed Brown for (i = 0, next = com->next; next; next = next->next, i++) { 5886f3c3dcfSJed Brown if (lvecs[i]) { 5896f3c3dcfSJed Brown Vec global; 590c5d31e75SLisandro Dalcin const PetscScalar *array; 5916f3c3dcfSJed Brown PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3); 5929566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 5939566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 5949566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart)); 5959566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i])); 5969566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i])); 5979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 5989566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 5999566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 6006f3c3dcfSJed Brown } 6016f3c3dcfSJed Brown } 6023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6036f3c3dcfSJed Brown } 6046f3c3dcfSJed Brown 60547c6ae99SBarry Smith /*@C 60647c6ae99SBarry Smith DMCompositeGather - Gathers into a global packed vector from its individual local vectors 60747c6ae99SBarry Smith 60820f4b53cSBarry Smith Collective 60947c6ae99SBarry Smith 610d8d19677SJose E. Roman Input Parameters: 611dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 612dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES` 613a4e35b19SJacob Faibussowitsch . gvec - the global vector 614a4e35b19SJacob Faibussowitsch - ... - the individual sequential vectors, `NULL` for any that are not needed 61547c6ae99SBarry Smith 61647c6ae99SBarry Smith Level: advanced 61747c6ae99SBarry Smith 61860225df5SJacob Faibussowitsch Fortran Notes: 619dce8aebaSBarry Smith Fortran users should use `DMCompositeGatherArray()` 620f5f57ec0SBarry Smith 621dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 622db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 623db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 62447c6ae99SBarry Smith @*/ 625d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...) 626d71ae5a4SJacob Faibussowitsch { 62747c6ae99SBarry Smith va_list Argp; 62847c6ae99SBarry Smith struct DMCompositeLink *next; 62947c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 630c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt; 63171b14b3eSStefano Zampini PetscBool flg; 63247c6ae99SBarry Smith 63347c6ae99SBarry Smith PetscFunctionBegin; 63447c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 635064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3); 6369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 6377a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 63848a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 63947c6ae99SBarry Smith 64015229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 6411dac896bSSatish Balay va_start(Argp, gvec); 6428fd8f222SJed Brown for (cnt = 3, next = com->next; next; cnt++, next = next->next) { 6439ae5db72SJed Brown Vec local; 6449ae5db72SJed Brown local = va_arg(Argp, Vec); 6459ae5db72SJed Brown if (local) { 64647c6ae99SBarry Smith PetscScalar *array; 6479ae5db72SJed Brown Vec global; 64863a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt)); 6499566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 6509566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 6519566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 6529566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global)); 6539566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global)); 6549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 6559566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 6569566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 65747c6ae99SBarry Smith } 65847c6ae99SBarry Smith } 65947c6ae99SBarry Smith va_end(Argp); 6603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66147c6ae99SBarry Smith } 66247c6ae99SBarry Smith 6636f3c3dcfSJed Brown /*@ 6646f3c3dcfSJed Brown DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors 6656f3c3dcfSJed Brown 66620f4b53cSBarry Smith Collective 6676f3c3dcfSJed Brown 668d8d19677SJose E. Roman Input Parameters: 669dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 6706f3c3dcfSJed Brown . gvec - the global vector 671dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES` 6726f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed 6736f3c3dcfSJed Brown 6746f3c3dcfSJed Brown Level: advanced 6756f3c3dcfSJed Brown 676dce8aebaSBarry Smith Note: 677dce8aebaSBarry Smith This is a non-variadic alternative to `DMCompositeGather()`. 6786f3c3dcfSJed Brown 679dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 680db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 681db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`, 6826f3c3dcfSJed Brown @*/ 683d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs) 684d71ae5a4SJacob Faibussowitsch { 6856f3c3dcfSJed Brown struct DMCompositeLink *next; 6866f3c3dcfSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 6876f3c3dcfSJed Brown PetscInt i; 68871b14b3eSStefano Zampini PetscBool flg; 6896f3c3dcfSJed Brown 6906f3c3dcfSJed Brown PetscFunctionBegin; 6916f3c3dcfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 692064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3); 6939566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 6947a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 69548a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 6966f3c3dcfSJed Brown 69715229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 6986f3c3dcfSJed Brown for (next = com->next, i = 0; next; next = next->next, i++) { 6996f3c3dcfSJed Brown if (lvecs[i]) { 7006f3c3dcfSJed Brown PetscScalar *array; 7016f3c3dcfSJed Brown Vec global; 702064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4); 7039566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 7049566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 7059566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 7069566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global)); 7079566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global)); 7089566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 7099566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 7109566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 7116f3c3dcfSJed Brown } 7126f3c3dcfSJed Brown } 7133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7146f3c3dcfSJed Brown } 7156f3c3dcfSJed Brown 716f5f57ec0SBarry Smith /*@ 717dce8aebaSBarry Smith DMCompositeAddDM - adds a `DM` vector to a `DMCOMPOSITE` 71847c6ae99SBarry Smith 71920f4b53cSBarry Smith Collective 72047c6ae99SBarry Smith 721d8d19677SJose E. Roman Input Parameters: 722dce8aebaSBarry Smith + dmc - the `DMCOMPOSITE` object 723dce8aebaSBarry Smith - dm - the `DM` object 72447c6ae99SBarry Smith 72547c6ae99SBarry Smith Level: advanced 72647c6ae99SBarry Smith 72742747ad1SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeGather()`, `DMCreateGlobalVector()`, 728db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 729db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 73047c6ae99SBarry Smith @*/ 731d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeAddDM(DM dmc, DM dm) 732d71ae5a4SJacob Faibussowitsch { 73306ebdd98SJed Brown PetscInt n, nlocal; 73447c6ae99SBarry Smith struct DMCompositeLink *mine, *next; 73506ebdd98SJed Brown Vec global, local; 73647c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dmc->data; 73771b14b3eSStefano Zampini PetscBool flg; 73847c6ae99SBarry Smith 73947c6ae99SBarry Smith PetscFunctionBegin; 74047c6ae99SBarry Smith PetscValidHeaderSpecific(dmc, DM_CLASSID, 1); 74147c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 7429566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg)); 7437a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 74447c6ae99SBarry Smith next = com->next; 74528b400f6SJacob Faibussowitsch PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite"); 74647c6ae99SBarry Smith 74747c6ae99SBarry Smith /* create new link */ 7489566063dSJacob Faibussowitsch PetscCall(PetscNew(&mine)); 7499566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 750cac35740SJunchao Zhang PetscCall(DMCreateGlobalVector(dm, &global)); // Not using DMGetGlobalVector() since it will stash the vector with a type decided by dm, 751cac35740SJunchao Zhang PetscCall(VecGetLocalSize(global, &n)); // while we may want to set dmc's vectype later, say via DMSetFromOptions(dmc), and we 752cac35740SJunchao Zhang PetscCall(VecDestroy(&global)); // want to propagate the type to dm. 753cac35740SJunchao Zhang PetscCall(DMCreateLocalVector(dm, &local)); // Not using DMGetLocalVector(), same reason as above. 7549566063dSJacob Faibussowitsch PetscCall(VecGetSize(local, &nlocal)); 755cac35740SJunchao Zhang PetscCall(VecDestroy(&local)); 7568865f1eaSKarl Rupp 75747c6ae99SBarry Smith mine->n = n; 75806ebdd98SJed Brown mine->nlocal = nlocal; 75947c6ae99SBarry Smith mine->dm = dm; 7600298fd71SBarry Smith mine->next = NULL; 76147c6ae99SBarry Smith com->n += n; 7627ac2b803SAlex Fikl com->nghost += nlocal; 76347c6ae99SBarry Smith 76447c6ae99SBarry Smith /* add to end of list */ 7658865f1eaSKarl Rupp if (!next) com->next = mine; 7668865f1eaSKarl Rupp else { 76747c6ae99SBarry Smith while (next->next) next = next->next; 76847c6ae99SBarry Smith next->next = mine; 76947c6ae99SBarry Smith } 77047c6ae99SBarry Smith com->nDM++; 77147c6ae99SBarry Smith com->nmine++; 7723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 77347c6ae99SBarry Smith } 77447c6ae99SBarry Smith 7759804daf3SBarry Smith #include <petscdraw.h> 776d6acfc2dSPierre Jolivet PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecView_MPI(Vec, PetscViewer); 777d6acfc2dSPierre Jolivet 77866976f2fSJacob Faibussowitsch static PetscErrorCode VecView_DMComposite(Vec gvec, PetscViewer viewer) 779d71ae5a4SJacob Faibussowitsch { 78047c6ae99SBarry Smith DM dm; 78147c6ae99SBarry Smith struct DMCompositeLink *next; 78247c6ae99SBarry Smith PetscBool isdraw; 783cef07954SSatish Balay DM_Composite *com; 78447c6ae99SBarry Smith 78547c6ae99SBarry Smith PetscFunctionBegin; 7869566063dSJacob Faibussowitsch PetscCall(VecGetDM(gvec, &dm)); 7877a8be351SBarry Smith PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite"); 78847c6ae99SBarry Smith com = (DM_Composite *)dm->data; 78947c6ae99SBarry Smith next = com->next; 79047c6ae99SBarry Smith 7919566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 79247c6ae99SBarry Smith if (!isdraw) { 79347c6ae99SBarry Smith /* do I really want to call this? */ 7949566063dSJacob Faibussowitsch PetscCall(VecView_MPI(gvec, viewer)); 79547c6ae99SBarry Smith } else { 79647c6ae99SBarry Smith PetscInt cnt = 0; 79747c6ae99SBarry Smith 79815229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 79947c6ae99SBarry Smith while (next) { 80047c6ae99SBarry Smith Vec vec; 801ca4278abSLisandro Dalcin const PetscScalar *array; 80247c6ae99SBarry Smith PetscInt bs; 80347c6ae99SBarry Smith 8049ae5db72SJed Brown /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */ 8059566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &vec)); 8069566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 8079566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart)); 8089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 8099566063dSJacob Faibussowitsch PetscCall(VecView(vec, viewer)); 8109566063dSJacob Faibussowitsch PetscCall(VecResetArray(vec)); 8119566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vec, &bs)); 8129566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &vec)); 8139566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer, bs)); 81447c6ae99SBarry Smith cnt += bs; 81547c6ae99SBarry Smith next = next->next; 81647c6ae99SBarry Smith } 8179566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt)); 81847c6ae99SBarry Smith } 8193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82047c6ae99SBarry Smith } 82147c6ae99SBarry Smith 82266976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec) 823d71ae5a4SJacob Faibussowitsch { 82447c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 82547c6ae99SBarry Smith 82647c6ae99SBarry Smith PetscFunctionBegin; 82747c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8289566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 8299566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 8309566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec)); 8319566063dSJacob Faibussowitsch PetscCall(VecSetType(*gvec, dm->vectype)); 8329566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*gvec, com->n, com->N)); 8339566063dSJacob Faibussowitsch PetscCall(VecSetDM(*gvec, dm)); 83457d50842SBarry Smith PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (PetscErrorCodeFn *)VecView_DMComposite)); 8353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 83647c6ae99SBarry Smith } 83747c6ae99SBarry Smith 83866976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec) 839d71ae5a4SJacob Faibussowitsch { 84047c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 84147c6ae99SBarry Smith 84247c6ae99SBarry Smith PetscFunctionBegin; 84347c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 84447c6ae99SBarry Smith if (!com->setup) { 8459566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 8469566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 84747c6ae99SBarry Smith } 8489566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, lvec)); 8499566063dSJacob Faibussowitsch PetscCall(VecSetType(*lvec, dm->vectype)); 8509566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE)); 8519566063dSJacob Faibussowitsch PetscCall(VecSetDM(*lvec, dm)); 8523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 85347c6ae99SBarry Smith } 85447c6ae99SBarry Smith 85547c6ae99SBarry Smith /*@C 856dce8aebaSBarry Smith DMCompositeGetISLocalToGlobalMappings - gets an `ISLocalToGlobalMapping` for each `DM` in the `DMCOMPOSITE`, maps to the composite global space 85747c6ae99SBarry Smith 85820f4b53cSBarry Smith Collective; No Fortran Support 85947c6ae99SBarry Smith 86047c6ae99SBarry Smith Input Parameter: 861dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 86247c6ae99SBarry Smith 8632fe279fdSBarry Smith Output Parameter: 8649ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes 865dce8aebaSBarry Smith all the ghost points that individual ghosted `DMDA` may have. 86647c6ae99SBarry Smith 86747c6ae99SBarry Smith Level: advanced 86847c6ae99SBarry Smith 869dce8aebaSBarry Smith Note: 870f13dfd9eSBarry Smith Each entry of `ltogs` should be destroyed with `ISLocalToGlobalMappingDestroy()`, `ltogs` should be freed with `PetscFree()`. 87147c6ae99SBarry Smith 872dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 873db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 874c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 87547c6ae99SBarry Smith @*/ 876cc4c1da9SBarry Smith PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping *ltogs[]) 877d71ae5a4SJacob Faibussowitsch { 87847c6ae99SBarry Smith PetscInt i, *idx, n, cnt; 87947c6ae99SBarry Smith struct DMCompositeLink *next; 88047c6ae99SBarry Smith PetscMPIInt rank; 88147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 88271b14b3eSStefano Zampini PetscBool flg; 88347c6ae99SBarry Smith 88447c6ae99SBarry Smith PetscFunctionBegin; 88547c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8869566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 8877a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 8889566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 8899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM, ltogs)); 89047c6ae99SBarry Smith next = com->next; 8919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 89247c6ae99SBarry Smith 89315229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 89447c6ae99SBarry Smith cnt = 0; 89547c6ae99SBarry Smith while (next) { 8966eb61c8cSJed Brown ISLocalToGlobalMapping ltog; 8976eb61c8cSJed Brown PetscMPIInt size; 89886994e45SJed Brown const PetscInt *suboff, *indices; 8996eb61c8cSJed Brown Vec global; 90047c6ae99SBarry Smith 9016eb61c8cSJed Brown /* Get sub-DM global indices for each local dof */ 9029566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(next->dm, <og)); 9039566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n)); 9049566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices)); 9059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idx)); 90647c6ae99SBarry Smith 9076eb61c8cSJed Brown /* Get the offsets for the sub-DM global vector */ 9089566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 9099566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(global, &suboff)); 9109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size)); 9116eb61c8cSJed Brown 9126eb61c8cSJed Brown /* Shift the sub-DM definition of the global space to the composite global space */ 9136eb61c8cSJed Brown for (i = 0; i < n; i++) { 91486994e45SJed Brown PetscInt subi = indices[i], lo = 0, hi = size, t; 915e333b035SStefano Zampini /* There's no consensus on what a negative index means, 916e333b035SStefano Zampini except for skipping when setting the values in vectors and matrices */ 9179371c9d4SSatish Balay if (subi < 0) { 9189371c9d4SSatish Balay idx[i] = subi - next->grstarts[rank]; 9199371c9d4SSatish Balay continue; 9209371c9d4SSatish Balay } 9216eb61c8cSJed Brown /* Binary search to find which rank owns subi */ 9226eb61c8cSJed Brown while (hi - lo > 1) { 9236eb61c8cSJed Brown t = lo + (hi - lo) / 2; 9246eb61c8cSJed Brown if (suboff[t] > subi) hi = t; 9256eb61c8cSJed Brown else lo = t; 9266eb61c8cSJed Brown } 9276eb61c8cSJed Brown idx[i] = subi - suboff[lo] + next->grstarts[lo]; 9286eb61c8cSJed Brown } 9299566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices)); 9309566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt])); 9319566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 93247c6ae99SBarry Smith next = next->next; 93347c6ae99SBarry Smith cnt++; 93447c6ae99SBarry Smith } 9353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 93647c6ae99SBarry Smith } 93747c6ae99SBarry Smith 93887c85e80SJed Brown /*@C 9399ae5db72SJed Brown DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector 94087c85e80SJed Brown 94120f4b53cSBarry Smith Not Collective; No Fortran Support 94287c85e80SJed Brown 9434165533cSJose E. Roman Input Parameter: 944dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` 94587c85e80SJed Brown 9464165533cSJose E. Roman Output Parameter: 94715229ffcSPierre Jolivet . is - array of serial index sets for each component of the `DMCOMPOSITE` 94887c85e80SJed Brown 94987c85e80SJed Brown Level: intermediate 95087c85e80SJed Brown 95187c85e80SJed Brown Notes: 95287c85e80SJed Brown At present, a composite local vector does not normally exist. This function is used to provide index sets for 95315229ffcSPierre Jolivet `MatGetLocalSubMatrix()`. In the future, the scatters for each entry in the `DMCOMPOSITE` may be merged into a single 9549ae5db72SJed Brown scatter to a composite local vector. The user should not typically need to know which is being done. 95587c85e80SJed Brown 956dce8aebaSBarry Smith To get the composite global indices at all local points (including ghosts), use `DMCompositeGetISLocalToGlobalMappings()`. 95787c85e80SJed Brown 958dce8aebaSBarry Smith To get index sets for pieces of the composite global vector, use `DMCompositeGetGlobalISs()`. 95987c85e80SJed Brown 960dce8aebaSBarry Smith Each returned `IS` should be destroyed with `ISDestroy()`, the array should be freed with `PetscFree()`. 96187c85e80SJed Brown 962f13dfd9eSBarry Smith Fortran Note: 963ce78bad3SBarry Smith Use `DMCompositeRestoreLocalISs()` to release the `is`. 964f13dfd9eSBarry Smith 965f13dfd9eSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`, 966f13dfd9eSBarry Smith `MatCreateLocalRef()`, `DMCompositeGetNumberDM()` 96787c85e80SJed Brown @*/ 968cc4c1da9SBarry Smith PetscErrorCode DMCompositeGetLocalISs(DM dm, IS *is[]) 969d71ae5a4SJacob Faibussowitsch { 97087c85e80SJed Brown DM_Composite *com = (DM_Composite *)dm->data; 97187c85e80SJed Brown struct DMCompositeLink *link; 97287c85e80SJed Brown PetscInt cnt, start; 97371b14b3eSStefano Zampini PetscBool flg; 97487c85e80SJed Brown 97587c85e80SJed Brown PetscFunctionBegin; 97687c85e80SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 9774f572ea9SToby Isaac PetscAssertPointer(is, 2); 9789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 9797a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 9809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nmine, is)); 98106ebdd98SJed Brown for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) { 982520db06cSJed Brown PetscInt bs; 9839566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt])); 9849566063dSJacob Faibussowitsch PetscCall(DMGetBlockSize(link->dm, &bs)); 9859566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize((*is)[cnt], bs)); 986520db06cSJed Brown } 9873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 98887c85e80SJed Brown } 98987c85e80SJed Brown 99047c6ae99SBarry Smith /*@C 991dce8aebaSBarry Smith DMCompositeGetGlobalISs - Gets the index sets for each composed object in a `DMCOMPOSITE` 99247c6ae99SBarry Smith 99320f4b53cSBarry Smith Collective 99447c6ae99SBarry Smith 99547c6ae99SBarry Smith Input Parameter: 996dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 99747c6ae99SBarry Smith 9982fe279fdSBarry Smith Output Parameter: 99947c6ae99SBarry Smith . is - the array of index sets 100047c6ae99SBarry Smith 100147c6ae99SBarry Smith Level: advanced 100247c6ae99SBarry Smith 100347c6ae99SBarry Smith Notes: 1004f13dfd9eSBarry Smith The `is` entries should be destroyed with `ISDestroy()`, `is` should be freed with `PetscFree()` 100547c6ae99SBarry Smith 100647c6ae99SBarry Smith These could be used to extract a subset of vector entries for a "multi-physics" preconditioner 100747c6ae99SBarry Smith 1008dce8aebaSBarry Smith Use `DMCompositeGetLocalISs()` for index sets in the packed local numbering, and 1009dce8aebaSBarry Smith `DMCompositeGetISLocalToGlobalMappings()` for to map local sub-`DM` (including ghost) indices to packed global 10106eb61c8cSJed Brown indices. 101147c6ae99SBarry Smith 1012ce78bad3SBarry Smith Fortran Note: 1013ce78bad3SBarry Smith Use `DMCompositeRestoreGlobalISs()` to release the `is`. 1014f3cb0f7eSJed Brown 1015dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1016db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 1017c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 101847c6ae99SBarry Smith @*/ 1019d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[]) 1020d71ae5a4SJacob Faibussowitsch { 102166bb578eSMark Adams PetscInt cnt = 0; 102247c6ae99SBarry Smith struct DMCompositeLink *next; 102347c6ae99SBarry Smith PetscMPIInt rank; 102447c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 102571b14b3eSStefano Zampini PetscBool flg; 102647c6ae99SBarry Smith 102747c6ae99SBarry Smith PetscFunctionBegin; 102847c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 10307a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 10317a8be351SBarry Smith PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before"); 10329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM, is)); 103347c6ae99SBarry Smith next = com->next; 10349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 103547c6ae99SBarry Smith 103615229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 103747c6ae99SBarry Smith while (next) { 1038e5e52638SMatthew G. Knepley PetscDS prob; 1039e5e52638SMatthew G. Knepley 10409566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt])); 10419566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 1042e5e52638SMatthew G. Knepley if (prob) { 104365c226d8SMatthew G. Knepley MatNullSpace space; 104465c226d8SMatthew G. Knepley Mat pmat; 10450f21e855SMatthew G. Knepley PetscObject disc; 10460f21e855SMatthew G. Knepley PetscInt Nf; 104765c226d8SMatthew G. Knepley 10489566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 1049f24dd8d2SMatthew G. Knepley if (cnt < Nf) { 10509566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(prob, cnt, &disc)); 10519566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space)); 10529566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space)); 10539566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space)); 10549566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space)); 10559566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat)); 10569566063dSJacob Faibussowitsch if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat)); 105765c226d8SMatthew G. Knepley } 1058f24dd8d2SMatthew G. Knepley } 105947c6ae99SBarry Smith cnt++; 106047c6ae99SBarry Smith next = next->next; 106147c6ae99SBarry Smith } 10623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106347c6ae99SBarry Smith } 106447c6ae99SBarry Smith 106566976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) 1066d71ae5a4SJacob Faibussowitsch { 10674d343eeaSMatthew G Knepley PetscInt nDM; 10684d343eeaSMatthew G Knepley DM *dms; 10694d343eeaSMatthew G Knepley PetscInt i; 10704d343eeaSMatthew G Knepley 10714d343eeaSMatthew G Knepley PetscFunctionBegin; 10729566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 10738865f1eaSKarl Rupp if (numFields) *numFields = nDM; 10749566063dSJacob Faibussowitsch PetscCall(DMCompositeGetGlobalISs(dm, fields)); 10754d343eeaSMatthew G Knepley if (fieldNames) { 10769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, &dms)); 10779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, fieldNames)); 10789566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms)); 10794d343eeaSMatthew G Knepley for (i = 0; i < nDM; i++) { 10804d343eeaSMatthew G Knepley char buf[256]; 10814d343eeaSMatthew G Knepley const char *splitname; 10824d343eeaSMatthew G Knepley 10834d343eeaSMatthew G Knepley /* Split naming precedence: object name, prefix, number */ 10844d343eeaSMatthew G Knepley splitname = ((PetscObject)dm)->name; 10854d343eeaSMatthew G Knepley if (!splitname) { 10869566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname)); 10874d343eeaSMatthew G Knepley if (splitname) { 10884d343eeaSMatthew G Knepley size_t len; 10899566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(buf, splitname, sizeof(buf))); 10908caf3d72SBarry Smith buf[sizeof(buf) - 1] = 0; 10919566063dSJacob Faibussowitsch PetscCall(PetscStrlen(buf, &len)); 10924d343eeaSMatthew G Knepley if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */ 10934d343eeaSMatthew G Knepley splitname = buf; 10944d343eeaSMatthew G Knepley } 10954d343eeaSMatthew G Knepley } 10964d343eeaSMatthew G Knepley if (!splitname) { 109763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i)); 10984d343eeaSMatthew G Knepley splitname = buf; 10994d343eeaSMatthew G Knepley } 11009566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i])); 11014d343eeaSMatthew G Knepley } 11029566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 11034d343eeaSMatthew G Knepley } 11043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11054d343eeaSMatthew G Knepley } 11064d343eeaSMatthew G Knepley 1107e7c4fc90SDmitry Karpeev /* 1108e7c4fc90SDmitry Karpeev This could take over from DMCreateFieldIS(), as it is more general, 11090298fd71SBarry Smith making DMCreateFieldIS() a special case -- calling with dmlist == NULL; 1110e7c4fc90SDmitry Karpeev At this point it's probably best to be less intrusive, however. 1111e7c4fc90SDmitry Karpeev */ 111266976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 1113d71ae5a4SJacob Faibussowitsch { 1114e7c4fc90SDmitry Karpeev PetscInt nDM; 1115e7c4fc90SDmitry Karpeev PetscInt i; 1116e7c4fc90SDmitry Karpeev 1117e7c4fc90SDmitry Karpeev PetscFunctionBegin; 11189566063dSJacob Faibussowitsch PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist)); 1119e7c4fc90SDmitry Karpeev if (dmlist) { 11209566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 11219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, dmlist)); 11229566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, *dmlist)); 112348a46eb9SPierre Jolivet for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i]))); 1124e7c4fc90SDmitry Karpeev } 11253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1126e7c4fc90SDmitry Karpeev } 1127e7c4fc90SDmitry Karpeev 112847c6ae99SBarry Smith /*@C 1129dce8aebaSBarry Smith DMCompositeGetLocalVectors - Gets local vectors for each part of a `DMCOMPOSITE` 1130dce8aebaSBarry Smith Use `DMCompositeRestoreLocalVectors()` to return them. 113147c6ae99SBarry Smith 113220f4b53cSBarry Smith Not Collective; No Fortran Support 113347c6ae99SBarry Smith 113447c6ae99SBarry Smith Input Parameter: 1135dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 113647c6ae99SBarry Smith 113747c6ae99SBarry Smith Output Parameter: 1138a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`s 113947c6ae99SBarry Smith 114047c6ae99SBarry Smith Level: advanced 114147c6ae99SBarry Smith 1142dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1143db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1144db781477SPatrick Sanan `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 114547c6ae99SBarry Smith @*/ 1146d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...) 1147d71ae5a4SJacob Faibussowitsch { 114847c6ae99SBarry Smith va_list Argp; 114947c6ae99SBarry Smith struct DMCompositeLink *next; 115047c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 115171b14b3eSStefano Zampini PetscBool flg; 115247c6ae99SBarry Smith 115347c6ae99SBarry Smith PetscFunctionBegin; 115447c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11559566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 11567a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 115747c6ae99SBarry Smith next = com->next; 115815229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 115947c6ae99SBarry Smith va_start(Argp, dm); 116047c6ae99SBarry Smith while (next) { 116147c6ae99SBarry Smith Vec *vec; 116247c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 11639566063dSJacob Faibussowitsch if (vec) PetscCall(DMGetLocalVector(next->dm, vec)); 116447c6ae99SBarry Smith next = next->next; 116547c6ae99SBarry Smith } 116647c6ae99SBarry Smith va_end(Argp); 11673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116847c6ae99SBarry Smith } 116947c6ae99SBarry Smith 117047c6ae99SBarry Smith /*@C 1171dce8aebaSBarry Smith DMCompositeRestoreLocalVectors - Restores local vectors for each part of a `DMCOMPOSITE` 117247c6ae99SBarry Smith 117320f4b53cSBarry Smith Not Collective; No Fortran Support 117447c6ae99SBarry Smith 117547c6ae99SBarry Smith Input Parameter: 1176dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 117747c6ae99SBarry Smith 117847c6ae99SBarry Smith Output Parameter: 1179a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec` 118047c6ae99SBarry Smith 118147c6ae99SBarry Smith Level: advanced 118247c6ae99SBarry Smith 1183dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1184db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1185db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 118647c6ae99SBarry Smith @*/ 1187d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...) 1188d71ae5a4SJacob Faibussowitsch { 118947c6ae99SBarry Smith va_list Argp; 119047c6ae99SBarry Smith struct DMCompositeLink *next; 119147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 119271b14b3eSStefano Zampini PetscBool flg; 119347c6ae99SBarry Smith 119447c6ae99SBarry Smith PetscFunctionBegin; 119547c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 11977a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 119847c6ae99SBarry Smith next = com->next; 119915229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 120047c6ae99SBarry Smith va_start(Argp, dm); 120147c6ae99SBarry Smith while (next) { 120247c6ae99SBarry Smith Vec *vec; 120347c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 12049566063dSJacob Faibussowitsch if (vec) PetscCall(DMRestoreLocalVector(next->dm, vec)); 120547c6ae99SBarry Smith next = next->next; 120647c6ae99SBarry Smith } 120747c6ae99SBarry Smith va_end(Argp); 12083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120947c6ae99SBarry Smith } 121047c6ae99SBarry Smith 121147c6ae99SBarry Smith /*@C 1212dce8aebaSBarry Smith DMCompositeGetEntries - Gets the `DM` for each entry in a `DMCOMPOSITE`. 121347c6ae99SBarry Smith 121447c6ae99SBarry Smith Not Collective 121547c6ae99SBarry Smith 121647c6ae99SBarry Smith Input Parameter: 1217dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 121847c6ae99SBarry Smith 121947c6ae99SBarry Smith Output Parameter: 1220a4e35b19SJacob Faibussowitsch . ... - the individual entries `DM` 122147c6ae99SBarry Smith 122247c6ae99SBarry Smith Level: advanced 122347c6ae99SBarry Smith 122460225df5SJacob Faibussowitsch Fortran Notes: 12255d83a8b1SBarry Smith Use `DMCompositeGetEntriesArray()` 1226f5f57ec0SBarry Smith 1227dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()` 1228db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 122960225df5SJacob Faibussowitsch `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()` 123047c6ae99SBarry Smith @*/ 1231d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntries(DM dm, ...) 1232d71ae5a4SJacob Faibussowitsch { 123347c6ae99SBarry Smith va_list Argp; 123447c6ae99SBarry Smith struct DMCompositeLink *next; 123547c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 123671b14b3eSStefano Zampini PetscBool flg; 123747c6ae99SBarry Smith 123847c6ae99SBarry Smith PetscFunctionBegin; 123947c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12409566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 12417a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 124247c6ae99SBarry Smith next = com->next; 124315229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 124447c6ae99SBarry Smith va_start(Argp, dm); 124547c6ae99SBarry Smith while (next) { 124647c6ae99SBarry Smith DM *dmn; 124747c6ae99SBarry Smith dmn = va_arg(Argp, DM *); 12489ae5db72SJed Brown if (dmn) *dmn = next->dm; 124947c6ae99SBarry Smith next = next->next; 125047c6ae99SBarry Smith } 125147c6ae99SBarry Smith va_end(Argp); 12523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125347c6ae99SBarry Smith } 125447c6ae99SBarry Smith 1255cc4c1da9SBarry Smith /*@ 1256dce8aebaSBarry Smith DMCompositeGetEntriesArray - Gets the DM for each entry in a `DMCOMPOSITE` 12572fa5ba8aSJed Brown 12582fa5ba8aSJed Brown Not Collective 12592fa5ba8aSJed Brown 12602fa5ba8aSJed Brown Input Parameter: 1261dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 1262907376e6SBarry Smith 1263907376e6SBarry Smith Output Parameter: 1264dce8aebaSBarry Smith . dms - array of sufficient length (see `DMCompositeGetNumberDM()`) to hold the individual `DM` 12652fa5ba8aSJed Brown 12662fa5ba8aSJed Brown Level: advanced 12672fa5ba8aSJed Brown 1268dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()` 1269db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 127060225df5SJacob Faibussowitsch `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()` 12712fa5ba8aSJed Brown @*/ 1272d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntriesArray(DM dm, DM dms[]) 1273d71ae5a4SJacob Faibussowitsch { 12742fa5ba8aSJed Brown struct DMCompositeLink *next; 12752fa5ba8aSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 12762fa5ba8aSJed Brown PetscInt i; 127771b14b3eSStefano Zampini PetscBool flg; 12782fa5ba8aSJed Brown 12792fa5ba8aSJed Brown PetscFunctionBegin; 12802fa5ba8aSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 12827a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 128315229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 12842fa5ba8aSJed Brown for (next = com->next, i = 0; next; next = next->next, i++) dms[i] = next->dm; 12853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12862fa5ba8aSJed Brown } 12872fa5ba8aSJed Brown 1288e10fd815SStefano Zampini typedef struct { 1289e10fd815SStefano Zampini DM dm; 1290e10fd815SStefano Zampini PetscViewer *subv; 1291e10fd815SStefano Zampini Vec *vecs; 1292e10fd815SStefano Zampini } GLVisViewerCtx; 1293e10fd815SStefano Zampini 1294e6aa7a3bSBarry Smith static PetscErrorCode DestroyGLVisViewerCtx_Private(void **vctx) 1295d71ae5a4SJacob Faibussowitsch { 1296e6aa7a3bSBarry Smith GLVisViewerCtx *ctx = *(GLVisViewerCtx **)vctx; 1297e10fd815SStefano Zampini PetscInt i, n; 1298e10fd815SStefano Zampini 1299e10fd815SStefano Zampini PetscFunctionBegin; 13009566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm, &n)); 130148a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(PetscViewerDestroy(&ctx->subv[i])); 13029566063dSJacob Faibussowitsch PetscCall(PetscFree2(ctx->subv, ctx->vecs)); 13039566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx->dm)); 13049566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 13053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1306e10fd815SStefano Zampini } 1307e10fd815SStefano Zampini 1308d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 1309d71ae5a4SJacob Faibussowitsch { 1310e10fd815SStefano Zampini Vec X = (Vec)oX; 1311e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 1312e10fd815SStefano Zampini PetscInt i, n, cumf; 1313e10fd815SStefano Zampini 1314e10fd815SStefano Zampini PetscFunctionBegin; 13159566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm, &n)); 13169566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(ctx->dm, X, n, NULL, ctx->vecs)); 1317e10fd815SStefano Zampini for (i = 0, cumf = 0; i < n; i++) { 1318e10fd815SStefano Zampini PetscErrorCode (*g2l)(PetscObject, PetscInt, PetscObject[], void *); 1319e10fd815SStefano Zampini void *fctx; 1320e10fd815SStefano Zampini PetscInt nfi; 1321e10fd815SStefano Zampini 132234e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nfi, NULL, NULL, &g2l, NULL, &fctx)); 1323e10fd815SStefano Zampini if (!nfi) continue; 13241baa6e33SBarry Smith if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i], nfi, oXfield + cumf, fctx)); 1325f4f49eeaSPierre Jolivet else PetscCall(VecCopy(ctx->vecs[i], (Vec)oXfield[cumf])); 1326e10fd815SStefano Zampini cumf += nfi; 1327e10fd815SStefano Zampini } 13289566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(ctx->dm, X, n, NULL, ctx->vecs)); 13293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1330e10fd815SStefano Zampini } 1331e10fd815SStefano Zampini 1332d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer) 1333d71ae5a4SJacob Faibussowitsch { 1334e10fd815SStefano Zampini DM dm = (DM)odm, *dms; 1335e10fd815SStefano Zampini Vec *Ufds; 1336e10fd815SStefano Zampini GLVisViewerCtx *ctx; 1337e10fd815SStefano Zampini PetscInt i, n, tnf, *sdim; 1338e10fd815SStefano Zampini char **fecs; 1339e10fd815SStefano Zampini 1340e10fd815SStefano Zampini PetscFunctionBegin; 13419566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 13429566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 1343e10fd815SStefano Zampini ctx->dm = dm; 13449566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &n)); 13459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dms)); 13469566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms)); 13479566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &ctx->subv, n, &ctx->vecs)); 1348e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1349e10fd815SStefano Zampini PetscInt nf; 1350e10fd815SStefano Zampini 13519566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PetscObjectComm(odm), &ctx->subv[i])); 13529566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(ctx->subv[i], PETSCVIEWERGLVIS)); 135334e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisSetDM_Internal(ctx->subv[i], (PetscObject)dms[i])); 135434e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, NULL, NULL, NULL, NULL, NULL)); 1355e10fd815SStefano Zampini tnf += nf; 1356e10fd815SStefano Zampini } 13579566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 13589566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(tnf, &fecs, tnf, &sdim, tnf, &Ufds)); 1359e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1360e10fd815SStefano Zampini PetscInt *sd, nf, f; 1361e10fd815SStefano Zampini const char **fec; 1362e10fd815SStefano Zampini Vec *Uf; 1363e10fd815SStefano Zampini 136434e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, &fec, &sd, NULL, (PetscObject **)&Uf, NULL)); 1365e10fd815SStefano Zampini for (f = 0; f < nf; f++) { 13669566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec[f], &fecs[tnf + f])); 1367e10fd815SStefano Zampini Ufds[tnf + f] = Uf[f]; 1368e10fd815SStefano Zampini sdim[tnf + f] = sd[f]; 1369e10fd815SStefano Zampini } 1370e10fd815SStefano Zampini tnf += nf; 1371e10fd815SStefano Zampini } 13729566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetFields(viewer, tnf, (const char **)fecs, sdim, DMCompositeSampleGLVisFields_Private, (PetscObject *)Ufds, ctx, DestroyGLVisViewerCtx_Private)); 137348a46eb9SPierre Jolivet for (i = 0; i < tnf; i++) PetscCall(PetscFree(fecs[i])); 13749566063dSJacob Faibussowitsch PetscCall(PetscFree3(fecs, sdim, Ufds)); 13753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1376e10fd815SStefano Zampini } 1377e10fd815SStefano Zampini 137866976f2fSJacob Faibussowitsch static PetscErrorCode DMRefine_Composite(DM dmi, MPI_Comm comm, DM *fine) 1379d71ae5a4SJacob Faibussowitsch { 138047c6ae99SBarry Smith struct DMCompositeLink *next; 138147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dmi->data; 138247c6ae99SBarry Smith DM dm; 138347c6ae99SBarry Smith 138447c6ae99SBarry Smith PetscFunctionBegin; 138547c6ae99SBarry Smith PetscValidHeaderSpecific(dmi, DM_CLASSID, 1); 138648a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm)); 13879566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi)); 138847c6ae99SBarry Smith next = com->next; 13899566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm, fine)); 139047c6ae99SBarry Smith 139115229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 139247c6ae99SBarry Smith while (next) { 13939566063dSJacob Faibussowitsch PetscCall(DMRefine(next->dm, comm, &dm)); 13949566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine, dm)); 13959566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 139647c6ae99SBarry Smith next = next->next; 139747c6ae99SBarry Smith } 13983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 139947c6ae99SBarry Smith } 140047c6ae99SBarry Smith 140166976f2fSJacob Faibussowitsch static PetscErrorCode DMCoarsen_Composite(DM dmi, MPI_Comm comm, DM *fine) 1402d71ae5a4SJacob Faibussowitsch { 140314354c39SJed Brown struct DMCompositeLink *next; 140414354c39SJed Brown DM_Composite *com = (DM_Composite *)dmi->data; 140514354c39SJed Brown DM dm; 140614354c39SJed Brown 140714354c39SJed Brown PetscFunctionBegin; 140814354c39SJed Brown PetscValidHeaderSpecific(dmi, DM_CLASSID, 1); 14099566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi)); 141048a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm)); 141114354c39SJed Brown next = com->next; 14129566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm, fine)); 141314354c39SJed Brown 141415229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 141514354c39SJed Brown while (next) { 14169566063dSJacob Faibussowitsch PetscCall(DMCoarsen(next->dm, comm, &dm)); 14179566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine, dm)); 14189566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 141914354c39SJed Brown next = next->next; 142014354c39SJed Brown } 14213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 142214354c39SJed Brown } 142347c6ae99SBarry Smith 142466976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_Composite(DM coarse, DM fine, Mat *A, Vec *v) 1425d71ae5a4SJacob Faibussowitsch { 14269ae5db72SJed Brown PetscInt m, n, M, N, nDM, i; 142747c6ae99SBarry Smith struct DMCompositeLink *nextc; 142847c6ae99SBarry Smith struct DMCompositeLink *nextf; 142925296bd5SBarry Smith Vec gcoarse, gfine, *vecs; 143047c6ae99SBarry Smith DM_Composite *comcoarse = (DM_Composite *)coarse->data; 143147c6ae99SBarry Smith DM_Composite *comfine = (DM_Composite *)fine->data; 14329ae5db72SJed Brown Mat *mats; 143347c6ae99SBarry Smith 143447c6ae99SBarry Smith PetscFunctionBegin; 143547c6ae99SBarry Smith PetscValidHeaderSpecific(coarse, DM_CLASSID, 1); 143647c6ae99SBarry Smith PetscValidHeaderSpecific(fine, DM_CLASSID, 2); 14379566063dSJacob Faibussowitsch PetscCall(DMSetUp(coarse)); 14389566063dSJacob Faibussowitsch PetscCall(DMSetUp(fine)); 143947c6ae99SBarry Smith /* use global vectors only for determining matrix layout */ 14409566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(coarse, &gcoarse)); 14419566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fine, &gfine)); 14429566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gcoarse, &n)); 14439566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gfine, &m)); 14449566063dSJacob Faibussowitsch PetscCall(VecGetSize(gcoarse, &N)); 14459566063dSJacob Faibussowitsch PetscCall(VecGetSize(gfine, &M)); 14469566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(coarse, &gcoarse)); 14479566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fine, &gfine)); 144847c6ae99SBarry Smith 14499ae5db72SJed Brown nDM = comfine->nDM; 145063a3b9bcSJacob 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); 14519566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nDM * nDM, &mats)); 145248a46eb9SPierre Jolivet if (v) PetscCall(PetscCalloc1(nDM, &vecs)); 145347c6ae99SBarry Smith 145415229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 14559ae5db72SJed Brown for (nextc = comcoarse->next, nextf = comfine->next, i = 0; nextc; nextc = nextc->next, nextf = nextf->next, i++) { 14561baa6e33SBarry Smith if (!v) PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], NULL)); 14571baa6e33SBarry Smith else PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], &vecs[i])); 145847c6ae99SBarry Smith } 14599566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, nDM, NULL, mats, A)); 14601baa6e33SBarry Smith if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, vecs, v)); 14619566063dSJacob Faibussowitsch for (i = 0; i < nDM * nDM; i++) PetscCall(MatDestroy(&mats[i])); 14629566063dSJacob Faibussowitsch PetscCall(PetscFree(mats)); 146325296bd5SBarry Smith if (v) { 14649566063dSJacob Faibussowitsch for (i = 0; i < nDM; i++) PetscCall(VecDestroy(&vecs[i])); 14659566063dSJacob Faibussowitsch PetscCall(PetscFree(vecs)); 146625296bd5SBarry Smith } 14673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 146847c6ae99SBarry Smith } 146947c6ae99SBarry Smith 1470d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm) 1471d71ae5a4SJacob Faibussowitsch { 14721411c6eeSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 14731411c6eeSJed Brown ISLocalToGlobalMapping *ltogs; 1474f7efa3c7SJed Brown PetscInt i; 14751411c6eeSJed Brown 14761411c6eeSJed Brown PetscFunctionBegin; 14771411c6eeSJed Brown /* Set the ISLocalToGlobalMapping on the new matrix */ 14789566063dSJacob Faibussowitsch PetscCall(DMCompositeGetISLocalToGlobalMappings(dm, <ogs)); 14799566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm), com->nDM, ltogs, &dm->ltogmap)); 14809566063dSJacob Faibussowitsch for (i = 0; i < com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(<ogs[i])); 14819566063dSJacob Faibussowitsch PetscCall(PetscFree(ltogs)); 14823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14831411c6eeSJed Brown } 14841411c6eeSJed Brown 148566976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateColoring_Composite(DM dm, ISColoringType ctype, ISColoring *coloring) 1486d71ae5a4SJacob Faibussowitsch { 148747c6ae99SBarry Smith PetscInt n, i, cnt; 148847c6ae99SBarry Smith ISColoringValue *colors; 148947c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 149047c6ae99SBarry Smith ISColoringValue maxcol = 0; 149147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 149247c6ae99SBarry Smith 149347c6ae99SBarry Smith PetscFunctionBegin; 149447c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 149508401ef6SPierre Jolivet PetscCheck(ctype != IS_COLORING_LOCAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only global coloring supported"); 1496966bd95aSPierre Jolivet PetscCheck(ctype == IS_COLORING_GLOBAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unknown ISColoringType"); 149747c6ae99SBarry Smith n = com->n; 14989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &colors)); /* freed in ISColoringDestroy() */ 149947c6ae99SBarry Smith 15009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL)); 150147c6ae99SBarry Smith if (dense) { 15026497c311SBarry Smith PetscCall(ISColoringValueCast(com->N, &maxcol)); 15036497c311SBarry Smith for (i = 0; i < n; i++) PetscCall(ISColoringValueCast(com->rstart + i, colors + i)); 150447c6ae99SBarry Smith } else { 150547c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 150647c6ae99SBarry Smith PetscMPIInt rank; 150747c6ae99SBarry Smith 15089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 150947c6ae99SBarry Smith cnt = 0; 151047c6ae99SBarry Smith while (next) { 151147c6ae99SBarry Smith ISColoring lcoloring; 151247c6ae99SBarry Smith 15139566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(next->dm, IS_COLORING_GLOBAL, &lcoloring)); 15146497c311SBarry Smith for (i = 0; i < lcoloring->N; i++) PetscCall(ISColoringValueCast(maxcol + lcoloring->colors[i], colors + cnt++)); 151547c6ae99SBarry Smith maxcol += lcoloring->n; 15169566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&lcoloring)); 151747c6ae99SBarry Smith next = next->next; 151847c6ae99SBarry Smith } 151947c6ae99SBarry Smith } 15209566063dSJacob Faibussowitsch PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), maxcol, n, colors, PETSC_OWN_POINTER, coloring)); 15213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152247c6ae99SBarry Smith } 152347c6ae99SBarry Smith 152466976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) 1525d71ae5a4SJacob Faibussowitsch { 152647c6ae99SBarry Smith struct DMCompositeLink *next; 1527*51a5d4e4SZach Atkins const PetscScalar *garray, *garrayhead; 1528*51a5d4e4SZach Atkins PetscScalar *larray, *larrayhead; 152947c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 153047c6ae99SBarry Smith 153147c6ae99SBarry Smith PetscFunctionBegin; 153247c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 153347c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 153439d35262SVincent Le Chenadec 153548a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 153639d35262SVincent Le Chenadec 1537*51a5d4e4SZach Atkins PetscCall(VecGetArrayRead(gvec, &garray)); 1538*51a5d4e4SZach Atkins garrayhead = garray; 15399566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec, &larray)); 1540*51a5d4e4SZach Atkins larrayhead = larray; 154147c6ae99SBarry Smith 154215229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 154339d35262SVincent Le Chenadec next = com->next; 154447c6ae99SBarry Smith while (next) { 154547c6ae99SBarry Smith Vec local, global; 154647c6ae99SBarry Smith PetscInt N; 154747c6ae99SBarry Smith 15489566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 15499566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global, &N)); 1550*51a5d4e4SZach Atkins PetscCall(VecPlaceArray(global, garrayhead)); 15519566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local)); 1552*51a5d4e4SZach Atkins PetscCall(VecPlaceArray(local, larrayhead)); 15539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, mode, local)); 15549566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, mode, local)); 15559566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 15569566063dSJacob Faibussowitsch PetscCall(VecResetArray(local)); 15579566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 15589566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local)); 155939d35262SVincent Le Chenadec 1560*51a5d4e4SZach Atkins larrayhead += next->nlocal; 1561*51a5d4e4SZach Atkins garrayhead += next->n; 156247c6ae99SBarry Smith next = next->next; 156347c6ae99SBarry Smith } 1564*51a5d4e4SZach Atkins PetscCall(VecRestoreArrayRead(gvec, &garray)); 1565*51a5d4e4SZach Atkins PetscCall(VecRestoreArray(lvec, &larray)); 15663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 156747c6ae99SBarry Smith } 156847c6ae99SBarry Smith 156966976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) 1570d71ae5a4SJacob Faibussowitsch { 15710c010503SBarry Smith PetscFunctionBegin; 157239d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 157339d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 157439d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 4); 15753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 157639d35262SVincent Le Chenadec } 157739d35262SVincent Le Chenadec 157866976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1579d71ae5a4SJacob Faibussowitsch { 158039d35262SVincent Le Chenadec struct DMCompositeLink *next; 1581*51a5d4e4SZach Atkins const PetscScalar *larray, *larrayhead; 1582*51a5d4e4SZach Atkins PetscScalar *garray, *garrayhead; 158339d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite *)dm->data; 158439d35262SVincent Le Chenadec 158539d35262SVincent Le Chenadec PetscFunctionBegin; 158639d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 158739d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 158839d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 158939d35262SVincent Le Chenadec 159048a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 159139d35262SVincent Le Chenadec 1592*51a5d4e4SZach Atkins PetscCall(VecGetArrayRead(lvec, &larray)); 1593*51a5d4e4SZach Atkins larrayhead = larray; 15949566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &garray)); 1595*51a5d4e4SZach Atkins garrayhead = garray; 159639d35262SVincent Le Chenadec 159715229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 159839d35262SVincent Le Chenadec next = com->next; 159939d35262SVincent Le Chenadec while (next) { 160039d35262SVincent Le Chenadec Vec global, local; 160139d35262SVincent Le Chenadec 16029566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local)); 1603*51a5d4e4SZach Atkins PetscCall(VecPlaceArray(local, larrayhead)); 16049566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 1605*51a5d4e4SZach Atkins PetscCall(VecPlaceArray(global, garrayhead)); 16069566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, local, mode, global)); 16079566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, local, mode, global)); 16089566063dSJacob Faibussowitsch PetscCall(VecResetArray(local)); 16099566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 16109566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 16119566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local)); 161239d35262SVincent Le Chenadec 1613*51a5d4e4SZach Atkins garrayhead += next->n; 1614*51a5d4e4SZach Atkins larrayhead += next->nlocal; 161539d35262SVincent Le Chenadec next = next->next; 161639d35262SVincent Le Chenadec } 161739d35262SVincent Le Chenadec 1618*51a5d4e4SZach Atkins PetscCall(VecRestoreArray(gvec, &garray)); 1619*51a5d4e4SZach Atkins PetscCall(VecRestoreArrayRead(lvec, &larray)); 16203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 162139d35262SVincent Le Chenadec } 162239d35262SVincent Le Chenadec 162366976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1624d71ae5a4SJacob Faibussowitsch { 162539d35262SVincent Le Chenadec PetscFunctionBegin; 162639d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 162739d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 162839d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 16293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 163039d35262SVincent Le Chenadec } 163139d35262SVincent Le Chenadec 163266976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2) 1633d71ae5a4SJacob Faibussowitsch { 163439d35262SVincent Le Chenadec struct DMCompositeLink *next; 1635*51a5d4e4SZach Atkins const PetscScalar *array1, *array1head; 1636*51a5d4e4SZach Atkins PetscScalar *array2, *array2head; 163739d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite *)dm->data; 163839d35262SVincent Le Chenadec 163939d35262SVincent Le Chenadec PetscFunctionBegin; 164039d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 164139d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2); 164239d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4); 164339d35262SVincent Le Chenadec 164448a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 164539d35262SVincent Le Chenadec 1646*51a5d4e4SZach Atkins PetscCall(VecGetArrayRead(vec1, &array1)); 1647*51a5d4e4SZach Atkins array1head = array1; 16489566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec2, &array2)); 1649*51a5d4e4SZach Atkins array2head = array2; 165039d35262SVincent Le Chenadec 165115229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 165239d35262SVincent Le Chenadec next = com->next; 165339d35262SVincent Le Chenadec while (next) { 165439d35262SVincent Le Chenadec Vec local1, local2; 165539d35262SVincent Le Chenadec 16569566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local1)); 1657*51a5d4e4SZach Atkins PetscCall(VecPlaceArray(local1, array1head)); 16589566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local2)); 1659*51a5d4e4SZach Atkins PetscCall(VecPlaceArray(local2, array2head)); 16609566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2)); 16619566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2)); 16629566063dSJacob Faibussowitsch PetscCall(VecResetArray(local2)); 16639566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local2)); 16649566063dSJacob Faibussowitsch PetscCall(VecResetArray(local1)); 16659566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local1)); 166639d35262SVincent Le Chenadec 1667*51a5d4e4SZach Atkins array1head += next->nlocal; 1668*51a5d4e4SZach Atkins array2head += next->nlocal; 166939d35262SVincent Le Chenadec next = next->next; 167039d35262SVincent Le Chenadec } 167139d35262SVincent Le Chenadec 1672*51a5d4e4SZach Atkins PetscCall(VecRestoreArrayRead(vec1, &array1)); 1673*51a5d4e4SZach Atkins PetscCall(VecRestoreArray(vec2, &array2)); 16743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 167539d35262SVincent Le Chenadec } 167639d35262SVincent Le Chenadec 167766976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1678d71ae5a4SJacob Faibussowitsch { 167939d35262SVincent Le Chenadec PetscFunctionBegin; 168039d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 168139d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 168239d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 16833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16840c010503SBarry Smith } 168547c6ae99SBarry Smith 16866ae3a549SBarry Smith /*MC 1687dce8aebaSBarry Smith DMCOMPOSITE = "composite" - A `DM` object that is used to manage data for a collection of `DM` 16886ae3a549SBarry Smith 16896ae3a549SBarry Smith Level: intermediate 16906ae3a549SBarry Smith 1691db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()` 16926ae3a549SBarry Smith M*/ 16936ae3a549SBarry Smith 1694d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p) 1695d71ae5a4SJacob Faibussowitsch { 1696a4121054SBarry Smith DM_Composite *com; 1697a4121054SBarry Smith 1698a4121054SBarry Smith PetscFunctionBegin; 16994dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&com)); 1700a4121054SBarry Smith p->data = com; 1701a4121054SBarry Smith com->n = 0; 17027ac2b803SAlex Fikl com->nghost = 0; 17030298fd71SBarry Smith com->next = NULL; 1704a4121054SBarry Smith com->nDM = 0; 1705a4121054SBarry Smith 1706a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Composite; 1707a4121054SBarry Smith p->ops->createlocalvector = DMCreateLocalVector_Composite; 1708184d77edSJed Brown p->ops->getlocaltoglobalmapping = DMGetLocalToGlobalMapping_Composite; 17094d343eeaSMatthew G Knepley p->ops->createfieldis = DMCreateFieldIS_Composite; 171016621825SDmitry Karpeev p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite; 1711a4121054SBarry Smith p->ops->refine = DMRefine_Composite; 171214354c39SJed Brown p->ops->coarsen = DMCoarsen_Composite; 171325296bd5SBarry Smith p->ops->createinterpolation = DMCreateInterpolation_Composite; 171425296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Composite; 1715e727c939SJed Brown p->ops->getcoloring = DMCreateColoring_Composite; 1716a4121054SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite; 1717a4121054SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite; 171839d35262SVincent Le Chenadec p->ops->localtoglobalbegin = DMLocalToGlobalBegin_Composite; 171939d35262SVincent Le Chenadec p->ops->localtoglobalend = DMLocalToGlobalEnd_Composite; 172039d35262SVincent Le Chenadec p->ops->localtolocalbegin = DMLocalToLocalBegin_Composite; 172139d35262SVincent Le Chenadec p->ops->localtolocalend = DMLocalToLocalEnd_Composite; 1722a4121054SBarry Smith p->ops->destroy = DMDestroy_Composite; 1723a4121054SBarry Smith p->ops->view = DMView_Composite; 1724a4121054SBarry Smith p->ops->setup = DMSetUp_Composite; 1725e10fd815SStefano Zampini 17269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite)); 17273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1728a4121054SBarry Smith } 1729a4121054SBarry Smith 17301fd49c25SBarry Smith /*@ 1731dce8aebaSBarry Smith DMCompositeCreate - Creates a `DMCOMPOSITE`, used to generate "composite" 17320c010503SBarry Smith vectors made up of several subvectors. 17330c010503SBarry Smith 1734d083f849SBarry Smith Collective 173547c6ae99SBarry Smith 173647c6ae99SBarry Smith Input Parameter: 17370c010503SBarry Smith . comm - the processors that will share the global vector 17380c010503SBarry Smith 17392fe279fdSBarry Smith Output Parameter: 1740dce8aebaSBarry Smith . packer - the `DMCOMPOSITE` object 174147c6ae99SBarry Smith 174247c6ae99SBarry Smith Level: advanced 174347c6ae99SBarry Smith 174460225df5SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCreate()` 1745db781477SPatrick Sanan `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()` 1746db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 174747c6ae99SBarry Smith @*/ 1748d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer) 1749d71ae5a4SJacob Faibussowitsch { 175047c6ae99SBarry Smith PetscFunctionBegin; 17514f572ea9SToby Isaac PetscAssertPointer(packer, 2); 17529566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, packer)); 17539566063dSJacob Faibussowitsch PetscCall(DMSetType(*packer, DMCOMPOSITE)); 17543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 175547c6ae99SBarry Smith } 1756