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 36d71ae5a4SJacob Faibussowitsch 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 56d71ae5a4SJacob Faibussowitsch 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 /* --------------------------------------------------------------------------------------*/ 82d71ae5a4SJacob Faibussowitsch 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: 158*a4e35b19SJacob 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)); 18847c6ae99SBarry Smith /* loop over packed objects, handling one at at 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 } 3377ac2b803SAlex Fikl 3383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3397ac2b803SAlex Fikl } 3407ac2b803SAlex Fikl 34147c6ae99SBarry Smith /*@C 342dce8aebaSBarry Smith DMCompositeRestoreAccess - Returns the vectors obtained with `DMCompositeGetAccess()` 34347c6ae99SBarry Smith representation. 34447c6ae99SBarry Smith 34520f4b53cSBarry Smith Collective 34647c6ae99SBarry Smith 3479ae5db72SJed Brown Input Parameters: 348dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 34947c6ae99SBarry Smith . gvec - the global vector 350*a4e35b19SJacob Faibussowitsch - ... - the individual parallel vectors, `NULL` for those that are not needed 35147c6ae99SBarry Smith 35247c6ae99SBarry Smith Level: advanced 35347c6ae99SBarry Smith 354dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 355db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`, 35642747ad1SJacob Faibussowitsch `DMCompositeGetAccess()` 35747c6ae99SBarry Smith @*/ 358d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...) 359d71ae5a4SJacob Faibussowitsch { 36047c6ae99SBarry Smith va_list Argp; 36147c6ae99SBarry Smith struct DMCompositeLink *next; 36247c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 3635edff71fSBarry Smith PetscInt readonly; 36471b14b3eSStefano Zampini PetscBool flg; 36547c6ae99SBarry Smith 36647c6ae99SBarry Smith PetscFunctionBegin; 36747c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 36847c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 3699566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 3707a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 37147c6ae99SBarry Smith next = com->next; 37248a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 37347c6ae99SBarry Smith 3749566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec, &readonly)); 37547c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 37647c6ae99SBarry Smith va_start(Argp, gvec); 37747c6ae99SBarry Smith while (next) { 37847c6ae99SBarry Smith Vec *vec; 37947c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 3809ae5db72SJed Brown if (vec) { 3819566063dSJacob Faibussowitsch PetscCall(VecResetArray(*vec)); 3821baa6e33SBarry Smith if (readonly) PetscCall(VecLockReadPop(*vec)); 3839566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, vec)); 38447c6ae99SBarry Smith } 38547c6ae99SBarry Smith next = next->next; 38647c6ae99SBarry Smith } 38747c6ae99SBarry Smith va_end(Argp); 3883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38947c6ae99SBarry Smith } 39047c6ae99SBarry Smith 391f73e5cebSJed Brown /*@C 392dce8aebaSBarry Smith DMCompositeRestoreAccessArray - Returns the vectors obtained with `DMCompositeGetAccessArray()` 393f73e5cebSJed Brown 39420f4b53cSBarry Smith Collective 395f73e5cebSJed Brown 396f73e5cebSJed Brown Input Parameters: 397dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 398f73e5cebSJed Brown . pvec - packed vector 399f73e5cebSJed Brown . nwanted - number of vectors wanted 4000298fd71SBarry Smith . wanted - sorted array of vectors wanted, or NULL to get all vectors 401f73e5cebSJed Brown - vecs - array of global vectors to return 402f73e5cebSJed Brown 403f73e5cebSJed Brown Level: advanced 404f73e5cebSJed Brown 405dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 406f73e5cebSJed Brown @*/ 407d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 408d71ae5a4SJacob Faibussowitsch { 409f73e5cebSJed Brown struct DMCompositeLink *link; 410f73e5cebSJed Brown PetscInt i, wnum; 411f73e5cebSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 412bee642f7SBarry Smith PetscInt readonly; 41371b14b3eSStefano Zampini PetscBool flg; 414f73e5cebSJed Brown 415f73e5cebSJed Brown PetscFunctionBegin; 416f73e5cebSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 417f73e5cebSJed Brown PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 4189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 4197a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 42048a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 421f73e5cebSJed Brown 4229566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 423f73e5cebSJed Brown for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 424f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 4259566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum])); 42648a46eb9SPierre Jolivet if (readonly) PetscCall(VecLockReadPop(vecs[wnum])); 4279566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum])); 428f73e5cebSJed Brown wnum++; 429f73e5cebSJed Brown } 430f73e5cebSJed Brown } 4313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 432f73e5cebSJed Brown } 433f73e5cebSJed Brown 4347ac2b803SAlex Fikl /*@C 435dce8aebaSBarry Smith DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with `DMCompositeGetLocalAccessArray()`. 4367ac2b803SAlex Fikl 43720f4b53cSBarry Smith Collective 4387ac2b803SAlex Fikl 4397ac2b803SAlex Fikl Input Parameters: 440dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 4417ac2b803SAlex Fikl . pvec - packed vector 4427ac2b803SAlex Fikl . nwanted - number of vectors wanted 4437ac2b803SAlex Fikl . wanted - sorted array of vectors wanted, or NULL to restore all vectors 4447ac2b803SAlex Fikl - vecs - array of local vectors to return 4457ac2b803SAlex Fikl 4467ac2b803SAlex Fikl Level: advanced 4477ac2b803SAlex Fikl 448dce8aebaSBarry Smith Note: 449dce8aebaSBarry Smith nwanted and wanted must match the values given to `DMCompositeGetLocalAccessArray()` 4507ac2b803SAlex Fikl otherwise the call will fail. 4517ac2b803SAlex Fikl 452dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`, 453db781477SPatrick Sanan `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, 454db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeGather()` 4557ac2b803SAlex Fikl @*/ 456d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 457d71ae5a4SJacob Faibussowitsch { 4587ac2b803SAlex Fikl struct DMCompositeLink *link; 4597ac2b803SAlex Fikl PetscInt i, wnum; 4607ac2b803SAlex Fikl DM_Composite *com = (DM_Composite *)dm->data; 4617ac2b803SAlex Fikl PetscInt readonly; 46271b14b3eSStefano Zampini PetscBool flg; 4637ac2b803SAlex Fikl 4647ac2b803SAlex Fikl PetscFunctionBegin; 4657ac2b803SAlex Fikl PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4667ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 4679566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 4687a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 46948a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 4707ac2b803SAlex Fikl 4719566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 4727ac2b803SAlex Fikl for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 4737ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 4749566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum])); 47548a46eb9SPierre Jolivet if (readonly) PetscCall(VecLockReadPop(vecs[wnum])); 4769566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum])); 4777ac2b803SAlex Fikl wnum++; 4787ac2b803SAlex Fikl } 4797ac2b803SAlex Fikl } 4803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4817ac2b803SAlex Fikl } 4827ac2b803SAlex Fikl 48347c6ae99SBarry Smith /*@C 48447c6ae99SBarry Smith DMCompositeScatter - Scatters from a global packed vector into its individual local vectors 48547c6ae99SBarry Smith 48620f4b53cSBarry Smith Collective 48747c6ae99SBarry Smith 4889ae5db72SJed Brown Input Parameters: 489dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 49047c6ae99SBarry Smith . gvec - the global vector 491*a4e35b19SJacob Faibussowitsch - ... - the individual sequential vectors, `NULL` for those that are not needed 49247c6ae99SBarry Smith 49347c6ae99SBarry Smith Level: advanced 49447c6ae99SBarry Smith 495dce8aebaSBarry Smith Note: 496dce8aebaSBarry Smith `DMCompositeScatterArray()` is a non-variadic alternative that is often more convenient for library callers and is 4976f3c3dcfSJed Brown accessible from Fortran. 4986f3c3dcfSJed Brown 499dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 500db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 501db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 502db781477SPatrick Sanan `DMCompositeScatterArray()` 50347c6ae99SBarry Smith @*/ 504d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...) 505d71ae5a4SJacob Faibussowitsch { 50647c6ae99SBarry Smith va_list Argp; 50747c6ae99SBarry Smith struct DMCompositeLink *next; 508c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt; 50947c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 51071b14b3eSStefano Zampini PetscBool flg; 51147c6ae99SBarry Smith 51247c6ae99SBarry Smith PetscFunctionBegin; 51347c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 51447c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 5167a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 51748a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 51847c6ae99SBarry Smith 51947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 52047c6ae99SBarry Smith va_start(Argp, gvec); 5218fd8f222SJed Brown for (cnt = 3, next = com->next; next; cnt++, next = next->next) { 5229ae5db72SJed Brown Vec local; 5239ae5db72SJed Brown local = va_arg(Argp, Vec); 5249ae5db72SJed Brown if (local) { 5259ae5db72SJed Brown Vec global; 5265edff71fSBarry Smith const PetscScalar *array; 52763a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt)); 5289566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 5299566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 5309566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 5319566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local)); 5329566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local)); 5339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 5349566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 5359566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 53647c6ae99SBarry Smith } 53747c6ae99SBarry Smith } 53847c6ae99SBarry Smith va_end(Argp); 5393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54047c6ae99SBarry Smith } 54147c6ae99SBarry Smith 5426f3c3dcfSJed Brown /*@ 5436f3c3dcfSJed Brown DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors 5446f3c3dcfSJed Brown 54520f4b53cSBarry Smith Collective 5466f3c3dcfSJed Brown 5476f3c3dcfSJed Brown Input Parameters: 548dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 5496f3c3dcfSJed Brown . gvec - the global vector 550a2b725a8SWilliam Gropp - lvecs - array of local vectors, NULL for any that are not needed 5516f3c3dcfSJed Brown 5526f3c3dcfSJed Brown Level: advanced 5536f3c3dcfSJed Brown 5546f3c3dcfSJed Brown Note: 555dce8aebaSBarry Smith This is a non-variadic alternative to `DMCompositeScatter()` 5566f3c3dcfSJed Brown 557dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()` 558db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 559db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 5606f3c3dcfSJed Brown @*/ 561d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs) 562d71ae5a4SJacob Faibussowitsch { 5636f3c3dcfSJed Brown struct DMCompositeLink *next; 5646f3c3dcfSJed Brown PetscInt i; 5656f3c3dcfSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 56671b14b3eSStefano Zampini PetscBool flg; 5676f3c3dcfSJed Brown 5686f3c3dcfSJed Brown PetscFunctionBegin; 5696f3c3dcfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5706f3c3dcfSJed Brown PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 5719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 5727a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 57348a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 5746f3c3dcfSJed Brown 5756f3c3dcfSJed Brown /* loop over packed objects, handling one at at time */ 5766f3c3dcfSJed Brown for (i = 0, next = com->next; next; next = next->next, i++) { 5776f3c3dcfSJed Brown if (lvecs[i]) { 5786f3c3dcfSJed Brown Vec global; 579c5d31e75SLisandro Dalcin const PetscScalar *array; 5806f3c3dcfSJed Brown PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3); 5819566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 5829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 5839566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart)); 5849566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i])); 5859566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i])); 5869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 5879566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 5889566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 5896f3c3dcfSJed Brown } 5906f3c3dcfSJed Brown } 5913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5926f3c3dcfSJed Brown } 5936f3c3dcfSJed Brown 59447c6ae99SBarry Smith /*@C 59547c6ae99SBarry Smith DMCompositeGather - Gathers into a global packed vector from its individual local vectors 59647c6ae99SBarry Smith 59720f4b53cSBarry Smith Collective 59847c6ae99SBarry Smith 599d8d19677SJose E. Roman Input Parameters: 600dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 601dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES` 602*a4e35b19SJacob Faibussowitsch . gvec - the global vector 603*a4e35b19SJacob Faibussowitsch - ... - the individual sequential vectors, `NULL` for any that are not needed 60447c6ae99SBarry Smith 60547c6ae99SBarry Smith Level: advanced 60647c6ae99SBarry Smith 60760225df5SJacob Faibussowitsch Fortran Notes: 608dce8aebaSBarry Smith Fortran users should use `DMCompositeGatherArray()` 609f5f57ec0SBarry Smith 610dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 611db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 612db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 61347c6ae99SBarry Smith @*/ 614d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...) 615d71ae5a4SJacob Faibussowitsch { 61647c6ae99SBarry Smith va_list Argp; 61747c6ae99SBarry Smith struct DMCompositeLink *next; 61847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 619c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt; 62071b14b3eSStefano Zampini PetscBool flg; 62147c6ae99SBarry Smith 62247c6ae99SBarry Smith PetscFunctionBegin; 62347c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 624064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3); 6259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 6267a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 62748a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 62847c6ae99SBarry Smith 62947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 6301dac896bSSatish Balay va_start(Argp, gvec); 6318fd8f222SJed Brown for (cnt = 3, next = com->next; next; cnt++, next = next->next) { 6329ae5db72SJed Brown Vec local; 6339ae5db72SJed Brown local = va_arg(Argp, Vec); 6349ae5db72SJed Brown if (local) { 63547c6ae99SBarry Smith PetscScalar *array; 6369ae5db72SJed Brown Vec global; 63763a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt)); 6389566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 6399566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 6409566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 6419566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global)); 6429566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global)); 6439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 6449566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 6459566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 64647c6ae99SBarry Smith } 64747c6ae99SBarry Smith } 64847c6ae99SBarry Smith va_end(Argp); 6493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65047c6ae99SBarry Smith } 65147c6ae99SBarry Smith 6526f3c3dcfSJed Brown /*@ 6536f3c3dcfSJed Brown DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors 6546f3c3dcfSJed Brown 65520f4b53cSBarry Smith Collective 6566f3c3dcfSJed Brown 657d8d19677SJose E. Roman Input Parameters: 658dce8aebaSBarry Smith + dm - the `DMCOMPOSITE` object 6596f3c3dcfSJed Brown . gvec - the global vector 660dce8aebaSBarry Smith . imode - `INSERT_VALUES` or `ADD_VALUES` 6616f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed 6626f3c3dcfSJed Brown 6636f3c3dcfSJed Brown Level: advanced 6646f3c3dcfSJed Brown 665dce8aebaSBarry Smith Note: 666dce8aebaSBarry Smith This is a non-variadic alternative to `DMCompositeGather()`. 6676f3c3dcfSJed Brown 668dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 669db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 670db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`, 6716f3c3dcfSJed Brown @*/ 672d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs) 673d71ae5a4SJacob Faibussowitsch { 6746f3c3dcfSJed Brown struct DMCompositeLink *next; 6756f3c3dcfSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 6766f3c3dcfSJed Brown PetscInt i; 67771b14b3eSStefano Zampini PetscBool flg; 6786f3c3dcfSJed Brown 6796f3c3dcfSJed Brown PetscFunctionBegin; 6806f3c3dcfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 681064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3); 6829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 6837a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 68448a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 6856f3c3dcfSJed Brown 6866f3c3dcfSJed Brown /* loop over packed objects, handling one at at time */ 6876f3c3dcfSJed Brown for (next = com->next, i = 0; next; next = next->next, i++) { 6886f3c3dcfSJed Brown if (lvecs[i]) { 6896f3c3dcfSJed Brown PetscScalar *array; 6906f3c3dcfSJed Brown Vec global; 691064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4); 6929566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 6939566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 6949566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 6959566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global)); 6969566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global)); 6979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 6989566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 6999566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 7006f3c3dcfSJed Brown } 7016f3c3dcfSJed Brown } 7023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7036f3c3dcfSJed Brown } 7046f3c3dcfSJed Brown 705f5f57ec0SBarry Smith /*@ 706dce8aebaSBarry Smith DMCompositeAddDM - adds a `DM` vector to a `DMCOMPOSITE` 70747c6ae99SBarry Smith 70820f4b53cSBarry Smith Collective 70947c6ae99SBarry Smith 710d8d19677SJose E. Roman Input Parameters: 711dce8aebaSBarry Smith + dmc - the `DMCOMPOSITE` object 712dce8aebaSBarry Smith - dm - the `DM` object 71347c6ae99SBarry Smith 71447c6ae99SBarry Smith Level: advanced 71547c6ae99SBarry Smith 71642747ad1SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeGather()`, `DMCreateGlobalVector()`, 717db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 718db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 71947c6ae99SBarry Smith @*/ 720d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeAddDM(DM dmc, DM dm) 721d71ae5a4SJacob Faibussowitsch { 72206ebdd98SJed Brown PetscInt n, nlocal; 72347c6ae99SBarry Smith struct DMCompositeLink *mine, *next; 72406ebdd98SJed Brown Vec global, local; 72547c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dmc->data; 72671b14b3eSStefano Zampini PetscBool flg; 72747c6ae99SBarry Smith 72847c6ae99SBarry Smith PetscFunctionBegin; 72947c6ae99SBarry Smith PetscValidHeaderSpecific(dmc, DM_CLASSID, 1); 73047c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 7319566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg)); 7327a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 73347c6ae99SBarry Smith next = com->next; 73428b400f6SJacob Faibussowitsch PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite"); 73547c6ae99SBarry Smith 73647c6ae99SBarry Smith /* create new link */ 7379566063dSJacob Faibussowitsch PetscCall(PetscNew(&mine)); 7389566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 7399566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &global)); 7409566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global, &n)); 7419566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &global)); 7429566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &local)); 7439566063dSJacob Faibussowitsch PetscCall(VecGetSize(local, &nlocal)); 7449566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &local)); 7458865f1eaSKarl Rupp 74647c6ae99SBarry Smith mine->n = n; 74706ebdd98SJed Brown mine->nlocal = nlocal; 74847c6ae99SBarry Smith mine->dm = dm; 7490298fd71SBarry Smith mine->next = NULL; 75047c6ae99SBarry Smith com->n += n; 7517ac2b803SAlex Fikl com->nghost += nlocal; 75247c6ae99SBarry Smith 75347c6ae99SBarry Smith /* add to end of list */ 7548865f1eaSKarl Rupp if (!next) com->next = mine; 7558865f1eaSKarl Rupp else { 75647c6ae99SBarry Smith while (next->next) next = next->next; 75747c6ae99SBarry Smith next->next = mine; 75847c6ae99SBarry Smith } 75947c6ae99SBarry Smith com->nDM++; 76047c6ae99SBarry Smith com->nmine++; 7613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 76247c6ae99SBarry Smith } 76347c6ae99SBarry Smith 7649804daf3SBarry Smith #include <petscdraw.h> 76526887b52SJed Brown PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer); 766d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_DMComposite(Vec gvec, PetscViewer viewer) 767d71ae5a4SJacob Faibussowitsch { 76847c6ae99SBarry Smith DM dm; 76947c6ae99SBarry Smith struct DMCompositeLink *next; 77047c6ae99SBarry Smith PetscBool isdraw; 771cef07954SSatish Balay DM_Composite *com; 77247c6ae99SBarry Smith 77347c6ae99SBarry Smith PetscFunctionBegin; 7749566063dSJacob Faibussowitsch PetscCall(VecGetDM(gvec, &dm)); 7757a8be351SBarry Smith PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite"); 77647c6ae99SBarry Smith com = (DM_Composite *)dm->data; 77747c6ae99SBarry Smith next = com->next; 77847c6ae99SBarry Smith 7799566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 78047c6ae99SBarry Smith if (!isdraw) { 78147c6ae99SBarry Smith /* do I really want to call this? */ 7829566063dSJacob Faibussowitsch PetscCall(VecView_MPI(gvec, viewer)); 78347c6ae99SBarry Smith } else { 78447c6ae99SBarry Smith PetscInt cnt = 0; 78547c6ae99SBarry Smith 78647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 78747c6ae99SBarry Smith while (next) { 78847c6ae99SBarry Smith Vec vec; 789ca4278abSLisandro Dalcin const PetscScalar *array; 79047c6ae99SBarry Smith PetscInt bs; 79147c6ae99SBarry Smith 7929ae5db72SJed Brown /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */ 7939566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &vec)); 7949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 7959566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart)); 7969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 7979566063dSJacob Faibussowitsch PetscCall(VecView(vec, viewer)); 7989566063dSJacob Faibussowitsch PetscCall(VecResetArray(vec)); 7999566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vec, &bs)); 8009566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &vec)); 8019566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer, bs)); 80247c6ae99SBarry Smith cnt += bs; 80347c6ae99SBarry Smith next = next->next; 80447c6ae99SBarry Smith } 8059566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt)); 80647c6ae99SBarry Smith } 8073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80847c6ae99SBarry Smith } 80947c6ae99SBarry Smith 810d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec) 811d71ae5a4SJacob Faibussowitsch { 81247c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 81347c6ae99SBarry Smith 81447c6ae99SBarry Smith PetscFunctionBegin; 81547c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8169566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 8179566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 8189566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec)); 8199566063dSJacob Faibussowitsch PetscCall(VecSetType(*gvec, dm->vectype)); 8209566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*gvec, com->n, com->N)); 8219566063dSJacob Faibussowitsch PetscCall(VecSetDM(*gvec, dm)); 8229566063dSJacob Faibussowitsch PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (void (*)(void))VecView_DMComposite)); 8233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82447c6ae99SBarry Smith } 82547c6ae99SBarry Smith 826d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec) 827d71ae5a4SJacob Faibussowitsch { 82847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 82947c6ae99SBarry Smith 83047c6ae99SBarry Smith PetscFunctionBegin; 83147c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 83247c6ae99SBarry Smith if (!com->setup) { 8339566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 8349566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 83547c6ae99SBarry Smith } 8369566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, lvec)); 8379566063dSJacob Faibussowitsch PetscCall(VecSetType(*lvec, dm->vectype)); 8389566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE)); 8399566063dSJacob Faibussowitsch PetscCall(VecSetDM(*lvec, dm)); 8403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84147c6ae99SBarry Smith } 84247c6ae99SBarry Smith 84347c6ae99SBarry Smith /*@C 844dce8aebaSBarry Smith DMCompositeGetISLocalToGlobalMappings - gets an `ISLocalToGlobalMapping` for each `DM` in the `DMCOMPOSITE`, maps to the composite global space 84547c6ae99SBarry Smith 84620f4b53cSBarry Smith Collective; No Fortran Support 84747c6ae99SBarry Smith 84847c6ae99SBarry Smith Input Parameter: 849dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 85047c6ae99SBarry Smith 8512fe279fdSBarry Smith Output Parameter: 8529ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes 853dce8aebaSBarry Smith all the ghost points that individual ghosted `DMDA` may have. 85447c6ae99SBarry Smith 85547c6ae99SBarry Smith Level: advanced 85647c6ae99SBarry Smith 857dce8aebaSBarry Smith Note: 858dce8aebaSBarry Smith Each entry of ltogs should be destroyed with `ISLocalToGlobalMappingDestroy()`, the ltogs array should be freed with `PetscFree()`. 85947c6ae99SBarry Smith 860dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 861db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 862c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 86347c6ae99SBarry Smith @*/ 864d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping **ltogs) 865d71ae5a4SJacob Faibussowitsch { 86647c6ae99SBarry Smith PetscInt i, *idx, n, cnt; 86747c6ae99SBarry Smith struct DMCompositeLink *next; 86847c6ae99SBarry Smith PetscMPIInt rank; 86947c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 87071b14b3eSStefano Zampini PetscBool flg; 87147c6ae99SBarry Smith 87247c6ae99SBarry Smith PetscFunctionBegin; 87347c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8749566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 8757a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 8769566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 8779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM, ltogs)); 87847c6ae99SBarry Smith next = com->next; 8799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 88047c6ae99SBarry Smith 88147c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 88247c6ae99SBarry Smith cnt = 0; 88347c6ae99SBarry Smith while (next) { 8846eb61c8cSJed Brown ISLocalToGlobalMapping ltog; 8856eb61c8cSJed Brown PetscMPIInt size; 88686994e45SJed Brown const PetscInt *suboff, *indices; 8876eb61c8cSJed Brown Vec global; 88847c6ae99SBarry Smith 8896eb61c8cSJed Brown /* Get sub-DM global indices for each local dof */ 8909566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(next->dm, <og)); 8919566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n)); 8929566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices)); 8939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idx)); 89447c6ae99SBarry Smith 8956eb61c8cSJed Brown /* Get the offsets for the sub-DM global vector */ 8969566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 8979566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(global, &suboff)); 8989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size)); 8996eb61c8cSJed Brown 9006eb61c8cSJed Brown /* Shift the sub-DM definition of the global space to the composite global space */ 9016eb61c8cSJed Brown for (i = 0; i < n; i++) { 90286994e45SJed Brown PetscInt subi = indices[i], lo = 0, hi = size, t; 903e333b035SStefano Zampini /* There's no consensus on what a negative index means, 904e333b035SStefano Zampini except for skipping when setting the values in vectors and matrices */ 9059371c9d4SSatish Balay if (subi < 0) { 9069371c9d4SSatish Balay idx[i] = subi - next->grstarts[rank]; 9079371c9d4SSatish Balay continue; 9089371c9d4SSatish Balay } 9096eb61c8cSJed Brown /* Binary search to find which rank owns subi */ 9106eb61c8cSJed Brown while (hi - lo > 1) { 9116eb61c8cSJed Brown t = lo + (hi - lo) / 2; 9126eb61c8cSJed Brown if (suboff[t] > subi) hi = t; 9136eb61c8cSJed Brown else lo = t; 9146eb61c8cSJed Brown } 9156eb61c8cSJed Brown idx[i] = subi - suboff[lo] + next->grstarts[lo]; 9166eb61c8cSJed Brown } 9179566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices)); 9189566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt])); 9199566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 92047c6ae99SBarry Smith next = next->next; 92147c6ae99SBarry Smith cnt++; 92247c6ae99SBarry Smith } 9233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 92447c6ae99SBarry Smith } 92547c6ae99SBarry Smith 92687c85e80SJed Brown /*@C 9279ae5db72SJed Brown DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector 92887c85e80SJed Brown 92920f4b53cSBarry Smith Not Collective; No Fortran Support 93087c85e80SJed Brown 9314165533cSJose E. Roman Input Parameter: 932dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` 93387c85e80SJed Brown 9344165533cSJose E. Roman Output Parameter: 935dce8aebaSBarry Smith . is - array of serial index sets for each each component of the `DMCOMPOSITE` 93687c85e80SJed Brown 93787c85e80SJed Brown Level: intermediate 93887c85e80SJed Brown 93987c85e80SJed Brown Notes: 94087c85e80SJed Brown At present, a composite local vector does not normally exist. This function is used to provide index sets for 941dce8aebaSBarry Smith `MatGetLocalSubMatrix()`. In the future, the scatters for each entry in the `DMCOMPOSITE` may be be merged into a single 9429ae5db72SJed Brown scatter to a composite local vector. The user should not typically need to know which is being done. 94387c85e80SJed Brown 944dce8aebaSBarry Smith To get the composite global indices at all local points (including ghosts), use `DMCompositeGetISLocalToGlobalMappings()`. 94587c85e80SJed Brown 946dce8aebaSBarry Smith To get index sets for pieces of the composite global vector, use `DMCompositeGetGlobalISs()`. 94787c85e80SJed Brown 948dce8aebaSBarry Smith Each returned `IS` should be destroyed with `ISDestroy()`, the array should be freed with `PetscFree()`. 94987c85e80SJed Brown 950dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`, `MatCreateLocalRef()` 95187c85e80SJed Brown @*/ 952d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalISs(DM dm, IS **is) 953d71ae5a4SJacob Faibussowitsch { 95487c85e80SJed Brown DM_Composite *com = (DM_Composite *)dm->data; 95587c85e80SJed Brown struct DMCompositeLink *link; 95687c85e80SJed Brown PetscInt cnt, start; 95771b14b3eSStefano Zampini PetscBool flg; 95887c85e80SJed Brown 95987c85e80SJed Brown PetscFunctionBegin; 96087c85e80SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 9614f572ea9SToby Isaac PetscAssertPointer(is, 2); 9629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 9637a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 9649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nmine, is)); 96506ebdd98SJed Brown for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) { 966520db06cSJed Brown PetscInt bs; 9679566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt])); 9689566063dSJacob Faibussowitsch PetscCall(DMGetBlockSize(link->dm, &bs)); 9699566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize((*is)[cnt], bs)); 970520db06cSJed Brown } 9713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97287c85e80SJed Brown } 97387c85e80SJed Brown 97447c6ae99SBarry Smith /*@C 975dce8aebaSBarry Smith DMCompositeGetGlobalISs - Gets the index sets for each composed object in a `DMCOMPOSITE` 97647c6ae99SBarry Smith 97720f4b53cSBarry Smith Collective 97847c6ae99SBarry Smith 97947c6ae99SBarry Smith Input Parameter: 980dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 98147c6ae99SBarry Smith 9822fe279fdSBarry Smith Output Parameter: 98347c6ae99SBarry Smith . is - the array of index sets 98447c6ae99SBarry Smith 98547c6ae99SBarry Smith Level: advanced 98647c6ae99SBarry Smith 98747c6ae99SBarry Smith Notes: 988dce8aebaSBarry Smith The is entries should be destroyed with `ISDestroy()`, the is array should be freed with `PetscFree()` 98947c6ae99SBarry Smith 99047c6ae99SBarry Smith These could be used to extract a subset of vector entries for a "multi-physics" preconditioner 99147c6ae99SBarry Smith 992dce8aebaSBarry Smith Use `DMCompositeGetLocalISs()` for index sets in the packed local numbering, and 993dce8aebaSBarry Smith `DMCompositeGetISLocalToGlobalMappings()` for to map local sub-`DM` (including ghost) indices to packed global 9946eb61c8cSJed Brown indices. 99547c6ae99SBarry Smith 99660225df5SJacob Faibussowitsch Fortran Notes: 997dce8aebaSBarry Smith The output argument 'is' must be an allocated array of sufficient length, which can be learned using `DMCompositeGetNumberDM()`. 998f3cb0f7eSJed Brown 999dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1000db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 1001c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 100247c6ae99SBarry Smith @*/ 1003d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[]) 1004d71ae5a4SJacob Faibussowitsch { 100566bb578eSMark Adams PetscInt cnt = 0; 100647c6ae99SBarry Smith struct DMCompositeLink *next; 100747c6ae99SBarry Smith PetscMPIInt rank; 100847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 100971b14b3eSStefano Zampini PetscBool flg; 101047c6ae99SBarry Smith 101147c6ae99SBarry Smith PetscFunctionBegin; 101247c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10139566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 10147a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 10157a8be351SBarry Smith PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before"); 10169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM, is)); 101747c6ae99SBarry Smith next = com->next; 10189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 101947c6ae99SBarry Smith 102047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 102147c6ae99SBarry Smith while (next) { 1022e5e52638SMatthew G. Knepley PetscDS prob; 1023e5e52638SMatthew G. Knepley 10249566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt])); 10259566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 1026e5e52638SMatthew G. Knepley if (prob) { 102765c226d8SMatthew G. Knepley MatNullSpace space; 102865c226d8SMatthew G. Knepley Mat pmat; 10290f21e855SMatthew G. Knepley PetscObject disc; 10300f21e855SMatthew G. Knepley PetscInt Nf; 103165c226d8SMatthew G. Knepley 10329566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 1033f24dd8d2SMatthew G. Knepley if (cnt < Nf) { 10349566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(prob, cnt, &disc)); 10359566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space)); 10369566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space)); 10379566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space)); 10389566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space)); 10399566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat)); 10409566063dSJacob Faibussowitsch if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat)); 104165c226d8SMatthew G. Knepley } 1042f24dd8d2SMatthew G. Knepley } 104347c6ae99SBarry Smith cnt++; 104447c6ae99SBarry Smith next = next->next; 104547c6ae99SBarry Smith } 10463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 104747c6ae99SBarry Smith } 104847c6ae99SBarry Smith 1049d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) 1050d71ae5a4SJacob Faibussowitsch { 10514d343eeaSMatthew G Knepley PetscInt nDM; 10524d343eeaSMatthew G Knepley DM *dms; 10534d343eeaSMatthew G Knepley PetscInt i; 10544d343eeaSMatthew G Knepley 10554d343eeaSMatthew G Knepley PetscFunctionBegin; 10569566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 10578865f1eaSKarl Rupp if (numFields) *numFields = nDM; 10589566063dSJacob Faibussowitsch PetscCall(DMCompositeGetGlobalISs(dm, fields)); 10594d343eeaSMatthew G Knepley if (fieldNames) { 10609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, &dms)); 10619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, fieldNames)); 10629566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms)); 10634d343eeaSMatthew G Knepley for (i = 0; i < nDM; i++) { 10644d343eeaSMatthew G Knepley char buf[256]; 10654d343eeaSMatthew G Knepley const char *splitname; 10664d343eeaSMatthew G Knepley 10674d343eeaSMatthew G Knepley /* Split naming precedence: object name, prefix, number */ 10684d343eeaSMatthew G Knepley splitname = ((PetscObject)dm)->name; 10694d343eeaSMatthew G Knepley if (!splitname) { 10709566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname)); 10714d343eeaSMatthew G Knepley if (splitname) { 10724d343eeaSMatthew G Knepley size_t len; 10739566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(buf, splitname, sizeof(buf))); 10748caf3d72SBarry Smith buf[sizeof(buf) - 1] = 0; 10759566063dSJacob Faibussowitsch PetscCall(PetscStrlen(buf, &len)); 10764d343eeaSMatthew G Knepley if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */ 10774d343eeaSMatthew G Knepley splitname = buf; 10784d343eeaSMatthew G Knepley } 10794d343eeaSMatthew G Knepley } 10804d343eeaSMatthew G Knepley if (!splitname) { 108163a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i)); 10824d343eeaSMatthew G Knepley splitname = buf; 10834d343eeaSMatthew G Knepley } 10849566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i])); 10854d343eeaSMatthew G Knepley } 10869566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 10874d343eeaSMatthew G Knepley } 10883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10894d343eeaSMatthew G Knepley } 10904d343eeaSMatthew G Knepley 1091e7c4fc90SDmitry Karpeev /* 1092e7c4fc90SDmitry Karpeev This could take over from DMCreateFieldIS(), as it is more general, 10930298fd71SBarry Smith making DMCreateFieldIS() a special case -- calling with dmlist == NULL; 1094e7c4fc90SDmitry Karpeev At this point it's probably best to be less intrusive, however. 1095e7c4fc90SDmitry Karpeev */ 1096d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 1097d71ae5a4SJacob Faibussowitsch { 1098e7c4fc90SDmitry Karpeev PetscInt nDM; 1099e7c4fc90SDmitry Karpeev PetscInt i; 1100e7c4fc90SDmitry Karpeev 1101e7c4fc90SDmitry Karpeev PetscFunctionBegin; 11029566063dSJacob Faibussowitsch PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist)); 1103e7c4fc90SDmitry Karpeev if (dmlist) { 11049566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 11059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, dmlist)); 11069566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, *dmlist)); 110748a46eb9SPierre Jolivet for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i]))); 1108e7c4fc90SDmitry Karpeev } 11093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1110e7c4fc90SDmitry Karpeev } 1111e7c4fc90SDmitry Karpeev 111247c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 111347c6ae99SBarry Smith /*@C 1114dce8aebaSBarry Smith DMCompositeGetLocalVectors - Gets local vectors for each part of a `DMCOMPOSITE` 1115dce8aebaSBarry Smith Use `DMCompositeRestoreLocalVectors()` to return them. 111647c6ae99SBarry Smith 111720f4b53cSBarry Smith Not Collective; No Fortran Support 111847c6ae99SBarry Smith 111947c6ae99SBarry Smith Input Parameter: 1120dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 112147c6ae99SBarry Smith 112247c6ae99SBarry Smith Output Parameter: 1123*a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec`s 112447c6ae99SBarry Smith 112547c6ae99SBarry Smith Level: advanced 112647c6ae99SBarry Smith 1127dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1128db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1129db781477SPatrick Sanan `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 113047c6ae99SBarry Smith @*/ 1131d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...) 1132d71ae5a4SJacob Faibussowitsch { 113347c6ae99SBarry Smith va_list Argp; 113447c6ae99SBarry Smith struct DMCompositeLink *next; 113547c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 113671b14b3eSStefano Zampini PetscBool flg; 113747c6ae99SBarry Smith 113847c6ae99SBarry Smith PetscFunctionBegin; 113947c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11409566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 11417a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 114247c6ae99SBarry Smith next = com->next; 114347c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 114447c6ae99SBarry Smith va_start(Argp, dm); 114547c6ae99SBarry Smith while (next) { 114647c6ae99SBarry Smith Vec *vec; 114747c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 11489566063dSJacob Faibussowitsch if (vec) PetscCall(DMGetLocalVector(next->dm, vec)); 114947c6ae99SBarry Smith next = next->next; 115047c6ae99SBarry Smith } 115147c6ae99SBarry Smith va_end(Argp); 11523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115347c6ae99SBarry Smith } 115447c6ae99SBarry Smith 115547c6ae99SBarry Smith /*@C 1156dce8aebaSBarry Smith DMCompositeRestoreLocalVectors - Restores local vectors for each part of a `DMCOMPOSITE` 115747c6ae99SBarry Smith 115820f4b53cSBarry Smith Not Collective; No Fortran Support 115947c6ae99SBarry Smith 116047c6ae99SBarry Smith Input Parameter: 1161dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 116247c6ae99SBarry Smith 116347c6ae99SBarry Smith Output Parameter: 1164*a4e35b19SJacob Faibussowitsch . ... - the individual sequential `Vec` 116547c6ae99SBarry Smith 116647c6ae99SBarry Smith Level: advanced 116747c6ae99SBarry Smith 1168dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1169db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1170db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 117147c6ae99SBarry Smith @*/ 1172d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...) 1173d71ae5a4SJacob Faibussowitsch { 117447c6ae99SBarry Smith va_list Argp; 117547c6ae99SBarry Smith struct DMCompositeLink *next; 117647c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 117771b14b3eSStefano Zampini PetscBool flg; 117847c6ae99SBarry Smith 117947c6ae99SBarry Smith PetscFunctionBegin; 118047c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 11827a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 118347c6ae99SBarry Smith next = com->next; 118447c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 118547c6ae99SBarry Smith va_start(Argp, dm); 118647c6ae99SBarry Smith while (next) { 118747c6ae99SBarry Smith Vec *vec; 118847c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 11899566063dSJacob Faibussowitsch if (vec) PetscCall(DMRestoreLocalVector(next->dm, vec)); 119047c6ae99SBarry Smith next = next->next; 119147c6ae99SBarry Smith } 119247c6ae99SBarry Smith va_end(Argp); 11933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119447c6ae99SBarry Smith } 119547c6ae99SBarry Smith 119647c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 119747c6ae99SBarry Smith /*@C 1198dce8aebaSBarry Smith DMCompositeGetEntries - Gets the `DM` for each entry in a `DMCOMPOSITE`. 119947c6ae99SBarry Smith 120047c6ae99SBarry Smith Not Collective 120147c6ae99SBarry Smith 120247c6ae99SBarry Smith Input Parameter: 1203dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 120447c6ae99SBarry Smith 120547c6ae99SBarry Smith Output Parameter: 1206*a4e35b19SJacob Faibussowitsch . ... - the individual entries `DM` 120747c6ae99SBarry Smith 120847c6ae99SBarry Smith Level: advanced 120947c6ae99SBarry Smith 121060225df5SJacob Faibussowitsch Fortran Notes: 1211dce8aebaSBarry Smith Available as `DMCompositeGetEntries()` for one output `DM`, DMCompositeGetEntries2() for 2, etc 1212f5f57ec0SBarry Smith 1213dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()` 1214db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 121560225df5SJacob Faibussowitsch `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()` 121647c6ae99SBarry Smith @*/ 1217d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntries(DM dm, ...) 1218d71ae5a4SJacob Faibussowitsch { 121947c6ae99SBarry Smith va_list Argp; 122047c6ae99SBarry Smith struct DMCompositeLink *next; 122147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 122271b14b3eSStefano Zampini PetscBool flg; 122347c6ae99SBarry Smith 122447c6ae99SBarry Smith PetscFunctionBegin; 122547c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 12277a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 122847c6ae99SBarry Smith next = com->next; 122947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 123047c6ae99SBarry Smith va_start(Argp, dm); 123147c6ae99SBarry Smith while (next) { 123247c6ae99SBarry Smith DM *dmn; 123347c6ae99SBarry Smith dmn = va_arg(Argp, DM *); 12349ae5db72SJed Brown if (dmn) *dmn = next->dm; 123547c6ae99SBarry Smith next = next->next; 123647c6ae99SBarry Smith } 123747c6ae99SBarry Smith va_end(Argp); 12383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123947c6ae99SBarry Smith } 124047c6ae99SBarry Smith 1241dbab29e1SMark F. Adams /*@C 1242dce8aebaSBarry Smith DMCompositeGetEntriesArray - Gets the DM for each entry in a `DMCOMPOSITE` 12432fa5ba8aSJed Brown 12442fa5ba8aSJed Brown Not Collective 12452fa5ba8aSJed Brown 12462fa5ba8aSJed Brown Input Parameter: 1247dce8aebaSBarry Smith . dm - the `DMCOMPOSITE` object 1248907376e6SBarry Smith 1249907376e6SBarry Smith Output Parameter: 1250dce8aebaSBarry Smith . dms - array of sufficient length (see `DMCompositeGetNumberDM()`) to hold the individual `DM` 12512fa5ba8aSJed Brown 12522fa5ba8aSJed Brown Level: advanced 12532fa5ba8aSJed Brown 1254dce8aebaSBarry Smith .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()` 1255db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 125660225df5SJacob Faibussowitsch `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()` 12572fa5ba8aSJed Brown @*/ 1258d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntriesArray(DM dm, DM dms[]) 1259d71ae5a4SJacob Faibussowitsch { 12602fa5ba8aSJed Brown struct DMCompositeLink *next; 12612fa5ba8aSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 12622fa5ba8aSJed Brown PetscInt i; 126371b14b3eSStefano Zampini PetscBool flg; 12642fa5ba8aSJed Brown 12652fa5ba8aSJed Brown PetscFunctionBegin; 12662fa5ba8aSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12679566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 12687a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 12692fa5ba8aSJed Brown /* loop over packed objects, handling one at at time */ 12702fa5ba8aSJed Brown for (next = com->next, i = 0; next; next = next->next, i++) dms[i] = next->dm; 12713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12722fa5ba8aSJed Brown } 12732fa5ba8aSJed Brown 1274e10fd815SStefano Zampini typedef struct { 1275e10fd815SStefano Zampini DM dm; 1276e10fd815SStefano Zampini PetscViewer *subv; 1277e10fd815SStefano Zampini Vec *vecs; 1278e10fd815SStefano Zampini } GLVisViewerCtx; 1279e10fd815SStefano Zampini 1280d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 1281d71ae5a4SJacob Faibussowitsch { 1282e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 1283e10fd815SStefano Zampini PetscInt i, n; 1284e10fd815SStefano Zampini 1285e10fd815SStefano Zampini PetscFunctionBegin; 12869566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm, &n)); 128748a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(PetscViewerDestroy(&ctx->subv[i])); 12889566063dSJacob Faibussowitsch PetscCall(PetscFree2(ctx->subv, ctx->vecs)); 12899566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx->dm)); 12909566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 12913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1292e10fd815SStefano Zampini } 1293e10fd815SStefano Zampini 1294d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 1295d71ae5a4SJacob Faibussowitsch { 1296e10fd815SStefano Zampini Vec X = (Vec)oX; 1297e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 1298e10fd815SStefano Zampini PetscInt i, n, cumf; 1299e10fd815SStefano Zampini 1300e10fd815SStefano Zampini PetscFunctionBegin; 13019566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm, &n)); 13029566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(ctx->dm, X, n, NULL, ctx->vecs)); 1303e10fd815SStefano Zampini for (i = 0, cumf = 0; i < n; i++) { 1304e10fd815SStefano Zampini PetscErrorCode (*g2l)(PetscObject, PetscInt, PetscObject[], void *); 1305e10fd815SStefano Zampini void *fctx; 1306e10fd815SStefano Zampini PetscInt nfi; 1307e10fd815SStefano Zampini 13089566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nfi, NULL, NULL, &g2l, NULL, &fctx)); 1309e10fd815SStefano Zampini if (!nfi) continue; 13101baa6e33SBarry Smith if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i], nfi, oXfield + cumf, fctx)); 13111baa6e33SBarry Smith else PetscCall(VecCopy(ctx->vecs[i], (Vec)(oXfield[cumf]))); 1312e10fd815SStefano Zampini cumf += nfi; 1313e10fd815SStefano Zampini } 13149566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(ctx->dm, X, n, NULL, ctx->vecs)); 13153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1316e10fd815SStefano Zampini } 1317e10fd815SStefano Zampini 1318d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer) 1319d71ae5a4SJacob Faibussowitsch { 1320e10fd815SStefano Zampini DM dm = (DM)odm, *dms; 1321e10fd815SStefano Zampini Vec *Ufds; 1322e10fd815SStefano Zampini GLVisViewerCtx *ctx; 1323e10fd815SStefano Zampini PetscInt i, n, tnf, *sdim; 1324e10fd815SStefano Zampini char **fecs; 1325e10fd815SStefano Zampini 1326e10fd815SStefano Zampini PetscFunctionBegin; 13279566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 13289566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 1329e10fd815SStefano Zampini ctx->dm = dm; 13309566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &n)); 13319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dms)); 13329566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms)); 13339566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &ctx->subv, n, &ctx->vecs)); 1334e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1335e10fd815SStefano Zampini PetscInt nf; 1336e10fd815SStefano Zampini 13379566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PetscObjectComm(odm), &ctx->subv[i])); 13389566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(ctx->subv[i], PETSCVIEWERGLVIS)); 13399566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetDM_Private(ctx->subv[i], (PetscObject)dms[i])); 13409566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nf, NULL, NULL, NULL, NULL, NULL)); 1341e10fd815SStefano Zampini tnf += nf; 1342e10fd815SStefano Zampini } 13439566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 13449566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(tnf, &fecs, tnf, &sdim, tnf, &Ufds)); 1345e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1346e10fd815SStefano Zampini PetscInt *sd, nf, f; 1347e10fd815SStefano Zampini const char **fec; 1348e10fd815SStefano Zampini Vec *Uf; 1349e10fd815SStefano Zampini 13509566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nf, &fec, &sd, NULL, (PetscObject **)&Uf, NULL)); 1351e10fd815SStefano Zampini for (f = 0; f < nf; f++) { 13529566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec[f], &fecs[tnf + f])); 1353e10fd815SStefano Zampini Ufds[tnf + f] = Uf[f]; 1354e10fd815SStefano Zampini sdim[tnf + f] = sd[f]; 1355e10fd815SStefano Zampini } 1356e10fd815SStefano Zampini tnf += nf; 1357e10fd815SStefano Zampini } 13589566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetFields(viewer, tnf, (const char **)fecs, sdim, DMCompositeSampleGLVisFields_Private, (PetscObject *)Ufds, ctx, DestroyGLVisViewerCtx_Private)); 135948a46eb9SPierre Jolivet for (i = 0; i < tnf; i++) PetscCall(PetscFree(fecs[i])); 13609566063dSJacob Faibussowitsch PetscCall(PetscFree3(fecs, sdim, Ufds)); 13613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1362e10fd815SStefano Zampini } 1363e10fd815SStefano Zampini 1364d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_Composite(DM dmi, MPI_Comm comm, DM *fine) 1365d71ae5a4SJacob Faibussowitsch { 136647c6ae99SBarry Smith struct DMCompositeLink *next; 136747c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dmi->data; 136847c6ae99SBarry Smith DM dm; 136947c6ae99SBarry Smith 137047c6ae99SBarry Smith PetscFunctionBegin; 137147c6ae99SBarry Smith PetscValidHeaderSpecific(dmi, DM_CLASSID, 1); 137248a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm)); 13739566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi)); 137447c6ae99SBarry Smith next = com->next; 13759566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm, fine)); 137647c6ae99SBarry Smith 137747c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 137847c6ae99SBarry Smith while (next) { 13799566063dSJacob Faibussowitsch PetscCall(DMRefine(next->dm, comm, &dm)); 13809566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine, dm)); 13819566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 138247c6ae99SBarry Smith next = next->next; 138347c6ae99SBarry Smith } 13843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 138547c6ae99SBarry Smith } 138647c6ae99SBarry Smith 1387d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsen_Composite(DM dmi, MPI_Comm comm, DM *fine) 1388d71ae5a4SJacob Faibussowitsch { 138914354c39SJed Brown struct DMCompositeLink *next; 139014354c39SJed Brown DM_Composite *com = (DM_Composite *)dmi->data; 139114354c39SJed Brown DM dm; 139214354c39SJed Brown 139314354c39SJed Brown PetscFunctionBegin; 139414354c39SJed Brown PetscValidHeaderSpecific(dmi, DM_CLASSID, 1); 13959566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi)); 139648a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm)); 139714354c39SJed Brown next = com->next; 13989566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm, fine)); 139914354c39SJed Brown 140014354c39SJed Brown /* loop over packed objects, handling one at at time */ 140114354c39SJed Brown while (next) { 14029566063dSJacob Faibussowitsch PetscCall(DMCoarsen(next->dm, comm, &dm)); 14039566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine, dm)); 14049566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 140514354c39SJed Brown next = next->next; 140614354c39SJed Brown } 14073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 140814354c39SJed Brown } 140947c6ae99SBarry Smith 1410d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInterpolation_Composite(DM coarse, DM fine, Mat *A, Vec *v) 1411d71ae5a4SJacob Faibussowitsch { 14129ae5db72SJed Brown PetscInt m, n, M, N, nDM, i; 141347c6ae99SBarry Smith struct DMCompositeLink *nextc; 141447c6ae99SBarry Smith struct DMCompositeLink *nextf; 141525296bd5SBarry Smith Vec gcoarse, gfine, *vecs; 141647c6ae99SBarry Smith DM_Composite *comcoarse = (DM_Composite *)coarse->data; 141747c6ae99SBarry Smith DM_Composite *comfine = (DM_Composite *)fine->data; 14189ae5db72SJed Brown Mat *mats; 141947c6ae99SBarry Smith 142047c6ae99SBarry Smith PetscFunctionBegin; 142147c6ae99SBarry Smith PetscValidHeaderSpecific(coarse, DM_CLASSID, 1); 142247c6ae99SBarry Smith PetscValidHeaderSpecific(fine, DM_CLASSID, 2); 14239566063dSJacob Faibussowitsch PetscCall(DMSetUp(coarse)); 14249566063dSJacob Faibussowitsch PetscCall(DMSetUp(fine)); 142547c6ae99SBarry Smith /* use global vectors only for determining matrix layout */ 14269566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(coarse, &gcoarse)); 14279566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fine, &gfine)); 14289566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gcoarse, &n)); 14299566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gfine, &m)); 14309566063dSJacob Faibussowitsch PetscCall(VecGetSize(gcoarse, &N)); 14319566063dSJacob Faibussowitsch PetscCall(VecGetSize(gfine, &M)); 14329566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(coarse, &gcoarse)); 14339566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fine, &gfine)); 143447c6ae99SBarry Smith 14359ae5db72SJed Brown nDM = comfine->nDM; 143663a3b9bcSJacob 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); 14379566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nDM * nDM, &mats)); 143848a46eb9SPierre Jolivet if (v) PetscCall(PetscCalloc1(nDM, &vecs)); 143947c6ae99SBarry Smith 144047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 14419ae5db72SJed Brown for (nextc = comcoarse->next, nextf = comfine->next, i = 0; nextc; nextc = nextc->next, nextf = nextf->next, i++) { 14421baa6e33SBarry Smith if (!v) PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], NULL)); 14431baa6e33SBarry Smith else PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], &vecs[i])); 144447c6ae99SBarry Smith } 14459566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, nDM, NULL, mats, A)); 14461baa6e33SBarry Smith if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, vecs, v)); 14479566063dSJacob Faibussowitsch for (i = 0; i < nDM * nDM; i++) PetscCall(MatDestroy(&mats[i])); 14489566063dSJacob Faibussowitsch PetscCall(PetscFree(mats)); 144925296bd5SBarry Smith if (v) { 14509566063dSJacob Faibussowitsch for (i = 0; i < nDM; i++) PetscCall(VecDestroy(&vecs[i])); 14519566063dSJacob Faibussowitsch PetscCall(PetscFree(vecs)); 145225296bd5SBarry Smith } 14533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 145447c6ae99SBarry Smith } 145547c6ae99SBarry Smith 1456d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm) 1457d71ae5a4SJacob Faibussowitsch { 14581411c6eeSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 14591411c6eeSJed Brown ISLocalToGlobalMapping *ltogs; 1460f7efa3c7SJed Brown PetscInt i; 14611411c6eeSJed Brown 14621411c6eeSJed Brown PetscFunctionBegin; 14631411c6eeSJed Brown /* Set the ISLocalToGlobalMapping on the new matrix */ 14649566063dSJacob Faibussowitsch PetscCall(DMCompositeGetISLocalToGlobalMappings(dm, <ogs)); 14659566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm), com->nDM, ltogs, &dm->ltogmap)); 14669566063dSJacob Faibussowitsch for (i = 0; i < com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(<ogs[i])); 14679566063dSJacob Faibussowitsch PetscCall(PetscFree(ltogs)); 14683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14691411c6eeSJed Brown } 14701411c6eeSJed Brown 1471d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateColoring_Composite(DM dm, ISColoringType ctype, ISColoring *coloring) 1472d71ae5a4SJacob Faibussowitsch { 147347c6ae99SBarry Smith PetscInt n, i, cnt; 147447c6ae99SBarry Smith ISColoringValue *colors; 147547c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 147647c6ae99SBarry Smith ISColoringValue maxcol = 0; 147747c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 147847c6ae99SBarry Smith 147947c6ae99SBarry Smith PetscFunctionBegin; 148047c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 148108401ef6SPierre Jolivet PetscCheck(ctype != IS_COLORING_LOCAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only global coloring supported"); 1482f7d195e4SLawrence Mitchell if (ctype == IS_COLORING_GLOBAL) { 148347c6ae99SBarry Smith n = com->n; 1484ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unknown ISColoringType"); 14859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &colors)); /* freed in ISColoringDestroy() */ 148647c6ae99SBarry Smith 14879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL)); 148847c6ae99SBarry Smith if (dense) { 1489ad540459SPierre Jolivet for (i = 0; i < n; i++) colors[i] = (ISColoringValue)(com->rstart + i); 149047c6ae99SBarry Smith maxcol = com->N; 149147c6ae99SBarry Smith } else { 149247c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 149347c6ae99SBarry Smith PetscMPIInt rank; 149447c6ae99SBarry Smith 14959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 149647c6ae99SBarry Smith cnt = 0; 149747c6ae99SBarry Smith while (next) { 149847c6ae99SBarry Smith ISColoring lcoloring; 149947c6ae99SBarry Smith 15009566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(next->dm, IS_COLORING_GLOBAL, &lcoloring)); 1501ad540459SPierre Jolivet for (i = 0; i < lcoloring->N; i++) colors[cnt++] = maxcol + lcoloring->colors[i]; 150247c6ae99SBarry Smith maxcol += lcoloring->n; 15039566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&lcoloring)); 150447c6ae99SBarry Smith next = next->next; 150547c6ae99SBarry Smith } 150647c6ae99SBarry Smith } 15079566063dSJacob Faibussowitsch PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), maxcol, n, colors, PETSC_OWN_POINTER, coloring)); 15083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150947c6ae99SBarry Smith } 151047c6ae99SBarry Smith 1511d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) 1512d71ae5a4SJacob Faibussowitsch { 151347c6ae99SBarry Smith struct DMCompositeLink *next; 151447c6ae99SBarry Smith PetscScalar *garray, *larray; 151547c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 151647c6ae99SBarry Smith 151747c6ae99SBarry Smith PetscFunctionBegin; 151847c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 151947c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 152039d35262SVincent Le Chenadec 152148a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 152239d35262SVincent Le Chenadec 15239566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &garray)); 15249566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec, &larray)); 152547c6ae99SBarry Smith 152647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 152739d35262SVincent Le Chenadec next = com->next; 152847c6ae99SBarry Smith while (next) { 152947c6ae99SBarry Smith Vec local, global; 153047c6ae99SBarry Smith PetscInt N; 153147c6ae99SBarry Smith 15329566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 15339566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global, &N)); 15349566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, garray)); 15359566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local)); 15369566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local, larray)); 15379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, mode, local)); 15389566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, mode, local)); 15399566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 15409566063dSJacob Faibussowitsch PetscCall(VecResetArray(local)); 15419566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 15429566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local)); 154339d35262SVincent Le Chenadec 154406ebdd98SJed Brown larray += next->nlocal; 154539d35262SVincent Le Chenadec garray += next->n; 154647c6ae99SBarry Smith next = next->next; 154747c6ae99SBarry Smith } 154847c6ae99SBarry Smith 15499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, NULL)); 15509566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lvec, NULL)); 15513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 155247c6ae99SBarry Smith } 155347c6ae99SBarry Smith 1554d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) 1555d71ae5a4SJacob Faibussowitsch { 15560c010503SBarry Smith PetscFunctionBegin; 155739d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 155839d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 155939d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 4); 15603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 156139d35262SVincent Le Chenadec } 156239d35262SVincent Le Chenadec 1563d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1564d71ae5a4SJacob Faibussowitsch { 156539d35262SVincent Le Chenadec struct DMCompositeLink *next; 156639d35262SVincent Le Chenadec PetscScalar *larray, *garray; 156739d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite *)dm->data; 156839d35262SVincent Le Chenadec 156939d35262SVincent Le Chenadec PetscFunctionBegin; 157039d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 157139d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 157239d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 157339d35262SVincent Le Chenadec 157448a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 157539d35262SVincent Le Chenadec 15769566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec, &larray)); 15779566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &garray)); 157839d35262SVincent Le Chenadec 157939d35262SVincent Le Chenadec /* loop over packed objects, handling one at at time */ 158039d35262SVincent Le Chenadec next = com->next; 158139d35262SVincent Le Chenadec while (next) { 158239d35262SVincent Le Chenadec Vec global, local; 158339d35262SVincent Le Chenadec 15849566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local)); 15859566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local, larray)); 15869566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 15879566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, garray)); 15889566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, local, mode, global)); 15899566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, local, mode, global)); 15909566063dSJacob Faibussowitsch PetscCall(VecResetArray(local)); 15919566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 15929566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 15939566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local)); 159439d35262SVincent Le Chenadec 159539d35262SVincent Le Chenadec garray += next->n; 159639d35262SVincent Le Chenadec larray += next->nlocal; 159739d35262SVincent Le Chenadec next = next->next; 159839d35262SVincent Le Chenadec } 159939d35262SVincent Le Chenadec 16009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, NULL)); 16019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lvec, NULL)); 160239d35262SVincent Le Chenadec 16033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 160439d35262SVincent Le Chenadec } 160539d35262SVincent Le Chenadec 1606d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1607d71ae5a4SJacob Faibussowitsch { 160839d35262SVincent Le Chenadec PetscFunctionBegin; 160939d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 161039d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 161139d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 16123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 161339d35262SVincent Le Chenadec } 161439d35262SVincent Le Chenadec 1615d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2) 1616d71ae5a4SJacob Faibussowitsch { 161739d35262SVincent Le Chenadec struct DMCompositeLink *next; 161839d35262SVincent Le Chenadec PetscScalar *array1, *array2; 161939d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite *)dm->data; 162039d35262SVincent Le Chenadec 162139d35262SVincent Le Chenadec PetscFunctionBegin; 162239d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 162339d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2); 162439d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4); 162539d35262SVincent Le Chenadec 162648a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 162739d35262SVincent Le Chenadec 16289566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec1, &array1)); 16299566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec2, &array2)); 163039d35262SVincent Le Chenadec 163139d35262SVincent Le Chenadec /* loop over packed objects, handling one at at time */ 163239d35262SVincent Le Chenadec next = com->next; 163339d35262SVincent Le Chenadec while (next) { 163439d35262SVincent Le Chenadec Vec local1, local2; 163539d35262SVincent Le Chenadec 16369566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local1)); 16379566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local1, array1)); 16389566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local2)); 16399566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local2, array2)); 16409566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2)); 16419566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2)); 16429566063dSJacob Faibussowitsch PetscCall(VecResetArray(local2)); 16439566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local2)); 16449566063dSJacob Faibussowitsch PetscCall(VecResetArray(local1)); 16459566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local1)); 164639d35262SVincent Le Chenadec 164739d35262SVincent Le Chenadec array1 += next->nlocal; 164839d35262SVincent Le Chenadec array2 += next->nlocal; 164939d35262SVincent Le Chenadec next = next->next; 165039d35262SVincent Le Chenadec } 165139d35262SVincent Le Chenadec 16529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec1, NULL)); 16539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec2, NULL)); 165439d35262SVincent Le Chenadec 16553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 165639d35262SVincent Le Chenadec } 165739d35262SVincent Le Chenadec 1658d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1659d71ae5a4SJacob Faibussowitsch { 166039d35262SVincent Le Chenadec PetscFunctionBegin; 166139d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 166239d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 166339d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 16643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16650c010503SBarry Smith } 166647c6ae99SBarry Smith 16676ae3a549SBarry Smith /*MC 1668dce8aebaSBarry Smith DMCOMPOSITE = "composite" - A `DM` object that is used to manage data for a collection of `DM` 16696ae3a549SBarry Smith 16706ae3a549SBarry Smith Level: intermediate 16716ae3a549SBarry Smith 1672db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()` 16736ae3a549SBarry Smith M*/ 16746ae3a549SBarry Smith 1675d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p) 1676d71ae5a4SJacob Faibussowitsch { 1677a4121054SBarry Smith DM_Composite *com; 1678a4121054SBarry Smith 1679a4121054SBarry Smith PetscFunctionBegin; 16804dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&com)); 1681a4121054SBarry Smith p->data = com; 1682a4121054SBarry Smith com->n = 0; 16837ac2b803SAlex Fikl com->nghost = 0; 16840298fd71SBarry Smith com->next = NULL; 1685a4121054SBarry Smith com->nDM = 0; 1686a4121054SBarry Smith 1687a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Composite; 1688a4121054SBarry Smith p->ops->createlocalvector = DMCreateLocalVector_Composite; 1689184d77edSJed Brown p->ops->getlocaltoglobalmapping = DMGetLocalToGlobalMapping_Composite; 16904d343eeaSMatthew G Knepley p->ops->createfieldis = DMCreateFieldIS_Composite; 169116621825SDmitry Karpeev p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite; 1692a4121054SBarry Smith p->ops->refine = DMRefine_Composite; 169314354c39SJed Brown p->ops->coarsen = DMCoarsen_Composite; 169425296bd5SBarry Smith p->ops->createinterpolation = DMCreateInterpolation_Composite; 169525296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Composite; 1696e727c939SJed Brown p->ops->getcoloring = DMCreateColoring_Composite; 1697a4121054SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite; 1698a4121054SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite; 169939d35262SVincent Le Chenadec p->ops->localtoglobalbegin = DMLocalToGlobalBegin_Composite; 170039d35262SVincent Le Chenadec p->ops->localtoglobalend = DMLocalToGlobalEnd_Composite; 170139d35262SVincent Le Chenadec p->ops->localtolocalbegin = DMLocalToLocalBegin_Composite; 170239d35262SVincent Le Chenadec p->ops->localtolocalend = DMLocalToLocalEnd_Composite; 1703a4121054SBarry Smith p->ops->destroy = DMDestroy_Composite; 1704a4121054SBarry Smith p->ops->view = DMView_Composite; 1705a4121054SBarry Smith p->ops->setup = DMSetUp_Composite; 1706e10fd815SStefano Zampini 17079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite)); 17083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1709a4121054SBarry Smith } 1710a4121054SBarry Smith 17111fd49c25SBarry Smith /*@ 1712dce8aebaSBarry Smith DMCompositeCreate - Creates a `DMCOMPOSITE`, used to generate "composite" 17130c010503SBarry Smith vectors made up of several subvectors. 17140c010503SBarry Smith 1715d083f849SBarry Smith Collective 171647c6ae99SBarry Smith 171747c6ae99SBarry Smith Input Parameter: 17180c010503SBarry Smith . comm - the processors that will share the global vector 17190c010503SBarry Smith 17202fe279fdSBarry Smith Output Parameter: 1721dce8aebaSBarry Smith . packer - the `DMCOMPOSITE` object 172247c6ae99SBarry Smith 172347c6ae99SBarry Smith Level: advanced 172447c6ae99SBarry Smith 172560225df5SJacob Faibussowitsch .seealso: `DMCOMPOSITE`, `DM`, `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCreate()` 1726db781477SPatrick Sanan `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()` 1727db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 172847c6ae99SBarry Smith @*/ 1729d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer) 1730d71ae5a4SJacob Faibussowitsch { 173147c6ae99SBarry Smith PetscFunctionBegin; 17324f572ea9SToby Isaac PetscAssertPointer(packer, 2); 17339566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, packer)); 17349566063dSJacob Faibussowitsch PetscCall(DMSetType(*packer, DMCOMPOSITE)); 17353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173647c6ae99SBarry Smith } 1737