147c6ae99SBarry Smith 2ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h> /*I "petscdmcomposite.h" I*/ 3af0996ceSBarry Smith #include <petsc/private/isimpl.h> 4e10fd815SStefano Zampini #include <petsc/private/glvisviewerimpl.h> 52764a2aaSMatthew G. Knepley #include <petscds.h> 647c6ae99SBarry Smith 747c6ae99SBarry Smith /*@C 847c6ae99SBarry Smith DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the 9bebe2cf6SSatish Balay separate components (DMs) in a DMto build the correct matrix nonzero structure. 1047c6ae99SBarry Smith 11d083f849SBarry Smith Logically Collective 1247c6ae99SBarry Smith 13d8d19677SJose E. Roman Input Parameters: 1447c6ae99SBarry Smith + dm - the composite object 1547c6ae99SBarry Smith - formcouplelocations - routine to set the nonzero locations in the matrix 1647c6ae99SBarry Smith 1747c6ae99SBarry Smith Level: advanced 1847c6ae99SBarry Smith 19f5f57ec0SBarry Smith Not available from Fortran 20f5f57ec0SBarry Smith 2195452b02SPatrick Sanan Notes: 2295452b02SPatrick Sanan See DMSetApplicationContext() and DMGetApplicationContext() for how to get user information into 2347c6ae99SBarry Smith this routine 2447c6ae99SBarry Smith 2547c6ae99SBarry Smith @*/ 26*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeSetCoupling(DM dm, PetscErrorCode (*FormCoupleLocations)(DM, Mat, PetscInt *, PetscInt *, PetscInt, PetscInt, PetscInt, PetscInt)) 27*d71ae5a4SJacob Faibussowitsch { 2847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 2971b14b3eSStefano Zampini PetscBool flg; 3047c6ae99SBarry Smith 3147c6ae99SBarry Smith PetscFunctionBegin; 329566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 337a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 3447c6ae99SBarry Smith com->FormCoupleLocations = FormCoupleLocations; 3547c6ae99SBarry Smith PetscFunctionReturn(0); 3647c6ae99SBarry Smith } 3747c6ae99SBarry Smith 38*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDestroy_Composite(DM dm) 39*d71ae5a4SJacob Faibussowitsch { 4047c6ae99SBarry Smith struct DMCompositeLink *next, *prev; 4147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 4247c6ae99SBarry Smith 4347c6ae99SBarry Smith PetscFunctionBegin; 4447c6ae99SBarry Smith next = com->next; 4547c6ae99SBarry Smith while (next) { 4647c6ae99SBarry Smith prev = next; 4747c6ae99SBarry Smith next = next->next; 489566063dSJacob Faibussowitsch PetscCall(DMDestroy(&prev->dm)); 499566063dSJacob Faibussowitsch PetscCall(PetscFree(prev->grstarts)); 509566063dSJacob Faibussowitsch PetscCall(PetscFree(prev)); 5147c6ae99SBarry Smith } 529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL)); 53435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 549566063dSJacob Faibussowitsch PetscCall(PetscFree(com)); 5547c6ae99SBarry Smith PetscFunctionReturn(0); 5647c6ae99SBarry Smith } 5747c6ae99SBarry Smith 58*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_Composite(DM dm, PetscViewer v) 59*d71ae5a4SJacob Faibussowitsch { 6047c6ae99SBarry Smith PetscBool iascii; 6147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 6247c6ae99SBarry Smith 6347c6ae99SBarry Smith PetscFunctionBegin; 649566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii)); 6547c6ae99SBarry Smith if (iascii) { 6647c6ae99SBarry Smith struct DMCompositeLink *lnk = com->next; 6747c6ae99SBarry Smith PetscInt i; 6847c6ae99SBarry Smith 699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "DM (%s)\n", ((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix")); 7063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, " contains %" PetscInt_FMT " DMs\n", com->nDM)); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v)); 7247c6ae99SBarry Smith for (i = 0; lnk; lnk = lnk->next, i++) { 7363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "Link %" PetscInt_FMT ": DM of type %s\n", i, ((PetscObject)lnk->dm)->type_name)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v)); 759566063dSJacob Faibussowitsch PetscCall(DMView(lnk->dm, v)); 769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v)); 7747c6ae99SBarry Smith } 789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v)); 7947c6ae99SBarry Smith } 8047c6ae99SBarry Smith PetscFunctionReturn(0); 8147c6ae99SBarry Smith } 8247c6ae99SBarry Smith 8347c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/ 84*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUp_Composite(DM dm) 85*d71ae5a4SJacob Faibussowitsch { 8647c6ae99SBarry Smith PetscInt nprev = 0; 8747c6ae99SBarry Smith PetscMPIInt rank, size; 8847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 8947c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 9047c6ae99SBarry Smith PetscLayout map; 9147c6ae99SBarry Smith 9247c6ae99SBarry Smith PetscFunctionBegin; 9328b400f6SJacob Faibussowitsch PetscCheck(!com->setup, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Packer has already been setup"); 949566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)dm), &map)); 959566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(map, com->n)); 969566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetSize(map, PETSC_DETERMINE)); 979566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(map, 1)); 989566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(map)); 999566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &com->N)); 1009566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(map, &com->rstart, NULL)); 1019566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&map)); 10247c6ae99SBarry Smith 1039ae5db72SJed Brown /* now set the rstart for each linked vector */ 1049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 10647c6ae99SBarry Smith while (next) { 10747c6ae99SBarry Smith next->rstart = nprev; 10806ebdd98SJed Brown nprev += next->n; 10947c6ae99SBarry Smith next->grstart = com->rstart + next->rstart; 1109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &next->grstarts)); 1119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&next->grstart, 1, MPIU_INT, next->grstarts, 1, MPIU_INT, PetscObjectComm((PetscObject)dm))); 11247c6ae99SBarry Smith next = next->next; 11347c6ae99SBarry Smith } 11447c6ae99SBarry Smith com->setup = PETSC_TRUE; 11547c6ae99SBarry Smith PetscFunctionReturn(0); 11647c6ae99SBarry Smith } 11747c6ae99SBarry Smith 11847c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/ 11947c6ae99SBarry Smith 12073e31fe2SJed Brown /*@ 12147c6ae99SBarry Smith DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite 12247c6ae99SBarry Smith representation. 12347c6ae99SBarry Smith 12447c6ae99SBarry Smith Not Collective 12547c6ae99SBarry Smith 12647c6ae99SBarry Smith Input Parameter: 12747c6ae99SBarry Smith . dm - the packer object 12847c6ae99SBarry Smith 12947c6ae99SBarry Smith Output Parameter: 13047c6ae99SBarry Smith . nDM - the number of DMs 13147c6ae99SBarry Smith 13247c6ae99SBarry Smith Level: beginner 13347c6ae99SBarry Smith 13447c6ae99SBarry Smith @*/ 135*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetNumberDM(DM dm, PetscInt *nDM) 136*d71ae5a4SJacob Faibussowitsch { 13747c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 13871b14b3eSStefano Zampini PetscBool flg; 1395fd66863SKarl Rupp 14047c6ae99SBarry Smith PetscFunctionBegin; 14147c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1429566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1437a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 14447c6ae99SBarry Smith *nDM = com->nDM; 14547c6ae99SBarry Smith PetscFunctionReturn(0); 14647c6ae99SBarry Smith } 14747c6ae99SBarry Smith 14847c6ae99SBarry Smith /*@C 14947c6ae99SBarry Smith DMCompositeGetAccess - Allows one to access the individual packed vectors in their global 15047c6ae99SBarry Smith representation. 15147c6ae99SBarry Smith 152d083f849SBarry Smith Collective on dm 15347c6ae99SBarry Smith 1549ae5db72SJed Brown Input Parameters: 15547c6ae99SBarry Smith + dm - the packer object 1569ae5db72SJed Brown - gvec - the global vector 1579ae5db72SJed Brown 1589ae5db72SJed Brown Output Parameters: 1590298fd71SBarry Smith . Vec* ... - the packed parallel vectors, NULL for those that are not needed 16047c6ae99SBarry Smith 16195452b02SPatrick Sanan Notes: 16295452b02SPatrick Sanan Use DMCompositeRestoreAccess() to return the vectors when you no longer need them 16347c6ae99SBarry Smith 164f73e5cebSJed Brown Fortran Notes: 165f73e5cebSJed Brown 166f73e5cebSJed Brown Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4) 167f73e5cebSJed Brown or use the alternative interface DMCompositeGetAccessArray(). 168f73e5cebSJed Brown 16947c6ae99SBarry Smith Level: advanced 17047c6ae99SBarry Smith 171db781477SPatrick Sanan .seealso: `DMCompositeGetEntries()`, `DMCompositeScatter()` 17247c6ae99SBarry Smith @*/ 173*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccess(DM dm, Vec gvec, ...) 174*d71ae5a4SJacob Faibussowitsch { 17547c6ae99SBarry Smith va_list Argp; 17647c6ae99SBarry Smith struct DMCompositeLink *next; 17747c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 1785edff71fSBarry Smith PetscInt readonly; 17971b14b3eSStefano Zampini PetscBool flg; 18047c6ae99SBarry Smith 18147c6ae99SBarry Smith PetscFunctionBegin; 18247c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 18347c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 1849566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 1857a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 18647c6ae99SBarry Smith next = com->next; 18748a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 18847c6ae99SBarry Smith 1899566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec, &readonly)); 19047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 19147c6ae99SBarry Smith va_start(Argp, gvec); 19247c6ae99SBarry Smith while (next) { 19347c6ae99SBarry Smith Vec *vec; 19447c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 1959ae5db72SJed Brown if (vec) { 1969566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, vec)); 1975edff71fSBarry Smith if (readonly) { 1985edff71fSBarry Smith const PetscScalar *array; 1999566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 2009566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(*vec, array + next->rstart)); 2019566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(*vec)); 2029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 2035edff71fSBarry Smith } else { 2045edff71fSBarry Smith PetscScalar *array; 2059566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 2069566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(*vec, array + next->rstart)); 2079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 20847c6ae99SBarry Smith } 2095edff71fSBarry Smith } 21047c6ae99SBarry Smith next = next->next; 21147c6ae99SBarry Smith } 21247c6ae99SBarry Smith va_end(Argp); 21347c6ae99SBarry Smith PetscFunctionReturn(0); 21447c6ae99SBarry Smith } 21547c6ae99SBarry Smith 216f73e5cebSJed Brown /*@C 217f73e5cebSJed Brown DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global 218f73e5cebSJed Brown representation. 219f73e5cebSJed Brown 220d083f849SBarry Smith Collective on dm 221f73e5cebSJed Brown 222f73e5cebSJed Brown Input Parameters: 223f73e5cebSJed Brown + dm - the packer object 224f73e5cebSJed Brown . pvec - packed vector 225f73e5cebSJed Brown . nwanted - number of vectors wanted 2260298fd71SBarry Smith - wanted - sorted array of vectors wanted, or NULL to get all vectors 227f73e5cebSJed Brown 228f73e5cebSJed Brown Output Parameters: 229f73e5cebSJed Brown . vecs - array of requested global vectors (must be allocated) 230f73e5cebSJed Brown 23195452b02SPatrick Sanan Notes: 23295452b02SPatrick Sanan Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them 233f73e5cebSJed Brown 234f73e5cebSJed Brown Level: advanced 235f73e5cebSJed Brown 236db781477SPatrick Sanan .seealso: `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 237f73e5cebSJed Brown @*/ 238*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 239*d71ae5a4SJacob Faibussowitsch { 240f73e5cebSJed Brown struct DMCompositeLink *link; 241f73e5cebSJed Brown PetscInt i, wnum; 242f73e5cebSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 243bee642f7SBarry Smith PetscInt readonly; 24471b14b3eSStefano Zampini PetscBool flg; 245f73e5cebSJed Brown 246f73e5cebSJed Brown PetscFunctionBegin; 247f73e5cebSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 248f73e5cebSJed Brown PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 2507a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 25148a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 252f73e5cebSJed Brown 2539566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 254f73e5cebSJed Brown for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 255f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 256f73e5cebSJed Brown Vec v; 2579566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(link->dm, &v)); 258bee642f7SBarry Smith if (readonly) { 259bee642f7SBarry Smith const PetscScalar *array; 2609566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvec, &array)); 2619566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + link->rstart)); 2629566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(v)); 2639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvec, &array)); 264bee642f7SBarry Smith } else { 265bee642f7SBarry Smith PetscScalar *array; 2669566063dSJacob Faibussowitsch PetscCall(VecGetArray(pvec, &array)); 2679566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + link->rstart)); 2689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pvec, &array)); 269bee642f7SBarry Smith } 270f73e5cebSJed Brown vecs[wnum++] = v; 271f73e5cebSJed Brown } 272f73e5cebSJed Brown } 273f73e5cebSJed Brown PetscFunctionReturn(0); 274f73e5cebSJed Brown } 275f73e5cebSJed Brown 2767ac2b803SAlex Fikl /*@C 2777ac2b803SAlex Fikl DMCompositeGetLocalAccessArray - Allows one to access the individual 2787ac2b803SAlex Fikl packed vectors in their local representation. 2797ac2b803SAlex Fikl 280d083f849SBarry Smith Collective on dm. 2817ac2b803SAlex Fikl 2827ac2b803SAlex Fikl Input Parameters: 2837ac2b803SAlex Fikl + dm - the packer object 2847ac2b803SAlex Fikl . pvec - packed vector 2857ac2b803SAlex Fikl . nwanted - number of vectors wanted 2867ac2b803SAlex Fikl - wanted - sorted array of vectors wanted, or NULL to get all vectors 2877ac2b803SAlex Fikl 2887ac2b803SAlex Fikl Output Parameters: 2897ac2b803SAlex Fikl . vecs - array of requested local vectors (must be allocated) 2907ac2b803SAlex Fikl 29195452b02SPatrick Sanan Notes: 29295452b02SPatrick Sanan Use DMCompositeRestoreLocalAccessArray() to return the vectors 2937ac2b803SAlex Fikl when you no longer need them. 2947ac2b803SAlex Fikl 2957ac2b803SAlex Fikl Level: advanced 2967ac2b803SAlex Fikl 297db781477SPatrick Sanan .seealso: `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`, 298db781477SPatrick Sanan `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 2997ac2b803SAlex Fikl @*/ 300*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 301*d71ae5a4SJacob Faibussowitsch { 3027ac2b803SAlex Fikl struct DMCompositeLink *link; 3037ac2b803SAlex Fikl PetscInt i, wnum; 3047ac2b803SAlex Fikl DM_Composite *com = (DM_Composite *)dm->data; 3057ac2b803SAlex Fikl PetscInt readonly; 3067ac2b803SAlex Fikl PetscInt nlocal = 0; 30771b14b3eSStefano Zampini PetscBool flg; 3087ac2b803SAlex Fikl 3097ac2b803SAlex Fikl PetscFunctionBegin; 3107ac2b803SAlex Fikl PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3117ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 3129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 3137a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 31448a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 3157ac2b803SAlex Fikl 3169566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 3177ac2b803SAlex Fikl for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 3187ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 3197ac2b803SAlex Fikl Vec v; 3209566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(link->dm, &v)); 3217ac2b803SAlex Fikl if (readonly) { 3227ac2b803SAlex Fikl const PetscScalar *array; 3239566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvec, &array)); 3249566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + nlocal)); 325b1c3483dSMark 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 3269566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(v)); 3279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvec, &array)); 3287ac2b803SAlex Fikl } else { 3297ac2b803SAlex Fikl PetscScalar *array; 3309566063dSJacob Faibussowitsch PetscCall(VecGetArray(pvec, &array)); 3319566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v, array + nlocal)); 3329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pvec, &array)); 3337ac2b803SAlex Fikl } 3347ac2b803SAlex Fikl vecs[wnum++] = v; 3357ac2b803SAlex Fikl } 3367ac2b803SAlex Fikl 3377ac2b803SAlex Fikl nlocal += link->nlocal; 3387ac2b803SAlex Fikl } 3397ac2b803SAlex Fikl 3407ac2b803SAlex Fikl PetscFunctionReturn(0); 3417ac2b803SAlex Fikl } 3427ac2b803SAlex Fikl 34347c6ae99SBarry Smith /*@C 344aa219208SBarry Smith DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess() 34547c6ae99SBarry Smith representation. 34647c6ae99SBarry Smith 347d083f849SBarry Smith Collective on dm 34847c6ae99SBarry Smith 3499ae5db72SJed Brown Input Parameters: 35047c6ae99SBarry Smith + dm - the packer object 35147c6ae99SBarry Smith . gvec - the global vector 3520298fd71SBarry Smith - Vec* ... - the individual parallel vectors, NULL for those that are not needed 35347c6ae99SBarry Smith 35447c6ae99SBarry Smith Level: advanced 35547c6ae99SBarry Smith 356db781477SPatrick Sanan .seealso `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 357db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`, 358db781477SPatrick Sanan `DMCompositeRestoreAccess()`, `DMCompositeGetAccess()` 35947c6ae99SBarry Smith 36047c6ae99SBarry Smith @*/ 361*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccess(DM dm, Vec gvec, ...) 362*d71ae5a4SJacob Faibussowitsch { 36347c6ae99SBarry Smith va_list Argp; 36447c6ae99SBarry Smith struct DMCompositeLink *next; 36547c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 3665edff71fSBarry Smith PetscInt readonly; 36771b14b3eSStefano Zampini PetscBool flg; 36847c6ae99SBarry Smith 36947c6ae99SBarry Smith PetscFunctionBegin; 37047c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 37147c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 3729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 3737a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 37447c6ae99SBarry Smith next = com->next; 37548a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 37647c6ae99SBarry Smith 3779566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec, &readonly)); 37847c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 37947c6ae99SBarry Smith va_start(Argp, gvec); 38047c6ae99SBarry Smith while (next) { 38147c6ae99SBarry Smith Vec *vec; 38247c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 3839ae5db72SJed Brown if (vec) { 3849566063dSJacob Faibussowitsch PetscCall(VecResetArray(*vec)); 3851baa6e33SBarry Smith if (readonly) PetscCall(VecLockReadPop(*vec)); 3869566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, vec)); 38747c6ae99SBarry Smith } 38847c6ae99SBarry Smith next = next->next; 38947c6ae99SBarry Smith } 39047c6ae99SBarry Smith va_end(Argp); 39147c6ae99SBarry Smith PetscFunctionReturn(0); 39247c6ae99SBarry Smith } 39347c6ae99SBarry Smith 394f73e5cebSJed Brown /*@C 395f73e5cebSJed Brown DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray() 396f73e5cebSJed Brown 397d083f849SBarry Smith Collective on dm 398f73e5cebSJed Brown 399f73e5cebSJed Brown Input Parameters: 400f73e5cebSJed Brown + dm - the packer object 401f73e5cebSJed Brown . pvec - packed vector 402f73e5cebSJed Brown . nwanted - number of vectors wanted 4030298fd71SBarry Smith . wanted - sorted array of vectors wanted, or NULL to get all vectors 404f73e5cebSJed Brown - vecs - array of global vectors to return 405f73e5cebSJed Brown 406f73e5cebSJed Brown Level: advanced 407f73e5cebSJed Brown 408db781477SPatrick Sanan .seealso: `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 409f73e5cebSJed Brown @*/ 410*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 411*d71ae5a4SJacob Faibussowitsch { 412f73e5cebSJed Brown struct DMCompositeLink *link; 413f73e5cebSJed Brown PetscInt i, wnum; 414f73e5cebSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 415bee642f7SBarry Smith PetscInt readonly; 41671b14b3eSStefano Zampini PetscBool flg; 417f73e5cebSJed Brown 418f73e5cebSJed Brown PetscFunctionBegin; 419f73e5cebSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 420f73e5cebSJed Brown PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 4219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 4227a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 42348a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 424f73e5cebSJed Brown 4259566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 426f73e5cebSJed Brown for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 427f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 4289566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum])); 42948a46eb9SPierre Jolivet if (readonly) PetscCall(VecLockReadPop(vecs[wnum])); 4309566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(link->dm, &vecs[wnum])); 431f73e5cebSJed Brown wnum++; 432f73e5cebSJed Brown } 433f73e5cebSJed Brown } 434f73e5cebSJed Brown PetscFunctionReturn(0); 435f73e5cebSJed Brown } 436f73e5cebSJed Brown 4377ac2b803SAlex Fikl /*@C 4387ac2b803SAlex Fikl DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray(). 4397ac2b803SAlex Fikl 440d083f849SBarry Smith Collective on dm. 4417ac2b803SAlex Fikl 4427ac2b803SAlex Fikl Input Parameters: 4437ac2b803SAlex Fikl + dm - the packer object 4447ac2b803SAlex Fikl . pvec - packed vector 4457ac2b803SAlex Fikl . nwanted - number of vectors wanted 4467ac2b803SAlex Fikl . wanted - sorted array of vectors wanted, or NULL to restore all vectors 4477ac2b803SAlex Fikl - vecs - array of local vectors to return 4487ac2b803SAlex Fikl 4497ac2b803SAlex Fikl Level: advanced 4507ac2b803SAlex Fikl 4517ac2b803SAlex Fikl Notes: 4527ac2b803SAlex Fikl nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray() 4537ac2b803SAlex Fikl otherwise the call will fail. 4547ac2b803SAlex Fikl 455db781477SPatrick Sanan .seealso: `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`, 456db781477SPatrick Sanan `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, 457db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeGather()` 4587ac2b803SAlex Fikl @*/ 459*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm, Vec pvec, PetscInt nwanted, const PetscInt *wanted, Vec *vecs) 460*d71ae5a4SJacob Faibussowitsch { 4617ac2b803SAlex Fikl struct DMCompositeLink *link; 4627ac2b803SAlex Fikl PetscInt i, wnum; 4637ac2b803SAlex Fikl DM_Composite *com = (DM_Composite *)dm->data; 4647ac2b803SAlex Fikl PetscInt readonly; 46571b14b3eSStefano Zampini PetscBool flg; 4667ac2b803SAlex Fikl 4677ac2b803SAlex Fikl PetscFunctionBegin; 4687ac2b803SAlex Fikl PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4697ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec, VEC_CLASSID, 2); 4709566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 4717a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 47248a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 4737ac2b803SAlex Fikl 4749566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec, &readonly)); 4757ac2b803SAlex Fikl for (i = 0, wnum = 0, link = com->next; link && wnum < nwanted; i++, link = link->next) { 4767ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 4779566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum])); 47848a46eb9SPierre Jolivet if (readonly) PetscCall(VecLockReadPop(vecs[wnum])); 4799566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(link->dm, &vecs[wnum])); 4807ac2b803SAlex Fikl wnum++; 4817ac2b803SAlex Fikl } 4827ac2b803SAlex Fikl } 4837ac2b803SAlex Fikl PetscFunctionReturn(0); 4847ac2b803SAlex Fikl } 4857ac2b803SAlex Fikl 48647c6ae99SBarry Smith /*@C 48747c6ae99SBarry Smith DMCompositeScatter - Scatters from a global packed vector into its individual local vectors 48847c6ae99SBarry Smith 489d083f849SBarry Smith Collective on dm 49047c6ae99SBarry Smith 4919ae5db72SJed Brown Input Parameters: 49247c6ae99SBarry Smith + dm - the packer object 49347c6ae99SBarry Smith . gvec - the global vector 4940298fd71SBarry Smith - Vec ... - the individual sequential vectors, NULL for those that are not needed 49547c6ae99SBarry Smith 49647c6ae99SBarry Smith Level: advanced 49747c6ae99SBarry Smith 4986f3c3dcfSJed Brown Notes: 4996f3c3dcfSJed Brown DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is 5006f3c3dcfSJed Brown accessible from Fortran. 5016f3c3dcfSJed Brown 502db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 503db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 504db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 505db781477SPatrick Sanan `DMCompositeScatterArray()` 50647c6ae99SBarry Smith 50747c6ae99SBarry Smith @*/ 508*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatter(DM dm, Vec gvec, ...) 509*d71ae5a4SJacob Faibussowitsch { 51047c6ae99SBarry Smith va_list Argp; 51147c6ae99SBarry Smith struct DMCompositeLink *next; 512c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt; 51347c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 51471b14b3eSStefano Zampini PetscBool flg; 51547c6ae99SBarry Smith 51647c6ae99SBarry Smith PetscFunctionBegin; 51747c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 51847c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 5199566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 5207a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 52148a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 52247c6ae99SBarry Smith 52347c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 52447c6ae99SBarry Smith va_start(Argp, gvec); 5258fd8f222SJed Brown for (cnt = 3, next = com->next; next; cnt++, next = next->next) { 5269ae5db72SJed Brown Vec local; 5279ae5db72SJed Brown local = va_arg(Argp, Vec); 5289ae5db72SJed Brown if (local) { 5299ae5db72SJed Brown Vec global; 5305edff71fSBarry Smith const PetscScalar *array; 53163a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt)); 5329566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 5339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 5349566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 5359566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, local)); 5369566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, local)); 5379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 5389566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 5399566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 54047c6ae99SBarry Smith } 54147c6ae99SBarry Smith } 54247c6ae99SBarry Smith va_end(Argp); 54347c6ae99SBarry Smith PetscFunctionReturn(0); 54447c6ae99SBarry Smith } 54547c6ae99SBarry Smith 5466f3c3dcfSJed Brown /*@ 5476f3c3dcfSJed Brown DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors 5486f3c3dcfSJed Brown 549d083f849SBarry Smith Collective on dm 5506f3c3dcfSJed Brown 5516f3c3dcfSJed Brown Input Parameters: 5526f3c3dcfSJed Brown + dm - the packer object 5536f3c3dcfSJed Brown . gvec - the global vector 554a2b725a8SWilliam Gropp - lvecs - array of local vectors, NULL for any that are not needed 5556f3c3dcfSJed Brown 5566f3c3dcfSJed Brown Level: advanced 5576f3c3dcfSJed Brown 5586f3c3dcfSJed Brown Note: 559907376e6SBarry Smith This is a non-variadic alternative to DMCompositeScatter() 5606f3c3dcfSJed Brown 561db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()` 562db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 563db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 5646f3c3dcfSJed Brown 5656f3c3dcfSJed Brown @*/ 566*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeScatterArray(DM dm, Vec gvec, Vec *lvecs) 567*d71ae5a4SJacob Faibussowitsch { 5686f3c3dcfSJed Brown struct DMCompositeLink *next; 5696f3c3dcfSJed Brown PetscInt i; 5706f3c3dcfSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 57171b14b3eSStefano Zampini PetscBool flg; 5726f3c3dcfSJed Brown 5736f3c3dcfSJed Brown PetscFunctionBegin; 5746f3c3dcfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5756f3c3dcfSJed Brown PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 5769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 5777a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 57848a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 5796f3c3dcfSJed Brown 5806f3c3dcfSJed Brown /* loop over packed objects, handling one at at time */ 5816f3c3dcfSJed Brown for (i = 0, next = com->next; next; next = next->next, i++) { 5826f3c3dcfSJed Brown if (lvecs[i]) { 5836f3c3dcfSJed Brown Vec global; 584c5d31e75SLisandro Dalcin const PetscScalar *array; 5856f3c3dcfSJed Brown PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 3); 5869566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 5879566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 5889566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, (PetscScalar *)array + next->rstart)); 5899566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, INSERT_VALUES, lvecs[i])); 5909566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, INSERT_VALUES, lvecs[i])); 5919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 5929566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 5939566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 5946f3c3dcfSJed Brown } 5956f3c3dcfSJed Brown } 5966f3c3dcfSJed Brown PetscFunctionReturn(0); 5976f3c3dcfSJed Brown } 5986f3c3dcfSJed Brown 59947c6ae99SBarry Smith /*@C 60047c6ae99SBarry Smith DMCompositeGather - Gathers into a global packed vector from its individual local vectors 60147c6ae99SBarry Smith 602d083f849SBarry Smith Collective on dm 60347c6ae99SBarry Smith 604d8d19677SJose E. Roman Input Parameters: 60547c6ae99SBarry Smith + dm - the packer object 60647c6ae99SBarry Smith . gvec - the global vector 607907376e6SBarry Smith . imode - INSERT_VALUES or ADD_VALUES 6080298fd71SBarry Smith - Vec ... - the individual sequential vectors, NULL for any that are not needed 60947c6ae99SBarry Smith 61047c6ae99SBarry Smith Level: advanced 61147c6ae99SBarry Smith 612f5f57ec0SBarry Smith Not available from Fortran, Fortran users can use DMCompositeGatherArray() 613f5f57ec0SBarry Smith 614db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 615db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 616db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 61747c6ae99SBarry Smith 61847c6ae99SBarry Smith @*/ 619*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGather(DM dm, InsertMode imode, Vec gvec, ...) 620*d71ae5a4SJacob Faibussowitsch { 62147c6ae99SBarry Smith va_list Argp; 62247c6ae99SBarry Smith struct DMCompositeLink *next; 62347c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 624c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt; 62571b14b3eSStefano Zampini PetscBool flg; 62647c6ae99SBarry Smith 62747c6ae99SBarry Smith PetscFunctionBegin; 62847c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 629064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3); 6309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 6317a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 63248a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 63347c6ae99SBarry Smith 63447c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 6351dac896bSSatish Balay va_start(Argp, gvec); 6368fd8f222SJed Brown for (cnt = 3, next = com->next; next; cnt++, next = next->next) { 6379ae5db72SJed Brown Vec local; 6389ae5db72SJed Brown local = va_arg(Argp, Vec); 6399ae5db72SJed Brown if (local) { 64047c6ae99SBarry Smith PetscScalar *array; 6419ae5db72SJed Brown Vec global; 64263a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local, VEC_CLASSID, (int)cnt)); 6439566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 6449566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 6459566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 6469566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, local, imode, global)); 6479566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, local, imode, global)); 6489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 6499566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 6509566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 65147c6ae99SBarry Smith } 65247c6ae99SBarry Smith } 65347c6ae99SBarry Smith va_end(Argp); 65447c6ae99SBarry Smith PetscFunctionReturn(0); 65547c6ae99SBarry Smith } 65647c6ae99SBarry Smith 6576f3c3dcfSJed Brown /*@ 6586f3c3dcfSJed Brown DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors 6596f3c3dcfSJed Brown 660d083f849SBarry Smith Collective on dm 6616f3c3dcfSJed Brown 662d8d19677SJose E. Roman Input Parameters: 6636f3c3dcfSJed Brown + dm - the packer object 6646f3c3dcfSJed Brown . gvec - the global vector 665907376e6SBarry Smith . imode - INSERT_VALUES or ADD_VALUES 6666f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed 6676f3c3dcfSJed Brown 6686f3c3dcfSJed Brown Level: advanced 6696f3c3dcfSJed Brown 6706f3c3dcfSJed Brown Notes: 6716f3c3dcfSJed Brown This is a non-variadic alternative to DMCompositeGather(). 6726f3c3dcfSJed Brown 673db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 674db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 675db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`, 6766f3c3dcfSJed Brown @*/ 677*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGatherArray(DM dm, InsertMode imode, Vec gvec, Vec *lvecs) 678*d71ae5a4SJacob Faibussowitsch { 6796f3c3dcfSJed Brown struct DMCompositeLink *next; 6806f3c3dcfSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 6816f3c3dcfSJed Brown PetscInt i; 68271b14b3eSStefano Zampini PetscBool flg; 6836f3c3dcfSJed Brown 6846f3c3dcfSJed Brown PetscFunctionBegin; 6856f3c3dcfSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 686064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec, VEC_CLASSID, 3); 6879566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 6887a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 68948a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 6906f3c3dcfSJed Brown 6916f3c3dcfSJed Brown /* loop over packed objects, handling one at at time */ 6926f3c3dcfSJed Brown for (next = com->next, i = 0; next; next = next->next, i++) { 6936f3c3dcfSJed Brown if (lvecs[i]) { 6946f3c3dcfSJed Brown PetscScalar *array; 6956f3c3dcfSJed Brown Vec global; 696064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(lvecs[i], VEC_CLASSID, 4); 6979566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 6989566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &array)); 6999566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, array + next->rstart)); 7009566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, lvecs[i], imode, global)); 7019566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, lvecs[i], imode, global)); 7029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, &array)); 7039566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 7049566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 7056f3c3dcfSJed Brown } 7066f3c3dcfSJed Brown } 7076f3c3dcfSJed Brown PetscFunctionReturn(0); 7086f3c3dcfSJed Brown } 7096f3c3dcfSJed Brown 710f5f57ec0SBarry Smith /*@ 711aa219208SBarry Smith DMCompositeAddDM - adds a DM vector to a DMComposite 71247c6ae99SBarry Smith 713d083f849SBarry Smith Collective on dm 71447c6ae99SBarry Smith 715d8d19677SJose E. Roman Input Parameters: 7169b52a9b5SPatrick Sanan + dmc - the DMComposite (packer) object 717f5f57ec0SBarry Smith - dm - the DM object 71847c6ae99SBarry Smith 71947c6ae99SBarry Smith Level: advanced 72047c6ae99SBarry Smith 721db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeGather()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 722db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 723db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 72447c6ae99SBarry Smith 72547c6ae99SBarry Smith @*/ 726*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeAddDM(DM dmc, DM dm) 727*d71ae5a4SJacob Faibussowitsch { 72806ebdd98SJed Brown PetscInt n, nlocal; 72947c6ae99SBarry Smith struct DMCompositeLink *mine, *next; 73006ebdd98SJed Brown Vec global, local; 73147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dmc->data; 73271b14b3eSStefano Zampini PetscBool flg; 73347c6ae99SBarry Smith 73447c6ae99SBarry Smith PetscFunctionBegin; 73547c6ae99SBarry Smith PetscValidHeaderSpecific(dmc, DM_CLASSID, 1); 73647c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 7379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dmc, DMCOMPOSITE, &flg)); 7387a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 73947c6ae99SBarry Smith next = com->next; 74028b400f6SJacob Faibussowitsch PetscCheck(!com->setup, PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONGSTATE, "Cannot add a DM once you have used the DMComposite"); 74147c6ae99SBarry Smith 74247c6ae99SBarry Smith /* create new link */ 7439566063dSJacob Faibussowitsch PetscCall(PetscNew(&mine)); 7449566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 7459566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &global)); 7469566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global, &n)); 7479566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &global)); 7489566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &local)); 7499566063dSJacob Faibussowitsch PetscCall(VecGetSize(local, &nlocal)); 7509566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &local)); 7518865f1eaSKarl Rupp 75247c6ae99SBarry Smith mine->n = n; 75306ebdd98SJed Brown mine->nlocal = nlocal; 75447c6ae99SBarry Smith mine->dm = dm; 7550298fd71SBarry Smith mine->next = NULL; 75647c6ae99SBarry Smith com->n += n; 7577ac2b803SAlex Fikl com->nghost += nlocal; 75847c6ae99SBarry Smith 75947c6ae99SBarry Smith /* add to end of list */ 7608865f1eaSKarl Rupp if (!next) com->next = mine; 7618865f1eaSKarl Rupp else { 76247c6ae99SBarry Smith while (next->next) next = next->next; 76347c6ae99SBarry Smith next->next = mine; 76447c6ae99SBarry Smith } 76547c6ae99SBarry Smith com->nDM++; 76647c6ae99SBarry Smith com->nmine++; 76747c6ae99SBarry Smith PetscFunctionReturn(0); 76847c6ae99SBarry Smith } 76947c6ae99SBarry Smith 7709804daf3SBarry Smith #include <petscdraw.h> 77126887b52SJed Brown PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer); 772*d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_DMComposite(Vec gvec, PetscViewer viewer) 773*d71ae5a4SJacob Faibussowitsch { 77447c6ae99SBarry Smith DM dm; 77547c6ae99SBarry Smith struct DMCompositeLink *next; 77647c6ae99SBarry Smith PetscBool isdraw; 777cef07954SSatish Balay DM_Composite *com; 77847c6ae99SBarry Smith 77947c6ae99SBarry Smith PetscFunctionBegin; 7809566063dSJacob Faibussowitsch PetscCall(VecGetDM(gvec, &dm)); 7817a8be351SBarry Smith PetscCheck(dm, PetscObjectComm((PetscObject)gvec), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMComposite"); 78247c6ae99SBarry Smith com = (DM_Composite *)dm->data; 78347c6ae99SBarry Smith next = com->next; 78447c6ae99SBarry Smith 7859566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 78647c6ae99SBarry Smith if (!isdraw) { 78747c6ae99SBarry Smith /* do I really want to call this? */ 7889566063dSJacob Faibussowitsch PetscCall(VecView_MPI(gvec, viewer)); 78947c6ae99SBarry Smith } else { 79047c6ae99SBarry Smith PetscInt cnt = 0; 79147c6ae99SBarry Smith 79247c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 79347c6ae99SBarry Smith while (next) { 79447c6ae99SBarry Smith Vec vec; 795ca4278abSLisandro Dalcin const PetscScalar *array; 79647c6ae99SBarry Smith PetscInt bs; 79747c6ae99SBarry Smith 7989ae5db72SJed Brown /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */ 7999566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &vec)); 8009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec, &array)); 8019566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(vec, (PetscScalar *)array + next->rstart)); 8029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec, &array)); 8039566063dSJacob Faibussowitsch PetscCall(VecView(vec, viewer)); 8049566063dSJacob Faibussowitsch PetscCall(VecResetArray(vec)); 8059566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vec, &bs)); 8069566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &vec)); 8079566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer, bs)); 80847c6ae99SBarry Smith cnt += bs; 80947c6ae99SBarry Smith next = next->next; 81047c6ae99SBarry Smith } 8119566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer, -cnt)); 81247c6ae99SBarry Smith } 81347c6ae99SBarry Smith PetscFunctionReturn(0); 81447c6ae99SBarry Smith } 81547c6ae99SBarry Smith 816*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateGlobalVector_Composite(DM dm, Vec *gvec) 817*d71ae5a4SJacob Faibussowitsch { 81847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 81947c6ae99SBarry Smith 82047c6ae99SBarry Smith PetscFunctionBegin; 82147c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8229566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 8239566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 8249566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec)); 8259566063dSJacob Faibussowitsch PetscCall(VecSetType(*gvec, dm->vectype)); 8269566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*gvec, com->n, com->N)); 8279566063dSJacob Faibussowitsch PetscCall(VecSetDM(*gvec, dm)); 8289566063dSJacob Faibussowitsch PetscCall(VecSetOperation(*gvec, VECOP_VIEW, (void (*)(void))VecView_DMComposite)); 82947c6ae99SBarry Smith PetscFunctionReturn(0); 83047c6ae99SBarry Smith } 83147c6ae99SBarry Smith 832*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalVector_Composite(DM dm, Vec *lvec) 833*d71ae5a4SJacob Faibussowitsch { 83447c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 83547c6ae99SBarry Smith 83647c6ae99SBarry Smith PetscFunctionBegin; 83747c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 83847c6ae99SBarry Smith if (!com->setup) { 8399566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 8409566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 84147c6ae99SBarry Smith } 8429566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, lvec)); 8439566063dSJacob Faibussowitsch PetscCall(VecSetType(*lvec, dm->vectype)); 8449566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*lvec, com->nghost, PETSC_DECIDE)); 8459566063dSJacob Faibussowitsch PetscCall(VecSetDM(*lvec, dm)); 84647c6ae99SBarry Smith PetscFunctionReturn(0); 84747c6ae99SBarry Smith } 84847c6ae99SBarry Smith 84947c6ae99SBarry Smith /*@C 8509ae5db72SJed Brown DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space 85147c6ae99SBarry Smith 85206ebdd98SJed Brown Collective on DM 85347c6ae99SBarry Smith 85447c6ae99SBarry Smith Input Parameter: 85547c6ae99SBarry Smith . dm - the packer object 85647c6ae99SBarry Smith 85747c6ae99SBarry Smith Output Parameters: 8589ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes 8599ae5db72SJed Brown all the ghost points that individual ghosted DMDA's may have. 86047c6ae99SBarry Smith 86147c6ae99SBarry Smith Level: advanced 86247c6ae99SBarry Smith 86347c6ae99SBarry Smith Notes: 8646eb61c8cSJed Brown Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree(). 86547c6ae99SBarry Smith 866f5f57ec0SBarry Smith Not available from Fortran 867f5f57ec0SBarry Smith 868db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 869db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 870c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 87147c6ae99SBarry Smith 87247c6ae99SBarry Smith @*/ 873*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm, ISLocalToGlobalMapping **ltogs) 874*d71ae5a4SJacob Faibussowitsch { 87547c6ae99SBarry Smith PetscInt i, *idx, n, cnt; 87647c6ae99SBarry Smith struct DMCompositeLink *next; 87747c6ae99SBarry Smith PetscMPIInt rank; 87847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 87971b14b3eSStefano Zampini PetscBool flg; 88047c6ae99SBarry Smith 88147c6ae99SBarry Smith PetscFunctionBegin; 88247c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 8847a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 8859566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 8869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM, ltogs)); 88747c6ae99SBarry Smith next = com->next; 8889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 88947c6ae99SBarry Smith 89047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 89147c6ae99SBarry Smith cnt = 0; 89247c6ae99SBarry Smith while (next) { 8936eb61c8cSJed Brown ISLocalToGlobalMapping ltog; 8946eb61c8cSJed Brown PetscMPIInt size; 89586994e45SJed Brown const PetscInt *suboff, *indices; 8966eb61c8cSJed Brown Vec global; 89747c6ae99SBarry Smith 8986eb61c8cSJed Brown /* Get sub-DM global indices for each local dof */ 8999566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(next->dm, <og)); 9009566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltog, &n)); 9019566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltog, &indices)); 9029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idx)); 90347c6ae99SBarry Smith 9046eb61c8cSJed Brown /* Get the offsets for the sub-DM global vector */ 9059566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 9069566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(global, &suboff)); 9079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global), &size)); 9086eb61c8cSJed Brown 9096eb61c8cSJed Brown /* Shift the sub-DM definition of the global space to the composite global space */ 9106eb61c8cSJed Brown for (i = 0; i < n; i++) { 91186994e45SJed Brown PetscInt subi = indices[i], lo = 0, hi = size, t; 912e333b035SStefano Zampini /* There's no consensus on what a negative index means, 913e333b035SStefano Zampini except for skipping when setting the values in vectors and matrices */ 9149371c9d4SSatish Balay if (subi < 0) { 9159371c9d4SSatish Balay idx[i] = subi - next->grstarts[rank]; 9169371c9d4SSatish Balay continue; 9179371c9d4SSatish Balay } 9186eb61c8cSJed Brown /* Binary search to find which rank owns subi */ 9196eb61c8cSJed Brown while (hi - lo > 1) { 9206eb61c8cSJed Brown t = lo + (hi - lo) / 2; 9216eb61c8cSJed Brown if (suboff[t] > subi) hi = t; 9226eb61c8cSJed Brown else lo = t; 9236eb61c8cSJed Brown } 9246eb61c8cSJed Brown idx[i] = subi - suboff[lo] + next->grstarts[lo]; 9256eb61c8cSJed Brown } 9269566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog, &indices)); 9279566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, n, idx, PETSC_OWN_POINTER, &(*ltogs)[cnt])); 9289566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 92947c6ae99SBarry Smith next = next->next; 93047c6ae99SBarry Smith cnt++; 93147c6ae99SBarry Smith } 93247c6ae99SBarry Smith PetscFunctionReturn(0); 93347c6ae99SBarry Smith } 93447c6ae99SBarry Smith 93587c85e80SJed Brown /*@C 9369ae5db72SJed Brown DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector 93787c85e80SJed Brown 93887c85e80SJed Brown Not Collective 93987c85e80SJed Brown 9404165533cSJose E. Roman Input Parameter: 94187c85e80SJed Brown . dm - composite DM 94287c85e80SJed Brown 9434165533cSJose E. Roman Output Parameter: 94487c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite 94587c85e80SJed Brown 94687c85e80SJed Brown Level: intermediate 94787c85e80SJed Brown 94887c85e80SJed Brown Notes: 94987c85e80SJed Brown At present, a composite local vector does not normally exist. This function is used to provide index sets for 95087c85e80SJed Brown MatGetLocalSubMatrix(). In the future, the scatters for each entry in the DMComposite may be be merged into a single 9519ae5db72SJed Brown scatter to a composite local vector. The user should not typically need to know which is being done. 95287c85e80SJed Brown 95387c85e80SJed Brown To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings(). 95487c85e80SJed Brown 95587c85e80SJed Brown To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs(). 95687c85e80SJed Brown 95787c85e80SJed Brown Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree(). 95887c85e80SJed Brown 959f5f57ec0SBarry Smith Not available from Fortran 960f5f57ec0SBarry Smith 961db781477SPatrick Sanan .seealso: `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`, `MatCreateLocalRef()` 96287c85e80SJed Brown @*/ 963*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalISs(DM dm, IS **is) 964*d71ae5a4SJacob Faibussowitsch { 96587c85e80SJed Brown DM_Composite *com = (DM_Composite *)dm->data; 96687c85e80SJed Brown struct DMCompositeLink *link; 96787c85e80SJed Brown PetscInt cnt, start; 96871b14b3eSStefano Zampini PetscBool flg; 96987c85e80SJed Brown 97087c85e80SJed Brown PetscFunctionBegin; 97187c85e80SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 97287c85e80SJed Brown PetscValidPointer(is, 2); 9739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 9747a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 9759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nmine, is)); 97606ebdd98SJed Brown for (cnt = 0, start = 0, link = com->next; link; start += link->nlocal, cnt++, link = link->next) { 977520db06cSJed Brown PetscInt bs; 9789566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, link->nlocal, start, 1, &(*is)[cnt])); 9799566063dSJacob Faibussowitsch PetscCall(DMGetBlockSize(link->dm, &bs)); 9809566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize((*is)[cnt], bs)); 981520db06cSJed Brown } 98287c85e80SJed Brown PetscFunctionReturn(0); 98387c85e80SJed Brown } 98487c85e80SJed Brown 98547c6ae99SBarry Smith /*@C 98647c6ae99SBarry Smith DMCompositeGetGlobalISs - Gets the index sets for each composed object 98747c6ae99SBarry Smith 988d083f849SBarry Smith Collective on dm 98947c6ae99SBarry Smith 99047c6ae99SBarry Smith Input Parameter: 99147c6ae99SBarry Smith . dm - the packer object 99247c6ae99SBarry Smith 99347c6ae99SBarry Smith Output Parameters: 99447c6ae99SBarry Smith . is - the array of index sets 99547c6ae99SBarry Smith 99647c6ae99SBarry Smith Level: advanced 99747c6ae99SBarry Smith 99847c6ae99SBarry Smith Notes: 99947c6ae99SBarry Smith The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree() 100047c6ae99SBarry Smith 100147c6ae99SBarry Smith These could be used to extract a subset of vector entries for a "multi-physics" preconditioner 100247c6ae99SBarry Smith 10036eb61c8cSJed Brown Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and 10046eb61c8cSJed Brown DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global 10056eb61c8cSJed Brown indices. 100647c6ae99SBarry Smith 1007f3cb0f7eSJed Brown Fortran Notes: 1008f3cb0f7eSJed Brown 1009f3cb0f7eSJed Brown The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM(). 1010f3cb0f7eSJed Brown 1011db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1012db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 1013c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 101447c6ae99SBarry Smith 101547c6ae99SBarry Smith @*/ 1016*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetGlobalISs(DM dm, IS *is[]) 1017*d71ae5a4SJacob Faibussowitsch { 101866bb578eSMark Adams PetscInt cnt = 0; 101947c6ae99SBarry Smith struct DMCompositeLink *next; 102047c6ae99SBarry Smith PetscMPIInt rank; 102147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 102271b14b3eSStefano Zampini PetscBool flg; 102347c6ae99SBarry Smith 102447c6ae99SBarry Smith PetscFunctionBegin; 102547c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 10277a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 10287a8be351SBarry Smith PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMSetUp() before"); 10299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM, is)); 103047c6ae99SBarry Smith next = com->next; 10319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 103247c6ae99SBarry Smith 103347c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 103447c6ae99SBarry Smith while (next) { 1035e5e52638SMatthew G. Knepley PetscDS prob; 1036e5e52638SMatthew G. Knepley 10379566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm), next->n, next->grstart, 1, &(*is)[cnt])); 10389566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 1039e5e52638SMatthew G. Knepley if (prob) { 104065c226d8SMatthew G. Knepley MatNullSpace space; 104165c226d8SMatthew G. Knepley Mat pmat; 10420f21e855SMatthew G. Knepley PetscObject disc; 10430f21e855SMatthew G. Knepley PetscInt Nf; 104465c226d8SMatthew G. Knepley 10459566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 1046f24dd8d2SMatthew G. Knepley if (cnt < Nf) { 10479566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(prob, cnt, &disc)); 10489566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject *)&space)); 10499566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nullspace", (PetscObject)space)); 10509566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject *)&space)); 10519566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "nearnullspace", (PetscObject)space)); 10529566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject *)&pmat)); 10539566063dSJacob Faibussowitsch if (pmat) PetscCall(PetscObjectCompose((PetscObject)(*is)[cnt], "pmat", (PetscObject)pmat)); 105465c226d8SMatthew G. Knepley } 1055f24dd8d2SMatthew G. Knepley } 105647c6ae99SBarry Smith cnt++; 105747c6ae99SBarry Smith next = next->next; 105847c6ae99SBarry Smith } 105947c6ae99SBarry Smith PetscFunctionReturn(0); 106047c6ae99SBarry Smith } 106147c6ae99SBarry Smith 1062*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) 1063*d71ae5a4SJacob Faibussowitsch { 10644d343eeaSMatthew G Knepley PetscInt nDM; 10654d343eeaSMatthew G Knepley DM *dms; 10664d343eeaSMatthew G Knepley PetscInt i; 10674d343eeaSMatthew G Knepley 10684d343eeaSMatthew G Knepley PetscFunctionBegin; 10699566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 10708865f1eaSKarl Rupp if (numFields) *numFields = nDM; 10719566063dSJacob Faibussowitsch PetscCall(DMCompositeGetGlobalISs(dm, fields)); 10724d343eeaSMatthew G Knepley if (fieldNames) { 10739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, &dms)); 10749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, fieldNames)); 10759566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms)); 10764d343eeaSMatthew G Knepley for (i = 0; i < nDM; i++) { 10774d343eeaSMatthew G Knepley char buf[256]; 10784d343eeaSMatthew G Knepley const char *splitname; 10794d343eeaSMatthew G Knepley 10804d343eeaSMatthew G Knepley /* Split naming precedence: object name, prefix, number */ 10814d343eeaSMatthew G Knepley splitname = ((PetscObject)dm)->name; 10824d343eeaSMatthew G Knepley if (!splitname) { 10839566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i], &splitname)); 10844d343eeaSMatthew G Knepley if (splitname) { 10854d343eeaSMatthew G Knepley size_t len; 10869566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(buf, splitname, sizeof(buf))); 10878caf3d72SBarry Smith buf[sizeof(buf) - 1] = 0; 10889566063dSJacob Faibussowitsch PetscCall(PetscStrlen(buf, &len)); 10894d343eeaSMatthew G Knepley if (buf[len - 1] == '_') buf[len - 1] = 0; /* Remove trailing underscore if it was used */ 10904d343eeaSMatthew G Knepley splitname = buf; 10914d343eeaSMatthew G Knepley } 10924d343eeaSMatthew G Knepley } 10934d343eeaSMatthew G Knepley if (!splitname) { 109463a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, i)); 10954d343eeaSMatthew G Knepley splitname = buf; 10964d343eeaSMatthew G Knepley } 10979566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(splitname, &(*fieldNames)[i])); 10984d343eeaSMatthew G Knepley } 10999566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 11004d343eeaSMatthew G Knepley } 11014d343eeaSMatthew G Knepley PetscFunctionReturn(0); 11024d343eeaSMatthew G Knepley } 11034d343eeaSMatthew G Knepley 1104e7c4fc90SDmitry Karpeev /* 1105e7c4fc90SDmitry Karpeev This could take over from DMCreateFieldIS(), as it is more general, 11060298fd71SBarry Smith making DMCreateFieldIS() a special case -- calling with dmlist == NULL; 1107e7c4fc90SDmitry Karpeev At this point it's probably best to be less intrusive, however. 1108e7c4fc90SDmitry Karpeev */ 1109*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 1110*d71ae5a4SJacob Faibussowitsch { 1111e7c4fc90SDmitry Karpeev PetscInt nDM; 1112e7c4fc90SDmitry Karpeev PetscInt i; 1113e7c4fc90SDmitry Karpeev 1114e7c4fc90SDmitry Karpeev PetscFunctionBegin; 11159566063dSJacob Faibussowitsch PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist)); 1116e7c4fc90SDmitry Karpeev if (dmlist) { 11179566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 11189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, dmlist)); 11199566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, *dmlist)); 112048a46eb9SPierre Jolivet for (i = 0; i < nDM; i++) PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i]))); 1121e7c4fc90SDmitry Karpeev } 1122e7c4fc90SDmitry Karpeev PetscFunctionReturn(0); 1123e7c4fc90SDmitry Karpeev } 1124e7c4fc90SDmitry Karpeev 112547c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 112647c6ae99SBarry Smith /*@C 11279ae5db72SJed Brown DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite. 112847c6ae99SBarry Smith Use DMCompositeRestoreLocalVectors() to return them. 112947c6ae99SBarry Smith 113047c6ae99SBarry Smith Not Collective 113147c6ae99SBarry Smith 113247c6ae99SBarry Smith Input Parameter: 113347c6ae99SBarry Smith . dm - the packer object 113447c6ae99SBarry Smith 113547c6ae99SBarry Smith Output Parameter: 11369ae5db72SJed Brown . Vec ... - the individual sequential Vecs 113747c6ae99SBarry Smith 113847c6ae99SBarry Smith Level: advanced 113947c6ae99SBarry Smith 1140f5f57ec0SBarry Smith Not available from Fortran 1141f5f57ec0SBarry Smith 1142db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1143db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1144db781477SPatrick Sanan `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 114547c6ae99SBarry Smith 114647c6ae99SBarry Smith @*/ 1147*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetLocalVectors(DM dm, ...) 1148*d71ae5a4SJacob Faibussowitsch { 114947c6ae99SBarry Smith va_list Argp; 115047c6ae99SBarry Smith struct DMCompositeLink *next; 115147c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 115271b14b3eSStefano Zampini PetscBool flg; 115347c6ae99SBarry Smith 115447c6ae99SBarry Smith PetscFunctionBegin; 115547c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11569566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 11577a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 115847c6ae99SBarry Smith next = com->next; 115947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 116047c6ae99SBarry Smith va_start(Argp, dm); 116147c6ae99SBarry Smith while (next) { 116247c6ae99SBarry Smith Vec *vec; 116347c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 11649566063dSJacob Faibussowitsch if (vec) PetscCall(DMGetLocalVector(next->dm, vec)); 116547c6ae99SBarry Smith next = next->next; 116647c6ae99SBarry Smith } 116747c6ae99SBarry Smith va_end(Argp); 116847c6ae99SBarry Smith PetscFunctionReturn(0); 116947c6ae99SBarry Smith } 117047c6ae99SBarry Smith 117147c6ae99SBarry Smith /*@C 11729ae5db72SJed Brown DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite. 117347c6ae99SBarry Smith 117447c6ae99SBarry Smith Not Collective 117547c6ae99SBarry Smith 117647c6ae99SBarry Smith Input Parameter: 117747c6ae99SBarry Smith . dm - the packer object 117847c6ae99SBarry Smith 117947c6ae99SBarry Smith Output Parameter: 11809ae5db72SJed Brown . Vec ... - the individual sequential Vecs 118147c6ae99SBarry Smith 118247c6ae99SBarry Smith Level: advanced 118347c6ae99SBarry Smith 1184f5f57ec0SBarry Smith Not available from Fortran 1185f5f57ec0SBarry Smith 1186db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1187db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1188db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 118947c6ae99SBarry Smith 119047c6ae99SBarry Smith @*/ 1191*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeRestoreLocalVectors(DM dm, ...) 1192*d71ae5a4SJacob Faibussowitsch { 119347c6ae99SBarry Smith va_list Argp; 119447c6ae99SBarry Smith struct DMCompositeLink *next; 119547c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 119671b14b3eSStefano Zampini PetscBool flg; 119747c6ae99SBarry Smith 119847c6ae99SBarry Smith PetscFunctionBegin; 119947c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 12017a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 120247c6ae99SBarry Smith next = com->next; 120347c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 120447c6ae99SBarry Smith va_start(Argp, dm); 120547c6ae99SBarry Smith while (next) { 120647c6ae99SBarry Smith Vec *vec; 120747c6ae99SBarry Smith vec = va_arg(Argp, Vec *); 12089566063dSJacob Faibussowitsch if (vec) PetscCall(DMRestoreLocalVector(next->dm, vec)); 120947c6ae99SBarry Smith next = next->next; 121047c6ae99SBarry Smith } 121147c6ae99SBarry Smith va_end(Argp); 121247c6ae99SBarry Smith PetscFunctionReturn(0); 121347c6ae99SBarry Smith } 121447c6ae99SBarry Smith 121547c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 121647c6ae99SBarry Smith /*@C 12179ae5db72SJed Brown DMCompositeGetEntries - Gets the DM for each entry in a DMComposite. 121847c6ae99SBarry Smith 121947c6ae99SBarry Smith Not Collective 122047c6ae99SBarry Smith 122147c6ae99SBarry Smith Input Parameter: 122247c6ae99SBarry Smith . dm - the packer object 122347c6ae99SBarry Smith 122447c6ae99SBarry Smith Output Parameter: 12259ae5db72SJed Brown . DM ... - the individual entries (DMs) 122647c6ae99SBarry Smith 122747c6ae99SBarry Smith Level: advanced 122847c6ae99SBarry Smith 122995452b02SPatrick Sanan Fortran Notes: 123095452b02SPatrick Sanan Available as DMCompositeGetEntries() for one output DM, DMCompositeGetEntries2() for 2, etc 1231f5f57ec0SBarry Smith 1232db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()` 1233db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1234db781477SPatrick Sanan `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, 1235db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()` 123647c6ae99SBarry Smith 123747c6ae99SBarry Smith @*/ 1238*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntries(DM dm, ...) 1239*d71ae5a4SJacob Faibussowitsch { 124047c6ae99SBarry Smith va_list Argp; 124147c6ae99SBarry Smith struct DMCompositeLink *next; 124247c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 124371b14b3eSStefano Zampini PetscBool flg; 124447c6ae99SBarry Smith 124547c6ae99SBarry Smith PetscFunctionBegin; 124647c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 12487a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 124947c6ae99SBarry Smith next = com->next; 125047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 125147c6ae99SBarry Smith va_start(Argp, dm); 125247c6ae99SBarry Smith while (next) { 125347c6ae99SBarry Smith DM *dmn; 125447c6ae99SBarry Smith dmn = va_arg(Argp, DM *); 12559ae5db72SJed Brown if (dmn) *dmn = next->dm; 125647c6ae99SBarry Smith next = next->next; 125747c6ae99SBarry Smith } 125847c6ae99SBarry Smith va_end(Argp); 125947c6ae99SBarry Smith PetscFunctionReturn(0); 126047c6ae99SBarry Smith } 126147c6ae99SBarry Smith 1262dbab29e1SMark F. Adams /*@C 12632fa5ba8aSJed Brown DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite. 12642fa5ba8aSJed Brown 12652fa5ba8aSJed Brown Not Collective 12662fa5ba8aSJed Brown 12672fa5ba8aSJed Brown Input Parameter: 1268907376e6SBarry Smith . dm - the packer object 1269907376e6SBarry Smith 1270907376e6SBarry Smith Output Parameter: 1271907376e6SBarry Smith . dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs 12722fa5ba8aSJed Brown 12732fa5ba8aSJed Brown Level: advanced 12742fa5ba8aSJed Brown 1275db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()` 1276db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1277db781477SPatrick Sanan `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, 1278db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()` 12792fa5ba8aSJed Brown 12802fa5ba8aSJed Brown @*/ 1281*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeGetEntriesArray(DM dm, DM dms[]) 1282*d71ae5a4SJacob Faibussowitsch { 12832fa5ba8aSJed Brown struct DMCompositeLink *next; 12842fa5ba8aSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 12852fa5ba8aSJed Brown PetscInt i; 128671b14b3eSStefano Zampini PetscBool flg; 12872fa5ba8aSJed Brown 12882fa5ba8aSJed Brown PetscFunctionBegin; 12892fa5ba8aSJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12909566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMCOMPOSITE, &flg)); 12917a8be351SBarry Smith PetscCheck(flg, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Not for type %s", ((PetscObject)dm)->type_name); 12922fa5ba8aSJed Brown /* loop over packed objects, handling one at at time */ 12932fa5ba8aSJed Brown for (next = com->next, i = 0; next; next = next->next, i++) dms[i] = next->dm; 12942fa5ba8aSJed Brown PetscFunctionReturn(0); 12952fa5ba8aSJed Brown } 12962fa5ba8aSJed Brown 1297e10fd815SStefano Zampini typedef struct { 1298e10fd815SStefano Zampini DM dm; 1299e10fd815SStefano Zampini PetscViewer *subv; 1300e10fd815SStefano Zampini Vec *vecs; 1301e10fd815SStefano Zampini } GLVisViewerCtx; 1302e10fd815SStefano Zampini 1303*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 1304*d71ae5a4SJacob Faibussowitsch { 1305e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 1306e10fd815SStefano Zampini PetscInt i, n; 1307e10fd815SStefano Zampini 1308e10fd815SStefano Zampini PetscFunctionBegin; 13099566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm, &n)); 131048a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(PetscViewerDestroy(&ctx->subv[i])); 13119566063dSJacob Faibussowitsch PetscCall(PetscFree2(ctx->subv, ctx->vecs)); 13129566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx->dm)); 13139566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 1314e10fd815SStefano Zampini PetscFunctionReturn(0); 1315e10fd815SStefano Zampini } 1316e10fd815SStefano Zampini 1317*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 1318*d71ae5a4SJacob Faibussowitsch { 1319e10fd815SStefano Zampini Vec X = (Vec)oX; 1320e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx *)vctx; 1321e10fd815SStefano Zampini PetscInt i, n, cumf; 1322e10fd815SStefano Zampini 1323e10fd815SStefano Zampini PetscFunctionBegin; 13249566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm, &n)); 13259566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(ctx->dm, X, n, NULL, ctx->vecs)); 1326e10fd815SStefano Zampini for (i = 0, cumf = 0; i < n; i++) { 1327e10fd815SStefano Zampini PetscErrorCode (*g2l)(PetscObject, PetscInt, PetscObject[], void *); 1328e10fd815SStefano Zampini void *fctx; 1329e10fd815SStefano Zampini PetscInt nfi; 1330e10fd815SStefano Zampini 13319566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nfi, NULL, NULL, &g2l, NULL, &fctx)); 1332e10fd815SStefano Zampini if (!nfi) continue; 13331baa6e33SBarry Smith if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i], nfi, oXfield + cumf, fctx)); 13341baa6e33SBarry Smith else PetscCall(VecCopy(ctx->vecs[i], (Vec)(oXfield[cumf]))); 1335e10fd815SStefano Zampini cumf += nfi; 1336e10fd815SStefano Zampini } 13379566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(ctx->dm, X, n, NULL, ctx->vecs)); 1338e10fd815SStefano Zampini PetscFunctionReturn(0); 1339e10fd815SStefano Zampini } 1340e10fd815SStefano Zampini 1341*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer) 1342*d71ae5a4SJacob Faibussowitsch { 1343e10fd815SStefano Zampini DM dm = (DM)odm, *dms; 1344e10fd815SStefano Zampini Vec *Ufds; 1345e10fd815SStefano Zampini GLVisViewerCtx *ctx; 1346e10fd815SStefano Zampini PetscInt i, n, tnf, *sdim; 1347e10fd815SStefano Zampini char **fecs; 1348e10fd815SStefano Zampini 1349e10fd815SStefano Zampini PetscFunctionBegin; 13509566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 13519566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 1352e10fd815SStefano Zampini ctx->dm = dm; 13539566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &n)); 13549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &dms)); 13559566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms)); 13569566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &ctx->subv, n, &ctx->vecs)); 1357e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1358e10fd815SStefano Zampini PetscInt nf; 1359e10fd815SStefano Zampini 13609566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PetscObjectComm(odm), &ctx->subv[i])); 13619566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(ctx->subv[i], PETSCVIEWERGLVIS)); 13629566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetDM_Private(ctx->subv[i], (PetscObject)dms[i])); 13639566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nf, NULL, NULL, NULL, NULL, NULL)); 1364e10fd815SStefano Zampini tnf += nf; 1365e10fd815SStefano Zampini } 13669566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 13679566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(tnf, &fecs, tnf, &sdim, tnf, &Ufds)); 1368e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1369e10fd815SStefano Zampini PetscInt *sd, nf, f; 1370e10fd815SStefano Zampini const char **fec; 1371e10fd815SStefano Zampini Vec *Uf; 1372e10fd815SStefano Zampini 13739566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i], &nf, &fec, &sd, NULL, (PetscObject **)&Uf, NULL)); 1374e10fd815SStefano Zampini for (f = 0; f < nf; f++) { 13759566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec[f], &fecs[tnf + f])); 1376e10fd815SStefano Zampini Ufds[tnf + f] = Uf[f]; 1377e10fd815SStefano Zampini sdim[tnf + f] = sd[f]; 1378e10fd815SStefano Zampini } 1379e10fd815SStefano Zampini tnf += nf; 1380e10fd815SStefano Zampini } 13819566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetFields(viewer, tnf, (const char **)fecs, sdim, DMCompositeSampleGLVisFields_Private, (PetscObject *)Ufds, ctx, DestroyGLVisViewerCtx_Private)); 138248a46eb9SPierre Jolivet for (i = 0; i < tnf; i++) PetscCall(PetscFree(fecs[i])); 13839566063dSJacob Faibussowitsch PetscCall(PetscFree3(fecs, sdim, Ufds)); 1384e10fd815SStefano Zampini PetscFunctionReturn(0); 1385e10fd815SStefano Zampini } 1386e10fd815SStefano Zampini 1387*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_Composite(DM dmi, MPI_Comm comm, DM *fine) 1388*d71ae5a4SJacob Faibussowitsch { 138947c6ae99SBarry Smith struct DMCompositeLink *next; 139047c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dmi->data; 139147c6ae99SBarry Smith DM dm; 139247c6ae99SBarry Smith 139347c6ae99SBarry Smith PetscFunctionBegin; 139447c6ae99SBarry Smith PetscValidHeaderSpecific(dmi, DM_CLASSID, 1); 139548a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm)); 13969566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi)); 139747c6ae99SBarry Smith next = com->next; 13989566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm, fine)); 139947c6ae99SBarry Smith 140047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 140147c6ae99SBarry Smith while (next) { 14029566063dSJacob Faibussowitsch PetscCall(DMRefine(next->dm, comm, &dm)); 14039566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine, dm)); 14049566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 140547c6ae99SBarry Smith next = next->next; 140647c6ae99SBarry Smith } 140747c6ae99SBarry Smith PetscFunctionReturn(0); 140847c6ae99SBarry Smith } 140947c6ae99SBarry Smith 1410*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsen_Composite(DM dmi, MPI_Comm comm, DM *fine) 1411*d71ae5a4SJacob Faibussowitsch { 141214354c39SJed Brown struct DMCompositeLink *next; 141314354c39SJed Brown DM_Composite *com = (DM_Composite *)dmi->data; 141414354c39SJed Brown DM dm; 141514354c39SJed Brown 141614354c39SJed Brown PetscFunctionBegin; 141714354c39SJed Brown PetscValidHeaderSpecific(dmi, DM_CLASSID, 1); 14189566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi)); 141948a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmi, &comm)); 142014354c39SJed Brown next = com->next; 14219566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm, fine)); 142214354c39SJed Brown 142314354c39SJed Brown /* loop over packed objects, handling one at at time */ 142414354c39SJed Brown while (next) { 14259566063dSJacob Faibussowitsch PetscCall(DMCoarsen(next->dm, comm, &dm)); 14269566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine, dm)); 14279566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 142814354c39SJed Brown next = next->next; 142914354c39SJed Brown } 143014354c39SJed Brown PetscFunctionReturn(0); 143114354c39SJed Brown } 143247c6ae99SBarry Smith 1433*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateInterpolation_Composite(DM coarse, DM fine, Mat *A, Vec *v) 1434*d71ae5a4SJacob Faibussowitsch { 14359ae5db72SJed Brown PetscInt m, n, M, N, nDM, i; 143647c6ae99SBarry Smith struct DMCompositeLink *nextc; 143747c6ae99SBarry Smith struct DMCompositeLink *nextf; 143825296bd5SBarry Smith Vec gcoarse, gfine, *vecs; 143947c6ae99SBarry Smith DM_Composite *comcoarse = (DM_Composite *)coarse->data; 144047c6ae99SBarry Smith DM_Composite *comfine = (DM_Composite *)fine->data; 14419ae5db72SJed Brown Mat *mats; 144247c6ae99SBarry Smith 144347c6ae99SBarry Smith PetscFunctionBegin; 144447c6ae99SBarry Smith PetscValidHeaderSpecific(coarse, DM_CLASSID, 1); 144547c6ae99SBarry Smith PetscValidHeaderSpecific(fine, DM_CLASSID, 2); 14469566063dSJacob Faibussowitsch PetscCall(DMSetUp(coarse)); 14479566063dSJacob Faibussowitsch PetscCall(DMSetUp(fine)); 144847c6ae99SBarry Smith /* use global vectors only for determining matrix layout */ 14499566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(coarse, &gcoarse)); 14509566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fine, &gfine)); 14519566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gcoarse, &n)); 14529566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gfine, &m)); 14539566063dSJacob Faibussowitsch PetscCall(VecGetSize(gcoarse, &N)); 14549566063dSJacob Faibussowitsch PetscCall(VecGetSize(gfine, &M)); 14559566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(coarse, &gcoarse)); 14569566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fine, &gfine)); 145747c6ae99SBarry Smith 14589ae5db72SJed Brown nDM = comfine->nDM; 145963a3b9bcSJacob 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); 14609566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nDM * nDM, &mats)); 146148a46eb9SPierre Jolivet if (v) PetscCall(PetscCalloc1(nDM, &vecs)); 146247c6ae99SBarry Smith 146347c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 14649ae5db72SJed Brown for (nextc = comcoarse->next, nextf = comfine->next, i = 0; nextc; nextc = nextc->next, nextf = nextf->next, i++) { 14651baa6e33SBarry Smith if (!v) PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], NULL)); 14661baa6e33SBarry Smith else PetscCall(DMCreateInterpolation(nextc->dm, nextf->dm, &mats[i * nDM + i], &vecs[i])); 146747c6ae99SBarry Smith } 14689566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, nDM, NULL, mats, A)); 14691baa6e33SBarry Smith if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine), nDM, NULL, vecs, v)); 14709566063dSJacob Faibussowitsch for (i = 0; i < nDM * nDM; i++) PetscCall(MatDestroy(&mats[i])); 14719566063dSJacob Faibussowitsch PetscCall(PetscFree(mats)); 147225296bd5SBarry Smith if (v) { 14739566063dSJacob Faibussowitsch for (i = 0; i < nDM; i++) PetscCall(VecDestroy(&vecs[i])); 14749566063dSJacob Faibussowitsch PetscCall(PetscFree(vecs)); 147525296bd5SBarry Smith } 147647c6ae99SBarry Smith PetscFunctionReturn(0); 147747c6ae99SBarry Smith } 147847c6ae99SBarry Smith 1479*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm) 1480*d71ae5a4SJacob Faibussowitsch { 14811411c6eeSJed Brown DM_Composite *com = (DM_Composite *)dm->data; 14821411c6eeSJed Brown ISLocalToGlobalMapping *ltogs; 1483f7efa3c7SJed Brown PetscInt i; 14841411c6eeSJed Brown 14851411c6eeSJed Brown PetscFunctionBegin; 14861411c6eeSJed Brown /* Set the ISLocalToGlobalMapping on the new matrix */ 14879566063dSJacob Faibussowitsch PetscCall(DMCompositeGetISLocalToGlobalMappings(dm, <ogs)); 14889566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm), com->nDM, ltogs, &dm->ltogmap)); 14899566063dSJacob Faibussowitsch for (i = 0; i < com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(<ogs[i])); 14909566063dSJacob Faibussowitsch PetscCall(PetscFree(ltogs)); 14911411c6eeSJed Brown PetscFunctionReturn(0); 14921411c6eeSJed Brown } 14931411c6eeSJed Brown 1494*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateColoring_Composite(DM dm, ISColoringType ctype, ISColoring *coloring) 1495*d71ae5a4SJacob Faibussowitsch { 149647c6ae99SBarry Smith PetscInt n, i, cnt; 149747c6ae99SBarry Smith ISColoringValue *colors; 149847c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 149947c6ae99SBarry Smith ISColoringValue maxcol = 0; 150047c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 150147c6ae99SBarry Smith 150247c6ae99SBarry Smith PetscFunctionBegin; 150347c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 150408401ef6SPierre Jolivet PetscCheck(ctype != IS_COLORING_LOCAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only global coloring supported"); 1505f7d195e4SLawrence Mitchell if (ctype == IS_COLORING_GLOBAL) { 150647c6ae99SBarry Smith n = com->n; 1507ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Unknown ISColoringType"); 15089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &colors)); /* freed in ISColoringDestroy() */ 150947c6ae99SBarry Smith 15109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dmcomposite_dense_jacobian", &dense, NULL)); 151147c6ae99SBarry Smith if (dense) { 1512ad540459SPierre Jolivet for (i = 0; i < n; i++) colors[i] = (ISColoringValue)(com->rstart + i); 151347c6ae99SBarry Smith maxcol = com->N; 151447c6ae99SBarry Smith } else { 151547c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 151647c6ae99SBarry Smith PetscMPIInt rank; 151747c6ae99SBarry Smith 15189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 151947c6ae99SBarry Smith cnt = 0; 152047c6ae99SBarry Smith while (next) { 152147c6ae99SBarry Smith ISColoring lcoloring; 152247c6ae99SBarry Smith 15239566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(next->dm, IS_COLORING_GLOBAL, &lcoloring)); 1524ad540459SPierre Jolivet for (i = 0; i < lcoloring->N; i++) colors[cnt++] = maxcol + lcoloring->colors[i]; 152547c6ae99SBarry Smith maxcol += lcoloring->n; 15269566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&lcoloring)); 152747c6ae99SBarry Smith next = next->next; 152847c6ae99SBarry Smith } 152947c6ae99SBarry Smith } 15309566063dSJacob Faibussowitsch PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), maxcol, n, colors, PETSC_OWN_POINTER, coloring)); 153147c6ae99SBarry Smith PetscFunctionReturn(0); 153247c6ae99SBarry Smith } 153347c6ae99SBarry Smith 1534*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) 1535*d71ae5a4SJacob Faibussowitsch { 153647c6ae99SBarry Smith struct DMCompositeLink *next; 153747c6ae99SBarry Smith PetscScalar *garray, *larray; 153847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 153947c6ae99SBarry Smith 154047c6ae99SBarry Smith PetscFunctionBegin; 154147c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 154247c6ae99SBarry Smith PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 154339d35262SVincent Le Chenadec 154448a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 154539d35262SVincent Le Chenadec 15469566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &garray)); 15479566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec, &larray)); 154847c6ae99SBarry Smith 154947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 155039d35262SVincent Le Chenadec next = com->next; 155147c6ae99SBarry Smith while (next) { 155247c6ae99SBarry Smith Vec local, global; 155347c6ae99SBarry Smith PetscInt N; 155447c6ae99SBarry Smith 15559566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 15569566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global, &N)); 15579566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, garray)); 15589566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local)); 15599566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local, larray)); 15609566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm, global, mode, local)); 15619566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm, global, mode, local)); 15629566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 15639566063dSJacob Faibussowitsch PetscCall(VecResetArray(local)); 15649566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 15659566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local)); 156639d35262SVincent Le Chenadec 156706ebdd98SJed Brown larray += next->nlocal; 156839d35262SVincent Le Chenadec garray += next->n; 156947c6ae99SBarry Smith next = next->next; 157047c6ae99SBarry Smith } 157147c6ae99SBarry Smith 15729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, NULL)); 15739566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lvec, NULL)); 157447c6ae99SBarry Smith PetscFunctionReturn(0); 157547c6ae99SBarry Smith } 157647c6ae99SBarry Smith 1577*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm, Vec gvec, InsertMode mode, Vec lvec) 1578*d71ae5a4SJacob Faibussowitsch { 15790c010503SBarry Smith PetscFunctionBegin; 158039d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 158139d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 2); 158239d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 4); 158339d35262SVincent Le Chenadec PetscFunctionReturn(0); 158439d35262SVincent Le Chenadec } 158539d35262SVincent Le Chenadec 1586*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1587*d71ae5a4SJacob Faibussowitsch { 158839d35262SVincent Le Chenadec struct DMCompositeLink *next; 158939d35262SVincent Le Chenadec PetscScalar *larray, *garray; 159039d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite *)dm->data; 159139d35262SVincent Le Chenadec 159239d35262SVincent Le Chenadec PetscFunctionBegin; 159339d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 159439d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 159539d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 159639d35262SVincent Le Chenadec 159748a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 159839d35262SVincent Le Chenadec 15999566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec, &larray)); 16009566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec, &garray)); 160139d35262SVincent Le Chenadec 160239d35262SVincent Le Chenadec /* loop over packed objects, handling one at at time */ 160339d35262SVincent Le Chenadec next = com->next; 160439d35262SVincent Le Chenadec while (next) { 160539d35262SVincent Le Chenadec Vec global, local; 160639d35262SVincent Le Chenadec 16079566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local)); 16089566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local, larray)); 16099566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm, &global)); 16109566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global, garray)); 16119566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm, local, mode, global)); 16129566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm, local, mode, global)); 16139566063dSJacob Faibussowitsch PetscCall(VecResetArray(local)); 16149566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 16159566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm, &global)); 16169566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local)); 161739d35262SVincent Le Chenadec 161839d35262SVincent Le Chenadec garray += next->n; 161939d35262SVincent Le Chenadec larray += next->nlocal; 162039d35262SVincent Le Chenadec next = next->next; 162139d35262SVincent Le Chenadec } 162239d35262SVincent Le Chenadec 16239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec, NULL)); 16249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lvec, NULL)); 162539d35262SVincent Le Chenadec 162639d35262SVincent Le Chenadec PetscFunctionReturn(0); 162739d35262SVincent Le Chenadec } 162839d35262SVincent Le Chenadec 1629*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1630*d71ae5a4SJacob Faibussowitsch { 163139d35262SVincent Le Chenadec PetscFunctionBegin; 163239d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 163339d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 163439d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 163539d35262SVincent Le Chenadec PetscFunctionReturn(0); 163639d35262SVincent Le Chenadec } 163739d35262SVincent Le Chenadec 1638*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalBegin_Composite(DM dm, Vec vec1, InsertMode mode, Vec vec2) 1639*d71ae5a4SJacob Faibussowitsch { 164039d35262SVincent Le Chenadec struct DMCompositeLink *next; 164139d35262SVincent Le Chenadec PetscScalar *array1, *array2; 164239d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite *)dm->data; 164339d35262SVincent Le Chenadec 164439d35262SVincent Le Chenadec PetscFunctionBegin; 164539d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 164639d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec1, VEC_CLASSID, 2); 164739d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec2, VEC_CLASSID, 4); 164839d35262SVincent Le Chenadec 164948a46eb9SPierre Jolivet if (!com->setup) PetscCall(DMSetUp(dm)); 165039d35262SVincent Le Chenadec 16519566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec1, &array1)); 16529566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec2, &array2)); 165339d35262SVincent Le Chenadec 165439d35262SVincent Le Chenadec /* loop over packed objects, handling one at at time */ 165539d35262SVincent Le Chenadec next = com->next; 165639d35262SVincent Le Chenadec while (next) { 165739d35262SVincent Le Chenadec Vec local1, local2; 165839d35262SVincent Le Chenadec 16599566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local1)); 16609566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local1, array1)); 16619566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm, &local2)); 16629566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local2, array2)); 16639566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalBegin(next->dm, local1, mode, local2)); 16649566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalEnd(next->dm, local1, mode, local2)); 16659566063dSJacob Faibussowitsch PetscCall(VecResetArray(local2)); 16669566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local2)); 16679566063dSJacob Faibussowitsch PetscCall(VecResetArray(local1)); 16689566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm, &local1)); 166939d35262SVincent Le Chenadec 167039d35262SVincent Le Chenadec array1 += next->nlocal; 167139d35262SVincent Le Chenadec array2 += next->nlocal; 167239d35262SVincent Le Chenadec next = next->next; 167339d35262SVincent Le Chenadec } 167439d35262SVincent Le Chenadec 16759566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec1, NULL)); 16769566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec2, NULL)); 167739d35262SVincent Le Chenadec 167839d35262SVincent Le Chenadec PetscFunctionReturn(0); 167939d35262SVincent Le Chenadec } 168039d35262SVincent Le Chenadec 1681*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalEnd_Composite(DM dm, Vec lvec, InsertMode mode, Vec gvec) 1682*d71ae5a4SJacob Faibussowitsch { 168339d35262SVincent Le Chenadec PetscFunctionBegin; 168439d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 168539d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec, VEC_CLASSID, 2); 168639d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec, VEC_CLASSID, 4); 16870c010503SBarry Smith PetscFunctionReturn(0); 16880c010503SBarry Smith } 168947c6ae99SBarry Smith 16906ae3a549SBarry Smith /*MC 16916ae3a549SBarry Smith DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs 16926ae3a549SBarry Smith 16936ae3a549SBarry Smith Level: intermediate 16946ae3a549SBarry Smith 1695db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()` 16966ae3a549SBarry Smith M*/ 16976ae3a549SBarry Smith 1698*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p) 1699*d71ae5a4SJacob Faibussowitsch { 1700a4121054SBarry Smith DM_Composite *com; 1701a4121054SBarry Smith 1702a4121054SBarry Smith PetscFunctionBegin; 17034dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&com)); 1704a4121054SBarry Smith p->data = com; 1705a4121054SBarry Smith com->n = 0; 17067ac2b803SAlex Fikl com->nghost = 0; 17070298fd71SBarry Smith com->next = NULL; 1708a4121054SBarry Smith com->nDM = 0; 1709a4121054SBarry Smith 1710a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Composite; 1711a4121054SBarry Smith p->ops->createlocalvector = DMCreateLocalVector_Composite; 1712184d77edSJed Brown p->ops->getlocaltoglobalmapping = DMGetLocalToGlobalMapping_Composite; 17134d343eeaSMatthew G Knepley p->ops->createfieldis = DMCreateFieldIS_Composite; 171416621825SDmitry Karpeev p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite; 1715a4121054SBarry Smith p->ops->refine = DMRefine_Composite; 171614354c39SJed Brown p->ops->coarsen = DMCoarsen_Composite; 171725296bd5SBarry Smith p->ops->createinterpolation = DMCreateInterpolation_Composite; 171825296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Composite; 1719e727c939SJed Brown p->ops->getcoloring = DMCreateColoring_Composite; 1720a4121054SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite; 1721a4121054SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite; 172239d35262SVincent Le Chenadec p->ops->localtoglobalbegin = DMLocalToGlobalBegin_Composite; 172339d35262SVincent Le Chenadec p->ops->localtoglobalend = DMLocalToGlobalEnd_Composite; 172439d35262SVincent Le Chenadec p->ops->localtolocalbegin = DMLocalToLocalBegin_Composite; 172539d35262SVincent Le Chenadec p->ops->localtolocalend = DMLocalToLocalEnd_Composite; 1726a4121054SBarry Smith p->ops->destroy = DMDestroy_Composite; 1727a4121054SBarry Smith p->ops->view = DMView_Composite; 1728a4121054SBarry Smith p->ops->setup = DMSetUp_Composite; 1729e10fd815SStefano Zampini 17309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)p, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Composite)); 1731a4121054SBarry Smith PetscFunctionReturn(0); 1732a4121054SBarry Smith } 1733a4121054SBarry Smith 17341fd49c25SBarry Smith /*@ 17350c010503SBarry Smith DMCompositeCreate - Creates a vector packer, used to generate "composite" 17360c010503SBarry Smith vectors made up of several subvectors. 17370c010503SBarry Smith 1738d083f849SBarry Smith Collective 173947c6ae99SBarry Smith 174047c6ae99SBarry Smith Input Parameter: 17410c010503SBarry Smith . comm - the processors that will share the global vector 17420c010503SBarry Smith 17430c010503SBarry Smith Output Parameters: 17440c010503SBarry Smith . packer - the packer object 174547c6ae99SBarry Smith 174647c6ae99SBarry Smith Level: advanced 174747c6ae99SBarry Smith 1748c2e3fba1SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCOMPOSITE`, `DMCreate()` 1749db781477SPatrick Sanan `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()` 1750db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 175147c6ae99SBarry Smith 175247c6ae99SBarry Smith @*/ 1753*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCompositeCreate(MPI_Comm comm, DM *packer) 1754*d71ae5a4SJacob Faibussowitsch { 175547c6ae99SBarry Smith PetscFunctionBegin; 17560c010503SBarry Smith PetscValidPointer(packer, 2); 17579566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, packer)); 17589566063dSJacob Faibussowitsch PetscCall(DMSetType(*packer, DMCOMPOSITE)); 175947c6ae99SBarry Smith PetscFunctionReturn(0); 176047c6ae99SBarry Smith } 1761