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 1147c6ae99SBarry Smith 1247c6ae99SBarry Smith Logically Collective on MPI_Comm 1347c6ae99SBarry Smith 1447c6ae99SBarry Smith Input Parameter: 1547c6ae99SBarry Smith + dm - the composite object 1647c6ae99SBarry Smith - formcouplelocations - routine to set the nonzero locations in the matrix 1747c6ae99SBarry Smith 1847c6ae99SBarry Smith Level: advanced 1947c6ae99SBarry Smith 20f5f57ec0SBarry Smith Not available from Fortran 21f5f57ec0SBarry Smith 221b2093e4SBarry Smith Notes: 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; 29*71b14b3eSStefano Zampini PetscBool flg; 30*71b14b3eSStefano Zampini PetscErrorCode ierr; 3147c6ae99SBarry Smith 3247c6ae99SBarry Smith PetscFunctionBegin; 33*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 34*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 3547c6ae99SBarry Smith com->FormCoupleLocations = FormCoupleLocations; 3647c6ae99SBarry Smith PetscFunctionReturn(0); 3747c6ae99SBarry Smith } 3847c6ae99SBarry Smith 396bf464f9SBarry Smith PetscErrorCode DMDestroy_Composite(DM dm) 4047c6ae99SBarry Smith { 4147c6ae99SBarry Smith PetscErrorCode ierr; 4247c6ae99SBarry Smith struct DMCompositeLink *next, *prev; 4347c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 4447c6ae99SBarry Smith 4547c6ae99SBarry Smith PetscFunctionBegin; 4647c6ae99SBarry Smith next = com->next; 4747c6ae99SBarry Smith while (next) { 4847c6ae99SBarry Smith prev = next; 4947c6ae99SBarry Smith next = next->next; 50fcfd50ebSBarry Smith ierr = DMDestroy(&prev->dm);CHKERRQ(ierr); 5147c6ae99SBarry Smith ierr = PetscFree(prev->grstarts);CHKERRQ(ierr); 5247c6ae99SBarry Smith ierr = PetscFree(prev);CHKERRQ(ierr); 5347c6ae99SBarry Smith } 54e10fd815SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL);CHKERRQ(ierr); 55435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 56435a35e8SMatthew G Knepley ierr = PetscFree(com);CHKERRQ(ierr); 5747c6ae99SBarry Smith PetscFunctionReturn(0); 5847c6ae99SBarry Smith } 5947c6ae99SBarry Smith 607087cfbeSBarry Smith PetscErrorCode DMView_Composite(DM dm,PetscViewer v) 6147c6ae99SBarry Smith { 6247c6ae99SBarry Smith PetscErrorCode ierr; 6347c6ae99SBarry Smith PetscBool iascii; 6447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 6547c6ae99SBarry Smith 6647c6ae99SBarry Smith PetscFunctionBegin; 67251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 6847c6ae99SBarry Smith if (iascii) { 6947c6ae99SBarry Smith struct DMCompositeLink *lnk = com->next; 7047c6ae99SBarry Smith PetscInt i; 7147c6ae99SBarry Smith 7247c6ae99SBarry Smith ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix");CHKERRQ(ierr); 739ae5db72SJed Brown ierr = PetscViewerASCIIPrintf(v," contains %D DMs\n",com->nDM);CHKERRQ(ierr); 7447c6ae99SBarry Smith ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 7547c6ae99SBarry Smith for (i=0; lnk; lnk=lnk->next,i++) { 769ae5db72SJed Brown ierr = PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr); 7747c6ae99SBarry Smith ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 7847c6ae99SBarry Smith ierr = DMView(lnk->dm,v);CHKERRQ(ierr); 7947c6ae99SBarry Smith ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 8047c6ae99SBarry Smith } 8147c6ae99SBarry Smith ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 8247c6ae99SBarry Smith } 8347c6ae99SBarry Smith PetscFunctionReturn(0); 8447c6ae99SBarry Smith } 8547c6ae99SBarry Smith 8647c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/ 877087cfbeSBarry Smith PetscErrorCode DMSetUp_Composite(DM dm) 8847c6ae99SBarry Smith { 8947c6ae99SBarry Smith PetscErrorCode ierr; 9047c6ae99SBarry Smith PetscInt nprev = 0; 9147c6ae99SBarry Smith PetscMPIInt rank,size; 9247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 9347c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 9447c6ae99SBarry Smith PetscLayout map; 9547c6ae99SBarry Smith 9647c6ae99SBarry Smith PetscFunctionBegin; 97ce94432eSBarry Smith if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup"); 98ce94432eSBarry Smith ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm),&map);CHKERRQ(ierr); 9947c6ae99SBarry Smith ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr); 10047c6ae99SBarry Smith ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr); 10147c6ae99SBarry Smith ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr); 10247c6ae99SBarry Smith ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 10347c6ae99SBarry Smith ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr); 1040298fd71SBarry Smith ierr = PetscLayoutGetRange(map,&com->rstart,NULL);CHKERRQ(ierr); 105fcfd50ebSBarry Smith ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 10647c6ae99SBarry Smith 1079ae5db72SJed Brown /* now set the rstart for each linked vector */ 108ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 109ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 11047c6ae99SBarry Smith while (next) { 11147c6ae99SBarry Smith next->rstart = nprev; 11206ebdd98SJed Brown nprev += next->n; 11347c6ae99SBarry Smith next->grstart = com->rstart + next->rstart; 114785e854fSJed Brown ierr = PetscMalloc1(size,&next->grstarts);CHKERRQ(ierr); 115ce94432eSBarry Smith ierr = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 11647c6ae99SBarry Smith next = next->next; 11747c6ae99SBarry Smith } 11847c6ae99SBarry Smith com->setup = PETSC_TRUE; 11947c6ae99SBarry Smith PetscFunctionReturn(0); 12047c6ae99SBarry Smith } 12147c6ae99SBarry Smith 12247c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/ 12347c6ae99SBarry Smith 12473e31fe2SJed Brown /*@ 12547c6ae99SBarry Smith DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite 12647c6ae99SBarry Smith representation. 12747c6ae99SBarry Smith 12847c6ae99SBarry Smith Not Collective 12947c6ae99SBarry Smith 13047c6ae99SBarry Smith Input Parameter: 13147c6ae99SBarry Smith . dm - the packer object 13247c6ae99SBarry Smith 13347c6ae99SBarry Smith Output Parameter: 13447c6ae99SBarry Smith . nDM - the number of DMs 13547c6ae99SBarry Smith 13647c6ae99SBarry Smith Level: beginner 13747c6ae99SBarry Smith 13847c6ae99SBarry Smith @*/ 1397087cfbeSBarry Smith PetscErrorCode DMCompositeGetNumberDM(DM dm,PetscInt *nDM) 14047c6ae99SBarry Smith { 14147c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 142*71b14b3eSStefano Zampini PetscBool flg; 143*71b14b3eSStefano Zampini PetscErrorCode ierr; 1445fd66863SKarl Rupp 14547c6ae99SBarry Smith PetscFunctionBegin; 14647c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 147*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 148*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 14947c6ae99SBarry Smith *nDM = com->nDM; 15047c6ae99SBarry Smith PetscFunctionReturn(0); 15147c6ae99SBarry Smith } 15247c6ae99SBarry Smith 15347c6ae99SBarry Smith 15447c6ae99SBarry Smith /*@C 15547c6ae99SBarry Smith DMCompositeGetAccess - Allows one to access the individual packed vectors in their global 15647c6ae99SBarry Smith representation. 15747c6ae99SBarry Smith 15847c6ae99SBarry Smith Collective on DMComposite 15947c6ae99SBarry Smith 1609ae5db72SJed Brown Input Parameters: 16147c6ae99SBarry Smith + dm - the packer object 1629ae5db72SJed Brown - gvec - the global vector 1639ae5db72SJed Brown 1649ae5db72SJed Brown Output Parameters: 1650298fd71SBarry Smith . Vec* ... - the packed parallel vectors, NULL for those that are not needed 16647c6ae99SBarry Smith 16747c6ae99SBarry Smith Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them 16847c6ae99SBarry Smith 169f73e5cebSJed Brown Fortran Notes: 170f73e5cebSJed Brown 171f73e5cebSJed Brown Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4) 172f73e5cebSJed Brown or use the alternative interface DMCompositeGetAccessArray(). 173f73e5cebSJed Brown 17447c6ae99SBarry Smith Level: advanced 17547c6ae99SBarry Smith 176f73e5cebSJed Brown .seealso: DMCompositeGetEntries(), DMCompositeScatter() 17747c6ae99SBarry Smith @*/ 1787087cfbeSBarry Smith PetscErrorCode DMCompositeGetAccess(DM dm,Vec gvec,...) 17947c6ae99SBarry Smith { 18047c6ae99SBarry Smith va_list Argp; 18147c6ae99SBarry Smith PetscErrorCode ierr; 18247c6ae99SBarry Smith struct DMCompositeLink *next; 18347c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 1845edff71fSBarry Smith PetscInt readonly; 185*71b14b3eSStefano Zampini PetscBool flg; 18647c6ae99SBarry Smith 18747c6ae99SBarry Smith PetscFunctionBegin; 18847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 18947c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 190*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 191*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 19247c6ae99SBarry Smith next = com->next; 19347c6ae99SBarry Smith if (!com->setup) { 194d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 19547c6ae99SBarry Smith } 19647c6ae99SBarry Smith 1975edff71fSBarry Smith ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr); 19847c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 19947c6ae99SBarry Smith va_start(Argp,gvec); 20047c6ae99SBarry Smith while (next) { 20147c6ae99SBarry Smith Vec *vec; 20247c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 2039ae5db72SJed Brown if (vec) { 2049ae5db72SJed Brown ierr = DMGetGlobalVector(next->dm,vec);CHKERRQ(ierr); 2055edff71fSBarry Smith if (readonly) { 2065edff71fSBarry Smith const PetscScalar *array; 2075edff71fSBarry Smith ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr); 2085edff71fSBarry Smith ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr); 2095edff71fSBarry Smith ierr = VecLockPush(*vec);CHKERRQ(ierr); 2105edff71fSBarry Smith ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr); 2115edff71fSBarry Smith } else { 2125edff71fSBarry Smith PetscScalar *array; 2139ae5db72SJed Brown ierr = VecGetArray(gvec,&array);CHKERRQ(ierr); 2149ae5db72SJed Brown ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr); 2159ae5db72SJed Brown ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr); 21647c6ae99SBarry Smith } 2175edff71fSBarry Smith } 21847c6ae99SBarry Smith next = next->next; 21947c6ae99SBarry Smith } 22047c6ae99SBarry Smith va_end(Argp); 22147c6ae99SBarry Smith PetscFunctionReturn(0); 22247c6ae99SBarry Smith } 22347c6ae99SBarry Smith 224f73e5cebSJed Brown /*@C 225f73e5cebSJed Brown DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global 226f73e5cebSJed Brown representation. 227f73e5cebSJed Brown 228f73e5cebSJed Brown Collective on DMComposite 229f73e5cebSJed Brown 230f73e5cebSJed Brown Input Parameters: 231f73e5cebSJed Brown + dm - the packer object 232f73e5cebSJed Brown . pvec - packed vector 233f73e5cebSJed Brown . nwanted - number of vectors wanted 2340298fd71SBarry Smith - wanted - sorted array of vectors wanted, or NULL to get all vectors 235f73e5cebSJed Brown 236f73e5cebSJed Brown Output Parameters: 237f73e5cebSJed Brown . vecs - array of requested global vectors (must be allocated) 238f73e5cebSJed Brown 239f73e5cebSJed Brown Notes: Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them 240f73e5cebSJed Brown 241f73e5cebSJed Brown Level: advanced 242f73e5cebSJed Brown 243f73e5cebSJed Brown .seealso: DMCompositeGetAccess(), DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather() 244f73e5cebSJed Brown @*/ 245f73e5cebSJed Brown PetscErrorCode DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 246f73e5cebSJed Brown { 247f73e5cebSJed Brown PetscErrorCode ierr; 248f73e5cebSJed Brown struct DMCompositeLink *link; 249f73e5cebSJed Brown PetscInt i,wnum; 250f73e5cebSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 251bee642f7SBarry Smith PetscInt readonly; 252*71b14b3eSStefano Zampini PetscBool flg; 253f73e5cebSJed Brown 254f73e5cebSJed Brown PetscFunctionBegin; 255f73e5cebSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 256f73e5cebSJed Brown PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 257*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 258*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 259f73e5cebSJed Brown if (!com->setup) { 260f73e5cebSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 261f73e5cebSJed Brown } 262f73e5cebSJed Brown 263bee642f7SBarry Smith ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr); 264f73e5cebSJed Brown for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 265f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 266f73e5cebSJed Brown Vec v; 267f73e5cebSJed Brown ierr = DMGetGlobalVector(link->dm,&v);CHKERRQ(ierr); 268bee642f7SBarry Smith if (readonly) { 269bee642f7SBarry Smith const PetscScalar *array; 270bee642f7SBarry Smith ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr); 271bee642f7SBarry Smith ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr); 272bee642f7SBarry Smith ierr = VecLockPush(v);CHKERRQ(ierr); 273bee642f7SBarry Smith ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr); 274bee642f7SBarry Smith } else { 275bee642f7SBarry Smith PetscScalar *array; 276f73e5cebSJed Brown ierr = VecGetArray(pvec,&array);CHKERRQ(ierr); 277f73e5cebSJed Brown ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr); 278f73e5cebSJed Brown ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr); 279bee642f7SBarry Smith } 280f73e5cebSJed Brown vecs[wnum++] = v; 281f73e5cebSJed Brown } 282f73e5cebSJed Brown } 283f73e5cebSJed Brown PetscFunctionReturn(0); 284f73e5cebSJed Brown } 285f73e5cebSJed Brown 2867ac2b803SAlex Fikl /*@C 2877ac2b803SAlex Fikl DMCompositeGetLocalAccessArray - Allows one to access the individual 2887ac2b803SAlex Fikl packed vectors in their local representation. 2897ac2b803SAlex Fikl 2907ac2b803SAlex Fikl Collective on DMComposite. 2917ac2b803SAlex Fikl 2927ac2b803SAlex Fikl Input Parameters: 2937ac2b803SAlex Fikl + dm - the packer object 2947ac2b803SAlex Fikl . pvec - packed vector 2957ac2b803SAlex Fikl . nwanted - number of vectors wanted 2967ac2b803SAlex Fikl - wanted - sorted array of vectors wanted, or NULL to get all vectors 2977ac2b803SAlex Fikl 2987ac2b803SAlex Fikl Output Parameters: 2997ac2b803SAlex Fikl . vecs - array of requested local vectors (must be allocated) 3007ac2b803SAlex Fikl 3017ac2b803SAlex Fikl Notes: Use DMCompositeRestoreLocalAccessArray() to return the vectors 3027ac2b803SAlex Fikl when you no longer need them. 3037ac2b803SAlex Fikl 3047ac2b803SAlex Fikl Level: advanced 3057ac2b803SAlex Fikl 3067ac2b803SAlex Fikl .seealso: DMCompositeRestoreLocalAccessArray(), DMCompositeGetAccess(), 3077ac2b803SAlex Fikl DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather() 3087ac2b803SAlex Fikl @*/ 3097ac2b803SAlex Fikl PetscErrorCode DMCompositeGetLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 3107ac2b803SAlex Fikl { 3117ac2b803SAlex Fikl PetscErrorCode ierr; 3127ac2b803SAlex Fikl struct DMCompositeLink *link; 3137ac2b803SAlex Fikl PetscInt i,wnum; 3147ac2b803SAlex Fikl DM_Composite *com = (DM_Composite*)dm->data; 3157ac2b803SAlex Fikl PetscInt readonly; 3167ac2b803SAlex Fikl PetscInt nlocal = 0; 317*71b14b3eSStefano Zampini PetscBool flg; 3187ac2b803SAlex Fikl 3197ac2b803SAlex Fikl PetscFunctionBegin; 3207ac2b803SAlex Fikl PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3217ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 322*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 323*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 3247ac2b803SAlex Fikl if (!com->setup) { 3257ac2b803SAlex Fikl ierr = DMSetUp(dm);CHKERRQ(ierr); 3267ac2b803SAlex Fikl } 3277ac2b803SAlex Fikl 3287ac2b803SAlex Fikl ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr); 3297ac2b803SAlex Fikl for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 3307ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 3317ac2b803SAlex Fikl Vec v; 3327ac2b803SAlex Fikl ierr = DMGetLocalVector(link->dm,&v);CHKERRQ(ierr); 3337ac2b803SAlex Fikl if (readonly) { 3347ac2b803SAlex Fikl const PetscScalar *array; 3357ac2b803SAlex Fikl ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr); 3367ac2b803SAlex Fikl ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr); 3377ac2b803SAlex Fikl ierr = VecLockPush(v);CHKERRQ(ierr); 3387ac2b803SAlex Fikl ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr); 3397ac2b803SAlex Fikl } else { 3407ac2b803SAlex Fikl PetscScalar *array; 3417ac2b803SAlex Fikl ierr = VecGetArray(pvec,&array);CHKERRQ(ierr); 3427ac2b803SAlex Fikl ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr); 3437ac2b803SAlex Fikl ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr); 3447ac2b803SAlex Fikl } 3457ac2b803SAlex Fikl vecs[wnum++] = v; 3467ac2b803SAlex Fikl } 3477ac2b803SAlex Fikl 3487ac2b803SAlex Fikl nlocal += link->nlocal; 3497ac2b803SAlex Fikl } 3507ac2b803SAlex Fikl 3517ac2b803SAlex Fikl PetscFunctionReturn(0); 3527ac2b803SAlex Fikl } 3537ac2b803SAlex Fikl 35447c6ae99SBarry Smith /*@C 355aa219208SBarry Smith DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess() 35647c6ae99SBarry Smith representation. 35747c6ae99SBarry Smith 35847c6ae99SBarry Smith Collective on DMComposite 35947c6ae99SBarry Smith 3609ae5db72SJed Brown Input Parameters: 36147c6ae99SBarry Smith + dm - the packer object 36247c6ae99SBarry Smith . gvec - the global vector 3630298fd71SBarry Smith - Vec* ... - the individual parallel vectors, NULL for those that are not needed 36447c6ae99SBarry Smith 36547c6ae99SBarry Smith Level: advanced 36647c6ae99SBarry Smith 3679ae5db72SJed Brown .seealso DMCompositeAddDM(), DMCreateGlobalVector(), 3686eb61c8cSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(), 369aa219208SBarry Smith DMCompositeRestoreAccess(), DMCompositeGetAccess() 37047c6ae99SBarry Smith 37147c6ae99SBarry Smith @*/ 3727087cfbeSBarry Smith PetscErrorCode DMCompositeRestoreAccess(DM dm,Vec gvec,...) 37347c6ae99SBarry Smith { 37447c6ae99SBarry Smith va_list Argp; 37547c6ae99SBarry Smith PetscErrorCode ierr; 37647c6ae99SBarry Smith struct DMCompositeLink *next; 37747c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 3785edff71fSBarry Smith PetscInt readonly; 379*71b14b3eSStefano Zampini PetscBool flg; 38047c6ae99SBarry Smith 38147c6ae99SBarry Smith PetscFunctionBegin; 38247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 38347c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 384*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 385*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 38647c6ae99SBarry Smith next = com->next; 38747c6ae99SBarry Smith if (!com->setup) { 388d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 38947c6ae99SBarry Smith } 39047c6ae99SBarry Smith 3915edff71fSBarry Smith ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr); 39247c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 39347c6ae99SBarry Smith va_start(Argp,gvec); 39447c6ae99SBarry Smith while (next) { 39547c6ae99SBarry Smith Vec *vec; 39647c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 3979ae5db72SJed Brown if (vec) { 3989ae5db72SJed Brown ierr = VecResetArray(*vec);CHKERRQ(ierr); 3995edff71fSBarry Smith if (readonly) { 4005edff71fSBarry Smith ierr = VecLockPop(*vec);CHKERRQ(ierr); 4015edff71fSBarry Smith } 402bee642f7SBarry Smith ierr = DMRestoreGlobalVector(next->dm,vec);CHKERRQ(ierr); 40347c6ae99SBarry Smith } 40447c6ae99SBarry Smith next = next->next; 40547c6ae99SBarry Smith } 40647c6ae99SBarry Smith va_end(Argp); 40747c6ae99SBarry Smith PetscFunctionReturn(0); 40847c6ae99SBarry Smith } 40947c6ae99SBarry Smith 410f73e5cebSJed Brown /*@C 411f73e5cebSJed Brown DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray() 412f73e5cebSJed Brown 413f73e5cebSJed Brown Collective on DMComposite 414f73e5cebSJed Brown 415f73e5cebSJed Brown Input Parameters: 416f73e5cebSJed Brown + dm - the packer object 417f73e5cebSJed Brown . pvec - packed vector 418f73e5cebSJed Brown . nwanted - number of vectors wanted 4190298fd71SBarry Smith . wanted - sorted array of vectors wanted, or NULL to get all vectors 420f73e5cebSJed Brown - vecs - array of global vectors to return 421f73e5cebSJed Brown 422f73e5cebSJed Brown Level: advanced 423f73e5cebSJed Brown 424f73e5cebSJed Brown .seealso: DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), DMCompositeScatter(), DMCompositeGather() 425f73e5cebSJed Brown @*/ 426f73e5cebSJed Brown PetscErrorCode DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 427f73e5cebSJed Brown { 428f73e5cebSJed Brown PetscErrorCode ierr; 429f73e5cebSJed Brown struct DMCompositeLink *link; 430f73e5cebSJed Brown PetscInt i,wnum; 431f73e5cebSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 432bee642f7SBarry Smith PetscInt readonly; 433*71b14b3eSStefano Zampini PetscBool flg; 434f73e5cebSJed Brown 435f73e5cebSJed Brown PetscFunctionBegin; 436f73e5cebSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 437f73e5cebSJed Brown PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 438*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 439*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 440f73e5cebSJed Brown if (!com->setup) { 441f73e5cebSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 442f73e5cebSJed Brown } 443f73e5cebSJed Brown 444bee642f7SBarry Smith ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr); 445f73e5cebSJed Brown for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 446f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 447f73e5cebSJed Brown ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr); 448bee642f7SBarry Smith if (readonly) { 449bee642f7SBarry Smith ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr); 450bee642f7SBarry Smith } 451f73e5cebSJed Brown ierr = DMRestoreGlobalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr); 452f73e5cebSJed Brown wnum++; 453f73e5cebSJed Brown } 454f73e5cebSJed Brown } 455f73e5cebSJed Brown PetscFunctionReturn(0); 456f73e5cebSJed Brown } 457f73e5cebSJed Brown 4587ac2b803SAlex Fikl /*@C 4597ac2b803SAlex Fikl DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray(). 4607ac2b803SAlex Fikl 4617ac2b803SAlex Fikl Collective on DMComposite. 4627ac2b803SAlex Fikl 4637ac2b803SAlex Fikl Input Parameters: 4647ac2b803SAlex Fikl + dm - the packer object 4657ac2b803SAlex Fikl . pvec - packed vector 4667ac2b803SAlex Fikl . nwanted - number of vectors wanted 4677ac2b803SAlex Fikl . wanted - sorted array of vectors wanted, or NULL to restore all vectors 4687ac2b803SAlex Fikl - vecs - array of local vectors to return 4697ac2b803SAlex Fikl 4707ac2b803SAlex Fikl Level: advanced 4717ac2b803SAlex Fikl 4727ac2b803SAlex Fikl Notes: 4737ac2b803SAlex Fikl nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray() 4747ac2b803SAlex Fikl otherwise the call will fail. 4757ac2b803SAlex Fikl 4767ac2b803SAlex Fikl .seealso: DMCompositeGetLocalAccessArray(), DMCompositeRestoreAccessArray(), 4777ac2b803SAlex Fikl DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), 4787ac2b803SAlex Fikl DMCompositeScatter(), DMCompositeGather() 4797ac2b803SAlex Fikl @*/ 4807ac2b803SAlex Fikl PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 4817ac2b803SAlex Fikl { 4827ac2b803SAlex Fikl PetscErrorCode ierr; 4837ac2b803SAlex Fikl struct DMCompositeLink *link; 4847ac2b803SAlex Fikl PetscInt i,wnum; 4857ac2b803SAlex Fikl DM_Composite *com = (DM_Composite*)dm->data; 4867ac2b803SAlex Fikl PetscInt readonly; 487*71b14b3eSStefano Zampini PetscBool flg; 4887ac2b803SAlex Fikl 4897ac2b803SAlex Fikl PetscFunctionBegin; 4907ac2b803SAlex Fikl PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4917ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 492*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 493*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 4947ac2b803SAlex Fikl if (!com->setup) { 4957ac2b803SAlex Fikl ierr = DMSetUp(dm);CHKERRQ(ierr); 4967ac2b803SAlex Fikl } 4977ac2b803SAlex Fikl 4987ac2b803SAlex Fikl ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr); 4997ac2b803SAlex Fikl for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 5007ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 5017ac2b803SAlex Fikl ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr); 5027ac2b803SAlex Fikl if (readonly) { 5037ac2b803SAlex Fikl ierr = VecLockPop(vecs[wnum]);CHKERRQ(ierr); 5047ac2b803SAlex Fikl } 5057ac2b803SAlex Fikl ierr = DMRestoreLocalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr); 5067ac2b803SAlex Fikl wnum++; 5077ac2b803SAlex Fikl } 5087ac2b803SAlex Fikl } 5097ac2b803SAlex Fikl PetscFunctionReturn(0); 5107ac2b803SAlex Fikl } 5117ac2b803SAlex Fikl 51247c6ae99SBarry Smith /*@C 51347c6ae99SBarry Smith DMCompositeScatter - Scatters from a global packed vector into its individual local vectors 51447c6ae99SBarry Smith 51547c6ae99SBarry Smith Collective on DMComposite 51647c6ae99SBarry Smith 5179ae5db72SJed Brown Input Parameters: 51847c6ae99SBarry Smith + dm - the packer object 51947c6ae99SBarry Smith . gvec - the global vector 5200298fd71SBarry Smith - Vec ... - the individual sequential vectors, NULL for those that are not needed 52147c6ae99SBarry Smith 52247c6ae99SBarry Smith Level: advanced 52347c6ae99SBarry Smith 5246f3c3dcfSJed Brown Notes: 5256f3c3dcfSJed Brown DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is 5266f3c3dcfSJed Brown accessible from Fortran. 5276f3c3dcfSJed Brown 5289ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 5296eb61c8cSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 53047c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 5316f3c3dcfSJed Brown DMCompositeScatterArray() 53247c6ae99SBarry Smith 53347c6ae99SBarry Smith @*/ 5347087cfbeSBarry Smith PetscErrorCode DMCompositeScatter(DM dm,Vec gvec,...) 53547c6ae99SBarry Smith { 53647c6ae99SBarry Smith va_list Argp; 53747c6ae99SBarry Smith PetscErrorCode ierr; 53847c6ae99SBarry Smith struct DMCompositeLink *next; 5398fd8f222SJed Brown PetscInt cnt; 54047c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 541*71b14b3eSStefano Zampini PetscBool flg; 54247c6ae99SBarry Smith 54347c6ae99SBarry Smith PetscFunctionBegin; 54447c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 54547c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 546*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 547*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 54847c6ae99SBarry Smith if (!com->setup) { 549d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 55047c6ae99SBarry Smith } 55147c6ae99SBarry Smith 55247c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 55347c6ae99SBarry Smith va_start(Argp,gvec); 5548fd8f222SJed Brown for (cnt=3,next=com->next; next; cnt++,next=next->next) { 5559ae5db72SJed Brown Vec local; 5569ae5db72SJed Brown local = va_arg(Argp, Vec); 5579ae5db72SJed Brown if (local) { 5589ae5db72SJed Brown Vec global; 5595edff71fSBarry Smith const PetscScalar *array; 5609ae5db72SJed Brown PetscValidHeaderSpecific(local,VEC_CLASSID,cnt); 5619ae5db72SJed Brown ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 5625edff71fSBarry Smith ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr); 5639ae5db72SJed Brown ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr); 5649ae5db72SJed Brown ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr); 5659ae5db72SJed Brown ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr); 5665edff71fSBarry Smith ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr); 5679ae5db72SJed Brown ierr = VecResetArray(global);CHKERRQ(ierr); 5689ae5db72SJed Brown ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 56947c6ae99SBarry Smith } 57047c6ae99SBarry Smith } 57147c6ae99SBarry Smith va_end(Argp); 57247c6ae99SBarry Smith PetscFunctionReturn(0); 57347c6ae99SBarry Smith } 57447c6ae99SBarry Smith 5756f3c3dcfSJed Brown /*@ 5766f3c3dcfSJed Brown DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors 5776f3c3dcfSJed Brown 5786f3c3dcfSJed Brown Collective on DMComposite 5796f3c3dcfSJed Brown 5806f3c3dcfSJed Brown Input Parameters: 5816f3c3dcfSJed Brown + dm - the packer object 5826f3c3dcfSJed Brown . gvec - the global vector 5836f3c3dcfSJed Brown . lvecs - array of local vectors, NULL for any that are not needed 5846f3c3dcfSJed Brown 5856f3c3dcfSJed Brown Level: advanced 5866f3c3dcfSJed Brown 5876f3c3dcfSJed Brown Note: 588907376e6SBarry Smith This is a non-variadic alternative to DMCompositeScatter() 5896f3c3dcfSJed Brown 5906f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector() 5916f3c3dcfSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 5926f3c3dcfSJed Brown DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 5936f3c3dcfSJed Brown 5946f3c3dcfSJed Brown @*/ 5956f3c3dcfSJed Brown PetscErrorCode DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs) 5966f3c3dcfSJed Brown { 5976f3c3dcfSJed Brown PetscErrorCode ierr; 5986f3c3dcfSJed Brown struct DMCompositeLink *next; 5996f3c3dcfSJed Brown PetscInt i; 6006f3c3dcfSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 601*71b14b3eSStefano Zampini PetscBool flg; 6026f3c3dcfSJed Brown 6036f3c3dcfSJed Brown PetscFunctionBegin; 6046f3c3dcfSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6056f3c3dcfSJed Brown PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 606*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 607*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 6086f3c3dcfSJed Brown if (!com->setup) { 6096f3c3dcfSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 6106f3c3dcfSJed Brown } 6116f3c3dcfSJed Brown 6126f3c3dcfSJed Brown /* loop over packed objects, handling one at at time */ 6136f3c3dcfSJed Brown for (i=0,next=com->next; next; next=next->next,i++) { 6146f3c3dcfSJed Brown if (lvecs[i]) { 6156f3c3dcfSJed Brown Vec global; 616c5d31e75SLisandro Dalcin const PetscScalar *array; 6176f3c3dcfSJed Brown PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3); 6186f3c3dcfSJed Brown ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 619c5d31e75SLisandro Dalcin ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr); 620c5d31e75SLisandro Dalcin ierr = VecPlaceArray(global,(PetscScalar*)array+next->rstart);CHKERRQ(ierr); 6216f3c3dcfSJed Brown ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr); 6226f3c3dcfSJed Brown ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr); 623c5d31e75SLisandro Dalcin ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr); 6246f3c3dcfSJed Brown ierr = VecResetArray(global);CHKERRQ(ierr); 6256f3c3dcfSJed Brown ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 6266f3c3dcfSJed Brown } 6276f3c3dcfSJed Brown } 6286f3c3dcfSJed Brown PetscFunctionReturn(0); 6296f3c3dcfSJed Brown } 6306f3c3dcfSJed Brown 63147c6ae99SBarry Smith /*@C 63247c6ae99SBarry Smith DMCompositeGather - Gathers into a global packed vector from its individual local vectors 63347c6ae99SBarry Smith 63447c6ae99SBarry Smith Collective on DMComposite 63547c6ae99SBarry Smith 63647c6ae99SBarry Smith Input Parameter: 63747c6ae99SBarry Smith + dm - the packer object 63847c6ae99SBarry Smith . gvec - the global vector 639907376e6SBarry Smith . imode - INSERT_VALUES or ADD_VALUES 6400298fd71SBarry Smith - Vec ... - the individual sequential vectors, NULL for any that are not needed 64147c6ae99SBarry Smith 64247c6ae99SBarry Smith Level: advanced 64347c6ae99SBarry Smith 644f5f57ec0SBarry Smith Not available from Fortran, Fortran users can use DMCompositeGatherArray() 645f5f57ec0SBarry Smith 6469ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 6476eb61c8cSJed Brown DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 64847c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 64947c6ae99SBarry Smith 65047c6ae99SBarry Smith @*/ 6511dac896bSSatish Balay PetscErrorCode DMCompositeGather(DM dm,InsertMode imode,Vec gvec,...) 65247c6ae99SBarry Smith { 65347c6ae99SBarry Smith va_list Argp; 65447c6ae99SBarry Smith PetscErrorCode ierr; 65547c6ae99SBarry Smith struct DMCompositeLink *next; 65647c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 6578fd8f222SJed Brown PetscInt cnt; 658*71b14b3eSStefano Zampini PetscBool flg; 65947c6ae99SBarry Smith 66047c6ae99SBarry Smith PetscFunctionBegin; 66147c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 66247c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 663*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 664*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 66547c6ae99SBarry Smith if (!com->setup) { 666d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 66747c6ae99SBarry Smith } 66847c6ae99SBarry Smith 66947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 6701dac896bSSatish Balay va_start(Argp,gvec); 6718fd8f222SJed Brown for (cnt=3,next=com->next; next; cnt++,next=next->next) { 6729ae5db72SJed Brown Vec local; 6739ae5db72SJed Brown local = va_arg(Argp, Vec); 6749ae5db72SJed Brown if (local) { 67547c6ae99SBarry Smith PetscScalar *array; 6769ae5db72SJed Brown Vec global; 6779ae5db72SJed Brown PetscValidHeaderSpecific(local,VEC_CLASSID,cnt); 6789ae5db72SJed Brown ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 6799ae5db72SJed Brown ierr = VecGetArray(gvec,&array);CHKERRQ(ierr); 6809ae5db72SJed Brown ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr); 6819ae5db72SJed Brown ierr = DMLocalToGlobalBegin(next->dm,local,imode,global);CHKERRQ(ierr); 6829ae5db72SJed Brown ierr = DMLocalToGlobalEnd(next->dm,local,imode,global);CHKERRQ(ierr); 6839ae5db72SJed Brown ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr); 6849ae5db72SJed Brown ierr = VecResetArray(global);CHKERRQ(ierr); 6859ae5db72SJed Brown ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 68647c6ae99SBarry Smith } 68747c6ae99SBarry Smith } 68847c6ae99SBarry Smith va_end(Argp); 68947c6ae99SBarry Smith PetscFunctionReturn(0); 69047c6ae99SBarry Smith } 69147c6ae99SBarry Smith 6926f3c3dcfSJed Brown /*@ 6936f3c3dcfSJed Brown DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors 6946f3c3dcfSJed Brown 6956f3c3dcfSJed Brown Collective on DMComposite 6966f3c3dcfSJed Brown 6976f3c3dcfSJed Brown Input Parameter: 6986f3c3dcfSJed Brown + dm - the packer object 6996f3c3dcfSJed Brown . gvec - the global vector 700907376e6SBarry Smith . imode - INSERT_VALUES or ADD_VALUES 7016f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed 7026f3c3dcfSJed Brown 7036f3c3dcfSJed Brown Level: advanced 7046f3c3dcfSJed Brown 7056f3c3dcfSJed Brown Notes: 7066f3c3dcfSJed Brown This is a non-variadic alternative to DMCompositeGather(). 7076f3c3dcfSJed Brown 7086f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 7096f3c3dcfSJed Brown DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 7106f3c3dcfSJed Brown DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries(), 7116f3c3dcfSJed Brown @*/ 7121dac896bSSatish Balay PetscErrorCode DMCompositeGatherArray(DM dm,InsertMode imode,Vec gvec,Vec *lvecs) 7136f3c3dcfSJed Brown { 7146f3c3dcfSJed Brown PetscErrorCode ierr; 7156f3c3dcfSJed Brown struct DMCompositeLink *next; 7166f3c3dcfSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 7176f3c3dcfSJed Brown PetscInt i; 718*71b14b3eSStefano Zampini PetscBool flg; 7196f3c3dcfSJed Brown 7206f3c3dcfSJed Brown PetscFunctionBegin; 7216f3c3dcfSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 7226f3c3dcfSJed Brown PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 723*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 724*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 7256f3c3dcfSJed Brown if (!com->setup) { 7266f3c3dcfSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 7276f3c3dcfSJed Brown } 7286f3c3dcfSJed Brown 7296f3c3dcfSJed Brown /* loop over packed objects, handling one at at time */ 7306f3c3dcfSJed Brown for (next=com->next,i=0; next; next=next->next,i++) { 7316f3c3dcfSJed Brown if (lvecs[i]) { 7326f3c3dcfSJed Brown PetscScalar *array; 7336f3c3dcfSJed Brown Vec global; 7346f3c3dcfSJed Brown PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3); 7356f3c3dcfSJed Brown ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 7366f3c3dcfSJed Brown ierr = VecGetArray(gvec,&array);CHKERRQ(ierr); 7376f3c3dcfSJed Brown ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr); 7386f3c3dcfSJed Brown ierr = DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global);CHKERRQ(ierr); 7396f3c3dcfSJed Brown ierr = DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global);CHKERRQ(ierr); 7406f3c3dcfSJed Brown ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr); 7416f3c3dcfSJed Brown ierr = VecResetArray(global);CHKERRQ(ierr); 7426f3c3dcfSJed Brown ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 7436f3c3dcfSJed Brown } 7446f3c3dcfSJed Brown } 7456f3c3dcfSJed Brown PetscFunctionReturn(0); 7466f3c3dcfSJed Brown } 7476f3c3dcfSJed Brown 748f5f57ec0SBarry Smith /*@ 749aa219208SBarry Smith DMCompositeAddDM - adds a DM vector to a DMComposite 75047c6ae99SBarry Smith 75147c6ae99SBarry Smith Collective on DMComposite 75247c6ae99SBarry Smith 75347c6ae99SBarry Smith Input Parameter: 7549b52a9b5SPatrick Sanan + dmc - the DMComposite (packer) object 755f5f57ec0SBarry Smith - dm - the DM object 75647c6ae99SBarry Smith 75747c6ae99SBarry Smith Level: advanced 75847c6ae99SBarry Smith 7590c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(), 7606eb61c8cSJed Brown DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 76147c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 76247c6ae99SBarry Smith 76347c6ae99SBarry Smith @*/ 7647087cfbeSBarry Smith PetscErrorCode DMCompositeAddDM(DM dmc,DM dm) 76547c6ae99SBarry Smith { 76647c6ae99SBarry Smith PetscErrorCode ierr; 76706ebdd98SJed Brown PetscInt n,nlocal; 76847c6ae99SBarry Smith struct DMCompositeLink *mine,*next; 76906ebdd98SJed Brown Vec global,local; 77047c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dmc->data; 771*71b14b3eSStefano Zampini PetscBool flg; 77247c6ae99SBarry Smith 77347c6ae99SBarry Smith PetscFunctionBegin; 77447c6ae99SBarry Smith PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 77547c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,2); 776*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dmc,DMCOMPOSITE,&flg);CHKERRQ(ierr); 777*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 77847c6ae99SBarry Smith next = com->next; 779ce94432eSBarry Smith if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite"); 78047c6ae99SBarry Smith 78147c6ae99SBarry Smith /* create new link */ 782b00a9115SJed Brown ierr = PetscNew(&mine);CHKERRQ(ierr); 78347c6ae99SBarry Smith ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 78447c6ae99SBarry Smith ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr); 78547c6ae99SBarry Smith ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr); 78647c6ae99SBarry Smith ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr); 78706ebdd98SJed Brown ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr); 78806ebdd98SJed Brown ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr); 78906ebdd98SJed Brown ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr); 7908865f1eaSKarl Rupp 79147c6ae99SBarry Smith mine->n = n; 79206ebdd98SJed Brown mine->nlocal = nlocal; 79347c6ae99SBarry Smith mine->dm = dm; 7940298fd71SBarry Smith mine->next = NULL; 79547c6ae99SBarry Smith com->n += n; 7967ac2b803SAlex Fikl com->nghost += nlocal; 79747c6ae99SBarry Smith 79847c6ae99SBarry Smith /* add to end of list */ 7998865f1eaSKarl Rupp if (!next) com->next = mine; 8008865f1eaSKarl Rupp else { 80147c6ae99SBarry Smith while (next->next) next = next->next; 80247c6ae99SBarry Smith next->next = mine; 80347c6ae99SBarry Smith } 80447c6ae99SBarry Smith com->nDM++; 80547c6ae99SBarry Smith com->nmine++; 80647c6ae99SBarry Smith PetscFunctionReturn(0); 80747c6ae99SBarry Smith } 80847c6ae99SBarry Smith 8099804daf3SBarry Smith #include <petscdraw.h> 81026887b52SJed Brown PETSC_EXTERN PetscErrorCode VecView_MPI(Vec,PetscViewer); 8117087cfbeSBarry Smith PetscErrorCode VecView_DMComposite(Vec gvec,PetscViewer viewer) 81247c6ae99SBarry Smith { 81347c6ae99SBarry Smith DM dm; 81447c6ae99SBarry Smith PetscErrorCode ierr; 81547c6ae99SBarry Smith struct DMCompositeLink *next; 81647c6ae99SBarry Smith PetscBool isdraw; 817cef07954SSatish Balay DM_Composite *com; 81847c6ae99SBarry Smith 81947c6ae99SBarry Smith PetscFunctionBegin; 820c688c046SMatthew G Knepley ierr = VecGetDM(gvec, &dm);CHKERRQ(ierr); 821ce94432eSBarry Smith if (!dm) SETERRQ(PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite"); 82247c6ae99SBarry Smith com = (DM_Composite*)dm->data; 82347c6ae99SBarry Smith next = com->next; 82447c6ae99SBarry Smith 825251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 82647c6ae99SBarry Smith if (!isdraw) { 82747c6ae99SBarry Smith /* do I really want to call this? */ 82847c6ae99SBarry Smith ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr); 82947c6ae99SBarry Smith } else { 83047c6ae99SBarry Smith PetscInt cnt = 0; 83147c6ae99SBarry Smith 83247c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 83347c6ae99SBarry Smith while (next) { 83447c6ae99SBarry Smith Vec vec; 835ca4278abSLisandro Dalcin const PetscScalar *array; 83647c6ae99SBarry Smith PetscInt bs; 83747c6ae99SBarry Smith 8389ae5db72SJed Brown /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */ 8399ae5db72SJed Brown ierr = DMGetGlobalVector(next->dm,&vec);CHKERRQ(ierr); 840ca4278abSLisandro Dalcin ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr); 841ca4278abSLisandro Dalcin ierr = VecPlaceArray(vec,(PetscScalar*)array+next->rstart);CHKERRQ(ierr); 842ca4278abSLisandro Dalcin ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr); 84347c6ae99SBarry Smith ierr = VecView(vec,viewer);CHKERRQ(ierr); 8449ae5db72SJed Brown ierr = VecResetArray(vec);CHKERRQ(ierr); 845ca4278abSLisandro Dalcin ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr); 8469ae5db72SJed Brown ierr = DMRestoreGlobalVector(next->dm,&vec);CHKERRQ(ierr); 84747c6ae99SBarry Smith ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr); 84847c6ae99SBarry Smith cnt += bs; 84947c6ae99SBarry Smith next = next->next; 85047c6ae99SBarry Smith } 85147c6ae99SBarry Smith ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr); 85247c6ae99SBarry Smith } 85347c6ae99SBarry Smith PetscFunctionReturn(0); 85447c6ae99SBarry Smith } 85547c6ae99SBarry Smith 8567087cfbeSBarry Smith PetscErrorCode DMCreateGlobalVector_Composite(DM dm,Vec *gvec) 85747c6ae99SBarry Smith { 85847c6ae99SBarry Smith PetscErrorCode ierr; 85947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 86047c6ae99SBarry Smith 86147c6ae99SBarry Smith PetscFunctionBegin; 86247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 863d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 864ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),com->n,com->N,gvec);CHKERRQ(ierr); 865c688c046SMatthew G Knepley ierr = VecSetDM(*gvec, dm);CHKERRQ(ierr); 86647c6ae99SBarry Smith ierr = VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);CHKERRQ(ierr); 86747c6ae99SBarry Smith PetscFunctionReturn(0); 86847c6ae99SBarry Smith } 86947c6ae99SBarry Smith 8707087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_Composite(DM dm,Vec *lvec) 87147c6ae99SBarry Smith { 87247c6ae99SBarry Smith PetscErrorCode ierr; 87347c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 87447c6ae99SBarry Smith 87547c6ae99SBarry Smith PetscFunctionBegin; 87647c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 87747c6ae99SBarry Smith if (!com->setup) { 878d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 87947c6ae99SBarry Smith } 880f0e01b1fSVincent Le Chenadec ierr = VecCreateSeq(PETSC_COMM_SELF,com->nghost,lvec);CHKERRQ(ierr); 881c688c046SMatthew G Knepley ierr = VecSetDM(*lvec, dm);CHKERRQ(ierr); 88247c6ae99SBarry Smith PetscFunctionReturn(0); 88347c6ae99SBarry Smith } 88447c6ae99SBarry Smith 88547c6ae99SBarry Smith /*@C 8869ae5db72SJed Brown DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space 88747c6ae99SBarry Smith 88806ebdd98SJed Brown Collective on DM 88947c6ae99SBarry Smith 89047c6ae99SBarry Smith Input Parameter: 89147c6ae99SBarry Smith . dm - the packer object 89247c6ae99SBarry Smith 89347c6ae99SBarry Smith Output Parameters: 8949ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes 8959ae5db72SJed Brown all the ghost points that individual ghosted DMDA's may have. 89647c6ae99SBarry Smith 89747c6ae99SBarry Smith Level: advanced 89847c6ae99SBarry Smith 89947c6ae99SBarry Smith Notes: 9006eb61c8cSJed Brown Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree(). 90147c6ae99SBarry Smith 902f5f57ec0SBarry Smith Not available from Fortran 903f5f57ec0SBarry Smith 9049ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 90547c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(), 90647c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries() 90747c6ae99SBarry Smith 90847c6ae99SBarry Smith @*/ 9097087cfbeSBarry Smith PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs) 91047c6ae99SBarry Smith { 91147c6ae99SBarry Smith PetscErrorCode ierr; 91247c6ae99SBarry Smith PetscInt i,*idx,n,cnt; 91347c6ae99SBarry Smith struct DMCompositeLink *next; 91447c6ae99SBarry Smith PetscMPIInt rank; 91547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 916*71b14b3eSStefano Zampini PetscBool flg; 91747c6ae99SBarry Smith 91847c6ae99SBarry Smith PetscFunctionBegin; 91947c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 920*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 921*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 922728e99d6SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 923854ce69bSBarry Smith ierr = PetscMalloc1(com->nDM,ltogs);CHKERRQ(ierr); 92447c6ae99SBarry Smith next = com->next; 925ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 92647c6ae99SBarry Smith 92747c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 92847c6ae99SBarry Smith cnt = 0; 92947c6ae99SBarry Smith while (next) { 9306eb61c8cSJed Brown ISLocalToGlobalMapping ltog; 9316eb61c8cSJed Brown PetscMPIInt size; 93286994e45SJed Brown const PetscInt *suboff,*indices; 9336eb61c8cSJed Brown Vec global; 93447c6ae99SBarry Smith 9356eb61c8cSJed Brown /* Get sub-DM global indices for each local dof */ 9361411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(next->dm,<og);CHKERRQ(ierr); 9376eb61c8cSJed Brown ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr); 93886994e45SJed Brown ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr); 939785e854fSJed Brown ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr); 94047c6ae99SBarry Smith 9416eb61c8cSJed Brown /* Get the offsets for the sub-DM global vector */ 9426eb61c8cSJed Brown ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 9436eb61c8cSJed Brown ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr); 944ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)global),&size);CHKERRQ(ierr); 9456eb61c8cSJed Brown 9466eb61c8cSJed Brown /* Shift the sub-DM definition of the global space to the composite global space */ 9476eb61c8cSJed Brown for (i=0; i<n; i++) { 94886994e45SJed Brown PetscInt subi = indices[i],lo = 0,hi = size,t; 9496eb61c8cSJed Brown /* Binary search to find which rank owns subi */ 9506eb61c8cSJed Brown while (hi-lo > 1) { 9516eb61c8cSJed Brown t = lo + (hi-lo)/2; 9526eb61c8cSJed Brown if (suboff[t] > subi) hi = t; 9536eb61c8cSJed Brown else lo = t; 9546eb61c8cSJed Brown } 9556eb61c8cSJed Brown idx[i] = subi - suboff[lo] + next->grstarts[lo]; 9566eb61c8cSJed Brown } 95786994e45SJed Brown ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr); 958f0413b6fSBarry Smith ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr); 9596eb61c8cSJed Brown ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 96047c6ae99SBarry Smith next = next->next; 96147c6ae99SBarry Smith cnt++; 96247c6ae99SBarry Smith } 96347c6ae99SBarry Smith PetscFunctionReturn(0); 96447c6ae99SBarry Smith } 96547c6ae99SBarry Smith 96687c85e80SJed Brown /*@C 9679ae5db72SJed Brown DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector 96887c85e80SJed Brown 96987c85e80SJed Brown Not Collective 97087c85e80SJed Brown 97187c85e80SJed Brown Input Arguments: 97287c85e80SJed Brown . dm - composite DM 97387c85e80SJed Brown 97487c85e80SJed Brown Output Arguments: 97587c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite 97687c85e80SJed Brown 97787c85e80SJed Brown Level: intermediate 97887c85e80SJed Brown 97987c85e80SJed Brown Notes: 98087c85e80SJed Brown At present, a composite local vector does not normally exist. This function is used to provide index sets for 98187c85e80SJed Brown MatGetLocalSubMatrix(). In the future, the scatters for each entry in the DMComposite may be be merged into a single 9829ae5db72SJed Brown scatter to a composite local vector. The user should not typically need to know which is being done. 98387c85e80SJed Brown 98487c85e80SJed Brown To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings(). 98587c85e80SJed Brown 98687c85e80SJed Brown To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs(). 98787c85e80SJed Brown 98887c85e80SJed Brown Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree(). 98987c85e80SJed Brown 990f5f57ec0SBarry Smith Not available from Fortran 991f5f57ec0SBarry Smith 99287c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef() 99387c85e80SJed Brown @*/ 9947087cfbeSBarry Smith PetscErrorCode DMCompositeGetLocalISs(DM dm,IS **is) 99587c85e80SJed Brown { 99687c85e80SJed Brown PetscErrorCode ierr; 99787c85e80SJed Brown DM_Composite *com = (DM_Composite*)dm->data; 99887c85e80SJed Brown struct DMCompositeLink *link; 99987c85e80SJed Brown PetscInt cnt,start; 1000*71b14b3eSStefano Zampini PetscBool flg; 100187c85e80SJed Brown 100287c85e80SJed Brown PetscFunctionBegin; 100387c85e80SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 100487c85e80SJed Brown PetscValidPointer(is,2); 1005*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 1006*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 1007785e854fSJed Brown ierr = PetscMalloc1(com->nmine,is);CHKERRQ(ierr); 100806ebdd98SJed Brown for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) { 1009520db06cSJed Brown PetscInt bs; 10109ae5db72SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr); 10111411c6eeSJed Brown ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr); 1012520db06cSJed Brown ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr); 1013520db06cSJed Brown } 101487c85e80SJed Brown PetscFunctionReturn(0); 101587c85e80SJed Brown } 101687c85e80SJed Brown 101747c6ae99SBarry Smith /*@C 101847c6ae99SBarry Smith DMCompositeGetGlobalISs - Gets the index sets for each composed object 101947c6ae99SBarry Smith 102047c6ae99SBarry Smith Collective on DMComposite 102147c6ae99SBarry Smith 102247c6ae99SBarry Smith Input Parameter: 102347c6ae99SBarry Smith . dm - the packer object 102447c6ae99SBarry Smith 102547c6ae99SBarry Smith Output Parameters: 102647c6ae99SBarry Smith . is - the array of index sets 102747c6ae99SBarry Smith 102847c6ae99SBarry Smith Level: advanced 102947c6ae99SBarry Smith 103047c6ae99SBarry Smith Notes: 103147c6ae99SBarry Smith The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree() 103247c6ae99SBarry Smith 103347c6ae99SBarry Smith These could be used to extract a subset of vector entries for a "multi-physics" preconditioner 103447c6ae99SBarry Smith 10356eb61c8cSJed Brown Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and 10366eb61c8cSJed Brown DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global 10376eb61c8cSJed Brown indices. 103847c6ae99SBarry Smith 1039f3cb0f7eSJed Brown Fortran Notes: 1040f3cb0f7eSJed Brown 1041f3cb0f7eSJed Brown The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM(). 1042f3cb0f7eSJed Brown 10439ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 104447c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(), 104547c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries() 104647c6ae99SBarry Smith 104747c6ae99SBarry Smith @*/ 10487087cfbeSBarry Smith PetscErrorCode DMCompositeGetGlobalISs(DM dm,IS *is[]) 104947c6ae99SBarry Smith { 105047c6ae99SBarry Smith PetscErrorCode ierr; 105166bb578eSMark Adams PetscInt cnt = 0; 105247c6ae99SBarry Smith struct DMCompositeLink *next; 105347c6ae99SBarry Smith PetscMPIInt rank; 105447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 1055*71b14b3eSStefano Zampini PetscBool flg; 105647c6ae99SBarry Smith 105747c6ae99SBarry Smith PetscFunctionBegin; 105847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1059*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 1060*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 1061*71b14b3eSStefano Zampini if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() before"); 1062854ce69bSBarry Smith ierr = PetscMalloc1(com->nDM,is);CHKERRQ(ierr); 106347c6ae99SBarry Smith next = com->next; 1064ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 106547c6ae99SBarry Smith 106647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 106747c6ae99SBarry Smith while (next) { 106866bb578eSMark Adams ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt]);CHKERRQ(ierr); 10690f21e855SMatthew G. Knepley if (dm->prob) { 107065c226d8SMatthew G. Knepley MatNullSpace space; 107165c226d8SMatthew G. Knepley Mat pmat; 10720f21e855SMatthew G. Knepley PetscObject disc; 10730f21e855SMatthew G. Knepley PetscInt Nf; 107465c226d8SMatthew G. Knepley 10752764a2aaSMatthew G. Knepley ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr); 1076f24dd8d2SMatthew G. Knepley if (cnt < Nf) { 10772764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(dm->prob, cnt, &disc);CHKERRQ(ierr); 10780f21e855SMatthew G. Knepley ierr = PetscObjectQuery(disc, "nullspace", (PetscObject*) &space);CHKERRQ(ierr); 1079aac2dd2dSMatthew G. Knepley if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space);CHKERRQ(ierr);} 10800f21e855SMatthew G. Knepley ierr = PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space);CHKERRQ(ierr); 1081aac2dd2dSMatthew G. Knepley if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space);CHKERRQ(ierr);} 10820f21e855SMatthew G. Knepley ierr = PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr); 1083aac2dd2dSMatthew G. Knepley if (pmat) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat);CHKERRQ(ierr);} 108465c226d8SMatthew G. Knepley } 1085f24dd8d2SMatthew G. Knepley } 108647c6ae99SBarry Smith cnt++; 108747c6ae99SBarry Smith next = next->next; 108847c6ae99SBarry Smith } 108947c6ae99SBarry Smith PetscFunctionReturn(0); 109047c6ae99SBarry Smith } 109147c6ae99SBarry Smith 109221c9b008SJed Brown PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields) 10934d343eeaSMatthew G Knepley { 10944d343eeaSMatthew G Knepley PetscInt nDM; 10954d343eeaSMatthew G Knepley DM *dms; 10964d343eeaSMatthew G Knepley PetscInt i; 10974d343eeaSMatthew G Knepley PetscErrorCode ierr; 10984d343eeaSMatthew G Knepley 10994d343eeaSMatthew G Knepley PetscFunctionBegin; 11004d343eeaSMatthew G Knepley ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr); 11018865f1eaSKarl Rupp if (numFields) *numFields = nDM; 11024d343eeaSMatthew G Knepley ierr = DMCompositeGetGlobalISs(dm, fields);CHKERRQ(ierr); 11034d343eeaSMatthew G Knepley if (fieldNames) { 1104785e854fSJed Brown ierr = PetscMalloc1(nDM, &dms);CHKERRQ(ierr); 1105785e854fSJed Brown ierr = PetscMalloc1(nDM, fieldNames);CHKERRQ(ierr); 11064d343eeaSMatthew G Knepley ierr = DMCompositeGetEntriesArray(dm, dms);CHKERRQ(ierr); 11074d343eeaSMatthew G Knepley for (i=0; i<nDM; i++) { 11084d343eeaSMatthew G Knepley char buf[256]; 11094d343eeaSMatthew G Knepley const char *splitname; 11104d343eeaSMatthew G Knepley 11114d343eeaSMatthew G Knepley /* Split naming precedence: object name, prefix, number */ 11124d343eeaSMatthew G Knepley splitname = ((PetscObject) dm)->name; 11134d343eeaSMatthew G Knepley if (!splitname) { 11144d343eeaSMatthew G Knepley ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr); 11154d343eeaSMatthew G Knepley if (splitname) { 11164d343eeaSMatthew G Knepley size_t len; 11178caf3d72SBarry Smith ierr = PetscStrncpy(buf,splitname,sizeof(buf));CHKERRQ(ierr); 11188caf3d72SBarry Smith buf[sizeof(buf) - 1] = 0; 11194d343eeaSMatthew G Knepley ierr = PetscStrlen(buf,&len);CHKERRQ(ierr); 11204d343eeaSMatthew G Knepley if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */ 11214d343eeaSMatthew G Knepley splitname = buf; 11224d343eeaSMatthew G Knepley } 11234d343eeaSMatthew G Knepley } 11244d343eeaSMatthew G Knepley if (!splitname) { 11258caf3d72SBarry Smith ierr = PetscSNPrintf(buf,sizeof(buf),"%D",i);CHKERRQ(ierr); 11264d343eeaSMatthew G Knepley splitname = buf; 11274d343eeaSMatthew G Knepley } 112821c9b008SJed Brown ierr = PetscStrallocpy(splitname,&(*fieldNames)[i]);CHKERRQ(ierr); 11294d343eeaSMatthew G Knepley } 11304d343eeaSMatthew G Knepley ierr = PetscFree(dms);CHKERRQ(ierr); 11314d343eeaSMatthew G Knepley } 11324d343eeaSMatthew G Knepley PetscFunctionReturn(0); 11334d343eeaSMatthew G Knepley } 11344d343eeaSMatthew G Knepley 1135e7c4fc90SDmitry Karpeev /* 1136e7c4fc90SDmitry Karpeev This could take over from DMCreateFieldIS(), as it is more general, 11370298fd71SBarry Smith making DMCreateFieldIS() a special case -- calling with dmlist == NULL; 1138e7c4fc90SDmitry Karpeev At this point it's probably best to be less intrusive, however. 1139e7c4fc90SDmitry Karpeev */ 114016621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist) 1141e7c4fc90SDmitry Karpeev { 1142e7c4fc90SDmitry Karpeev PetscInt nDM; 1143e7c4fc90SDmitry Karpeev PetscInt i; 1144e7c4fc90SDmitry Karpeev PetscErrorCode ierr; 1145e7c4fc90SDmitry Karpeev 1146e7c4fc90SDmitry Karpeev PetscFunctionBegin; 1147e7c4fc90SDmitry Karpeev ierr = DMCreateFieldIS_Composite(dm, len, namelist, islist);CHKERRQ(ierr); 1148e7c4fc90SDmitry Karpeev if (dmlist) { 1149e7c4fc90SDmitry Karpeev ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr); 1150785e854fSJed Brown ierr = PetscMalloc1(nDM, dmlist);CHKERRQ(ierr); 1151e7c4fc90SDmitry Karpeev ierr = DMCompositeGetEntriesArray(dm, *dmlist);CHKERRQ(ierr); 1152e7c4fc90SDmitry Karpeev for (i=0; i<nDM; i++) { 1153e7c4fc90SDmitry Karpeev ierr = PetscObjectReference((PetscObject)((*dmlist)[i]));CHKERRQ(ierr); 1154e7c4fc90SDmitry Karpeev } 1155e7c4fc90SDmitry Karpeev } 1156e7c4fc90SDmitry Karpeev PetscFunctionReturn(0); 1157e7c4fc90SDmitry Karpeev } 1158e7c4fc90SDmitry Karpeev 1159e7c4fc90SDmitry Karpeev 1160e7c4fc90SDmitry Karpeev 116147c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 116247c6ae99SBarry Smith /*@C 11639ae5db72SJed Brown DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite. 116447c6ae99SBarry Smith Use DMCompositeRestoreLocalVectors() to return them. 116547c6ae99SBarry Smith 116647c6ae99SBarry Smith Not Collective 116747c6ae99SBarry Smith 116847c6ae99SBarry Smith Input Parameter: 116947c6ae99SBarry Smith . dm - the packer object 117047c6ae99SBarry Smith 117147c6ae99SBarry Smith Output Parameter: 11729ae5db72SJed Brown . Vec ... - the individual sequential Vecs 117347c6ae99SBarry Smith 117447c6ae99SBarry Smith Level: advanced 117547c6ae99SBarry Smith 1176f5f57ec0SBarry Smith Not available from Fortran 1177f5f57ec0SBarry Smith 11789ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 11796eb61c8cSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 118047c6ae99SBarry Smith DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries() 118147c6ae99SBarry Smith 118247c6ae99SBarry Smith @*/ 11837087cfbeSBarry Smith PetscErrorCode DMCompositeGetLocalVectors(DM dm,...) 118447c6ae99SBarry Smith { 118547c6ae99SBarry Smith va_list Argp; 118647c6ae99SBarry Smith PetscErrorCode ierr; 118747c6ae99SBarry Smith struct DMCompositeLink *next; 118847c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 1189*71b14b3eSStefano Zampini PetscBool flg; 119047c6ae99SBarry Smith 119147c6ae99SBarry Smith PetscFunctionBegin; 119247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1193*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 1194*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 119547c6ae99SBarry Smith next = com->next; 119647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 119747c6ae99SBarry Smith va_start(Argp,dm); 119847c6ae99SBarry Smith while (next) { 119947c6ae99SBarry Smith Vec *vec; 120047c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 120106930112SJed Brown if (vec) {ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);} 120247c6ae99SBarry Smith next = next->next; 120347c6ae99SBarry Smith } 120447c6ae99SBarry Smith va_end(Argp); 120547c6ae99SBarry Smith PetscFunctionReturn(0); 120647c6ae99SBarry Smith } 120747c6ae99SBarry Smith 120847c6ae99SBarry Smith /*@C 12099ae5db72SJed Brown DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite. 121047c6ae99SBarry Smith 121147c6ae99SBarry Smith Not Collective 121247c6ae99SBarry Smith 121347c6ae99SBarry Smith Input Parameter: 121447c6ae99SBarry Smith . dm - the packer object 121547c6ae99SBarry Smith 121647c6ae99SBarry Smith Output Parameter: 12179ae5db72SJed Brown . Vec ... - the individual sequential Vecs 121847c6ae99SBarry Smith 121947c6ae99SBarry Smith Level: advanced 122047c6ae99SBarry Smith 1221f5f57ec0SBarry Smith Not available from Fortran 1222f5f57ec0SBarry Smith 12239ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 12246eb61c8cSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 122547c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries() 122647c6ae99SBarry Smith 122747c6ae99SBarry Smith @*/ 12287087cfbeSBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors(DM dm,...) 122947c6ae99SBarry Smith { 123047c6ae99SBarry Smith va_list Argp; 123147c6ae99SBarry Smith PetscErrorCode ierr; 123247c6ae99SBarry Smith struct DMCompositeLink *next; 123347c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 1234*71b14b3eSStefano Zampini PetscBool flg; 123547c6ae99SBarry Smith 123647c6ae99SBarry Smith PetscFunctionBegin; 123747c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1238*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 1239*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 124047c6ae99SBarry Smith next = com->next; 124147c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 124247c6ae99SBarry Smith va_start(Argp,dm); 124347c6ae99SBarry Smith while (next) { 124447c6ae99SBarry Smith Vec *vec; 124547c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 124606930112SJed Brown if (vec) {ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);} 124747c6ae99SBarry Smith next = next->next; 124847c6ae99SBarry Smith } 124947c6ae99SBarry Smith va_end(Argp); 125047c6ae99SBarry Smith PetscFunctionReturn(0); 125147c6ae99SBarry Smith } 125247c6ae99SBarry Smith 125347c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 125447c6ae99SBarry Smith /*@C 12559ae5db72SJed Brown DMCompositeGetEntries - Gets the DM for each entry in a DMComposite. 125647c6ae99SBarry Smith 125747c6ae99SBarry Smith Not Collective 125847c6ae99SBarry Smith 125947c6ae99SBarry Smith Input Parameter: 126047c6ae99SBarry Smith . dm - the packer object 126147c6ae99SBarry Smith 126247c6ae99SBarry Smith Output Parameter: 12639ae5db72SJed Brown . DM ... - the individual entries (DMs) 126447c6ae99SBarry Smith 126547c6ae99SBarry Smith Level: advanced 126647c6ae99SBarry Smith 1267f5f57ec0SBarry Smith Fortran Notes: Available as DMCompositeGetEntries() for one output DM, DMCompositeGetEntries2() for 2, etc 1268f5f57ec0SBarry Smith 12692fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray() 12706eb61c8cSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 127147c6ae99SBarry Smith DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(), DMCompositeScatter(), 127247c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors() 127347c6ae99SBarry Smith 127447c6ae99SBarry Smith @*/ 12757087cfbeSBarry Smith PetscErrorCode DMCompositeGetEntries(DM dm,...) 127647c6ae99SBarry Smith { 127747c6ae99SBarry Smith va_list Argp; 127847c6ae99SBarry Smith struct DMCompositeLink *next; 127947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 1280*71b14b3eSStefano Zampini PetscBool flg; 1281*71b14b3eSStefano Zampini PetscErrorCode ierr; 128247c6ae99SBarry Smith 128347c6ae99SBarry Smith PetscFunctionBegin; 128447c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1285*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 1286*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 128747c6ae99SBarry Smith next = com->next; 128847c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 128947c6ae99SBarry Smith va_start(Argp,dm); 129047c6ae99SBarry Smith while (next) { 129147c6ae99SBarry Smith DM *dmn; 129247c6ae99SBarry Smith dmn = va_arg(Argp, DM*); 12939ae5db72SJed Brown if (dmn) *dmn = next->dm; 129447c6ae99SBarry Smith next = next->next; 129547c6ae99SBarry Smith } 129647c6ae99SBarry Smith va_end(Argp); 129747c6ae99SBarry Smith PetscFunctionReturn(0); 129847c6ae99SBarry Smith } 129947c6ae99SBarry Smith 1300dbab29e1SMark F. Adams /*@C 13012fa5ba8aSJed Brown DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite. 13022fa5ba8aSJed Brown 13032fa5ba8aSJed Brown Not Collective 13042fa5ba8aSJed Brown 13052fa5ba8aSJed Brown Input Parameter: 1306907376e6SBarry Smith . dm - the packer object 1307907376e6SBarry Smith 1308907376e6SBarry Smith Output Parameter: 1309907376e6SBarry Smith . dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs 13102fa5ba8aSJed Brown 13112fa5ba8aSJed Brown Level: advanced 13122fa5ba8aSJed Brown 13132fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries() 13142fa5ba8aSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 13152fa5ba8aSJed Brown DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(), DMCompositeScatter(), 13162fa5ba8aSJed Brown DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors() 13172fa5ba8aSJed Brown 13182fa5ba8aSJed Brown @*/ 13192fa5ba8aSJed Brown PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[]) 13202fa5ba8aSJed Brown { 13212fa5ba8aSJed Brown struct DMCompositeLink *next; 13222fa5ba8aSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 13232fa5ba8aSJed Brown PetscInt i; 1324*71b14b3eSStefano Zampini PetscBool flg; 1325*71b14b3eSStefano Zampini PetscErrorCode ierr; 13262fa5ba8aSJed Brown 13272fa5ba8aSJed Brown PetscFunctionBegin; 13282fa5ba8aSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1329*71b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 1330*71b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 13312fa5ba8aSJed Brown /* loop over packed objects, handling one at at time */ 13322fa5ba8aSJed Brown for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm; 13332fa5ba8aSJed Brown PetscFunctionReturn(0); 13342fa5ba8aSJed Brown } 13352fa5ba8aSJed Brown 1336e10fd815SStefano Zampini typedef struct { 1337e10fd815SStefano Zampini DM dm; 1338e10fd815SStefano Zampini PetscViewer *subv; 1339e10fd815SStefano Zampini Vec *vecs; 1340e10fd815SStefano Zampini } GLVisViewerCtx; 1341e10fd815SStefano Zampini 1342e10fd815SStefano Zampini static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 1343e10fd815SStefano Zampini { 1344e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 1345e10fd815SStefano Zampini PetscInt i,n; 1346e10fd815SStefano Zampini PetscErrorCode ierr; 1347e10fd815SStefano Zampini 1348e10fd815SStefano Zampini PetscFunctionBegin; 1349e10fd815SStefano Zampini ierr = DMCompositeGetNumberDM(ctx->dm,&n);CHKERRQ(ierr); 1350e10fd815SStefano Zampini for (i = 0; i < n; i++) { 1351e10fd815SStefano Zampini ierr = PetscViewerDestroy(&ctx->subv[i]);CHKERRQ(ierr); 1352e10fd815SStefano Zampini } 1353e10fd815SStefano Zampini ierr = PetscFree2(ctx->subv,ctx->vecs);CHKERRQ(ierr); 1354e10fd815SStefano Zampini ierr = DMDestroy(&ctx->dm);CHKERRQ(ierr); 1355e10fd815SStefano Zampini ierr = PetscFree(ctx);CHKERRQ(ierr); 1356e10fd815SStefano Zampini PetscFunctionReturn(0); 1357e10fd815SStefano Zampini } 1358e10fd815SStefano Zampini 1359e10fd815SStefano Zampini static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 1360e10fd815SStefano Zampini { 1361e10fd815SStefano Zampini Vec X = (Vec)oX; 1362e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 1363e10fd815SStefano Zampini PetscInt i,n,cumf; 1364e10fd815SStefano Zampini PetscErrorCode ierr; 1365e10fd815SStefano Zampini 1366e10fd815SStefano Zampini PetscFunctionBegin; 1367e10fd815SStefano Zampini ierr = DMCompositeGetNumberDM(ctx->dm,&n);CHKERRQ(ierr); 1368e10fd815SStefano Zampini ierr = DMCompositeGetAccessArray(ctx->dm,X,n,NULL,ctx->vecs);CHKERRQ(ierr); 1369e10fd815SStefano Zampini for (i = 0, cumf = 0; i < n; i++) { 1370e10fd815SStefano Zampini PetscErrorCode (*g2l)(PetscObject,PetscInt,PetscObject[],void*); 1371e10fd815SStefano Zampini void *fctx; 1372e10fd815SStefano Zampini PetscInt nfi; 1373e10fd815SStefano Zampini 1374e10fd815SStefano Zampini ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nfi,NULL,NULL,&g2l,NULL,&fctx);CHKERRQ(ierr); 1375e10fd815SStefano Zampini if (!nfi) continue; 1376e10fd815SStefano Zampini if (g2l) { 1377e10fd815SStefano Zampini ierr = (*g2l)((PetscObject)ctx->vecs[i],nfi,oXfield+cumf,fctx);CHKERRQ(ierr); 1378e10fd815SStefano Zampini } else { 1379e10fd815SStefano Zampini ierr = VecCopy(ctx->vecs[i],(Vec)(oXfield[cumf]));CHKERRQ(ierr); 1380e10fd815SStefano Zampini } 1381e10fd815SStefano Zampini cumf += nfi; 1382e10fd815SStefano Zampini } 1383e10fd815SStefano Zampini ierr = DMCompositeRestoreAccessArray(ctx->dm,X,n,NULL,ctx->vecs);CHKERRQ(ierr); 1384e10fd815SStefano Zampini PetscFunctionReturn(0); 1385e10fd815SStefano Zampini } 1386e10fd815SStefano Zampini 1387e10fd815SStefano Zampini static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer) 1388e10fd815SStefano Zampini { 1389e10fd815SStefano Zampini DM dm = (DM)odm, *dms; 1390e10fd815SStefano Zampini Vec *Ufds; 1391e10fd815SStefano Zampini GLVisViewerCtx *ctx; 1392e10fd815SStefano Zampini PetscInt i,n,tnf,*sdim; 1393e10fd815SStefano Zampini char **fecs; 1394e10fd815SStefano Zampini PetscErrorCode ierr; 1395e10fd815SStefano Zampini 1396e10fd815SStefano Zampini PetscFunctionBegin; 1397e10fd815SStefano Zampini ierr = PetscNew(&ctx);CHKERRQ(ierr); 1398e10fd815SStefano Zampini ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 1399e10fd815SStefano Zampini ctx->dm = dm; 1400e10fd815SStefano Zampini ierr = DMCompositeGetNumberDM(dm,&n);CHKERRQ(ierr); 1401e10fd815SStefano Zampini ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 1402e10fd815SStefano Zampini ierr = DMCompositeGetEntriesArray(dm,dms);CHKERRQ(ierr); 1403e10fd815SStefano Zampini ierr = PetscMalloc2(n,&ctx->subv,n,&ctx->vecs);CHKERRQ(ierr); 1404e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1405e10fd815SStefano Zampini PetscInt nf; 1406e10fd815SStefano Zampini 1407e10fd815SStefano Zampini ierr = PetscViewerCreate(PetscObjectComm(odm),&ctx->subv[i]);CHKERRQ(ierr); 1408e10fd815SStefano Zampini ierr = PetscViewerSetType(ctx->subv[i],PETSCVIEWERGLVIS);CHKERRQ(ierr); 1409e10fd815SStefano Zampini ierr = PetscViewerGLVisSetDM_Private(ctx->subv[i],(PetscObject)dms[i]);CHKERRQ(ierr); 1410e10fd815SStefano Zampini ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1411e10fd815SStefano Zampini tnf += nf; 1412e10fd815SStefano Zampini } 1413e10fd815SStefano Zampini ierr = PetscFree(dms);CHKERRQ(ierr); 1414e10fd815SStefano Zampini ierr = PetscMalloc3(tnf,&fecs,tnf,&sdim,tnf,&Ufds);CHKERRQ(ierr); 1415e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1416e10fd815SStefano Zampini PetscInt *sd,nf,f; 1417e10fd815SStefano Zampini const char **fec; 1418e10fd815SStefano Zampini Vec *Uf; 1419e10fd815SStefano Zampini 1420e10fd815SStefano Zampini ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,&fec,&sd,NULL,(PetscObject**)&Uf,NULL);CHKERRQ(ierr); 1421e10fd815SStefano Zampini for (f = 0; f < nf; f++) { 1422e10fd815SStefano Zampini ierr = PetscStrallocpy(fec[f],&fecs[tnf+f]);CHKERRQ(ierr); 1423e10fd815SStefano Zampini Ufds[tnf+f] = Uf[f]; 1424e10fd815SStefano Zampini sdim[tnf+f] = sd[f]; 1425e10fd815SStefano Zampini } 1426e10fd815SStefano Zampini tnf += nf; 1427e10fd815SStefano Zampini } 1428e10fd815SStefano Zampini ierr = PetscViewerGLVisSetFields(viewer,tnf,(const char**)fecs,sdim,DMCompositeSampleGLVisFields_Private,(PetscObject*)Ufds,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr); 1429e10fd815SStefano Zampini for (i = 0; i < tnf; i++) { 1430e10fd815SStefano Zampini ierr = PetscFree(fecs[i]);CHKERRQ(ierr); 1431e10fd815SStefano Zampini } 1432e10fd815SStefano Zampini ierr = PetscFree3(fecs,sdim,Ufds);CHKERRQ(ierr); 1433e10fd815SStefano Zampini PetscFunctionReturn(0); 1434e10fd815SStefano Zampini } 1435e10fd815SStefano Zampini 14367087cfbeSBarry Smith PetscErrorCode DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine) 143747c6ae99SBarry Smith { 143847c6ae99SBarry Smith PetscErrorCode ierr; 143947c6ae99SBarry Smith struct DMCompositeLink *next; 144047c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dmi->data; 144147c6ae99SBarry Smith DM dm; 144247c6ae99SBarry Smith 144347c6ae99SBarry Smith PetscFunctionBegin; 144447c6ae99SBarry Smith PetscValidHeaderSpecific(dmi,DM_CLASSID,1); 1445ce94432eSBarry Smith if (comm == MPI_COMM_NULL) { 1446ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr); 1447ce94432eSBarry Smith } 14482ce3a92bSJed Brown ierr = DMSetUp(dmi);CHKERRQ(ierr); 144947c6ae99SBarry Smith next = com->next; 145047c6ae99SBarry Smith ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr); 145147c6ae99SBarry Smith 145247c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 145347c6ae99SBarry Smith while (next) { 145447c6ae99SBarry Smith ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr); 145547c6ae99SBarry Smith ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr); 145647c6ae99SBarry Smith ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 145747c6ae99SBarry Smith next = next->next; 145847c6ae99SBarry Smith } 145947c6ae99SBarry Smith PetscFunctionReturn(0); 146047c6ae99SBarry Smith } 146147c6ae99SBarry Smith 146214354c39SJed Brown PetscErrorCode DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine) 146314354c39SJed Brown { 146414354c39SJed Brown PetscErrorCode ierr; 146514354c39SJed Brown struct DMCompositeLink *next; 146614354c39SJed Brown DM_Composite *com = (DM_Composite*)dmi->data; 146714354c39SJed Brown DM dm; 146814354c39SJed Brown 146914354c39SJed Brown PetscFunctionBegin; 147014354c39SJed Brown PetscValidHeaderSpecific(dmi,DM_CLASSID,1); 14712ce3a92bSJed Brown ierr = DMSetUp(dmi);CHKERRQ(ierr); 14722ee06e3bSJed Brown if (comm == MPI_COMM_NULL) { 147325296bd5SBarry Smith ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr); 147425296bd5SBarry Smith } 147514354c39SJed Brown next = com->next; 147614354c39SJed Brown ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr); 147714354c39SJed Brown 147814354c39SJed Brown /* loop over packed objects, handling one at at time */ 147914354c39SJed Brown while (next) { 148014354c39SJed Brown ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr); 148114354c39SJed Brown ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr); 148214354c39SJed Brown ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 148314354c39SJed Brown next = next->next; 148414354c39SJed Brown } 148514354c39SJed Brown PetscFunctionReturn(0); 148614354c39SJed Brown } 148747c6ae99SBarry Smith 1488e727c939SJed Brown PetscErrorCode DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v) 148947c6ae99SBarry Smith { 149047c6ae99SBarry Smith PetscErrorCode ierr; 14919ae5db72SJed Brown PetscInt m,n,M,N,nDM,i; 149247c6ae99SBarry Smith struct DMCompositeLink *nextc; 149347c6ae99SBarry Smith struct DMCompositeLink *nextf; 149425296bd5SBarry Smith Vec gcoarse,gfine,*vecs; 149547c6ae99SBarry Smith DM_Composite *comcoarse = (DM_Composite*)coarse->data; 149647c6ae99SBarry Smith DM_Composite *comfine = (DM_Composite*)fine->data; 14979ae5db72SJed Brown Mat *mats; 149847c6ae99SBarry Smith 149947c6ae99SBarry Smith PetscFunctionBegin; 150047c6ae99SBarry Smith PetscValidHeaderSpecific(coarse,DM_CLASSID,1); 150147c6ae99SBarry Smith PetscValidHeaderSpecific(fine,DM_CLASSID,2); 1502f692024eSJed Brown ierr = DMSetUp(coarse);CHKERRQ(ierr); 1503f692024eSJed Brown ierr = DMSetUp(fine);CHKERRQ(ierr); 150447c6ae99SBarry Smith /* use global vectors only for determining matrix layout */ 15059ae5db72SJed Brown ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr); 15069ae5db72SJed Brown ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr); 150747c6ae99SBarry Smith ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr); 150847c6ae99SBarry Smith ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr); 150947c6ae99SBarry Smith ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr); 151047c6ae99SBarry Smith ierr = VecGetSize(gfine,&M);CHKERRQ(ierr); 15119ae5db72SJed Brown ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr); 15129ae5db72SJed Brown ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr); 151347c6ae99SBarry Smith 15149ae5db72SJed Brown nDM = comfine->nDM; 1515ce94432eSBarry Smith if (nDM != comcoarse->nDM) SETERRQ2(PetscObjectComm((PetscObject)fine),PETSC_ERR_ARG_INCOMP,"Fine DMComposite has %D entries, but coarse has %D",nDM,comcoarse->nDM); 15161795a4d1SJed Brown ierr = PetscCalloc1(nDM*nDM,&mats);CHKERRQ(ierr); 151725296bd5SBarry Smith if (v) { 15181795a4d1SJed Brown ierr = PetscCalloc1(nDM,&vecs);CHKERRQ(ierr); 151925296bd5SBarry Smith } 152047c6ae99SBarry Smith 152147c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 15229ae5db72SJed Brown for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) { 152325296bd5SBarry Smith if (!v) { 15240298fd71SBarry Smith ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL);CHKERRQ(ierr); 152525296bd5SBarry Smith } else { 152625296bd5SBarry Smith ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr); 152725296bd5SBarry Smith } 152847c6ae99SBarry Smith } 1529ce94432eSBarry Smith ierr = MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A);CHKERRQ(ierr); 153025296bd5SBarry Smith if (v) { 1531ce94432eSBarry Smith ierr = VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v);CHKERRQ(ierr); 153225296bd5SBarry Smith } 15339ae5db72SJed Brown for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);} 15349ae5db72SJed Brown ierr = PetscFree(mats);CHKERRQ(ierr); 153525296bd5SBarry Smith if (v) { 153625296bd5SBarry Smith for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);} 153725296bd5SBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 153825296bd5SBarry Smith } 153947c6ae99SBarry Smith PetscFunctionReturn(0); 154047c6ae99SBarry Smith } 154147c6ae99SBarry Smith 1542184d77edSJed Brown static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm) 15431411c6eeSJed Brown { 15441411c6eeSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 15451411c6eeSJed Brown ISLocalToGlobalMapping *ltogs; 1546f7efa3c7SJed Brown PetscInt i; 15471411c6eeSJed Brown PetscErrorCode ierr; 15481411c6eeSJed Brown 15491411c6eeSJed Brown PetscFunctionBegin; 15501411c6eeSJed Brown /* Set the ISLocalToGlobalMapping on the new matrix */ 15511411c6eeSJed Brown ierr = DMCompositeGetISLocalToGlobalMappings(dm,<ogs);CHKERRQ(ierr); 1552ce94432eSBarry Smith ierr = ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr); 15539ae5db72SJed Brown for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(<ogs[i]);CHKERRQ(ierr);} 15541411c6eeSJed Brown ierr = PetscFree(ltogs);CHKERRQ(ierr); 15551411c6eeSJed Brown PetscFunctionReturn(0); 15561411c6eeSJed Brown } 15571411c6eeSJed Brown 15581411c6eeSJed Brown 1559b412c318SBarry Smith PetscErrorCode DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring) 156047c6ae99SBarry Smith { 156147c6ae99SBarry Smith PetscErrorCode ierr; 156247c6ae99SBarry Smith PetscInt n,i,cnt; 156347c6ae99SBarry Smith ISColoringValue *colors; 156447c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 156547c6ae99SBarry Smith ISColoringValue maxcol = 0; 156647c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 156747c6ae99SBarry Smith 156847c6ae99SBarry Smith PetscFunctionBegin; 156947c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 15705bdb020cSBarry Smith if (ctype == IS_COLORING_LOCAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported"); 1571e3247f34SBarry Smith else if (ctype == IS_COLORING_GLOBAL) { 157247c6ae99SBarry Smith n = com->n; 1573ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType"); 1574785e854fSJed Brown ierr = PetscMalloc1(n,&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */ 157547c6ae99SBarry Smith 1576c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL);CHKERRQ(ierr); 157747c6ae99SBarry Smith if (dense) { 157847c6ae99SBarry Smith for (i=0; i<n; i++) { 157947c6ae99SBarry Smith colors[i] = (ISColoringValue)(com->rstart + i); 158047c6ae99SBarry Smith } 158147c6ae99SBarry Smith maxcol = com->N; 158247c6ae99SBarry Smith } else { 158347c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 158447c6ae99SBarry Smith PetscMPIInt rank; 158547c6ae99SBarry Smith 1586ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 158747c6ae99SBarry Smith cnt = 0; 158847c6ae99SBarry Smith while (next) { 158947c6ae99SBarry Smith ISColoring lcoloring; 159047c6ae99SBarry Smith 1591b412c318SBarry Smith ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring);CHKERRQ(ierr); 159247c6ae99SBarry Smith for (i=0; i<lcoloring->N; i++) { 159347c6ae99SBarry Smith colors[cnt++] = maxcol + lcoloring->colors[i]; 159447c6ae99SBarry Smith } 159547c6ae99SBarry Smith maxcol += lcoloring->n; 1596fcfd50ebSBarry Smith ierr = ISColoringDestroy(&lcoloring);CHKERRQ(ierr); 159747c6ae99SBarry Smith next = next->next; 159847c6ae99SBarry Smith } 159947c6ae99SBarry Smith } 1600aaf3ff59SMatthew G. Knepley ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,PETSC_OWN_POINTER,coloring);CHKERRQ(ierr); 160147c6ae99SBarry Smith PetscFunctionReturn(0); 160247c6ae99SBarry Smith } 160347c6ae99SBarry Smith 16047087cfbeSBarry Smith PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec) 160547c6ae99SBarry Smith { 160647c6ae99SBarry Smith PetscErrorCode ierr; 160747c6ae99SBarry Smith struct DMCompositeLink *next; 160847c6ae99SBarry Smith PetscScalar *garray,*larray; 160947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 161047c6ae99SBarry Smith 161147c6ae99SBarry Smith PetscFunctionBegin; 161247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 161347c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 161439d35262SVincent Le Chenadec 161547c6ae99SBarry Smith if (!com->setup) { 1616d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 161747c6ae99SBarry Smith } 161839d35262SVincent Le Chenadec 161947c6ae99SBarry Smith ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr); 162047c6ae99SBarry Smith ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr); 162147c6ae99SBarry Smith 162247c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 162339d35262SVincent Le Chenadec next = com->next; 162447c6ae99SBarry Smith while (next) { 162547c6ae99SBarry Smith Vec local,global; 162647c6ae99SBarry Smith PetscInt N; 162747c6ae99SBarry Smith 162847c6ae99SBarry Smith ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 162947c6ae99SBarry Smith ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr); 163047c6ae99SBarry Smith ierr = VecPlaceArray(global,garray);CHKERRQ(ierr); 163147c6ae99SBarry Smith ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr); 163247c6ae99SBarry Smith ierr = VecPlaceArray(local,larray);CHKERRQ(ierr); 163347c6ae99SBarry Smith ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr); 163447c6ae99SBarry Smith ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr); 163547c6ae99SBarry Smith ierr = VecResetArray(global);CHKERRQ(ierr); 163647c6ae99SBarry Smith ierr = VecResetArray(local);CHKERRQ(ierr); 163747c6ae99SBarry Smith ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 163839d35262SVincent Le Chenadec ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr); 163939d35262SVincent Le Chenadec 164006ebdd98SJed Brown larray += next->nlocal; 164139d35262SVincent Le Chenadec garray += next->n; 164247c6ae99SBarry Smith next = next->next; 164347c6ae99SBarry Smith } 164447c6ae99SBarry Smith 16450298fd71SBarry Smith ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr); 16460298fd71SBarry Smith ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr); 164747c6ae99SBarry Smith PetscFunctionReturn(0); 164847c6ae99SBarry Smith } 164947c6ae99SBarry Smith 16507087cfbeSBarry Smith PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec) 16510c010503SBarry Smith { 16520c010503SBarry Smith PetscFunctionBegin; 165339d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 165439d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 165539d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,4); 165639d35262SVincent Le Chenadec PetscFunctionReturn(0); 165739d35262SVincent Le Chenadec } 165839d35262SVincent Le Chenadec 165939d35262SVincent Le Chenadec PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec) 166039d35262SVincent Le Chenadec { 166139d35262SVincent Le Chenadec PetscErrorCode ierr; 166239d35262SVincent Le Chenadec struct DMCompositeLink *next; 166339d35262SVincent Le Chenadec PetscScalar *larray,*garray; 166439d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite*)dm->data; 166539d35262SVincent Le Chenadec 166639d35262SVincent Le Chenadec PetscFunctionBegin; 166739d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 166839d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,2); 166939d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,4); 167039d35262SVincent Le Chenadec 167139d35262SVincent Le Chenadec if (!com->setup) { 167239d35262SVincent Le Chenadec ierr = DMSetUp(dm);CHKERRQ(ierr); 167339d35262SVincent Le Chenadec } 167439d35262SVincent Le Chenadec 167539d35262SVincent Le Chenadec ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr); 167639d35262SVincent Le Chenadec ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr); 167739d35262SVincent Le Chenadec 167839d35262SVincent Le Chenadec /* loop over packed objects, handling one at at time */ 167939d35262SVincent Le Chenadec next = com->next; 168039d35262SVincent Le Chenadec while (next) { 168139d35262SVincent Le Chenadec Vec global,local; 168239d35262SVincent Le Chenadec 168339d35262SVincent Le Chenadec ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr); 168439d35262SVincent Le Chenadec ierr = VecPlaceArray(local,larray);CHKERRQ(ierr); 168539d35262SVincent Le Chenadec ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 168639d35262SVincent Le Chenadec ierr = VecPlaceArray(global,garray);CHKERRQ(ierr); 168739d35262SVincent Le Chenadec ierr = DMLocalToGlobalBegin(next->dm,local,mode,global);CHKERRQ(ierr); 168839d35262SVincent Le Chenadec ierr = DMLocalToGlobalEnd(next->dm,local,mode,global);CHKERRQ(ierr); 168939d35262SVincent Le Chenadec ierr = VecResetArray(local);CHKERRQ(ierr); 169039d35262SVincent Le Chenadec ierr = VecResetArray(global);CHKERRQ(ierr); 169139d35262SVincent Le Chenadec ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 169239d35262SVincent Le Chenadec ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr); 169339d35262SVincent Le Chenadec 169439d35262SVincent Le Chenadec garray += next->n; 169539d35262SVincent Le Chenadec larray += next->nlocal; 169639d35262SVincent Le Chenadec next = next->next; 169739d35262SVincent Le Chenadec } 169839d35262SVincent Le Chenadec 169939d35262SVincent Le Chenadec ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr); 170039d35262SVincent Le Chenadec ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr); 170139d35262SVincent Le Chenadec 170239d35262SVincent Le Chenadec PetscFunctionReturn(0); 170339d35262SVincent Le Chenadec } 170439d35262SVincent Le Chenadec 170539d35262SVincent Le Chenadec PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec) 170639d35262SVincent Le Chenadec { 170739d35262SVincent Le Chenadec PetscFunctionBegin; 170839d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 170939d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,2); 171039d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,4); 171139d35262SVincent Le Chenadec PetscFunctionReturn(0); 171239d35262SVincent Le Chenadec } 171339d35262SVincent Le Chenadec 171439d35262SVincent Le Chenadec PetscErrorCode DMLocalToLocalBegin_Composite(DM dm,Vec vec1,InsertMode mode,Vec vec2) 171539d35262SVincent Le Chenadec { 171639d35262SVincent Le Chenadec PetscErrorCode ierr; 171739d35262SVincent Le Chenadec struct DMCompositeLink *next; 171839d35262SVincent Le Chenadec PetscScalar *array1,*array2; 171939d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite*)dm->data; 172039d35262SVincent Le Chenadec 172139d35262SVincent Le Chenadec PetscFunctionBegin; 172239d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 172339d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec1,VEC_CLASSID,2); 172439d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec2,VEC_CLASSID,4); 172539d35262SVincent Le Chenadec 172639d35262SVincent Le Chenadec if (!com->setup) { 172739d35262SVincent Le Chenadec ierr = DMSetUp(dm);CHKERRQ(ierr); 172839d35262SVincent Le Chenadec } 172939d35262SVincent Le Chenadec 173039d35262SVincent Le Chenadec ierr = VecGetArray(vec1,&array1);CHKERRQ(ierr); 173139d35262SVincent Le Chenadec ierr = VecGetArray(vec2,&array2);CHKERRQ(ierr); 173239d35262SVincent Le Chenadec 173339d35262SVincent Le Chenadec /* loop over packed objects, handling one at at time */ 173439d35262SVincent Le Chenadec next = com->next; 173539d35262SVincent Le Chenadec while (next) { 173639d35262SVincent Le Chenadec Vec local1,local2; 173739d35262SVincent Le Chenadec 173839d35262SVincent Le Chenadec ierr = DMGetLocalVector(next->dm,&local1);CHKERRQ(ierr); 173939d35262SVincent Le Chenadec ierr = VecPlaceArray(local1,array1);CHKERRQ(ierr); 174039d35262SVincent Le Chenadec ierr = DMGetLocalVector(next->dm,&local2);CHKERRQ(ierr); 174139d35262SVincent Le Chenadec ierr = VecPlaceArray(local2,array2);CHKERRQ(ierr); 174239d35262SVincent Le Chenadec ierr = DMLocalToLocalBegin(next->dm,local1,mode,local2);CHKERRQ(ierr); 174339d35262SVincent Le Chenadec ierr = DMLocalToLocalEnd(next->dm,local1,mode,local2);CHKERRQ(ierr); 174439d35262SVincent Le Chenadec ierr = VecResetArray(local2);CHKERRQ(ierr); 174539d35262SVincent Le Chenadec ierr = DMRestoreLocalVector(next->dm,&local2);CHKERRQ(ierr); 174639d35262SVincent Le Chenadec ierr = VecResetArray(local1);CHKERRQ(ierr); 174739d35262SVincent Le Chenadec ierr = DMRestoreLocalVector(next->dm,&local1);CHKERRQ(ierr); 174839d35262SVincent Le Chenadec 174939d35262SVincent Le Chenadec array1 += next->nlocal; 175039d35262SVincent Le Chenadec array2 += next->nlocal; 175139d35262SVincent Le Chenadec next = next->next; 175239d35262SVincent Le Chenadec } 175339d35262SVincent Le Chenadec 175439d35262SVincent Le Chenadec ierr = VecRestoreArray(vec1,NULL);CHKERRQ(ierr); 175539d35262SVincent Le Chenadec ierr = VecRestoreArray(vec2,NULL);CHKERRQ(ierr); 175639d35262SVincent Le Chenadec 175739d35262SVincent Le Chenadec PetscFunctionReturn(0); 175839d35262SVincent Le Chenadec } 175939d35262SVincent Le Chenadec 176039d35262SVincent Le Chenadec PetscErrorCode DMLocalToLocalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec) 176139d35262SVincent Le Chenadec { 176239d35262SVincent Le Chenadec PetscFunctionBegin; 176339d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 176439d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,2); 176539d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,4); 17660c010503SBarry Smith PetscFunctionReturn(0); 17670c010503SBarry Smith } 176847c6ae99SBarry Smith 17696ae3a549SBarry Smith /*MC 17706ae3a549SBarry Smith DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs 17716ae3a549SBarry Smith 17726ae3a549SBarry Smith Level: intermediate 17736ae3a549SBarry Smith 17741abcec8cSBarry Smith .seealso: DMType, DM, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate() 17756ae3a549SBarry Smith M*/ 17766ae3a549SBarry Smith 17776ae3a549SBarry Smith 17788cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p) 1779a4121054SBarry Smith { 1780a4121054SBarry Smith PetscErrorCode ierr; 1781a4121054SBarry Smith DM_Composite *com; 1782a4121054SBarry Smith 1783a4121054SBarry Smith PetscFunctionBegin; 1784b00a9115SJed Brown ierr = PetscNewLog(p,&com);CHKERRQ(ierr); 1785a4121054SBarry Smith p->data = com; 1786a4121054SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr); 1787a4121054SBarry Smith com->n = 0; 17887ac2b803SAlex Fikl com->nghost = 0; 17890298fd71SBarry Smith com->next = NULL; 1790a4121054SBarry Smith com->nDM = 0; 1791a4121054SBarry Smith 1792a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Composite; 1793a4121054SBarry Smith p->ops->createlocalvector = DMCreateLocalVector_Composite; 1794184d77edSJed Brown p->ops->getlocaltoglobalmapping = DMGetLocalToGlobalMapping_Composite; 17954d343eeaSMatthew G Knepley p->ops->createfieldis = DMCreateFieldIS_Composite; 179616621825SDmitry Karpeev p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite; 1797a4121054SBarry Smith p->ops->refine = DMRefine_Composite; 179814354c39SJed Brown p->ops->coarsen = DMCoarsen_Composite; 179925296bd5SBarry Smith p->ops->createinterpolation = DMCreateInterpolation_Composite; 180025296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Composite; 1801e727c939SJed Brown p->ops->getcoloring = DMCreateColoring_Composite; 1802a4121054SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite; 1803a4121054SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite; 180439d35262SVincent Le Chenadec p->ops->localtoglobalbegin = DMLocalToGlobalBegin_Composite; 180539d35262SVincent Le Chenadec p->ops->localtoglobalend = DMLocalToGlobalEnd_Composite; 180639d35262SVincent Le Chenadec p->ops->localtolocalbegin = DMLocalToLocalBegin_Composite; 180739d35262SVincent Le Chenadec p->ops->localtolocalend = DMLocalToLocalEnd_Composite; 1808a4121054SBarry Smith p->ops->destroy = DMDestroy_Composite; 1809a4121054SBarry Smith p->ops->view = DMView_Composite; 1810a4121054SBarry Smith p->ops->setup = DMSetUp_Composite; 1811e10fd815SStefano Zampini 1812e10fd815SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)p,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Composite);CHKERRQ(ierr); 1813a4121054SBarry Smith PetscFunctionReturn(0); 1814a4121054SBarry Smith } 1815a4121054SBarry Smith 18161fd49c25SBarry Smith /*@ 18170c010503SBarry Smith DMCompositeCreate - Creates a vector packer, used to generate "composite" 18180c010503SBarry Smith vectors made up of several subvectors. 18190c010503SBarry Smith 18200c010503SBarry Smith Collective on MPI_Comm 182147c6ae99SBarry Smith 182247c6ae99SBarry Smith Input Parameter: 18230c010503SBarry Smith . comm - the processors that will share the global vector 18240c010503SBarry Smith 18250c010503SBarry Smith Output Parameters: 18260c010503SBarry Smith . packer - the packer object 182747c6ae99SBarry Smith 182847c6ae99SBarry Smith Level: advanced 182947c6ae99SBarry Smith 18301abcec8cSBarry Smith .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(), DMCOMPOSITE,DMCreate() 18316eb61c8cSJed Brown DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess() 183247c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 183347c6ae99SBarry Smith 183447c6ae99SBarry Smith @*/ 18357087cfbeSBarry Smith PetscErrorCode DMCompositeCreate(MPI_Comm comm,DM *packer) 183647c6ae99SBarry Smith { 18370c010503SBarry Smith PetscErrorCode ierr; 18380c010503SBarry Smith 183947c6ae99SBarry Smith PetscFunctionBegin; 18400c010503SBarry Smith PetscValidPointer(packer,2); 1841a4121054SBarry Smith ierr = DMCreate(comm,packer);CHKERRQ(ierr); 1842a4121054SBarry Smith ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr); 184347c6ae99SBarry Smith PetscFunctionReturn(0); 184447c6ae99SBarry Smith } 1845