147c6ae99SBarry Smith #define PETSCDM_DLL 247c6ae99SBarry Smith 3e1589f56SBarry Smith #include "petscdm.h" /*I "petscdm.h" I*/ 447c6ae99SBarry Smith #include "private/dmimpl.h" 547c6ae99SBarry Smith #include "petscmat.h" 647c6ae99SBarry Smith 747c6ae99SBarry Smith /* 847c6ae99SBarry Smith rstart is where an array/subvector starts in the global parallel vector, so arrays 947c6ae99SBarry Smith rstarts are meaningless (and set to the previous one) except on the processor where the array lives 1047c6ae99SBarry Smith */ 1147c6ae99SBarry Smith 1247c6ae99SBarry Smith typedef enum {DMCOMPOSITE_ARRAY, DMCOMPOSITE_DM} DMCompositeLinkType; 1347c6ae99SBarry Smith 1447c6ae99SBarry Smith struct DMCompositeLink { 1547c6ae99SBarry Smith DMCompositeLinkType type; 1647c6ae99SBarry Smith struct DMCompositeLink *next; 1747c6ae99SBarry Smith PetscInt n,rstart; /* rstart is relative to this process */ 1847c6ae99SBarry Smith PetscInt grstart; /* grstart is relative to all processes */ 1947c6ae99SBarry Smith 2047c6ae99SBarry Smith /* only used for DMCOMPOSITE_DM */ 2147c6ae99SBarry Smith PetscInt *grstarts; /* global row for first unknown of this DM on each process */ 2247c6ae99SBarry Smith DM dm; 2347c6ae99SBarry Smith 2447c6ae99SBarry Smith /* only used for DMCOMPOSITE_ARRAY */ 2547c6ae99SBarry Smith PetscMPIInt rank; /* process where array unknowns live */ 2647c6ae99SBarry Smith }; 2747c6ae99SBarry Smith 2847c6ae99SBarry Smith typedef struct { 2947c6ae99SBarry Smith PetscInt n,N,rstart; /* rstart is relative to all processors, n unknowns owned by this process, N is total unknowns */ 30aa219208SBarry Smith PetscInt nghost; /* number of all local entries include DMDA ghost points and any shared redundant arrays */ 3147c6ae99SBarry Smith PetscInt nDM,nredundant,nmine; /* how many DM's and seperate redundant arrays used to build DM(nmine is ones on this process) */ 3247c6ae99SBarry Smith PetscBool setup; /* after this is set, cannot add new links to the DM*/ 3347c6ae99SBarry Smith struct DMCompositeLink *next; 3447c6ae99SBarry Smith 3547c6ae99SBarry Smith PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt); 3647c6ae99SBarry Smith } DM_Composite; 3747c6ae99SBarry Smith 3847c6ae99SBarry Smith #undef __FUNCT__ 3947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeSetCoupling" 4047c6ae99SBarry Smith /*@C 4147c6ae99SBarry Smith DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the 42aa219208SBarry Smith seperate components (DMDA's and arrays) in a DMto build the correct matrix nonzero structure. 4347c6ae99SBarry Smith 4447c6ae99SBarry Smith 4547c6ae99SBarry Smith Logically Collective on MPI_Comm 4647c6ae99SBarry Smith 4747c6ae99SBarry Smith Input Parameter: 4847c6ae99SBarry Smith + dm - the composite object 4947c6ae99SBarry Smith - formcouplelocations - routine to set the nonzero locations in the matrix 5047c6ae99SBarry Smith 5147c6ae99SBarry Smith Level: advanced 5247c6ae99SBarry Smith 5347c6ae99SBarry Smith Notes: See DMCompositeSetContext() and DMCompositeGetContext() for how to get user information into 5447c6ae99SBarry Smith this routine 5547c6ae99SBarry Smith 5647c6ae99SBarry Smith @*/ 5747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt)) 5847c6ae99SBarry Smith { 5947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 6047c6ae99SBarry Smith 6147c6ae99SBarry Smith PetscFunctionBegin; 6247c6ae99SBarry Smith com->FormCoupleLocations = FormCoupleLocations; 6347c6ae99SBarry Smith PetscFunctionReturn(0); 6447c6ae99SBarry Smith } 6547c6ae99SBarry Smith 6647c6ae99SBarry Smith #undef __FUNCT__ 6747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeSetContext" 6847c6ae99SBarry Smith /*@ 6947c6ae99SBarry Smith DMCompositeSetContext - Allows user to stash data they may need within the form coupling routine they 7047c6ae99SBarry Smith set with DMCompositeSetCoupling() 7147c6ae99SBarry Smith 7247c6ae99SBarry Smith 7347c6ae99SBarry Smith Not Collective 7447c6ae99SBarry Smith 7547c6ae99SBarry Smith Input Parameter: 7647c6ae99SBarry Smith + dm - the composite object 7747c6ae99SBarry Smith - ctx - the user supplied context 7847c6ae99SBarry Smith 7947c6ae99SBarry Smith Level: advanced 8047c6ae99SBarry Smith 8147c6ae99SBarry Smith Notes: Use DMCompositeGetContext() to retrieve the context when needed. 8247c6ae99SBarry Smith 8347c6ae99SBarry Smith @*/ 8447c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeSetContext(DM dm,void *ctx) 8547c6ae99SBarry Smith { 8647c6ae99SBarry Smith PetscFunctionBegin; 8747c6ae99SBarry Smith dm->ctx = ctx; 8847c6ae99SBarry Smith PetscFunctionReturn(0); 8947c6ae99SBarry Smith } 9047c6ae99SBarry Smith 9147c6ae99SBarry Smith #undef __FUNCT__ 9247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetContext" 9347c6ae99SBarry Smith /*@ 9447c6ae99SBarry Smith DMCompositeGetContext - Access the context set with DMCompositeSetContext() 9547c6ae99SBarry Smith 9647c6ae99SBarry Smith 9747c6ae99SBarry Smith Not Collective 9847c6ae99SBarry Smith 9947c6ae99SBarry Smith Input Parameter: 10047c6ae99SBarry Smith . dm - the composite object 10147c6ae99SBarry Smith 10247c6ae99SBarry Smith Output Parameter: 10347c6ae99SBarry Smith . ctx - the user supplied context 10447c6ae99SBarry Smith 10547c6ae99SBarry Smith Level: advanced 10647c6ae99SBarry Smith 10747c6ae99SBarry Smith Notes: Use DMCompositeGetContext() to retrieve the context when needed. 10847c6ae99SBarry Smith 10947c6ae99SBarry Smith @*/ 11047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetContext(DM dm,void **ctx) 11147c6ae99SBarry Smith { 11247c6ae99SBarry Smith PetscFunctionBegin; 11347c6ae99SBarry Smith *ctx = dm->ctx; 11447c6ae99SBarry Smith PetscFunctionReturn(0); 11547c6ae99SBarry Smith } 11647c6ae99SBarry Smith 11747c6ae99SBarry Smith 11847c6ae99SBarry Smith 11947c6ae99SBarry Smith extern PetscErrorCode DMDestroy_Private(DM,PetscBool *); 12047c6ae99SBarry Smith 12147c6ae99SBarry Smith #undef __FUNCT__ 1220c010503SBarry Smith #define __FUNCT__ "DMDestroy_Composite" 1230c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDestroy_Composite(DM dm) 12447c6ae99SBarry Smith { 12547c6ae99SBarry Smith PetscErrorCode ierr; 12647c6ae99SBarry Smith struct DMCompositeLink *next, *prev; 12747c6ae99SBarry Smith PetscBool done; 12847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 12947c6ae99SBarry Smith 13047c6ae99SBarry Smith PetscFunctionBegin; 13147c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 13247c6ae99SBarry Smith ierr = DMDestroy_Private((DM)dm,&done);CHKERRQ(ierr); 13347c6ae99SBarry Smith if (!done) PetscFunctionReturn(0); 13447c6ae99SBarry Smith 13547c6ae99SBarry Smith next = com->next; 13647c6ae99SBarry Smith while (next) { 13747c6ae99SBarry Smith prev = next; 13847c6ae99SBarry Smith next = next->next; 13947c6ae99SBarry Smith if (prev->type == DMCOMPOSITE_DM) { 14047c6ae99SBarry Smith ierr = DMDestroy(prev->dm);CHKERRQ(ierr); 14147c6ae99SBarry Smith } 14247c6ae99SBarry Smith if (prev->grstarts) { 14347c6ae99SBarry Smith ierr = PetscFree(prev->grstarts);CHKERRQ(ierr); 14447c6ae99SBarry Smith } 14547c6ae99SBarry Smith ierr = PetscFree(prev);CHKERRQ(ierr); 14647c6ae99SBarry Smith } 14747c6ae99SBarry Smith ierr = PetscFree(com);CHKERRQ(ierr); 14847c6ae99SBarry Smith ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr); 14947c6ae99SBarry Smith PetscFunctionReturn(0); 15047c6ae99SBarry Smith } 15147c6ae99SBarry Smith 15247c6ae99SBarry Smith #undef __FUNCT__ 1530c010503SBarry Smith #define __FUNCT__ "DMView_Composite" 1540c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMView_Composite(DM dm,PetscViewer v) 15547c6ae99SBarry Smith { 15647c6ae99SBarry Smith PetscErrorCode ierr; 15747c6ae99SBarry Smith PetscBool iascii; 15847c6ae99SBarry Smith DM_Composite *com = (DM_Composite *)dm->data; 15947c6ae99SBarry Smith 16047c6ae99SBarry Smith PetscFunctionBegin; 16147c6ae99SBarry Smith ierr = PetscTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 16247c6ae99SBarry Smith if (iascii) { 16347c6ae99SBarry Smith struct DMCompositeLink *lnk = com->next; 16447c6ae99SBarry Smith PetscInt i; 16547c6ae99SBarry Smith 16647c6ae99SBarry Smith ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix?((PetscObject)dm)->prefix:"no prefix");CHKERRQ(ierr); 16747c6ae99SBarry Smith ierr = PetscViewerASCIIPrintf(v," contains %d DMs and %d redundant arrays\n",com->nDM,com->nredundant);CHKERRQ(ierr); 16847c6ae99SBarry Smith ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 16947c6ae99SBarry Smith for (i=0; lnk; lnk=lnk->next,i++) { 17047c6ae99SBarry Smith if (lnk->dm) { 17147c6ae99SBarry Smith ierr = PetscViewerASCIIPrintf(v,"Link %d: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr); 17247c6ae99SBarry Smith ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 17347c6ae99SBarry Smith ierr = DMView(lnk->dm,v);CHKERRQ(ierr); 17447c6ae99SBarry Smith ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 17547c6ae99SBarry Smith } else { 17647c6ae99SBarry Smith ierr = PetscViewerASCIIPrintf(v,"Link %d: Redundant array of size %d owned by rank %d\n",i,lnk->n,lnk->rank);CHKERRQ(ierr); 17747c6ae99SBarry Smith } 17847c6ae99SBarry Smith } 17947c6ae99SBarry Smith ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 18047c6ae99SBarry Smith } 18147c6ae99SBarry Smith PetscFunctionReturn(0); 18247c6ae99SBarry Smith } 18347c6ae99SBarry Smith 18447c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/ 18547c6ae99SBarry Smith #undef __FUNCT__ 186d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp_Composite" 187d7bf68aeSBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSetUp_Composite(DM dm) 18847c6ae99SBarry Smith { 18947c6ae99SBarry Smith PetscErrorCode ierr; 19047c6ae99SBarry Smith PetscInt nprev = 0; 19147c6ae99SBarry Smith PetscMPIInt rank,size; 19247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 19347c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 19447c6ae99SBarry Smith PetscLayout map; 19547c6ae99SBarry Smith 19647c6ae99SBarry Smith PetscFunctionBegin; 19747c6ae99SBarry Smith if (com->setup) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup"); 19847c6ae99SBarry Smith ierr = PetscLayoutCreate(((PetscObject)dm)->comm,&map);CHKERRQ(ierr); 19947c6ae99SBarry Smith ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr); 20047c6ae99SBarry Smith ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr); 20147c6ae99SBarry Smith ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr); 20247c6ae99SBarry Smith ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 20347c6ae99SBarry Smith ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr); 20447c6ae99SBarry Smith ierr = PetscLayoutGetRange(map,&com->rstart,PETSC_NULL);CHKERRQ(ierr); 20547c6ae99SBarry Smith ierr = PetscLayoutDestroy(map);CHKERRQ(ierr); 20647c6ae99SBarry Smith 20747c6ae99SBarry Smith /* now set the rstart for each linked array/vector */ 20847c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 20947c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)dm)->comm,&size);CHKERRQ(ierr); 21047c6ae99SBarry Smith while (next) { 21147c6ae99SBarry Smith next->rstart = nprev; 21247c6ae99SBarry Smith if ((rank == next->rank) || next->type != DMCOMPOSITE_ARRAY) nprev += next->n; 21347c6ae99SBarry Smith next->grstart = com->rstart + next->rstart; 21447c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 21547c6ae99SBarry Smith ierr = MPI_Bcast(&next->grstart,1,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr); 21647c6ae99SBarry Smith } else { 21747c6ae99SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&next->grstarts);CHKERRQ(ierr); 21847c6ae99SBarry Smith ierr = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr); 21947c6ae99SBarry Smith } 22047c6ae99SBarry Smith next = next->next; 22147c6ae99SBarry Smith } 22247c6ae99SBarry Smith com->setup = PETSC_TRUE; 22347c6ae99SBarry Smith PetscFunctionReturn(0); 22447c6ae99SBarry Smith } 22547c6ae99SBarry Smith 22647c6ae99SBarry Smith 22747c6ae99SBarry Smith #undef __FUNCT__ 22847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess_Array" 22947c6ae99SBarry Smith PetscErrorCode DMCompositeGetAccess_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar **array) 23047c6ae99SBarry Smith { 23147c6ae99SBarry Smith PetscErrorCode ierr; 23247c6ae99SBarry Smith PetscScalar *varray; 23347c6ae99SBarry Smith PetscMPIInt rank; 23447c6ae99SBarry Smith 23547c6ae99SBarry Smith PetscFunctionBegin; 23647c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 23747c6ae99SBarry Smith if (array) { 23847c6ae99SBarry Smith if (rank == mine->rank) { 23947c6ae99SBarry Smith ierr = VecGetArray(vec,&varray);CHKERRQ(ierr); 24047c6ae99SBarry Smith *array = varray + mine->rstart; 24147c6ae99SBarry Smith ierr = VecRestoreArray(vec,&varray);CHKERRQ(ierr); 24247c6ae99SBarry Smith } else { 24347c6ae99SBarry Smith *array = 0; 24447c6ae99SBarry Smith } 24547c6ae99SBarry Smith } 24647c6ae99SBarry Smith PetscFunctionReturn(0); 24747c6ae99SBarry Smith } 24847c6ae99SBarry Smith 24947c6ae99SBarry Smith #undef __FUNCT__ 25047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess_DM" 25147c6ae99SBarry Smith PetscErrorCode DMCompositeGetAccess_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec *global) 25247c6ae99SBarry Smith { 25347c6ae99SBarry Smith PetscErrorCode ierr; 25447c6ae99SBarry Smith PetscScalar *array; 25547c6ae99SBarry Smith 25647c6ae99SBarry Smith PetscFunctionBegin; 25747c6ae99SBarry Smith if (global) { 25847c6ae99SBarry Smith ierr = DMGetGlobalVector(mine->dm,global);CHKERRQ(ierr); 25947c6ae99SBarry Smith ierr = VecGetArray(vec,&array);CHKERRQ(ierr); 26047c6ae99SBarry Smith ierr = VecPlaceArray(*global,array+mine->rstart);CHKERRQ(ierr); 26147c6ae99SBarry Smith ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr); 26247c6ae99SBarry Smith } 26347c6ae99SBarry Smith PetscFunctionReturn(0); 26447c6ae99SBarry Smith } 26547c6ae99SBarry Smith 26647c6ae99SBarry Smith #undef __FUNCT__ 26747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess_Array" 26847c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreAccess_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar **array) 26947c6ae99SBarry Smith { 27047c6ae99SBarry Smith PetscFunctionBegin; 27147c6ae99SBarry Smith PetscFunctionReturn(0); 27247c6ae99SBarry Smith } 27347c6ae99SBarry Smith 27447c6ae99SBarry Smith #undef __FUNCT__ 27547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess_DM" 27647c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreAccess_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec *global) 27747c6ae99SBarry Smith { 27847c6ae99SBarry Smith PetscErrorCode ierr; 27947c6ae99SBarry Smith 28047c6ae99SBarry Smith PetscFunctionBegin; 28147c6ae99SBarry Smith if (global) { 28247c6ae99SBarry Smith ierr = VecResetArray(*global);CHKERRQ(ierr); 28347c6ae99SBarry Smith ierr = DMRestoreGlobalVector(mine->dm,global);CHKERRQ(ierr); 28447c6ae99SBarry Smith } 28547c6ae99SBarry Smith PetscFunctionReturn(0); 28647c6ae99SBarry Smith } 28747c6ae99SBarry Smith 28847c6ae99SBarry Smith #undef __FUNCT__ 28947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter_Array" 29047c6ae99SBarry Smith PetscErrorCode DMCompositeScatter_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar *array) 29147c6ae99SBarry Smith { 29247c6ae99SBarry Smith PetscErrorCode ierr; 29347c6ae99SBarry Smith PetscScalar *varray; 29447c6ae99SBarry Smith PetscMPIInt rank; 29547c6ae99SBarry Smith 29647c6ae99SBarry Smith PetscFunctionBegin; 29747c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 29847c6ae99SBarry Smith if (rank == mine->rank) { 29947c6ae99SBarry Smith ierr = VecGetArray(vec,&varray);CHKERRQ(ierr); 30047c6ae99SBarry Smith ierr = PetscMemcpy(array,varray+mine->rstart,mine->n*sizeof(PetscScalar));CHKERRQ(ierr); 30147c6ae99SBarry Smith ierr = VecRestoreArray(vec,&varray);CHKERRQ(ierr); 30247c6ae99SBarry Smith } 30347c6ae99SBarry Smith ierr = MPI_Bcast(array,mine->n,MPIU_SCALAR,mine->rank,((PetscObject)dm)->comm);CHKERRQ(ierr); 30447c6ae99SBarry Smith PetscFunctionReturn(0); 30547c6ae99SBarry Smith } 30647c6ae99SBarry Smith 30747c6ae99SBarry Smith #undef __FUNCT__ 30847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter_DM" 30947c6ae99SBarry Smith PetscErrorCode DMCompositeScatter_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec local) 31047c6ae99SBarry Smith { 31147c6ae99SBarry Smith PetscErrorCode ierr; 31247c6ae99SBarry Smith PetscScalar *array; 31347c6ae99SBarry Smith Vec global; 31447c6ae99SBarry Smith 31547c6ae99SBarry Smith PetscFunctionBegin; 31647c6ae99SBarry Smith ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr); 31747c6ae99SBarry Smith ierr = VecGetArray(vec,&array);CHKERRQ(ierr); 31847c6ae99SBarry Smith ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr); 31947c6ae99SBarry Smith ierr = DMGlobalToLocalBegin(mine->dm,global,INSERT_VALUES,local);CHKERRQ(ierr); 32047c6ae99SBarry Smith ierr = DMGlobalToLocalEnd(mine->dm,global,INSERT_VALUES,local);CHKERRQ(ierr); 32147c6ae99SBarry Smith ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr); 32247c6ae99SBarry Smith ierr = VecResetArray(global);CHKERRQ(ierr); 32347c6ae99SBarry Smith ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr); 32447c6ae99SBarry Smith PetscFunctionReturn(0); 32547c6ae99SBarry Smith } 32647c6ae99SBarry Smith 32747c6ae99SBarry Smith #undef __FUNCT__ 32847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_Array" 329*df0c820aSJed Brown PetscErrorCode DMCompositeGather_Array(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,PetscScalar *array) 33047c6ae99SBarry Smith { 33147c6ae99SBarry Smith PetscErrorCode ierr; 33247c6ae99SBarry Smith PetscScalar *varray; 33347c6ae99SBarry Smith PetscMPIInt rank; 33447c6ae99SBarry Smith 33547c6ae99SBarry Smith PetscFunctionBegin; 33647c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 33747c6ae99SBarry Smith if (rank == mine->rank) { 33847c6ae99SBarry Smith ierr = VecGetArray(vec,&varray);CHKERRQ(ierr); 33947c6ae99SBarry Smith if (varray+mine->rstart == array) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"You need not DMCompositeGather() into objects obtained via DMCompositeGetAccess()"); 340*df0c820aSJed Brown switch (imode) { 341*df0c820aSJed Brown case INSERT_VALUES: 34247c6ae99SBarry Smith ierr = PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(PetscScalar));CHKERRQ(ierr); 343*df0c820aSJed Brown break; 344*df0c820aSJed Brown case ADD_VALUES: { 345*df0c820aSJed Brown PetscInt i; 346*df0c820aSJed Brown for (i=0; i<mine->n; i++) varray[mine->rstart+i] += array[i]; 347*df0c820aSJed Brown } break; 348*df0c820aSJed Brown default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"imode"); 349*df0c820aSJed Brown } 35047c6ae99SBarry Smith ierr = VecRestoreArray(vec,&varray);CHKERRQ(ierr); 35147c6ae99SBarry Smith } 35247c6ae99SBarry Smith PetscFunctionReturn(0); 35347c6ae99SBarry Smith } 35447c6ae99SBarry Smith 35547c6ae99SBarry Smith #undef __FUNCT__ 35647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_DM" 357*df0c820aSJed Brown PetscErrorCode DMCompositeGather_DM(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,Vec local) 35847c6ae99SBarry Smith { 35947c6ae99SBarry Smith PetscErrorCode ierr; 36047c6ae99SBarry Smith PetscScalar *array; 36147c6ae99SBarry Smith Vec global; 36247c6ae99SBarry Smith 36347c6ae99SBarry Smith PetscFunctionBegin; 36447c6ae99SBarry Smith ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr); 36547c6ae99SBarry Smith ierr = VecGetArray(vec,&array);CHKERRQ(ierr); 36647c6ae99SBarry Smith ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr); 367*df0c820aSJed Brown ierr = DMLocalToGlobalBegin(mine->dm,local,imode,global);CHKERRQ(ierr); 368*df0c820aSJed Brown ierr = DMLocalToGlobalEnd(mine->dm,local,imode,global);CHKERRQ(ierr); 36947c6ae99SBarry Smith ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr); 37047c6ae99SBarry Smith ierr = VecResetArray(global);CHKERRQ(ierr); 37147c6ae99SBarry Smith ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr); 37247c6ae99SBarry Smith PetscFunctionReturn(0); 37347c6ae99SBarry Smith } 37447c6ae99SBarry Smith 37547c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/ 37647c6ae99SBarry Smith 37747c6ae99SBarry Smith #include <stdarg.h> 37847c6ae99SBarry Smith 37947c6ae99SBarry Smith #undef __FUNCT__ 38047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetNumberDM" 38147c6ae99SBarry Smith /*@C 38247c6ae99SBarry Smith DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite 38347c6ae99SBarry Smith representation. 38447c6ae99SBarry Smith 38547c6ae99SBarry Smith Not Collective 38647c6ae99SBarry Smith 38747c6ae99SBarry Smith Input Parameter: 38847c6ae99SBarry Smith . dm - the packer object 38947c6ae99SBarry Smith 39047c6ae99SBarry Smith Output Parameter: 39147c6ae99SBarry Smith . nDM - the number of DMs 39247c6ae99SBarry Smith 39347c6ae99SBarry Smith Level: beginner 39447c6ae99SBarry Smith 39547c6ae99SBarry Smith @*/ 39647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetNumberDM(DM dm,PetscInt *nDM) 39747c6ae99SBarry Smith { 39847c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 39947c6ae99SBarry Smith PetscFunctionBegin; 40047c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 40147c6ae99SBarry Smith *nDM = com->nDM; 40247c6ae99SBarry Smith PetscFunctionReturn(0); 40347c6ae99SBarry Smith } 40447c6ae99SBarry Smith 40547c6ae99SBarry Smith 40647c6ae99SBarry Smith #undef __FUNCT__ 40747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess" 40847c6ae99SBarry Smith /*@C 40947c6ae99SBarry Smith DMCompositeGetAccess - Allows one to access the individual packed vectors in their global 41047c6ae99SBarry Smith representation. 41147c6ae99SBarry Smith 41247c6ae99SBarry Smith Collective on DMComposite 41347c6ae99SBarry Smith 41447c6ae99SBarry Smith Input Parameter: 41547c6ae99SBarry Smith + dm - the packer object 41647c6ae99SBarry Smith . gvec - the global vector 41747c6ae99SBarry Smith - ... - the individual sequential or parallel objects (arrays or vectors) 41847c6ae99SBarry Smith 41947c6ae99SBarry Smith Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them 42047c6ae99SBarry Smith 42147c6ae99SBarry Smith Level: advanced 42247c6ae99SBarry Smith 42347c6ae99SBarry Smith @*/ 42447c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetAccess(DM dm,Vec gvec,...) 42547c6ae99SBarry Smith { 42647c6ae99SBarry Smith va_list Argp; 42747c6ae99SBarry Smith PetscErrorCode ierr; 42847c6ae99SBarry Smith struct DMCompositeLink *next; 42947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 43047c6ae99SBarry Smith 43147c6ae99SBarry Smith PetscFunctionBegin; 43247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 43347c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 43447c6ae99SBarry Smith next = com->next; 43547c6ae99SBarry Smith if (!com->setup) { 436d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 43747c6ae99SBarry Smith } 43847c6ae99SBarry Smith 43947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 44047c6ae99SBarry Smith va_start(Argp,gvec); 44147c6ae99SBarry Smith while (next) { 44247c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 44347c6ae99SBarry Smith PetscScalar **array; 44447c6ae99SBarry Smith array = va_arg(Argp, PetscScalar**); 44547c6ae99SBarry Smith ierr = DMCompositeGetAccess_Array(dm,next,gvec,array);CHKERRQ(ierr); 44647c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 44747c6ae99SBarry Smith Vec *vec; 44847c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 44947c6ae99SBarry Smith ierr = DMCompositeGetAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr); 45047c6ae99SBarry Smith } else { 45147c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 45247c6ae99SBarry Smith } 45347c6ae99SBarry Smith next = next->next; 45447c6ae99SBarry Smith } 45547c6ae99SBarry Smith va_end(Argp); 45647c6ae99SBarry Smith PetscFunctionReturn(0); 45747c6ae99SBarry Smith } 45847c6ae99SBarry Smith 45947c6ae99SBarry Smith #undef __FUNCT__ 46047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess" 46147c6ae99SBarry Smith /*@C 462aa219208SBarry Smith DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess() 46347c6ae99SBarry Smith representation. 46447c6ae99SBarry Smith 46547c6ae99SBarry Smith Collective on DMComposite 46647c6ae99SBarry Smith 46747c6ae99SBarry Smith Input Parameter: 46847c6ae99SBarry Smith + dm - the packer object 46947c6ae99SBarry Smith . gvec - the global vector 47047c6ae99SBarry Smith - ... - the individual sequential or parallel objects (arrays or vectors) 47147c6ae99SBarry Smith 47247c6ae99SBarry Smith Level: advanced 47347c6ae99SBarry Smith 4740c010503SBarry Smith .seealso DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 47547c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeScatter(), 476aa219208SBarry Smith DMCompositeRestoreAccess(), DMCompositeGetAccess() 47747c6ae99SBarry Smith 47847c6ae99SBarry Smith @*/ 47947c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreAccess(DM dm,Vec gvec,...) 48047c6ae99SBarry Smith { 48147c6ae99SBarry Smith va_list Argp; 48247c6ae99SBarry Smith PetscErrorCode ierr; 48347c6ae99SBarry Smith struct DMCompositeLink *next; 48447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 48547c6ae99SBarry Smith 48647c6ae99SBarry Smith PetscFunctionBegin; 48747c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 48847c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 48947c6ae99SBarry Smith next = com->next; 49047c6ae99SBarry Smith if (!com->setup) { 491d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 49247c6ae99SBarry Smith } 49347c6ae99SBarry Smith 49447c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 49547c6ae99SBarry Smith va_start(Argp,gvec); 49647c6ae99SBarry Smith while (next) { 49747c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 49847c6ae99SBarry Smith PetscScalar **array; 49947c6ae99SBarry Smith array = va_arg(Argp, PetscScalar**); 50047c6ae99SBarry Smith ierr = DMCompositeRestoreAccess_Array(dm,next,gvec,array);CHKERRQ(ierr); 50147c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 50247c6ae99SBarry Smith Vec *vec; 50347c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 50447c6ae99SBarry Smith ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr); 50547c6ae99SBarry Smith } else { 50647c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 50747c6ae99SBarry Smith } 50847c6ae99SBarry Smith next = next->next; 50947c6ae99SBarry Smith } 51047c6ae99SBarry Smith va_end(Argp); 51147c6ae99SBarry Smith PetscFunctionReturn(0); 51247c6ae99SBarry Smith } 51347c6ae99SBarry Smith 51447c6ae99SBarry Smith #undef __FUNCT__ 51547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter" 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 52147c6ae99SBarry Smith Input Parameter: 52247c6ae99SBarry Smith + dm - the packer object 52347c6ae99SBarry Smith . gvec - the global vector 52447c6ae99SBarry Smith - ... - the individual sequential objects (arrays or vectors) 52547c6ae99SBarry Smith 52647c6ae99SBarry Smith Level: advanced 52747c6ae99SBarry Smith 5280c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 52947c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 53047c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 53147c6ae99SBarry Smith 53247c6ae99SBarry Smith @*/ 53347c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeScatter(DM dm,Vec gvec,...) 53447c6ae99SBarry Smith { 53547c6ae99SBarry Smith va_list Argp; 53647c6ae99SBarry Smith PetscErrorCode ierr; 53747c6ae99SBarry Smith struct DMCompositeLink *next; 53847c6ae99SBarry Smith PetscInt cnt = 3; 53947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 54047c6ae99SBarry Smith 54147c6ae99SBarry Smith PetscFunctionBegin; 54247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 54347c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 54447c6ae99SBarry Smith next = com->next; 54547c6ae99SBarry Smith if (!com->setup) { 546d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 54747c6ae99SBarry Smith } 54847c6ae99SBarry Smith 54947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 55047c6ae99SBarry Smith va_start(Argp,gvec); 55147c6ae99SBarry Smith while (next) { 55247c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 55347c6ae99SBarry Smith PetscScalar *array; 55447c6ae99SBarry Smith array = va_arg(Argp, PetscScalar*); 55547c6ae99SBarry Smith ierr = DMCompositeScatter_Array(dm,next,gvec,array);CHKERRQ(ierr); 55647c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 55747c6ae99SBarry Smith Vec vec; 55847c6ae99SBarry Smith vec = va_arg(Argp, Vec); 55947c6ae99SBarry Smith PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt); 56047c6ae99SBarry Smith ierr = DMCompositeScatter_DM(dm,next,gvec,vec);CHKERRQ(ierr); 56147c6ae99SBarry Smith } else { 56247c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 56347c6ae99SBarry Smith } 56447c6ae99SBarry Smith cnt++; 56547c6ae99SBarry Smith next = next->next; 56647c6ae99SBarry Smith } 56747c6ae99SBarry Smith va_end(Argp); 56847c6ae99SBarry Smith PetscFunctionReturn(0); 56947c6ae99SBarry Smith } 57047c6ae99SBarry Smith 57147c6ae99SBarry Smith #undef __FUNCT__ 57247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather" 57347c6ae99SBarry Smith /*@C 57447c6ae99SBarry Smith DMCompositeGather - Gathers into a global packed vector from its individual local vectors 57547c6ae99SBarry Smith 57647c6ae99SBarry Smith Collective on DMComposite 57747c6ae99SBarry Smith 57847c6ae99SBarry Smith Input Parameter: 57947c6ae99SBarry Smith + dm - the packer object 58047c6ae99SBarry Smith . gvec - the global vector 58147c6ae99SBarry Smith - ... - the individual sequential objects (arrays or vectors) 58247c6ae99SBarry Smith 58347c6ae99SBarry Smith Level: advanced 58447c6ae99SBarry Smith 5850c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 58647c6ae99SBarry Smith DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 58747c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 58847c6ae99SBarry Smith 58947c6ae99SBarry Smith @*/ 590*df0c820aSJed Brown PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...) 59147c6ae99SBarry Smith { 59247c6ae99SBarry Smith va_list Argp; 59347c6ae99SBarry Smith PetscErrorCode ierr; 59447c6ae99SBarry Smith struct DMCompositeLink *next; 59547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 59647c6ae99SBarry Smith 59747c6ae99SBarry Smith PetscFunctionBegin; 59847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 59947c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 60047c6ae99SBarry Smith next = com->next; 60147c6ae99SBarry Smith if (!com->setup) { 602d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 60347c6ae99SBarry Smith } 60447c6ae99SBarry Smith 60547c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 606*df0c820aSJed Brown va_start(Argp,imode); 60747c6ae99SBarry Smith while (next) { 60847c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 60947c6ae99SBarry Smith PetscScalar *array; 61047c6ae99SBarry Smith array = va_arg(Argp, PetscScalar*); 611*df0c820aSJed Brown ierr = DMCompositeGather_Array(dm,next,gvec,imode,array);CHKERRQ(ierr); 61247c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 61347c6ae99SBarry Smith Vec vec; 61447c6ae99SBarry Smith vec = va_arg(Argp, Vec); 61547c6ae99SBarry Smith PetscValidHeaderSpecific(vec,VEC_CLASSID,3); 616*df0c820aSJed Brown ierr = DMCompositeGather_DM(dm,next,gvec,imode,vec);CHKERRQ(ierr); 61747c6ae99SBarry Smith } else { 61847c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 61947c6ae99SBarry Smith } 62047c6ae99SBarry Smith next = next->next; 62147c6ae99SBarry Smith } 62247c6ae99SBarry Smith va_end(Argp); 62347c6ae99SBarry Smith PetscFunctionReturn(0); 62447c6ae99SBarry Smith } 62547c6ae99SBarry Smith 62647c6ae99SBarry Smith #undef __FUNCT__ 62747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddArray" 62847c6ae99SBarry Smith /*@C 62947c6ae99SBarry Smith DMCompositeAddArray - adds an "redundant" array to a DMComposite. The array values will 63047c6ae99SBarry Smith be stored in part of the array on process orank. 63147c6ae99SBarry Smith 63247c6ae99SBarry Smith Collective on DMComposite 63347c6ae99SBarry Smith 63447c6ae99SBarry Smith Input Parameter: 63547c6ae99SBarry Smith + dm - the packer object 63647c6ae99SBarry Smith . orank - the process on which the array entries officially live, this number must be 63747c6ae99SBarry Smith the same on all processes. 63847c6ae99SBarry Smith - n - the length of the array 63947c6ae99SBarry Smith 64047c6ae99SBarry Smith Level: advanced 64147c6ae99SBarry Smith 6420c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(), 64347c6ae99SBarry Smith DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 64447c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 64547c6ae99SBarry Smith 64647c6ae99SBarry Smith @*/ 64747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeAddArray(DM dm,PetscMPIInt orank,PetscInt n) 64847c6ae99SBarry Smith { 64947c6ae99SBarry Smith struct DMCompositeLink *mine,*next; 65047c6ae99SBarry Smith PetscErrorCode ierr; 65147c6ae99SBarry Smith PetscMPIInt rank; 65247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 65347c6ae99SBarry Smith 65447c6ae99SBarry Smith PetscFunctionBegin; 65547c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 65647c6ae99SBarry Smith next = com->next; 65747c6ae99SBarry Smith if (com->setup) { 65847c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add an array once you have used the DMComposite"); 65947c6ae99SBarry Smith } 66047c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG) 66147c6ae99SBarry Smith { 66247c6ae99SBarry Smith PetscMPIInt orankmax; 66347c6ae99SBarry Smith ierr = MPI_Allreduce(&orank,&orankmax,1,MPI_INT,MPI_MAX,((PetscObject)dm)->comm);CHKERRQ(ierr); 66447c6ae99SBarry Smith if (orank != orankmax) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"orank %d must be equal on all processes, another process has value %d",orank,orankmax); 66547c6ae99SBarry Smith } 66647c6ae99SBarry Smith #endif 66747c6ae99SBarry Smith 66847c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 66947c6ae99SBarry Smith /* create new link */ 67047c6ae99SBarry Smith ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr); 67147c6ae99SBarry Smith mine->n = n; 67247c6ae99SBarry Smith mine->rank = orank; 67347c6ae99SBarry Smith mine->dm = PETSC_NULL; 67447c6ae99SBarry Smith mine->type = DMCOMPOSITE_ARRAY; 67547c6ae99SBarry Smith mine->next = PETSC_NULL; 67647c6ae99SBarry Smith if (rank == mine->rank) {com->n += n;com->nmine++;} 67747c6ae99SBarry Smith 67847c6ae99SBarry Smith /* add to end of list */ 67947c6ae99SBarry Smith if (!next) { 68047c6ae99SBarry Smith com->next = mine; 68147c6ae99SBarry Smith } else { 68247c6ae99SBarry Smith while (next->next) next = next->next; 68347c6ae99SBarry Smith next->next = mine; 68447c6ae99SBarry Smith } 68547c6ae99SBarry Smith com->nredundant++; 68647c6ae99SBarry Smith PetscFunctionReturn(0); 68747c6ae99SBarry Smith } 68847c6ae99SBarry Smith 68947c6ae99SBarry Smith #undef __FUNCT__ 69047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddDM" 69147c6ae99SBarry Smith /*@C 692aa219208SBarry Smith DMCompositeAddDM - adds a DM vector to a DMComposite 69347c6ae99SBarry Smith 69447c6ae99SBarry Smith Collective on DMComposite 69547c6ae99SBarry Smith 69647c6ae99SBarry Smith Input Parameter: 69747c6ae99SBarry Smith + dm - the packer object 69847c6ae99SBarry Smith - dm - the DM object, if the DM is a da you will need to caste it with a (DM) 69947c6ae99SBarry Smith 70047c6ae99SBarry Smith Level: advanced 70147c6ae99SBarry Smith 7020c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(), 70347c6ae99SBarry Smith DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 70447c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 70547c6ae99SBarry Smith 70647c6ae99SBarry Smith @*/ 70747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeAddDM(DM dmc,DM dm) 70847c6ae99SBarry Smith { 70947c6ae99SBarry Smith PetscErrorCode ierr; 71047c6ae99SBarry Smith PetscInt n; 71147c6ae99SBarry Smith struct DMCompositeLink *mine,*next; 71247c6ae99SBarry Smith Vec global; 71347c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dmc->data; 71447c6ae99SBarry Smith 71547c6ae99SBarry Smith PetscFunctionBegin; 71647c6ae99SBarry Smith PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 71747c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,2); 71847c6ae99SBarry Smith next = com->next; 719aa219208SBarry Smith if (com->setup) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite"); 72047c6ae99SBarry Smith 72147c6ae99SBarry Smith /* create new link */ 72247c6ae99SBarry Smith ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr); 72347c6ae99SBarry Smith ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 72447c6ae99SBarry Smith ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr); 72547c6ae99SBarry Smith ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr); 72647c6ae99SBarry Smith ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr); 72747c6ae99SBarry Smith mine->n = n; 72847c6ae99SBarry Smith mine->dm = dm; 72947c6ae99SBarry Smith mine->type = DMCOMPOSITE_DM; 73047c6ae99SBarry Smith mine->next = PETSC_NULL; 73147c6ae99SBarry Smith com->n += n; 73247c6ae99SBarry Smith 73347c6ae99SBarry Smith /* add to end of list */ 73447c6ae99SBarry Smith if (!next) { 73547c6ae99SBarry Smith com->next = mine; 73647c6ae99SBarry Smith } else { 73747c6ae99SBarry Smith while (next->next) next = next->next; 73847c6ae99SBarry Smith next->next = mine; 73947c6ae99SBarry Smith } 74047c6ae99SBarry Smith com->nDM++; 74147c6ae99SBarry Smith com->nmine++; 74247c6ae99SBarry Smith PetscFunctionReturn(0); 74347c6ae99SBarry Smith } 74447c6ae99SBarry Smith 74547c6ae99SBarry Smith extern PetscErrorCode PETSCDM_DLLEXPORT VecView_MPI(Vec,PetscViewer); 74647c6ae99SBarry Smith EXTERN_C_BEGIN 74747c6ae99SBarry Smith #undef __FUNCT__ 74847c6ae99SBarry Smith #define __FUNCT__ "VecView_DMComposite" 74947c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT VecView_DMComposite(Vec gvec,PetscViewer viewer) 75047c6ae99SBarry Smith { 75147c6ae99SBarry Smith DM dm; 75247c6ae99SBarry Smith PetscErrorCode ierr; 75347c6ae99SBarry Smith struct DMCompositeLink *next; 75447c6ae99SBarry Smith PetscBool isdraw; 755cef07954SSatish Balay DM_Composite *com; 75647c6ae99SBarry Smith 75747c6ae99SBarry Smith PetscFunctionBegin; 75847c6ae99SBarry Smith ierr = PetscObjectQuery((PetscObject)gvec,"DMComposite",(PetscObject*)&dm);CHKERRQ(ierr); 75947c6ae99SBarry Smith if (!dm) SETERRQ(((PetscObject)gvec)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite"); 76047c6ae99SBarry Smith com = (DM_Composite*)dm->data; 76147c6ae99SBarry Smith next = com->next; 76247c6ae99SBarry Smith 76347c6ae99SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 76447c6ae99SBarry Smith if (!isdraw) { 76547c6ae99SBarry Smith /* do I really want to call this? */ 76647c6ae99SBarry Smith ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr); 76747c6ae99SBarry Smith } else { 76847c6ae99SBarry Smith PetscInt cnt = 0; 76947c6ae99SBarry Smith 77047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 77147c6ae99SBarry Smith while (next) { 77247c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 77347c6ae99SBarry Smith PetscScalar *array; 77447c6ae99SBarry Smith ierr = DMCompositeGetAccess_Array(dm,next,gvec,&array);CHKERRQ(ierr); 77547c6ae99SBarry Smith 77647c6ae99SBarry Smith /*skip it for now */ 77747c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 77847c6ae99SBarry Smith Vec vec; 77947c6ae99SBarry Smith PetscInt bs; 78047c6ae99SBarry Smith 78147c6ae99SBarry Smith ierr = DMCompositeGetAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr); 78247c6ae99SBarry Smith ierr = VecView(vec,viewer);CHKERRQ(ierr); 78347c6ae99SBarry Smith ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr); 78447c6ae99SBarry Smith ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr); 78547c6ae99SBarry Smith ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr); 78647c6ae99SBarry Smith cnt += bs; 78747c6ae99SBarry Smith } else { 78847c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 78947c6ae99SBarry Smith } 79047c6ae99SBarry Smith next = next->next; 79147c6ae99SBarry Smith } 79247c6ae99SBarry Smith ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr); 79347c6ae99SBarry Smith } 79447c6ae99SBarry Smith PetscFunctionReturn(0); 79547c6ae99SBarry Smith } 79647c6ae99SBarry Smith EXTERN_C_END 79747c6ae99SBarry Smith 79847c6ae99SBarry Smith 79947c6ae99SBarry Smith #undef __FUNCT__ 8000c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Composite" 8010c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateGlobalVector_Composite(DM dm,Vec *gvec) 80247c6ae99SBarry Smith { 80347c6ae99SBarry Smith PetscErrorCode ierr; 80447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 80547c6ae99SBarry Smith 80647c6ae99SBarry Smith PetscFunctionBegin; 80747c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 80847c6ae99SBarry Smith if (!com->setup) { 809d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 81047c6ae99SBarry Smith } 81147c6ae99SBarry Smith ierr = VecCreateMPI(((PetscObject)dm)->comm,com->n,com->N,gvec);CHKERRQ(ierr); 81247c6ae99SBarry Smith ierr = PetscObjectCompose((PetscObject)*gvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr); 81347c6ae99SBarry Smith ierr = VecSetOperation(*gvec,VECOP_VIEW,(void(*)(void))VecView_DMComposite);CHKERRQ(ierr); 81447c6ae99SBarry Smith PetscFunctionReturn(0); 81547c6ae99SBarry Smith } 81647c6ae99SBarry Smith 81747c6ae99SBarry Smith #undef __FUNCT__ 8180c010503SBarry Smith #define __FUNCT__ "DMCreateLocalVector_Composite" 8190c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateLocalVector_Composite(DM dm,Vec *lvec) 82047c6ae99SBarry Smith { 82147c6ae99SBarry Smith PetscErrorCode ierr; 82247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 82347c6ae99SBarry Smith 82447c6ae99SBarry Smith PetscFunctionBegin; 82547c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 82647c6ae99SBarry Smith if (!com->setup) { 827d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 82847c6ae99SBarry Smith } 82947c6ae99SBarry Smith ierr = VecCreateSeq(((PetscObject)dm)->comm,com->nghost,lvec);CHKERRQ(ierr); 83047c6ae99SBarry Smith ierr = PetscObjectCompose((PetscObject)*lvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr); 83147c6ae99SBarry Smith PetscFunctionReturn(0); 83247c6ae99SBarry Smith } 83347c6ae99SBarry Smith 83447c6ae99SBarry Smith #undef __FUNCT__ 83547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalISs" 83647c6ae99SBarry Smith /*@C 83747c6ae99SBarry Smith DMCompositeGetLocalISs - gets an IS for each DM/array in the DMComposite, include ghost points 83847c6ae99SBarry Smith 83947c6ae99SBarry Smith Collective on DMComposite 84047c6ae99SBarry Smith 84147c6ae99SBarry Smith Input Parameter: 84247c6ae99SBarry Smith . dm - the packer object 84347c6ae99SBarry Smith 84447c6ae99SBarry Smith Output Parameters: 84547c6ae99SBarry Smith . is - the individual indices for each packed vector/array. Note that this includes 846aa219208SBarry Smith all the ghost points that individual ghosted DMDA's may have. Also each process has an 84747c6ae99SBarry Smith is for EACH redundant array (not just the local redundant arrays). 84847c6ae99SBarry Smith 84947c6ae99SBarry Smith Level: advanced 85047c6ae99SBarry Smith 85147c6ae99SBarry Smith Notes: 85247c6ae99SBarry Smith The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree() 85347c6ae99SBarry Smith 85447c6ae99SBarry Smith Use DMCompositeGetGlobalISs() for non-ghosted ISs. 85547c6ae99SBarry Smith 8560c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 85747c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(), 85847c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries() 85947c6ae99SBarry Smith 86047c6ae99SBarry Smith @*/ 86147c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetGlobalIndices(DM dm,IS *is[]) 86247c6ae99SBarry Smith { 86347c6ae99SBarry Smith PetscErrorCode ierr; 86447c6ae99SBarry Smith PetscInt i,*idx,n,cnt; 86547c6ae99SBarry Smith struct DMCompositeLink *next; 86647c6ae99SBarry Smith Vec global,dglobal; 86747c6ae99SBarry Smith PF pf; 86847c6ae99SBarry Smith PetscScalar *array; 86947c6ae99SBarry Smith PetscMPIInt rank; 87047c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 87147c6ae99SBarry Smith 87247c6ae99SBarry Smith PetscFunctionBegin; 87347c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 87447c6ae99SBarry Smith ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr); 87547c6ae99SBarry Smith next = com->next; 8760c010503SBarry Smith ierr = DMCreateGlobalVector(dm,&global);CHKERRQ(ierr); 87747c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 87847c6ae99SBarry Smith 87947c6ae99SBarry Smith /* put 0 to N-1 into the global vector */ 88047c6ae99SBarry Smith ierr = PFCreate(PETSC_COMM_WORLD,1,1,&pf);CHKERRQ(ierr); 88147c6ae99SBarry Smith ierr = PFSetType(pf,PFIDENTITY,PETSC_NULL);CHKERRQ(ierr); 88247c6ae99SBarry Smith ierr = PFApplyVec(pf,PETSC_NULL,global);CHKERRQ(ierr); 88347c6ae99SBarry Smith ierr = PFDestroy(pf);CHKERRQ(ierr); 88447c6ae99SBarry Smith 88547c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 88647c6ae99SBarry Smith cnt = 0; 88747c6ae99SBarry Smith while (next) { 88847c6ae99SBarry Smith 88947c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 89047c6ae99SBarry Smith 89147c6ae99SBarry Smith ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr); 89247c6ae99SBarry Smith if (rank == next->rank) { 89347c6ae99SBarry Smith ierr = VecGetArray(global,&array);CHKERRQ(ierr); 89447c6ae99SBarry Smith array += next->rstart; 89547c6ae99SBarry Smith for (i=0; i<next->n; i++) idx[i] = (PetscInt)PetscRealPart(array[i]); 89647c6ae99SBarry Smith array -= next->rstart; 89747c6ae99SBarry Smith ierr = VecRestoreArray(global,&array);CHKERRQ(ierr); 89847c6ae99SBarry Smith } 89947c6ae99SBarry Smith ierr = MPI_Bcast(idx,next->n,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr); 90047c6ae99SBarry Smith ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr); 90147c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 90247c6ae99SBarry Smith Vec local; 90347c6ae99SBarry Smith 90447c6ae99SBarry Smith ierr = DMCreateLocalVector(next->dm,&local);CHKERRQ(ierr); 90547c6ae99SBarry Smith ierr = VecGetArray(global,&array);CHKERRQ(ierr); 90647c6ae99SBarry Smith array += next->rstart; 90747c6ae99SBarry Smith ierr = DMGetGlobalVector(next->dm,&dglobal);CHKERRQ(ierr); 90847c6ae99SBarry Smith ierr = VecPlaceArray(dglobal,array);CHKERRQ(ierr); 90947c6ae99SBarry Smith ierr = DMGlobalToLocalBegin(next->dm,dglobal,INSERT_VALUES,local);CHKERRQ(ierr); 91047c6ae99SBarry Smith ierr = DMGlobalToLocalEnd(next->dm,dglobal,INSERT_VALUES,local);CHKERRQ(ierr); 91147c6ae99SBarry Smith array -= next->rstart; 91247c6ae99SBarry Smith ierr = VecRestoreArray(global,&array);CHKERRQ(ierr); 91347c6ae99SBarry Smith ierr = VecResetArray(dglobal);CHKERRQ(ierr); 91447c6ae99SBarry Smith ierr = DMRestoreGlobalVector(next->dm,&dglobal);CHKERRQ(ierr); 91547c6ae99SBarry Smith 91647c6ae99SBarry Smith ierr = VecGetArray(local,&array);CHKERRQ(ierr); 91747c6ae99SBarry Smith ierr = VecGetSize(local,&n);CHKERRQ(ierr); 91847c6ae99SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr); 91947c6ae99SBarry Smith for (i=0; i<n; i++) idx[i] = (PetscInt)PetscRealPart(array[i]); 92047c6ae99SBarry Smith ierr = VecRestoreArray(local,&array);CHKERRQ(ierr); 92147c6ae99SBarry Smith ierr = VecDestroy(local);CHKERRQ(ierr); 92247c6ae99SBarry Smith ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr); 92347c6ae99SBarry Smith 92447c6ae99SBarry Smith } else SETERRQ(((PetscObject)global)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 92547c6ae99SBarry Smith next = next->next; 92647c6ae99SBarry Smith cnt++; 92747c6ae99SBarry Smith } 92847c6ae99SBarry Smith ierr = VecDestroy(global);CHKERRQ(ierr); 92947c6ae99SBarry Smith PetscFunctionReturn(0); 93047c6ae99SBarry Smith } 93147c6ae99SBarry Smith 93247c6ae99SBarry Smith #undef __FUNCT__ 93347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetGlobalISs" 93447c6ae99SBarry Smith /*@C 93547c6ae99SBarry Smith DMCompositeGetGlobalISs - Gets the index sets for each composed object 93647c6ae99SBarry Smith 93747c6ae99SBarry Smith Collective on DMComposite 93847c6ae99SBarry Smith 93947c6ae99SBarry Smith Input Parameter: 94047c6ae99SBarry Smith . dm - the packer object 94147c6ae99SBarry Smith 94247c6ae99SBarry Smith Output Parameters: 94347c6ae99SBarry Smith . is - the array of index sets 94447c6ae99SBarry Smith 94547c6ae99SBarry Smith Level: advanced 94647c6ae99SBarry Smith 94747c6ae99SBarry Smith Notes: 94847c6ae99SBarry Smith The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree() 94947c6ae99SBarry Smith 95047c6ae99SBarry Smith The number of IS on each process will/may be different when redundant arrays are used 95147c6ae99SBarry Smith 95247c6ae99SBarry Smith These could be used to extract a subset of vector entries for a "multi-physics" preconditioner 95347c6ae99SBarry Smith 95447c6ae99SBarry Smith Use DMCompositeGetLocalISs() for index sets that include ghost points 95547c6ae99SBarry Smith 9560c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 95747c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(), 95847c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries() 95947c6ae99SBarry Smith 96047c6ae99SBarry Smith @*/ 96147c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetGlobalISs(DM dm,IS *is[]) 96247c6ae99SBarry Smith { 96347c6ae99SBarry Smith PetscErrorCode ierr; 96447c6ae99SBarry Smith PetscInt cnt = 0,*idx,i; 96547c6ae99SBarry Smith struct DMCompositeLink *next; 96647c6ae99SBarry Smith PetscMPIInt rank; 96747c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 96847c6ae99SBarry Smith 96947c6ae99SBarry Smith PetscFunctionBegin; 97047c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 97147c6ae99SBarry Smith ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr); 97247c6ae99SBarry Smith next = com->next; 97347c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 97447c6ae99SBarry Smith 97547c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 97647c6ae99SBarry Smith while (next) { 97747c6ae99SBarry Smith 97847c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 97947c6ae99SBarry Smith 98047c6ae99SBarry Smith if (rank == next->rank) { 98147c6ae99SBarry Smith ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr); 98247c6ae99SBarry Smith for (i=0; i<next->n; i++) idx[i] = next->grstart + i; 98347c6ae99SBarry Smith ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr); 98447c6ae99SBarry Smith cnt++; 98547c6ae99SBarry Smith } 98647c6ae99SBarry Smith 98747c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 98847c6ae99SBarry Smith ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr); 98947c6ae99SBarry Smith for (i=0; i<next->n; i++) idx[i] = next->grstart + i; 99047c6ae99SBarry Smith ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr); 99147c6ae99SBarry Smith cnt++; 99247c6ae99SBarry Smith } else { 99347c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 99447c6ae99SBarry Smith } 99547c6ae99SBarry Smith next = next->next; 99647c6ae99SBarry Smith } 99747c6ae99SBarry Smith PetscFunctionReturn(0); 99847c6ae99SBarry Smith } 99947c6ae99SBarry Smith 100047c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 100147c6ae99SBarry Smith #undef __FUNCT__ 100247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_Array" 100347c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array) 100447c6ae99SBarry Smith { 100547c6ae99SBarry Smith PetscErrorCode ierr; 100647c6ae99SBarry Smith PetscFunctionBegin; 100747c6ae99SBarry Smith if (array) { 100847c6ae99SBarry Smith ierr = PetscMalloc(mine->n*sizeof(PetscScalar),array);CHKERRQ(ierr); 100947c6ae99SBarry Smith } 101047c6ae99SBarry Smith PetscFunctionReturn(0); 101147c6ae99SBarry Smith } 101247c6ae99SBarry Smith 101347c6ae99SBarry Smith #undef __FUNCT__ 101447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_DM" 101547c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local) 101647c6ae99SBarry Smith { 101747c6ae99SBarry Smith PetscErrorCode ierr; 101847c6ae99SBarry Smith PetscFunctionBegin; 101947c6ae99SBarry Smith if (local) { 102047c6ae99SBarry Smith ierr = DMGetLocalVector(mine->dm,local);CHKERRQ(ierr); 102147c6ae99SBarry Smith } 102247c6ae99SBarry Smith PetscFunctionReturn(0); 102347c6ae99SBarry Smith } 102447c6ae99SBarry Smith 102547c6ae99SBarry Smith #undef __FUNCT__ 102647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_Array" 102747c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array) 102847c6ae99SBarry Smith { 102947c6ae99SBarry Smith PetscErrorCode ierr; 103047c6ae99SBarry Smith PetscFunctionBegin; 103147c6ae99SBarry Smith if (array) { 103247c6ae99SBarry Smith ierr = PetscFree(*array);CHKERRQ(ierr); 103347c6ae99SBarry Smith } 103447c6ae99SBarry Smith PetscFunctionReturn(0); 103547c6ae99SBarry Smith } 103647c6ae99SBarry Smith 103747c6ae99SBarry Smith #undef __FUNCT__ 103847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_DM" 103947c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local) 104047c6ae99SBarry Smith { 104147c6ae99SBarry Smith PetscErrorCode ierr; 104247c6ae99SBarry Smith PetscFunctionBegin; 104347c6ae99SBarry Smith if (local) { 104447c6ae99SBarry Smith ierr = DMRestoreLocalVector(mine->dm,local);CHKERRQ(ierr); 104547c6ae99SBarry Smith } 104647c6ae99SBarry Smith PetscFunctionReturn(0); 104747c6ae99SBarry Smith } 104847c6ae99SBarry Smith 104947c6ae99SBarry Smith #undef __FUNCT__ 105047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors" 105147c6ae99SBarry Smith /*@C 105247c6ae99SBarry Smith DMCompositeGetLocalVectors - Gets local vectors and arrays for each part of a DMComposite. 105347c6ae99SBarry Smith Use DMCompositeRestoreLocalVectors() to return them. 105447c6ae99SBarry Smith 105547c6ae99SBarry Smith Not Collective 105647c6ae99SBarry Smith 105747c6ae99SBarry Smith Input Parameter: 105847c6ae99SBarry Smith . dm - the packer object 105947c6ae99SBarry Smith 106047c6ae99SBarry Smith Output Parameter: 106147c6ae99SBarry Smith . ... - the individual sequential objects (arrays or vectors) 106247c6ae99SBarry Smith 106347c6ae99SBarry Smith Level: advanced 106447c6ae99SBarry Smith 10650c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 106647c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 106747c6ae99SBarry Smith DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries() 106847c6ae99SBarry Smith 106947c6ae99SBarry Smith @*/ 107047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetLocalVectors(DM dm,...) 107147c6ae99SBarry Smith { 107247c6ae99SBarry Smith va_list Argp; 107347c6ae99SBarry Smith PetscErrorCode ierr; 107447c6ae99SBarry Smith struct DMCompositeLink *next; 107547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 107647c6ae99SBarry Smith 107747c6ae99SBarry Smith PetscFunctionBegin; 107847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 107947c6ae99SBarry Smith next = com->next; 108047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 108147c6ae99SBarry Smith va_start(Argp,dm); 108247c6ae99SBarry Smith while (next) { 108347c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 108447c6ae99SBarry Smith PetscScalar **array; 108547c6ae99SBarry Smith array = va_arg(Argp, PetscScalar**); 108647c6ae99SBarry Smith ierr = DMCompositeGetLocalVectors_Array(dm,next,array);CHKERRQ(ierr); 108747c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 108847c6ae99SBarry Smith Vec *vec; 108947c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 109047c6ae99SBarry Smith ierr = DMCompositeGetLocalVectors_DM(dm,next,vec);CHKERRQ(ierr); 109147c6ae99SBarry Smith } else { 109247c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 109347c6ae99SBarry Smith } 109447c6ae99SBarry Smith next = next->next; 109547c6ae99SBarry Smith } 109647c6ae99SBarry Smith va_end(Argp); 109747c6ae99SBarry Smith PetscFunctionReturn(0); 109847c6ae99SBarry Smith } 109947c6ae99SBarry Smith 110047c6ae99SBarry Smith #undef __FUNCT__ 110147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors" 110247c6ae99SBarry Smith /*@C 110347c6ae99SBarry Smith DMCompositeRestoreLocalVectors - Restores local vectors and arrays for each part of a DMComposite. 110447c6ae99SBarry Smith 110547c6ae99SBarry Smith Not Collective 110647c6ae99SBarry Smith 110747c6ae99SBarry Smith Input Parameter: 110847c6ae99SBarry Smith . dm - the packer object 110947c6ae99SBarry Smith 111047c6ae99SBarry Smith Output Parameter: 111147c6ae99SBarry Smith . ... - the individual sequential objects (arrays or vectors) 111247c6ae99SBarry Smith 111347c6ae99SBarry Smith Level: advanced 111447c6ae99SBarry Smith 11150c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 111647c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 111747c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries() 111847c6ae99SBarry Smith 111947c6ae99SBarry Smith @*/ 112047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreLocalVectors(DM dm,...) 112147c6ae99SBarry Smith { 112247c6ae99SBarry Smith va_list Argp; 112347c6ae99SBarry Smith PetscErrorCode ierr; 112447c6ae99SBarry Smith struct DMCompositeLink *next; 112547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 112647c6ae99SBarry Smith 112747c6ae99SBarry Smith PetscFunctionBegin; 112847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 112947c6ae99SBarry Smith next = com->next; 113047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 113147c6ae99SBarry Smith va_start(Argp,dm); 113247c6ae99SBarry Smith while (next) { 113347c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 113447c6ae99SBarry Smith PetscScalar **array; 113547c6ae99SBarry Smith array = va_arg(Argp, PetscScalar**); 113647c6ae99SBarry Smith ierr = DMCompositeRestoreLocalVectors_Array(dm,next,array);CHKERRQ(ierr); 113747c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 113847c6ae99SBarry Smith Vec *vec; 113947c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 114047c6ae99SBarry Smith ierr = DMCompositeRestoreLocalVectors_DM(dm,next,vec);CHKERRQ(ierr); 114147c6ae99SBarry Smith } else { 114247c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 114347c6ae99SBarry Smith } 114447c6ae99SBarry Smith next = next->next; 114547c6ae99SBarry Smith } 114647c6ae99SBarry Smith va_end(Argp); 114747c6ae99SBarry Smith PetscFunctionReturn(0); 114847c6ae99SBarry Smith } 114947c6ae99SBarry Smith 115047c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 115147c6ae99SBarry Smith #undef __FUNCT__ 115247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_Array" 115347c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_Array(DM dm,struct DMCompositeLink *mine,PetscInt *n) 115447c6ae99SBarry Smith { 115547c6ae99SBarry Smith PetscFunctionBegin; 115647c6ae99SBarry Smith if (n) *n = mine->n; 115747c6ae99SBarry Smith PetscFunctionReturn(0); 115847c6ae99SBarry Smith } 115947c6ae99SBarry Smith 116047c6ae99SBarry Smith #undef __FUNCT__ 116147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_DM" 116247c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_DM(DM dmi,struct DMCompositeLink *mine,DM *dm) 116347c6ae99SBarry Smith { 116447c6ae99SBarry Smith PetscFunctionBegin; 116547c6ae99SBarry Smith if (dm) *dm = mine->dm; 116647c6ae99SBarry Smith PetscFunctionReturn(0); 116747c6ae99SBarry Smith } 116847c6ae99SBarry Smith 116947c6ae99SBarry Smith #undef __FUNCT__ 117047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries" 117147c6ae99SBarry Smith /*@C 1172aa219208SBarry Smith DMCompositeGetEntries - Gets the DM, redundant size, etc for each entry in a DMComposite. 117347c6ae99SBarry Smith 117447c6ae99SBarry Smith Not Collective 117547c6ae99SBarry Smith 117647c6ae99SBarry Smith Input Parameter: 117747c6ae99SBarry Smith . dm - the packer object 117847c6ae99SBarry Smith 117947c6ae99SBarry Smith Output Parameter: 1180aa219208SBarry Smith . ... - the individual entries, DMs or integer sizes) 118147c6ae99SBarry Smith 118247c6ae99SBarry Smith Level: advanced 118347c6ae99SBarry Smith 11840c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 118547c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 118647c6ae99SBarry Smith DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(), DMCompositeScatter(), 118747c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors() 118847c6ae99SBarry Smith 118947c6ae99SBarry Smith @*/ 119047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetEntries(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; 119647c6ae99SBarry Smith 119747c6ae99SBarry Smith PetscFunctionBegin; 119847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 119947c6ae99SBarry Smith next = com->next; 120047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 120147c6ae99SBarry Smith va_start(Argp,dm); 120247c6ae99SBarry Smith while (next) { 120347c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 120447c6ae99SBarry Smith PetscInt *n; 120547c6ae99SBarry Smith n = va_arg(Argp, PetscInt*); 120647c6ae99SBarry Smith ierr = DMCompositeGetEntries_Array(dm,next,n);CHKERRQ(ierr); 120747c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 120847c6ae99SBarry Smith DM *dmn; 120947c6ae99SBarry Smith dmn = va_arg(Argp, DM*); 121047c6ae99SBarry Smith ierr = DMCompositeGetEntries_DM(dm,next,dmn);CHKERRQ(ierr); 121147c6ae99SBarry Smith } else { 121247c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 121347c6ae99SBarry Smith } 121447c6ae99SBarry Smith next = next->next; 121547c6ae99SBarry Smith } 121647c6ae99SBarry Smith va_end(Argp); 121747c6ae99SBarry Smith PetscFunctionReturn(0); 121847c6ae99SBarry Smith } 121947c6ae99SBarry Smith 122047c6ae99SBarry Smith #undef __FUNCT__ 12210c010503SBarry Smith #define __FUNCT__ "DMRefine_Composite" 12220c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine) 122347c6ae99SBarry Smith { 122447c6ae99SBarry Smith PetscErrorCode ierr; 122547c6ae99SBarry Smith struct DMCompositeLink *next; 122647c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dmi->data; 122747c6ae99SBarry Smith DM dm; 122847c6ae99SBarry Smith 122947c6ae99SBarry Smith PetscFunctionBegin; 123047c6ae99SBarry Smith PetscValidHeaderSpecific(dmi,DM_CLASSID,1); 123147c6ae99SBarry Smith next = com->next; 123247c6ae99SBarry Smith ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr); 123347c6ae99SBarry Smith 123447c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 123547c6ae99SBarry Smith while (next) { 123647c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 123747c6ae99SBarry Smith ierr = DMCompositeAddArray(*fine,next->rank,next->n);CHKERRQ(ierr); 123847c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 123947c6ae99SBarry Smith ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr); 124047c6ae99SBarry Smith ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr); 124147c6ae99SBarry Smith ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 124247c6ae99SBarry Smith } else { 124347c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 124447c6ae99SBarry Smith } 124547c6ae99SBarry Smith next = next->next; 124647c6ae99SBarry Smith } 124747c6ae99SBarry Smith PetscFunctionReturn(0); 124847c6ae99SBarry Smith } 124947c6ae99SBarry Smith 125047c6ae99SBarry Smith #include "petscmat.h" 125147c6ae99SBarry Smith 125247c6ae99SBarry Smith struct MatPackLink { 125347c6ae99SBarry Smith Mat A; 125447c6ae99SBarry Smith struct MatPackLink *next; 125547c6ae99SBarry Smith }; 125647c6ae99SBarry Smith 125747c6ae99SBarry Smith struct MatPack { 125847c6ae99SBarry Smith DM right,left; 125947c6ae99SBarry Smith struct MatPackLink *next; 126047c6ae99SBarry Smith }; 126147c6ae99SBarry Smith 126247c6ae99SBarry Smith #undef __FUNCT__ 126347c6ae99SBarry Smith #define __FUNCT__ "MatMultBoth_Shell_Pack" 126447c6ae99SBarry Smith PetscErrorCode MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscBool add) 126547c6ae99SBarry Smith { 126647c6ae99SBarry Smith struct MatPack *mpack; 126747c6ae99SBarry Smith struct DMCompositeLink *xnext,*ynext; 126847c6ae99SBarry Smith struct MatPackLink *anext; 126947c6ae99SBarry Smith PetscScalar *xarray,*yarray; 127047c6ae99SBarry Smith PetscErrorCode ierr; 127147c6ae99SBarry Smith PetscInt i; 127247c6ae99SBarry Smith Vec xglobal,yglobal; 127347c6ae99SBarry Smith PetscMPIInt rank; 127447c6ae99SBarry Smith DM_Composite *comright; 127547c6ae99SBarry Smith DM_Composite *comleft; 127647c6ae99SBarry Smith 127747c6ae99SBarry Smith PetscFunctionBegin; 127847c6ae99SBarry Smith ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr); 127947c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr); 128047c6ae99SBarry Smith comright = (DM_Composite*)mpack->right->data; 128147c6ae99SBarry Smith comleft = (DM_Composite*)mpack->left->data; 128247c6ae99SBarry Smith xnext = comright->next; 128347c6ae99SBarry Smith ynext = comleft->next; 128447c6ae99SBarry Smith anext = mpack->next; 128547c6ae99SBarry Smith 128647c6ae99SBarry Smith while (xnext) { 128747c6ae99SBarry Smith if (xnext->type == DMCOMPOSITE_ARRAY) { 128847c6ae99SBarry Smith if (rank == xnext->rank) { 128947c6ae99SBarry Smith ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 129047c6ae99SBarry Smith ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 129147c6ae99SBarry Smith if (add) { 129247c6ae99SBarry Smith for (i=0; i<xnext->n; i++) { 129347c6ae99SBarry Smith yarray[ynext->rstart+i] += xarray[xnext->rstart+i]; 129447c6ae99SBarry Smith } 129547c6ae99SBarry Smith } else { 129647c6ae99SBarry Smith ierr = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr); 129747c6ae99SBarry Smith } 129847c6ae99SBarry Smith ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 129947c6ae99SBarry Smith ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 130047c6ae99SBarry Smith } 130147c6ae99SBarry Smith } else if (xnext->type == DMCOMPOSITE_DM) { 130247c6ae99SBarry Smith ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 130347c6ae99SBarry Smith ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 130447c6ae99SBarry Smith ierr = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr); 130547c6ae99SBarry Smith ierr = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr); 130647c6ae99SBarry Smith ierr = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr); 130747c6ae99SBarry Smith ierr = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr); 130847c6ae99SBarry Smith if (add) { 130947c6ae99SBarry Smith ierr = MatMultAdd(anext->A,xglobal,yglobal,yglobal);CHKERRQ(ierr); 131047c6ae99SBarry Smith } else { 131147c6ae99SBarry Smith ierr = MatMult(anext->A,xglobal,yglobal);CHKERRQ(ierr); 131247c6ae99SBarry Smith } 131347c6ae99SBarry Smith ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 131447c6ae99SBarry Smith ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 131547c6ae99SBarry Smith ierr = VecResetArray(xglobal);CHKERRQ(ierr); 131647c6ae99SBarry Smith ierr = VecResetArray(yglobal);CHKERRQ(ierr); 131747c6ae99SBarry Smith ierr = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr); 131847c6ae99SBarry Smith ierr = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr); 131947c6ae99SBarry Smith anext = anext->next; 132047c6ae99SBarry Smith } else { 132147c6ae99SBarry Smith SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 132247c6ae99SBarry Smith } 132347c6ae99SBarry Smith xnext = xnext->next; 132447c6ae99SBarry Smith ynext = ynext->next; 132547c6ae99SBarry Smith } 132647c6ae99SBarry Smith PetscFunctionReturn(0); 132747c6ae99SBarry Smith } 132847c6ae99SBarry Smith 132947c6ae99SBarry Smith #undef __FUNCT__ 133047c6ae99SBarry Smith #define __FUNCT__ "MatMultAdd_Shell_Pack" 133147c6ae99SBarry Smith PetscErrorCode MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z) 133247c6ae99SBarry Smith { 133347c6ae99SBarry Smith PetscErrorCode ierr; 133447c6ae99SBarry Smith PetscFunctionBegin; 133547c6ae99SBarry Smith if (z != y) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Handles y == z only"); 133647c6ae99SBarry Smith ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);CHKERRQ(ierr); 133747c6ae99SBarry Smith PetscFunctionReturn(0); 133847c6ae99SBarry Smith } 133947c6ae99SBarry Smith 134047c6ae99SBarry Smith #undef __FUNCT__ 134147c6ae99SBarry Smith #define __FUNCT__ "MatMult_Shell_Pack" 134247c6ae99SBarry Smith PetscErrorCode MatMult_Shell_Pack(Mat A,Vec x,Vec y) 134347c6ae99SBarry Smith { 134447c6ae99SBarry Smith PetscErrorCode ierr; 134547c6ae99SBarry Smith PetscFunctionBegin; 134647c6ae99SBarry Smith ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);CHKERRQ(ierr); 134747c6ae99SBarry Smith PetscFunctionReturn(0); 134847c6ae99SBarry Smith } 134947c6ae99SBarry Smith 135047c6ae99SBarry Smith #undef __FUNCT__ 135147c6ae99SBarry Smith #define __FUNCT__ "MatMultTranspose_Shell_Pack" 135247c6ae99SBarry Smith PetscErrorCode MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y) 135347c6ae99SBarry Smith { 135447c6ae99SBarry Smith struct MatPack *mpack; 135547c6ae99SBarry Smith struct DMCompositeLink *xnext,*ynext; 135647c6ae99SBarry Smith struct MatPackLink *anext; 135747c6ae99SBarry Smith PetscScalar *xarray,*yarray; 135847c6ae99SBarry Smith PetscErrorCode ierr; 135947c6ae99SBarry Smith Vec xglobal,yglobal; 136047c6ae99SBarry Smith PetscMPIInt rank; 136147c6ae99SBarry Smith DM_Composite *comright; 136247c6ae99SBarry Smith DM_Composite *comleft; 136347c6ae99SBarry Smith 136447c6ae99SBarry Smith PetscFunctionBegin; 136547c6ae99SBarry Smith ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr); 136647c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr); 136747c6ae99SBarry Smith comright = (DM_Composite*)mpack->right->data; 136847c6ae99SBarry Smith comleft = (DM_Composite*)mpack->left->data; 136947c6ae99SBarry Smith ynext = comright->next; 137047c6ae99SBarry Smith xnext = comleft->next; 137147c6ae99SBarry Smith anext = mpack->next; 137247c6ae99SBarry Smith 137347c6ae99SBarry Smith while (xnext) { 137447c6ae99SBarry Smith if (xnext->type == DMCOMPOSITE_ARRAY) { 137547c6ae99SBarry Smith if (rank == ynext->rank) { 137647c6ae99SBarry Smith ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 137747c6ae99SBarry Smith ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 137847c6ae99SBarry Smith ierr = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr); 137947c6ae99SBarry Smith ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 138047c6ae99SBarry Smith ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 138147c6ae99SBarry Smith } 138247c6ae99SBarry Smith } else if (xnext->type == DMCOMPOSITE_DM) { 138347c6ae99SBarry Smith ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 138447c6ae99SBarry Smith ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 138547c6ae99SBarry Smith ierr = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr); 138647c6ae99SBarry Smith ierr = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr); 138747c6ae99SBarry Smith ierr = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr); 138847c6ae99SBarry Smith ierr = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr); 138947c6ae99SBarry Smith ierr = MatMultTranspose(anext->A,xglobal,yglobal);CHKERRQ(ierr); 139047c6ae99SBarry Smith ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 139147c6ae99SBarry Smith ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 139247c6ae99SBarry Smith ierr = VecResetArray(xglobal);CHKERRQ(ierr); 139347c6ae99SBarry Smith ierr = VecResetArray(yglobal);CHKERRQ(ierr); 139447c6ae99SBarry Smith ierr = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr); 139547c6ae99SBarry Smith ierr = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr); 139647c6ae99SBarry Smith anext = anext->next; 139747c6ae99SBarry Smith } else { 139847c6ae99SBarry Smith SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 139947c6ae99SBarry Smith } 140047c6ae99SBarry Smith xnext = xnext->next; 140147c6ae99SBarry Smith ynext = ynext->next; 140247c6ae99SBarry Smith } 140347c6ae99SBarry Smith PetscFunctionReturn(0); 140447c6ae99SBarry Smith } 140547c6ae99SBarry Smith 140647c6ae99SBarry Smith #undef __FUNCT__ 140747c6ae99SBarry Smith #define __FUNCT__ "MatDestroy_Shell_Pack" 140847c6ae99SBarry Smith PetscErrorCode MatDestroy_Shell_Pack(Mat A) 140947c6ae99SBarry Smith { 141047c6ae99SBarry Smith struct MatPack *mpack; 141147c6ae99SBarry Smith struct MatPackLink *anext,*oldanext; 141247c6ae99SBarry Smith PetscErrorCode ierr; 141347c6ae99SBarry Smith 141447c6ae99SBarry Smith PetscFunctionBegin; 141547c6ae99SBarry Smith ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr); 141647c6ae99SBarry Smith anext = mpack->next; 141747c6ae99SBarry Smith 141847c6ae99SBarry Smith while (anext) { 141947c6ae99SBarry Smith ierr = MatDestroy(anext->A);CHKERRQ(ierr); 142047c6ae99SBarry Smith oldanext = anext; 142147c6ae99SBarry Smith anext = anext->next; 142247c6ae99SBarry Smith ierr = PetscFree(oldanext);CHKERRQ(ierr); 142347c6ae99SBarry Smith } 142447c6ae99SBarry Smith ierr = PetscFree(mpack);CHKERRQ(ierr); 142547c6ae99SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 142647c6ae99SBarry Smith PetscFunctionReturn(0); 142747c6ae99SBarry Smith } 142847c6ae99SBarry Smith 142947c6ae99SBarry Smith #undef __FUNCT__ 14300c010503SBarry Smith #define __FUNCT__ "DMGetInterpolation_Composite" 14310c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v) 143247c6ae99SBarry Smith { 143347c6ae99SBarry Smith PetscErrorCode ierr; 143447c6ae99SBarry Smith PetscInt m,n,M,N; 143547c6ae99SBarry Smith struct DMCompositeLink *nextc; 143647c6ae99SBarry Smith struct DMCompositeLink *nextf; 143747c6ae99SBarry Smith struct MatPackLink *nextmat,*pnextmat = 0; 143847c6ae99SBarry Smith struct MatPack *mpack; 143947c6ae99SBarry Smith Vec gcoarse,gfine; 144047c6ae99SBarry Smith DM_Composite *comcoarse = (DM_Composite*)coarse->data; 144147c6ae99SBarry Smith DM_Composite *comfine = (DM_Composite*)fine->data; 144247c6ae99SBarry Smith 144347c6ae99SBarry Smith PetscFunctionBegin; 144447c6ae99SBarry Smith PetscValidHeaderSpecific(coarse,DM_CLASSID,1); 144547c6ae99SBarry Smith PetscValidHeaderSpecific(fine,DM_CLASSID,2); 144647c6ae99SBarry Smith nextc = comcoarse->next; 144747c6ae99SBarry Smith nextf = comfine->next; 144847c6ae99SBarry Smith /* use global vectors only for determining matrix layout */ 14490c010503SBarry Smith ierr = DMCreateGlobalVector(coarse,&gcoarse);CHKERRQ(ierr); 14500c010503SBarry Smith ierr = DMCreateGlobalVector(fine,&gfine);CHKERRQ(ierr); 145147c6ae99SBarry Smith ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr); 145247c6ae99SBarry Smith ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr); 145347c6ae99SBarry Smith ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr); 145447c6ae99SBarry Smith ierr = VecGetSize(gfine,&M);CHKERRQ(ierr); 145547c6ae99SBarry Smith ierr = VecDestroy(gcoarse);CHKERRQ(ierr); 145647c6ae99SBarry Smith ierr = VecDestroy(gfine);CHKERRQ(ierr); 145747c6ae99SBarry Smith 145847c6ae99SBarry Smith ierr = PetscNew(struct MatPack,&mpack);CHKERRQ(ierr); 145947c6ae99SBarry Smith mpack->right = coarse; 146047c6ae99SBarry Smith mpack->left = fine; 146147c6ae99SBarry Smith ierr = MatCreate(((PetscObject)fine)->comm,A);CHKERRQ(ierr); 146247c6ae99SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 146347c6ae99SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 146447c6ae99SBarry Smith ierr = MatShellSetContext(*A,mpack);CHKERRQ(ierr); 146547c6ae99SBarry Smith ierr = MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);CHKERRQ(ierr); 146647c6ae99SBarry Smith ierr = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);CHKERRQ(ierr); 146747c6ae99SBarry Smith ierr = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);CHKERRQ(ierr); 146847c6ae99SBarry Smith ierr = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);CHKERRQ(ierr); 146947c6ae99SBarry Smith 147047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 147147c6ae99SBarry Smith while (nextc) { 147247c6ae99SBarry Smith if (nextc->type != nextf->type) SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Two DM have different layout"); 147347c6ae99SBarry Smith 147447c6ae99SBarry Smith if (nextc->type == DMCOMPOSITE_ARRAY) { 147547c6ae99SBarry Smith ; 147647c6ae99SBarry Smith } else if (nextc->type == DMCOMPOSITE_DM) { 147747c6ae99SBarry Smith ierr = PetscNew(struct MatPackLink,&nextmat);CHKERRQ(ierr); 147847c6ae99SBarry Smith nextmat->next = 0; 147947c6ae99SBarry Smith if (pnextmat) { 148047c6ae99SBarry Smith pnextmat->next = nextmat; 148147c6ae99SBarry Smith pnextmat = nextmat; 148247c6ae99SBarry Smith } else { 148347c6ae99SBarry Smith pnextmat = nextmat; 148447c6ae99SBarry Smith mpack->next = nextmat; 148547c6ae99SBarry Smith } 148647c6ae99SBarry Smith ierr = DMGetInterpolation(nextc->dm,nextf->dm,&nextmat->A,PETSC_NULL);CHKERRQ(ierr); 148747c6ae99SBarry Smith } else { 148847c6ae99SBarry Smith SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 148947c6ae99SBarry Smith } 149047c6ae99SBarry Smith nextc = nextc->next; 149147c6ae99SBarry Smith nextf = nextf->next; 149247c6ae99SBarry Smith } 149347c6ae99SBarry Smith PetscFunctionReturn(0); 149447c6ae99SBarry Smith } 149547c6ae99SBarry Smith 149647c6ae99SBarry Smith #undef __FUNCT__ 14970c010503SBarry Smith #define __FUNCT__ "DMGetMatrix_Composite" 14980c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix_Composite(DM dm, const MatType mtype,Mat *J) 149947c6ae99SBarry Smith { 150047c6ae99SBarry Smith PetscErrorCode ierr; 150147c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 150247c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 150347c6ae99SBarry Smith PetscInt m,*dnz,*onz,i,j,mA; 150447c6ae99SBarry Smith Mat Atmp; 150547c6ae99SBarry Smith PetscMPIInt rank; 150647c6ae99SBarry Smith PetscScalar zero = 0.0; 150747c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 150847c6ae99SBarry Smith 150947c6ae99SBarry Smith PetscFunctionBegin; 151047c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 151147c6ae99SBarry Smith 151247c6ae99SBarry Smith /* use global vector to determine layout needed for matrix */ 151347c6ae99SBarry Smith m = com->n; 151447c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 151547c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr); 151647c6ae99SBarry Smith ierr = MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 151747c6ae99SBarry Smith ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 151847c6ae99SBarry Smith 151947c6ae99SBarry Smith /* 152047c6ae99SBarry Smith Extremely inefficient but will compute entire Jacobian for testing 152147c6ae99SBarry Smith */ 1522671f6225SBarry Smith ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr); 152347c6ae99SBarry Smith if (dense) { 152447c6ae99SBarry Smith PetscInt rstart,rend,*indices; 152547c6ae99SBarry Smith PetscScalar *values; 152647c6ae99SBarry Smith 152747c6ae99SBarry Smith mA = com->N; 152847c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,mA,PETSC_NULL,mA-m,PETSC_NULL);CHKERRQ(ierr); 152947c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,mA,PETSC_NULL);CHKERRQ(ierr); 153047c6ae99SBarry Smith 153147c6ae99SBarry Smith ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr); 153247c6ae99SBarry Smith ierr = PetscMalloc2(mA,PetscScalar,&values,mA,PetscInt,&indices);CHKERRQ(ierr); 153347c6ae99SBarry Smith ierr = PetscMemzero(values,mA*sizeof(PetscScalar));CHKERRQ(ierr); 153447c6ae99SBarry Smith for (i=0; i<mA; i++) indices[i] = i; 153547c6ae99SBarry Smith for (i=rstart; i<rend; i++) { 153647c6ae99SBarry Smith ierr = MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES);CHKERRQ(ierr); 153747c6ae99SBarry Smith } 153847c6ae99SBarry Smith ierr = PetscFree2(values,indices);CHKERRQ(ierr); 153947c6ae99SBarry Smith ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 154047c6ae99SBarry Smith ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 154147c6ae99SBarry Smith PetscFunctionReturn(0); 154247c6ae99SBarry Smith } 154347c6ae99SBarry Smith 154447c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)dm)->comm,m,m,dnz,onz);CHKERRQ(ierr); 154547c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 154647c6ae99SBarry Smith next = com->next; 154747c6ae99SBarry Smith while (next) { 154847c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 154947c6ae99SBarry Smith if (rank == next->rank) { /* zero the "little" block */ 155047c6ae99SBarry Smith for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) { 155147c6ae99SBarry Smith for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) { 155247c6ae99SBarry Smith ierr = MatPreallocateSet(j,1,&i,dnz,onz);CHKERRQ(ierr); 155347c6ae99SBarry Smith } 155447c6ae99SBarry Smith } 155547c6ae99SBarry Smith } 155647c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 155747c6ae99SBarry Smith PetscInt nc,rstart,*ccols,maxnc; 155847c6ae99SBarry Smith const PetscInt *cols,*rstarts; 155947c6ae99SBarry Smith PetscMPIInt proc; 156047c6ae99SBarry Smith 156147c6ae99SBarry Smith ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr); 156247c6ae99SBarry Smith ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr); 156347c6ae99SBarry Smith ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr); 156447c6ae99SBarry Smith ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr); 156547c6ae99SBarry Smith 156647c6ae99SBarry Smith maxnc = 0; 156747c6ae99SBarry Smith for (i=0; i<mA; i++) { 156847c6ae99SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 156947c6ae99SBarry Smith ierr = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 157047c6ae99SBarry Smith maxnc = PetscMax(nc,maxnc); 157147c6ae99SBarry Smith } 157247c6ae99SBarry Smith ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr); 157347c6ae99SBarry Smith for (i=0; i<mA; i++) { 157447c6ae99SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr); 157547c6ae99SBarry Smith /* remap the columns taking into how much they are shifted on each process */ 157647c6ae99SBarry Smith for (j=0; j<nc; j++) { 157747c6ae99SBarry Smith proc = 0; 157847c6ae99SBarry Smith while (cols[j] >= rstarts[proc+1]) proc++; 157947c6ae99SBarry Smith ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 158047c6ae99SBarry Smith } 158147c6ae99SBarry Smith ierr = MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz);CHKERRQ(ierr); 158247c6ae99SBarry Smith ierr = MatRestoreRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr); 158347c6ae99SBarry Smith } 158447c6ae99SBarry Smith ierr = PetscFree(ccols);CHKERRQ(ierr); 158547c6ae99SBarry Smith ierr = MatDestroy(Atmp);CHKERRQ(ierr); 158647c6ae99SBarry Smith } else { 158747c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 158847c6ae99SBarry Smith } 158947c6ae99SBarry Smith next = next->next; 159047c6ae99SBarry Smith } 159147c6ae99SBarry Smith if (com->FormCoupleLocations) { 159247c6ae99SBarry Smith ierr = (*com->FormCoupleLocations)(dm,PETSC_NULL,dnz,onz,__rstart,__nrows,__start,__end);CHKERRQ(ierr); 159347c6ae99SBarry Smith } 159447c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);CHKERRQ(ierr); 159547c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,0,dnz);CHKERRQ(ierr); 159647c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 159747c6ae99SBarry Smith 159847c6ae99SBarry Smith next = com->next; 159947c6ae99SBarry Smith while (next) { 160047c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 160147c6ae99SBarry Smith if (rank == next->rank) { 160247c6ae99SBarry Smith for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) { 160347c6ae99SBarry Smith for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) { 160447c6ae99SBarry Smith ierr = MatSetValues(*J,1,&j,1,&i,&zero,INSERT_VALUES);CHKERRQ(ierr); 160547c6ae99SBarry Smith } 160647c6ae99SBarry Smith } 160747c6ae99SBarry Smith } 160847c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 160947c6ae99SBarry Smith PetscInt nc,rstart,row,maxnc,*ccols; 161047c6ae99SBarry Smith const PetscInt *cols,*rstarts; 161147c6ae99SBarry Smith const PetscScalar *values; 161247c6ae99SBarry Smith PetscMPIInt proc; 161347c6ae99SBarry Smith 161447c6ae99SBarry Smith ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr); 161547c6ae99SBarry Smith ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr); 161647c6ae99SBarry Smith ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr); 161747c6ae99SBarry Smith ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr); 161847c6ae99SBarry Smith maxnc = 0; 161947c6ae99SBarry Smith for (i=0; i<mA; i++) { 162047c6ae99SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 162147c6ae99SBarry Smith ierr = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 162247c6ae99SBarry Smith maxnc = PetscMax(nc,maxnc); 162347c6ae99SBarry Smith } 162447c6ae99SBarry Smith ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr); 162547c6ae99SBarry Smith for (i=0; i<mA; i++) { 162647c6ae99SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr); 162747c6ae99SBarry Smith for (j=0; j<nc; j++) { 162847c6ae99SBarry Smith proc = 0; 162947c6ae99SBarry Smith while (cols[j] >= rstarts[proc+1]) proc++; 163047c6ae99SBarry Smith ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 163147c6ae99SBarry Smith } 163247c6ae99SBarry Smith row = com->rstart+next->rstart+i; 163347c6ae99SBarry Smith ierr = MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES);CHKERRQ(ierr); 163447c6ae99SBarry Smith ierr = MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr); 163547c6ae99SBarry Smith } 163647c6ae99SBarry Smith ierr = PetscFree(ccols);CHKERRQ(ierr); 163747c6ae99SBarry Smith ierr = MatDestroy(Atmp);CHKERRQ(ierr); 163847c6ae99SBarry Smith } else { 163947c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 164047c6ae99SBarry Smith } 164147c6ae99SBarry Smith next = next->next; 164247c6ae99SBarry Smith } 164347c6ae99SBarry Smith if (com->FormCoupleLocations) { 164447c6ae99SBarry Smith PetscInt __rstart; 164547c6ae99SBarry Smith ierr = MatGetOwnershipRange(*J,&__rstart,PETSC_NULL);CHKERRQ(ierr); 164647c6ae99SBarry Smith ierr = (*com->FormCoupleLocations)(dm,*J,PETSC_NULL,PETSC_NULL,__rstart,0,0,0);CHKERRQ(ierr); 164747c6ae99SBarry Smith } 164847c6ae99SBarry Smith ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164947c6ae99SBarry Smith ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165047c6ae99SBarry Smith PetscFunctionReturn(0); 165147c6ae99SBarry Smith } 165247c6ae99SBarry Smith 165347c6ae99SBarry Smith #undef __FUNCT__ 16540c010503SBarry Smith #define __FUNCT__ "DMGetColoring_Composite" 16550c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring) 165647c6ae99SBarry Smith { 165747c6ae99SBarry Smith PetscErrorCode ierr; 165847c6ae99SBarry Smith PetscInt n,i,cnt; 165947c6ae99SBarry Smith ISColoringValue *colors; 166047c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 166147c6ae99SBarry Smith ISColoringValue maxcol = 0; 166247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 166347c6ae99SBarry Smith 166447c6ae99SBarry Smith PetscFunctionBegin; 166547c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 166647c6ae99SBarry Smith if (ctype == IS_COLORING_GHOSTED) { 166747c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" ); 166847c6ae99SBarry Smith } else if (ctype == IS_COLORING_GLOBAL) { 166947c6ae99SBarry Smith n = com->n; 167047c6ae99SBarry Smith } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType"); 167147c6ae99SBarry Smith ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */ 167247c6ae99SBarry Smith 1673671f6225SBarry Smith ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr); 167447c6ae99SBarry Smith if (dense) { 167547c6ae99SBarry Smith for (i=0; i<n; i++) { 167647c6ae99SBarry Smith colors[i] = (ISColoringValue)(com->rstart + i); 167747c6ae99SBarry Smith } 167847c6ae99SBarry Smith maxcol = com->N; 167947c6ae99SBarry Smith } else { 168047c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 168147c6ae99SBarry Smith PetscMPIInt rank; 168247c6ae99SBarry Smith 168347c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 168447c6ae99SBarry Smith cnt = 0; 168547c6ae99SBarry Smith while (next) { 168647c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 168747c6ae99SBarry Smith if (rank == next->rank) { /* each column gets is own color */ 168847c6ae99SBarry Smith for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) { 168947c6ae99SBarry Smith colors[cnt++] = maxcol++; 169047c6ae99SBarry Smith } 169147c6ae99SBarry Smith } 169247c6ae99SBarry Smith ierr = MPI_Bcast(&maxcol,1,MPIU_COLORING_VALUE,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr); 169347c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 169447c6ae99SBarry Smith ISColoring lcoloring; 169547c6ae99SBarry Smith 169647c6ae99SBarry Smith ierr = DMGetColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr); 169747c6ae99SBarry Smith for (i=0; i<lcoloring->N; i++) { 169847c6ae99SBarry Smith colors[cnt++] = maxcol + lcoloring->colors[i]; 169947c6ae99SBarry Smith } 170047c6ae99SBarry Smith maxcol += lcoloring->n; 170147c6ae99SBarry Smith ierr = ISColoringDestroy(lcoloring);CHKERRQ(ierr); 170247c6ae99SBarry Smith } else { 170347c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 170447c6ae99SBarry Smith } 170547c6ae99SBarry Smith next = next->next; 170647c6ae99SBarry Smith } 170747c6ae99SBarry Smith } 170847c6ae99SBarry Smith ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr); 170947c6ae99SBarry Smith PetscFunctionReturn(0); 171047c6ae99SBarry Smith } 171147c6ae99SBarry Smith 171247c6ae99SBarry Smith #undef __FUNCT__ 17130c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin_Composite" 17140c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec) 171547c6ae99SBarry Smith { 171647c6ae99SBarry Smith PetscErrorCode ierr; 171747c6ae99SBarry Smith struct DMCompositeLink *next; 171847c6ae99SBarry Smith PetscInt cnt = 3; 171947c6ae99SBarry Smith PetscMPIInt rank; 172047c6ae99SBarry Smith PetscScalar *garray,*larray; 172147c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 172247c6ae99SBarry Smith 172347c6ae99SBarry Smith PetscFunctionBegin; 172447c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 172547c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 172647c6ae99SBarry Smith next = com->next; 172747c6ae99SBarry Smith if (!com->setup) { 1728d7bf68aeSBarry Smith ierr = DMSetUp(dm);CHKERRQ(ierr); 172947c6ae99SBarry Smith } 173047c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 173147c6ae99SBarry Smith ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr); 173247c6ae99SBarry Smith ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr); 173347c6ae99SBarry Smith 173447c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 173547c6ae99SBarry Smith while (next) { 173647c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 173747c6ae99SBarry Smith if (rank == next->rank) { 173847c6ae99SBarry Smith ierr = PetscMemcpy(larray,garray,next->n*sizeof(PetscScalar));CHKERRQ(ierr); 173947c6ae99SBarry Smith garray += next->n; 174047c6ae99SBarry Smith } 174147c6ae99SBarry Smith /* does not handle ADD_VALUES */ 174247c6ae99SBarry Smith ierr = MPI_Bcast(larray,next->n,MPIU_SCALAR,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr); 174347c6ae99SBarry Smith larray += next->n; 174447c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 174547c6ae99SBarry Smith Vec local,global; 174647c6ae99SBarry Smith PetscInt N; 174747c6ae99SBarry Smith 174847c6ae99SBarry Smith ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 174947c6ae99SBarry Smith ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr); 175047c6ae99SBarry Smith ierr = VecPlaceArray(global,garray);CHKERRQ(ierr); 175147c6ae99SBarry Smith ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr); 175247c6ae99SBarry Smith ierr = VecPlaceArray(local,larray);CHKERRQ(ierr); 175347c6ae99SBarry Smith ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr); 175447c6ae99SBarry Smith ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr); 175547c6ae99SBarry Smith ierr = VecResetArray(global);CHKERRQ(ierr); 175647c6ae99SBarry Smith ierr = VecResetArray(local);CHKERRQ(ierr); 175747c6ae99SBarry Smith ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 175847c6ae99SBarry Smith ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr); 175947c6ae99SBarry Smith larray += next->n; 176047c6ae99SBarry Smith } else { 176147c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 176247c6ae99SBarry Smith } 176347c6ae99SBarry Smith cnt++; 176447c6ae99SBarry Smith next = next->next; 176547c6ae99SBarry Smith } 176647c6ae99SBarry Smith 176747c6ae99SBarry Smith ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr); 176847c6ae99SBarry Smith ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr); 176947c6ae99SBarry Smith PetscFunctionReturn(0); 177047c6ae99SBarry Smith } 177147c6ae99SBarry Smith 177247c6ae99SBarry Smith #undef __FUNCT__ 17730c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd_Composite" 17740c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec) 17750c010503SBarry Smith { 17760c010503SBarry Smith PetscFunctionBegin; 17770c010503SBarry Smith PetscFunctionReturn(0); 17780c010503SBarry Smith } 177947c6ae99SBarry Smith 1780a4121054SBarry Smith EXTERN_C_BEGIN 1781a4121054SBarry Smith #undef __FUNCT__ 1782a4121054SBarry Smith #define __FUNCT__ "DMCreate_Composite" 1783a4121054SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreate_Composite(DM p) 1784a4121054SBarry Smith { 1785a4121054SBarry Smith PetscErrorCode ierr; 1786a4121054SBarry Smith DM_Composite *com; 1787a4121054SBarry Smith 1788a4121054SBarry Smith PetscFunctionBegin; 1789a4121054SBarry Smith ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr); 1790a4121054SBarry Smith p->data = com; 1791a4121054SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr); 1792a4121054SBarry Smith com->n = 0; 1793a4121054SBarry Smith com->next = PETSC_NULL; 1794a4121054SBarry Smith com->nredundant = 0; 1795a4121054SBarry Smith com->nDM = 0; 1796a4121054SBarry Smith 1797a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Composite; 1798a4121054SBarry Smith p->ops->createlocalvector = DMCreateLocalVector_Composite; 1799a4121054SBarry Smith p->ops->refine = DMRefine_Composite; 1800a4121054SBarry Smith p->ops->getinterpolation = DMGetInterpolation_Composite; 1801a4121054SBarry Smith p->ops->getmatrix = DMGetMatrix_Composite; 1802a4121054SBarry Smith p->ops->getcoloring = DMGetColoring_Composite; 1803a4121054SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite; 1804a4121054SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite; 1805a4121054SBarry Smith p->ops->destroy = DMDestroy_Composite; 1806a4121054SBarry Smith p->ops->view = DMView_Composite; 1807a4121054SBarry Smith p->ops->setup = DMSetUp_Composite; 1808a4121054SBarry Smith PetscFunctionReturn(0); 1809a4121054SBarry Smith } 1810a4121054SBarry Smith EXTERN_C_END 1811a4121054SBarry Smith 18120c010503SBarry Smith #undef __FUNCT__ 18130c010503SBarry Smith #define __FUNCT__ "DMCompositeCreate" 18140c010503SBarry Smith /*@C 18150c010503SBarry Smith DMCompositeCreate - Creates a vector packer, used to generate "composite" 18160c010503SBarry Smith vectors made up of several subvectors. 18170c010503SBarry Smith 18180c010503SBarry Smith Collective on MPI_Comm 181947c6ae99SBarry Smith 182047c6ae99SBarry Smith Input Parameter: 18210c010503SBarry Smith . comm - the processors that will share the global vector 18220c010503SBarry Smith 18230c010503SBarry Smith Output Parameters: 18240c010503SBarry Smith . packer - the packer object 182547c6ae99SBarry Smith 182647c6ae99SBarry Smith Level: advanced 182747c6ae99SBarry Smith 18280c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCompositeScatter(), 18290c010503SBarry Smith DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess() 183047c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 183147c6ae99SBarry Smith 183247c6ae99SBarry Smith @*/ 18330c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeCreate(MPI_Comm comm,DM *packer) 183447c6ae99SBarry Smith { 18350c010503SBarry Smith PetscErrorCode ierr; 18360c010503SBarry Smith 183747c6ae99SBarry Smith PetscFunctionBegin; 18380c010503SBarry Smith PetscValidPointer(packer,2); 1839a4121054SBarry Smith ierr = DMCreate(comm,packer);CHKERRQ(ierr); 1840a4121054SBarry Smith ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr); 184147c6ae99SBarry Smith PetscFunctionReturn(0); 184247c6ae99SBarry Smith } 1843