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