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 @*/ 267087cfbeSBarry Smith PetscErrorCode DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt)) 2747c6ae99SBarry Smith { 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 386bf464f9SBarry Smith PetscErrorCode DMDestroy_Composite(DM dm) 3947c6ae99SBarry Smith { 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 587087cfbeSBarry Smith PetscErrorCode DMView_Composite(DM dm,PetscViewer v) 5947c6ae99SBarry Smith { 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 /* --------------------------------------------------------------------------------------*/ 847087cfbeSBarry Smith PetscErrorCode DMSetUp_Composite(DM dm) 8547c6ae99SBarry Smith { 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 @*/ 1357087cfbeSBarry Smith PetscErrorCode DMCompositeGetNumberDM(DM dm,PetscInt *nDM) 13647c6ae99SBarry Smith { 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 @*/ 1737087cfbeSBarry Smith PetscErrorCode DMCompositeGetAccess(DM dm,Vec gvec,...) 17447c6ae99SBarry Smith { 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; 18747c6ae99SBarry Smith if (!com->setup) { 1889566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 18947c6ae99SBarry Smith } 19047c6ae99SBarry Smith 1919566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec,&readonly)); 19247c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 19347c6ae99SBarry Smith va_start(Argp,gvec); 19447c6ae99SBarry Smith while (next) { 19547c6ae99SBarry Smith Vec *vec; 19647c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 1979ae5db72SJed Brown if (vec) { 1989566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm,vec)); 1995edff71fSBarry Smith if (readonly) { 2005edff71fSBarry Smith const PetscScalar *array; 2019566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec,&array)); 2029566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(*vec,array+next->rstart)); 2039566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(*vec)); 2049566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec,&array)); 2055edff71fSBarry Smith } else { 2065edff71fSBarry Smith PetscScalar *array; 2079566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec,&array)); 2089566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(*vec,array+next->rstart)); 2099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec,&array)); 21047c6ae99SBarry Smith } 2115edff71fSBarry Smith } 21247c6ae99SBarry Smith next = next->next; 21347c6ae99SBarry Smith } 21447c6ae99SBarry Smith va_end(Argp); 21547c6ae99SBarry Smith PetscFunctionReturn(0); 21647c6ae99SBarry Smith } 21747c6ae99SBarry Smith 218f73e5cebSJed Brown /*@C 219f73e5cebSJed Brown DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global 220f73e5cebSJed Brown representation. 221f73e5cebSJed Brown 222d083f849SBarry Smith Collective on dm 223f73e5cebSJed Brown 224f73e5cebSJed Brown Input Parameters: 225f73e5cebSJed Brown + dm - the packer object 226f73e5cebSJed Brown . pvec - packed vector 227f73e5cebSJed Brown . nwanted - number of vectors wanted 2280298fd71SBarry Smith - wanted - sorted array of vectors wanted, or NULL to get all vectors 229f73e5cebSJed Brown 230f73e5cebSJed Brown Output Parameters: 231f73e5cebSJed Brown . vecs - array of requested global vectors (must be allocated) 232f73e5cebSJed Brown 23395452b02SPatrick Sanan Notes: 23495452b02SPatrick Sanan Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them 235f73e5cebSJed Brown 236f73e5cebSJed Brown Level: advanced 237f73e5cebSJed Brown 238db781477SPatrick Sanan .seealso: `DMCompositeGetAccess()`, `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 239f73e5cebSJed Brown @*/ 240f73e5cebSJed Brown PetscErrorCode DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 241f73e5cebSJed Brown { 242f73e5cebSJed Brown struct DMCompositeLink *link; 243f73e5cebSJed Brown PetscInt i,wnum; 244f73e5cebSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 245bee642f7SBarry Smith PetscInt readonly; 24671b14b3eSStefano Zampini PetscBool flg; 247f73e5cebSJed Brown 248f73e5cebSJed Brown PetscFunctionBegin; 249f73e5cebSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 250f73e5cebSJed Brown PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 2527a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 253f73e5cebSJed Brown if (!com->setup) { 2549566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 255f73e5cebSJed Brown } 256f73e5cebSJed Brown 2579566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec,&readonly)); 258f73e5cebSJed Brown for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 259f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 260f73e5cebSJed Brown Vec v; 2619566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(link->dm,&v)); 262bee642f7SBarry Smith if (readonly) { 263bee642f7SBarry Smith const PetscScalar *array; 2649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvec,&array)); 2659566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v,array+link->rstart)); 2669566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(v)); 2679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvec,&array)); 268bee642f7SBarry Smith } else { 269bee642f7SBarry Smith PetscScalar *array; 2709566063dSJacob Faibussowitsch PetscCall(VecGetArray(pvec,&array)); 2719566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v,array+link->rstart)); 2729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pvec,&array)); 273bee642f7SBarry Smith } 274f73e5cebSJed Brown vecs[wnum++] = v; 275f73e5cebSJed Brown } 276f73e5cebSJed Brown } 277f73e5cebSJed Brown PetscFunctionReturn(0); 278f73e5cebSJed Brown } 279f73e5cebSJed Brown 2807ac2b803SAlex Fikl /*@C 2817ac2b803SAlex Fikl DMCompositeGetLocalAccessArray - Allows one to access the individual 2827ac2b803SAlex Fikl packed vectors in their local representation. 2837ac2b803SAlex Fikl 284d083f849SBarry Smith Collective on dm. 2857ac2b803SAlex Fikl 2867ac2b803SAlex Fikl Input Parameters: 2877ac2b803SAlex Fikl + dm - the packer object 2887ac2b803SAlex Fikl . pvec - packed vector 2897ac2b803SAlex Fikl . nwanted - number of vectors wanted 2907ac2b803SAlex Fikl - wanted - sorted array of vectors wanted, or NULL to get all vectors 2917ac2b803SAlex Fikl 2927ac2b803SAlex Fikl Output Parameters: 2937ac2b803SAlex Fikl . vecs - array of requested local vectors (must be allocated) 2947ac2b803SAlex Fikl 29595452b02SPatrick Sanan Notes: 29695452b02SPatrick Sanan Use DMCompositeRestoreLocalAccessArray() to return the vectors 2977ac2b803SAlex Fikl when you no longer need them. 2987ac2b803SAlex Fikl 2997ac2b803SAlex Fikl Level: advanced 3007ac2b803SAlex Fikl 301db781477SPatrick Sanan .seealso: `DMCompositeRestoreLocalAccessArray()`, `DMCompositeGetAccess()`, 302db781477SPatrick Sanan `DMCompositeGetEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 3037ac2b803SAlex Fikl @*/ 3047ac2b803SAlex Fikl PetscErrorCode DMCompositeGetLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 3057ac2b803SAlex Fikl { 3067ac2b803SAlex Fikl struct DMCompositeLink *link; 3077ac2b803SAlex Fikl PetscInt i,wnum; 3087ac2b803SAlex Fikl DM_Composite *com = (DM_Composite*)dm->data; 3097ac2b803SAlex Fikl PetscInt readonly; 3107ac2b803SAlex Fikl PetscInt nlocal = 0; 31171b14b3eSStefano Zampini PetscBool flg; 3127ac2b803SAlex Fikl 3137ac2b803SAlex Fikl PetscFunctionBegin; 3147ac2b803SAlex Fikl PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3157ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 3169566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 3177a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 3187ac2b803SAlex Fikl if (!com->setup) { 3199566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 3207ac2b803SAlex Fikl } 3217ac2b803SAlex Fikl 3229566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec,&readonly)); 3237ac2b803SAlex Fikl for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 3247ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 3257ac2b803SAlex Fikl Vec v; 3269566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(link->dm,&v)); 3277ac2b803SAlex Fikl if (readonly) { 3287ac2b803SAlex Fikl const PetscScalar *array; 3299566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pvec,&array)); 3309566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v,array+nlocal)); 331b1c3483dSMark 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 3329566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(v)); 3339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pvec,&array)); 3347ac2b803SAlex Fikl } else { 3357ac2b803SAlex Fikl PetscScalar *array; 3369566063dSJacob Faibussowitsch PetscCall(VecGetArray(pvec,&array)); 3379566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(v,array+nlocal)); 3389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pvec,&array)); 3397ac2b803SAlex Fikl } 3407ac2b803SAlex Fikl vecs[wnum++] = v; 3417ac2b803SAlex Fikl } 3427ac2b803SAlex Fikl 3437ac2b803SAlex Fikl nlocal += link->nlocal; 3447ac2b803SAlex Fikl } 3457ac2b803SAlex Fikl 3467ac2b803SAlex Fikl PetscFunctionReturn(0); 3477ac2b803SAlex Fikl } 3487ac2b803SAlex Fikl 34947c6ae99SBarry Smith /*@C 350aa219208SBarry Smith DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess() 35147c6ae99SBarry Smith representation. 35247c6ae99SBarry Smith 353d083f849SBarry Smith Collective on dm 35447c6ae99SBarry Smith 3559ae5db72SJed Brown Input Parameters: 35647c6ae99SBarry Smith + dm - the packer object 35747c6ae99SBarry Smith . gvec - the global vector 3580298fd71SBarry Smith - Vec* ... - the individual parallel vectors, NULL for those that are not needed 35947c6ae99SBarry Smith 36047c6ae99SBarry Smith Level: advanced 36147c6ae99SBarry Smith 362db781477SPatrick Sanan .seealso `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 363db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeScatter()`, 364db781477SPatrick Sanan `DMCompositeRestoreAccess()`, `DMCompositeGetAccess()` 36547c6ae99SBarry Smith 36647c6ae99SBarry Smith @*/ 3677087cfbeSBarry Smith PetscErrorCode DMCompositeRestoreAccess(DM dm,Vec gvec,...) 36847c6ae99SBarry Smith { 36947c6ae99SBarry Smith va_list Argp; 37047c6ae99SBarry Smith struct DMCompositeLink *next; 37147c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 3725edff71fSBarry Smith PetscInt readonly; 37371b14b3eSStefano Zampini PetscBool flg; 37447c6ae99SBarry Smith 37547c6ae99SBarry Smith PetscFunctionBegin; 37647c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 37747c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 3789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 3797a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 38047c6ae99SBarry Smith next = com->next; 38147c6ae99SBarry Smith if (!com->setup) { 3829566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 38347c6ae99SBarry Smith } 38447c6ae99SBarry Smith 3859566063dSJacob Faibussowitsch PetscCall(VecLockGet(gvec,&readonly)); 38647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 38747c6ae99SBarry Smith va_start(Argp,gvec); 38847c6ae99SBarry Smith while (next) { 38947c6ae99SBarry Smith Vec *vec; 39047c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 3919ae5db72SJed Brown if (vec) { 3929566063dSJacob Faibussowitsch PetscCall(VecResetArray(*vec)); 3931baa6e33SBarry Smith if (readonly) PetscCall(VecLockReadPop(*vec)); 3949566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm,vec)); 39547c6ae99SBarry Smith } 39647c6ae99SBarry Smith next = next->next; 39747c6ae99SBarry Smith } 39847c6ae99SBarry Smith va_end(Argp); 39947c6ae99SBarry Smith PetscFunctionReturn(0); 40047c6ae99SBarry Smith } 40147c6ae99SBarry Smith 402f73e5cebSJed Brown /*@C 403f73e5cebSJed Brown DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray() 404f73e5cebSJed Brown 405d083f849SBarry Smith Collective on dm 406f73e5cebSJed Brown 407f73e5cebSJed Brown Input Parameters: 408f73e5cebSJed Brown + dm - the packer object 409f73e5cebSJed Brown . pvec - packed vector 410f73e5cebSJed Brown . nwanted - number of vectors wanted 4110298fd71SBarry Smith . wanted - sorted array of vectors wanted, or NULL to get all vectors 412f73e5cebSJed Brown - vecs - array of global vectors to return 413f73e5cebSJed Brown 414f73e5cebSJed Brown Level: advanced 415f73e5cebSJed Brown 416db781477SPatrick Sanan .seealso: `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, `DMCompositeScatter()`, `DMCompositeGather()` 417f73e5cebSJed Brown @*/ 418f73e5cebSJed Brown PetscErrorCode DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 419f73e5cebSJed Brown { 420f73e5cebSJed Brown struct DMCompositeLink *link; 421f73e5cebSJed Brown PetscInt i,wnum; 422f73e5cebSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 423bee642f7SBarry Smith PetscInt readonly; 42471b14b3eSStefano Zampini PetscBool flg; 425f73e5cebSJed Brown 426f73e5cebSJed Brown PetscFunctionBegin; 427f73e5cebSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 428f73e5cebSJed Brown PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 4299566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 4307a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 431f73e5cebSJed Brown if (!com->setup) { 4329566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 433f73e5cebSJed Brown } 434f73e5cebSJed Brown 4359566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec,&readonly)); 436f73e5cebSJed Brown for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 437f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 4389566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum])); 439bee642f7SBarry Smith if (readonly) { 4409566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(vecs[wnum])); 441bee642f7SBarry Smith } 4429566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(link->dm,&vecs[wnum])); 443f73e5cebSJed Brown wnum++; 444f73e5cebSJed Brown } 445f73e5cebSJed Brown } 446f73e5cebSJed Brown PetscFunctionReturn(0); 447f73e5cebSJed Brown } 448f73e5cebSJed Brown 4497ac2b803SAlex Fikl /*@C 4507ac2b803SAlex Fikl DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray(). 4517ac2b803SAlex Fikl 452d083f849SBarry Smith Collective on dm. 4537ac2b803SAlex Fikl 4547ac2b803SAlex Fikl Input Parameters: 4557ac2b803SAlex Fikl + dm - the packer object 4567ac2b803SAlex Fikl . pvec - packed vector 4577ac2b803SAlex Fikl . nwanted - number of vectors wanted 4587ac2b803SAlex Fikl . wanted - sorted array of vectors wanted, or NULL to restore all vectors 4597ac2b803SAlex Fikl - vecs - array of local vectors to return 4607ac2b803SAlex Fikl 4617ac2b803SAlex Fikl Level: advanced 4627ac2b803SAlex Fikl 4637ac2b803SAlex Fikl Notes: 4647ac2b803SAlex Fikl nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray() 4657ac2b803SAlex Fikl otherwise the call will fail. 4667ac2b803SAlex Fikl 467db781477SPatrick Sanan .seealso: `DMCompositeGetLocalAccessArray()`, `DMCompositeRestoreAccessArray()`, 468db781477SPatrick Sanan `DMCompositeRestoreAccess()`, `DMCompositeRestoreEntries()`, 469db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeGather()` 4707ac2b803SAlex Fikl @*/ 4717ac2b803SAlex Fikl PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 4727ac2b803SAlex Fikl { 4737ac2b803SAlex Fikl struct DMCompositeLink *link; 4747ac2b803SAlex Fikl PetscInt i,wnum; 4757ac2b803SAlex Fikl DM_Composite *com = (DM_Composite*)dm->data; 4767ac2b803SAlex Fikl PetscInt readonly; 47771b14b3eSStefano Zampini PetscBool flg; 4787ac2b803SAlex Fikl 4797ac2b803SAlex Fikl PetscFunctionBegin; 4807ac2b803SAlex Fikl PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4817ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 4829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 4837a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 4847ac2b803SAlex Fikl if (!com->setup) { 4859566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 4867ac2b803SAlex Fikl } 4877ac2b803SAlex Fikl 4889566063dSJacob Faibussowitsch PetscCall(VecLockGet(pvec,&readonly)); 4897ac2b803SAlex Fikl for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 4907ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 4919566063dSJacob Faibussowitsch PetscCall(VecResetArray(vecs[wnum])); 4927ac2b803SAlex Fikl if (readonly) { 4939566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(vecs[wnum])); 4947ac2b803SAlex Fikl } 4959566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(link->dm,&vecs[wnum])); 4967ac2b803SAlex Fikl wnum++; 4977ac2b803SAlex Fikl } 4987ac2b803SAlex Fikl } 4997ac2b803SAlex Fikl PetscFunctionReturn(0); 5007ac2b803SAlex Fikl } 5017ac2b803SAlex Fikl 50247c6ae99SBarry Smith /*@C 50347c6ae99SBarry Smith DMCompositeScatter - Scatters from a global packed vector into its individual local vectors 50447c6ae99SBarry Smith 505d083f849SBarry Smith Collective on dm 50647c6ae99SBarry Smith 5079ae5db72SJed Brown Input Parameters: 50847c6ae99SBarry Smith + dm - the packer object 50947c6ae99SBarry Smith . gvec - the global vector 5100298fd71SBarry Smith - Vec ... - the individual sequential vectors, NULL for those that are not needed 51147c6ae99SBarry Smith 51247c6ae99SBarry Smith Level: advanced 51347c6ae99SBarry Smith 5146f3c3dcfSJed Brown Notes: 5156f3c3dcfSJed Brown DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is 5166f3c3dcfSJed Brown accessible from Fortran. 5176f3c3dcfSJed Brown 518db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 519db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 520db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 521db781477SPatrick Sanan `DMCompositeScatterArray()` 52247c6ae99SBarry Smith 52347c6ae99SBarry Smith @*/ 5247087cfbeSBarry Smith PetscErrorCode DMCompositeScatter(DM dm,Vec gvec,...) 52547c6ae99SBarry Smith { 52647c6ae99SBarry Smith va_list Argp; 52747c6ae99SBarry Smith struct DMCompositeLink *next; 528*c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt; 52947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 53071b14b3eSStefano Zampini PetscBool flg; 53147c6ae99SBarry Smith 53247c6ae99SBarry Smith PetscFunctionBegin; 53347c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 53447c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 5359566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 5367a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 53747c6ae99SBarry Smith if (!com->setup) { 5389566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 53947c6ae99SBarry Smith } 54047c6ae99SBarry Smith 54147c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 54247c6ae99SBarry Smith va_start(Argp,gvec); 5438fd8f222SJed Brown for (cnt=3,next=com->next; next; cnt++,next=next->next) { 5449ae5db72SJed Brown Vec local; 5459ae5db72SJed Brown local = va_arg(Argp, Vec); 5469ae5db72SJed Brown if (local) { 5479ae5db72SJed Brown Vec global; 5485edff71fSBarry Smith const PetscScalar *array; 54963a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local,VEC_CLASSID,(int)cnt)); 5509566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm,&global)); 5519566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec,&array)); 5529566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global,array+next->rstart)); 5539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local)); 5549566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local)); 5559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec,&array)); 5569566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 5579566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm,&global)); 55847c6ae99SBarry Smith } 55947c6ae99SBarry Smith } 56047c6ae99SBarry Smith va_end(Argp); 56147c6ae99SBarry Smith PetscFunctionReturn(0); 56247c6ae99SBarry Smith } 56347c6ae99SBarry Smith 5646f3c3dcfSJed Brown /*@ 5656f3c3dcfSJed Brown DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors 5666f3c3dcfSJed Brown 567d083f849SBarry Smith Collective on dm 5686f3c3dcfSJed Brown 5696f3c3dcfSJed Brown Input Parameters: 5706f3c3dcfSJed Brown + dm - the packer object 5716f3c3dcfSJed Brown . gvec - the global vector 572a2b725a8SWilliam Gropp - lvecs - array of local vectors, NULL for any that are not needed 5736f3c3dcfSJed Brown 5746f3c3dcfSJed Brown Level: advanced 5756f3c3dcfSJed Brown 5766f3c3dcfSJed Brown Note: 577907376e6SBarry Smith This is a non-variadic alternative to DMCompositeScatter() 5786f3c3dcfSJed Brown 579db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()` 580db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 581db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 5826f3c3dcfSJed Brown 5836f3c3dcfSJed Brown @*/ 5846f3c3dcfSJed Brown PetscErrorCode DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs) 5856f3c3dcfSJed Brown { 5866f3c3dcfSJed Brown struct DMCompositeLink *next; 5876f3c3dcfSJed Brown PetscInt i; 5886f3c3dcfSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 58971b14b3eSStefano Zampini PetscBool flg; 5906f3c3dcfSJed Brown 5916f3c3dcfSJed Brown PetscFunctionBegin; 5926f3c3dcfSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5936f3c3dcfSJed Brown PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 5949566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 5957a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 5966f3c3dcfSJed Brown if (!com->setup) { 5979566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 5986f3c3dcfSJed Brown } 5996f3c3dcfSJed Brown 6006f3c3dcfSJed Brown /* loop over packed objects, handling one at at time */ 6016f3c3dcfSJed Brown for (i=0,next=com->next; next; next=next->next,i++) { 6026f3c3dcfSJed Brown if (lvecs[i]) { 6036f3c3dcfSJed Brown Vec global; 604c5d31e75SLisandro Dalcin const PetscScalar *array; 6056f3c3dcfSJed Brown PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3); 6069566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm,&global)); 6079566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec,&array)); 6089566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global,(PetscScalar*)array+next->rstart)); 6099566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i])); 6109566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i])); 6119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec,&array)); 6129566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 6139566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm,&global)); 6146f3c3dcfSJed Brown } 6156f3c3dcfSJed Brown } 6166f3c3dcfSJed Brown PetscFunctionReturn(0); 6176f3c3dcfSJed Brown } 6186f3c3dcfSJed Brown 61947c6ae99SBarry Smith /*@C 62047c6ae99SBarry Smith DMCompositeGather - Gathers into a global packed vector from its individual local vectors 62147c6ae99SBarry Smith 622d083f849SBarry Smith Collective on dm 62347c6ae99SBarry Smith 624d8d19677SJose E. Roman Input Parameters: 62547c6ae99SBarry Smith + dm - the packer object 62647c6ae99SBarry Smith . gvec - the global vector 627907376e6SBarry Smith . imode - INSERT_VALUES or ADD_VALUES 6280298fd71SBarry Smith - Vec ... - the individual sequential vectors, NULL for any that are not needed 62947c6ae99SBarry Smith 63047c6ae99SBarry Smith Level: advanced 63147c6ae99SBarry Smith 632f5f57ec0SBarry Smith Not available from Fortran, Fortran users can use DMCompositeGatherArray() 633f5f57ec0SBarry Smith 634db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 635db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 636db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 63747c6ae99SBarry Smith 63847c6ae99SBarry Smith @*/ 6391dac896bSSatish Balay PetscErrorCode DMCompositeGather(DM dm,InsertMode imode,Vec gvec,...) 64047c6ae99SBarry Smith { 64147c6ae99SBarry Smith va_list Argp; 64247c6ae99SBarry Smith struct DMCompositeLink *next; 64347c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 644*c599c493SJunchao Zhang PETSC_UNUSED PetscInt cnt; 64571b14b3eSStefano Zampini PetscBool flg; 64647c6ae99SBarry Smith 64747c6ae99SBarry Smith PetscFunctionBegin; 64847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 649064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec,VEC_CLASSID,3); 6509566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 6517a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 65247c6ae99SBarry Smith if (!com->setup) { 6539566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 65447c6ae99SBarry Smith } 65547c6ae99SBarry Smith 65647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 6571dac896bSSatish Balay va_start(Argp,gvec); 6588fd8f222SJed Brown for (cnt=3,next=com->next; next; cnt++,next=next->next) { 6599ae5db72SJed Brown Vec local; 6609ae5db72SJed Brown local = va_arg(Argp, Vec); 6619ae5db72SJed Brown if (local) { 66247c6ae99SBarry Smith PetscScalar *array; 6639ae5db72SJed Brown Vec global; 66463a3b9bcSJacob Faibussowitsch PetscDisableStaticAnalyzerForExpressionUnderstandingThatThisIsDangerousAndBugprone(PetscValidHeaderSpecific(local,VEC_CLASSID,(int)cnt)); 6659566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm,&global)); 6669566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec,&array)); 6679566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global,array+next->rstart)); 6689566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm,local,imode,global)); 6699566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm,local,imode,global)); 6709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec,&array)); 6719566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 6729566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm,&global)); 67347c6ae99SBarry Smith } 67447c6ae99SBarry Smith } 67547c6ae99SBarry Smith va_end(Argp); 67647c6ae99SBarry Smith PetscFunctionReturn(0); 67747c6ae99SBarry Smith } 67847c6ae99SBarry Smith 6796f3c3dcfSJed Brown /*@ 6806f3c3dcfSJed Brown DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors 6816f3c3dcfSJed Brown 682d083f849SBarry Smith Collective on dm 6836f3c3dcfSJed Brown 684d8d19677SJose E. Roman Input Parameters: 6856f3c3dcfSJed Brown + dm - the packer object 6866f3c3dcfSJed Brown . gvec - the global vector 687907376e6SBarry Smith . imode - INSERT_VALUES or ADD_VALUES 6886f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed 6896f3c3dcfSJed Brown 6906f3c3dcfSJed Brown Level: advanced 6916f3c3dcfSJed Brown 6926f3c3dcfSJed Brown Notes: 6936f3c3dcfSJed Brown This is a non-variadic alternative to DMCompositeGather(). 6946f3c3dcfSJed Brown 695db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 696db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 697db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()`, 6986f3c3dcfSJed Brown @*/ 6991dac896bSSatish Balay PetscErrorCode DMCompositeGatherArray(DM dm,InsertMode imode,Vec gvec,Vec *lvecs) 7006f3c3dcfSJed Brown { 7016f3c3dcfSJed Brown struct DMCompositeLink *next; 7026f3c3dcfSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 7036f3c3dcfSJed Brown PetscInt i; 70471b14b3eSStefano Zampini PetscBool flg; 7056f3c3dcfSJed Brown 7066f3c3dcfSJed Brown PetscFunctionBegin; 7076f3c3dcfSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 708064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(gvec,VEC_CLASSID,3); 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 7107a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 7116f3c3dcfSJed Brown if (!com->setup) { 7129566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 7136f3c3dcfSJed Brown } 7146f3c3dcfSJed Brown 7156f3c3dcfSJed Brown /* loop over packed objects, handling one at at time */ 7166f3c3dcfSJed Brown for (next=com->next,i=0; next; next=next->next,i++) { 7176f3c3dcfSJed Brown if (lvecs[i]) { 7186f3c3dcfSJed Brown PetscScalar *array; 7196f3c3dcfSJed Brown Vec global; 720064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,4); 7219566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm,&global)); 7229566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec,&array)); 7239566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global,array+next->rstart)); 7249566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global)); 7259566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global)); 7269566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec,&array)); 7279566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 7289566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm,&global)); 7296f3c3dcfSJed Brown } 7306f3c3dcfSJed Brown } 7316f3c3dcfSJed Brown PetscFunctionReturn(0); 7326f3c3dcfSJed Brown } 7336f3c3dcfSJed Brown 734f5f57ec0SBarry Smith /*@ 735aa219208SBarry Smith DMCompositeAddDM - adds a DM vector to a DMComposite 73647c6ae99SBarry Smith 737d083f849SBarry Smith Collective on dm 73847c6ae99SBarry Smith 739d8d19677SJose E. Roman Input Parameters: 7409b52a9b5SPatrick Sanan + dmc - the DMComposite (packer) object 741f5f57ec0SBarry Smith - dm - the DM object 74247c6ae99SBarry Smith 74347c6ae99SBarry Smith Level: advanced 74447c6ae99SBarry Smith 745db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeGather()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 746db781477SPatrick Sanan `DMCompositeScatter()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 747db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 74847c6ae99SBarry Smith 74947c6ae99SBarry Smith @*/ 7507087cfbeSBarry Smith PetscErrorCode DMCompositeAddDM(DM dmc,DM dm) 75147c6ae99SBarry Smith { 75206ebdd98SJed Brown PetscInt n,nlocal; 75347c6ae99SBarry Smith struct DMCompositeLink *mine,*next; 75406ebdd98SJed Brown Vec global,local; 75547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dmc->data; 75671b14b3eSStefano Zampini PetscBool flg; 75747c6ae99SBarry Smith 75847c6ae99SBarry Smith PetscFunctionBegin; 75947c6ae99SBarry Smith PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 76047c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,2); 7619566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dmc,DMCOMPOSITE,&flg)); 7627a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 76347c6ae99SBarry Smith next = com->next; 76428b400f6SJacob Faibussowitsch PetscCheck(!com->setup,PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite"); 76547c6ae99SBarry Smith 76647c6ae99SBarry Smith /* create new link */ 7679566063dSJacob Faibussowitsch PetscCall(PetscNew(&mine)); 7689566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 7699566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm,&global)); 7709566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global,&n)); 7719566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm,&global)); 7729566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&local)); 7739566063dSJacob Faibussowitsch PetscCall(VecGetSize(local,&nlocal)); 7749566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&local)); 7758865f1eaSKarl Rupp 77647c6ae99SBarry Smith mine->n = n; 77706ebdd98SJed Brown mine->nlocal = nlocal; 77847c6ae99SBarry Smith mine->dm = dm; 7790298fd71SBarry Smith mine->next = NULL; 78047c6ae99SBarry Smith com->n += n; 7817ac2b803SAlex Fikl com->nghost += nlocal; 78247c6ae99SBarry Smith 78347c6ae99SBarry Smith /* add to end of list */ 7848865f1eaSKarl Rupp if (!next) com->next = mine; 7858865f1eaSKarl Rupp else { 78647c6ae99SBarry Smith while (next->next) next = next->next; 78747c6ae99SBarry Smith next->next = mine; 78847c6ae99SBarry Smith } 78947c6ae99SBarry Smith com->nDM++; 79047c6ae99SBarry Smith com->nmine++; 79147c6ae99SBarry Smith PetscFunctionReturn(0); 79247c6ae99SBarry Smith } 79347c6ae99SBarry Smith 7949804daf3SBarry Smith #include <petscdraw.h> 79526887b52SJed Brown PETSC_EXTERN PetscErrorCode VecView_MPI(Vec,PetscViewer); 7967087cfbeSBarry Smith PetscErrorCode VecView_DMComposite(Vec gvec,PetscViewer viewer) 79747c6ae99SBarry Smith { 79847c6ae99SBarry Smith DM dm; 79947c6ae99SBarry Smith struct DMCompositeLink *next; 80047c6ae99SBarry Smith PetscBool isdraw; 801cef07954SSatish Balay DM_Composite *com; 80247c6ae99SBarry Smith 80347c6ae99SBarry Smith PetscFunctionBegin; 8049566063dSJacob Faibussowitsch PetscCall(VecGetDM(gvec, &dm)); 8057a8be351SBarry Smith PetscCheck(dm,PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite"); 80647c6ae99SBarry Smith com = (DM_Composite*)dm->data; 80747c6ae99SBarry Smith next = com->next; 80847c6ae99SBarry Smith 8099566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 81047c6ae99SBarry Smith if (!isdraw) { 81147c6ae99SBarry Smith /* do I really want to call this? */ 8129566063dSJacob Faibussowitsch PetscCall(VecView_MPI(gvec,viewer)); 81347c6ae99SBarry Smith } else { 81447c6ae99SBarry Smith PetscInt cnt = 0; 81547c6ae99SBarry Smith 81647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 81747c6ae99SBarry Smith while (next) { 81847c6ae99SBarry Smith Vec vec; 819ca4278abSLisandro Dalcin const PetscScalar *array; 82047c6ae99SBarry Smith PetscInt bs; 82147c6ae99SBarry Smith 8229ae5db72SJed Brown /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */ 8239566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm,&vec)); 8249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gvec,&array)); 8259566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(vec,(PetscScalar*)array+next->rstart)); 8269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gvec,&array)); 8279566063dSJacob Faibussowitsch PetscCall(VecView(vec,viewer)); 8289566063dSJacob Faibussowitsch PetscCall(VecResetArray(vec)); 8299566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vec,&bs)); 8309566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm,&vec)); 8319566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer,bs)); 83247c6ae99SBarry Smith cnt += bs; 83347c6ae99SBarry Smith next = next->next; 83447c6ae99SBarry Smith } 8359566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawBaseAdd(viewer,-cnt)); 83647c6ae99SBarry Smith } 83747c6ae99SBarry Smith PetscFunctionReturn(0); 83847c6ae99SBarry Smith } 83947c6ae99SBarry Smith 8407087cfbeSBarry Smith PetscErrorCode DMCreateGlobalVector_Composite(DM dm,Vec *gvec) 84147c6ae99SBarry Smith { 84247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 84347c6ae99SBarry Smith 84447c6ae99SBarry Smith PetscFunctionBegin; 84547c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 8469566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 8479566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 8489566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm),gvec)); 8499566063dSJacob Faibussowitsch PetscCall(VecSetType(*gvec,dm->vectype)); 8509566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*gvec,com->n,com->N)); 8519566063dSJacob Faibussowitsch PetscCall(VecSetDM(*gvec, dm)); 8529566063dSJacob Faibussowitsch PetscCall(VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite)); 85347c6ae99SBarry Smith PetscFunctionReturn(0); 85447c6ae99SBarry Smith } 85547c6ae99SBarry Smith 8567087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_Composite(DM dm,Vec *lvec) 85747c6ae99SBarry Smith { 85847c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 85947c6ae99SBarry Smith 86047c6ae99SBarry Smith PetscFunctionBegin; 86147c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 86247c6ae99SBarry Smith if (!com->setup) { 8639566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 8649566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 86547c6ae99SBarry Smith } 8669566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,lvec)); 8679566063dSJacob Faibussowitsch PetscCall(VecSetType(*lvec,dm->vectype)); 8689566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*lvec,com->nghost,PETSC_DECIDE)); 8699566063dSJacob Faibussowitsch PetscCall(VecSetDM(*lvec, dm)); 87047c6ae99SBarry Smith PetscFunctionReturn(0); 87147c6ae99SBarry Smith } 87247c6ae99SBarry Smith 87347c6ae99SBarry Smith /*@C 8749ae5db72SJed Brown DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space 87547c6ae99SBarry Smith 87606ebdd98SJed Brown Collective on DM 87747c6ae99SBarry Smith 87847c6ae99SBarry Smith Input Parameter: 87947c6ae99SBarry Smith . dm - the packer object 88047c6ae99SBarry Smith 88147c6ae99SBarry Smith Output Parameters: 8829ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes 8839ae5db72SJed Brown all the ghost points that individual ghosted DMDA's may have. 88447c6ae99SBarry Smith 88547c6ae99SBarry Smith Level: advanced 88647c6ae99SBarry Smith 88747c6ae99SBarry Smith Notes: 8886eb61c8cSJed Brown Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree(). 88947c6ae99SBarry Smith 890f5f57ec0SBarry Smith Not available from Fortran 891f5f57ec0SBarry Smith 892db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 893db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 894c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 89547c6ae99SBarry Smith 89647c6ae99SBarry Smith @*/ 8977087cfbeSBarry Smith PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs) 89847c6ae99SBarry Smith { 89947c6ae99SBarry Smith PetscInt i,*idx,n,cnt; 90047c6ae99SBarry Smith struct DMCompositeLink *next; 90147c6ae99SBarry Smith PetscMPIInt rank; 90247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 90371b14b3eSStefano Zampini PetscBool flg; 90447c6ae99SBarry Smith 90547c6ae99SBarry Smith PetscFunctionBegin; 90647c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 9079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 9087a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 9099566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 9109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM,ltogs)); 91147c6ae99SBarry Smith next = com->next; 9129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 91347c6ae99SBarry Smith 91447c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 91547c6ae99SBarry Smith cnt = 0; 91647c6ae99SBarry Smith while (next) { 9176eb61c8cSJed Brown ISLocalToGlobalMapping ltog; 9186eb61c8cSJed Brown PetscMPIInt size; 91986994e45SJed Brown const PetscInt *suboff,*indices; 9206eb61c8cSJed Brown Vec global; 92147c6ae99SBarry Smith 9226eb61c8cSJed Brown /* Get sub-DM global indices for each local dof */ 9239566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(next->dm,<og)); 9249566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltog,&n)); 9259566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltog,&indices)); 9269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&idx)); 92747c6ae99SBarry Smith 9286eb61c8cSJed Brown /* Get the offsets for the sub-DM global vector */ 9299566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm,&global)); 9309566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(global,&suboff)); 9319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)global),&size)); 9326eb61c8cSJed Brown 9336eb61c8cSJed Brown /* Shift the sub-DM definition of the global space to the composite global space */ 9346eb61c8cSJed Brown for (i=0; i<n; i++) { 93586994e45SJed Brown PetscInt subi = indices[i],lo = 0,hi = size,t; 936e333b035SStefano Zampini /* There's no consensus on what a negative index means, 937e333b035SStefano Zampini except for skipping when setting the values in vectors and matrices */ 938e333b035SStefano Zampini if (subi < 0) { idx[i] = subi - next->grstarts[rank]; continue; } 9396eb61c8cSJed Brown /* Binary search to find which rank owns subi */ 9406eb61c8cSJed Brown while (hi-lo > 1) { 9416eb61c8cSJed Brown t = lo + (hi-lo)/2; 9426eb61c8cSJed Brown if (suboff[t] > subi) hi = t; 9436eb61c8cSJed Brown else lo = t; 9446eb61c8cSJed Brown } 9456eb61c8cSJed Brown idx[i] = subi - suboff[lo] + next->grstarts[lo]; 9466eb61c8cSJed Brown } 9479566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltog,&indices)); 9489566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt])); 9499566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm,&global)); 95047c6ae99SBarry Smith next = next->next; 95147c6ae99SBarry Smith cnt++; 95247c6ae99SBarry Smith } 95347c6ae99SBarry Smith PetscFunctionReturn(0); 95447c6ae99SBarry Smith } 95547c6ae99SBarry Smith 95687c85e80SJed Brown /*@C 9579ae5db72SJed Brown DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector 95887c85e80SJed Brown 95987c85e80SJed Brown Not Collective 96087c85e80SJed Brown 9614165533cSJose E. Roman Input Parameter: 96287c85e80SJed Brown . dm - composite DM 96387c85e80SJed Brown 9644165533cSJose E. Roman Output Parameter: 96587c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite 96687c85e80SJed Brown 96787c85e80SJed Brown Level: intermediate 96887c85e80SJed Brown 96987c85e80SJed Brown Notes: 97087c85e80SJed Brown At present, a composite local vector does not normally exist. This function is used to provide index sets for 97187c85e80SJed Brown MatGetLocalSubMatrix(). In the future, the scatters for each entry in the DMComposite may be be merged into a single 9729ae5db72SJed Brown scatter to a composite local vector. The user should not typically need to know which is being done. 97387c85e80SJed Brown 97487c85e80SJed Brown To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings(). 97587c85e80SJed Brown 97687c85e80SJed Brown To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs(). 97787c85e80SJed Brown 97887c85e80SJed Brown Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree(). 97987c85e80SJed Brown 980f5f57ec0SBarry Smith Not available from Fortran 981f5f57ec0SBarry Smith 982db781477SPatrick Sanan .seealso: `DMCompositeGetGlobalISs()`, `DMCompositeGetISLocalToGlobalMappings()`, `MatGetLocalSubMatrix()`, `MatCreateLocalRef()` 98387c85e80SJed Brown @*/ 9847087cfbeSBarry Smith PetscErrorCode DMCompositeGetLocalISs(DM dm,IS **is) 98587c85e80SJed Brown { 98687c85e80SJed Brown DM_Composite *com = (DM_Composite*)dm->data; 98787c85e80SJed Brown struct DMCompositeLink *link; 98887c85e80SJed Brown PetscInt cnt,start; 98971b14b3eSStefano Zampini PetscBool flg; 99087c85e80SJed Brown 99187c85e80SJed Brown PetscFunctionBegin; 99287c85e80SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 99387c85e80SJed Brown PetscValidPointer(is,2); 9949566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 9957a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 9969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nmine,is)); 99706ebdd98SJed Brown for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) { 998520db06cSJed Brown PetscInt bs; 9999566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt])); 10009566063dSJacob Faibussowitsch PetscCall(DMGetBlockSize(link->dm,&bs)); 10019566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize((*is)[cnt],bs)); 1002520db06cSJed Brown } 100387c85e80SJed Brown PetscFunctionReturn(0); 100487c85e80SJed Brown } 100587c85e80SJed Brown 100647c6ae99SBarry Smith /*@C 100747c6ae99SBarry Smith DMCompositeGetGlobalISs - Gets the index sets for each composed object 100847c6ae99SBarry Smith 1009d083f849SBarry Smith Collective on dm 101047c6ae99SBarry Smith 101147c6ae99SBarry Smith Input Parameter: 101247c6ae99SBarry Smith . dm - the packer object 101347c6ae99SBarry Smith 101447c6ae99SBarry Smith Output Parameters: 101547c6ae99SBarry Smith . is - the array of index sets 101647c6ae99SBarry Smith 101747c6ae99SBarry Smith Level: advanced 101847c6ae99SBarry Smith 101947c6ae99SBarry Smith Notes: 102047c6ae99SBarry Smith The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree() 102147c6ae99SBarry Smith 102247c6ae99SBarry Smith These could be used to extract a subset of vector entries for a "multi-physics" preconditioner 102347c6ae99SBarry Smith 10246eb61c8cSJed Brown Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and 10256eb61c8cSJed Brown DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global 10266eb61c8cSJed Brown indices. 102747c6ae99SBarry Smith 1028f3cb0f7eSJed Brown Fortran Notes: 1029f3cb0f7eSJed Brown 1030f3cb0f7eSJed Brown The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM(). 1031f3cb0f7eSJed Brown 1032db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1033db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetAccess()`, `DMCompositeScatter()`, 1034c2e3fba1SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 103547c6ae99SBarry Smith 103647c6ae99SBarry Smith @*/ 10377087cfbeSBarry Smith PetscErrorCode DMCompositeGetGlobalISs(DM dm,IS *is[]) 103847c6ae99SBarry Smith { 103966bb578eSMark Adams PetscInt cnt = 0; 104047c6ae99SBarry Smith struct DMCompositeLink *next; 104147c6ae99SBarry Smith PetscMPIInt rank; 104247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 104371b14b3eSStefano Zampini PetscBool flg; 104447c6ae99SBarry Smith 104547c6ae99SBarry Smith PetscFunctionBegin; 104647c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 10479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 10487a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 10497a8be351SBarry Smith PetscCheck(dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() before"); 10509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(com->nDM,is)); 105147c6ae99SBarry Smith next = com->next; 10529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 105347c6ae99SBarry Smith 105447c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 105547c6ae99SBarry Smith while (next) { 1056e5e52638SMatthew G. Knepley PetscDS prob; 1057e5e52638SMatthew G. Knepley 10589566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt])); 10599566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 1060e5e52638SMatthew G. Knepley if (prob) { 106165c226d8SMatthew G. Knepley MatNullSpace space; 106265c226d8SMatthew G. Knepley Mat pmat; 10630f21e855SMatthew G. Knepley PetscObject disc; 10640f21e855SMatthew G. Knepley PetscInt Nf; 106565c226d8SMatthew G. Knepley 10669566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(prob, &Nf)); 1067f24dd8d2SMatthew G. Knepley if (cnt < Nf) { 10689566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(prob, cnt, &disc)); 10699566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nullspace", (PetscObject*) &space)); 10709566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space)); 10719566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space)); 10729566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space)); 10739566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat)); 10749566063dSJacob Faibussowitsch if (pmat) PetscCall(PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat)); 107565c226d8SMatthew G. Knepley } 1076f24dd8d2SMatthew G. Knepley } 107747c6ae99SBarry Smith cnt++; 107847c6ae99SBarry Smith next = next->next; 107947c6ae99SBarry Smith } 108047c6ae99SBarry Smith PetscFunctionReturn(0); 108147c6ae99SBarry Smith } 108247c6ae99SBarry Smith 108321c9b008SJed Brown PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields) 10844d343eeaSMatthew G Knepley { 10854d343eeaSMatthew G Knepley PetscInt nDM; 10864d343eeaSMatthew G Knepley DM *dms; 10874d343eeaSMatthew G Knepley PetscInt i; 10884d343eeaSMatthew G Knepley 10894d343eeaSMatthew G Knepley PetscFunctionBegin; 10909566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 10918865f1eaSKarl Rupp if (numFields) *numFields = nDM; 10929566063dSJacob Faibussowitsch PetscCall(DMCompositeGetGlobalISs(dm, fields)); 10934d343eeaSMatthew G Knepley if (fieldNames) { 10949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, &dms)); 10959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, fieldNames)); 10969566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, dms)); 10974d343eeaSMatthew G Knepley for (i=0; i<nDM; i++) { 10984d343eeaSMatthew G Knepley char buf[256]; 10994d343eeaSMatthew G Knepley const char *splitname; 11004d343eeaSMatthew G Knepley 11014d343eeaSMatthew G Knepley /* Split naming precedence: object name, prefix, number */ 11024d343eeaSMatthew G Knepley splitname = ((PetscObject) dm)->name; 11034d343eeaSMatthew G Knepley if (!splitname) { 11049566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname)); 11054d343eeaSMatthew G Knepley if (splitname) { 11064d343eeaSMatthew G Knepley size_t len; 11079566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(buf,splitname,sizeof(buf))); 11088caf3d72SBarry Smith buf[sizeof(buf) - 1] = 0; 11099566063dSJacob Faibussowitsch PetscCall(PetscStrlen(buf,&len)); 11104d343eeaSMatthew G Knepley if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */ 11114d343eeaSMatthew G Knepley splitname = buf; 11124d343eeaSMatthew G Knepley } 11134d343eeaSMatthew G Knepley } 11144d343eeaSMatthew G Knepley if (!splitname) { 111563a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf,sizeof(buf),"%" PetscInt_FMT,i)); 11164d343eeaSMatthew G Knepley splitname = buf; 11174d343eeaSMatthew G Knepley } 11189566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(splitname,&(*fieldNames)[i])); 11194d343eeaSMatthew G Knepley } 11209566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 11214d343eeaSMatthew G Knepley } 11224d343eeaSMatthew G Knepley PetscFunctionReturn(0); 11234d343eeaSMatthew G Knepley } 11244d343eeaSMatthew G Knepley 1125e7c4fc90SDmitry Karpeev /* 1126e7c4fc90SDmitry Karpeev This could take over from DMCreateFieldIS(), as it is more general, 11270298fd71SBarry Smith making DMCreateFieldIS() a special case -- calling with dmlist == NULL; 1128e7c4fc90SDmitry Karpeev At this point it's probably best to be less intrusive, however. 1129e7c4fc90SDmitry Karpeev */ 113016621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist) 1131e7c4fc90SDmitry Karpeev { 1132e7c4fc90SDmitry Karpeev PetscInt nDM; 1133e7c4fc90SDmitry Karpeev PetscInt i; 1134e7c4fc90SDmitry Karpeev 1135e7c4fc90SDmitry Karpeev PetscFunctionBegin; 11369566063dSJacob Faibussowitsch PetscCall(DMCreateFieldIS_Composite(dm, len, namelist, islist)); 1137e7c4fc90SDmitry Karpeev if (dmlist) { 11389566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm, &nDM)); 11399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDM, dmlist)); 11409566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm, *dmlist)); 1141e7c4fc90SDmitry Karpeev for (i=0; i<nDM; i++) { 11429566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)((*dmlist)[i]))); 1143e7c4fc90SDmitry Karpeev } 1144e7c4fc90SDmitry Karpeev } 1145e7c4fc90SDmitry Karpeev PetscFunctionReturn(0); 1146e7c4fc90SDmitry Karpeev } 1147e7c4fc90SDmitry Karpeev 114847c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 114947c6ae99SBarry Smith /*@C 11509ae5db72SJed Brown DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite. 115147c6ae99SBarry Smith Use DMCompositeRestoreLocalVectors() to return them. 115247c6ae99SBarry Smith 115347c6ae99SBarry Smith Not Collective 115447c6ae99SBarry Smith 115547c6ae99SBarry Smith Input Parameter: 115647c6ae99SBarry Smith . dm - the packer object 115747c6ae99SBarry Smith 115847c6ae99SBarry Smith Output Parameter: 11599ae5db72SJed Brown . Vec ... - the individual sequential Vecs 116047c6ae99SBarry Smith 116147c6ae99SBarry Smith Level: advanced 116247c6ae99SBarry Smith 1163f5f57ec0SBarry Smith Not available from Fortran 1164f5f57ec0SBarry Smith 1165db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1166db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1167db781477SPatrick Sanan `DMCompositeRestoreLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 116847c6ae99SBarry Smith 116947c6ae99SBarry Smith @*/ 11707087cfbeSBarry Smith PetscErrorCode DMCompositeGetLocalVectors(DM dm,...) 117147c6ae99SBarry Smith { 117247c6ae99SBarry Smith va_list Argp; 117347c6ae99SBarry Smith struct DMCompositeLink *next; 117447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 117571b14b3eSStefano Zampini PetscBool flg; 117647c6ae99SBarry Smith 117747c6ae99SBarry Smith PetscFunctionBegin; 117847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 11799566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 11807a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 118147c6ae99SBarry Smith next = com->next; 118247c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 118347c6ae99SBarry Smith va_start(Argp,dm); 118447c6ae99SBarry Smith while (next) { 118547c6ae99SBarry Smith Vec *vec; 118647c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 11879566063dSJacob Faibussowitsch if (vec) PetscCall(DMGetLocalVector(next->dm,vec)); 118847c6ae99SBarry Smith next = next->next; 118947c6ae99SBarry Smith } 119047c6ae99SBarry Smith va_end(Argp); 119147c6ae99SBarry Smith PetscFunctionReturn(0); 119247c6ae99SBarry Smith } 119347c6ae99SBarry Smith 119447c6ae99SBarry Smith /*@C 11959ae5db72SJed Brown DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite. 119647c6ae99SBarry Smith 119747c6ae99SBarry Smith Not Collective 119847c6ae99SBarry Smith 119947c6ae99SBarry Smith Input Parameter: 120047c6ae99SBarry Smith . dm - the packer object 120147c6ae99SBarry Smith 120247c6ae99SBarry Smith Output Parameter: 12039ae5db72SJed Brown . Vec ... - the individual sequential Vecs 120447c6ae99SBarry Smith 120547c6ae99SBarry Smith Level: advanced 120647c6ae99SBarry Smith 1207f5f57ec0SBarry Smith Not available from Fortran 1208f5f57ec0SBarry Smith 1209db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, 1210db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1211db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, `DMCompositeGetEntries()` 121247c6ae99SBarry Smith 121347c6ae99SBarry Smith @*/ 12147087cfbeSBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors(DM dm,...) 121547c6ae99SBarry Smith { 121647c6ae99SBarry Smith va_list Argp; 121747c6ae99SBarry Smith struct DMCompositeLink *next; 121847c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 121971b14b3eSStefano Zampini PetscBool flg; 122047c6ae99SBarry Smith 122147c6ae99SBarry Smith PetscFunctionBegin; 122247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 12239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 12247a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 122547c6ae99SBarry Smith next = com->next; 122647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 122747c6ae99SBarry Smith va_start(Argp,dm); 122847c6ae99SBarry Smith while (next) { 122947c6ae99SBarry Smith Vec *vec; 123047c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 12319566063dSJacob Faibussowitsch if (vec) PetscCall(DMRestoreLocalVector(next->dm,vec)); 123247c6ae99SBarry Smith next = next->next; 123347c6ae99SBarry Smith } 123447c6ae99SBarry Smith va_end(Argp); 123547c6ae99SBarry Smith PetscFunctionReturn(0); 123647c6ae99SBarry Smith } 123747c6ae99SBarry Smith 123847c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 123947c6ae99SBarry Smith /*@C 12409ae5db72SJed Brown DMCompositeGetEntries - Gets the DM for each entry in a DMComposite. 124147c6ae99SBarry Smith 124247c6ae99SBarry Smith Not Collective 124347c6ae99SBarry Smith 124447c6ae99SBarry Smith Input Parameter: 124547c6ae99SBarry Smith . dm - the packer object 124647c6ae99SBarry Smith 124747c6ae99SBarry Smith Output Parameter: 12489ae5db72SJed Brown . DM ... - the individual entries (DMs) 124947c6ae99SBarry Smith 125047c6ae99SBarry Smith Level: advanced 125147c6ae99SBarry Smith 125295452b02SPatrick Sanan Fortran Notes: 125395452b02SPatrick Sanan Available as DMCompositeGetEntries() for one output DM, DMCompositeGetEntries2() for 2, etc 1254f5f57ec0SBarry Smith 1255db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntriesArray()` 1256db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1257db781477SPatrick Sanan `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, 1258db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()` 125947c6ae99SBarry Smith 126047c6ae99SBarry Smith @*/ 12617087cfbeSBarry Smith PetscErrorCode DMCompositeGetEntries(DM dm,...) 126247c6ae99SBarry Smith { 126347c6ae99SBarry Smith va_list Argp; 126447c6ae99SBarry Smith struct DMCompositeLink *next; 126547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 126671b14b3eSStefano Zampini PetscBool flg; 126747c6ae99SBarry Smith 126847c6ae99SBarry Smith PetscFunctionBegin; 126947c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 12709566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 12717a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 127247c6ae99SBarry Smith next = com->next; 127347c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 127447c6ae99SBarry Smith va_start(Argp,dm); 127547c6ae99SBarry Smith while (next) { 127647c6ae99SBarry Smith DM *dmn; 127747c6ae99SBarry Smith dmn = va_arg(Argp, DM*); 12789ae5db72SJed Brown if (dmn) *dmn = next->dm; 127947c6ae99SBarry Smith next = next->next; 128047c6ae99SBarry Smith } 128147c6ae99SBarry Smith va_end(Argp); 128247c6ae99SBarry Smith PetscFunctionReturn(0); 128347c6ae99SBarry Smith } 128447c6ae99SBarry Smith 1285dbab29e1SMark F. Adams /*@C 12862fa5ba8aSJed Brown DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite. 12872fa5ba8aSJed Brown 12882fa5ba8aSJed Brown Not Collective 12892fa5ba8aSJed Brown 12902fa5ba8aSJed Brown Input Parameter: 1291907376e6SBarry Smith . dm - the packer object 1292907376e6SBarry Smith 1293907376e6SBarry Smith Output Parameter: 1294907376e6SBarry Smith . dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs 12952fa5ba8aSJed Brown 12962fa5ba8aSJed Brown Level: advanced 12972fa5ba8aSJed Brown 1298db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCreateGlobalVector()`, `DMCompositeGetEntries()` 1299db781477SPatrick Sanan `DMCompositeGather()`, `DMCompositeCreate()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()`, 1300db781477SPatrick Sanan `DMCompositeRestoreLocalVectors()`, `DMCompositeGetLocalVectors()`, `DMCompositeScatter()`, 1301db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()` 13022fa5ba8aSJed Brown 13032fa5ba8aSJed Brown @*/ 13042fa5ba8aSJed Brown PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[]) 13052fa5ba8aSJed Brown { 13062fa5ba8aSJed Brown struct DMCompositeLink *next; 13072fa5ba8aSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 13082fa5ba8aSJed Brown PetscInt i; 130971b14b3eSStefano Zampini PetscBool flg; 13102fa5ba8aSJed Brown 13112fa5ba8aSJed Brown PetscFunctionBegin; 13122fa5ba8aSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 13139566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg)); 13147a8be351SBarry Smith PetscCheck(flg,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 13152fa5ba8aSJed Brown /* loop over packed objects, handling one at at time */ 13162fa5ba8aSJed Brown for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm; 13172fa5ba8aSJed Brown PetscFunctionReturn(0); 13182fa5ba8aSJed Brown } 13192fa5ba8aSJed Brown 1320e10fd815SStefano Zampini typedef struct { 1321e10fd815SStefano Zampini DM dm; 1322e10fd815SStefano Zampini PetscViewer *subv; 1323e10fd815SStefano Zampini Vec *vecs; 1324e10fd815SStefano Zampini } GLVisViewerCtx; 1325e10fd815SStefano Zampini 1326e10fd815SStefano Zampini static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 1327e10fd815SStefano Zampini { 1328e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 1329e10fd815SStefano Zampini PetscInt i,n; 1330e10fd815SStefano Zampini 1331e10fd815SStefano Zampini PetscFunctionBegin; 13329566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm,&n)); 1333e10fd815SStefano Zampini for (i = 0; i < n; i++) { 13349566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&ctx->subv[i])); 1335e10fd815SStefano Zampini } 13369566063dSJacob Faibussowitsch PetscCall(PetscFree2(ctx->subv,ctx->vecs)); 13379566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx->dm)); 13389566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 1339e10fd815SStefano Zampini PetscFunctionReturn(0); 1340e10fd815SStefano Zampini } 1341e10fd815SStefano Zampini 1342e10fd815SStefano Zampini static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 1343e10fd815SStefano Zampini { 1344e10fd815SStefano Zampini Vec X = (Vec)oX; 1345e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 1346e10fd815SStefano Zampini PetscInt i,n,cumf; 1347e10fd815SStefano Zampini 1348e10fd815SStefano Zampini PetscFunctionBegin; 13499566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(ctx->dm,&n)); 13509566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(ctx->dm,X,n,NULL,ctx->vecs)); 1351e10fd815SStefano Zampini for (i = 0, cumf = 0; i < n; i++) { 1352e10fd815SStefano Zampini PetscErrorCode (*g2l)(PetscObject,PetscInt,PetscObject[],void*); 1353e10fd815SStefano Zampini void *fctx; 1354e10fd815SStefano Zampini PetscInt nfi; 1355e10fd815SStefano Zampini 13569566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i],&nfi,NULL,NULL,&g2l,NULL,&fctx)); 1357e10fd815SStefano Zampini if (!nfi) continue; 13581baa6e33SBarry Smith if (g2l) PetscCall((*g2l)((PetscObject)ctx->vecs[i],nfi,oXfield+cumf,fctx)); 13591baa6e33SBarry Smith else PetscCall(VecCopy(ctx->vecs[i],(Vec)(oXfield[cumf]))); 1360e10fd815SStefano Zampini cumf += nfi; 1361e10fd815SStefano Zampini } 13629566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(ctx->dm,X,n,NULL,ctx->vecs)); 1363e10fd815SStefano Zampini PetscFunctionReturn(0); 1364e10fd815SStefano Zampini } 1365e10fd815SStefano Zampini 1366e10fd815SStefano Zampini static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer) 1367e10fd815SStefano Zampini { 1368e10fd815SStefano Zampini DM dm = (DM)odm, *dms; 1369e10fd815SStefano Zampini Vec *Ufds; 1370e10fd815SStefano Zampini GLVisViewerCtx *ctx; 1371e10fd815SStefano Zampini PetscInt i,n,tnf,*sdim; 1372e10fd815SStefano Zampini char **fecs; 1373e10fd815SStefano Zampini 1374e10fd815SStefano Zampini PetscFunctionBegin; 13759566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 13769566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 1377e10fd815SStefano Zampini ctx->dm = dm; 13789566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(dm,&n)); 13799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&dms)); 13809566063dSJacob Faibussowitsch PetscCall(DMCompositeGetEntriesArray(dm,dms)); 13819566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n,&ctx->subv,n,&ctx->vecs)); 1382e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1383e10fd815SStefano Zampini PetscInt nf; 1384e10fd815SStefano Zampini 13859566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PetscObjectComm(odm),&ctx->subv[i])); 13869566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(ctx->subv[i],PETSCVIEWERGLVIS)); 13879566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetDM_Private(ctx->subv[i],(PetscObject)dms[i])); 13889566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,NULL,NULL,NULL,NULL,NULL)); 1389e10fd815SStefano Zampini tnf += nf; 1390e10fd815SStefano Zampini } 13919566063dSJacob Faibussowitsch PetscCall(PetscFree(dms)); 13929566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(tnf,&fecs,tnf,&sdim,tnf,&Ufds)); 1393e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1394e10fd815SStefano Zampini PetscInt *sd,nf,f; 1395e10fd815SStefano Zampini const char **fec; 1396e10fd815SStefano Zampini Vec *Uf; 1397e10fd815SStefano Zampini 13989566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,&fec,&sd,NULL,(PetscObject**)&Uf,NULL)); 1399e10fd815SStefano Zampini for (f = 0; f < nf; f++) { 14009566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec[f],&fecs[tnf+f])); 1401e10fd815SStefano Zampini Ufds[tnf+f] = Uf[f]; 1402e10fd815SStefano Zampini sdim[tnf+f] = sd[f]; 1403e10fd815SStefano Zampini } 1404e10fd815SStefano Zampini tnf += nf; 1405e10fd815SStefano Zampini } 14069566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetFields(viewer,tnf,(const char**)fecs,sdim,DMCompositeSampleGLVisFields_Private,(PetscObject*)Ufds,ctx,DestroyGLVisViewerCtx_Private)); 1407e10fd815SStefano Zampini for (i = 0; i < tnf; i++) { 14089566063dSJacob Faibussowitsch PetscCall(PetscFree(fecs[i])); 1409e10fd815SStefano Zampini } 14109566063dSJacob Faibussowitsch PetscCall(PetscFree3(fecs,sdim,Ufds)); 1411e10fd815SStefano Zampini PetscFunctionReturn(0); 1412e10fd815SStefano Zampini } 1413e10fd815SStefano Zampini 14147087cfbeSBarry Smith PetscErrorCode DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine) 141547c6ae99SBarry Smith { 141647c6ae99SBarry Smith struct DMCompositeLink *next; 141747c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dmi->data; 141847c6ae99SBarry Smith DM dm; 141947c6ae99SBarry Smith 142047c6ae99SBarry Smith PetscFunctionBegin; 142147c6ae99SBarry Smith PetscValidHeaderSpecific(dmi,DM_CLASSID,1); 1422ce94432eSBarry Smith if (comm == MPI_COMM_NULL) { 14239566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dmi,&comm)); 1424ce94432eSBarry Smith } 14259566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi)); 142647c6ae99SBarry Smith next = com->next; 14279566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm,fine)); 142847c6ae99SBarry Smith 142947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 143047c6ae99SBarry Smith while (next) { 14319566063dSJacob Faibussowitsch PetscCall(DMRefine(next->dm,comm,&dm)); 14329566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine,dm)); 14339566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 143447c6ae99SBarry Smith next = next->next; 143547c6ae99SBarry Smith } 143647c6ae99SBarry Smith PetscFunctionReturn(0); 143747c6ae99SBarry Smith } 143847c6ae99SBarry Smith 143914354c39SJed Brown PetscErrorCode DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine) 144014354c39SJed Brown { 144114354c39SJed Brown struct DMCompositeLink *next; 144214354c39SJed Brown DM_Composite *com = (DM_Composite*)dmi->data; 144314354c39SJed Brown DM dm; 144414354c39SJed Brown 144514354c39SJed Brown PetscFunctionBegin; 144614354c39SJed Brown PetscValidHeaderSpecific(dmi,DM_CLASSID,1); 14479566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmi)); 14482ee06e3bSJed Brown if (comm == MPI_COMM_NULL) { 14499566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dmi,&comm)); 145025296bd5SBarry Smith } 145114354c39SJed Brown next = com->next; 14529566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(comm,fine)); 145314354c39SJed Brown 145414354c39SJed Brown /* loop over packed objects, handling one at at time */ 145514354c39SJed Brown while (next) { 14569566063dSJacob Faibussowitsch PetscCall(DMCoarsen(next->dm,comm,&dm)); 14579566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(*fine,dm)); 14589566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 145914354c39SJed Brown next = next->next; 146014354c39SJed Brown } 146114354c39SJed Brown PetscFunctionReturn(0); 146214354c39SJed Brown } 146347c6ae99SBarry Smith 1464e727c939SJed Brown PetscErrorCode DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v) 146547c6ae99SBarry Smith { 14669ae5db72SJed Brown PetscInt m,n,M,N,nDM,i; 146747c6ae99SBarry Smith struct DMCompositeLink *nextc; 146847c6ae99SBarry Smith struct DMCompositeLink *nextf; 146925296bd5SBarry Smith Vec gcoarse,gfine,*vecs; 147047c6ae99SBarry Smith DM_Composite *comcoarse = (DM_Composite*)coarse->data; 147147c6ae99SBarry Smith DM_Composite *comfine = (DM_Composite*)fine->data; 14729ae5db72SJed Brown Mat *mats; 147347c6ae99SBarry Smith 147447c6ae99SBarry Smith PetscFunctionBegin; 147547c6ae99SBarry Smith PetscValidHeaderSpecific(coarse,DM_CLASSID,1); 147647c6ae99SBarry Smith PetscValidHeaderSpecific(fine,DM_CLASSID,2); 14779566063dSJacob Faibussowitsch PetscCall(DMSetUp(coarse)); 14789566063dSJacob Faibussowitsch PetscCall(DMSetUp(fine)); 147947c6ae99SBarry Smith /* use global vectors only for determining matrix layout */ 14809566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(coarse,&gcoarse)); 14819566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fine,&gfine)); 14829566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gcoarse,&n)); 14839566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(gfine,&m)); 14849566063dSJacob Faibussowitsch PetscCall(VecGetSize(gcoarse,&N)); 14859566063dSJacob Faibussowitsch PetscCall(VecGetSize(gfine,&M)); 14869566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(coarse,&gcoarse)); 14879566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fine,&gfine)); 148847c6ae99SBarry Smith 14899ae5db72SJed Brown nDM = comfine->nDM; 149063a3b9bcSJacob 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); 14919566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nDM*nDM,&mats)); 149225296bd5SBarry Smith if (v) { 14939566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nDM,&vecs)); 149425296bd5SBarry Smith } 149547c6ae99SBarry Smith 149647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 14979ae5db72SJed Brown for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) { 14981baa6e33SBarry Smith if (!v) PetscCall(DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL)); 14991baa6e33SBarry Smith else PetscCall(DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i])); 150047c6ae99SBarry Smith } 15019566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A)); 15021baa6e33SBarry Smith if (v) PetscCall(VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v)); 15039566063dSJacob Faibussowitsch for (i=0; i<nDM*nDM; i++) PetscCall(MatDestroy(&mats[i])); 15049566063dSJacob Faibussowitsch PetscCall(PetscFree(mats)); 150525296bd5SBarry Smith if (v) { 15069566063dSJacob Faibussowitsch for (i=0; i<nDM; i++) PetscCall(VecDestroy(&vecs[i])); 15079566063dSJacob Faibussowitsch PetscCall(PetscFree(vecs)); 150825296bd5SBarry Smith } 150947c6ae99SBarry Smith PetscFunctionReturn(0); 151047c6ae99SBarry Smith } 151147c6ae99SBarry Smith 1512184d77edSJed Brown static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm) 15131411c6eeSJed Brown { 15141411c6eeSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 15151411c6eeSJed Brown ISLocalToGlobalMapping *ltogs; 1516f7efa3c7SJed Brown PetscInt i; 15171411c6eeSJed Brown 15181411c6eeSJed Brown PetscFunctionBegin; 15191411c6eeSJed Brown /* Set the ISLocalToGlobalMapping on the new matrix */ 15209566063dSJacob Faibussowitsch PetscCall(DMCompositeGetISLocalToGlobalMappings(dm,<ogs)); 15219566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap)); 15229566063dSJacob Faibussowitsch for (i=0; i<com->nDM; i++) PetscCall(ISLocalToGlobalMappingDestroy(<ogs[i])); 15239566063dSJacob Faibussowitsch PetscCall(PetscFree(ltogs)); 15241411c6eeSJed Brown PetscFunctionReturn(0); 15251411c6eeSJed Brown } 15261411c6eeSJed Brown 1527b412c318SBarry Smith PetscErrorCode DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring) 152847c6ae99SBarry Smith { 152947c6ae99SBarry Smith PetscInt n,i,cnt; 153047c6ae99SBarry Smith ISColoringValue *colors; 153147c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 153247c6ae99SBarry Smith ISColoringValue maxcol = 0; 153347c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 153447c6ae99SBarry Smith 153547c6ae99SBarry Smith PetscFunctionBegin; 153647c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 153708401ef6SPierre Jolivet PetscCheck(ctype != IS_COLORING_LOCAL,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported"); 1538e3247f34SBarry Smith else if (ctype == IS_COLORING_GLOBAL) { 153947c6ae99SBarry Smith n = com->n; 1540ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType"); 15419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&colors)); /* freed in ISColoringDestroy() */ 154247c6ae99SBarry Smith 15439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL)); 154447c6ae99SBarry Smith if (dense) { 154547c6ae99SBarry Smith for (i=0; i<n; i++) { 154647c6ae99SBarry Smith colors[i] = (ISColoringValue)(com->rstart + i); 154747c6ae99SBarry Smith } 154847c6ae99SBarry Smith maxcol = com->N; 154947c6ae99SBarry Smith } else { 155047c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 155147c6ae99SBarry Smith PetscMPIInt rank; 155247c6ae99SBarry Smith 15539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 155447c6ae99SBarry Smith cnt = 0; 155547c6ae99SBarry Smith while (next) { 155647c6ae99SBarry Smith ISColoring lcoloring; 155747c6ae99SBarry Smith 15589566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring)); 155947c6ae99SBarry Smith for (i=0; i<lcoloring->N; i++) { 156047c6ae99SBarry Smith colors[cnt++] = maxcol + lcoloring->colors[i]; 156147c6ae99SBarry Smith } 156247c6ae99SBarry Smith maxcol += lcoloring->n; 15639566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&lcoloring)); 156447c6ae99SBarry Smith next = next->next; 156547c6ae99SBarry Smith } 156647c6ae99SBarry Smith } 15679566063dSJacob Faibussowitsch PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,PETSC_OWN_POINTER,coloring)); 156847c6ae99SBarry Smith PetscFunctionReturn(0); 156947c6ae99SBarry Smith } 157047c6ae99SBarry Smith 15717087cfbeSBarry Smith PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec) 157247c6ae99SBarry Smith { 157347c6ae99SBarry Smith struct DMCompositeLink *next; 157447c6ae99SBarry Smith PetscScalar *garray,*larray; 157547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 157647c6ae99SBarry Smith 157747c6ae99SBarry Smith PetscFunctionBegin; 157847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 157947c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 158039d35262SVincent Le Chenadec 158147c6ae99SBarry Smith if (!com->setup) { 15829566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 158347c6ae99SBarry Smith } 158439d35262SVincent Le Chenadec 15859566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec,&garray)); 15869566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec,&larray)); 158747c6ae99SBarry Smith 158847c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 158939d35262SVincent Le Chenadec next = com->next; 159047c6ae99SBarry Smith while (next) { 159147c6ae99SBarry Smith Vec local,global; 159247c6ae99SBarry Smith PetscInt N; 159347c6ae99SBarry Smith 15949566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm,&global)); 15959566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global,&N)); 15969566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global,garray)); 15979566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm,&local)); 15989566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local,larray)); 15999566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(next->dm,global,mode,local)); 16009566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(next->dm,global,mode,local)); 16019566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 16029566063dSJacob Faibussowitsch PetscCall(VecResetArray(local)); 16039566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm,&global)); 16049566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm,&local)); 160539d35262SVincent Le Chenadec 160606ebdd98SJed Brown larray += next->nlocal; 160739d35262SVincent Le Chenadec garray += next->n; 160847c6ae99SBarry Smith next = next->next; 160947c6ae99SBarry Smith } 161047c6ae99SBarry Smith 16119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec,NULL)); 16129566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lvec,NULL)); 161347c6ae99SBarry Smith PetscFunctionReturn(0); 161447c6ae99SBarry Smith } 161547c6ae99SBarry Smith 16167087cfbeSBarry Smith PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec) 16170c010503SBarry Smith { 16180c010503SBarry Smith PetscFunctionBegin; 161939d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 162039d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 162139d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,4); 162239d35262SVincent Le Chenadec PetscFunctionReturn(0); 162339d35262SVincent Le Chenadec } 162439d35262SVincent Le Chenadec 162539d35262SVincent Le Chenadec PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec) 162639d35262SVincent Le Chenadec { 162739d35262SVincent Le Chenadec struct DMCompositeLink *next; 162839d35262SVincent Le Chenadec PetscScalar *larray,*garray; 162939d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite*)dm->data; 163039d35262SVincent Le Chenadec 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 163639d35262SVincent Le Chenadec if (!com->setup) { 16379566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 163839d35262SVincent Le Chenadec } 163939d35262SVincent Le Chenadec 16409566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec,&larray)); 16419566063dSJacob Faibussowitsch PetscCall(VecGetArray(gvec,&garray)); 164239d35262SVincent Le Chenadec 164339d35262SVincent Le Chenadec /* loop over packed objects, handling one at at time */ 164439d35262SVincent Le Chenadec next = com->next; 164539d35262SVincent Le Chenadec while (next) { 164639d35262SVincent Le Chenadec Vec global,local; 164739d35262SVincent Le Chenadec 16489566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm,&local)); 16499566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local,larray)); 16509566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(next->dm,&global)); 16519566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(global,garray)); 16529566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(next->dm,local,mode,global)); 16539566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(next->dm,local,mode,global)); 16549566063dSJacob Faibussowitsch PetscCall(VecResetArray(local)); 16559566063dSJacob Faibussowitsch PetscCall(VecResetArray(global)); 16569566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(next->dm,&global)); 16579566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm,&local)); 165839d35262SVincent Le Chenadec 165939d35262SVincent Le Chenadec garray += next->n; 166039d35262SVincent Le Chenadec larray += next->nlocal; 166139d35262SVincent Le Chenadec next = next->next; 166239d35262SVincent Le Chenadec } 166339d35262SVincent Le Chenadec 16649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gvec,NULL)); 16659566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lvec,NULL)); 166639d35262SVincent Le Chenadec 166739d35262SVincent Le Chenadec PetscFunctionReturn(0); 166839d35262SVincent Le Chenadec } 166939d35262SVincent Le Chenadec 167039d35262SVincent Le Chenadec PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec) 167139d35262SVincent Le Chenadec { 167239d35262SVincent Le Chenadec PetscFunctionBegin; 167339d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 167439d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,2); 167539d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,4); 167639d35262SVincent Le Chenadec PetscFunctionReturn(0); 167739d35262SVincent Le Chenadec } 167839d35262SVincent Le Chenadec 167939d35262SVincent Le Chenadec PetscErrorCode DMLocalToLocalBegin_Composite(DM dm,Vec vec1,InsertMode mode,Vec vec2) 168039d35262SVincent Le Chenadec { 168139d35262SVincent Le Chenadec struct DMCompositeLink *next; 168239d35262SVincent Le Chenadec PetscScalar *array1,*array2; 168339d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite*)dm->data; 168439d35262SVincent Le Chenadec 168539d35262SVincent Le Chenadec PetscFunctionBegin; 168639d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 168739d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec1,VEC_CLASSID,2); 168839d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec2,VEC_CLASSID,4); 168939d35262SVincent Le Chenadec 169039d35262SVincent Le Chenadec if (!com->setup) { 16919566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 169239d35262SVincent Le Chenadec } 169339d35262SVincent Le Chenadec 16949566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec1,&array1)); 16959566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec2,&array2)); 169639d35262SVincent Le Chenadec 169739d35262SVincent Le Chenadec /* loop over packed objects, handling one at at time */ 169839d35262SVincent Le Chenadec next = com->next; 169939d35262SVincent Le Chenadec while (next) { 170039d35262SVincent Le Chenadec Vec local1,local2; 170139d35262SVincent Le Chenadec 17029566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm,&local1)); 17039566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local1,array1)); 17049566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(next->dm,&local2)); 17059566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(local2,array2)); 17069566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalBegin(next->dm,local1,mode,local2)); 17079566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalEnd(next->dm,local1,mode,local2)); 17089566063dSJacob Faibussowitsch PetscCall(VecResetArray(local2)); 17099566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm,&local2)); 17109566063dSJacob Faibussowitsch PetscCall(VecResetArray(local1)); 17119566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(next->dm,&local1)); 171239d35262SVincent Le Chenadec 171339d35262SVincent Le Chenadec array1 += next->nlocal; 171439d35262SVincent Le Chenadec array2 += next->nlocal; 171539d35262SVincent Le Chenadec next = next->next; 171639d35262SVincent Le Chenadec } 171739d35262SVincent Le Chenadec 17189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec1,NULL)); 17199566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec2,NULL)); 172039d35262SVincent Le Chenadec 172139d35262SVincent Le Chenadec PetscFunctionReturn(0); 172239d35262SVincent Le Chenadec } 172339d35262SVincent Le Chenadec 172439d35262SVincent Le Chenadec PetscErrorCode DMLocalToLocalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec) 172539d35262SVincent Le Chenadec { 172639d35262SVincent Le Chenadec PetscFunctionBegin; 172739d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 172839d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,2); 172939d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,4); 17300c010503SBarry Smith PetscFunctionReturn(0); 17310c010503SBarry Smith } 173247c6ae99SBarry Smith 17336ae3a549SBarry Smith /*MC 17346ae3a549SBarry Smith DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs 17356ae3a549SBarry Smith 17366ae3a549SBarry Smith Level: intermediate 17376ae3a549SBarry Smith 1738db781477SPatrick Sanan .seealso: `DMType`, `DM`, `DMDACreate()`, `DMCreate()`, `DMSetType()`, `DMCompositeCreate()` 17396ae3a549SBarry Smith M*/ 17406ae3a549SBarry Smith 17418cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p) 1742a4121054SBarry Smith { 1743a4121054SBarry Smith DM_Composite *com; 1744a4121054SBarry Smith 1745a4121054SBarry Smith PetscFunctionBegin; 17469566063dSJacob Faibussowitsch PetscCall(PetscNewLog(p,&com)); 1747a4121054SBarry Smith p->data = com; 1748a4121054SBarry Smith com->n = 0; 17497ac2b803SAlex Fikl com->nghost = 0; 17500298fd71SBarry Smith com->next = NULL; 1751a4121054SBarry Smith com->nDM = 0; 1752a4121054SBarry Smith 1753a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Composite; 1754a4121054SBarry Smith p->ops->createlocalvector = DMCreateLocalVector_Composite; 1755184d77edSJed Brown p->ops->getlocaltoglobalmapping = DMGetLocalToGlobalMapping_Composite; 17564d343eeaSMatthew G Knepley p->ops->createfieldis = DMCreateFieldIS_Composite; 175716621825SDmitry Karpeev p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite; 1758a4121054SBarry Smith p->ops->refine = DMRefine_Composite; 175914354c39SJed Brown p->ops->coarsen = DMCoarsen_Composite; 176025296bd5SBarry Smith p->ops->createinterpolation = DMCreateInterpolation_Composite; 176125296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Composite; 1762e727c939SJed Brown p->ops->getcoloring = DMCreateColoring_Composite; 1763a4121054SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite; 1764a4121054SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite; 176539d35262SVincent Le Chenadec p->ops->localtoglobalbegin = DMLocalToGlobalBegin_Composite; 176639d35262SVincent Le Chenadec p->ops->localtoglobalend = DMLocalToGlobalEnd_Composite; 176739d35262SVincent Le Chenadec p->ops->localtolocalbegin = DMLocalToLocalBegin_Composite; 176839d35262SVincent Le Chenadec p->ops->localtolocalend = DMLocalToLocalEnd_Composite; 1769a4121054SBarry Smith p->ops->destroy = DMDestroy_Composite; 1770a4121054SBarry Smith p->ops->view = DMView_Composite; 1771a4121054SBarry Smith p->ops->setup = DMSetUp_Composite; 1772e10fd815SStefano Zampini 17739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)p,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Composite)); 1774a4121054SBarry Smith PetscFunctionReturn(0); 1775a4121054SBarry Smith } 1776a4121054SBarry Smith 17771fd49c25SBarry Smith /*@ 17780c010503SBarry Smith DMCompositeCreate - Creates a vector packer, used to generate "composite" 17790c010503SBarry Smith vectors made up of several subvectors. 17800c010503SBarry Smith 1781d083f849SBarry Smith Collective 178247c6ae99SBarry Smith 178347c6ae99SBarry Smith Input Parameter: 17840c010503SBarry Smith . comm - the processors that will share the global vector 17850c010503SBarry Smith 17860c010503SBarry Smith Output Parameters: 17870c010503SBarry Smith . packer - the packer object 178847c6ae99SBarry Smith 178947c6ae99SBarry Smith Level: advanced 179047c6ae99SBarry Smith 1791c2e3fba1SPatrick Sanan .seealso `DMDestroy()`, `DMCompositeAddDM()`, `DMCompositeScatter()`, `DMCOMPOSITE`, `DMCreate()` 1792db781477SPatrick Sanan `DMCompositeGather()`, `DMCreateGlobalVector()`, `DMCompositeGetISLocalToGlobalMappings()`, `DMCompositeGetAccess()` 1793db781477SPatrick Sanan `DMCompositeGetLocalVectors()`, `DMCompositeRestoreLocalVectors()`, `DMCompositeGetEntries()` 179447c6ae99SBarry Smith 179547c6ae99SBarry Smith @*/ 17967087cfbeSBarry Smith PetscErrorCode DMCompositeCreate(MPI_Comm comm,DM *packer) 179747c6ae99SBarry Smith { 179847c6ae99SBarry Smith PetscFunctionBegin; 17990c010503SBarry Smith PetscValidPointer(packer,2); 18009566063dSJacob Faibussowitsch PetscCall(DMCreate(comm,packer)); 18019566063dSJacob Faibussowitsch PetscCall(DMSetType(*packer,DMCOMPOSITE)); 180247c6ae99SBarry Smith PetscFunctionReturn(0); 180347c6ae99SBarry Smith } 1804