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 2295452b02SPatrick Sanan Notes: 2395452b02SPatrick Sanan See DMSetApplicationContext() and DMGetApplicationContext() for how to get user information into 2447c6ae99SBarry Smith this routine 2547c6ae99SBarry Smith 2647c6ae99SBarry Smith @*/ 277087cfbeSBarry Smith PetscErrorCode DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt)) 2847c6ae99SBarry Smith { 2947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 3071b14b3eSStefano Zampini PetscBool flg; 3171b14b3eSStefano Zampini PetscErrorCode ierr; 3247c6ae99SBarry Smith 3347c6ae99SBarry Smith PetscFunctionBegin; 3471b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 3571b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 3647c6ae99SBarry Smith com->FormCoupleLocations = FormCoupleLocations; 3747c6ae99SBarry Smith PetscFunctionReturn(0); 3847c6ae99SBarry Smith } 3947c6ae99SBarry Smith 406bf464f9SBarry Smith PetscErrorCode DMDestroy_Composite(DM dm) 4147c6ae99SBarry Smith { 4247c6ae99SBarry Smith PetscErrorCode ierr; 4347c6ae99SBarry Smith struct DMCompositeLink *next, *prev; 4447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 4547c6ae99SBarry Smith 4647c6ae99SBarry Smith PetscFunctionBegin; 4747c6ae99SBarry Smith next = com->next; 4847c6ae99SBarry Smith while (next) { 4947c6ae99SBarry Smith prev = next; 5047c6ae99SBarry Smith next = next->next; 51fcfd50ebSBarry Smith ierr = DMDestroy(&prev->dm);CHKERRQ(ierr); 5247c6ae99SBarry Smith ierr = PetscFree(prev->grstarts);CHKERRQ(ierr); 5347c6ae99SBarry Smith ierr = PetscFree(prev);CHKERRQ(ierr); 5447c6ae99SBarry Smith } 55e10fd815SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL);CHKERRQ(ierr); 56435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 57435a35e8SMatthew G Knepley ierr = PetscFree(com);CHKERRQ(ierr); 5847c6ae99SBarry Smith PetscFunctionReturn(0); 5947c6ae99SBarry Smith } 6047c6ae99SBarry Smith 617087cfbeSBarry Smith PetscErrorCode DMView_Composite(DM dm,PetscViewer v) 6247c6ae99SBarry Smith { 6347c6ae99SBarry Smith PetscErrorCode ierr; 6447c6ae99SBarry Smith PetscBool iascii; 6547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 6647c6ae99SBarry Smith 6747c6ae99SBarry Smith PetscFunctionBegin; 68251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 6947c6ae99SBarry Smith if (iascii) { 7047c6ae99SBarry Smith struct DMCompositeLink *lnk = com->next; 7147c6ae99SBarry Smith PetscInt i; 7247c6ae99SBarry Smith 7347c6ae99SBarry Smith ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix ? ((PetscObject)dm)->prefix : "no prefix");CHKERRQ(ierr); 749ae5db72SJed Brown ierr = PetscViewerASCIIPrintf(v," contains %D DMs\n",com->nDM);CHKERRQ(ierr); 7547c6ae99SBarry Smith ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 7647c6ae99SBarry Smith for (i=0; lnk; lnk=lnk->next,i++) { 779ae5db72SJed Brown ierr = PetscViewerASCIIPrintf(v,"Link %D: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr); 7847c6ae99SBarry Smith ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 7947c6ae99SBarry Smith ierr = DMView(lnk->dm,v);CHKERRQ(ierr); 8047c6ae99SBarry Smith ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 8147c6ae99SBarry Smith } 8247c6ae99SBarry Smith ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 8347c6ae99SBarry Smith } 8447c6ae99SBarry Smith PetscFunctionReturn(0); 8547c6ae99SBarry Smith } 8647c6ae99SBarry Smith 8747c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/ 887087cfbeSBarry Smith PetscErrorCode DMSetUp_Composite(DM dm) 8947c6ae99SBarry Smith { 9047c6ae99SBarry Smith PetscErrorCode ierr; 9147c6ae99SBarry Smith PetscInt nprev = 0; 9247c6ae99SBarry Smith PetscMPIInt rank,size; 9347c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 9447c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 9547c6ae99SBarry Smith PetscLayout map; 9647c6ae99SBarry Smith 9747c6ae99SBarry Smith PetscFunctionBegin; 98ce94432eSBarry Smith if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup"); 99ce94432eSBarry Smith ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)dm),&map);CHKERRQ(ierr); 10047c6ae99SBarry Smith ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr); 10147c6ae99SBarry Smith ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr); 10247c6ae99SBarry Smith ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr); 10347c6ae99SBarry Smith ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 10447c6ae99SBarry Smith ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr); 1050298fd71SBarry Smith ierr = PetscLayoutGetRange(map,&com->rstart,NULL);CHKERRQ(ierr); 106fcfd50ebSBarry Smith ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 10747c6ae99SBarry Smith 1089ae5db72SJed Brown /* now set the rstart for each linked vector */ 109ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 110ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 11147c6ae99SBarry Smith while (next) { 11247c6ae99SBarry Smith next->rstart = nprev; 11306ebdd98SJed Brown nprev += next->n; 11447c6ae99SBarry Smith next->grstart = com->rstart + next->rstart; 115785e854fSJed Brown ierr = PetscMalloc1(size,&next->grstarts);CHKERRQ(ierr); 116ce94432eSBarry Smith ierr = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 11747c6ae99SBarry Smith next = next->next; 11847c6ae99SBarry Smith } 11947c6ae99SBarry Smith com->setup = PETSC_TRUE; 12047c6ae99SBarry Smith PetscFunctionReturn(0); 12147c6ae99SBarry Smith } 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/ 12447c6ae99SBarry Smith 12573e31fe2SJed Brown /*@ 12647c6ae99SBarry Smith DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite 12747c6ae99SBarry Smith representation. 12847c6ae99SBarry Smith 12947c6ae99SBarry Smith Not Collective 13047c6ae99SBarry Smith 13147c6ae99SBarry Smith Input Parameter: 13247c6ae99SBarry Smith . dm - the packer object 13347c6ae99SBarry Smith 13447c6ae99SBarry Smith Output Parameter: 13547c6ae99SBarry Smith . nDM - the number of DMs 13647c6ae99SBarry Smith 13747c6ae99SBarry Smith Level: beginner 13847c6ae99SBarry Smith 13947c6ae99SBarry Smith @*/ 1407087cfbeSBarry Smith PetscErrorCode DMCompositeGetNumberDM(DM dm,PetscInt *nDM) 14147c6ae99SBarry Smith { 14247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 14371b14b3eSStefano Zampini PetscBool flg; 14471b14b3eSStefano Zampini PetscErrorCode ierr; 1455fd66863SKarl Rupp 14647c6ae99SBarry Smith PetscFunctionBegin; 14747c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 14871b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 14971b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 15047c6ae99SBarry Smith *nDM = com->nDM; 15147c6ae99SBarry Smith PetscFunctionReturn(0); 15247c6ae99SBarry Smith } 15347c6ae99SBarry Smith 15447c6ae99SBarry Smith 15547c6ae99SBarry Smith /*@C 15647c6ae99SBarry Smith DMCompositeGetAccess - Allows one to access the individual packed vectors in their global 15747c6ae99SBarry Smith representation. 15847c6ae99SBarry Smith 15947c6ae99SBarry Smith Collective on DMComposite 16047c6ae99SBarry Smith 1619ae5db72SJed Brown Input Parameters: 16247c6ae99SBarry Smith + dm - the packer object 1639ae5db72SJed Brown - gvec - the global vector 1649ae5db72SJed Brown 1659ae5db72SJed Brown Output Parameters: 1660298fd71SBarry Smith . Vec* ... - the packed parallel vectors, NULL for those that are not needed 16747c6ae99SBarry Smith 16895452b02SPatrick Sanan Notes: 16995452b02SPatrick Sanan Use DMCompositeRestoreAccess() to return the vectors when you no longer need them 17047c6ae99SBarry Smith 171f73e5cebSJed Brown Fortran Notes: 172f73e5cebSJed Brown 173f73e5cebSJed Brown Fortran callers must use numbered versions of this routine, e.g., DMCompositeGetAccess4(dm,gvec,vec1,vec2,vec3,vec4) 174f73e5cebSJed Brown or use the alternative interface DMCompositeGetAccessArray(). 175f73e5cebSJed Brown 17647c6ae99SBarry Smith Level: advanced 17747c6ae99SBarry Smith 178f73e5cebSJed Brown .seealso: DMCompositeGetEntries(), DMCompositeScatter() 17947c6ae99SBarry Smith @*/ 1807087cfbeSBarry Smith PetscErrorCode DMCompositeGetAccess(DM dm,Vec gvec,...) 18147c6ae99SBarry Smith { 18247c6ae99SBarry Smith va_list Argp; 18347c6ae99SBarry Smith PetscErrorCode ierr; 18447c6ae99SBarry Smith struct DMCompositeLink *next; 18547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 1865edff71fSBarry Smith PetscInt readonly; 18771b14b3eSStefano Zampini PetscBool flg; 18847c6ae99SBarry Smith 18947c6ae99SBarry Smith PetscFunctionBegin; 19047c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 19147c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 19271b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 19371b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 19447c6ae99SBarry Smith next = com->next; 19547c6ae99SBarry Smith if (!com->setup) { 196d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 19747c6ae99SBarry Smith } 19847c6ae99SBarry Smith 1995edff71fSBarry Smith ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr); 20047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 20147c6ae99SBarry Smith va_start(Argp,gvec); 20247c6ae99SBarry Smith while (next) { 20347c6ae99SBarry Smith Vec *vec; 20447c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 2059ae5db72SJed Brown if (vec) { 2069ae5db72SJed Brown ierr = DMGetGlobalVector(next->dm,vec);CHKERRQ(ierr); 2075edff71fSBarry Smith if (readonly) { 2085edff71fSBarry Smith const PetscScalar *array; 2095edff71fSBarry Smith ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr); 2105edff71fSBarry Smith ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr); 2118860a134SJunchao Zhang ierr = VecLockReadPush(*vec);CHKERRQ(ierr); 2125edff71fSBarry Smith ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr); 2135edff71fSBarry Smith } else { 2145edff71fSBarry Smith PetscScalar *array; 2159ae5db72SJed Brown ierr = VecGetArray(gvec,&array);CHKERRQ(ierr); 2169ae5db72SJed Brown ierr = VecPlaceArray(*vec,array+next->rstart);CHKERRQ(ierr); 2179ae5db72SJed Brown ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr); 21847c6ae99SBarry Smith } 2195edff71fSBarry Smith } 22047c6ae99SBarry Smith next = next->next; 22147c6ae99SBarry Smith } 22247c6ae99SBarry Smith va_end(Argp); 22347c6ae99SBarry Smith PetscFunctionReturn(0); 22447c6ae99SBarry Smith } 22547c6ae99SBarry Smith 226f73e5cebSJed Brown /*@C 227f73e5cebSJed Brown DMCompositeGetAccessArray - Allows one to access the individual packed vectors in their global 228f73e5cebSJed Brown representation. 229f73e5cebSJed Brown 230f73e5cebSJed Brown Collective on DMComposite 231f73e5cebSJed Brown 232f73e5cebSJed Brown Input Parameters: 233f73e5cebSJed Brown + dm - the packer object 234f73e5cebSJed Brown . pvec - packed vector 235f73e5cebSJed Brown . nwanted - number of vectors wanted 2360298fd71SBarry Smith - wanted - sorted array of vectors wanted, or NULL to get all vectors 237f73e5cebSJed Brown 238f73e5cebSJed Brown Output Parameters: 239f73e5cebSJed Brown . vecs - array of requested global vectors (must be allocated) 240f73e5cebSJed Brown 24195452b02SPatrick Sanan Notes: 24295452b02SPatrick Sanan Use DMCompositeRestoreAccessArray() to return the vectors when you no longer need them 243f73e5cebSJed Brown 244f73e5cebSJed Brown Level: advanced 245f73e5cebSJed Brown 246f73e5cebSJed Brown .seealso: DMCompositeGetAccess(), DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather() 247f73e5cebSJed Brown @*/ 248f73e5cebSJed Brown PetscErrorCode DMCompositeGetAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 249f73e5cebSJed Brown { 250f73e5cebSJed Brown PetscErrorCode ierr; 251f73e5cebSJed Brown struct DMCompositeLink *link; 252f73e5cebSJed Brown PetscInt i,wnum; 253f73e5cebSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 254bee642f7SBarry Smith PetscInt readonly; 25571b14b3eSStefano Zampini PetscBool flg; 256f73e5cebSJed Brown 257f73e5cebSJed Brown PetscFunctionBegin; 258f73e5cebSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 259f73e5cebSJed Brown PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 26071b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 26171b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 262f73e5cebSJed Brown if (!com->setup) { 263f73e5cebSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 264f73e5cebSJed Brown } 265f73e5cebSJed Brown 266bee642f7SBarry Smith ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr); 267f73e5cebSJed Brown for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 268f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 269f73e5cebSJed Brown Vec v; 270f73e5cebSJed Brown ierr = DMGetGlobalVector(link->dm,&v);CHKERRQ(ierr); 271bee642f7SBarry Smith if (readonly) { 272bee642f7SBarry Smith const PetscScalar *array; 273bee642f7SBarry Smith ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr); 274bee642f7SBarry Smith ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr); 2758860a134SJunchao Zhang ierr = VecLockReadPush(v);CHKERRQ(ierr); 276bee642f7SBarry Smith ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr); 277bee642f7SBarry Smith } else { 278bee642f7SBarry Smith PetscScalar *array; 279f73e5cebSJed Brown ierr = VecGetArray(pvec,&array);CHKERRQ(ierr); 280f73e5cebSJed Brown ierr = VecPlaceArray(v,array+link->rstart);CHKERRQ(ierr); 281f73e5cebSJed Brown ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr); 282bee642f7SBarry Smith } 283f73e5cebSJed Brown vecs[wnum++] = v; 284f73e5cebSJed Brown } 285f73e5cebSJed Brown } 286f73e5cebSJed Brown PetscFunctionReturn(0); 287f73e5cebSJed Brown } 288f73e5cebSJed Brown 2897ac2b803SAlex Fikl /*@C 2907ac2b803SAlex Fikl DMCompositeGetLocalAccessArray - Allows one to access the individual 2917ac2b803SAlex Fikl packed vectors in their local representation. 2927ac2b803SAlex Fikl 2937ac2b803SAlex Fikl Collective on DMComposite. 2947ac2b803SAlex Fikl 2957ac2b803SAlex Fikl Input Parameters: 2967ac2b803SAlex Fikl + dm - the packer object 2977ac2b803SAlex Fikl . pvec - packed vector 2987ac2b803SAlex Fikl . nwanted - number of vectors wanted 2997ac2b803SAlex Fikl - wanted - sorted array of vectors wanted, or NULL to get all vectors 3007ac2b803SAlex Fikl 3017ac2b803SAlex Fikl Output Parameters: 3027ac2b803SAlex Fikl . vecs - array of requested local vectors (must be allocated) 3037ac2b803SAlex Fikl 30495452b02SPatrick Sanan Notes: 30595452b02SPatrick Sanan Use DMCompositeRestoreLocalAccessArray() to return the vectors 3067ac2b803SAlex Fikl when you no longer need them. 3077ac2b803SAlex Fikl 3087ac2b803SAlex Fikl Level: advanced 3097ac2b803SAlex Fikl 3107ac2b803SAlex Fikl .seealso: DMCompositeRestoreLocalAccessArray(), DMCompositeGetAccess(), 3117ac2b803SAlex Fikl DMCompositeGetEntries(), DMCompositeScatter(), DMCompositeGather() 3127ac2b803SAlex Fikl @*/ 3137ac2b803SAlex Fikl PetscErrorCode DMCompositeGetLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 3147ac2b803SAlex Fikl { 3157ac2b803SAlex Fikl PetscErrorCode ierr; 3167ac2b803SAlex Fikl struct DMCompositeLink *link; 3177ac2b803SAlex Fikl PetscInt i,wnum; 3187ac2b803SAlex Fikl DM_Composite *com = (DM_Composite*)dm->data; 3197ac2b803SAlex Fikl PetscInt readonly; 3207ac2b803SAlex Fikl PetscInt nlocal = 0; 32171b14b3eSStefano Zampini PetscBool flg; 3227ac2b803SAlex Fikl 3237ac2b803SAlex Fikl PetscFunctionBegin; 3247ac2b803SAlex Fikl PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3257ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 32671b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 32771b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 3287ac2b803SAlex Fikl if (!com->setup) { 3297ac2b803SAlex Fikl ierr = DMSetUp(dm);CHKERRQ(ierr); 3307ac2b803SAlex Fikl } 3317ac2b803SAlex Fikl 3327ac2b803SAlex Fikl ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr); 3337ac2b803SAlex Fikl for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 3347ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 3357ac2b803SAlex Fikl Vec v; 3367ac2b803SAlex Fikl ierr = DMGetLocalVector(link->dm,&v);CHKERRQ(ierr); 3377ac2b803SAlex Fikl if (readonly) { 3387ac2b803SAlex Fikl const PetscScalar *array; 3397ac2b803SAlex Fikl ierr = VecGetArrayRead(pvec,&array);CHKERRQ(ierr); 3407ac2b803SAlex Fikl ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr); 3418860a134SJunchao Zhang ierr = VecLockReadPush(v);CHKERRQ(ierr); 3427ac2b803SAlex Fikl ierr = VecRestoreArrayRead(pvec,&array);CHKERRQ(ierr); 3437ac2b803SAlex Fikl } else { 3447ac2b803SAlex Fikl PetscScalar *array; 3457ac2b803SAlex Fikl ierr = VecGetArray(pvec,&array);CHKERRQ(ierr); 3467ac2b803SAlex Fikl ierr = VecPlaceArray(v,array+nlocal);CHKERRQ(ierr); 3477ac2b803SAlex Fikl ierr = VecRestoreArray(pvec,&array);CHKERRQ(ierr); 3487ac2b803SAlex Fikl } 3497ac2b803SAlex Fikl vecs[wnum++] = v; 3507ac2b803SAlex Fikl } 3517ac2b803SAlex Fikl 3527ac2b803SAlex Fikl nlocal += link->nlocal; 3537ac2b803SAlex Fikl } 3547ac2b803SAlex Fikl 3557ac2b803SAlex Fikl PetscFunctionReturn(0); 3567ac2b803SAlex Fikl } 3577ac2b803SAlex Fikl 35847c6ae99SBarry Smith /*@C 359aa219208SBarry Smith DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess() 36047c6ae99SBarry Smith representation. 36147c6ae99SBarry Smith 36247c6ae99SBarry Smith Collective on DMComposite 36347c6ae99SBarry Smith 3649ae5db72SJed Brown Input Parameters: 36547c6ae99SBarry Smith + dm - the packer object 36647c6ae99SBarry Smith . gvec - the global vector 3670298fd71SBarry Smith - Vec* ... - the individual parallel vectors, NULL for those that are not needed 36847c6ae99SBarry Smith 36947c6ae99SBarry Smith Level: advanced 37047c6ae99SBarry Smith 3719ae5db72SJed Brown .seealso DMCompositeAddDM(), DMCreateGlobalVector(), 3726eb61c8cSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(), 373aa219208SBarry Smith DMCompositeRestoreAccess(), DMCompositeGetAccess() 37447c6ae99SBarry Smith 37547c6ae99SBarry Smith @*/ 3767087cfbeSBarry Smith PetscErrorCode DMCompositeRestoreAccess(DM dm,Vec gvec,...) 37747c6ae99SBarry Smith { 37847c6ae99SBarry Smith va_list Argp; 37947c6ae99SBarry Smith PetscErrorCode ierr; 38047c6ae99SBarry Smith struct DMCompositeLink *next; 38147c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 3825edff71fSBarry Smith PetscInt readonly; 38371b14b3eSStefano Zampini PetscBool flg; 38447c6ae99SBarry Smith 38547c6ae99SBarry Smith PetscFunctionBegin; 38647c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 38747c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 38871b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 38971b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 39047c6ae99SBarry Smith next = com->next; 39147c6ae99SBarry Smith if (!com->setup) { 392d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 39347c6ae99SBarry Smith } 39447c6ae99SBarry Smith 3955edff71fSBarry Smith ierr = VecLockGet(gvec,&readonly);CHKERRQ(ierr); 39647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 39747c6ae99SBarry Smith va_start(Argp,gvec); 39847c6ae99SBarry Smith while (next) { 39947c6ae99SBarry Smith Vec *vec; 40047c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 4019ae5db72SJed Brown if (vec) { 4029ae5db72SJed Brown ierr = VecResetArray(*vec);CHKERRQ(ierr); 4035edff71fSBarry Smith if (readonly) { 4048860a134SJunchao Zhang ierr = VecLockReadPop(*vec);CHKERRQ(ierr); 4055edff71fSBarry Smith } 406bee642f7SBarry Smith ierr = DMRestoreGlobalVector(next->dm,vec);CHKERRQ(ierr); 40747c6ae99SBarry Smith } 40847c6ae99SBarry Smith next = next->next; 40947c6ae99SBarry Smith } 41047c6ae99SBarry Smith va_end(Argp); 41147c6ae99SBarry Smith PetscFunctionReturn(0); 41247c6ae99SBarry Smith } 41347c6ae99SBarry Smith 414f73e5cebSJed Brown /*@C 415f73e5cebSJed Brown DMCompositeRestoreAccessArray - Returns the vectors obtained with DMCompositeGetAccessArray() 416f73e5cebSJed Brown 417f73e5cebSJed Brown Collective on DMComposite 418f73e5cebSJed Brown 419f73e5cebSJed Brown Input Parameters: 420f73e5cebSJed Brown + dm - the packer object 421f73e5cebSJed Brown . pvec - packed vector 422f73e5cebSJed Brown . nwanted - number of vectors wanted 4230298fd71SBarry Smith . wanted - sorted array of vectors wanted, or NULL to get all vectors 424f73e5cebSJed Brown - vecs - array of global vectors to return 425f73e5cebSJed Brown 426f73e5cebSJed Brown Level: advanced 427f73e5cebSJed Brown 428f73e5cebSJed Brown .seealso: DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), DMCompositeScatter(), DMCompositeGather() 429f73e5cebSJed Brown @*/ 430f73e5cebSJed Brown PetscErrorCode DMCompositeRestoreAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 431f73e5cebSJed Brown { 432f73e5cebSJed Brown PetscErrorCode ierr; 433f73e5cebSJed Brown struct DMCompositeLink *link; 434f73e5cebSJed Brown PetscInt i,wnum; 435f73e5cebSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 436bee642f7SBarry Smith PetscInt readonly; 43771b14b3eSStefano Zampini PetscBool flg; 438f73e5cebSJed Brown 439f73e5cebSJed Brown PetscFunctionBegin; 440f73e5cebSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 441f73e5cebSJed Brown PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 44271b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 44371b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 444f73e5cebSJed Brown if (!com->setup) { 445f73e5cebSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 446f73e5cebSJed Brown } 447f73e5cebSJed Brown 448bee642f7SBarry Smith ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr); 449f73e5cebSJed Brown for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 450f73e5cebSJed Brown if (!wanted || i == wanted[wnum]) { 451f73e5cebSJed Brown ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr); 452bee642f7SBarry Smith if (readonly) { 4538860a134SJunchao Zhang ierr = VecLockReadPop(vecs[wnum]);CHKERRQ(ierr); 454bee642f7SBarry Smith } 455f73e5cebSJed Brown ierr = DMRestoreGlobalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr); 456f73e5cebSJed Brown wnum++; 457f73e5cebSJed Brown } 458f73e5cebSJed Brown } 459f73e5cebSJed Brown PetscFunctionReturn(0); 460f73e5cebSJed Brown } 461f73e5cebSJed Brown 4627ac2b803SAlex Fikl /*@C 4637ac2b803SAlex Fikl DMCompositeRestoreLocalAccessArray - Returns the vectors obtained with DMCompositeGetLocalAccessArray(). 4647ac2b803SAlex Fikl 4657ac2b803SAlex Fikl Collective on DMComposite. 4667ac2b803SAlex Fikl 4677ac2b803SAlex Fikl Input Parameters: 4687ac2b803SAlex Fikl + dm - the packer object 4697ac2b803SAlex Fikl . pvec - packed vector 4707ac2b803SAlex Fikl . nwanted - number of vectors wanted 4717ac2b803SAlex Fikl . wanted - sorted array of vectors wanted, or NULL to restore all vectors 4727ac2b803SAlex Fikl - vecs - array of local vectors to return 4737ac2b803SAlex Fikl 4747ac2b803SAlex Fikl Level: advanced 4757ac2b803SAlex Fikl 4767ac2b803SAlex Fikl Notes: 4777ac2b803SAlex Fikl nwanted and wanted must match the values given to DMCompositeGetLocalAccessArray() 4787ac2b803SAlex Fikl otherwise the call will fail. 4797ac2b803SAlex Fikl 4807ac2b803SAlex Fikl .seealso: DMCompositeGetLocalAccessArray(), DMCompositeRestoreAccessArray(), 4817ac2b803SAlex Fikl DMCompositeRestoreAccess(), DMCompositeRestoreEntries(), 4827ac2b803SAlex Fikl DMCompositeScatter(), DMCompositeGather() 4837ac2b803SAlex Fikl @*/ 4847ac2b803SAlex Fikl PetscErrorCode DMCompositeRestoreLocalAccessArray(DM dm,Vec pvec,PetscInt nwanted,const PetscInt *wanted,Vec *vecs) 4857ac2b803SAlex Fikl { 4867ac2b803SAlex Fikl PetscErrorCode ierr; 4877ac2b803SAlex Fikl struct DMCompositeLink *link; 4887ac2b803SAlex Fikl PetscInt i,wnum; 4897ac2b803SAlex Fikl DM_Composite *com = (DM_Composite*)dm->data; 4907ac2b803SAlex Fikl PetscInt readonly; 49171b14b3eSStefano Zampini PetscBool flg; 4927ac2b803SAlex Fikl 4937ac2b803SAlex Fikl PetscFunctionBegin; 4947ac2b803SAlex Fikl PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4957ac2b803SAlex Fikl PetscValidHeaderSpecific(pvec,VEC_CLASSID,2); 49671b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 49771b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 4987ac2b803SAlex Fikl if (!com->setup) { 4997ac2b803SAlex Fikl ierr = DMSetUp(dm);CHKERRQ(ierr); 5007ac2b803SAlex Fikl } 5017ac2b803SAlex Fikl 5027ac2b803SAlex Fikl ierr = VecLockGet(pvec,&readonly);CHKERRQ(ierr); 5037ac2b803SAlex Fikl for (i=0,wnum=0,link=com->next; link && wnum<nwanted; i++,link=link->next) { 5047ac2b803SAlex Fikl if (!wanted || i == wanted[wnum]) { 5057ac2b803SAlex Fikl ierr = VecResetArray(vecs[wnum]);CHKERRQ(ierr); 5067ac2b803SAlex Fikl if (readonly) { 5078860a134SJunchao Zhang ierr = VecLockReadPop(vecs[wnum]);CHKERRQ(ierr); 5087ac2b803SAlex Fikl } 5097ac2b803SAlex Fikl ierr = DMRestoreLocalVector(link->dm,&vecs[wnum]);CHKERRQ(ierr); 5107ac2b803SAlex Fikl wnum++; 5117ac2b803SAlex Fikl } 5127ac2b803SAlex Fikl } 5137ac2b803SAlex Fikl PetscFunctionReturn(0); 5147ac2b803SAlex Fikl } 5157ac2b803SAlex Fikl 51647c6ae99SBarry Smith /*@C 51747c6ae99SBarry Smith DMCompositeScatter - Scatters from a global packed vector into its individual local vectors 51847c6ae99SBarry Smith 51947c6ae99SBarry Smith Collective on DMComposite 52047c6ae99SBarry Smith 5219ae5db72SJed Brown Input Parameters: 52247c6ae99SBarry Smith + dm - the packer object 52347c6ae99SBarry Smith . gvec - the global vector 5240298fd71SBarry Smith - Vec ... - the individual sequential vectors, NULL for those that are not needed 52547c6ae99SBarry Smith 52647c6ae99SBarry Smith Level: advanced 52747c6ae99SBarry Smith 5286f3c3dcfSJed Brown Notes: 5296f3c3dcfSJed Brown DMCompositeScatterArray() is a non-variadic alternative that is often more convenient for library callers and is 5306f3c3dcfSJed Brown accessible from Fortran. 5316f3c3dcfSJed Brown 5329ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 5336eb61c8cSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 53447c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 5356f3c3dcfSJed Brown DMCompositeScatterArray() 53647c6ae99SBarry Smith 53747c6ae99SBarry Smith @*/ 5387087cfbeSBarry Smith PetscErrorCode DMCompositeScatter(DM dm,Vec gvec,...) 53947c6ae99SBarry Smith { 54047c6ae99SBarry Smith va_list Argp; 54147c6ae99SBarry Smith PetscErrorCode ierr; 54247c6ae99SBarry Smith struct DMCompositeLink *next; 5438fd8f222SJed Brown PetscInt cnt; 54447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 54571b14b3eSStefano Zampini PetscBool flg; 54647c6ae99SBarry Smith 54747c6ae99SBarry Smith PetscFunctionBegin; 54847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 54947c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 55071b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 55171b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 55247c6ae99SBarry Smith if (!com->setup) { 553d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 55447c6ae99SBarry Smith } 55547c6ae99SBarry Smith 55647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 55747c6ae99SBarry Smith va_start(Argp,gvec); 5588fd8f222SJed Brown for (cnt=3,next=com->next; next; cnt++,next=next->next) { 5599ae5db72SJed Brown Vec local; 5609ae5db72SJed Brown local = va_arg(Argp, Vec); 5619ae5db72SJed Brown if (local) { 5629ae5db72SJed Brown Vec global; 5635edff71fSBarry Smith const PetscScalar *array; 5649ae5db72SJed Brown PetscValidHeaderSpecific(local,VEC_CLASSID,cnt); 5659ae5db72SJed Brown ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 5665edff71fSBarry Smith ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr); 5679ae5db72SJed Brown ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr); 5689ae5db72SJed Brown ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr); 5699ae5db72SJed Brown ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,local);CHKERRQ(ierr); 5705edff71fSBarry Smith ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr); 5719ae5db72SJed Brown ierr = VecResetArray(global);CHKERRQ(ierr); 5729ae5db72SJed Brown ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 57347c6ae99SBarry Smith } 57447c6ae99SBarry Smith } 57547c6ae99SBarry Smith va_end(Argp); 57647c6ae99SBarry Smith PetscFunctionReturn(0); 57747c6ae99SBarry Smith } 57847c6ae99SBarry Smith 5796f3c3dcfSJed Brown /*@ 5806f3c3dcfSJed Brown DMCompositeScatterArray - Scatters from a global packed vector into its individual local vectors 5816f3c3dcfSJed Brown 5826f3c3dcfSJed Brown Collective on DMComposite 5836f3c3dcfSJed Brown 5846f3c3dcfSJed Brown Input Parameters: 5856f3c3dcfSJed Brown + dm - the packer object 5866f3c3dcfSJed Brown . gvec - the global vector 587*a2b725a8SWilliam Gropp - lvecs - array of local vectors, NULL for any that are not needed 5886f3c3dcfSJed Brown 5896f3c3dcfSJed Brown Level: advanced 5906f3c3dcfSJed Brown 5916f3c3dcfSJed Brown Note: 592907376e6SBarry Smith This is a non-variadic alternative to DMCompositeScatter() 5936f3c3dcfSJed Brown 5946f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector() 5956f3c3dcfSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 5966f3c3dcfSJed Brown DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 5976f3c3dcfSJed Brown 5986f3c3dcfSJed Brown @*/ 5996f3c3dcfSJed Brown PetscErrorCode DMCompositeScatterArray(DM dm,Vec gvec,Vec *lvecs) 6006f3c3dcfSJed Brown { 6016f3c3dcfSJed Brown PetscErrorCode ierr; 6026f3c3dcfSJed Brown struct DMCompositeLink *next; 6036f3c3dcfSJed Brown PetscInt i; 6046f3c3dcfSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 60571b14b3eSStefano Zampini PetscBool flg; 6066f3c3dcfSJed Brown 6076f3c3dcfSJed Brown PetscFunctionBegin; 6086f3c3dcfSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6096f3c3dcfSJed Brown PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 61071b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 61171b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 6126f3c3dcfSJed Brown if (!com->setup) { 6136f3c3dcfSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 6146f3c3dcfSJed Brown } 6156f3c3dcfSJed Brown 6166f3c3dcfSJed Brown /* loop over packed objects, handling one at at time */ 6176f3c3dcfSJed Brown for (i=0,next=com->next; next; next=next->next,i++) { 6186f3c3dcfSJed Brown if (lvecs[i]) { 6196f3c3dcfSJed Brown Vec global; 620c5d31e75SLisandro Dalcin const PetscScalar *array; 6216f3c3dcfSJed Brown PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3); 6226f3c3dcfSJed Brown ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 623c5d31e75SLisandro Dalcin ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr); 624c5d31e75SLisandro Dalcin ierr = VecPlaceArray(global,(PetscScalar*)array+next->rstart);CHKERRQ(ierr); 6256f3c3dcfSJed Brown ierr = DMGlobalToLocalBegin(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr); 6266f3c3dcfSJed Brown ierr = DMGlobalToLocalEnd(next->dm,global,INSERT_VALUES,lvecs[i]);CHKERRQ(ierr); 627c5d31e75SLisandro Dalcin ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr); 6286f3c3dcfSJed Brown ierr = VecResetArray(global);CHKERRQ(ierr); 6296f3c3dcfSJed Brown ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 6306f3c3dcfSJed Brown } 6316f3c3dcfSJed Brown } 6326f3c3dcfSJed Brown PetscFunctionReturn(0); 6336f3c3dcfSJed Brown } 6346f3c3dcfSJed Brown 63547c6ae99SBarry Smith /*@C 63647c6ae99SBarry Smith DMCompositeGather - Gathers into a global packed vector from its individual local vectors 63747c6ae99SBarry Smith 63847c6ae99SBarry Smith Collective on DMComposite 63947c6ae99SBarry Smith 64047c6ae99SBarry Smith Input Parameter: 64147c6ae99SBarry Smith + dm - the packer object 64247c6ae99SBarry Smith . gvec - the global vector 643907376e6SBarry Smith . imode - INSERT_VALUES or ADD_VALUES 6440298fd71SBarry Smith - Vec ... - the individual sequential vectors, NULL for any that are not needed 64547c6ae99SBarry Smith 64647c6ae99SBarry Smith Level: advanced 64747c6ae99SBarry Smith 648f5f57ec0SBarry Smith Not available from Fortran, Fortran users can use DMCompositeGatherArray() 649f5f57ec0SBarry Smith 6509ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 6516eb61c8cSJed Brown DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 65247c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 65347c6ae99SBarry Smith 65447c6ae99SBarry Smith @*/ 6551dac896bSSatish Balay PetscErrorCode DMCompositeGather(DM dm,InsertMode imode,Vec gvec,...) 65647c6ae99SBarry Smith { 65747c6ae99SBarry Smith va_list Argp; 65847c6ae99SBarry Smith PetscErrorCode ierr; 65947c6ae99SBarry Smith struct DMCompositeLink *next; 66047c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 6618fd8f222SJed Brown PetscInt cnt; 66271b14b3eSStefano Zampini PetscBool flg; 66347c6ae99SBarry Smith 66447c6ae99SBarry Smith PetscFunctionBegin; 66547c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 66647c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 66771b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 66871b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 66947c6ae99SBarry Smith if (!com->setup) { 670d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 67147c6ae99SBarry Smith } 67247c6ae99SBarry Smith 67347c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 6741dac896bSSatish Balay va_start(Argp,gvec); 6758fd8f222SJed Brown for (cnt=3,next=com->next; next; cnt++,next=next->next) { 6769ae5db72SJed Brown Vec local; 6779ae5db72SJed Brown local = va_arg(Argp, Vec); 6789ae5db72SJed Brown if (local) { 67947c6ae99SBarry Smith PetscScalar *array; 6809ae5db72SJed Brown Vec global; 6819ae5db72SJed Brown PetscValidHeaderSpecific(local,VEC_CLASSID,cnt); 6829ae5db72SJed Brown ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 6839ae5db72SJed Brown ierr = VecGetArray(gvec,&array);CHKERRQ(ierr); 6849ae5db72SJed Brown ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr); 6859ae5db72SJed Brown ierr = DMLocalToGlobalBegin(next->dm,local,imode,global);CHKERRQ(ierr); 6869ae5db72SJed Brown ierr = DMLocalToGlobalEnd(next->dm,local,imode,global);CHKERRQ(ierr); 6879ae5db72SJed Brown ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr); 6889ae5db72SJed Brown ierr = VecResetArray(global);CHKERRQ(ierr); 6899ae5db72SJed Brown ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 69047c6ae99SBarry Smith } 69147c6ae99SBarry Smith } 69247c6ae99SBarry Smith va_end(Argp); 69347c6ae99SBarry Smith PetscFunctionReturn(0); 69447c6ae99SBarry Smith } 69547c6ae99SBarry Smith 6966f3c3dcfSJed Brown /*@ 6976f3c3dcfSJed Brown DMCompositeGatherArray - Gathers into a global packed vector from its individual local vectors 6986f3c3dcfSJed Brown 6996f3c3dcfSJed Brown Collective on DMComposite 7006f3c3dcfSJed Brown 7016f3c3dcfSJed Brown Input Parameter: 7026f3c3dcfSJed Brown + dm - the packer object 7036f3c3dcfSJed Brown . gvec - the global vector 704907376e6SBarry Smith . imode - INSERT_VALUES or ADD_VALUES 7056f3c3dcfSJed Brown - lvecs - the individual sequential vectors, NULL for any that are not needed 7066f3c3dcfSJed Brown 7076f3c3dcfSJed Brown Level: advanced 7086f3c3dcfSJed Brown 7096f3c3dcfSJed Brown Notes: 7106f3c3dcfSJed Brown This is a non-variadic alternative to DMCompositeGather(). 7116f3c3dcfSJed Brown 7126f3c3dcfSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 7136f3c3dcfSJed Brown DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 7146f3c3dcfSJed Brown DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries(), 7156f3c3dcfSJed Brown @*/ 7161dac896bSSatish Balay PetscErrorCode DMCompositeGatherArray(DM dm,InsertMode imode,Vec gvec,Vec *lvecs) 7176f3c3dcfSJed Brown { 7186f3c3dcfSJed Brown PetscErrorCode ierr; 7196f3c3dcfSJed Brown struct DMCompositeLink *next; 7206f3c3dcfSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 7216f3c3dcfSJed Brown PetscInt i; 72271b14b3eSStefano Zampini PetscBool flg; 7236f3c3dcfSJed Brown 7246f3c3dcfSJed Brown PetscFunctionBegin; 7256f3c3dcfSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 7266f3c3dcfSJed Brown PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 72771b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 72871b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 7296f3c3dcfSJed Brown if (!com->setup) { 7306f3c3dcfSJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 7316f3c3dcfSJed Brown } 7326f3c3dcfSJed Brown 7336f3c3dcfSJed Brown /* loop over packed objects, handling one at at time */ 7346f3c3dcfSJed Brown for (next=com->next,i=0; next; next=next->next,i++) { 7356f3c3dcfSJed Brown if (lvecs[i]) { 7366f3c3dcfSJed Brown PetscScalar *array; 7376f3c3dcfSJed Brown Vec global; 7386f3c3dcfSJed Brown PetscValidHeaderSpecific(lvecs[i],VEC_CLASSID,3); 7396f3c3dcfSJed Brown ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 7406f3c3dcfSJed Brown ierr = VecGetArray(gvec,&array);CHKERRQ(ierr); 7416f3c3dcfSJed Brown ierr = VecPlaceArray(global,array+next->rstart);CHKERRQ(ierr); 7426f3c3dcfSJed Brown ierr = DMLocalToGlobalBegin(next->dm,lvecs[i],imode,global);CHKERRQ(ierr); 7436f3c3dcfSJed Brown ierr = DMLocalToGlobalEnd(next->dm,lvecs[i],imode,global);CHKERRQ(ierr); 7446f3c3dcfSJed Brown ierr = VecRestoreArray(gvec,&array);CHKERRQ(ierr); 7456f3c3dcfSJed Brown ierr = VecResetArray(global);CHKERRQ(ierr); 7466f3c3dcfSJed Brown ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 7476f3c3dcfSJed Brown } 7486f3c3dcfSJed Brown } 7496f3c3dcfSJed Brown PetscFunctionReturn(0); 7506f3c3dcfSJed Brown } 7516f3c3dcfSJed Brown 752f5f57ec0SBarry Smith /*@ 753aa219208SBarry Smith DMCompositeAddDM - adds a DM vector to a DMComposite 75447c6ae99SBarry Smith 75547c6ae99SBarry Smith Collective on DMComposite 75647c6ae99SBarry Smith 75747c6ae99SBarry Smith Input Parameter: 7589b52a9b5SPatrick Sanan + dmc - the DMComposite (packer) object 759f5f57ec0SBarry Smith - dm - the DM object 76047c6ae99SBarry Smith 76147c6ae99SBarry Smith Level: advanced 76247c6ae99SBarry Smith 7630c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(), 7646eb61c8cSJed Brown DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 76547c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 76647c6ae99SBarry Smith 76747c6ae99SBarry Smith @*/ 7687087cfbeSBarry Smith PetscErrorCode DMCompositeAddDM(DM dmc,DM dm) 76947c6ae99SBarry Smith { 77047c6ae99SBarry Smith PetscErrorCode ierr; 77106ebdd98SJed Brown PetscInt n,nlocal; 77247c6ae99SBarry Smith struct DMCompositeLink *mine,*next; 77306ebdd98SJed Brown Vec global,local; 77447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dmc->data; 77571b14b3eSStefano Zampini PetscBool flg; 77647c6ae99SBarry Smith 77747c6ae99SBarry Smith PetscFunctionBegin; 77847c6ae99SBarry Smith PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 77947c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,2); 78071b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dmc,DMCOMPOSITE,&flg);CHKERRQ(ierr); 78171b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 78247c6ae99SBarry Smith next = com->next; 783ce94432eSBarry Smith if (com->setup) SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite"); 78447c6ae99SBarry Smith 78547c6ae99SBarry Smith /* create new link */ 786b00a9115SJed Brown ierr = PetscNew(&mine);CHKERRQ(ierr); 78747c6ae99SBarry Smith ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 78847c6ae99SBarry Smith ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr); 78947c6ae99SBarry Smith ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr); 79047c6ae99SBarry Smith ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr); 79106ebdd98SJed Brown ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr); 79206ebdd98SJed Brown ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr); 79306ebdd98SJed Brown ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr); 7948865f1eaSKarl Rupp 79547c6ae99SBarry Smith mine->n = n; 79606ebdd98SJed Brown mine->nlocal = nlocal; 79747c6ae99SBarry Smith mine->dm = dm; 7980298fd71SBarry Smith mine->next = NULL; 79947c6ae99SBarry Smith com->n += n; 8007ac2b803SAlex Fikl com->nghost += nlocal; 80147c6ae99SBarry Smith 80247c6ae99SBarry Smith /* add to end of list */ 8038865f1eaSKarl Rupp if (!next) com->next = mine; 8048865f1eaSKarl Rupp else { 80547c6ae99SBarry Smith while (next->next) next = next->next; 80647c6ae99SBarry Smith next->next = mine; 80747c6ae99SBarry Smith } 80847c6ae99SBarry Smith com->nDM++; 80947c6ae99SBarry Smith com->nmine++; 81047c6ae99SBarry Smith PetscFunctionReturn(0); 81147c6ae99SBarry Smith } 81247c6ae99SBarry Smith 8139804daf3SBarry Smith #include <petscdraw.h> 81426887b52SJed Brown PETSC_EXTERN PetscErrorCode VecView_MPI(Vec,PetscViewer); 8157087cfbeSBarry Smith PetscErrorCode VecView_DMComposite(Vec gvec,PetscViewer viewer) 81647c6ae99SBarry Smith { 81747c6ae99SBarry Smith DM dm; 81847c6ae99SBarry Smith PetscErrorCode ierr; 81947c6ae99SBarry Smith struct DMCompositeLink *next; 82047c6ae99SBarry Smith PetscBool isdraw; 821cef07954SSatish Balay DM_Composite *com; 82247c6ae99SBarry Smith 82347c6ae99SBarry Smith PetscFunctionBegin; 824c688c046SMatthew G Knepley ierr = VecGetDM(gvec, &dm);CHKERRQ(ierr); 825ce94432eSBarry Smith if (!dm) SETERRQ(PetscObjectComm((PetscObject)gvec),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite"); 82647c6ae99SBarry Smith com = (DM_Composite*)dm->data; 82747c6ae99SBarry Smith next = com->next; 82847c6ae99SBarry Smith 829251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 83047c6ae99SBarry Smith if (!isdraw) { 83147c6ae99SBarry Smith /* do I really want to call this? */ 83247c6ae99SBarry Smith ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr); 83347c6ae99SBarry Smith } else { 83447c6ae99SBarry Smith PetscInt cnt = 0; 83547c6ae99SBarry Smith 83647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 83747c6ae99SBarry Smith while (next) { 83847c6ae99SBarry Smith Vec vec; 839ca4278abSLisandro Dalcin const PetscScalar *array; 84047c6ae99SBarry Smith PetscInt bs; 84147c6ae99SBarry Smith 8429ae5db72SJed Brown /* Should use VecGetSubVector() eventually, but would need to forward the DM for that to work */ 8439ae5db72SJed Brown ierr = DMGetGlobalVector(next->dm,&vec);CHKERRQ(ierr); 844ca4278abSLisandro Dalcin ierr = VecGetArrayRead(gvec,&array);CHKERRQ(ierr); 845ca4278abSLisandro Dalcin ierr = VecPlaceArray(vec,(PetscScalar*)array+next->rstart);CHKERRQ(ierr); 846ca4278abSLisandro Dalcin ierr = VecRestoreArrayRead(gvec,&array);CHKERRQ(ierr); 84747c6ae99SBarry Smith ierr = VecView(vec,viewer);CHKERRQ(ierr); 8489ae5db72SJed Brown ierr = VecResetArray(vec);CHKERRQ(ierr); 849ca4278abSLisandro Dalcin ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr); 8509ae5db72SJed Brown ierr = DMRestoreGlobalVector(next->dm,&vec);CHKERRQ(ierr); 85147c6ae99SBarry Smith ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr); 85247c6ae99SBarry Smith cnt += bs; 85347c6ae99SBarry Smith next = next->next; 85447c6ae99SBarry Smith } 85547c6ae99SBarry Smith ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr); 85647c6ae99SBarry Smith } 85747c6ae99SBarry Smith PetscFunctionReturn(0); 85847c6ae99SBarry Smith } 85947c6ae99SBarry Smith 8607087cfbeSBarry Smith PetscErrorCode DMCreateGlobalVector_Composite(DM dm,Vec *gvec) 86147c6ae99SBarry Smith { 86247c6ae99SBarry Smith PetscErrorCode ierr; 86347c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 86447c6ae99SBarry Smith 86547c6ae99SBarry Smith PetscFunctionBegin; 86647c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 867d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 868ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),com->n,com->N,gvec);CHKERRQ(ierr); 869c688c046SMatthew G Knepley ierr = VecSetDM(*gvec, dm);CHKERRQ(ierr); 87047c6ae99SBarry Smith ierr = VecSetOperation(*gvec,VECOP_VIEW,(void (*)(void))VecView_DMComposite);CHKERRQ(ierr); 87147c6ae99SBarry Smith PetscFunctionReturn(0); 87247c6ae99SBarry Smith } 87347c6ae99SBarry Smith 8747087cfbeSBarry Smith PetscErrorCode DMCreateLocalVector_Composite(DM dm,Vec *lvec) 87547c6ae99SBarry Smith { 87647c6ae99SBarry Smith PetscErrorCode ierr; 87747c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 87847c6ae99SBarry Smith 87947c6ae99SBarry Smith PetscFunctionBegin; 88047c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 88147c6ae99SBarry Smith if (!com->setup) { 882d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 88347c6ae99SBarry Smith } 884f0e01b1fSVincent Le Chenadec ierr = VecCreateSeq(PETSC_COMM_SELF,com->nghost,lvec);CHKERRQ(ierr); 885c688c046SMatthew G Knepley ierr = VecSetDM(*lvec, dm);CHKERRQ(ierr); 88647c6ae99SBarry Smith PetscFunctionReturn(0); 88747c6ae99SBarry Smith } 88847c6ae99SBarry Smith 88947c6ae99SBarry Smith /*@C 8909ae5db72SJed Brown DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM in the DMComposite, maps to the composite global space 89147c6ae99SBarry Smith 89206ebdd98SJed Brown Collective on DM 89347c6ae99SBarry Smith 89447c6ae99SBarry Smith Input Parameter: 89547c6ae99SBarry Smith . dm - the packer object 89647c6ae99SBarry Smith 89747c6ae99SBarry Smith Output Parameters: 8989ae5db72SJed Brown . ltogs - the individual mappings for each packed vector. Note that this includes 8999ae5db72SJed Brown all the ghost points that individual ghosted DMDA's may have. 90047c6ae99SBarry Smith 90147c6ae99SBarry Smith Level: advanced 90247c6ae99SBarry Smith 90347c6ae99SBarry Smith Notes: 9046eb61c8cSJed Brown Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree(). 90547c6ae99SBarry Smith 906f5f57ec0SBarry Smith Not available from Fortran 907f5f57ec0SBarry Smith 9089ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 90947c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(), 91047c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries() 91147c6ae99SBarry Smith 91247c6ae99SBarry Smith @*/ 9137087cfbeSBarry Smith PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs) 91447c6ae99SBarry Smith { 91547c6ae99SBarry Smith PetscErrorCode ierr; 91647c6ae99SBarry Smith PetscInt i,*idx,n,cnt; 91747c6ae99SBarry Smith struct DMCompositeLink *next; 91847c6ae99SBarry Smith PetscMPIInt rank; 91947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 92071b14b3eSStefano Zampini PetscBool flg; 92147c6ae99SBarry Smith 92247c6ae99SBarry Smith PetscFunctionBegin; 92347c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 92471b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 92571b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 926728e99d6SJed Brown ierr = DMSetUp(dm);CHKERRQ(ierr); 927854ce69bSBarry Smith ierr = PetscMalloc1(com->nDM,ltogs);CHKERRQ(ierr); 92847c6ae99SBarry Smith next = com->next; 929ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 93047c6ae99SBarry Smith 93147c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 93247c6ae99SBarry Smith cnt = 0; 93347c6ae99SBarry Smith while (next) { 9346eb61c8cSJed Brown ISLocalToGlobalMapping ltog; 9356eb61c8cSJed Brown PetscMPIInt size; 93686994e45SJed Brown const PetscInt *suboff,*indices; 9376eb61c8cSJed Brown Vec global; 93847c6ae99SBarry Smith 9396eb61c8cSJed Brown /* Get sub-DM global indices for each local dof */ 9401411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(next->dm,<og);CHKERRQ(ierr); 9416eb61c8cSJed Brown ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr); 94286994e45SJed Brown ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr); 943785e854fSJed Brown ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr); 94447c6ae99SBarry Smith 9456eb61c8cSJed Brown /* Get the offsets for the sub-DM global vector */ 9466eb61c8cSJed Brown ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 9476eb61c8cSJed Brown ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr); 948ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)global),&size);CHKERRQ(ierr); 9496eb61c8cSJed Brown 9506eb61c8cSJed Brown /* Shift the sub-DM definition of the global space to the composite global space */ 9516eb61c8cSJed Brown for (i=0; i<n; i++) { 95286994e45SJed Brown PetscInt subi = indices[i],lo = 0,hi = size,t; 9536eb61c8cSJed Brown /* Binary search to find which rank owns subi */ 9546eb61c8cSJed Brown while (hi-lo > 1) { 9556eb61c8cSJed Brown t = lo + (hi-lo)/2; 9566eb61c8cSJed Brown if (suboff[t] > subi) hi = t; 9576eb61c8cSJed Brown else lo = t; 9586eb61c8cSJed Brown } 9596eb61c8cSJed Brown idx[i] = subi - suboff[lo] + next->grstarts[lo]; 9606eb61c8cSJed Brown } 96186994e45SJed Brown ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr); 962f0413b6fSBarry Smith ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm),1,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr); 9636eb61c8cSJed Brown ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 96447c6ae99SBarry Smith next = next->next; 96547c6ae99SBarry Smith cnt++; 96647c6ae99SBarry Smith } 96747c6ae99SBarry Smith PetscFunctionReturn(0); 96847c6ae99SBarry Smith } 96947c6ae99SBarry Smith 97087c85e80SJed Brown /*@C 9719ae5db72SJed Brown DMCompositeGetLocalISs - Gets index sets for each component of a composite local vector 97287c85e80SJed Brown 97387c85e80SJed Brown Not Collective 97487c85e80SJed Brown 97587c85e80SJed Brown Input Arguments: 97687c85e80SJed Brown . dm - composite DM 97787c85e80SJed Brown 97887c85e80SJed Brown Output Arguments: 97987c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite 98087c85e80SJed Brown 98187c85e80SJed Brown Level: intermediate 98287c85e80SJed Brown 98387c85e80SJed Brown Notes: 98487c85e80SJed Brown At present, a composite local vector does not normally exist. This function is used to provide index sets for 98587c85e80SJed Brown MatGetLocalSubMatrix(). In the future, the scatters for each entry in the DMComposite may be be merged into a single 9869ae5db72SJed Brown scatter to a composite local vector. The user should not typically need to know which is being done. 98787c85e80SJed Brown 98887c85e80SJed Brown To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings(). 98987c85e80SJed Brown 99087c85e80SJed Brown To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs(). 99187c85e80SJed Brown 99287c85e80SJed Brown Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree(). 99387c85e80SJed Brown 994f5f57ec0SBarry Smith Not available from Fortran 995f5f57ec0SBarry Smith 99687c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef() 99787c85e80SJed Brown @*/ 9987087cfbeSBarry Smith PetscErrorCode DMCompositeGetLocalISs(DM dm,IS **is) 99987c85e80SJed Brown { 100087c85e80SJed Brown PetscErrorCode ierr; 100187c85e80SJed Brown DM_Composite *com = (DM_Composite*)dm->data; 100287c85e80SJed Brown struct DMCompositeLink *link; 100387c85e80SJed Brown PetscInt cnt,start; 100471b14b3eSStefano Zampini PetscBool flg; 100587c85e80SJed Brown 100687c85e80SJed Brown PetscFunctionBegin; 100787c85e80SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 100887c85e80SJed Brown PetscValidPointer(is,2); 100971b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 101071b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 1011785e854fSJed Brown ierr = PetscMalloc1(com->nmine,is);CHKERRQ(ierr); 101206ebdd98SJed Brown for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) { 1013520db06cSJed Brown PetscInt bs; 10149ae5db72SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr); 10151411c6eeSJed Brown ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr); 1016520db06cSJed Brown ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr); 1017520db06cSJed Brown } 101887c85e80SJed Brown PetscFunctionReturn(0); 101987c85e80SJed Brown } 102087c85e80SJed Brown 102147c6ae99SBarry Smith /*@C 102247c6ae99SBarry Smith DMCompositeGetGlobalISs - Gets the index sets for each composed object 102347c6ae99SBarry Smith 102447c6ae99SBarry Smith Collective on DMComposite 102547c6ae99SBarry Smith 102647c6ae99SBarry Smith Input Parameter: 102747c6ae99SBarry Smith . dm - the packer object 102847c6ae99SBarry Smith 102947c6ae99SBarry Smith Output Parameters: 103047c6ae99SBarry Smith . is - the array of index sets 103147c6ae99SBarry Smith 103247c6ae99SBarry Smith Level: advanced 103347c6ae99SBarry Smith 103447c6ae99SBarry Smith Notes: 103547c6ae99SBarry Smith The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree() 103647c6ae99SBarry Smith 103747c6ae99SBarry Smith These could be used to extract a subset of vector entries for a "multi-physics" preconditioner 103847c6ae99SBarry Smith 10396eb61c8cSJed Brown Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and 10406eb61c8cSJed Brown DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global 10416eb61c8cSJed Brown indices. 104247c6ae99SBarry Smith 1043f3cb0f7eSJed Brown Fortran Notes: 1044f3cb0f7eSJed Brown 1045f3cb0f7eSJed Brown The output argument 'is' must be an allocated array of sufficient length, which can be learned using DMCompositeGetNumberDM(). 1046f3cb0f7eSJed Brown 10479ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 104847c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(), 104947c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries() 105047c6ae99SBarry Smith 105147c6ae99SBarry Smith @*/ 10527087cfbeSBarry Smith PetscErrorCode DMCompositeGetGlobalISs(DM dm,IS *is[]) 105347c6ae99SBarry Smith { 105447c6ae99SBarry Smith PetscErrorCode ierr; 105566bb578eSMark Adams PetscInt cnt = 0; 105647c6ae99SBarry Smith struct DMCompositeLink *next; 105747c6ae99SBarry Smith PetscMPIInt rank; 105847c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 105971b14b3eSStefano Zampini PetscBool flg; 106047c6ae99SBarry Smith 106147c6ae99SBarry Smith PetscFunctionBegin; 106247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 106371b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 106471b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 106571b14b3eSStefano Zampini if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() before"); 1066854ce69bSBarry Smith ierr = PetscMalloc1(com->nDM,is);CHKERRQ(ierr); 106747c6ae99SBarry Smith next = com->next; 1068ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 106947c6ae99SBarry Smith 107047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 107147c6ae99SBarry Smith while (next) { 1072e5e52638SMatthew G. Knepley PetscDS prob; 1073e5e52638SMatthew G. Knepley 107466bb578eSMark Adams ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),next->n,next->grstart,1,&(*is)[cnt]);CHKERRQ(ierr); 1075e5e52638SMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1076e5e52638SMatthew G. Knepley if (prob) { 107765c226d8SMatthew G. Knepley MatNullSpace space; 107865c226d8SMatthew G. Knepley Mat pmat; 10790f21e855SMatthew G. Knepley PetscObject disc; 10800f21e855SMatthew G. Knepley PetscInt Nf; 108165c226d8SMatthew G. Knepley 1082e5e52638SMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 1083f24dd8d2SMatthew G. Knepley if (cnt < Nf) { 1084e5e52638SMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, cnt, &disc);CHKERRQ(ierr); 10850f21e855SMatthew G. Knepley ierr = PetscObjectQuery(disc, "nullspace", (PetscObject*) &space);CHKERRQ(ierr); 1086aac2dd2dSMatthew G. Knepley if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nullspace", (PetscObject) space);CHKERRQ(ierr);} 10870f21e855SMatthew G. Knepley ierr = PetscObjectQuery(disc, "nearnullspace", (PetscObject*) &space);CHKERRQ(ierr); 1088aac2dd2dSMatthew G. Knepley if (space) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "nearnullspace", (PetscObject) space);CHKERRQ(ierr);} 10890f21e855SMatthew G. Knepley ierr = PetscObjectQuery(disc, "pmat", (PetscObject*) &pmat);CHKERRQ(ierr); 1090aac2dd2dSMatthew G. Knepley if (pmat) {ierr = PetscObjectCompose((PetscObject) (*is)[cnt], "pmat", (PetscObject) pmat);CHKERRQ(ierr);} 109165c226d8SMatthew G. Knepley } 1092f24dd8d2SMatthew G. Knepley } 109347c6ae99SBarry Smith cnt++; 109447c6ae99SBarry Smith next = next->next; 109547c6ae99SBarry Smith } 109647c6ae99SBarry Smith PetscFunctionReturn(0); 109747c6ae99SBarry Smith } 109847c6ae99SBarry Smith 109921c9b008SJed Brown PetscErrorCode DMCreateFieldIS_Composite(DM dm, PetscInt *numFields,char ***fieldNames, IS **fields) 11004d343eeaSMatthew G Knepley { 11014d343eeaSMatthew G Knepley PetscInt nDM; 11024d343eeaSMatthew G Knepley DM *dms; 11034d343eeaSMatthew G Knepley PetscInt i; 11044d343eeaSMatthew G Knepley PetscErrorCode ierr; 11054d343eeaSMatthew G Knepley 11064d343eeaSMatthew G Knepley PetscFunctionBegin; 11074d343eeaSMatthew G Knepley ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr); 11088865f1eaSKarl Rupp if (numFields) *numFields = nDM; 11094d343eeaSMatthew G Knepley ierr = DMCompositeGetGlobalISs(dm, fields);CHKERRQ(ierr); 11104d343eeaSMatthew G Knepley if (fieldNames) { 1111785e854fSJed Brown ierr = PetscMalloc1(nDM, &dms);CHKERRQ(ierr); 1112785e854fSJed Brown ierr = PetscMalloc1(nDM, fieldNames);CHKERRQ(ierr); 11134d343eeaSMatthew G Knepley ierr = DMCompositeGetEntriesArray(dm, dms);CHKERRQ(ierr); 11144d343eeaSMatthew G Knepley for (i=0; i<nDM; i++) { 11154d343eeaSMatthew G Knepley char buf[256]; 11164d343eeaSMatthew G Knepley const char *splitname; 11174d343eeaSMatthew G Knepley 11184d343eeaSMatthew G Knepley /* Split naming precedence: object name, prefix, number */ 11194d343eeaSMatthew G Knepley splitname = ((PetscObject) dm)->name; 11204d343eeaSMatthew G Knepley if (!splitname) { 11214d343eeaSMatthew G Knepley ierr = PetscObjectGetOptionsPrefix((PetscObject)dms[i],&splitname);CHKERRQ(ierr); 11224d343eeaSMatthew G Knepley if (splitname) { 11234d343eeaSMatthew G Knepley size_t len; 11248caf3d72SBarry Smith ierr = PetscStrncpy(buf,splitname,sizeof(buf));CHKERRQ(ierr); 11258caf3d72SBarry Smith buf[sizeof(buf) - 1] = 0; 11264d343eeaSMatthew G Knepley ierr = PetscStrlen(buf,&len);CHKERRQ(ierr); 11274d343eeaSMatthew G Knepley if (buf[len-1] == '_') buf[len-1] = 0; /* Remove trailing underscore if it was used */ 11284d343eeaSMatthew G Knepley splitname = buf; 11294d343eeaSMatthew G Knepley } 11304d343eeaSMatthew G Knepley } 11314d343eeaSMatthew G Knepley if (!splitname) { 11328caf3d72SBarry Smith ierr = PetscSNPrintf(buf,sizeof(buf),"%D",i);CHKERRQ(ierr); 11334d343eeaSMatthew G Knepley splitname = buf; 11344d343eeaSMatthew G Knepley } 113521c9b008SJed Brown ierr = PetscStrallocpy(splitname,&(*fieldNames)[i]);CHKERRQ(ierr); 11364d343eeaSMatthew G Knepley } 11374d343eeaSMatthew G Knepley ierr = PetscFree(dms);CHKERRQ(ierr); 11384d343eeaSMatthew G Knepley } 11394d343eeaSMatthew G Knepley PetscFunctionReturn(0); 11404d343eeaSMatthew G Knepley } 11414d343eeaSMatthew G Knepley 1142e7c4fc90SDmitry Karpeev /* 1143e7c4fc90SDmitry Karpeev This could take over from DMCreateFieldIS(), as it is more general, 11440298fd71SBarry Smith making DMCreateFieldIS() a special case -- calling with dmlist == NULL; 1145e7c4fc90SDmitry Karpeev At this point it's probably best to be less intrusive, however. 1146e7c4fc90SDmitry Karpeev */ 114716621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_Composite(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist) 1148e7c4fc90SDmitry Karpeev { 1149e7c4fc90SDmitry Karpeev PetscInt nDM; 1150e7c4fc90SDmitry Karpeev PetscInt i; 1151e7c4fc90SDmitry Karpeev PetscErrorCode ierr; 1152e7c4fc90SDmitry Karpeev 1153e7c4fc90SDmitry Karpeev PetscFunctionBegin; 1154e7c4fc90SDmitry Karpeev ierr = DMCreateFieldIS_Composite(dm, len, namelist, islist);CHKERRQ(ierr); 1155e7c4fc90SDmitry Karpeev if (dmlist) { 1156e7c4fc90SDmitry Karpeev ierr = DMCompositeGetNumberDM(dm, &nDM);CHKERRQ(ierr); 1157785e854fSJed Brown ierr = PetscMalloc1(nDM, dmlist);CHKERRQ(ierr); 1158e7c4fc90SDmitry Karpeev ierr = DMCompositeGetEntriesArray(dm, *dmlist);CHKERRQ(ierr); 1159e7c4fc90SDmitry Karpeev for (i=0; i<nDM; i++) { 1160e7c4fc90SDmitry Karpeev ierr = PetscObjectReference((PetscObject)((*dmlist)[i]));CHKERRQ(ierr); 1161e7c4fc90SDmitry Karpeev } 1162e7c4fc90SDmitry Karpeev } 1163e7c4fc90SDmitry Karpeev PetscFunctionReturn(0); 1164e7c4fc90SDmitry Karpeev } 1165e7c4fc90SDmitry Karpeev 1166e7c4fc90SDmitry Karpeev 1167e7c4fc90SDmitry Karpeev 116847c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 116947c6ae99SBarry Smith /*@C 11709ae5db72SJed Brown DMCompositeGetLocalVectors - Gets local vectors for each part of a DMComposite. 117147c6ae99SBarry Smith Use DMCompositeRestoreLocalVectors() to return them. 117247c6ae99SBarry Smith 117347c6ae99SBarry Smith Not Collective 117447c6ae99SBarry Smith 117547c6ae99SBarry Smith Input Parameter: 117647c6ae99SBarry Smith . dm - the packer object 117747c6ae99SBarry Smith 117847c6ae99SBarry Smith Output Parameter: 11799ae5db72SJed Brown . Vec ... - the individual sequential Vecs 118047c6ae99SBarry Smith 118147c6ae99SBarry Smith Level: advanced 118247c6ae99SBarry Smith 1183f5f57ec0SBarry Smith Not available from Fortran 1184f5f57ec0SBarry Smith 11859ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 11866eb61c8cSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 118747c6ae99SBarry Smith DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries() 118847c6ae99SBarry Smith 118947c6ae99SBarry Smith @*/ 11907087cfbeSBarry Smith PetscErrorCode DMCompositeGetLocalVectors(DM dm,...) 119147c6ae99SBarry Smith { 119247c6ae99SBarry Smith va_list Argp; 119347c6ae99SBarry Smith PetscErrorCode ierr; 119447c6ae99SBarry Smith struct DMCompositeLink *next; 119547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 119671b14b3eSStefano Zampini PetscBool flg; 119747c6ae99SBarry Smith 119847c6ae99SBarry Smith PetscFunctionBegin; 119947c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 120071b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 120171b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 120247c6ae99SBarry Smith next = com->next; 120347c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 120447c6ae99SBarry Smith va_start(Argp,dm); 120547c6ae99SBarry Smith while (next) { 120647c6ae99SBarry Smith Vec *vec; 120747c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 120806930112SJed Brown if (vec) {ierr = DMGetLocalVector(next->dm,vec);CHKERRQ(ierr);} 120947c6ae99SBarry Smith next = next->next; 121047c6ae99SBarry Smith } 121147c6ae99SBarry Smith va_end(Argp); 121247c6ae99SBarry Smith PetscFunctionReturn(0); 121347c6ae99SBarry Smith } 121447c6ae99SBarry Smith 121547c6ae99SBarry Smith /*@C 12169ae5db72SJed Brown DMCompositeRestoreLocalVectors - Restores local vectors for each part of a DMComposite. 121747c6ae99SBarry Smith 121847c6ae99SBarry Smith Not Collective 121947c6ae99SBarry Smith 122047c6ae99SBarry Smith Input Parameter: 122147c6ae99SBarry Smith . dm - the packer object 122247c6ae99SBarry Smith 122347c6ae99SBarry Smith Output Parameter: 12249ae5db72SJed Brown . Vec ... - the individual sequential Vecs 122547c6ae99SBarry Smith 122647c6ae99SBarry Smith Level: advanced 122747c6ae99SBarry Smith 1228f5f57ec0SBarry Smith Not available from Fortran 1229f5f57ec0SBarry Smith 12309ae5db72SJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), 12316eb61c8cSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 123247c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries() 123347c6ae99SBarry Smith 123447c6ae99SBarry Smith @*/ 12357087cfbeSBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors(DM dm,...) 123647c6ae99SBarry Smith { 123747c6ae99SBarry Smith va_list Argp; 123847c6ae99SBarry Smith PetscErrorCode ierr; 123947c6ae99SBarry Smith struct DMCompositeLink *next; 124047c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 124171b14b3eSStefano Zampini PetscBool flg; 124247c6ae99SBarry Smith 124347c6ae99SBarry Smith PetscFunctionBegin; 124447c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 124571b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 124671b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 124747c6ae99SBarry Smith next = com->next; 124847c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 124947c6ae99SBarry Smith va_start(Argp,dm); 125047c6ae99SBarry Smith while (next) { 125147c6ae99SBarry Smith Vec *vec; 125247c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 125306930112SJed Brown if (vec) {ierr = DMRestoreLocalVector(next->dm,vec);CHKERRQ(ierr);} 125447c6ae99SBarry Smith next = next->next; 125547c6ae99SBarry Smith } 125647c6ae99SBarry Smith va_end(Argp); 125747c6ae99SBarry Smith PetscFunctionReturn(0); 125847c6ae99SBarry Smith } 125947c6ae99SBarry Smith 126047c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 126147c6ae99SBarry Smith /*@C 12629ae5db72SJed Brown DMCompositeGetEntries - Gets the DM for each entry in a DMComposite. 126347c6ae99SBarry Smith 126447c6ae99SBarry Smith Not Collective 126547c6ae99SBarry Smith 126647c6ae99SBarry Smith Input Parameter: 126747c6ae99SBarry Smith . dm - the packer object 126847c6ae99SBarry Smith 126947c6ae99SBarry Smith Output Parameter: 12709ae5db72SJed Brown . DM ... - the individual entries (DMs) 127147c6ae99SBarry Smith 127247c6ae99SBarry Smith Level: advanced 127347c6ae99SBarry Smith 127495452b02SPatrick Sanan Fortran Notes: 127595452b02SPatrick Sanan Available as DMCompositeGetEntries() for one output DM, DMCompositeGetEntries2() for 2, etc 1276f5f57ec0SBarry Smith 12772fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntriesArray() 12786eb61c8cSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 127947c6ae99SBarry Smith DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(), DMCompositeScatter(), 128047c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors() 128147c6ae99SBarry Smith 128247c6ae99SBarry Smith @*/ 12837087cfbeSBarry Smith PetscErrorCode DMCompositeGetEntries(DM dm,...) 128447c6ae99SBarry Smith { 128547c6ae99SBarry Smith va_list Argp; 128647c6ae99SBarry Smith struct DMCompositeLink *next; 128747c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 128871b14b3eSStefano Zampini PetscBool flg; 128971b14b3eSStefano Zampini PetscErrorCode ierr; 129047c6ae99SBarry Smith 129147c6ae99SBarry Smith PetscFunctionBegin; 129247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 129371b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 129471b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 129547c6ae99SBarry Smith next = com->next; 129647c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 129747c6ae99SBarry Smith va_start(Argp,dm); 129847c6ae99SBarry Smith while (next) { 129947c6ae99SBarry Smith DM *dmn; 130047c6ae99SBarry Smith dmn = va_arg(Argp, DM*); 13019ae5db72SJed Brown if (dmn) *dmn = next->dm; 130247c6ae99SBarry Smith next = next->next; 130347c6ae99SBarry Smith } 130447c6ae99SBarry Smith va_end(Argp); 130547c6ae99SBarry Smith PetscFunctionReturn(0); 130647c6ae99SBarry Smith } 130747c6ae99SBarry Smith 1308dbab29e1SMark F. Adams /*@C 13092fa5ba8aSJed Brown DMCompositeGetEntriesArray - Gets the DM for each entry in a DMComposite. 13102fa5ba8aSJed Brown 13112fa5ba8aSJed Brown Not Collective 13122fa5ba8aSJed Brown 13132fa5ba8aSJed Brown Input Parameter: 1314907376e6SBarry Smith . dm - the packer object 1315907376e6SBarry Smith 1316907376e6SBarry Smith Output Parameter: 1317907376e6SBarry Smith . dms - array of sufficient length (see DMCompositeGetNumberDM()) to hold the individual DMs 13182fa5ba8aSJed Brown 13192fa5ba8aSJed Brown Level: advanced 13202fa5ba8aSJed Brown 13212fa5ba8aSJed Brown .seealso DMDestroy(), DMCompositeAddDM(), DMCreateGlobalVector(), DMCompositeGetEntries() 13222fa5ba8aSJed Brown DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(), 13232fa5ba8aSJed Brown DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(), DMCompositeScatter(), 13242fa5ba8aSJed Brown DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors() 13252fa5ba8aSJed Brown 13262fa5ba8aSJed Brown @*/ 13272fa5ba8aSJed Brown PetscErrorCode DMCompositeGetEntriesArray(DM dm,DM dms[]) 13282fa5ba8aSJed Brown { 13292fa5ba8aSJed Brown struct DMCompositeLink *next; 13302fa5ba8aSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 13312fa5ba8aSJed Brown PetscInt i; 133271b14b3eSStefano Zampini PetscBool flg; 133371b14b3eSStefano Zampini PetscErrorCode ierr; 13342fa5ba8aSJed Brown 13352fa5ba8aSJed Brown PetscFunctionBegin; 13362fa5ba8aSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 133771b14b3eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMCOMPOSITE,&flg);CHKERRQ(ierr); 133871b14b3eSStefano Zampini if (!flg) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for type %s",((PetscObject)dm)->type_name); 13392fa5ba8aSJed Brown /* loop over packed objects, handling one at at time */ 13402fa5ba8aSJed Brown for (next=com->next,i=0; next; next=next->next,i++) dms[i] = next->dm; 13412fa5ba8aSJed Brown PetscFunctionReturn(0); 13422fa5ba8aSJed Brown } 13432fa5ba8aSJed Brown 1344e10fd815SStefano Zampini typedef struct { 1345e10fd815SStefano Zampini DM dm; 1346e10fd815SStefano Zampini PetscViewer *subv; 1347e10fd815SStefano Zampini Vec *vecs; 1348e10fd815SStefano Zampini } GLVisViewerCtx; 1349e10fd815SStefano Zampini 1350e10fd815SStefano Zampini static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 1351e10fd815SStefano Zampini { 1352e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 1353e10fd815SStefano Zampini PetscInt i,n; 1354e10fd815SStefano Zampini PetscErrorCode ierr; 1355e10fd815SStefano Zampini 1356e10fd815SStefano Zampini PetscFunctionBegin; 1357e10fd815SStefano Zampini ierr = DMCompositeGetNumberDM(ctx->dm,&n);CHKERRQ(ierr); 1358e10fd815SStefano Zampini for (i = 0; i < n; i++) { 1359e10fd815SStefano Zampini ierr = PetscViewerDestroy(&ctx->subv[i]);CHKERRQ(ierr); 1360e10fd815SStefano Zampini } 1361e10fd815SStefano Zampini ierr = PetscFree2(ctx->subv,ctx->vecs);CHKERRQ(ierr); 1362e10fd815SStefano Zampini ierr = DMDestroy(&ctx->dm);CHKERRQ(ierr); 1363e10fd815SStefano Zampini ierr = PetscFree(ctx);CHKERRQ(ierr); 1364e10fd815SStefano Zampini PetscFunctionReturn(0); 1365e10fd815SStefano Zampini } 1366e10fd815SStefano Zampini 1367e10fd815SStefano Zampini static PetscErrorCode DMCompositeSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 1368e10fd815SStefano Zampini { 1369e10fd815SStefano Zampini Vec X = (Vec)oX; 1370e10fd815SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 1371e10fd815SStefano Zampini PetscInt i,n,cumf; 1372e10fd815SStefano Zampini PetscErrorCode ierr; 1373e10fd815SStefano Zampini 1374e10fd815SStefano Zampini PetscFunctionBegin; 1375e10fd815SStefano Zampini ierr = DMCompositeGetNumberDM(ctx->dm,&n);CHKERRQ(ierr); 1376e10fd815SStefano Zampini ierr = DMCompositeGetAccessArray(ctx->dm,X,n,NULL,ctx->vecs);CHKERRQ(ierr); 1377e10fd815SStefano Zampini for (i = 0, cumf = 0; i < n; i++) { 1378e10fd815SStefano Zampini PetscErrorCode (*g2l)(PetscObject,PetscInt,PetscObject[],void*); 1379e10fd815SStefano Zampini void *fctx; 1380e10fd815SStefano Zampini PetscInt nfi; 1381e10fd815SStefano Zampini 1382e10fd815SStefano Zampini ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nfi,NULL,NULL,&g2l,NULL,&fctx);CHKERRQ(ierr); 1383e10fd815SStefano Zampini if (!nfi) continue; 1384e10fd815SStefano Zampini if (g2l) { 1385e10fd815SStefano Zampini ierr = (*g2l)((PetscObject)ctx->vecs[i],nfi,oXfield+cumf,fctx);CHKERRQ(ierr); 1386e10fd815SStefano Zampini } else { 1387e10fd815SStefano Zampini ierr = VecCopy(ctx->vecs[i],(Vec)(oXfield[cumf]));CHKERRQ(ierr); 1388e10fd815SStefano Zampini } 1389e10fd815SStefano Zampini cumf += nfi; 1390e10fd815SStefano Zampini } 1391e10fd815SStefano Zampini ierr = DMCompositeRestoreAccessArray(ctx->dm,X,n,NULL,ctx->vecs);CHKERRQ(ierr); 1392e10fd815SStefano Zampini PetscFunctionReturn(0); 1393e10fd815SStefano Zampini } 1394e10fd815SStefano Zampini 1395e10fd815SStefano Zampini static PetscErrorCode DMSetUpGLVisViewer_Composite(PetscObject odm, PetscViewer viewer) 1396e10fd815SStefano Zampini { 1397e10fd815SStefano Zampini DM dm = (DM)odm, *dms; 1398e10fd815SStefano Zampini Vec *Ufds; 1399e10fd815SStefano Zampini GLVisViewerCtx *ctx; 1400e10fd815SStefano Zampini PetscInt i,n,tnf,*sdim; 1401e10fd815SStefano Zampini char **fecs; 1402e10fd815SStefano Zampini PetscErrorCode ierr; 1403e10fd815SStefano Zampini 1404e10fd815SStefano Zampini PetscFunctionBegin; 1405e10fd815SStefano Zampini ierr = PetscNew(&ctx);CHKERRQ(ierr); 1406e10fd815SStefano Zampini ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 1407e10fd815SStefano Zampini ctx->dm = dm; 1408e10fd815SStefano Zampini ierr = DMCompositeGetNumberDM(dm,&n);CHKERRQ(ierr); 1409e10fd815SStefano Zampini ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 1410e10fd815SStefano Zampini ierr = DMCompositeGetEntriesArray(dm,dms);CHKERRQ(ierr); 1411e10fd815SStefano Zampini ierr = PetscMalloc2(n,&ctx->subv,n,&ctx->vecs);CHKERRQ(ierr); 1412e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1413e10fd815SStefano Zampini PetscInt nf; 1414e10fd815SStefano Zampini 1415e10fd815SStefano Zampini ierr = PetscViewerCreate(PetscObjectComm(odm),&ctx->subv[i]);CHKERRQ(ierr); 1416e10fd815SStefano Zampini ierr = PetscViewerSetType(ctx->subv[i],PETSCVIEWERGLVIS);CHKERRQ(ierr); 1417e10fd815SStefano Zampini ierr = PetscViewerGLVisSetDM_Private(ctx->subv[i],(PetscObject)dms[i]);CHKERRQ(ierr); 1418e10fd815SStefano Zampini ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1419e10fd815SStefano Zampini tnf += nf; 1420e10fd815SStefano Zampini } 1421e10fd815SStefano Zampini ierr = PetscFree(dms);CHKERRQ(ierr); 1422e10fd815SStefano Zampini ierr = PetscMalloc3(tnf,&fecs,tnf,&sdim,tnf,&Ufds);CHKERRQ(ierr); 1423e10fd815SStefano Zampini for (i = 0, tnf = 0; i < n; i++) { 1424e10fd815SStefano Zampini PetscInt *sd,nf,f; 1425e10fd815SStefano Zampini const char **fec; 1426e10fd815SStefano Zampini Vec *Uf; 1427e10fd815SStefano Zampini 1428e10fd815SStefano Zampini ierr = PetscViewerGLVisGetFields_Private(ctx->subv[i],&nf,&fec,&sd,NULL,(PetscObject**)&Uf,NULL);CHKERRQ(ierr); 1429e10fd815SStefano Zampini for (f = 0; f < nf; f++) { 1430e10fd815SStefano Zampini ierr = PetscStrallocpy(fec[f],&fecs[tnf+f]);CHKERRQ(ierr); 1431e10fd815SStefano Zampini Ufds[tnf+f] = Uf[f]; 1432e10fd815SStefano Zampini sdim[tnf+f] = sd[f]; 1433e10fd815SStefano Zampini } 1434e10fd815SStefano Zampini tnf += nf; 1435e10fd815SStefano Zampini } 1436e10fd815SStefano Zampini ierr = PetscViewerGLVisSetFields(viewer,tnf,(const char**)fecs,sdim,DMCompositeSampleGLVisFields_Private,(PetscObject*)Ufds,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr); 1437e10fd815SStefano Zampini for (i = 0; i < tnf; i++) { 1438e10fd815SStefano Zampini ierr = PetscFree(fecs[i]);CHKERRQ(ierr); 1439e10fd815SStefano Zampini } 1440e10fd815SStefano Zampini ierr = PetscFree3(fecs,sdim,Ufds);CHKERRQ(ierr); 1441e10fd815SStefano Zampini PetscFunctionReturn(0); 1442e10fd815SStefano Zampini } 1443e10fd815SStefano Zampini 14447087cfbeSBarry Smith PetscErrorCode DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine) 144547c6ae99SBarry Smith { 144647c6ae99SBarry Smith PetscErrorCode ierr; 144747c6ae99SBarry Smith struct DMCompositeLink *next; 144847c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dmi->data; 144947c6ae99SBarry Smith DM dm; 145047c6ae99SBarry Smith 145147c6ae99SBarry Smith PetscFunctionBegin; 145247c6ae99SBarry Smith PetscValidHeaderSpecific(dmi,DM_CLASSID,1); 1453ce94432eSBarry Smith if (comm == MPI_COMM_NULL) { 1454ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr); 1455ce94432eSBarry Smith } 14562ce3a92bSJed Brown ierr = DMSetUp(dmi);CHKERRQ(ierr); 145747c6ae99SBarry Smith next = com->next; 145847c6ae99SBarry Smith ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr); 145947c6ae99SBarry Smith 146047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 146147c6ae99SBarry Smith while (next) { 146247c6ae99SBarry Smith ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr); 146347c6ae99SBarry Smith ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr); 146447c6ae99SBarry Smith ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 146547c6ae99SBarry Smith next = next->next; 146647c6ae99SBarry Smith } 146747c6ae99SBarry Smith PetscFunctionReturn(0); 146847c6ae99SBarry Smith } 146947c6ae99SBarry Smith 147014354c39SJed Brown PetscErrorCode DMCoarsen_Composite(DM dmi,MPI_Comm comm,DM *fine) 147114354c39SJed Brown { 147214354c39SJed Brown PetscErrorCode ierr; 147314354c39SJed Brown struct DMCompositeLink *next; 147414354c39SJed Brown DM_Composite *com = (DM_Composite*)dmi->data; 147514354c39SJed Brown DM dm; 147614354c39SJed Brown 147714354c39SJed Brown PetscFunctionBegin; 147814354c39SJed Brown PetscValidHeaderSpecific(dmi,DM_CLASSID,1); 14792ce3a92bSJed Brown ierr = DMSetUp(dmi);CHKERRQ(ierr); 14802ee06e3bSJed Brown if (comm == MPI_COMM_NULL) { 148125296bd5SBarry Smith ierr = PetscObjectGetComm((PetscObject)dmi,&comm);CHKERRQ(ierr); 148225296bd5SBarry Smith } 148314354c39SJed Brown next = com->next; 148414354c39SJed Brown ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr); 148514354c39SJed Brown 148614354c39SJed Brown /* loop over packed objects, handling one at at time */ 148714354c39SJed Brown while (next) { 148814354c39SJed Brown ierr = DMCoarsen(next->dm,comm,&dm);CHKERRQ(ierr); 148914354c39SJed Brown ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr); 149014354c39SJed Brown ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 149114354c39SJed Brown next = next->next; 149214354c39SJed Brown } 149314354c39SJed Brown PetscFunctionReturn(0); 149414354c39SJed Brown } 149547c6ae99SBarry Smith 1496e727c939SJed Brown PetscErrorCode DMCreateInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v) 149747c6ae99SBarry Smith { 149847c6ae99SBarry Smith PetscErrorCode ierr; 14999ae5db72SJed Brown PetscInt m,n,M,N,nDM,i; 150047c6ae99SBarry Smith struct DMCompositeLink *nextc; 150147c6ae99SBarry Smith struct DMCompositeLink *nextf; 150225296bd5SBarry Smith Vec gcoarse,gfine,*vecs; 150347c6ae99SBarry Smith DM_Composite *comcoarse = (DM_Composite*)coarse->data; 150447c6ae99SBarry Smith DM_Composite *comfine = (DM_Composite*)fine->data; 15059ae5db72SJed Brown Mat *mats; 150647c6ae99SBarry Smith 150747c6ae99SBarry Smith PetscFunctionBegin; 150847c6ae99SBarry Smith PetscValidHeaderSpecific(coarse,DM_CLASSID,1); 150947c6ae99SBarry Smith PetscValidHeaderSpecific(fine,DM_CLASSID,2); 1510f692024eSJed Brown ierr = DMSetUp(coarse);CHKERRQ(ierr); 1511f692024eSJed Brown ierr = DMSetUp(fine);CHKERRQ(ierr); 151247c6ae99SBarry Smith /* use global vectors only for determining matrix layout */ 15139ae5db72SJed Brown ierr = DMGetGlobalVector(coarse,&gcoarse);CHKERRQ(ierr); 15149ae5db72SJed Brown ierr = DMGetGlobalVector(fine,&gfine);CHKERRQ(ierr); 151547c6ae99SBarry Smith ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr); 151647c6ae99SBarry Smith ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr); 151747c6ae99SBarry Smith ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr); 151847c6ae99SBarry Smith ierr = VecGetSize(gfine,&M);CHKERRQ(ierr); 15199ae5db72SJed Brown ierr = DMRestoreGlobalVector(coarse,&gcoarse);CHKERRQ(ierr); 15209ae5db72SJed Brown ierr = DMRestoreGlobalVector(fine,&gfine);CHKERRQ(ierr); 152147c6ae99SBarry Smith 15229ae5db72SJed Brown nDM = comfine->nDM; 1523ce94432eSBarry 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); 15241795a4d1SJed Brown ierr = PetscCalloc1(nDM*nDM,&mats);CHKERRQ(ierr); 152525296bd5SBarry Smith if (v) { 15261795a4d1SJed Brown ierr = PetscCalloc1(nDM,&vecs);CHKERRQ(ierr); 152725296bd5SBarry Smith } 152847c6ae99SBarry Smith 152947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 15309ae5db72SJed Brown for (nextc=comcoarse->next,nextf=comfine->next,i=0; nextc; nextc=nextc->next,nextf=nextf->next,i++) { 153125296bd5SBarry Smith if (!v) { 15320298fd71SBarry Smith ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],NULL);CHKERRQ(ierr); 153325296bd5SBarry Smith } else { 153425296bd5SBarry Smith ierr = DMCreateInterpolation(nextc->dm,nextf->dm,&mats[i*nDM+i],&vecs[i]);CHKERRQ(ierr); 153525296bd5SBarry Smith } 153647c6ae99SBarry Smith } 1537ce94432eSBarry Smith ierr = MatCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,nDM,NULL,mats,A);CHKERRQ(ierr); 153825296bd5SBarry Smith if (v) { 1539ce94432eSBarry Smith ierr = VecCreateNest(PetscObjectComm((PetscObject)fine),nDM,NULL,vecs,v);CHKERRQ(ierr); 154025296bd5SBarry Smith } 15419ae5db72SJed Brown for (i=0; i<nDM*nDM; i++) {ierr = MatDestroy(&mats[i]);CHKERRQ(ierr);} 15429ae5db72SJed Brown ierr = PetscFree(mats);CHKERRQ(ierr); 154325296bd5SBarry Smith if (v) { 154425296bd5SBarry Smith for (i=0; i<nDM; i++) {ierr = VecDestroy(&vecs[i]);CHKERRQ(ierr);} 154525296bd5SBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 154625296bd5SBarry Smith } 154747c6ae99SBarry Smith PetscFunctionReturn(0); 154847c6ae99SBarry Smith } 154947c6ae99SBarry Smith 1550184d77edSJed Brown static PetscErrorCode DMGetLocalToGlobalMapping_Composite(DM dm) 15511411c6eeSJed Brown { 15521411c6eeSJed Brown DM_Composite *com = (DM_Composite*)dm->data; 15531411c6eeSJed Brown ISLocalToGlobalMapping *ltogs; 1554f7efa3c7SJed Brown PetscInt i; 15551411c6eeSJed Brown PetscErrorCode ierr; 15561411c6eeSJed Brown 15571411c6eeSJed Brown PetscFunctionBegin; 15581411c6eeSJed Brown /* Set the ISLocalToGlobalMapping on the new matrix */ 15591411c6eeSJed Brown ierr = DMCompositeGetISLocalToGlobalMappings(dm,<ogs);CHKERRQ(ierr); 1560ce94432eSBarry Smith ierr = ISLocalToGlobalMappingConcatenate(PetscObjectComm((PetscObject)dm),com->nDM,ltogs,&dm->ltogmap);CHKERRQ(ierr); 15619ae5db72SJed Brown for (i=0; i<com->nDM; i++) {ierr = ISLocalToGlobalMappingDestroy(<ogs[i]);CHKERRQ(ierr);} 15621411c6eeSJed Brown ierr = PetscFree(ltogs);CHKERRQ(ierr); 15631411c6eeSJed Brown PetscFunctionReturn(0); 15641411c6eeSJed Brown } 15651411c6eeSJed Brown 15661411c6eeSJed Brown 1567b412c318SBarry Smith PetscErrorCode DMCreateColoring_Composite(DM dm,ISColoringType ctype,ISColoring *coloring) 156847c6ae99SBarry Smith { 156947c6ae99SBarry Smith PetscErrorCode ierr; 157047c6ae99SBarry Smith PetscInt n,i,cnt; 157147c6ae99SBarry Smith ISColoringValue *colors; 157247c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 157347c6ae99SBarry Smith ISColoringValue maxcol = 0; 157447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 157547c6ae99SBarry Smith 157647c6ae99SBarry Smith PetscFunctionBegin; 157747c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 15785bdb020cSBarry Smith if (ctype == IS_COLORING_LOCAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only global coloring supported"); 1579e3247f34SBarry Smith else if (ctype == IS_COLORING_GLOBAL) { 158047c6ae99SBarry Smith n = com->n; 1581ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType"); 1582785e854fSJed Brown ierr = PetscMalloc1(n,&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */ 158347c6ae99SBarry Smith 1584c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL);CHKERRQ(ierr); 158547c6ae99SBarry Smith if (dense) { 158647c6ae99SBarry Smith for (i=0; i<n; i++) { 158747c6ae99SBarry Smith colors[i] = (ISColoringValue)(com->rstart + i); 158847c6ae99SBarry Smith } 158947c6ae99SBarry Smith maxcol = com->N; 159047c6ae99SBarry Smith } else { 159147c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 159247c6ae99SBarry Smith PetscMPIInt rank; 159347c6ae99SBarry Smith 1594ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 159547c6ae99SBarry Smith cnt = 0; 159647c6ae99SBarry Smith while (next) { 159747c6ae99SBarry Smith ISColoring lcoloring; 159847c6ae99SBarry Smith 1599b412c318SBarry Smith ierr = DMCreateColoring(next->dm,IS_COLORING_GLOBAL,&lcoloring);CHKERRQ(ierr); 160047c6ae99SBarry Smith for (i=0; i<lcoloring->N; i++) { 160147c6ae99SBarry Smith colors[cnt++] = maxcol + lcoloring->colors[i]; 160247c6ae99SBarry Smith } 160347c6ae99SBarry Smith maxcol += lcoloring->n; 1604fcfd50ebSBarry Smith ierr = ISColoringDestroy(&lcoloring);CHKERRQ(ierr); 160547c6ae99SBarry Smith next = next->next; 160647c6ae99SBarry Smith } 160747c6ae99SBarry Smith } 1608aaf3ff59SMatthew G. Knepley ierr = ISColoringCreate(PetscObjectComm((PetscObject)dm),maxcol,n,colors,PETSC_OWN_POINTER,coloring);CHKERRQ(ierr); 160947c6ae99SBarry Smith PetscFunctionReturn(0); 161047c6ae99SBarry Smith } 161147c6ae99SBarry Smith 16127087cfbeSBarry Smith PetscErrorCode DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec) 161347c6ae99SBarry Smith { 161447c6ae99SBarry Smith PetscErrorCode ierr; 161547c6ae99SBarry Smith struct DMCompositeLink *next; 161647c6ae99SBarry Smith PetscScalar *garray,*larray; 161747c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 161847c6ae99SBarry Smith 161947c6ae99SBarry Smith PetscFunctionBegin; 162047c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 162147c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 162239d35262SVincent Le Chenadec 162347c6ae99SBarry Smith if (!com->setup) { 1624d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 162547c6ae99SBarry Smith } 162639d35262SVincent Le Chenadec 162747c6ae99SBarry Smith ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr); 162847c6ae99SBarry Smith ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr); 162947c6ae99SBarry Smith 163047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 163139d35262SVincent Le Chenadec next = com->next; 163247c6ae99SBarry Smith while (next) { 163347c6ae99SBarry Smith Vec local,global; 163447c6ae99SBarry Smith PetscInt N; 163547c6ae99SBarry Smith 163647c6ae99SBarry Smith ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 163747c6ae99SBarry Smith ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr); 163847c6ae99SBarry Smith ierr = VecPlaceArray(global,garray);CHKERRQ(ierr); 163947c6ae99SBarry Smith ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr); 164047c6ae99SBarry Smith ierr = VecPlaceArray(local,larray);CHKERRQ(ierr); 164147c6ae99SBarry Smith ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr); 164247c6ae99SBarry Smith ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr); 164347c6ae99SBarry Smith ierr = VecResetArray(global);CHKERRQ(ierr); 164447c6ae99SBarry Smith ierr = VecResetArray(local);CHKERRQ(ierr); 164547c6ae99SBarry Smith ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 164639d35262SVincent Le Chenadec ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr); 164739d35262SVincent Le Chenadec 164806ebdd98SJed Brown larray += next->nlocal; 164939d35262SVincent Le Chenadec garray += next->n; 165047c6ae99SBarry Smith next = next->next; 165147c6ae99SBarry Smith } 165247c6ae99SBarry Smith 16530298fd71SBarry Smith ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr); 16540298fd71SBarry Smith ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr); 165547c6ae99SBarry Smith PetscFunctionReturn(0); 165647c6ae99SBarry Smith } 165747c6ae99SBarry Smith 16587087cfbeSBarry Smith PetscErrorCode DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec) 16590c010503SBarry Smith { 16600c010503SBarry Smith PetscFunctionBegin; 166139d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 166239d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 166339d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,4); 166439d35262SVincent Le Chenadec PetscFunctionReturn(0); 166539d35262SVincent Le Chenadec } 166639d35262SVincent Le Chenadec 166739d35262SVincent Le Chenadec PetscErrorCode DMLocalToGlobalBegin_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec) 166839d35262SVincent Le Chenadec { 166939d35262SVincent Le Chenadec PetscErrorCode ierr; 167039d35262SVincent Le Chenadec struct DMCompositeLink *next; 167139d35262SVincent Le Chenadec PetscScalar *larray,*garray; 167239d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite*)dm->data; 167339d35262SVincent Le Chenadec 167439d35262SVincent Le Chenadec PetscFunctionBegin; 167539d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 167639d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,2); 167739d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,4); 167839d35262SVincent Le Chenadec 167939d35262SVincent Le Chenadec if (!com->setup) { 168039d35262SVincent Le Chenadec ierr = DMSetUp(dm);CHKERRQ(ierr); 168139d35262SVincent Le Chenadec } 168239d35262SVincent Le Chenadec 168339d35262SVincent Le Chenadec ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr); 168439d35262SVincent Le Chenadec ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr); 168539d35262SVincent Le Chenadec 168639d35262SVincent Le Chenadec /* loop over packed objects, handling one at at time */ 168739d35262SVincent Le Chenadec next = com->next; 168839d35262SVincent Le Chenadec while (next) { 168939d35262SVincent Le Chenadec Vec global,local; 169039d35262SVincent Le Chenadec 169139d35262SVincent Le Chenadec ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr); 169239d35262SVincent Le Chenadec ierr = VecPlaceArray(local,larray);CHKERRQ(ierr); 169339d35262SVincent Le Chenadec ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 169439d35262SVincent Le Chenadec ierr = VecPlaceArray(global,garray);CHKERRQ(ierr); 169539d35262SVincent Le Chenadec ierr = DMLocalToGlobalBegin(next->dm,local,mode,global);CHKERRQ(ierr); 169639d35262SVincent Le Chenadec ierr = DMLocalToGlobalEnd(next->dm,local,mode,global);CHKERRQ(ierr); 169739d35262SVincent Le Chenadec ierr = VecResetArray(local);CHKERRQ(ierr); 169839d35262SVincent Le Chenadec ierr = VecResetArray(global);CHKERRQ(ierr); 169939d35262SVincent Le Chenadec ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 170039d35262SVincent Le Chenadec ierr = DMRestoreLocalVector(next->dm,&local);CHKERRQ(ierr); 170139d35262SVincent Le Chenadec 170239d35262SVincent Le Chenadec garray += next->n; 170339d35262SVincent Le Chenadec larray += next->nlocal; 170439d35262SVincent Le Chenadec next = next->next; 170539d35262SVincent Le Chenadec } 170639d35262SVincent Le Chenadec 170739d35262SVincent Le Chenadec ierr = VecRestoreArray(gvec,NULL);CHKERRQ(ierr); 170839d35262SVincent Le Chenadec ierr = VecRestoreArray(lvec,NULL);CHKERRQ(ierr); 170939d35262SVincent Le Chenadec 171039d35262SVincent Le Chenadec PetscFunctionReturn(0); 171139d35262SVincent Le Chenadec } 171239d35262SVincent Le Chenadec 171339d35262SVincent Le Chenadec PetscErrorCode DMLocalToGlobalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec) 171439d35262SVincent Le Chenadec { 171539d35262SVincent Le Chenadec PetscFunctionBegin; 171639d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 171739d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,2); 171839d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,4); 171939d35262SVincent Le Chenadec PetscFunctionReturn(0); 172039d35262SVincent Le Chenadec } 172139d35262SVincent Le Chenadec 172239d35262SVincent Le Chenadec PetscErrorCode DMLocalToLocalBegin_Composite(DM dm,Vec vec1,InsertMode mode,Vec vec2) 172339d35262SVincent Le Chenadec { 172439d35262SVincent Le Chenadec PetscErrorCode ierr; 172539d35262SVincent Le Chenadec struct DMCompositeLink *next; 172639d35262SVincent Le Chenadec PetscScalar *array1,*array2; 172739d35262SVincent Le Chenadec DM_Composite *com = (DM_Composite*)dm->data; 172839d35262SVincent Le Chenadec 172939d35262SVincent Le Chenadec PetscFunctionBegin; 173039d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 173139d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec1,VEC_CLASSID,2); 173239d35262SVincent Le Chenadec PetscValidHeaderSpecific(vec2,VEC_CLASSID,4); 173339d35262SVincent Le Chenadec 173439d35262SVincent Le Chenadec if (!com->setup) { 173539d35262SVincent Le Chenadec ierr = DMSetUp(dm);CHKERRQ(ierr); 173639d35262SVincent Le Chenadec } 173739d35262SVincent Le Chenadec 173839d35262SVincent Le Chenadec ierr = VecGetArray(vec1,&array1);CHKERRQ(ierr); 173939d35262SVincent Le Chenadec ierr = VecGetArray(vec2,&array2);CHKERRQ(ierr); 174039d35262SVincent Le Chenadec 174139d35262SVincent Le Chenadec /* loop over packed objects, handling one at at time */ 174239d35262SVincent Le Chenadec next = com->next; 174339d35262SVincent Le Chenadec while (next) { 174439d35262SVincent Le Chenadec Vec local1,local2; 174539d35262SVincent Le Chenadec 174639d35262SVincent Le Chenadec ierr = DMGetLocalVector(next->dm,&local1);CHKERRQ(ierr); 174739d35262SVincent Le Chenadec ierr = VecPlaceArray(local1,array1);CHKERRQ(ierr); 174839d35262SVincent Le Chenadec ierr = DMGetLocalVector(next->dm,&local2);CHKERRQ(ierr); 174939d35262SVincent Le Chenadec ierr = VecPlaceArray(local2,array2);CHKERRQ(ierr); 175039d35262SVincent Le Chenadec ierr = DMLocalToLocalBegin(next->dm,local1,mode,local2);CHKERRQ(ierr); 175139d35262SVincent Le Chenadec ierr = DMLocalToLocalEnd(next->dm,local1,mode,local2);CHKERRQ(ierr); 175239d35262SVincent Le Chenadec ierr = VecResetArray(local2);CHKERRQ(ierr); 175339d35262SVincent Le Chenadec ierr = DMRestoreLocalVector(next->dm,&local2);CHKERRQ(ierr); 175439d35262SVincent Le Chenadec ierr = VecResetArray(local1);CHKERRQ(ierr); 175539d35262SVincent Le Chenadec ierr = DMRestoreLocalVector(next->dm,&local1);CHKERRQ(ierr); 175639d35262SVincent Le Chenadec 175739d35262SVincent Le Chenadec array1 += next->nlocal; 175839d35262SVincent Le Chenadec array2 += next->nlocal; 175939d35262SVincent Le Chenadec next = next->next; 176039d35262SVincent Le Chenadec } 176139d35262SVincent Le Chenadec 176239d35262SVincent Le Chenadec ierr = VecRestoreArray(vec1,NULL);CHKERRQ(ierr); 176339d35262SVincent Le Chenadec ierr = VecRestoreArray(vec2,NULL);CHKERRQ(ierr); 176439d35262SVincent Le Chenadec 176539d35262SVincent Le Chenadec PetscFunctionReturn(0); 176639d35262SVincent Le Chenadec } 176739d35262SVincent Le Chenadec 176839d35262SVincent Le Chenadec PetscErrorCode DMLocalToLocalEnd_Composite(DM dm,Vec lvec,InsertMode mode,Vec gvec) 176939d35262SVincent Le Chenadec { 177039d35262SVincent Le Chenadec PetscFunctionBegin; 177139d35262SVincent Le Chenadec PetscValidHeaderSpecific(dm,DM_CLASSID,1); 177239d35262SVincent Le Chenadec PetscValidHeaderSpecific(lvec,VEC_CLASSID,2); 177339d35262SVincent Le Chenadec PetscValidHeaderSpecific(gvec,VEC_CLASSID,4); 17740c010503SBarry Smith PetscFunctionReturn(0); 17750c010503SBarry Smith } 177647c6ae99SBarry Smith 17776ae3a549SBarry Smith /*MC 17786ae3a549SBarry Smith DMCOMPOSITE = "composite" - A DM object that is used to manage data for a collection of DMs 17796ae3a549SBarry Smith 17806ae3a549SBarry Smith Level: intermediate 17816ae3a549SBarry Smith 17821abcec8cSBarry Smith .seealso: DMType, DM, DMDACreate(), DMCreate(), DMSetType(), DMCompositeCreate() 17836ae3a549SBarry Smith M*/ 17846ae3a549SBarry Smith 17856ae3a549SBarry Smith 17868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Composite(DM p) 1787a4121054SBarry Smith { 1788a4121054SBarry Smith PetscErrorCode ierr; 1789a4121054SBarry Smith DM_Composite *com; 1790a4121054SBarry Smith 1791a4121054SBarry Smith PetscFunctionBegin; 1792b00a9115SJed Brown ierr = PetscNewLog(p,&com);CHKERRQ(ierr); 1793a4121054SBarry Smith p->data = com; 1794a4121054SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr); 1795a4121054SBarry Smith com->n = 0; 17967ac2b803SAlex Fikl com->nghost = 0; 17970298fd71SBarry Smith com->next = NULL; 1798a4121054SBarry Smith com->nDM = 0; 1799a4121054SBarry Smith 1800a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Composite; 1801a4121054SBarry Smith p->ops->createlocalvector = DMCreateLocalVector_Composite; 1802184d77edSJed Brown p->ops->getlocaltoglobalmapping = DMGetLocalToGlobalMapping_Composite; 18034d343eeaSMatthew G Knepley p->ops->createfieldis = DMCreateFieldIS_Composite; 180416621825SDmitry Karpeev p->ops->createfielddecomposition = DMCreateFieldDecomposition_Composite; 1805a4121054SBarry Smith p->ops->refine = DMRefine_Composite; 180614354c39SJed Brown p->ops->coarsen = DMCoarsen_Composite; 180725296bd5SBarry Smith p->ops->createinterpolation = DMCreateInterpolation_Composite; 180825296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Composite; 1809e727c939SJed Brown p->ops->getcoloring = DMCreateColoring_Composite; 1810a4121054SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite; 1811a4121054SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite; 181239d35262SVincent Le Chenadec p->ops->localtoglobalbegin = DMLocalToGlobalBegin_Composite; 181339d35262SVincent Le Chenadec p->ops->localtoglobalend = DMLocalToGlobalEnd_Composite; 181439d35262SVincent Le Chenadec p->ops->localtolocalbegin = DMLocalToLocalBegin_Composite; 181539d35262SVincent Le Chenadec p->ops->localtolocalend = DMLocalToLocalEnd_Composite; 1816a4121054SBarry Smith p->ops->destroy = DMDestroy_Composite; 1817a4121054SBarry Smith p->ops->view = DMView_Composite; 1818a4121054SBarry Smith p->ops->setup = DMSetUp_Composite; 1819e10fd815SStefano Zampini 1820e10fd815SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)p,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Composite);CHKERRQ(ierr); 1821a4121054SBarry Smith PetscFunctionReturn(0); 1822a4121054SBarry Smith } 1823a4121054SBarry Smith 18241fd49c25SBarry Smith /*@ 18250c010503SBarry Smith DMCompositeCreate - Creates a vector packer, used to generate "composite" 18260c010503SBarry Smith vectors made up of several subvectors. 18270c010503SBarry Smith 18280c010503SBarry Smith Collective on MPI_Comm 182947c6ae99SBarry Smith 183047c6ae99SBarry Smith Input Parameter: 18310c010503SBarry Smith . comm - the processors that will share the global vector 18320c010503SBarry Smith 18330c010503SBarry Smith Output Parameters: 18340c010503SBarry Smith . packer - the packer object 183547c6ae99SBarry Smith 183647c6ae99SBarry Smith Level: advanced 183747c6ae99SBarry Smith 18381abcec8cSBarry Smith .seealso DMDestroy(), DMCompositeAddDM(), DMCompositeScatter(), DMCOMPOSITE,DMCreate() 18396eb61c8cSJed Brown DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess() 184047c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 184147c6ae99SBarry Smith 184247c6ae99SBarry Smith @*/ 18437087cfbeSBarry Smith PetscErrorCode DMCompositeCreate(MPI_Comm comm,DM *packer) 184447c6ae99SBarry Smith { 18450c010503SBarry Smith PetscErrorCode ierr; 18460c010503SBarry Smith 184747c6ae99SBarry Smith PetscFunctionBegin; 18480c010503SBarry Smith PetscValidPointer(packer,2); 1849a4121054SBarry Smith ierr = DMCreate(comm,packer);CHKERRQ(ierr); 1850a4121054SBarry Smith ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr); 185147c6ae99SBarry Smith PetscFunctionReturn(0); 185247c6ae99SBarry Smith } 1853