1ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h> /*I "petscdmcomposite.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/isimpl.h> 3e10fd815SStefano Zampini #include <petsc/private/glvisviewerimpl.h> 42764a2aaSMatthew G. Knepley #include <petscds.h> 547c6ae99SBarry Smith 647c6ae99SBarry Smith /*@C 747c6ae99SBarry Smith DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the 8dce8aebaSBarry Smith separate components `DM` in a `DMCOMPOSITE` to build the correct matrix nonzero structure. 947c6ae99SBarry Smith 1020f4b53cSBarry Smith Logically Collective; No Fortran Support 1147c6ae99SBarry Smith 12d8d19677SJose E. Roman Input Parameters: 1347c6ae99SBarry Smith + dm - the composite object 1460225df5SJacob Faibussowitsch - FormCoupleLocations - routine to set the nonzero locations in the matrix 1547c6ae99SBarry Smith 1647c6ae99SBarry Smith Level: advanced 1747c6ae99SBarry Smith 18dce8aebaSBarry Smith Note: 19dce8aebaSBarry Smith See `DMSetApplicationContext()` and `DMGetApplicationContext()` for how to get user information into 2047c6ae99SBarry Smith this routine 2147c6ae99SBarry Smith 22dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM` 2347c6ae99SBarry Smith @*/ 24d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeSetCoupling(DM dm, PetscErrorCode (*FormCoupleLocations)(DM, Mat, PetscInt *, PetscInt *, PetscInt, PetscInt, PetscInt, PetscInt)) 25d71ae5a4SJacob Faibussowitsch { 2647c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 2771b14b3eSStefano Zampini PetscBool flg; 2847c6ae99SBarry Smith 2947c6ae99SBarry Smith PetscFunctionBegin; 309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 317a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 3247c6ae99SBarry Smith com->FormCoupleLocations = FormCoupleLocations; 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3447c6ae99SBarry Smith } 3547c6ae99SBarry Smith 3666976f2fSJacob Faibussowitsch static PetscErrorCode DMDestroy_Composite(DM dm) 37d71ae5a4SJacob Faibussowitsch { 3847c6ae99SBarry Smith struct DMCompositeLink *next, *prev; 3947c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 4047c6ae99SBarry Smith 4147c6ae99SBarry Smith PetscFunctionBegin; 4247c6ae99SBarry Smith next = com->next; 4347c6ae99SBarry Smith while (next) { 4447c6ae99SBarry Smith prev = next; 4547c6ae99SBarry Smith next = next->next; 469566063dSJacob Faibussowitsch PetscCall(DMDestroy(&prev->dm)); 479566063dSJacob Faibussowitsch PetscCall(PetscFree(prev->grstarts)); 489566063dSJacob Faibussowitsch PetscCall(PetscFree(prev)); 4947c6ae99SBarry Smith } 509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL)); 51435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 529566063dSJacob Faibussowitsch PetscCall(PetscFree(com)); 533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5447c6ae99SBarry Smith } 5547c6ae99SBarry Smith 5666976f2fSJacob Faibussowitsch static PetscErrorCode DMView_Composite(DM dm, PetscViewer v) 57d71ae5a4SJacob Faibussowitsch { 5847c6ae99SBarry Smith PetscBool iascii; 5947c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 6047c6ae99SBarry Smith 6147c6ae99SBarry Smith PetscFunctionBegin; 629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii)); 6347c6ae99SBarry Smith if (iascii) { 6447c6ae99SBarry Smith struct DMCompositeLink *lnk = com->next; 6547c6ae99SBarry Smith PetscInt i; 6647c6ae99SBarry Smith 679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "DM (%s)\n", ((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix")); 6863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, " contains %" PetscInt_FMT " DMs\n", com->nDM)); 699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v)); 7047c6ae99SBarry Smith for (i = 0; lnk; lnk = lnk->next, i++) { 7163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "Link %" PetscInt_FMT ": DM of type %s\n", i, ((PetscObject)lnk->dm)->type_name)); 729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v)); 739566063dSJacob Faibussowitsch PetscCall(DMView(lnk->dm, v)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v)); 7547c6ae99SBarry Smith } 769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v)); 7747c6ae99SBarry Smith } 783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7947c6ae99SBarry Smith } 8047c6ae99SBarry Smith 8147c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/ 8266976f2fSJacob Faibussowitsch static PetscErrorCode DMSetUp_Composite(DM dm) 83d71ae5a4SJacob Faibussowitsch { 8447c6ae99SBarry Smith PetscInt nprev = 0; 8547c6ae99SBarry Smith PetscMPIInt rank, size; 8647c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 8747c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 8847c6ae99SBarry Smith PetscLayout map; 8947c6ae99SBarry Smith 9047c6ae99SBarry Smith PetscFunctionBegin; 9128b400f6SJacob Faibussowitsch PetscCheck(!com->setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Packer has already been setup"); 929566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &map)); 939566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(map, com->n)); 949566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetSize(map, PETSC_DETERMINE)); 959566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(map, 1)); 969566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(map)); 979566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &com->N)); 989566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(map, &com->rstart, NULL)); 999566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&map)); 10047c6ae99SBarry Smith 1019ae5db72SJed Brown /* now set the rstart for each linked vector */ 1029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 10447c6ae99SBarry Smith while (next) { 10547c6ae99SBarry Smith next->rstart = nprev; 10606ebdd98SJed Brown nprev += next->n; 10747c6ae99SBarry Smith next->grstart = com->rstart + next->rstart; 1089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &next->grstarts)); 1099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&next->grstart, 1, MPIU_INT, next->grstarts, 1, MPIU_INT, PetscObjectComm((PetscObject)dm))); 11047c6ae99SBarry Smith next = next->next; 11147c6ae99SBarry Smith } 11247c6ae99SBarry Smith com->setup = PETSC_TRUE; 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11447c6ae99SBarry Smith } 11547c6ae99SBarry Smith 11647c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/ 11747c6ae99SBarry Smith 11873e31fe2SJed Brown /*@ 119aaa8cc7dSPierre Jolivet DMCompositeGetNumberDM - Gets the number of `DM` objects in the `DMCOMPOSITE` 12047c6ae99SBarry Smith representation. 12147c6ae99SBarry Smith 12247c6ae99SBarry Smith Not Collective 12347c6ae99SBarry Smith 12447c6ae99SBarry Smith Input Parameter: 125dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 12647c6ae99SBarry Smith 12747c6ae99SBarry Smith Output Parameter: 128dce8aebaSBarry Smith . nDM - the number of `DM` 12947c6ae99SBarry Smith 13047c6ae99SBarry Smith Level: beginner 13147c6ae99SBarry Smith 132dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM` 13347c6ae99SBarry Smith @*/ 134d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetNumberDM(DM dm, PetscInt *nDM) 135d71ae5a4SJacob Faibussowitsch { 13647c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 13771b14b3eSStefano Zampini PetscBool flg; 1385fd66863SKarl Rupp 13947c6ae99SBarry Smith PetscFunctionBegin; 14047c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1419566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1427a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 14347c6ae99SBarry Smith *nDM = com->nDM; 1443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14547c6ae99SBarry Smith } 14647c6ae99SBarry Smith 14747c6ae99SBarry Smith /*@C 14847c6ae99SBarry Smith DMCompositeGetAccess - Allows one to access the individual packed vectors in their global 14947c6ae99SBarry Smith representation. 15047c6ae99SBarry Smith 15120f4b53cSBarry Smith Collective 15247c6ae99SBarry Smith 1539ae5db72SJed Brown Input Parameters: 154dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 1559ae5db72SJed Brown - gvec - the global vector 1569ae5db72SJed Brown 1572fe279fdSBarry Smith Output Parameter: 158a4e35b19SJacob Faibussowitsch . ... - the packed parallel vectors, `NULL` for those that are not needed 15947c6ae99SBarry Smith 16047c6ae99SBarry Smith Level: advanced 16147c6ae99SBarry Smith 162dce8aebaSBarry Smith Note: 163dce8aebaSBarry Smith Use `DMCompositeRestoreAccess()` to return the vectors when you no longer need them 164dce8aebaSBarry Smith 16560225df5SJacob Faibussowitsch Fortran Notes: 166dce8aebaSBarry Smith Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4) 167dce8aebaSBarry Smith or use the alternative interface `DMCompositeGetAccessArray()`. 168dce8aebaSBarry Smith 169dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetEntries()`, `DMCompositeScatter()` 17047c6ae99SBarry Smith @*/ 171d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccess(DM dm, Vec gvec, ...) 172d71ae5a4SJacob Faibussowitsch { 17347c6ae99SBarry Smith va_list Argp; 17447c6ae99SBarry Smith struct DMCompositeLink *next; 17547c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 1765edff71fSBarry Smith PetscInt readonly; 17771b14b3eSStefano Zampini PetscBool flg; 17847c6ae99SBarry Smith 17947c6ae99SBarry Smith PetscFunctionBegin; 18047c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 18147c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 1829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1837a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 18447c6ae99SBarry Smith next = com->next; 18548a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 18647c6ae99SBarry Smith 1879566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec, &readonly)); 18815229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 18947c6ae99SBarry Smith va_start(Argp, gvec); 19047c6ae99SBarry Smith while (next) { 19147c6ae99SBarry Smith Vec *vec; 19247c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 1939ae5db72SJed Brown if (vec) { 1949566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, vec)); 1955edff71fSBarry Smith if (readonly) { 1965edff71fSBarry Smith const PetscScalar *array; 1979566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 1989566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(*vec, array + next->rstart)); 1999566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(*vec)); 2009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 2015edff71fSBarry Smith } else { 2025edff71fSBarry Smith PetscScalar *array; 2039566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 2049566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(*vec, array + next->rstart)); 2059566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 20647c6ae99SBarry Smith } 2075edff71fSBarry Smith } 20847c6ae99SBarry Smith next = next->next; 20947c6ae99SBarry Smith } 21047c6ae99SBarry Smith va_end(Argp); 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21247c6ae99SBarry Smith } 21347c6ae99SBarry Smith 214f73e5cebSJed Brown /*@C 215f73e5cebSJed Brown DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global 216f73e5cebSJed Brown representation. 217f73e5cebSJed Brown 21820f4b53cSBarry Smith Collective 219f73e5cebSJed Brown 220f73e5cebSJed Brown Input Parameters: 221dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` 222f73e5cebSJed Brown . pvec - packed vector 223f73e5cebSJed Brown . nwanted - number of vectors wanted 2240298fd71SBarry Smith - wanted - sorted array of vectors wanted, or NULL to get all vectors 225f73e5cebSJed Brown 2262fe279fdSBarry Smith Output Parameter: 227f73e5cebSJed Brown . vecs - array of requested global vectors (must be allocated) 228f73e5cebSJed Brown 229f73e5cebSJed Brown Level: advanced 230f73e5cebSJed Brown 231dce8aebaSBarry Smith Note: 232dce8aebaSBarry Smith Use `DMCompositeRestoreAccessArray()` to return the vectors when you no longer need them 233dce8aebaSBarry Smith 234dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 235f73e5cebSJed Brown @*/ 236d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 237d71ae5a4SJacob Faibussowitsch { 238f73e5cebSJed Brown struct DMCompositeLink *link; 239f73e5cebSJed Brown PetscInt i, wnum; 240f73e5cebSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 241bee642f7SBarry Smith PetscInt readonly; 24271b14b3eSStefano Zampini PetscBool flg; 243f73e5cebSJed Brown 244f73e5cebSJed Brown PetscFunctionBegin; 245f73e5cebSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 246f73e5cebSJed Brown PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 2479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 2487a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 24948a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 250f73e5cebSJed Brown 2519566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 252f73e5cebSJed Brown for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 253f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 254f73e5cebSJed Brown Vec v; 2559566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(link->dm, &v)); 256bee642f7SBarry Smith if (readonly) { 257bee642f7SBarry Smith const PetscScalar *array; 2589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvec, &array)); 2599566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + link->rstart)); 2609566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(v)); 2619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvec, &array)); 262bee642f7SBarry Smith } else { 263bee642f7SBarry Smith PetscScalar *array; 2649566063dSJacob Faibussowitsch PetscCall(VecGetArray(pvec, &array)); 2659566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + link->rstart)); 2669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pvec, &array)); 267bee642f7SBarry Smith } 268f73e5cebSJed Brown vecs[wnum++] = v; 269f73e5cebSJed Brown } 270f73e5cebSJed Brown } 2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 272f73e5cebSJed Brown } 273f73e5cebSJed Brown 2747ac2b803SAlex Fikl /*@C 2757ac2b803SAlex Fikl DMCompositeGetLocalAccessArray - Allows one to access the individual 2767ac2b803SAlex Fikl packed vectors in their local representation. 2777ac2b803SAlex Fikl 27820f4b53cSBarry Smith Collective 2797ac2b803SAlex Fikl 2807ac2b803SAlex Fikl Input Parameters: 281dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` 2827ac2b803SAlex Fikl . pvec - packed vector 2837ac2b803SAlex Fikl . nwanted - number of vectors wanted 2847ac2b803SAlex Fikl - wanted - sorted array of vectors wanted, or NULL to get all vectors 2857ac2b803SAlex Fikl 2862fe279fdSBarry Smith Output Parameter: 2877ac2b803SAlex Fikl . vecs - array of requested local vectors (must be allocated) 2887ac2b803SAlex Fikl 2897ac2b803SAlex Fikl Level: advanced 2907ac2b803SAlex Fikl 291dce8aebaSBarry Smith Note: 292dce8aebaSBarry Smith Use `DMCompositeRestoreLocalAccessArray()` to return the vectors 293dce8aebaSBarry Smith when you no longer need them. 294dce8aebaSBarry Smith 295dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`, 296db781477SPatrick Sanan `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 2977ac2b803SAlex Fikl @*/ 298d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 299d71ae5a4SJacob Faibussowitsch { 3007ac2b803SAlex Fikl struct DMCompositeLink *link; 3017ac2b803SAlex Fikl PetscInt i, wnum; 3027ac2b803SAlex Fikl DM_Composite *com = (DM_Composite *)dm->data; 3037ac2b803SAlex Fikl PetscInt readonly; 3047ac2b803SAlex Fikl PetscInt nlocal = 0; 30571b14b3eSStefano Zampini PetscBool flg; 3067ac2b803SAlex Fikl 3077ac2b803SAlex Fikl PetscFunctionBegin; 3087ac2b803SAlex Fikl PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3097ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 3109566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 3117a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 31248a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 3137ac2b803SAlex Fikl 3149566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 3157ac2b803SAlex Fikl for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 3167ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 3177ac2b803SAlex Fikl Vec v; 3189566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(link->dm, &v)); 3197ac2b803SAlex Fikl if (readonly) { 3207ac2b803SAlex Fikl const PetscScalar *array; 3219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvec, &array)); 3229566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + nlocal)); 323b1c3483dSMark Adams // this method does not make sense. The local vectors are not updated with a global-to-local and the user can not do it because it is locked 3249566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(v)); 3259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvec, &array)); 3267ac2b803SAlex Fikl } else { 3277ac2b803SAlex Fikl PetscScalar *array; 3289566063dSJacob Faibussowitsch PetscCall(VecGetArray(pvec, &array)); 3299566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + nlocal)); 3309566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pvec, &array)); 3317ac2b803SAlex Fikl } 3327ac2b803SAlex Fikl vecs[wnum++] = v; 3337ac2b803SAlex Fikl } 3347ac2b803SAlex Fikl 3357ac2b803SAlex Fikl nlocal += link->nlocal; 3367ac2b803SAlex Fikl } 3373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3387ac2b803SAlex Fikl } 3397ac2b803SAlex Fikl 34047c6ae99SBarry Smith /*@C 341dce8aebaSBarry Smith DMCompositeRestoreAccess - Returns the vectors obtained with `DMCompositeGetAccess()` 34247c6ae99SBarry Smith representation. 34347c6ae99SBarry Smith 34420f4b53cSBarry Smith Collective 34547c6ae99SBarry Smith 3469ae5db72SJed Brown Input Parameters: 347dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 34847c6ae99SBarry Smith . gvec - the global vector 349a4e35b19SJacob Faibussowitsch - ... - the individual parallel vectors, `NULL` for those that are not needed 35047c6ae99SBarry Smith 35147c6ae99SBarry Smith Level: advanced 35247c6ae99SBarry Smith 353dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 354db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`, 35542747ad1SJacob Faibussowitsch `DMCompositeGetAccess()` 35647c6ae99SBarry Smith @*/ 357d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...) 358d71ae5a4SJacob Faibussowitsch { 35947c6ae99SBarry Smith va_list Argp; 36047c6ae99SBarry Smith struct DMCompositeLink *next; 36147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 3625edff71fSBarry Smith PetscInt readonly; 36371b14b3eSStefano Zampini PetscBool flg; 36447c6ae99SBarry Smith 36547c6ae99SBarry Smith PetscFunctionBegin; 36647c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 36747c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 3689566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 3697a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 37047c6ae99SBarry Smith next = com->next; 37148a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 37247c6ae99SBarry Smith 3739566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec, &readonly)); 37415229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 37547c6ae99SBarry Smith va_start(Argp, gvec); 37647c6ae99SBarry Smith while (next) { 37747c6ae99SBarry Smith Vec *vec; 37847c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 3799ae5db72SJed Brown if (vec) { 3809566063dSJacob Faibussowitsch PetscCall(VecResetArray(*vec)); 3811baa6e33SBarry Smith if (readonly) PetscCall(VecLockReadPop(*vec)); 3829566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, vec)); 38347c6ae99SBarry Smith } 38447c6ae99SBarry Smith next = next->next; 38547c6ae99SBarry Smith } 38647c6ae99SBarry Smith va_end(Argp); 3873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38847c6ae99SBarry Smith } 38947c6ae99SBarry Smith 390f73e5cebSJed Brown /*@C 391dce8aebaSBarry Smith DMCompositeRestoreAccessArray - Returns the vectors obtained with `DMCompositeGetAccessArray()` 392f73e5cebSJed Brown 39320f4b53cSBarry Smith Collective 394f73e5cebSJed Brown 395f73e5cebSJed Brown Input Parameters: 396dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 397f73e5cebSJed Brown . pvec - packed vector 398f73e5cebSJed Brown . nwanted - number of vectors wanted 3990298fd71SBarry Smith . wanted - sorted array of vectors wanted, or NULL to get all vectors 400f73e5cebSJed Brown - vecs - array of global vectors to return 401f73e5cebSJed Brown 402f73e5cebSJed Brown Level: advanced 403f73e5cebSJed Brown 404dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 405f73e5cebSJed Brown @*/ 406d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 407d71ae5a4SJacob Faibussowitsch { 408f73e5cebSJed Brown struct DMCompositeLink *link; 409f73e5cebSJed Brown PetscInt i, wnum; 410f73e5cebSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 411bee642f7SBarry Smith PetscInt readonly; 41271b14b3eSStefano Zampini PetscBool flg; 413f73e5cebSJed Brown 414f73e5cebSJed Brown PetscFunctionBegin; 415f73e5cebSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 416f73e5cebSJed Brown PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 4179566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 4187a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 41948a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 420f73e5cebSJed Brown 4219566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 422f73e5cebSJed Brown for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 423f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 4249566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum])); 42548a46eb9SPierre Jolivet if (readonly) PetscCall(VecLockReadPop(vecs[wnum])); 4269566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum])); 427f73e5cebSJed Brown wnum++; 428f73e5cebSJed Brown } 429f73e5cebSJed Brown } 4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 431f73e5cebSJed Brown } 432f73e5cebSJed Brown 4337ac2b803SAlex Fikl /*@C 434dce8aebaSBarry Smith DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with `DMCompositeGetLocalAccessArray()`. 4357ac2b803SAlex Fikl 43620f4b53cSBarry Smith Collective 4377ac2b803SAlex Fikl 4387ac2b803SAlex Fikl Input Parameters: 439dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 4407ac2b803SAlex Fikl . pvec - packed vector 4417ac2b803SAlex Fikl . nwanted - number of vectors wanted 4427ac2b803SAlex Fikl . wanted - sorted array of vectors wanted, or NULL to restore all vectors 4437ac2b803SAlex Fikl - vecs - array of local vectors to return 4447ac2b803SAlex Fikl 4457ac2b803SAlex Fikl Level: advanced 4467ac2b803SAlex Fikl 447dce8aebaSBarry Smith Note: 448dce8aebaSBarry Smith nwanted and wanted must match the values given to `DMCompositeGetLocalAccessArray()` 4497ac2b803SAlex Fikl otherwise the call will fail. 4507ac2b803SAlex Fikl 451dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`, 452db781477SPatrick Sanan `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, 453db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeGather()` 4547ac2b803SAlex Fikl @*/ 455d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 456d71ae5a4SJacob Faibussowitsch { 4577ac2b803SAlex Fikl struct DMCompositeLink *link; 4587ac2b803SAlex Fikl PetscInt i, wnum; 4597ac2b803SAlex Fikl DM_Composite *com = (DM_Composite *)dm->data; 4607ac2b803SAlex Fikl PetscInt readonly; 46171b14b3eSStefano Zampini PetscBool flg; 4627ac2b803SAlex Fikl 4637ac2b803SAlex Fikl PetscFunctionBegin; 4647ac2b803SAlex Fikl PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4657ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 4669566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 4677a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 46848a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 4697ac2b803SAlex Fikl 4709566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 4717ac2b803SAlex Fikl for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 4727ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 4739566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum])); 47448a46eb9SPierre Jolivet if (readonly) PetscCall(VecLockReadPop(vecs[wnum])); 4759566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum])); 4767ac2b803SAlex Fikl wnum++; 4777ac2b803SAlex Fikl } 4787ac2b803SAlex Fikl } 4793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4807ac2b803SAlex Fikl } 4817ac2b803SAlex Fikl 48247c6ae99SBarry Smith /*@C 48347c6ae99SBarry Smith DMCompositeScatter - Scatters from a global packed vector into its individual local vectors 48447c6ae99SBarry Smith 48520f4b53cSBarry Smith Collective 48647c6ae99SBarry Smith 4879ae5db72SJed Brown Input Parameters: 488dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 48947c6ae99SBarry Smith . gvec - the global vector 490a4e35b19SJacob Faibussowitsch - ... - the individual sequential vectors, `NULL` for those that are not needed 49147c6ae99SBarry Smith 49247c6ae99SBarry Smith Level: advanced 49347c6ae99SBarry Smith 494dce8aebaSBarry Smith Note: 495dce8aebaSBarry Smith `DMCompositeScatterArray()` is a non-variadic alternative that is often more convenient for library callers and is 4966f3c3dcfSJed Brown accessible from Fortran. 4976f3c3dcfSJed Brown 498dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 499db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 500db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 501db781477SPatrick Sanan `DMCompositeScatterArray()` 50247c6ae99SBarry Smith @*/ 503d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...) 504d71ae5a4SJacob Faibussowitsch { 50547c6ae99SBarry Smith va_list Argp; 50647c6ae99SBarry Smith struct DMCompositeLink *next; 507c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt; 50847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 50971b14b3eSStefano Zampini PetscBool flg; 51047c6ae99SBarry Smith 51147c6ae99SBarry Smith PetscFunctionBegin; 51247c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 51347c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 5149566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 5157a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 51648a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 51747c6ae99SBarry Smith 51815229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 51947c6ae99SBarry Smith va_start(Argp, gvec); 5208fd8f222SJed Brown for (cnt = 3, next = com->next; next; cnt++, next = next->next) { 5219ae5db72SJed Brown Vec local; 5229ae5db72SJed Brown local = va_arg(Argp, Vec); 5239ae5db72SJed Brown if (local) { 5249ae5db72SJed Brown Vec global; 5255edff71fSBarry Smith const PetscScalar *array; 52663a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt)); 5279566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 5289566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 5299566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 5309566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local)); 5319566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local)); 5329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 5339566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 5349566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 53547c6ae99SBarry Smith } 53647c6ae99SBarry Smith } 53747c6ae99SBarry Smith va_end(Argp); 5383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53947c6ae99SBarry Smith } 54047c6ae99SBarry Smith 5416f3c3dcfSJed Brown /*@ 5426f3c3dcfSJed Brown DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors 5436f3c3dcfSJed Brown 54420f4b53cSBarry Smith Collective 5456f3c3dcfSJed Brown 5466f3c3dcfSJed Brown Input Parameters: 547dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 5486f3c3dcfSJed Brown . gvec - the global vector 549a2b725a8SWilliam Gropp - lvecs - array of local vectors, NULL for any that are not needed 5506f3c3dcfSJed Brown 5516f3c3dcfSJed Brown Level: advanced 5526f3c3dcfSJed Brown 5536f3c3dcfSJed Brown Note: 554dce8aebaSBarry Smith This is a non-variadic alternative to `DMCompositeScatter()` 5556f3c3dcfSJed Brown 556dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()` 557db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 558db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 5596f3c3dcfSJed Brown @*/ 560d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs) 561d71ae5a4SJacob Faibussowitsch { 5626f3c3dcfSJed Brown struct DMCompositeLink *next; 5636f3c3dcfSJed Brown PetscInt i; 5646f3c3dcfSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 56571b14b3eSStefano Zampini PetscBool flg; 5666f3c3dcfSJed Brown 5676f3c3dcfSJed Brown PetscFunctionBegin; 5686f3c3dcfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5696f3c3dcfSJed Brown PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 5709566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 5717a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 57248a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 5736f3c3dcfSJed Brown 57415229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 5756f3c3dcfSJed Brown for (i = 0, next = com->next; next; next = next->next, i++) { 5766f3c3dcfSJed Brown if (lvecs[i]) { 5776f3c3dcfSJed Brown Vec global; 578c5d31e75SLisandro Dalcin const PetscScalar *array; 5796f3c3dcfSJed Brown PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3); 5809566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 5819566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 5829566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart)); 5839566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i])); 5849566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i])); 5859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 5869566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 5879566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 5886f3c3dcfSJed Brown } 5896f3c3dcfSJed Brown } 5903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5916f3c3dcfSJed Brown } 5926f3c3dcfSJed Brown 59347c6ae99SBarry Smith /*@C 59447c6ae99SBarry Smith DMCompositeGather - Gathers into a global packed vector from its individual local vectors 59547c6ae99SBarry Smith 59620f4b53cSBarry Smith Collective 59747c6ae99SBarry Smith 598d8d19677SJose E. Roman Input Parameters: 599dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 600dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES` 601a4e35b19SJacob Faibussowitsch . gvec - the global vector 602a4e35b19SJacob Faibussowitsch - ... - the individual sequential vectors, `NULL` for any that are not needed 60347c6ae99SBarry Smith 60447c6ae99SBarry Smith Level: advanced 60547c6ae99SBarry Smith 60660225df5SJacob Faibussowitsch Fortran Notes: 607dce8aebaSBarry Smith Fortran users should use `DMCompositeGatherArray()` 608f5f57ec0SBarry Smith 609dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 610db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 611db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 61247c6ae99SBarry Smith @*/ 613d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...) 614d71ae5a4SJacob Faibussowitsch { 61547c6ae99SBarry Smith va_list Argp; 61647c6ae99SBarry Smith struct DMCompositeLink *next; 61747c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 618c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt; 61971b14b3eSStefano Zampini PetscBool flg; 62047c6ae99SBarry Smith 62147c6ae99SBarry Smith PetscFunctionBegin; 62247c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 623064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3); 6249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 6257a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 62648a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 62747c6ae99SBarry Smith 62815229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 6291dac896bSSatish Balay va_start(Argp, gvec); 6308fd8f222SJed Brown for (cnt = 3, next = com->next; next; cnt++, next = next->next) { 6319ae5db72SJed Brown Vec local; 6329ae5db72SJed Brown local = va_arg(Argp, Vec); 6339ae5db72SJed Brown if (local) { 63447c6ae99SBarry Smith PetscScalar *array; 6359ae5db72SJed Brown Vec global; 63663a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt)); 6379566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 6389566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 6399566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 6409566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global)); 6419566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global)); 6429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 6439566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 6449566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 64547c6ae99SBarry Smith } 64647c6ae99SBarry Smith } 64747c6ae99SBarry Smith va_end(Argp); 6483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64947c6ae99SBarry Smith } 65047c6ae99SBarry Smith 6516f3c3dcfSJed Brown /*@ 6526f3c3dcfSJed Brown DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors 6536f3c3dcfSJed Brown 65420f4b53cSBarry Smith Collective 6556f3c3dcfSJed Brown 656d8d19677SJose E. Roman Input Parameters: 657dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 6586f3c3dcfSJed Brown . gvec - the global vector 659dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES` 6606f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed 6616f3c3dcfSJed Brown 6626f3c3dcfSJed Brown Level: advanced 6636f3c3dcfSJed Brown 664dce8aebaSBarry Smith Note: 665dce8aebaSBarry Smith This is a non-variadic alternative to `DMCompositeGather()`. 6666f3c3dcfSJed Brown 667dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 668db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 669db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`, 6706f3c3dcfSJed Brown @*/ 671d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs) 672d71ae5a4SJacob Faibussowitsch { 6736f3c3dcfSJed Brown struct DMCompositeLink *next; 6746f3c3dcfSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 6756f3c3dcfSJed Brown PetscInt i; 67671b14b3eSStefano Zampini PetscBool flg; 6776f3c3dcfSJed Brown 6786f3c3dcfSJed Brown PetscFunctionBegin; 6796f3c3dcfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 680064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3); 6819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 6827a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 68348a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 6846f3c3dcfSJed Brown 68515229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 6866f3c3dcfSJed Brown for (next = com->next, i = 0; next; next = next->next, i++) { 6876f3c3dcfSJed Brown if (lvecs[i]) { 6886f3c3dcfSJed Brown PetscScalar *array; 6896f3c3dcfSJed Brown Vec global; 690064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4); 6919566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 6929566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 6939566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 6949566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global)); 6959566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global)); 6969566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 6979566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 6989566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 6996f3c3dcfSJed Brown } 7006f3c3dcfSJed Brown } 7013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7026f3c3dcfSJed Brown } 7036f3c3dcfSJed Brown 704f5f57ec0SBarry Smith /*@ 705dce8aebaSBarry Smith DMCompositeAddDM - adds a `DM` vector to a `DMCOMPOSITE` 70647c6ae99SBarry Smith 70720f4b53cSBarry Smith Collective 70847c6ae99SBarry Smith 709d8d19677SJose E. Roman Input Parameters: 710dce8aebaSBarry Smith + dmc - the `DMCOMPOSITE` object 711dce8aebaSBarry Smith - dm - the `DM` object 71247c6ae99SBarry Smith 71347c6ae99SBarry Smith Level: advanced 71447c6ae99SBarry Smith 71542747ad1SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeGather()`, `DMCreateGlobalVector()`, 716db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 717db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 71847c6ae99SBarry Smith @*/ 719d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeAddDM(DM dmc, DM dm) 720d71ae5a4SJacob Faibussowitsch { 72106ebdd98SJed Brown PetscInt n, nlocal; 72247c6ae99SBarry Smith struct DMCompositeLink *mine, *next; 72306ebdd98SJed Brown Vec global, local; 72447c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dmc->data; 72571b14b3eSStefano Zampini PetscBool flg; 72647c6ae99SBarry Smith 72747c6ae99SBarry Smith PetscFunctionBegin; 72847c6ae99SBarry Smith PetscValidHeaderSpecific(dmc, DM_CLASSID, 1); 72947c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 7309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg)); 7317a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 73247c6ae99SBarry Smith next = com->next; 73328b400f6SJacob Faibussowitsch PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite"); 73447c6ae99SBarry Smith 73547c6ae99SBarry Smith /* create new link */ 7369566063dSJacob Faibussowitsch PetscCall(PetscNew(&mine)); 7379566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 7389566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &global)); 7399566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global, &n)); 7409566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &global)); 7419566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &local)); 7429566063dSJacob Faibussowitsch PetscCall(VecGetSize(local, &nlocal)); 7439566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &local)); 7448865f1eaSKarl Rupp 74547c6ae99SBarry Smith mine->n = n; 74606ebdd98SJed Brown mine->nlocal = nlocal; 74747c6ae99SBarry Smith mine->dm = dm; 7480298fd71SBarry Smith mine->next = NULL; 74947c6ae99SBarry Smith com->n += n; 7507ac2b803SAlex Fikl com->nghost += nlocal; 75147c6ae99SBarry Smith 75247c6ae99SBarry Smith /* add to end of list */ 7538865f1eaSKarl Rupp if (!next) com->next = mine; 7548865f1eaSKarl Rupp else { 75547c6ae99SBarry Smith while (next->next) next = next->next; 75647c6ae99SBarry Smith next->next = mine; 75747c6ae99SBarry Smith } 75847c6ae99SBarry Smith com->nDM++; 75947c6ae99SBarry Smith com->nmine++; 7603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 76147c6ae99SBarry Smith } 76247c6ae99SBarry Smith 7639804daf3SBarry Smith #include <petscdraw.h> 76426887b52SJed Brown PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer); 76566976f2fSJacob Faibussowitsch static PetscErrorCode VecView_DMComposite(Vec gvec, PetscViewer viewer) 766d71ae5a4SJacob Faibussowitsch { 76747c6ae99SBarry Smith DM dm; 76847c6ae99SBarry Smith struct DMCompositeLink *next; 76947c6ae99SBarry Smith PetscBool isdraw; 770cef07954SSatish Balay DM_Composite *com; 77147c6ae99SBarry Smith 77247c6ae99SBarry Smith PetscFunctionBegin; 7739566063dSJacob Faibussowitsch PetscCall(VecGetDM(gvec, &dm)); 7747a8be351SBarry Smith PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite"); 77547c6ae99SBarry Smith com = (DM_Composite *)dm->data; 77647c6ae99SBarry Smith next = com->next; 77747c6ae99SBarry Smith 7789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 77947c6ae99SBarry Smith if (!isdraw) { 78047c6ae99SBarry Smith /* do I really want to call this? */ 7819566063dSJacob Faibussowitsch PetscCall(VecView_MPI(gvec, viewer)); 78247c6ae99SBarry Smith } else { 78347c6ae99SBarry Smith PetscInt cnt = 0; 78447c6ae99SBarry Smith 78515229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 78647c6ae99SBarry Smith while (next) { 78747c6ae99SBarry Smith Vec vec; 788ca4278abSLisandro Dalcin const PetscScalar *array; 78947c6ae99SBarry Smith PetscInt bs; 79047c6ae99SBarry Smith 7919ae5db72SJed Brown /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */ 7929566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &vec)); 7939566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 7949566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart)); 7959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 7969566063dSJacob Faibussowitsch PetscCall(VecView(vec, viewer)); 7979566063dSJacob Faibussowitsch PetscCall(VecResetArray(vec)); 7989566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vec, &bs)); 7999566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &vec)); 8009566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer, bs)); 80147c6ae99SBarry Smith cnt += bs; 80247c6ae99SBarry Smith next = next->next; 80347c6ae99SBarry Smith } 8049566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt)); 80547c6ae99SBarry Smith } 8063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80747c6ae99SBarry Smith } 80847c6ae99SBarry Smith 80966976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec) 810d71ae5a4SJacob Faibussowitsch { 81147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 81247c6ae99SBarry Smith 81347c6ae99SBarry Smith PetscFunctionBegin; 81447c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8159566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 8169566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 8179566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec)); 8189566063dSJacob Faibussowitsch PetscCall(VecSetType(*gvec, dm->vectype)); 8199566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*gvec, com->n, com->N)); 8209566063dSJacob Faibussowitsch PetscCall(VecSetDM(*gvec, dm)); 8219566063dSJacob Faibussowitsch PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (void (*)(void))VecView_DMComposite)); 8223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82347c6ae99SBarry Smith } 82447c6ae99SBarry Smith 82566976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec) 826d71ae5a4SJacob Faibussowitsch { 82747c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 82847c6ae99SBarry Smith 82947c6ae99SBarry Smith PetscFunctionBegin; 83047c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 83147c6ae99SBarry Smith if (!com->setup) { 8329566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 8339566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 83447c6ae99SBarry Smith } 8359566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, lvec)); 8369566063dSJacob Faibussowitsch PetscCall(VecSetType(*lvec, dm->vectype)); 8379566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE)); 8389566063dSJacob Faibussowitsch PetscCall(VecSetDM(*lvec, dm)); 8393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84047c6ae99SBarry Smith } 84147c6ae99SBarry Smith 84247c6ae99SBarry Smith /*@C 843dce8aebaSBarry Smith DMCompositeGetISLocalToGlobalMappings - gets an `ISLocalToGlobalMapping` for each `DM` in the `DMCOMPOSITE`, maps to the composite global space 84447c6ae99SBarry Smith 84520f4b53cSBarry Smith Collective; No Fortran Support 84647c6ae99SBarry Smith 84747c6ae99SBarry Smith Input Parameter: 848dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 84947c6ae99SBarry Smith 8502fe279fdSBarry Smith Output Parameter: 8519ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes 852dce8aebaSBarry Smith all the ghost points that individual ghosted `DMDA` may have. 85347c6ae99SBarry Smith 85447c6ae99SBarry Smith Level: advanced 85547c6ae99SBarry Smith 856dce8aebaSBarry Smith Note: 857dce8aebaSBarry Smith Each entry of ltogs should be destroyed with `ISLocalToGlobalMappingDestroy()`, the ltogs array should be freed with `PetscFree()`. 85847c6ae99SBarry Smith 859dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 860db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 861c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 86247c6ae99SBarry Smith @*/ 863d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping **ltogs) 864d71ae5a4SJacob Faibussowitsch { 86547c6ae99SBarry Smith PetscInt i, *idx, n, cnt; 86647c6ae99SBarry Smith struct DMCompositeLink *next; 86747c6ae99SBarry Smith PetscMPIInt rank; 86847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 86971b14b3eSStefano Zampini PetscBool flg; 87047c6ae99SBarry Smith 87147c6ae99SBarry Smith PetscFunctionBegin; 87247c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 8747a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 8759566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 8769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM, ltogs)); 87747c6ae99SBarry Smith next = com->next; 8789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 87947c6ae99SBarry Smith 88015229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 88147c6ae99SBarry Smith cnt = 0; 88247c6ae99SBarry Smith while (next) { 8836eb61c8cSJed Brown ISLocalToGlobalMapping ltog; 8846eb61c8cSJed Brown PetscMPIInt size; 88586994e45SJed Brown const PetscInt *suboff, *indices; 8866eb61c8cSJed Brown Vec global; 88747c6ae99SBarry Smith 8886eb61c8cSJed Brown /* Get sub-DM global indices for each local dof */ 8899566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(next->dm, <og)); 8909566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n)); 8919566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices)); 8929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idx)); 89347c6ae99SBarry Smith 8946eb61c8cSJed Brown /* Get the offsets for the sub-DM global vector */ 8959566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 8969566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(global, &suboff)); 8979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size)); 8986eb61c8cSJed Brown 8996eb61c8cSJed Brown /* Shift the sub-DM definition of the global space to the composite global space */ 9006eb61c8cSJed Brown for (i = 0; i < n; i++) { 90186994e45SJed Brown PetscInt subi = indices[i], lo = 0, hi = size, t; 902e333b035SStefano Zampini /* There's no consensus on what a negative index means, 903e333b035SStefano Zampini except for skipping when setting the values in vectors and matrices */ 9049371c9d4SSatish Balay if (subi < 0) { 9059371c9d4SSatish Balay idx[i] = subi - next->grstarts[rank]; 9069371c9d4SSatish Balay continue; 9079371c9d4SSatish Balay } 9086eb61c8cSJed Brown /* Binary search to find which rank owns subi */ 9096eb61c8cSJed Brown while (hi - lo > 1) { 9106eb61c8cSJed Brown t = lo + (hi - lo) / 2; 9116eb61c8cSJed Brown if (suboff[t] > subi) hi = t; 9126eb61c8cSJed Brown else lo = t; 9136eb61c8cSJed Brown } 9146eb61c8cSJed Brown idx[i] = subi - suboff[lo] + next->grstarts[lo]; 9156eb61c8cSJed Brown } 9169566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices)); 9179566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt])); 9189566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 91947c6ae99SBarry Smith next = next->next; 92047c6ae99SBarry Smith cnt++; 92147c6ae99SBarry Smith } 9223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 92347c6ae99SBarry Smith } 92447c6ae99SBarry Smith 92587c85e80SJed Brown /*@C 9269ae5db72SJed Brown DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector 92787c85e80SJed Brown 92820f4b53cSBarry Smith Not Collective; No Fortran Support 92987c85e80SJed Brown 9304165533cSJose E. Roman Input Parameter: 931dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` 93287c85e80SJed Brown 9334165533cSJose E. Roman Output Parameter: 93415229ffcSPierre Jolivet . is - array of serial index sets for each component of the `DMCOMPOSITE` 93587c85e80SJed Brown 93687c85e80SJed Brown Level: intermediate 93787c85e80SJed Brown 93887c85e80SJed Brown Notes: 93987c85e80SJed Brown At present, a composite local vector does not normally exist. This function is used to provide index sets for 94015229ffcSPierre Jolivet `MatGetLocalSubMatrix()`. In the future, the scatters for each entry in the `DMCOMPOSITE` may be merged into a single 9419ae5db72SJed Brown scatter to a composite local vector. The user should not typically need to know which is being done. 94287c85e80SJed Brown 943dce8aebaSBarry Smith To get the composite global indices at all local points (including ghosts), use `DMCompositeGetISLocalToGlobalMappings()`. 94487c85e80SJed Brown 945dce8aebaSBarry Smith To get index sets for pieces of the composite global vector, use `DMCompositeGetGlobalISs()`. 94687c85e80SJed Brown 947dce8aebaSBarry Smith Each returned `IS` should be destroyed with `ISDestroy()`, the array should be freed with `PetscFree()`. 94887c85e80SJed Brown 949dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`, `MatCreateLocalRef()` 95087c85e80SJed Brown @*/ 951d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalISs(DM dm, IS **is) 952d71ae5a4SJacob Faibussowitsch { 95387c85e80SJed Brown DM_Composite *com = (DM_Composite *)dm->data; 95487c85e80SJed Brown struct DMCompositeLink *link; 95587c85e80SJed Brown PetscInt cnt, start; 95671b14b3eSStefano Zampini PetscBool flg; 95787c85e80SJed Brown 95887c85e80SJed Brown PetscFunctionBegin; 95987c85e80SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 9604f572ea9SToby Isaac PetscAssertPointer(is, 2); 9619566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 9627a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 9639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nmine, is)); 96406ebdd98SJed Brown for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) { 965520db06cSJed Brown PetscInt bs; 9669566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt])); 9679566063dSJacob Faibussowitsch PetscCall(DMGetBlockSize(link->dm, &bs)); 9689566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize((*is)[cnt], bs)); 969520db06cSJed Brown } 9703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97187c85e80SJed Brown } 97287c85e80SJed Brown 97347c6ae99SBarry Smith /*@C 974dce8aebaSBarry Smith DMCompositeGetGlobalISs - Gets the index sets for each composed object in a `DMCOMPOSITE` 97547c6ae99SBarry Smith 97620f4b53cSBarry Smith Collective 97747c6ae99SBarry Smith 97847c6ae99SBarry Smith Input Parameter: 979dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 98047c6ae99SBarry Smith 9812fe279fdSBarry Smith Output Parameter: 98247c6ae99SBarry Smith . is - the array of index sets 98347c6ae99SBarry Smith 98447c6ae99SBarry Smith Level: advanced 98547c6ae99SBarry Smith 98647c6ae99SBarry Smith Notes: 987dce8aebaSBarry Smith The is entries should be destroyed with `ISDestroy()`, the is array should be freed with `PetscFree()` 98847c6ae99SBarry Smith 98947c6ae99SBarry Smith These could be used to extract a subset of vector entries for a "multi-physics" preconditioner 99047c6ae99SBarry Smith 991dce8aebaSBarry Smith Use `DMCompositeGetLocalISs()` for index sets in the packed local numbering, and 992dce8aebaSBarry Smith `DMCompositeGetISLocalToGlobalMappings()` for to map local sub-`DM` (including ghost) indices to packed global 9936eb61c8cSJed Brown indices. 99447c6ae99SBarry Smith 99560225df5SJacob Faibussowitsch Fortran Notes: 996dce8aebaSBarry Smith The output argument 'is' must be an allocated array of sufficient length, which can be learned using `DMCompositeGetNumberDM()`. 997f3cb0f7eSJed Brown 998dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 999db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 1000c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 100147c6ae99SBarry Smith @*/ 1002d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[]) 1003d71ae5a4SJacob Faibussowitsch { 100466bb578eSMark Adams PetscInt cnt = 0; 100547c6ae99SBarry Smith struct DMCompositeLink *next; 100647c6ae99SBarry Smith PetscMPIInt rank; 100747c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 100871b14b3eSStefano Zampini PetscBool flg; 100947c6ae99SBarry Smith 101047c6ae99SBarry Smith PetscFunctionBegin; 101147c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 10137a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 10147a8be351SBarry Smith PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before"); 10159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM, is)); 101647c6ae99SBarry Smith next = com->next; 10179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 101847c6ae99SBarry Smith 101915229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 102047c6ae99SBarry Smith while (next) { 1021e5e52638SMatthew G. Knepley PetscDS prob; 1022e5e52638SMatthew G. Knepley 10239566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt])); 10249566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 1025e5e52638SMatthew G. Knepley if (prob) { 102665c226d8SMatthew G. Knepley MatNullSpace space; 102765c226d8SMatthew G. Knepley Mat pmat; 10280f21e855SMatthew G. Knepley PetscObject disc; 10290f21e855SMatthew G. Knepley PetscInt Nf; 103065c226d8SMatthew G. Knepley 10319566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 1032f24dd8d2SMatthew G. Knepley if (cnt < Nf) { 10339566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(prob, cnt, &disc)); 10349566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space)); 10359566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space)); 10369566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space)); 10379566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space)); 10389566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat)); 10399566063dSJacob Faibussowitsch if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat)); 104065c226d8SMatthew G. Knepley } 1041f24dd8d2SMatthew G. Knepley } 104247c6ae99SBarry Smith cnt++; 104347c6ae99SBarry Smith next = next->next; 104447c6ae99SBarry Smith } 10453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 104647c6ae99SBarry Smith } 104747c6ae99SBarry Smith 104866976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) 1049d71ae5a4SJacob Faibussowitsch { 10504d343eeaSMatthew G Knepley PetscInt nDM; 10514d343eeaSMatthew G Knepley DM *dms; 10524d343eeaSMatthew G Knepley PetscInt i; 10534d343eeaSMatthew G Knepley 10544d343eeaSMatthew G Knepley PetscFunctionBegin; 10559566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 10568865f1eaSKarl Rupp if (numFields) *numFields = nDM; 10579566063dSJacob Faibussowitsch PetscCall(DMCompositeGetGlobalISs(dm, fields)); 10584d343eeaSMatthew G Knepley if (fieldNames) { 10599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, &dms)); 10609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, fieldNames)); 10619566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms)); 10624d343eeaSMatthew G Knepley for (i = 0; i < nDM; i++) { 10634d343eeaSMatthew G Knepley char buf[256]; 10644d343eeaSMatthew G Knepley const char *splitname; 10654d343eeaSMatthew G Knepley 10664d343eeaSMatthew G Knepley /* Split naming precedence: object name, prefix, number */ 10674d343eeaSMatthew G Knepley splitname = ((PetscObject)dm)->name; 10684d343eeaSMatthew G Knepley if (!splitname) { 10699566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname)); 10704d343eeaSMatthew G Knepley if (splitname) { 10714d343eeaSMatthew G Knepley size_t len; 10729566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(buf, splitname, sizeof(buf))); 10738caf3d72SBarry Smith buf[sizeof(buf) - 1] = 0; 10749566063dSJacob Faibussowitsch PetscCall(PetscStrlen(buf, &len)); 10754d343eeaSMatthew G Knepley if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */ 10764d343eeaSMatthew G Knepley splitname = buf; 10774d343eeaSMatthew G Knepley } 10784d343eeaSMatthew G Knepley } 10794d343eeaSMatthew G Knepley if (!splitname) { 108063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i)); 10814d343eeaSMatthew G Knepley splitname = buf; 10824d343eeaSMatthew G Knepley } 10839566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i])); 10844d343eeaSMatthew G Knepley } 10859566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 10864d343eeaSMatthew G Knepley } 10873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10884d343eeaSMatthew G Knepley } 10894d343eeaSMatthew G Knepley 1090e7c4fc90SDmitry Karpeev /* 1091e7c4fc90SDmitry Karpeev This could take over from DMCreateFieldIS(), as it is more general, 10920298fd71SBarry Smith making DMCreateFieldIS() a special case -- calling with dmlist == NULL; 1093e7c4fc90SDmitry Karpeev At this point it's probably best to be less intrusive, however. 1094e7c4fc90SDmitry Karpeev */ 109566976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 1096d71ae5a4SJacob Faibussowitsch { 1097e7c4fc90SDmitry Karpeev PetscInt nDM; 1098e7c4fc90SDmitry Karpeev PetscInt i; 1099e7c4fc90SDmitry Karpeev 1100e7c4fc90SDmitry Karpeev PetscFunctionBegin; 11019566063dSJacob Faibussowitsch PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist)); 1102e7c4fc90SDmitry Karpeev if (dmlist) { 11039566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 11049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, dmlist)); 11059566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, *dmlist)); 110648a46eb9SPierre Jolivet for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i]))); 1107e7c4fc90SDmitry Karpeev } 11083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1109e7c4fc90SDmitry Karpeev } 1110e7c4fc90SDmitry Karpeev 111147c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 111247c6ae99SBarry Smith /*@C 1113dce8aebaSBarry Smith DMCompositeGetLocalVectors - Gets local vectors for each part of a `DMCOMPOSITE` 1114dce8aebaSBarry Smith Use `DMCompositeRestoreLocalVectors()` to return them. 111547c6ae99SBarry Smith 111620f4b53cSBarry Smith Not Collective; No Fortran Support 111747c6ae99SBarry Smith 111847c6ae99SBarry Smith Input Parameter: 1119dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 112047c6ae99SBarry Smith 112147c6ae99SBarry Smith Output Parameter: 1122a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`s 112347c6ae99SBarry Smith 112447c6ae99SBarry Smith Level: advanced 112547c6ae99SBarry Smith 1126dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1127db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1128db781477SPatrick Sanan `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 112947c6ae99SBarry Smith @*/ 1130d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...) 1131d71ae5a4SJacob Faibussowitsch { 113247c6ae99SBarry Smith va_list Argp; 113347c6ae99SBarry Smith struct DMCompositeLink *next; 113447c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 113571b14b3eSStefano Zampini PetscBool flg; 113647c6ae99SBarry Smith 113747c6ae99SBarry Smith PetscFunctionBegin; 113847c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11399566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 11407a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 114147c6ae99SBarry Smith next = com->next; 114215229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 114347c6ae99SBarry Smith va_start(Argp, dm); 114447c6ae99SBarry Smith while (next) { 114547c6ae99SBarry Smith Vec *vec; 114647c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 11479566063dSJacob Faibussowitsch if (vec) PetscCall(DMGetLocalVector(next->dm, vec)); 114847c6ae99SBarry Smith next = next->next; 114947c6ae99SBarry Smith } 115047c6ae99SBarry Smith va_end(Argp); 11513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115247c6ae99SBarry Smith } 115347c6ae99SBarry Smith 115447c6ae99SBarry Smith /*@C 1155dce8aebaSBarry Smith DMCompositeRestoreLocalVectors - Restores local vectors for each part of a `DMCOMPOSITE` 115647c6ae99SBarry Smith 115720f4b53cSBarry Smith Not Collective; No Fortran Support 115847c6ae99SBarry Smith 115947c6ae99SBarry Smith Input Parameter: 1160dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 116147c6ae99SBarry Smith 116247c6ae99SBarry Smith Output Parameter: 1163a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec` 116447c6ae99SBarry Smith 116547c6ae99SBarry Smith Level: advanced 116647c6ae99SBarry Smith 1167dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1168db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1169db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 117047c6ae99SBarry Smith @*/ 1171d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...) 1172d71ae5a4SJacob Faibussowitsch { 117347c6ae99SBarry Smith va_list Argp; 117447c6ae99SBarry Smith struct DMCompositeLink *next; 117547c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 117671b14b3eSStefano Zampini PetscBool flg; 117747c6ae99SBarry Smith 117847c6ae99SBarry Smith PetscFunctionBegin; 117947c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11809566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 11817a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 118247c6ae99SBarry Smith next = com->next; 118315229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 118447c6ae99SBarry Smith va_start(Argp, dm); 118547c6ae99SBarry Smith while (next) { 118647c6ae99SBarry Smith Vec *vec; 118747c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 11889566063dSJacob Faibussowitsch if (vec) PetscCall(DMRestoreLocalVector(next->dm, vec)); 118947c6ae99SBarry Smith next = next->next; 119047c6ae99SBarry Smith } 119147c6ae99SBarry Smith va_end(Argp); 11923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119347c6ae99SBarry Smith } 119447c6ae99SBarry Smith 119547c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 119647c6ae99SBarry Smith /*@C 1197dce8aebaSBarry Smith DMCompositeGetEntries - Gets the `DM` for each entry in a `DMCOMPOSITE`. 119847c6ae99SBarry Smith 119947c6ae99SBarry Smith Not Collective 120047c6ae99SBarry Smith 120147c6ae99SBarry Smith Input Parameter: 1202dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 120347c6ae99SBarry Smith 120447c6ae99SBarry Smith Output Parameter: 1205a4e35b19SJacob Faibussowitsch . ... - the individual entries `DM` 120647c6ae99SBarry Smith 120747c6ae99SBarry Smith Level: advanced 120847c6ae99SBarry Smith 120960225df5SJacob Faibussowitsch Fortran Notes: 1210dce8aebaSBarry Smith Available as `DMCompositeGetEntries()` for one output `DM`, DMCompositeGetEntries2() for 2, etc 1211f5f57ec0SBarry Smith 1212dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()` 1213db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 121460225df5SJacob Faibussowitsch `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()` 121547c6ae99SBarry Smith @*/ 1216d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntries(DM dm, ...) 1217d71ae5a4SJacob Faibussowitsch { 121847c6ae99SBarry Smith va_list Argp; 121947c6ae99SBarry Smith struct DMCompositeLink *next; 122047c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 122171b14b3eSStefano Zampini PetscBool flg; 122247c6ae99SBarry Smith 122347c6ae99SBarry Smith PetscFunctionBegin; 122447c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 12267a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 122747c6ae99SBarry Smith next = com->next; 122815229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 122947c6ae99SBarry Smith va_start(Argp, dm); 123047c6ae99SBarry Smith while (next) { 123147c6ae99SBarry Smith DM *dmn; 123247c6ae99SBarry Smith dmn = va_arg(Argp, DM *); 12339ae5db72SJed Brown if (dmn) *dmn = next->dm; 123447c6ae99SBarry Smith next = next->next; 123547c6ae99SBarry Smith } 123647c6ae99SBarry Smith va_end(Argp); 12373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123847c6ae99SBarry Smith } 123947c6ae99SBarry Smith 1240dbab29e1SMark F. Adams /*@C 1241dce8aebaSBarry Smith DMCompositeGetEntriesArray - Gets the DM for each entry in a `DMCOMPOSITE` 12422fa5ba8aSJed Brown 12432fa5ba8aSJed Brown Not Collective 12442fa5ba8aSJed Brown 12452fa5ba8aSJed Brown Input Parameter: 1246dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 1247907376e6SBarry Smith 1248907376e6SBarry Smith Output Parameter: 1249dce8aebaSBarry Smith . dms - array of sufficient length (see `DMCompositeGetNumberDM()`) to hold the individual `DM` 12502fa5ba8aSJed Brown 12512fa5ba8aSJed Brown Level: advanced 12522fa5ba8aSJed Brown 1253dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()` 1254db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 125560225df5SJacob Faibussowitsch `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()` 12562fa5ba8aSJed Brown @*/ 1257d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntriesArray(DM dm, DM dms[]) 1258d71ae5a4SJacob Faibussowitsch { 12592fa5ba8aSJed Brown struct DMCompositeLink *next; 12602fa5ba8aSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 12612fa5ba8aSJed Brown PetscInt i; 126271b14b3eSStefano Zampini PetscBool flg; 12632fa5ba8aSJed Brown 12642fa5ba8aSJed Brown PetscFunctionBegin; 12652fa5ba8aSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12669566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 12677a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 126815229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 12692fa5ba8aSJed Brown for (next = com->next, i = 0; next; next = next->next, i++) dms[i] = next->dm; 12703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12712fa5ba8aSJed Brown } 12722fa5ba8aSJed Brown 1273e10fd815SStefano Zampini typedef struct { 1274e10fd815SStefano Zampini DM dm; 1275e10fd815SStefano Zampini PetscViewer *subv; 1276e10fd815SStefano Zampini Vec *vecs; 1277e10fd815SStefano Zampini } GLVisViewerCtx; 1278e10fd815SStefano Zampini 1279d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 1280d71ae5a4SJacob Faibussowitsch { 1281e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 1282e10fd815SStefano Zampini PetscInt i, n; 1283e10fd815SStefano Zampini 1284e10fd815SStefano Zampini PetscFunctionBegin; 12859566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm, &n)); 128648a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(PetscViewerDestroy(&ctx->subv[i])); 12879566063dSJacob Faibussowitsch PetscCall(PetscFree2(ctx->subv, ctx->vecs)); 12889566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx->dm)); 12899566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 12903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1291e10fd815SStefano Zampini } 1292e10fd815SStefano Zampini 1293d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 1294d71ae5a4SJacob Faibussowitsch { 1295e10fd815SStefano Zampini Vec X = (Vec)oX; 1296e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 1297e10fd815SStefano Zampini PetscInt i, n, cumf; 1298e10fd815SStefano Zampini 1299e10fd815SStefano Zampini PetscFunctionBegin; 13009566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm, &n)); 13019566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(ctx->dm, X, n, NULL, ctx->vecs)); 1302e10fd815SStefano Zampini for (i = 0, cumf = 0; i < n; i++) { 1303e10fd815SStefano Zampini PetscErrorCode (*g2l)(PetscObject, PetscInt, PetscObject[], void *); 1304e10fd815SStefano Zampini void *fctx; 1305e10fd815SStefano Zampini PetscInt nfi; 1306e10fd815SStefano Zampini 130734e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nfi, NULL, NULL, &g2l, NULL, &fctx)); 1308e10fd815SStefano Zampini if (!nfi) continue; 13091baa6e33SBarry Smith if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i], nfi, oXfield + cumf, fctx)); 1310*f4f49eeaSPierre Jolivet else PetscCall(VecCopy(ctx->vecs[i], (Vec)oXfield[cumf])); 1311e10fd815SStefano Zampini cumf += nfi; 1312e10fd815SStefano Zampini } 13139566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(ctx->dm, X, n, NULL, ctx->vecs)); 13143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1315e10fd815SStefano Zampini } 1316e10fd815SStefano Zampini 1317d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer) 1318d71ae5a4SJacob Faibussowitsch { 1319e10fd815SStefano Zampini DM dm = (DM)odm, *dms; 1320e10fd815SStefano Zampini Vec *Ufds; 1321e10fd815SStefano Zampini GLVisViewerCtx *ctx; 1322e10fd815SStefano Zampini PetscInt i, n, tnf, *sdim; 1323e10fd815SStefano Zampini char **fecs; 1324e10fd815SStefano Zampini 1325e10fd815SStefano Zampini PetscFunctionBegin; 13269566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 13279566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 1328e10fd815SStefano Zampini ctx->dm = dm; 13299566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &n)); 13309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dms)); 13319566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms)); 13329566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &ctx->subv, n, &ctx->vecs)); 1333e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1334e10fd815SStefano Zampini PetscInt nf; 1335e10fd815SStefano Zampini 13369566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PetscObjectComm(odm), &ctx->subv[i])); 13379566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(ctx->subv[i], PETSCVIEWERGLVIS)); 133834e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisSetDM_Internal(ctx->subv[i], (PetscObject)dms[i])); 133934e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, NULL, NULL, NULL, NULL, NULL)); 1340e10fd815SStefano Zampini tnf += nf; 1341e10fd815SStefano Zampini } 13429566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 13439566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(tnf, &fecs, tnf, &sdim, tnf, &Ufds)); 1344e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1345e10fd815SStefano Zampini PetscInt *sd, nf, f; 1346e10fd815SStefano Zampini const char **fec; 1347e10fd815SStefano Zampini Vec *Uf; 1348e10fd815SStefano Zampini 134934e79e72SJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Internal(ctx->subv[i], &nf, &fec, &sd, NULL, (PetscObject **)&Uf, NULL)); 1350e10fd815SStefano Zampini for (f = 0; f < nf; f++) { 13519566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec[f], &fecs[tnf + f])); 1352e10fd815SStefano Zampini Ufds[tnf + f] = Uf[f]; 1353e10fd815SStefano Zampini sdim[tnf + f] = sd[f]; 1354e10fd815SStefano Zampini } 1355e10fd815SStefano Zampini tnf += nf; 1356e10fd815SStefano Zampini } 13579566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetFields(viewer, tnf, (const char **)fecs, sdim, DMCompositeSampleGLVisFields_Private, (PetscObject *)Ufds, ctx, DestroyGLVisViewerCtx_Private)); 135848a46eb9SPierre Jolivet for (i = 0; i < tnf; i++) PetscCall(PetscFree(fecs[i])); 13599566063dSJacob Faibussowitsch PetscCall(PetscFree3(fecs, sdim, Ufds)); 13603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1361e10fd815SStefano Zampini } 1362e10fd815SStefano Zampini 136366976f2fSJacob Faibussowitsch static PetscErrorCode DMRefine_Composite(DM dmi, MPI_Comm comm, DM *fine) 1364d71ae5a4SJacob Faibussowitsch { 136547c6ae99SBarry Smith struct DMCompositeLink *next; 136647c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dmi->data; 136747c6ae99SBarry Smith DM dm; 136847c6ae99SBarry Smith 136947c6ae99SBarry Smith PetscFunctionBegin; 137047c6ae99SBarry Smith PetscValidHeaderSpecific(dmi, DM_CLASSID, 1); 137148a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm)); 13729566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi)); 137347c6ae99SBarry Smith next = com->next; 13749566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm, fine)); 137547c6ae99SBarry Smith 137615229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 137747c6ae99SBarry Smith while (next) { 13789566063dSJacob Faibussowitsch PetscCall(DMRefine(next->dm, comm, &dm)); 13799566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine, dm)); 13809566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 138147c6ae99SBarry Smith next = next->next; 138247c6ae99SBarry Smith } 13833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 138447c6ae99SBarry Smith } 138547c6ae99SBarry Smith 138666976f2fSJacob Faibussowitsch static PetscErrorCode DMCoarsen_Composite(DM dmi, MPI_Comm comm, DM *fine) 1387d71ae5a4SJacob Faibussowitsch { 138814354c39SJed Brown struct DMCompositeLink *next; 138914354c39SJed Brown DM_Composite *com = (DM_Composite *)dmi->data; 139014354c39SJed Brown DM dm; 139114354c39SJed Brown 139214354c39SJed Brown PetscFunctionBegin; 139314354c39SJed Brown PetscValidHeaderSpecific(dmi, DM_CLASSID, 1); 13949566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi)); 139548a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm)); 139614354c39SJed Brown next = com->next; 13979566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm, fine)); 139814354c39SJed Brown 139915229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 140014354c39SJed Brown while (next) { 14019566063dSJacob Faibussowitsch PetscCall(DMCoarsen(next->dm, comm, &dm)); 14029566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine, dm)); 14039566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 140414354c39SJed Brown next = next->next; 140514354c39SJed Brown } 14063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 140714354c39SJed Brown } 140847c6ae99SBarry Smith 140966976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_Composite(DM coarse, DM fine, Mat *A, Vec *v) 1410d71ae5a4SJacob Faibussowitsch { 14119ae5db72SJed Brown PetscInt m, n, M, N, nDM, i; 141247c6ae99SBarry Smith struct DMCompositeLink *nextc; 141347c6ae99SBarry Smith struct DMCompositeLink *nextf; 141425296bd5SBarry Smith Vec gcoarse, gfine, *vecs; 141547c6ae99SBarry Smith DM_Composite *comcoarse = (DM_Composite *)coarse->data; 141647c6ae99SBarry Smith DM_Composite *comfine = (DM_Composite *)fine->data; 14179ae5db72SJed Brown Mat *mats; 141847c6ae99SBarry Smith 141947c6ae99SBarry Smith PetscFunctionBegin; 142047c6ae99SBarry Smith PetscValidHeaderSpecific(coarse, DM_CLASSID, 1); 142147c6ae99SBarry Smith PetscValidHeaderSpecific(fine, DM_CLASSID, 2); 14229566063dSJacob Faibussowitsch PetscCall(DMSetUp(coarse)); 14239566063dSJacob Faibussowitsch PetscCall(DMSetUp(fine)); 142447c6ae99SBarry Smith /* use global vectors only for determining matrix layout */ 14259566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(coarse, &gcoarse)); 14269566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fine, &gfine)); 14279566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gcoarse, &n)); 14289566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gfine, &m)); 14299566063dSJacob Faibussowitsch PetscCall(VecGetSize(gcoarse, &N)); 14309566063dSJacob Faibussowitsch PetscCall(VecGetSize(gfine, &M)); 14319566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(coarse, &gcoarse)); 14329566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fine, &gfine)); 143347c6ae99SBarry Smith 14349ae5db72SJed Brown nDM = comfine->nDM; 143563a3b9bcSJacob 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); 14369566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nDM * nDM, &mats)); 143748a46eb9SPierre Jolivet if (v) PetscCall(PetscCalloc1(nDM, &vecs)); 143847c6ae99SBarry Smith 143915229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 14409ae5db72SJed Brown for (nextc = comcoarse->next, nextf = comfine->next, i = 0; nextc; nextc = nextc->next, nextf = nextf->next, i++) { 14411baa6e33SBarry Smith if (!v) PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], NULL)); 14421baa6e33SBarry Smith else PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], &vecs[i])); 144347c6ae99SBarry Smith } 14449566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, nDM, NULL, mats, A)); 14451baa6e33SBarry Smith if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, vecs, v)); 14469566063dSJacob Faibussowitsch for (i = 0; i < nDM * nDM; i++) PetscCall(MatDestroy(&mats[i])); 14479566063dSJacob Faibussowitsch PetscCall(PetscFree(mats)); 144825296bd5SBarry Smith if (v) { 14499566063dSJacob Faibussowitsch for (i = 0; i < nDM; i++) PetscCall(VecDestroy(&vecs[i])); 14509566063dSJacob Faibussowitsch PetscCall(PetscFree(vecs)); 145125296bd5SBarry Smith } 14523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 145347c6ae99SBarry Smith } 145447c6ae99SBarry Smith 1455d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm) 1456d71ae5a4SJacob Faibussowitsch { 14571411c6eeSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 14581411c6eeSJed Brown ISLocalToGlobalMapping *ltogs; 1459f7efa3c7SJed Brown PetscInt i; 14601411c6eeSJed Brown 14611411c6eeSJed Brown PetscFunctionBegin; 14621411c6eeSJed Brown /* Set the ISLocalToGlobalMapping on the new matrix */ 14639566063dSJacob Faibussowitsch PetscCall(DMCompositeGetISLocalToGlobalMappings(dm, <ogs)); 14649566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm), com->nDM, ltogs, &dm->ltogmap)); 14659566063dSJacob Faibussowitsch for (i = 0; i < com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(<ogs[i])); 14669566063dSJacob Faibussowitsch PetscCall(PetscFree(ltogs)); 14673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14681411c6eeSJed Brown } 14691411c6eeSJed Brown 147066976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateColoring_Composite(DM dm, ISColoringType ctype, ISColoring *coloring) 1471d71ae5a4SJacob Faibussowitsch { 147247c6ae99SBarry Smith PetscInt n, i, cnt; 147347c6ae99SBarry Smith ISColoringValue *colors; 147447c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 147547c6ae99SBarry Smith ISColoringValue maxcol = 0; 147647c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 147747c6ae99SBarry Smith 147847c6ae99SBarry Smith PetscFunctionBegin; 147947c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 148008401ef6SPierre Jolivet PetscCheck(ctype != IS_COLORING_LOCAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only global coloring supported"); 1481f7d195e4SLawrence Mitchell if (ctype == IS_COLORING_GLOBAL) { 148247c6ae99SBarry Smith n = com->n; 1483ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unknown ISColoringType"); 14849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &colors)); /* freed in ISColoringDestroy() */ 148547c6ae99SBarry Smith 14869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL)); 148747c6ae99SBarry Smith if (dense) { 1488ad540459SPierre Jolivet for (i = 0; i < n; i++) colors[i] = (ISColoringValue)(com->rstart + i); 148947c6ae99SBarry Smith maxcol = com->N; 149047c6ae99SBarry Smith } else { 149147c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 149247c6ae99SBarry Smith PetscMPIInt rank; 149347c6ae99SBarry Smith 14949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 149547c6ae99SBarry Smith cnt = 0; 149647c6ae99SBarry Smith while (next) { 149747c6ae99SBarry Smith ISColoring lcoloring; 149847c6ae99SBarry Smith 14999566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(next->dm, IS_COLORING_GLOBAL, &lcoloring)); 1500ad540459SPierre Jolivet for (i = 0; i < lcoloring->N; i++) colors[cnt++] = maxcol + lcoloring->colors[i]; 150147c6ae99SBarry Smith maxcol += lcoloring->n; 15029566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&lcoloring)); 150347c6ae99SBarry Smith next = next->next; 150447c6ae99SBarry Smith } 150547c6ae99SBarry Smith } 15069566063dSJacob Faibussowitsch PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), maxcol, n, colors, PETSC_OWN_POINTER, coloring)); 15073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150847c6ae99SBarry Smith } 150947c6ae99SBarry Smith 151066976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) 1511d71ae5a4SJacob Faibussowitsch { 151247c6ae99SBarry Smith struct DMCompositeLink *next; 151347c6ae99SBarry Smith PetscScalar *garray, *larray; 151447c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 151547c6ae99SBarry Smith 151647c6ae99SBarry Smith PetscFunctionBegin; 151747c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 151847c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 151939d35262SVincent Le Chenadec 152048a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 152139d35262SVincent Le Chenadec 15229566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &garray)); 15239566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec, &larray)); 152447c6ae99SBarry Smith 152515229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 152639d35262SVincent Le Chenadec next = com->next; 152747c6ae99SBarry Smith while (next) { 152847c6ae99SBarry Smith Vec local, global; 152947c6ae99SBarry Smith PetscInt N; 153047c6ae99SBarry Smith 15319566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 15329566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global, &N)); 15339566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, garray)); 15349566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local)); 15359566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local, larray)); 15369566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, mode, local)); 15379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, mode, local)); 15389566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 15399566063dSJacob Faibussowitsch PetscCall(VecResetArray(local)); 15409566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 15419566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local)); 154239d35262SVincent Le Chenadec 154306ebdd98SJed Brown larray += next->nlocal; 154439d35262SVincent Le Chenadec garray += next->n; 154547c6ae99SBarry Smith next = next->next; 154647c6ae99SBarry Smith } 154747c6ae99SBarry Smith 15489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, NULL)); 15499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lvec, NULL)); 15503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 155147c6ae99SBarry Smith } 155247c6ae99SBarry Smith 155366976f2fSJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) 1554d71ae5a4SJacob Faibussowitsch { 15550c010503SBarry Smith PetscFunctionBegin; 155639d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 155739d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 155839d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 4); 15593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 156039d35262SVincent Le Chenadec } 156139d35262SVincent Le Chenadec 156266976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1563d71ae5a4SJacob Faibussowitsch { 156439d35262SVincent Le Chenadec struct DMCompositeLink *next; 156539d35262SVincent Le Chenadec PetscScalar *larray, *garray; 156639d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite *)dm->data; 156739d35262SVincent Le Chenadec 156839d35262SVincent Le Chenadec PetscFunctionBegin; 156939d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 157039d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 157139d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 157239d35262SVincent Le Chenadec 157348a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 157439d35262SVincent Le Chenadec 15759566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec, &larray)); 15769566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &garray)); 157739d35262SVincent Le Chenadec 157815229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 157939d35262SVincent Le Chenadec next = com->next; 158039d35262SVincent Le Chenadec while (next) { 158139d35262SVincent Le Chenadec Vec global, local; 158239d35262SVincent Le Chenadec 15839566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local)); 15849566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local, larray)); 15859566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 15869566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, garray)); 15879566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, local, mode, global)); 15889566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, local, mode, global)); 15899566063dSJacob Faibussowitsch PetscCall(VecResetArray(local)); 15909566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 15919566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 15929566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local)); 159339d35262SVincent Le Chenadec 159439d35262SVincent Le Chenadec garray += next->n; 159539d35262SVincent Le Chenadec larray += next->nlocal; 159639d35262SVincent Le Chenadec next = next->next; 159739d35262SVincent Le Chenadec } 159839d35262SVincent Le Chenadec 15999566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, NULL)); 16009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lvec, NULL)); 16013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 160239d35262SVincent Le Chenadec } 160339d35262SVincent Le Chenadec 160466976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1605d71ae5a4SJacob Faibussowitsch { 160639d35262SVincent Le Chenadec PetscFunctionBegin; 160739d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 160839d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 160939d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 16103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 161139d35262SVincent Le Chenadec } 161239d35262SVincent Le Chenadec 161366976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2) 1614d71ae5a4SJacob Faibussowitsch { 161539d35262SVincent Le Chenadec struct DMCompositeLink *next; 161639d35262SVincent Le Chenadec PetscScalar *array1, *array2; 161739d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite *)dm->data; 161839d35262SVincent Le Chenadec 161939d35262SVincent Le Chenadec PetscFunctionBegin; 162039d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 162139d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2); 162239d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4); 162339d35262SVincent Le Chenadec 162448a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 162539d35262SVincent Le Chenadec 16269566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec1, &array1)); 16279566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec2, &array2)); 162839d35262SVincent Le Chenadec 162915229ffcSPierre Jolivet /* loop over packed objects, handling one at a time */ 163039d35262SVincent Le Chenadec next = com->next; 163139d35262SVincent Le Chenadec while (next) { 163239d35262SVincent Le Chenadec Vec local1, local2; 163339d35262SVincent Le Chenadec 16349566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local1)); 16359566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local1, array1)); 16369566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local2)); 16379566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local2, array2)); 16389566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2)); 16399566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2)); 16409566063dSJacob Faibussowitsch PetscCall(VecResetArray(local2)); 16419566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local2)); 16429566063dSJacob Faibussowitsch PetscCall(VecResetArray(local1)); 16439566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local1)); 164439d35262SVincent Le Chenadec 164539d35262SVincent Le Chenadec array1 += next->nlocal; 164639d35262SVincent Le Chenadec array2 += next->nlocal; 164739d35262SVincent Le Chenadec next = next->next; 164839d35262SVincent Le Chenadec } 164939d35262SVincent Le Chenadec 16509566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec1, NULL)); 16519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec2, NULL)); 16523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 165339d35262SVincent Le Chenadec } 165439d35262SVincent Le Chenadec 165566976f2fSJacob Faibussowitsch static PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1656d71ae5a4SJacob Faibussowitsch { 165739d35262SVincent Le Chenadec PetscFunctionBegin; 165839d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 165939d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 166039d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 16613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16620c010503SBarry Smith } 166347c6ae99SBarry Smith 16646ae3a549SBarry Smith /*MC 1665dce8aebaSBarry Smith DMCOMPOSITE = "composite" - A `DM` object that is used to manage data for a collection of `DM` 16666ae3a549SBarry Smith 16676ae3a549SBarry Smith Level: intermediate 16686ae3a549SBarry Smith 1669db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()` 16706ae3a549SBarry Smith M*/ 16716ae3a549SBarry Smith 1672d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p) 1673d71ae5a4SJacob Faibussowitsch { 1674a4121054SBarry Smith DM_Composite *com; 1675a4121054SBarry Smith 1676a4121054SBarry Smith PetscFunctionBegin; 16774dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&com)); 1678a4121054SBarry Smith p->data = com; 1679a4121054SBarry Smith com->n = 0; 16807ac2b803SAlex Fikl com->nghost = 0; 16810298fd71SBarry Smith com->next = NULL; 1682a4121054SBarry Smith com->nDM = 0; 1683a4121054SBarry Smith 1684a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Composite; 1685a4121054SBarry Smith p->ops->createlocalvector = DMCreateLocalVector_Composite; 1686184d77edSJed Brown p->ops->getlocaltoglobalmapping = DMGetLocalToGlobalMapping_Composite; 16874d343eeaSMatthew G Knepley p->ops->createfieldis = DMCreateFieldIS_Composite; 168816621825SDmitry Karpeev p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite; 1689a4121054SBarry Smith p->ops->refine = DMRefine_Composite; 169014354c39SJed Brown p->ops->coarsen = DMCoarsen_Composite; 169125296bd5SBarry Smith p->ops->createinterpolation = DMCreateInterpolation_Composite; 169225296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Composite; 1693e727c939SJed Brown p->ops->getcoloring = DMCreateColoring_Composite; 1694a4121054SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite; 1695a4121054SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite; 169639d35262SVincent Le Chenadec p->ops->localtoglobalbegin = DMLocalToGlobalBegin_Composite; 169739d35262SVincent Le Chenadec p->ops->localtoglobalend = DMLocalToGlobalEnd_Composite; 169839d35262SVincent Le Chenadec p->ops->localtolocalbegin = DMLocalToLocalBegin_Composite; 169939d35262SVincent Le Chenadec p->ops->localtolocalend = DMLocalToLocalEnd_Composite; 1700a4121054SBarry Smith p->ops->destroy = DMDestroy_Composite; 1701a4121054SBarry Smith p->ops->view = DMView_Composite; 1702a4121054SBarry Smith p->ops->setup = DMSetUp_Composite; 1703e10fd815SStefano Zampini 17049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite)); 17053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1706a4121054SBarry Smith } 1707a4121054SBarry Smith 17081fd49c25SBarry Smith /*@ 1709dce8aebaSBarry Smith DMCompositeCreate - Creates a `DMCOMPOSITE`, used to generate "composite" 17100c010503SBarry Smith vectors made up of several subvectors. 17110c010503SBarry Smith 1712d083f849SBarry Smith Collective 171347c6ae99SBarry Smith 171447c6ae99SBarry Smith Input Parameter: 17150c010503SBarry Smith . comm - the processors that will share the global vector 17160c010503SBarry Smith 17172fe279fdSBarry Smith Output Parameter: 1718dce8aebaSBarry Smith . packer - the `DMCOMPOSITE` object 171947c6ae99SBarry Smith 172047c6ae99SBarry Smith Level: advanced 172147c6ae99SBarry Smith 172260225df5SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCreate()` 1723db781477SPatrick Sanan `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()` 1724db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 172547c6ae99SBarry Smith @*/ 1726d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer) 1727d71ae5a4SJacob Faibussowitsch { 172847c6ae99SBarry Smith PetscFunctionBegin; 17294f572ea9SToby Isaac PetscAssertPointer(packer, 2); 17309566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, packer)); 17319566063dSJacob Faibussowitsch PetscCall(DMSetType(*packer, DMCOMPOSITE)); 17323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173347c6ae99SBarry Smith } 1734