147c6ae99SBarry Smith #define PETSCDM_DLL 247c6ae99SBarry Smith 347c6ae99SBarry Smith #include "petscda.h" /*I "petscda.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 */ 3047c6ae99SBarry Smith PetscInt nghost; /* number of all local entries include DA 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 4247c6ae99SBarry Smith seperate components (DA'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__ 122*0c010503SBarry Smith #define __FUNCT__ "DMDestroy_Composite" 123*0c010503SBarry 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__ 153*0c010503SBarry Smith #define __FUNCT__ "DMView_Composite" 154*0c010503SBarry 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__ 18647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeSetUp" 18747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeSetUp(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" 32947c6ae99SBarry Smith PetscErrorCode DMCompositeGather_Array(DM dm,struct DMCompositeLink *mine,Vec vec,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()"); 34047c6ae99SBarry Smith ierr = PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(PetscScalar));CHKERRQ(ierr); 34147c6ae99SBarry Smith ierr = VecRestoreArray(vec,&varray);CHKERRQ(ierr); 34247c6ae99SBarry Smith } 34347c6ae99SBarry Smith PetscFunctionReturn(0); 34447c6ae99SBarry Smith } 34547c6ae99SBarry Smith 34647c6ae99SBarry Smith #undef __FUNCT__ 34747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_DM" 34847c6ae99SBarry Smith PetscErrorCode DMCompositeGather_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec local) 34947c6ae99SBarry Smith { 35047c6ae99SBarry Smith PetscErrorCode ierr; 35147c6ae99SBarry Smith PetscScalar *array; 35247c6ae99SBarry Smith Vec global; 35347c6ae99SBarry Smith 35447c6ae99SBarry Smith PetscFunctionBegin; 35547c6ae99SBarry Smith ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr); 35647c6ae99SBarry Smith ierr = VecGetArray(vec,&array);CHKERRQ(ierr); 35747c6ae99SBarry Smith ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr); 35847c6ae99SBarry Smith ierr = DMLocalToGlobal(mine->dm,local,INSERT_VALUES,global);CHKERRQ(ierr); 35947c6ae99SBarry Smith ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr); 36047c6ae99SBarry Smith ierr = VecResetArray(global);CHKERRQ(ierr); 36147c6ae99SBarry Smith ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr); 36247c6ae99SBarry Smith PetscFunctionReturn(0); 36347c6ae99SBarry Smith } 36447c6ae99SBarry Smith 36547c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/ 36647c6ae99SBarry Smith 36747c6ae99SBarry Smith #include <stdarg.h> 36847c6ae99SBarry Smith 36947c6ae99SBarry Smith #undef __FUNCT__ 37047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetNumberDM" 37147c6ae99SBarry Smith /*@C 37247c6ae99SBarry Smith DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite 37347c6ae99SBarry Smith representation. 37447c6ae99SBarry Smith 37547c6ae99SBarry Smith Not Collective 37647c6ae99SBarry Smith 37747c6ae99SBarry Smith Input Parameter: 37847c6ae99SBarry Smith . dm - the packer object 37947c6ae99SBarry Smith 38047c6ae99SBarry Smith Output Parameter: 38147c6ae99SBarry Smith . nDM - the number of DMs 38247c6ae99SBarry Smith 38347c6ae99SBarry Smith Level: beginner 38447c6ae99SBarry Smith 38547c6ae99SBarry Smith @*/ 38647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetNumberDM(DM dm,PetscInt *nDM) 38747c6ae99SBarry Smith { 38847c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 38947c6ae99SBarry Smith PetscFunctionBegin; 39047c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 39147c6ae99SBarry Smith *nDM = com->nDM; 39247c6ae99SBarry Smith PetscFunctionReturn(0); 39347c6ae99SBarry Smith } 39447c6ae99SBarry Smith 39547c6ae99SBarry Smith 39647c6ae99SBarry Smith #undef __FUNCT__ 39747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess" 39847c6ae99SBarry Smith /*@C 39947c6ae99SBarry Smith DMCompositeGetAccess - Allows one to access the individual packed vectors in their global 40047c6ae99SBarry Smith representation. 40147c6ae99SBarry Smith 40247c6ae99SBarry Smith Collective on DMComposite 40347c6ae99SBarry Smith 40447c6ae99SBarry Smith Input Parameter: 40547c6ae99SBarry Smith + dm - the packer object 40647c6ae99SBarry Smith . gvec - the global vector 40747c6ae99SBarry Smith - ... - the individual sequential or parallel objects (arrays or vectors) 40847c6ae99SBarry Smith 40947c6ae99SBarry Smith Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them 41047c6ae99SBarry Smith 41147c6ae99SBarry Smith Level: advanced 41247c6ae99SBarry Smith 41347c6ae99SBarry Smith @*/ 41447c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetAccess(DM dm,Vec gvec,...) 41547c6ae99SBarry Smith { 41647c6ae99SBarry Smith va_list Argp; 41747c6ae99SBarry Smith PetscErrorCode ierr; 41847c6ae99SBarry Smith struct DMCompositeLink *next; 41947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 42047c6ae99SBarry Smith 42147c6ae99SBarry Smith PetscFunctionBegin; 42247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 42347c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 42447c6ae99SBarry Smith next = com->next; 42547c6ae99SBarry Smith if (!com->setup) { 42647c6ae99SBarry Smith ierr = DMCompositeSetUp(dm);CHKERRQ(ierr); 42747c6ae99SBarry Smith } 42847c6ae99SBarry Smith 42947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 43047c6ae99SBarry Smith va_start(Argp,gvec); 43147c6ae99SBarry Smith while (next) { 43247c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 43347c6ae99SBarry Smith PetscScalar **array; 43447c6ae99SBarry Smith array = va_arg(Argp, PetscScalar**); 43547c6ae99SBarry Smith ierr = DMCompositeGetAccess_Array(dm,next,gvec,array);CHKERRQ(ierr); 43647c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 43747c6ae99SBarry Smith Vec *vec; 43847c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 43947c6ae99SBarry Smith ierr = DMCompositeGetAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr); 44047c6ae99SBarry Smith } else { 44147c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 44247c6ae99SBarry Smith } 44347c6ae99SBarry Smith next = next->next; 44447c6ae99SBarry Smith } 44547c6ae99SBarry Smith va_end(Argp); 44647c6ae99SBarry Smith PetscFunctionReturn(0); 44747c6ae99SBarry Smith } 44847c6ae99SBarry Smith 44947c6ae99SBarry Smith #undef __FUNCT__ 45047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess" 45147c6ae99SBarry Smith /*@C 45247c6ae99SBarry Smith DMCompositeRestoreAccess - Returns the vectors obtained with DACompositeGetAccess() 45347c6ae99SBarry Smith representation. 45447c6ae99SBarry Smith 45547c6ae99SBarry Smith Collective on DMComposite 45647c6ae99SBarry Smith 45747c6ae99SBarry Smith Input Parameter: 45847c6ae99SBarry Smith + dm - the packer object 45947c6ae99SBarry Smith . gvec - the global vector 46047c6ae99SBarry Smith - ... - the individual sequential or parallel objects (arrays or vectors) 46147c6ae99SBarry Smith 46247c6ae99SBarry Smith Level: advanced 46347c6ae99SBarry Smith 464*0c010503SBarry Smith .seealso DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 46547c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeScatter(), 46647c6ae99SBarry Smith DMCompositeRestoreAccess(), DACompositeGetAccess() 46747c6ae99SBarry Smith 46847c6ae99SBarry Smith @*/ 46947c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreAccess(DM dm,Vec gvec,...) 47047c6ae99SBarry Smith { 47147c6ae99SBarry Smith va_list Argp; 47247c6ae99SBarry Smith PetscErrorCode ierr; 47347c6ae99SBarry Smith struct DMCompositeLink *next; 47447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 47547c6ae99SBarry Smith 47647c6ae99SBarry Smith PetscFunctionBegin; 47747c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 47847c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 47947c6ae99SBarry Smith next = com->next; 48047c6ae99SBarry Smith if (!com->setup) { 48147c6ae99SBarry Smith ierr = DMCompositeSetUp(dm);CHKERRQ(ierr); 48247c6ae99SBarry Smith } 48347c6ae99SBarry Smith 48447c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 48547c6ae99SBarry Smith va_start(Argp,gvec); 48647c6ae99SBarry Smith while (next) { 48747c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 48847c6ae99SBarry Smith PetscScalar **array; 48947c6ae99SBarry Smith array = va_arg(Argp, PetscScalar**); 49047c6ae99SBarry Smith ierr = DMCompositeRestoreAccess_Array(dm,next,gvec,array);CHKERRQ(ierr); 49147c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 49247c6ae99SBarry Smith Vec *vec; 49347c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 49447c6ae99SBarry Smith ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr); 49547c6ae99SBarry Smith } else { 49647c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 49747c6ae99SBarry Smith } 49847c6ae99SBarry Smith next = next->next; 49947c6ae99SBarry Smith } 50047c6ae99SBarry Smith va_end(Argp); 50147c6ae99SBarry Smith PetscFunctionReturn(0); 50247c6ae99SBarry Smith } 50347c6ae99SBarry Smith 50447c6ae99SBarry Smith #undef __FUNCT__ 50547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter" 50647c6ae99SBarry Smith /*@C 50747c6ae99SBarry Smith DMCompositeScatter - Scatters from a global packed vector into its individual local vectors 50847c6ae99SBarry Smith 50947c6ae99SBarry Smith Collective on DMComposite 51047c6ae99SBarry Smith 51147c6ae99SBarry Smith Input Parameter: 51247c6ae99SBarry Smith + dm - the packer object 51347c6ae99SBarry Smith . gvec - the global vector 51447c6ae99SBarry Smith - ... - the individual sequential objects (arrays or vectors) 51547c6ae99SBarry Smith 51647c6ae99SBarry Smith Level: advanced 51747c6ae99SBarry Smith 518*0c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 51947c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 52047c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 52147c6ae99SBarry Smith 52247c6ae99SBarry Smith @*/ 52347c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeScatter(DM dm,Vec gvec,...) 52447c6ae99SBarry Smith { 52547c6ae99SBarry Smith va_list Argp; 52647c6ae99SBarry Smith PetscErrorCode ierr; 52747c6ae99SBarry Smith struct DMCompositeLink *next; 52847c6ae99SBarry Smith PetscInt cnt = 3; 52947c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 53047c6ae99SBarry Smith 53147c6ae99SBarry Smith PetscFunctionBegin; 53247c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 53347c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 53447c6ae99SBarry Smith next = com->next; 53547c6ae99SBarry Smith if (!com->setup) { 53647c6ae99SBarry Smith ierr = DMCompositeSetUp(dm);CHKERRQ(ierr); 53747c6ae99SBarry Smith } 53847c6ae99SBarry Smith 53947c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 54047c6ae99SBarry Smith va_start(Argp,gvec); 54147c6ae99SBarry Smith while (next) { 54247c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 54347c6ae99SBarry Smith PetscScalar *array; 54447c6ae99SBarry Smith array = va_arg(Argp, PetscScalar*); 54547c6ae99SBarry Smith ierr = DMCompositeScatter_Array(dm,next,gvec,array);CHKERRQ(ierr); 54647c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 54747c6ae99SBarry Smith Vec vec; 54847c6ae99SBarry Smith vec = va_arg(Argp, Vec); 54947c6ae99SBarry Smith PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt); 55047c6ae99SBarry Smith ierr = DMCompositeScatter_DM(dm,next,gvec,vec);CHKERRQ(ierr); 55147c6ae99SBarry Smith } else { 55247c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 55347c6ae99SBarry Smith } 55447c6ae99SBarry Smith cnt++; 55547c6ae99SBarry Smith next = next->next; 55647c6ae99SBarry Smith } 55747c6ae99SBarry Smith va_end(Argp); 55847c6ae99SBarry Smith PetscFunctionReturn(0); 55947c6ae99SBarry Smith } 56047c6ae99SBarry Smith 56147c6ae99SBarry Smith #undef __FUNCT__ 56247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather" 56347c6ae99SBarry Smith /*@C 56447c6ae99SBarry Smith DMCompositeGather - Gathers into a global packed vector from its individual local vectors 56547c6ae99SBarry Smith 56647c6ae99SBarry Smith Collective on DMComposite 56747c6ae99SBarry Smith 56847c6ae99SBarry Smith Input Parameter: 56947c6ae99SBarry Smith + dm - the packer object 57047c6ae99SBarry Smith . gvec - the global vector 57147c6ae99SBarry Smith - ... - the individual sequential objects (arrays or vectors) 57247c6ae99SBarry Smith 57347c6ae99SBarry Smith Level: advanced 57447c6ae99SBarry Smith 575*0c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 57647c6ae99SBarry Smith DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 57747c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 57847c6ae99SBarry Smith 57947c6ae99SBarry Smith @*/ 58047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGather(DM dm,Vec gvec,...) 58147c6ae99SBarry Smith { 58247c6ae99SBarry Smith va_list Argp; 58347c6ae99SBarry Smith PetscErrorCode ierr; 58447c6ae99SBarry Smith struct DMCompositeLink *next; 58547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 58647c6ae99SBarry Smith 58747c6ae99SBarry Smith PetscFunctionBegin; 58847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 58947c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 59047c6ae99SBarry Smith next = com->next; 59147c6ae99SBarry Smith if (!com->setup) { 59247c6ae99SBarry Smith ierr = DMCompositeSetUp(dm);CHKERRQ(ierr); 59347c6ae99SBarry Smith } 59447c6ae99SBarry Smith 59547c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 59647c6ae99SBarry Smith va_start(Argp,gvec); 59747c6ae99SBarry Smith while (next) { 59847c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 59947c6ae99SBarry Smith PetscScalar *array; 60047c6ae99SBarry Smith array = va_arg(Argp, PetscScalar*); 60147c6ae99SBarry Smith ierr = DMCompositeGather_Array(dm,next,gvec,array);CHKERRQ(ierr); 60247c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 60347c6ae99SBarry Smith Vec vec; 60447c6ae99SBarry Smith vec = va_arg(Argp, Vec); 60547c6ae99SBarry Smith PetscValidHeaderSpecific(vec,VEC_CLASSID,3); 60647c6ae99SBarry Smith ierr = DMCompositeGather_DM(dm,next,gvec,vec);CHKERRQ(ierr); 60747c6ae99SBarry Smith } else { 60847c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 60947c6ae99SBarry Smith } 61047c6ae99SBarry Smith next = next->next; 61147c6ae99SBarry Smith } 61247c6ae99SBarry Smith va_end(Argp); 61347c6ae99SBarry Smith PetscFunctionReturn(0); 61447c6ae99SBarry Smith } 61547c6ae99SBarry Smith 61647c6ae99SBarry Smith #undef __FUNCT__ 61747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddArray" 61847c6ae99SBarry Smith /*@C 61947c6ae99SBarry Smith DMCompositeAddArray - adds an "redundant" array to a DMComposite. The array values will 62047c6ae99SBarry Smith be stored in part of the array on process orank. 62147c6ae99SBarry Smith 62247c6ae99SBarry Smith Collective on DMComposite 62347c6ae99SBarry Smith 62447c6ae99SBarry Smith Input Parameter: 62547c6ae99SBarry Smith + dm - the packer object 62647c6ae99SBarry Smith . orank - the process on which the array entries officially live, this number must be 62747c6ae99SBarry Smith the same on all processes. 62847c6ae99SBarry Smith - n - the length of the array 62947c6ae99SBarry Smith 63047c6ae99SBarry Smith Level: advanced 63147c6ae99SBarry Smith 632*0c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(), 63347c6ae99SBarry Smith DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 63447c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 63547c6ae99SBarry Smith 63647c6ae99SBarry Smith @*/ 63747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeAddArray(DM dm,PetscMPIInt orank,PetscInt n) 63847c6ae99SBarry Smith { 63947c6ae99SBarry Smith struct DMCompositeLink *mine,*next; 64047c6ae99SBarry Smith PetscErrorCode ierr; 64147c6ae99SBarry Smith PetscMPIInt rank; 64247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 64347c6ae99SBarry Smith 64447c6ae99SBarry Smith PetscFunctionBegin; 64547c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 64647c6ae99SBarry Smith next = com->next; 64747c6ae99SBarry Smith if (com->setup) { 64847c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add an array once you have used the DMComposite"); 64947c6ae99SBarry Smith } 65047c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG) 65147c6ae99SBarry Smith { 65247c6ae99SBarry Smith PetscMPIInt orankmax; 65347c6ae99SBarry Smith ierr = MPI_Allreduce(&orank,&orankmax,1,MPI_INT,MPI_MAX,((PetscObject)dm)->comm);CHKERRQ(ierr); 65447c6ae99SBarry 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); 65547c6ae99SBarry Smith } 65647c6ae99SBarry Smith #endif 65747c6ae99SBarry Smith 65847c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 65947c6ae99SBarry Smith /* create new link */ 66047c6ae99SBarry Smith ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr); 66147c6ae99SBarry Smith mine->n = n; 66247c6ae99SBarry Smith mine->rank = orank; 66347c6ae99SBarry Smith mine->dm = PETSC_NULL; 66447c6ae99SBarry Smith mine->type = DMCOMPOSITE_ARRAY; 66547c6ae99SBarry Smith mine->next = PETSC_NULL; 66647c6ae99SBarry Smith if (rank == mine->rank) {com->n += n;com->nmine++;} 66747c6ae99SBarry Smith 66847c6ae99SBarry Smith /* add to end of list */ 66947c6ae99SBarry Smith if (!next) { 67047c6ae99SBarry Smith com->next = mine; 67147c6ae99SBarry Smith } else { 67247c6ae99SBarry Smith while (next->next) next = next->next; 67347c6ae99SBarry Smith next->next = mine; 67447c6ae99SBarry Smith } 67547c6ae99SBarry Smith com->nredundant++; 67647c6ae99SBarry Smith PetscFunctionReturn(0); 67747c6ae99SBarry Smith } 67847c6ae99SBarry Smith 67947c6ae99SBarry Smith #undef __FUNCT__ 68047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddDM" 68147c6ae99SBarry Smith /*@C 68247c6ae99SBarry Smith DMCompositeAddDM - adds a DM (includes DA) vector to a DMComposite 68347c6ae99SBarry Smith 68447c6ae99SBarry Smith Collective on DMComposite 68547c6ae99SBarry Smith 68647c6ae99SBarry Smith Input Parameter: 68747c6ae99SBarry Smith + dm - the packer object 68847c6ae99SBarry Smith - dm - the DM object, if the DM is a da you will need to caste it with a (DM) 68947c6ae99SBarry Smith 69047c6ae99SBarry Smith Level: advanced 69147c6ae99SBarry Smith 692*0c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(), 69347c6ae99SBarry Smith DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 69447c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 69547c6ae99SBarry Smith 69647c6ae99SBarry Smith @*/ 69747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeAddDM(DM dmc,DM dm) 69847c6ae99SBarry Smith { 69947c6ae99SBarry Smith PetscErrorCode ierr; 70047c6ae99SBarry Smith PetscInt n; 70147c6ae99SBarry Smith struct DMCompositeLink *mine,*next; 70247c6ae99SBarry Smith Vec global; 70347c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dmc->data; 70447c6ae99SBarry Smith 70547c6ae99SBarry Smith PetscFunctionBegin; 70647c6ae99SBarry Smith PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 70747c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,2); 70847c6ae99SBarry Smith next = com->next; 70947c6ae99SBarry Smith if (com->setup) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DA once you have used the DMComposite"); 71047c6ae99SBarry Smith 71147c6ae99SBarry Smith /* create new link */ 71247c6ae99SBarry Smith ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr); 71347c6ae99SBarry Smith ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 71447c6ae99SBarry Smith ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr); 71547c6ae99SBarry Smith ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr); 71647c6ae99SBarry Smith ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr); 71747c6ae99SBarry Smith mine->n = n; 71847c6ae99SBarry Smith mine->dm = dm; 71947c6ae99SBarry Smith mine->type = DMCOMPOSITE_DM; 72047c6ae99SBarry Smith mine->next = PETSC_NULL; 72147c6ae99SBarry Smith com->n += n; 72247c6ae99SBarry Smith 72347c6ae99SBarry Smith /* add to end of list */ 72447c6ae99SBarry Smith if (!next) { 72547c6ae99SBarry Smith com->next = mine; 72647c6ae99SBarry Smith } else { 72747c6ae99SBarry Smith while (next->next) next = next->next; 72847c6ae99SBarry Smith next->next = mine; 72947c6ae99SBarry Smith } 73047c6ae99SBarry Smith com->nDM++; 73147c6ae99SBarry Smith com->nmine++; 73247c6ae99SBarry Smith PetscFunctionReturn(0); 73347c6ae99SBarry Smith } 73447c6ae99SBarry Smith 73547c6ae99SBarry Smith extern PetscErrorCode PETSCDM_DLLEXPORT VecView_MPI(Vec,PetscViewer); 73647c6ae99SBarry Smith EXTERN_C_BEGIN 73747c6ae99SBarry Smith #undef __FUNCT__ 73847c6ae99SBarry Smith #define __FUNCT__ "VecView_DMComposite" 73947c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT VecView_DMComposite(Vec gvec,PetscViewer viewer) 74047c6ae99SBarry Smith { 74147c6ae99SBarry Smith DM dm; 74247c6ae99SBarry Smith PetscErrorCode ierr; 74347c6ae99SBarry Smith struct DMCompositeLink *next; 74447c6ae99SBarry Smith PetscBool isdraw; 74547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 74647c6ae99SBarry Smith 74747c6ae99SBarry Smith PetscFunctionBegin; 74847c6ae99SBarry Smith ierr = PetscObjectQuery((PetscObject)gvec,"DMComposite",(PetscObject*)&dm);CHKERRQ(ierr); 74947c6ae99SBarry Smith if (!dm) SETERRQ(((PetscObject)gvec)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite"); 75047c6ae99SBarry Smith com = (DM_Composite*)dm->data; 75147c6ae99SBarry Smith next = com->next; 75247c6ae99SBarry Smith 75347c6ae99SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 75447c6ae99SBarry Smith if (!isdraw) { 75547c6ae99SBarry Smith /* do I really want to call this? */ 75647c6ae99SBarry Smith ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr); 75747c6ae99SBarry Smith } else { 75847c6ae99SBarry Smith PetscInt cnt = 0; 75947c6ae99SBarry Smith 76047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 76147c6ae99SBarry Smith while (next) { 76247c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 76347c6ae99SBarry Smith PetscScalar *array; 76447c6ae99SBarry Smith ierr = DMCompositeGetAccess_Array(dm,next,gvec,&array);CHKERRQ(ierr); 76547c6ae99SBarry Smith 76647c6ae99SBarry Smith /*skip it for now */ 76747c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 76847c6ae99SBarry Smith Vec vec; 76947c6ae99SBarry Smith PetscInt bs; 77047c6ae99SBarry Smith 77147c6ae99SBarry Smith ierr = DMCompositeGetAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr); 77247c6ae99SBarry Smith ierr = VecView(vec,viewer);CHKERRQ(ierr); 77347c6ae99SBarry Smith ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr); 77447c6ae99SBarry Smith ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr); 77547c6ae99SBarry Smith ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr); 77647c6ae99SBarry Smith cnt += bs; 77747c6ae99SBarry Smith } else { 77847c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 77947c6ae99SBarry Smith } 78047c6ae99SBarry Smith next = next->next; 78147c6ae99SBarry Smith } 78247c6ae99SBarry Smith ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr); 78347c6ae99SBarry Smith } 78447c6ae99SBarry Smith PetscFunctionReturn(0); 78547c6ae99SBarry Smith } 78647c6ae99SBarry Smith EXTERN_C_END 78747c6ae99SBarry Smith 78847c6ae99SBarry Smith 78947c6ae99SBarry Smith #undef __FUNCT__ 790*0c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Composite" 791*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateGlobalVector_Composite(DM dm,Vec *gvec) 79247c6ae99SBarry Smith { 79347c6ae99SBarry Smith PetscErrorCode ierr; 79447c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 79547c6ae99SBarry Smith 79647c6ae99SBarry Smith PetscFunctionBegin; 79747c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 79847c6ae99SBarry Smith if (!com->setup) { 79947c6ae99SBarry Smith ierr = DMCompositeSetUp(dm);CHKERRQ(ierr); 80047c6ae99SBarry Smith } 80147c6ae99SBarry Smith ierr = VecCreateMPI(((PetscObject)dm)->comm,com->n,com->N,gvec);CHKERRQ(ierr); 80247c6ae99SBarry Smith ierr = PetscObjectCompose((PetscObject)*gvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr); 80347c6ae99SBarry Smith ierr = VecSetOperation(*gvec,VECOP_VIEW,(void(*)(void))VecView_DMComposite);CHKERRQ(ierr); 80447c6ae99SBarry Smith PetscFunctionReturn(0); 80547c6ae99SBarry Smith } 80647c6ae99SBarry Smith 80747c6ae99SBarry Smith #undef __FUNCT__ 808*0c010503SBarry Smith #define __FUNCT__ "DMCreateLocalVector_Composite" 809*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateLocalVector_Composite(DM dm,Vec *lvec) 81047c6ae99SBarry Smith { 81147c6ae99SBarry Smith PetscErrorCode ierr; 81247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 81347c6ae99SBarry Smith 81447c6ae99SBarry Smith PetscFunctionBegin; 81547c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 81647c6ae99SBarry Smith if (!com->setup) { 81747c6ae99SBarry Smith ierr = DMCompositeSetUp(dm);CHKERRQ(ierr); 81847c6ae99SBarry Smith } 81947c6ae99SBarry Smith ierr = VecCreateSeq(((PetscObject)dm)->comm,com->nghost,lvec);CHKERRQ(ierr); 82047c6ae99SBarry Smith ierr = PetscObjectCompose((PetscObject)*lvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr); 82147c6ae99SBarry Smith PetscFunctionReturn(0); 82247c6ae99SBarry Smith } 82347c6ae99SBarry Smith 82447c6ae99SBarry Smith #undef __FUNCT__ 82547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalISs" 82647c6ae99SBarry Smith /*@C 82747c6ae99SBarry Smith DMCompositeGetLocalISs - gets an IS for each DM/array in the DMComposite, include ghost points 82847c6ae99SBarry Smith 82947c6ae99SBarry Smith Collective on DMComposite 83047c6ae99SBarry Smith 83147c6ae99SBarry Smith Input Parameter: 83247c6ae99SBarry Smith . dm - the packer object 83347c6ae99SBarry Smith 83447c6ae99SBarry Smith Output Parameters: 83547c6ae99SBarry Smith . is - the individual indices for each packed vector/array. Note that this includes 83647c6ae99SBarry Smith all the ghost points that individual ghosted DA's may have. Also each process has an 83747c6ae99SBarry Smith is for EACH redundant array (not just the local redundant arrays). 83847c6ae99SBarry Smith 83947c6ae99SBarry Smith Level: advanced 84047c6ae99SBarry Smith 84147c6ae99SBarry Smith Notes: 84247c6ae99SBarry Smith The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree() 84347c6ae99SBarry Smith 84447c6ae99SBarry Smith Use DMCompositeGetGlobalISs() for non-ghosted ISs. 84547c6ae99SBarry Smith 846*0c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 84747c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(), 84847c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries() 84947c6ae99SBarry Smith 85047c6ae99SBarry Smith @*/ 85147c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetGlobalIndices(DM dm,IS *is[]) 85247c6ae99SBarry Smith { 85347c6ae99SBarry Smith PetscErrorCode ierr; 85447c6ae99SBarry Smith PetscInt i,*idx,n,cnt; 85547c6ae99SBarry Smith struct DMCompositeLink *next; 85647c6ae99SBarry Smith Vec global,dglobal; 85747c6ae99SBarry Smith PF pf; 85847c6ae99SBarry Smith PetscScalar *array; 85947c6ae99SBarry Smith PetscMPIInt rank; 86047c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 86147c6ae99SBarry Smith 86247c6ae99SBarry Smith PetscFunctionBegin; 86347c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 86447c6ae99SBarry Smith ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr); 86547c6ae99SBarry Smith next = com->next; 866*0c010503SBarry Smith ierr = DMCreateGlobalVector(dm,&global);CHKERRQ(ierr); 86747c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 86847c6ae99SBarry Smith 86947c6ae99SBarry Smith /* put 0 to N-1 into the global vector */ 87047c6ae99SBarry Smith ierr = PFCreate(PETSC_COMM_WORLD,1,1,&pf);CHKERRQ(ierr); 87147c6ae99SBarry Smith ierr = PFSetType(pf,PFIDENTITY,PETSC_NULL);CHKERRQ(ierr); 87247c6ae99SBarry Smith ierr = PFApplyVec(pf,PETSC_NULL,global);CHKERRQ(ierr); 87347c6ae99SBarry Smith ierr = PFDestroy(pf);CHKERRQ(ierr); 87447c6ae99SBarry Smith 87547c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 87647c6ae99SBarry Smith cnt = 0; 87747c6ae99SBarry Smith while (next) { 87847c6ae99SBarry Smith 87947c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 88047c6ae99SBarry Smith 88147c6ae99SBarry Smith ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr); 88247c6ae99SBarry Smith if (rank == next->rank) { 88347c6ae99SBarry Smith ierr = VecGetArray(global,&array);CHKERRQ(ierr); 88447c6ae99SBarry Smith array += next->rstart; 88547c6ae99SBarry Smith for (i=0; i<next->n; i++) idx[i] = (PetscInt)PetscRealPart(array[i]); 88647c6ae99SBarry Smith array -= next->rstart; 88747c6ae99SBarry Smith ierr = VecRestoreArray(global,&array);CHKERRQ(ierr); 88847c6ae99SBarry Smith } 88947c6ae99SBarry Smith ierr = MPI_Bcast(idx,next->n,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr); 89047c6ae99SBarry Smith ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr); 89147c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 89247c6ae99SBarry Smith Vec local; 89347c6ae99SBarry Smith 89447c6ae99SBarry Smith ierr = DMCreateLocalVector(next->dm,&local);CHKERRQ(ierr); 89547c6ae99SBarry Smith ierr = VecGetArray(global,&array);CHKERRQ(ierr); 89647c6ae99SBarry Smith array += next->rstart; 89747c6ae99SBarry Smith ierr = DMGetGlobalVector(next->dm,&dglobal);CHKERRQ(ierr); 89847c6ae99SBarry Smith ierr = VecPlaceArray(dglobal,array);CHKERRQ(ierr); 89947c6ae99SBarry Smith ierr = DMGlobalToLocalBegin(next->dm,dglobal,INSERT_VALUES,local);CHKERRQ(ierr); 90047c6ae99SBarry Smith ierr = DMGlobalToLocalEnd(next->dm,dglobal,INSERT_VALUES,local);CHKERRQ(ierr); 90147c6ae99SBarry Smith array -= next->rstart; 90247c6ae99SBarry Smith ierr = VecRestoreArray(global,&array);CHKERRQ(ierr); 90347c6ae99SBarry Smith ierr = VecResetArray(dglobal);CHKERRQ(ierr); 90447c6ae99SBarry Smith ierr = DMRestoreGlobalVector(next->dm,&dglobal);CHKERRQ(ierr); 90547c6ae99SBarry Smith 90647c6ae99SBarry Smith ierr = VecGetArray(local,&array);CHKERRQ(ierr); 90747c6ae99SBarry Smith ierr = VecGetSize(local,&n);CHKERRQ(ierr); 90847c6ae99SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr); 90947c6ae99SBarry Smith for (i=0; i<n; i++) idx[i] = (PetscInt)PetscRealPart(array[i]); 91047c6ae99SBarry Smith ierr = VecRestoreArray(local,&array);CHKERRQ(ierr); 91147c6ae99SBarry Smith ierr = VecDestroy(local);CHKERRQ(ierr); 91247c6ae99SBarry Smith ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr); 91347c6ae99SBarry Smith 91447c6ae99SBarry Smith } else SETERRQ(((PetscObject)global)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 91547c6ae99SBarry Smith next = next->next; 91647c6ae99SBarry Smith cnt++; 91747c6ae99SBarry Smith } 91847c6ae99SBarry Smith ierr = VecDestroy(global);CHKERRQ(ierr); 91947c6ae99SBarry Smith PetscFunctionReturn(0); 92047c6ae99SBarry Smith } 92147c6ae99SBarry Smith 92247c6ae99SBarry Smith #undef __FUNCT__ 92347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetGlobalISs" 92447c6ae99SBarry Smith /*@C 92547c6ae99SBarry Smith DMCompositeGetGlobalISs - Gets the index sets for each composed object 92647c6ae99SBarry Smith 92747c6ae99SBarry Smith Collective on DMComposite 92847c6ae99SBarry Smith 92947c6ae99SBarry Smith Input Parameter: 93047c6ae99SBarry Smith . dm - the packer object 93147c6ae99SBarry Smith 93247c6ae99SBarry Smith Output Parameters: 93347c6ae99SBarry Smith . is - the array of index sets 93447c6ae99SBarry Smith 93547c6ae99SBarry Smith Level: advanced 93647c6ae99SBarry Smith 93747c6ae99SBarry Smith Notes: 93847c6ae99SBarry Smith The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree() 93947c6ae99SBarry Smith 94047c6ae99SBarry Smith The number of IS on each process will/may be different when redundant arrays are used 94147c6ae99SBarry Smith 94247c6ae99SBarry Smith These could be used to extract a subset of vector entries for a "multi-physics" preconditioner 94347c6ae99SBarry Smith 94447c6ae99SBarry Smith Use DMCompositeGetLocalISs() for index sets that include ghost points 94547c6ae99SBarry Smith 946*0c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 94747c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(), 94847c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries() 94947c6ae99SBarry Smith 95047c6ae99SBarry Smith @*/ 95147c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetGlobalISs(DM dm,IS *is[]) 95247c6ae99SBarry Smith { 95347c6ae99SBarry Smith PetscErrorCode ierr; 95447c6ae99SBarry Smith PetscInt cnt = 0,*idx,i; 95547c6ae99SBarry Smith struct DMCompositeLink *next; 95647c6ae99SBarry Smith PetscMPIInt rank; 95747c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 95847c6ae99SBarry Smith 95947c6ae99SBarry Smith PetscFunctionBegin; 96047c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 96147c6ae99SBarry Smith ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr); 96247c6ae99SBarry Smith next = com->next; 96347c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 96447c6ae99SBarry Smith 96547c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 96647c6ae99SBarry Smith while (next) { 96747c6ae99SBarry Smith 96847c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 96947c6ae99SBarry Smith 97047c6ae99SBarry Smith if (rank == next->rank) { 97147c6ae99SBarry Smith ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr); 97247c6ae99SBarry Smith for (i=0; i<next->n; i++) idx[i] = next->grstart + i; 97347c6ae99SBarry Smith ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr); 97447c6ae99SBarry Smith cnt++; 97547c6ae99SBarry Smith } 97647c6ae99SBarry Smith 97747c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 97847c6ae99SBarry Smith ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr); 97947c6ae99SBarry Smith for (i=0; i<next->n; i++) idx[i] = next->grstart + i; 98047c6ae99SBarry Smith ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr); 98147c6ae99SBarry Smith cnt++; 98247c6ae99SBarry Smith } else { 98347c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 98447c6ae99SBarry Smith } 98547c6ae99SBarry Smith next = next->next; 98647c6ae99SBarry Smith } 98747c6ae99SBarry Smith PetscFunctionReturn(0); 98847c6ae99SBarry Smith } 98947c6ae99SBarry Smith 99047c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 99147c6ae99SBarry Smith #undef __FUNCT__ 99247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_Array" 99347c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array) 99447c6ae99SBarry Smith { 99547c6ae99SBarry Smith PetscErrorCode ierr; 99647c6ae99SBarry Smith PetscFunctionBegin; 99747c6ae99SBarry Smith if (array) { 99847c6ae99SBarry Smith ierr = PetscMalloc(mine->n*sizeof(PetscScalar),array);CHKERRQ(ierr); 99947c6ae99SBarry Smith } 100047c6ae99SBarry Smith PetscFunctionReturn(0); 100147c6ae99SBarry Smith } 100247c6ae99SBarry Smith 100347c6ae99SBarry Smith #undef __FUNCT__ 100447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_DM" 100547c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local) 100647c6ae99SBarry Smith { 100747c6ae99SBarry Smith PetscErrorCode ierr; 100847c6ae99SBarry Smith PetscFunctionBegin; 100947c6ae99SBarry Smith if (local) { 101047c6ae99SBarry Smith ierr = DMGetLocalVector(mine->dm,local);CHKERRQ(ierr); 101147c6ae99SBarry Smith } 101247c6ae99SBarry Smith PetscFunctionReturn(0); 101347c6ae99SBarry Smith } 101447c6ae99SBarry Smith 101547c6ae99SBarry Smith #undef __FUNCT__ 101647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_Array" 101747c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array) 101847c6ae99SBarry Smith { 101947c6ae99SBarry Smith PetscErrorCode ierr; 102047c6ae99SBarry Smith PetscFunctionBegin; 102147c6ae99SBarry Smith if (array) { 102247c6ae99SBarry Smith ierr = PetscFree(*array);CHKERRQ(ierr); 102347c6ae99SBarry Smith } 102447c6ae99SBarry Smith PetscFunctionReturn(0); 102547c6ae99SBarry Smith } 102647c6ae99SBarry Smith 102747c6ae99SBarry Smith #undef __FUNCT__ 102847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_DM" 102947c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local) 103047c6ae99SBarry Smith { 103147c6ae99SBarry Smith PetscErrorCode ierr; 103247c6ae99SBarry Smith PetscFunctionBegin; 103347c6ae99SBarry Smith if (local) { 103447c6ae99SBarry Smith ierr = DMRestoreLocalVector(mine->dm,local);CHKERRQ(ierr); 103547c6ae99SBarry Smith } 103647c6ae99SBarry Smith PetscFunctionReturn(0); 103747c6ae99SBarry Smith } 103847c6ae99SBarry Smith 103947c6ae99SBarry Smith #undef __FUNCT__ 104047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors" 104147c6ae99SBarry Smith /*@C 104247c6ae99SBarry Smith DMCompositeGetLocalVectors - Gets local vectors and arrays for each part of a DMComposite. 104347c6ae99SBarry Smith Use DMCompositeRestoreLocalVectors() to return them. 104447c6ae99SBarry Smith 104547c6ae99SBarry Smith Not Collective 104647c6ae99SBarry Smith 104747c6ae99SBarry Smith Input Parameter: 104847c6ae99SBarry Smith . dm - the packer object 104947c6ae99SBarry Smith 105047c6ae99SBarry Smith Output Parameter: 105147c6ae99SBarry Smith . ... - the individual sequential objects (arrays or vectors) 105247c6ae99SBarry Smith 105347c6ae99SBarry Smith Level: advanced 105447c6ae99SBarry Smith 1055*0c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 105647c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 105747c6ae99SBarry Smith DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries() 105847c6ae99SBarry Smith 105947c6ae99SBarry Smith @*/ 106047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetLocalVectors(DM dm,...) 106147c6ae99SBarry Smith { 106247c6ae99SBarry Smith va_list Argp; 106347c6ae99SBarry Smith PetscErrorCode ierr; 106447c6ae99SBarry Smith struct DMCompositeLink *next; 106547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 106647c6ae99SBarry Smith 106747c6ae99SBarry Smith PetscFunctionBegin; 106847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 106947c6ae99SBarry Smith next = com->next; 107047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 107147c6ae99SBarry Smith va_start(Argp,dm); 107247c6ae99SBarry Smith while (next) { 107347c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 107447c6ae99SBarry Smith PetscScalar **array; 107547c6ae99SBarry Smith array = va_arg(Argp, PetscScalar**); 107647c6ae99SBarry Smith ierr = DMCompositeGetLocalVectors_Array(dm,next,array);CHKERRQ(ierr); 107747c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 107847c6ae99SBarry Smith Vec *vec; 107947c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 108047c6ae99SBarry Smith ierr = DMCompositeGetLocalVectors_DM(dm,next,vec);CHKERRQ(ierr); 108147c6ae99SBarry Smith } else { 108247c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 108347c6ae99SBarry Smith } 108447c6ae99SBarry Smith next = next->next; 108547c6ae99SBarry Smith } 108647c6ae99SBarry Smith va_end(Argp); 108747c6ae99SBarry Smith PetscFunctionReturn(0); 108847c6ae99SBarry Smith } 108947c6ae99SBarry Smith 109047c6ae99SBarry Smith #undef __FUNCT__ 109147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors" 109247c6ae99SBarry Smith /*@C 109347c6ae99SBarry Smith DMCompositeRestoreLocalVectors - Restores local vectors and arrays for each part of a DMComposite. 109447c6ae99SBarry Smith 109547c6ae99SBarry Smith Not Collective 109647c6ae99SBarry Smith 109747c6ae99SBarry Smith Input Parameter: 109847c6ae99SBarry Smith . dm - the packer object 109947c6ae99SBarry Smith 110047c6ae99SBarry Smith Output Parameter: 110147c6ae99SBarry Smith . ... - the individual sequential objects (arrays or vectors) 110247c6ae99SBarry Smith 110347c6ae99SBarry Smith Level: advanced 110447c6ae99SBarry Smith 1105*0c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 110647c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 110747c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries() 110847c6ae99SBarry Smith 110947c6ae99SBarry Smith @*/ 111047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreLocalVectors(DM dm,...) 111147c6ae99SBarry Smith { 111247c6ae99SBarry Smith va_list Argp; 111347c6ae99SBarry Smith PetscErrorCode ierr; 111447c6ae99SBarry Smith struct DMCompositeLink *next; 111547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 111647c6ae99SBarry Smith 111747c6ae99SBarry Smith PetscFunctionBegin; 111847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 111947c6ae99SBarry Smith next = com->next; 112047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 112147c6ae99SBarry Smith va_start(Argp,dm); 112247c6ae99SBarry Smith while (next) { 112347c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 112447c6ae99SBarry Smith PetscScalar **array; 112547c6ae99SBarry Smith array = va_arg(Argp, PetscScalar**); 112647c6ae99SBarry Smith ierr = DMCompositeRestoreLocalVectors_Array(dm,next,array);CHKERRQ(ierr); 112747c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 112847c6ae99SBarry Smith Vec *vec; 112947c6ae99SBarry Smith vec = va_arg(Argp, Vec*); 113047c6ae99SBarry Smith ierr = DMCompositeRestoreLocalVectors_DM(dm,next,vec);CHKERRQ(ierr); 113147c6ae99SBarry Smith } else { 113247c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 113347c6ae99SBarry Smith } 113447c6ae99SBarry Smith next = next->next; 113547c6ae99SBarry Smith } 113647c6ae99SBarry Smith va_end(Argp); 113747c6ae99SBarry Smith PetscFunctionReturn(0); 113847c6ae99SBarry Smith } 113947c6ae99SBarry Smith 114047c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/ 114147c6ae99SBarry Smith #undef __FUNCT__ 114247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_Array" 114347c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_Array(DM dm,struct DMCompositeLink *mine,PetscInt *n) 114447c6ae99SBarry Smith { 114547c6ae99SBarry Smith PetscFunctionBegin; 114647c6ae99SBarry Smith if (n) *n = mine->n; 114747c6ae99SBarry Smith PetscFunctionReturn(0); 114847c6ae99SBarry Smith } 114947c6ae99SBarry Smith 115047c6ae99SBarry Smith #undef __FUNCT__ 115147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_DM" 115247c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_DM(DM dmi,struct DMCompositeLink *mine,DM *dm) 115347c6ae99SBarry Smith { 115447c6ae99SBarry Smith PetscFunctionBegin; 115547c6ae99SBarry Smith if (dm) *dm = mine->dm; 115647c6ae99SBarry Smith PetscFunctionReturn(0); 115747c6ae99SBarry Smith } 115847c6ae99SBarry Smith 115947c6ae99SBarry Smith #undef __FUNCT__ 116047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries" 116147c6ae99SBarry Smith /*@C 116247c6ae99SBarry Smith DMCompositeGetEntries - Gets the DA, redundant size, etc for each entry in a DMComposite. 116347c6ae99SBarry Smith 116447c6ae99SBarry Smith Not Collective 116547c6ae99SBarry Smith 116647c6ae99SBarry Smith Input Parameter: 116747c6ae99SBarry Smith . dm - the packer object 116847c6ae99SBarry Smith 116947c6ae99SBarry Smith Output Parameter: 117047c6ae99SBarry Smith . ... - the individual entries, DAs or integer sizes) 117147c6ae99SBarry Smith 117247c6ae99SBarry Smith Level: advanced 117347c6ae99SBarry Smith 1174*0c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(), 117547c6ae99SBarry Smith DMCompositeGather(), DMCompositeCreate(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess(), 117647c6ae99SBarry Smith DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(), DMCompositeScatter(), 117747c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors() 117847c6ae99SBarry Smith 117947c6ae99SBarry Smith @*/ 118047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetEntries(DM dm,...) 118147c6ae99SBarry Smith { 118247c6ae99SBarry Smith va_list Argp; 118347c6ae99SBarry Smith PetscErrorCode ierr; 118447c6ae99SBarry Smith struct DMCompositeLink *next; 118547c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 118647c6ae99SBarry Smith 118747c6ae99SBarry Smith PetscFunctionBegin; 118847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 118947c6ae99SBarry Smith next = com->next; 119047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 119147c6ae99SBarry Smith va_start(Argp,dm); 119247c6ae99SBarry Smith while (next) { 119347c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 119447c6ae99SBarry Smith PetscInt *n; 119547c6ae99SBarry Smith n = va_arg(Argp, PetscInt*); 119647c6ae99SBarry Smith ierr = DMCompositeGetEntries_Array(dm,next,n);CHKERRQ(ierr); 119747c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 119847c6ae99SBarry Smith DM *dmn; 119947c6ae99SBarry Smith dmn = va_arg(Argp, DM*); 120047c6ae99SBarry Smith ierr = DMCompositeGetEntries_DM(dm,next,dmn);CHKERRQ(ierr); 120147c6ae99SBarry Smith } else { 120247c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 120347c6ae99SBarry Smith } 120447c6ae99SBarry Smith next = next->next; 120547c6ae99SBarry Smith } 120647c6ae99SBarry Smith va_end(Argp); 120747c6ae99SBarry Smith PetscFunctionReturn(0); 120847c6ae99SBarry Smith } 120947c6ae99SBarry Smith 121047c6ae99SBarry Smith #undef __FUNCT__ 1211*0c010503SBarry Smith #define __FUNCT__ "DMRefine_Composite" 1212*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine) 121347c6ae99SBarry Smith { 121447c6ae99SBarry Smith PetscErrorCode ierr; 121547c6ae99SBarry Smith struct DMCompositeLink *next; 121647c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dmi->data; 121747c6ae99SBarry Smith DM dm; 121847c6ae99SBarry Smith 121947c6ae99SBarry Smith PetscFunctionBegin; 122047c6ae99SBarry Smith PetscValidHeaderSpecific(dmi,DM_CLASSID,1); 122147c6ae99SBarry Smith next = com->next; 122247c6ae99SBarry Smith ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr); 122347c6ae99SBarry Smith 122447c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 122547c6ae99SBarry Smith while (next) { 122647c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 122747c6ae99SBarry Smith ierr = DMCompositeAddArray(*fine,next->rank,next->n);CHKERRQ(ierr); 122847c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 122947c6ae99SBarry Smith ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr); 123047c6ae99SBarry Smith ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr); 123147c6ae99SBarry Smith ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 123247c6ae99SBarry Smith } else { 123347c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 123447c6ae99SBarry Smith } 123547c6ae99SBarry Smith next = next->next; 123647c6ae99SBarry Smith } 123747c6ae99SBarry Smith PetscFunctionReturn(0); 123847c6ae99SBarry Smith } 123947c6ae99SBarry Smith 124047c6ae99SBarry Smith #include "petscmat.h" 124147c6ae99SBarry Smith 124247c6ae99SBarry Smith struct MatPackLink { 124347c6ae99SBarry Smith Mat A; 124447c6ae99SBarry Smith struct MatPackLink *next; 124547c6ae99SBarry Smith }; 124647c6ae99SBarry Smith 124747c6ae99SBarry Smith struct MatPack { 124847c6ae99SBarry Smith DM right,left; 124947c6ae99SBarry Smith struct MatPackLink *next; 125047c6ae99SBarry Smith }; 125147c6ae99SBarry Smith 125247c6ae99SBarry Smith #undef __FUNCT__ 125347c6ae99SBarry Smith #define __FUNCT__ "MatMultBoth_Shell_Pack" 125447c6ae99SBarry Smith PetscErrorCode MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscBool add) 125547c6ae99SBarry Smith { 125647c6ae99SBarry Smith struct MatPack *mpack; 125747c6ae99SBarry Smith struct DMCompositeLink *xnext,*ynext; 125847c6ae99SBarry Smith struct MatPackLink *anext; 125947c6ae99SBarry Smith PetscScalar *xarray,*yarray; 126047c6ae99SBarry Smith PetscErrorCode ierr; 126147c6ae99SBarry Smith PetscInt i; 126247c6ae99SBarry Smith Vec xglobal,yglobal; 126347c6ae99SBarry Smith PetscMPIInt rank; 126447c6ae99SBarry Smith DM_Composite *comright; 126547c6ae99SBarry Smith DM_Composite *comleft; 126647c6ae99SBarry Smith 126747c6ae99SBarry Smith PetscFunctionBegin; 126847c6ae99SBarry Smith ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr); 126947c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr); 127047c6ae99SBarry Smith comright = (DM_Composite*)mpack->right->data; 127147c6ae99SBarry Smith comleft = (DM_Composite*)mpack->left->data; 127247c6ae99SBarry Smith xnext = comright->next; 127347c6ae99SBarry Smith ynext = comleft->next; 127447c6ae99SBarry Smith anext = mpack->next; 127547c6ae99SBarry Smith 127647c6ae99SBarry Smith while (xnext) { 127747c6ae99SBarry Smith if (xnext->type == DMCOMPOSITE_ARRAY) { 127847c6ae99SBarry Smith if (rank == xnext->rank) { 127947c6ae99SBarry Smith ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 128047c6ae99SBarry Smith ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 128147c6ae99SBarry Smith if (add) { 128247c6ae99SBarry Smith for (i=0; i<xnext->n; i++) { 128347c6ae99SBarry Smith yarray[ynext->rstart+i] += xarray[xnext->rstart+i]; 128447c6ae99SBarry Smith } 128547c6ae99SBarry Smith } else { 128647c6ae99SBarry Smith ierr = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr); 128747c6ae99SBarry Smith } 128847c6ae99SBarry Smith ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 128947c6ae99SBarry Smith ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 129047c6ae99SBarry Smith } 129147c6ae99SBarry Smith } else if (xnext->type == DMCOMPOSITE_DM) { 129247c6ae99SBarry Smith ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 129347c6ae99SBarry Smith ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 129447c6ae99SBarry Smith ierr = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr); 129547c6ae99SBarry Smith ierr = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr); 129647c6ae99SBarry Smith ierr = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr); 129747c6ae99SBarry Smith ierr = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr); 129847c6ae99SBarry Smith if (add) { 129947c6ae99SBarry Smith ierr = MatMultAdd(anext->A,xglobal,yglobal,yglobal);CHKERRQ(ierr); 130047c6ae99SBarry Smith } else { 130147c6ae99SBarry Smith ierr = MatMult(anext->A,xglobal,yglobal);CHKERRQ(ierr); 130247c6ae99SBarry Smith } 130347c6ae99SBarry Smith ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 130447c6ae99SBarry Smith ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 130547c6ae99SBarry Smith ierr = VecResetArray(xglobal);CHKERRQ(ierr); 130647c6ae99SBarry Smith ierr = VecResetArray(yglobal);CHKERRQ(ierr); 130747c6ae99SBarry Smith ierr = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr); 130847c6ae99SBarry Smith ierr = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr); 130947c6ae99SBarry Smith anext = anext->next; 131047c6ae99SBarry Smith } else { 131147c6ae99SBarry Smith SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 131247c6ae99SBarry Smith } 131347c6ae99SBarry Smith xnext = xnext->next; 131447c6ae99SBarry Smith ynext = ynext->next; 131547c6ae99SBarry Smith } 131647c6ae99SBarry Smith PetscFunctionReturn(0); 131747c6ae99SBarry Smith } 131847c6ae99SBarry Smith 131947c6ae99SBarry Smith #undef __FUNCT__ 132047c6ae99SBarry Smith #define __FUNCT__ "MatMultAdd_Shell_Pack" 132147c6ae99SBarry Smith PetscErrorCode MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z) 132247c6ae99SBarry Smith { 132347c6ae99SBarry Smith PetscErrorCode ierr; 132447c6ae99SBarry Smith PetscFunctionBegin; 132547c6ae99SBarry Smith if (z != y) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Handles y == z only"); 132647c6ae99SBarry Smith ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);CHKERRQ(ierr); 132747c6ae99SBarry Smith PetscFunctionReturn(0); 132847c6ae99SBarry Smith } 132947c6ae99SBarry Smith 133047c6ae99SBarry Smith #undef __FUNCT__ 133147c6ae99SBarry Smith #define __FUNCT__ "MatMult_Shell_Pack" 133247c6ae99SBarry Smith PetscErrorCode MatMult_Shell_Pack(Mat A,Vec x,Vec y) 133347c6ae99SBarry Smith { 133447c6ae99SBarry Smith PetscErrorCode ierr; 133547c6ae99SBarry Smith PetscFunctionBegin; 133647c6ae99SBarry Smith ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);CHKERRQ(ierr); 133747c6ae99SBarry Smith PetscFunctionReturn(0); 133847c6ae99SBarry Smith } 133947c6ae99SBarry Smith 134047c6ae99SBarry Smith #undef __FUNCT__ 134147c6ae99SBarry Smith #define __FUNCT__ "MatMultTranspose_Shell_Pack" 134247c6ae99SBarry Smith PetscErrorCode MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y) 134347c6ae99SBarry Smith { 134447c6ae99SBarry Smith struct MatPack *mpack; 134547c6ae99SBarry Smith struct DMCompositeLink *xnext,*ynext; 134647c6ae99SBarry Smith struct MatPackLink *anext; 134747c6ae99SBarry Smith PetscScalar *xarray,*yarray; 134847c6ae99SBarry Smith PetscErrorCode ierr; 134947c6ae99SBarry Smith Vec xglobal,yglobal; 135047c6ae99SBarry Smith PetscMPIInt rank; 135147c6ae99SBarry Smith DM_Composite *comright; 135247c6ae99SBarry Smith DM_Composite *comleft; 135347c6ae99SBarry Smith 135447c6ae99SBarry Smith PetscFunctionBegin; 135547c6ae99SBarry Smith ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr); 135647c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr); 135747c6ae99SBarry Smith comright = (DM_Composite*)mpack->right->data; 135847c6ae99SBarry Smith comleft = (DM_Composite*)mpack->left->data; 135947c6ae99SBarry Smith ynext = comright->next; 136047c6ae99SBarry Smith xnext = comleft->next; 136147c6ae99SBarry Smith anext = mpack->next; 136247c6ae99SBarry Smith 136347c6ae99SBarry Smith while (xnext) { 136447c6ae99SBarry Smith if (xnext->type == DMCOMPOSITE_ARRAY) { 136547c6ae99SBarry Smith if (rank == ynext->rank) { 136647c6ae99SBarry Smith ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 136747c6ae99SBarry Smith ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 136847c6ae99SBarry Smith ierr = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr); 136947c6ae99SBarry Smith ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 137047c6ae99SBarry Smith ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 137147c6ae99SBarry Smith } 137247c6ae99SBarry Smith } else if (xnext->type == DMCOMPOSITE_DM) { 137347c6ae99SBarry Smith ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 137447c6ae99SBarry Smith ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 137547c6ae99SBarry Smith ierr = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr); 137647c6ae99SBarry Smith ierr = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr); 137747c6ae99SBarry Smith ierr = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr); 137847c6ae99SBarry Smith ierr = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr); 137947c6ae99SBarry Smith ierr = MatMultTranspose(anext->A,xglobal,yglobal);CHKERRQ(ierr); 138047c6ae99SBarry Smith ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 138147c6ae99SBarry Smith ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 138247c6ae99SBarry Smith ierr = VecResetArray(xglobal);CHKERRQ(ierr); 138347c6ae99SBarry Smith ierr = VecResetArray(yglobal);CHKERRQ(ierr); 138447c6ae99SBarry Smith ierr = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr); 138547c6ae99SBarry Smith ierr = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr); 138647c6ae99SBarry Smith anext = anext->next; 138747c6ae99SBarry Smith } else { 138847c6ae99SBarry Smith SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 138947c6ae99SBarry Smith } 139047c6ae99SBarry Smith xnext = xnext->next; 139147c6ae99SBarry Smith ynext = ynext->next; 139247c6ae99SBarry Smith } 139347c6ae99SBarry Smith PetscFunctionReturn(0); 139447c6ae99SBarry Smith } 139547c6ae99SBarry Smith 139647c6ae99SBarry Smith #undef __FUNCT__ 139747c6ae99SBarry Smith #define __FUNCT__ "MatDestroy_Shell_Pack" 139847c6ae99SBarry Smith PetscErrorCode MatDestroy_Shell_Pack(Mat A) 139947c6ae99SBarry Smith { 140047c6ae99SBarry Smith struct MatPack *mpack; 140147c6ae99SBarry Smith struct MatPackLink *anext,*oldanext; 140247c6ae99SBarry Smith PetscErrorCode ierr; 140347c6ae99SBarry Smith 140447c6ae99SBarry Smith PetscFunctionBegin; 140547c6ae99SBarry Smith ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr); 140647c6ae99SBarry Smith anext = mpack->next; 140747c6ae99SBarry Smith 140847c6ae99SBarry Smith while (anext) { 140947c6ae99SBarry Smith ierr = MatDestroy(anext->A);CHKERRQ(ierr); 141047c6ae99SBarry Smith oldanext = anext; 141147c6ae99SBarry Smith anext = anext->next; 141247c6ae99SBarry Smith ierr = PetscFree(oldanext);CHKERRQ(ierr); 141347c6ae99SBarry Smith } 141447c6ae99SBarry Smith ierr = PetscFree(mpack);CHKERRQ(ierr); 141547c6ae99SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 141647c6ae99SBarry Smith PetscFunctionReturn(0); 141747c6ae99SBarry Smith } 141847c6ae99SBarry Smith 141947c6ae99SBarry Smith #undef __FUNCT__ 1420*0c010503SBarry Smith #define __FUNCT__ "DMGetInterpolation_Composite" 1421*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v) 142247c6ae99SBarry Smith { 142347c6ae99SBarry Smith PetscErrorCode ierr; 142447c6ae99SBarry Smith PetscInt m,n,M,N; 142547c6ae99SBarry Smith struct DMCompositeLink *nextc; 142647c6ae99SBarry Smith struct DMCompositeLink *nextf; 142747c6ae99SBarry Smith struct MatPackLink *nextmat,*pnextmat = 0; 142847c6ae99SBarry Smith struct MatPack *mpack; 142947c6ae99SBarry Smith Vec gcoarse,gfine; 143047c6ae99SBarry Smith DM_Composite *comcoarse = (DM_Composite*)coarse->data; 143147c6ae99SBarry Smith DM_Composite *comfine = (DM_Composite*)fine->data; 143247c6ae99SBarry Smith 143347c6ae99SBarry Smith PetscFunctionBegin; 143447c6ae99SBarry Smith PetscValidHeaderSpecific(coarse,DM_CLASSID,1); 143547c6ae99SBarry Smith PetscValidHeaderSpecific(fine,DM_CLASSID,2); 143647c6ae99SBarry Smith nextc = comcoarse->next; 143747c6ae99SBarry Smith nextf = comfine->next; 143847c6ae99SBarry Smith /* use global vectors only for determining matrix layout */ 1439*0c010503SBarry Smith ierr = DMCreateGlobalVector(coarse,&gcoarse);CHKERRQ(ierr); 1440*0c010503SBarry Smith ierr = DMCreateGlobalVector(fine,&gfine);CHKERRQ(ierr); 144147c6ae99SBarry Smith ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr); 144247c6ae99SBarry Smith ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr); 144347c6ae99SBarry Smith ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr); 144447c6ae99SBarry Smith ierr = VecGetSize(gfine,&M);CHKERRQ(ierr); 144547c6ae99SBarry Smith ierr = VecDestroy(gcoarse);CHKERRQ(ierr); 144647c6ae99SBarry Smith ierr = VecDestroy(gfine);CHKERRQ(ierr); 144747c6ae99SBarry Smith 144847c6ae99SBarry Smith ierr = PetscNew(struct MatPack,&mpack);CHKERRQ(ierr); 144947c6ae99SBarry Smith mpack->right = coarse; 145047c6ae99SBarry Smith mpack->left = fine; 145147c6ae99SBarry Smith ierr = MatCreate(((PetscObject)fine)->comm,A);CHKERRQ(ierr); 145247c6ae99SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 145347c6ae99SBarry Smith ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr); 145447c6ae99SBarry Smith ierr = MatShellSetContext(*A,mpack);CHKERRQ(ierr); 145547c6ae99SBarry Smith ierr = MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);CHKERRQ(ierr); 145647c6ae99SBarry Smith ierr = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);CHKERRQ(ierr); 145747c6ae99SBarry Smith ierr = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);CHKERRQ(ierr); 145847c6ae99SBarry Smith ierr = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);CHKERRQ(ierr); 145947c6ae99SBarry Smith 146047c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 146147c6ae99SBarry Smith while (nextc) { 146247c6ae99SBarry Smith if (nextc->type != nextf->type) SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Two DM have different layout"); 146347c6ae99SBarry Smith 146447c6ae99SBarry Smith if (nextc->type == DMCOMPOSITE_ARRAY) { 146547c6ae99SBarry Smith ; 146647c6ae99SBarry Smith } else if (nextc->type == DMCOMPOSITE_DM) { 146747c6ae99SBarry Smith ierr = PetscNew(struct MatPackLink,&nextmat);CHKERRQ(ierr); 146847c6ae99SBarry Smith nextmat->next = 0; 146947c6ae99SBarry Smith if (pnextmat) { 147047c6ae99SBarry Smith pnextmat->next = nextmat; 147147c6ae99SBarry Smith pnextmat = nextmat; 147247c6ae99SBarry Smith } else { 147347c6ae99SBarry Smith pnextmat = nextmat; 147447c6ae99SBarry Smith mpack->next = nextmat; 147547c6ae99SBarry Smith } 147647c6ae99SBarry Smith ierr = DMGetInterpolation(nextc->dm,nextf->dm,&nextmat->A,PETSC_NULL);CHKERRQ(ierr); 147747c6ae99SBarry Smith } else { 147847c6ae99SBarry Smith SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 147947c6ae99SBarry Smith } 148047c6ae99SBarry Smith nextc = nextc->next; 148147c6ae99SBarry Smith nextf = nextf->next; 148247c6ae99SBarry Smith } 148347c6ae99SBarry Smith PetscFunctionReturn(0); 148447c6ae99SBarry Smith } 148547c6ae99SBarry Smith 148647c6ae99SBarry Smith #undef __FUNCT__ 1487*0c010503SBarry Smith #define __FUNCT__ "DMGetMatrix_Composite" 1488*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix_Composite(DM dm, const MatType mtype,Mat *J) 148947c6ae99SBarry Smith { 149047c6ae99SBarry Smith PetscErrorCode ierr; 149147c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 149247c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 149347c6ae99SBarry Smith PetscInt m,*dnz,*onz,i,j,mA; 149447c6ae99SBarry Smith Mat Atmp; 149547c6ae99SBarry Smith PetscMPIInt rank; 149647c6ae99SBarry Smith PetscScalar zero = 0.0; 149747c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 149847c6ae99SBarry Smith 149947c6ae99SBarry Smith PetscFunctionBegin; 150047c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 150147c6ae99SBarry Smith 150247c6ae99SBarry Smith /* use global vector to determine layout needed for matrix */ 150347c6ae99SBarry Smith m = com->n; 150447c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 150547c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr); 150647c6ae99SBarry Smith ierr = MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 150747c6ae99SBarry Smith ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 150847c6ae99SBarry Smith 150947c6ae99SBarry Smith /* 151047c6ae99SBarry Smith Extremely inefficient but will compute entire Jacobian for testing 151147c6ae99SBarry Smith */ 151247c6ae99SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr); 151347c6ae99SBarry Smith if (dense) { 151447c6ae99SBarry Smith PetscInt rstart,rend,*indices; 151547c6ae99SBarry Smith PetscScalar *values; 151647c6ae99SBarry Smith 151747c6ae99SBarry Smith mA = com->N; 151847c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,mA,PETSC_NULL,mA-m,PETSC_NULL);CHKERRQ(ierr); 151947c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,mA,PETSC_NULL);CHKERRQ(ierr); 152047c6ae99SBarry Smith 152147c6ae99SBarry Smith ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr); 152247c6ae99SBarry Smith ierr = PetscMalloc2(mA,PetscScalar,&values,mA,PetscInt,&indices);CHKERRQ(ierr); 152347c6ae99SBarry Smith ierr = PetscMemzero(values,mA*sizeof(PetscScalar));CHKERRQ(ierr); 152447c6ae99SBarry Smith for (i=0; i<mA; i++) indices[i] = i; 152547c6ae99SBarry Smith for (i=rstart; i<rend; i++) { 152647c6ae99SBarry Smith ierr = MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES);CHKERRQ(ierr); 152747c6ae99SBarry Smith } 152847c6ae99SBarry Smith ierr = PetscFree2(values,indices);CHKERRQ(ierr); 152947c6ae99SBarry Smith ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 153047c6ae99SBarry Smith ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 153147c6ae99SBarry Smith PetscFunctionReturn(0); 153247c6ae99SBarry Smith } 153347c6ae99SBarry Smith 153447c6ae99SBarry Smith ierr = MatPreallocateInitialize(((PetscObject)dm)->comm,m,m,dnz,onz);CHKERRQ(ierr); 153547c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 153647c6ae99SBarry Smith next = com->next; 153747c6ae99SBarry Smith while (next) { 153847c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 153947c6ae99SBarry Smith if (rank == next->rank) { /* zero the "little" block */ 154047c6ae99SBarry Smith for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) { 154147c6ae99SBarry Smith for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) { 154247c6ae99SBarry Smith ierr = MatPreallocateSet(j,1,&i,dnz,onz);CHKERRQ(ierr); 154347c6ae99SBarry Smith } 154447c6ae99SBarry Smith } 154547c6ae99SBarry Smith } 154647c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 154747c6ae99SBarry Smith PetscInt nc,rstart,*ccols,maxnc; 154847c6ae99SBarry Smith const PetscInt *cols,*rstarts; 154947c6ae99SBarry Smith PetscMPIInt proc; 155047c6ae99SBarry Smith 155147c6ae99SBarry Smith ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr); 155247c6ae99SBarry Smith ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr); 155347c6ae99SBarry Smith ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr); 155447c6ae99SBarry Smith ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr); 155547c6ae99SBarry Smith 155647c6ae99SBarry Smith maxnc = 0; 155747c6ae99SBarry Smith for (i=0; i<mA; i++) { 155847c6ae99SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 155947c6ae99SBarry Smith ierr = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 156047c6ae99SBarry Smith maxnc = PetscMax(nc,maxnc); 156147c6ae99SBarry Smith } 156247c6ae99SBarry Smith ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr); 156347c6ae99SBarry Smith for (i=0; i<mA; i++) { 156447c6ae99SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr); 156547c6ae99SBarry Smith /* remap the columns taking into how much they are shifted on each process */ 156647c6ae99SBarry Smith for (j=0; j<nc; j++) { 156747c6ae99SBarry Smith proc = 0; 156847c6ae99SBarry Smith while (cols[j] >= rstarts[proc+1]) proc++; 156947c6ae99SBarry Smith ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 157047c6ae99SBarry Smith } 157147c6ae99SBarry Smith ierr = MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz);CHKERRQ(ierr); 157247c6ae99SBarry Smith ierr = MatRestoreRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr); 157347c6ae99SBarry Smith } 157447c6ae99SBarry Smith ierr = PetscFree(ccols);CHKERRQ(ierr); 157547c6ae99SBarry Smith ierr = MatDestroy(Atmp);CHKERRQ(ierr); 157647c6ae99SBarry Smith } else { 157747c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 157847c6ae99SBarry Smith } 157947c6ae99SBarry Smith next = next->next; 158047c6ae99SBarry Smith } 158147c6ae99SBarry Smith if (com->FormCoupleLocations) { 158247c6ae99SBarry Smith ierr = (*com->FormCoupleLocations)(dm,PETSC_NULL,dnz,onz,__rstart,__nrows,__start,__end);CHKERRQ(ierr); 158347c6ae99SBarry Smith } 158447c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);CHKERRQ(ierr); 158547c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,0,dnz);CHKERRQ(ierr); 158647c6ae99SBarry Smith ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 158747c6ae99SBarry Smith 158847c6ae99SBarry Smith next = com->next; 158947c6ae99SBarry Smith while (next) { 159047c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 159147c6ae99SBarry Smith if (rank == next->rank) { 159247c6ae99SBarry Smith for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) { 159347c6ae99SBarry Smith for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) { 159447c6ae99SBarry Smith ierr = MatSetValues(*J,1,&j,1,&i,&zero,INSERT_VALUES);CHKERRQ(ierr); 159547c6ae99SBarry Smith } 159647c6ae99SBarry Smith } 159747c6ae99SBarry Smith } 159847c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 159947c6ae99SBarry Smith PetscInt nc,rstart,row,maxnc,*ccols; 160047c6ae99SBarry Smith const PetscInt *cols,*rstarts; 160147c6ae99SBarry Smith const PetscScalar *values; 160247c6ae99SBarry Smith PetscMPIInt proc; 160347c6ae99SBarry Smith 160447c6ae99SBarry Smith ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr); 160547c6ae99SBarry Smith ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr); 160647c6ae99SBarry Smith ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr); 160747c6ae99SBarry Smith ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr); 160847c6ae99SBarry Smith maxnc = 0; 160947c6ae99SBarry Smith for (i=0; i<mA; i++) { 161047c6ae99SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 161147c6ae99SBarry Smith ierr = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 161247c6ae99SBarry Smith maxnc = PetscMax(nc,maxnc); 161347c6ae99SBarry Smith } 161447c6ae99SBarry Smith ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr); 161547c6ae99SBarry Smith for (i=0; i<mA; i++) { 161647c6ae99SBarry Smith ierr = MatGetRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr); 161747c6ae99SBarry Smith for (j=0; j<nc; j++) { 161847c6ae99SBarry Smith proc = 0; 161947c6ae99SBarry Smith while (cols[j] >= rstarts[proc+1]) proc++; 162047c6ae99SBarry Smith ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc]; 162147c6ae99SBarry Smith } 162247c6ae99SBarry Smith row = com->rstart+next->rstart+i; 162347c6ae99SBarry Smith ierr = MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES);CHKERRQ(ierr); 162447c6ae99SBarry Smith ierr = MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr); 162547c6ae99SBarry Smith } 162647c6ae99SBarry Smith ierr = PetscFree(ccols);CHKERRQ(ierr); 162747c6ae99SBarry Smith ierr = MatDestroy(Atmp);CHKERRQ(ierr); 162847c6ae99SBarry Smith } else { 162947c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 163047c6ae99SBarry Smith } 163147c6ae99SBarry Smith next = next->next; 163247c6ae99SBarry Smith } 163347c6ae99SBarry Smith if (com->FormCoupleLocations) { 163447c6ae99SBarry Smith PetscInt __rstart; 163547c6ae99SBarry Smith ierr = MatGetOwnershipRange(*J,&__rstart,PETSC_NULL);CHKERRQ(ierr); 163647c6ae99SBarry Smith ierr = (*com->FormCoupleLocations)(dm,*J,PETSC_NULL,PETSC_NULL,__rstart,0,0,0);CHKERRQ(ierr); 163747c6ae99SBarry Smith } 163847c6ae99SBarry Smith ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 163947c6ae99SBarry Smith ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164047c6ae99SBarry Smith PetscFunctionReturn(0); 164147c6ae99SBarry Smith } 164247c6ae99SBarry Smith 164347c6ae99SBarry Smith #undef __FUNCT__ 1644*0c010503SBarry Smith #define __FUNCT__ "DMGetColoring_Composite" 1645*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring) 164647c6ae99SBarry Smith { 164747c6ae99SBarry Smith PetscErrorCode ierr; 164847c6ae99SBarry Smith PetscInt n,i,cnt; 164947c6ae99SBarry Smith ISColoringValue *colors; 165047c6ae99SBarry Smith PetscBool dense = PETSC_FALSE; 165147c6ae99SBarry Smith ISColoringValue maxcol = 0; 165247c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 165347c6ae99SBarry Smith 165447c6ae99SBarry Smith PetscFunctionBegin; 165547c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 165647c6ae99SBarry Smith if (ctype == IS_COLORING_GHOSTED) { 165747c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" ); 165847c6ae99SBarry Smith } else if (ctype == IS_COLORING_GLOBAL) { 165947c6ae99SBarry Smith n = com->n; 166047c6ae99SBarry Smith } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType"); 166147c6ae99SBarry Smith ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */ 166247c6ae99SBarry Smith 166347c6ae99SBarry Smith ierr = PetscOptionsGetTruth(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr); 166447c6ae99SBarry Smith if (dense) { 166547c6ae99SBarry Smith for (i=0; i<n; i++) { 166647c6ae99SBarry Smith colors[i] = (ISColoringValue)(com->rstart + i); 166747c6ae99SBarry Smith } 166847c6ae99SBarry Smith maxcol = com->N; 166947c6ae99SBarry Smith } else { 167047c6ae99SBarry Smith struct DMCompositeLink *next = com->next; 167147c6ae99SBarry Smith PetscMPIInt rank; 167247c6ae99SBarry Smith 167347c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 167447c6ae99SBarry Smith cnt = 0; 167547c6ae99SBarry Smith while (next) { 167647c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 167747c6ae99SBarry Smith if (rank == next->rank) { /* each column gets is own color */ 167847c6ae99SBarry Smith for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) { 167947c6ae99SBarry Smith colors[cnt++] = maxcol++; 168047c6ae99SBarry Smith } 168147c6ae99SBarry Smith } 168247c6ae99SBarry Smith ierr = MPI_Bcast(&maxcol,1,MPIU_COLORING_VALUE,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr); 168347c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 168447c6ae99SBarry Smith ISColoring lcoloring; 168547c6ae99SBarry Smith 168647c6ae99SBarry Smith ierr = DMGetColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr); 168747c6ae99SBarry Smith for (i=0; i<lcoloring->N; i++) { 168847c6ae99SBarry Smith colors[cnt++] = maxcol + lcoloring->colors[i]; 168947c6ae99SBarry Smith } 169047c6ae99SBarry Smith maxcol += lcoloring->n; 169147c6ae99SBarry Smith ierr = ISColoringDestroy(lcoloring);CHKERRQ(ierr); 169247c6ae99SBarry Smith } else { 169347c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 169447c6ae99SBarry Smith } 169547c6ae99SBarry Smith next = next->next; 169647c6ae99SBarry Smith } 169747c6ae99SBarry Smith } 169847c6ae99SBarry Smith ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr); 169947c6ae99SBarry Smith PetscFunctionReturn(0); 170047c6ae99SBarry Smith } 170147c6ae99SBarry Smith 170247c6ae99SBarry Smith #undef __FUNCT__ 1703*0c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin_Composite" 1704*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec) 170547c6ae99SBarry Smith { 170647c6ae99SBarry Smith PetscErrorCode ierr; 170747c6ae99SBarry Smith struct DMCompositeLink *next; 170847c6ae99SBarry Smith PetscInt cnt = 3; 170947c6ae99SBarry Smith PetscMPIInt rank; 171047c6ae99SBarry Smith PetscScalar *garray,*larray; 171147c6ae99SBarry Smith DM_Composite *com = (DM_Composite*)dm->data; 171247c6ae99SBarry Smith 171347c6ae99SBarry Smith PetscFunctionBegin; 171447c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 171547c6ae99SBarry Smith PetscValidHeaderSpecific(gvec,VEC_CLASSID,2); 171647c6ae99SBarry Smith next = com->next; 171747c6ae99SBarry Smith if (!com->setup) { 171847c6ae99SBarry Smith ierr = DMCompositeSetUp(dm);CHKERRQ(ierr); 171947c6ae99SBarry Smith } 172047c6ae99SBarry Smith ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr); 172147c6ae99SBarry Smith ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr); 172247c6ae99SBarry Smith ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr); 172347c6ae99SBarry Smith 172447c6ae99SBarry Smith /* loop over packed objects, handling one at at time */ 172547c6ae99SBarry Smith while (next) { 172647c6ae99SBarry Smith if (next->type == DMCOMPOSITE_ARRAY) { 172747c6ae99SBarry Smith if (rank == next->rank) { 172847c6ae99SBarry Smith ierr = PetscMemcpy(larray,garray,next->n*sizeof(PetscScalar));CHKERRQ(ierr); 172947c6ae99SBarry Smith garray += next->n; 173047c6ae99SBarry Smith } 173147c6ae99SBarry Smith /* does not handle ADD_VALUES */ 173247c6ae99SBarry Smith ierr = MPI_Bcast(larray,next->n,MPIU_SCALAR,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr); 173347c6ae99SBarry Smith larray += next->n; 173447c6ae99SBarry Smith } else if (next->type == DMCOMPOSITE_DM) { 173547c6ae99SBarry Smith Vec local,global; 173647c6ae99SBarry Smith PetscInt N; 173747c6ae99SBarry Smith 173847c6ae99SBarry Smith ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr); 173947c6ae99SBarry Smith ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr); 174047c6ae99SBarry Smith ierr = VecPlaceArray(global,garray);CHKERRQ(ierr); 174147c6ae99SBarry Smith ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr); 174247c6ae99SBarry Smith ierr = VecPlaceArray(local,larray);CHKERRQ(ierr); 174347c6ae99SBarry Smith ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr); 174447c6ae99SBarry Smith ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr); 174547c6ae99SBarry Smith ierr = VecResetArray(global);CHKERRQ(ierr); 174647c6ae99SBarry Smith ierr = VecResetArray(local);CHKERRQ(ierr); 174747c6ae99SBarry Smith ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr); 174847c6ae99SBarry Smith ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr); 174947c6ae99SBarry Smith larray += next->n; 175047c6ae99SBarry Smith } else { 175147c6ae99SBarry Smith SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet"); 175247c6ae99SBarry Smith } 175347c6ae99SBarry Smith cnt++; 175447c6ae99SBarry Smith next = next->next; 175547c6ae99SBarry Smith } 175647c6ae99SBarry Smith 175747c6ae99SBarry Smith ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr); 175847c6ae99SBarry Smith ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr); 175947c6ae99SBarry Smith PetscFunctionReturn(0); 176047c6ae99SBarry Smith } 176147c6ae99SBarry Smith 176247c6ae99SBarry Smith #undef __FUNCT__ 1763*0c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd_Composite" 1764*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec) 1765*0c010503SBarry Smith { 1766*0c010503SBarry Smith PetscFunctionBegin; 1767*0c010503SBarry Smith PetscFunctionReturn(0); 1768*0c010503SBarry Smith } 176947c6ae99SBarry Smith 1770*0c010503SBarry Smith #undef __FUNCT__ 1771*0c010503SBarry Smith #define __FUNCT__ "DMCompositeCreate" 1772*0c010503SBarry Smith /*@C 1773*0c010503SBarry Smith DMCompositeCreate - Creates a vector packer, used to generate "composite" 1774*0c010503SBarry Smith vectors made up of several subvectors. 1775*0c010503SBarry Smith 1776*0c010503SBarry Smith Collective on MPI_Comm 177747c6ae99SBarry Smith 177847c6ae99SBarry Smith Input Parameter: 1779*0c010503SBarry Smith . comm - the processors that will share the global vector 1780*0c010503SBarry Smith 1781*0c010503SBarry Smith Output Parameters: 1782*0c010503SBarry Smith . packer - the packer object 178347c6ae99SBarry Smith 178447c6ae99SBarry Smith Level: advanced 178547c6ae99SBarry Smith 1786*0c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCompositeScatter(), 1787*0c010503SBarry Smith DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetGlobalIndices(), DMCompositeGetAccess() 178847c6ae99SBarry Smith DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries() 178947c6ae99SBarry Smith 179047c6ae99SBarry Smith @*/ 1791*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeCreate(MPI_Comm comm,DM *packer) 179247c6ae99SBarry Smith { 1793*0c010503SBarry Smith PetscErrorCode ierr; 1794*0c010503SBarry Smith DM p; 1795*0c010503SBarry Smith DM_Composite *com; 1796*0c010503SBarry Smith 179747c6ae99SBarry Smith PetscFunctionBegin; 1798*0c010503SBarry Smith PetscValidPointer(packer,2); 1799*0c010503SBarry Smith *packer = PETSC_NULL; 1800*0c010503SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES 1801*0c010503SBarry Smith ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr); 1802*0c010503SBarry Smith #endif 1803*0c010503SBarry Smith 1804*0c010503SBarry Smith ierr = PetscHeaderCreate(p,_p_DM,struct _DMOps,DM_CLASSID,0,"DM",comm,DMDestroy,DMView);CHKERRQ(ierr); 1805*0c010503SBarry Smith ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr); 1806*0c010503SBarry Smith p->data = com; 1807*0c010503SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr); 1808*0c010503SBarry Smith com->n = 0; 1809*0c010503SBarry Smith com->next = PETSC_NULL; 1810*0c010503SBarry Smith com->nredundant = 0; 1811*0c010503SBarry Smith com->nDM = 0; 1812*0c010503SBarry Smith 1813*0c010503SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Composite; 1814*0c010503SBarry Smith p->ops->createlocalvector = DMCreateLocalVector_Composite; 1815*0c010503SBarry Smith p->ops->refine = DMRefine_Composite; 1816*0c010503SBarry Smith p->ops->getinterpolation = DMGetInterpolation_Composite; 1817*0c010503SBarry Smith p->ops->getmatrix = DMGetMatrix_Composite; 1818*0c010503SBarry Smith p->ops->getcoloring = DMGetColoring_Composite; 1819*0c010503SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite; 1820*0c010503SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Composite; 1821*0c010503SBarry Smith p->ops->destroy = DMDestroy_Composite; 1822*0c010503SBarry Smith p->ops->view = DMView_Composite; 1823*0c010503SBarry Smith 1824*0c010503SBarry Smith *packer = p; 182547c6ae99SBarry Smith PetscFunctionReturn(0); 182647c6ae99SBarry Smith } 1827