xref: /petsc/src/dm/impls/composite/pack.c (revision 520db06cfec378d7c34fb86a68534e48a6477cdf)
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"
329acc1e270SJed Brown PetscErrorCode DMCompositeGather_Array(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,const 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()");
340acc1e270SJed Brown   }
341df0c820aSJed Brown   switch (imode) {
342df0c820aSJed Brown   case INSERT_VALUES:
343acc1e270SJed Brown     if (rank == mine->rank) {
34447c6ae99SBarry Smith       ierr = PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(PetscScalar));CHKERRQ(ierr);
345acc1e270SJed Brown     }
346df0c820aSJed Brown     break;
347df0c820aSJed Brown   case ADD_VALUES: {
348df0c820aSJed Brown     PetscInt          i;
349acc1e270SJed Brown     const PetscScalar *source;
350acc1e270SJed Brown     PetscScalar       *buffer,*dest;
351acc1e270SJed Brown     if (rank == mine->rank) {
352acc1e270SJed Brown       dest = &varray[mine->rstart];
353acc1e270SJed Brown #if defined(PETSC_HAVE_MPI_IN_PLACE)
354acc1e270SJed Brown       buffer = dest;
355acc1e270SJed Brown       source = MPI_IN_PLACE;
356acc1e270SJed Brown #else
357acc1e270SJed Brown       ierr = PetscMalloc(mine->n*sizeof(PetscScalar),&buffer);CHKERRQ(ierr);
358acc1e270SJed Brown       source = buffer;
359acc1e270SJed Brown #endif
360acc1e270SJed Brown       for (i=0; i<mine->n; i++) buffer[i] = varray[mine->rstart+i] + array[i];
361acc1e270SJed Brown     } else {
362acc1e270SJed Brown       source = array;
363acc1e270SJed Brown       dest   = PETSC_NULL;
364acc1e270SJed Brown     }
365acc1e270SJed Brown     ierr = MPI_Reduce((void*)source,dest,mine->n,MPIU_SCALAR,MPI_SUM,mine->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
366acc1e270SJed Brown #if !defined(PETSC_HAVE_MPI_IN_PLACE)
367acc1e270SJed Brown     if (rank == mine->rank) {ierr = PetscFree(source);CHKERRQ(ierr);}
368acc1e270SJed Brown #endif
369df0c820aSJed Brown   } break;
370df0c820aSJed Brown   default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"imode");
371df0c820aSJed Brown   }
372acc1e270SJed Brown   if (rank == mine->rank) {ierr = VecRestoreArray(vec,&varray);CHKERRQ(ierr);}
37347c6ae99SBarry Smith   PetscFunctionReturn(0);
37447c6ae99SBarry Smith }
37547c6ae99SBarry Smith 
37647c6ae99SBarry Smith #undef __FUNCT__
37747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_DM"
378df0c820aSJed Brown PetscErrorCode DMCompositeGather_DM(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,Vec local)
37947c6ae99SBarry Smith {
38047c6ae99SBarry Smith   PetscErrorCode ierr;
38147c6ae99SBarry Smith   PetscScalar    *array;
38247c6ae99SBarry Smith   Vec            global;
38347c6ae99SBarry Smith 
38447c6ae99SBarry Smith   PetscFunctionBegin;
38547c6ae99SBarry Smith   ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr);
38647c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
38747c6ae99SBarry Smith   ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr);
388df0c820aSJed Brown   ierr = DMLocalToGlobalBegin(mine->dm,local,imode,global);CHKERRQ(ierr);
389df0c820aSJed Brown   ierr = DMLocalToGlobalEnd(mine->dm,local,imode,global);CHKERRQ(ierr);
39047c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
39147c6ae99SBarry Smith   ierr = VecResetArray(global);CHKERRQ(ierr);
39247c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr);
39347c6ae99SBarry Smith   PetscFunctionReturn(0);
39447c6ae99SBarry Smith }
39547c6ae99SBarry Smith 
39647c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
39747c6ae99SBarry Smith 
39847c6ae99SBarry Smith #include <stdarg.h>
39947c6ae99SBarry Smith 
40047c6ae99SBarry Smith #undef __FUNCT__
40147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetNumberDM"
40247c6ae99SBarry Smith /*@C
40347c6ae99SBarry Smith     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
40447c6ae99SBarry Smith        representation.
40547c6ae99SBarry Smith 
40647c6ae99SBarry Smith     Not Collective
40747c6ae99SBarry Smith 
40847c6ae99SBarry Smith     Input Parameter:
40947c6ae99SBarry Smith .    dm - the packer object
41047c6ae99SBarry Smith 
41147c6ae99SBarry Smith     Output Parameter:
41247c6ae99SBarry Smith .     nDM - the number of DMs
41347c6ae99SBarry Smith 
41447c6ae99SBarry Smith     Level: beginner
41547c6ae99SBarry Smith 
41647c6ae99SBarry Smith @*/
41747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
41847c6ae99SBarry Smith {
41947c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
42047c6ae99SBarry Smith   PetscFunctionBegin;
42147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
42247c6ae99SBarry Smith   *nDM = com->nDM;
42347c6ae99SBarry Smith   PetscFunctionReturn(0);
42447c6ae99SBarry Smith }
42547c6ae99SBarry Smith 
42647c6ae99SBarry Smith 
42747c6ae99SBarry Smith #undef __FUNCT__
42847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess"
42947c6ae99SBarry Smith /*@C
43047c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
43147c6ae99SBarry Smith        representation.
43247c6ae99SBarry Smith 
43347c6ae99SBarry Smith     Collective on DMComposite
43447c6ae99SBarry Smith 
43547c6ae99SBarry Smith     Input Parameter:
43647c6ae99SBarry Smith +    dm - the packer object
43747c6ae99SBarry Smith .    gvec - the global vector
43847c6ae99SBarry Smith -    ... - the individual sequential or parallel objects (arrays or vectors)
43947c6ae99SBarry Smith 
44047c6ae99SBarry Smith     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
44147c6ae99SBarry Smith 
44247c6ae99SBarry Smith     Level: advanced
44347c6ae99SBarry Smith 
44447c6ae99SBarry Smith @*/
44547c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetAccess(DM dm,Vec gvec,...)
44647c6ae99SBarry Smith {
44747c6ae99SBarry Smith   va_list                Argp;
44847c6ae99SBarry Smith   PetscErrorCode         ierr;
44947c6ae99SBarry Smith   struct DMCompositeLink *next;
45047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
45147c6ae99SBarry Smith 
45247c6ae99SBarry Smith   PetscFunctionBegin;
45347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
45447c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
45547c6ae99SBarry Smith   next = com->next;
45647c6ae99SBarry Smith   if (!com->setup) {
457d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
45847c6ae99SBarry Smith   }
45947c6ae99SBarry Smith 
46047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
46147c6ae99SBarry Smith   va_start(Argp,gvec);
46247c6ae99SBarry Smith   while (next) {
46347c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
46447c6ae99SBarry Smith       PetscScalar **array;
46547c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
46647c6ae99SBarry Smith       ierr  = DMCompositeGetAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
46747c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
46847c6ae99SBarry Smith       Vec *vec;
46947c6ae99SBarry Smith       vec  = va_arg(Argp, Vec*);
47047c6ae99SBarry Smith       ierr = DMCompositeGetAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
47147c6ae99SBarry Smith     } else {
47247c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
47347c6ae99SBarry Smith     }
47447c6ae99SBarry Smith     next = next->next;
47547c6ae99SBarry Smith   }
47647c6ae99SBarry Smith   va_end(Argp);
47747c6ae99SBarry Smith   PetscFunctionReturn(0);
47847c6ae99SBarry Smith }
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith #undef __FUNCT__
48147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess"
48247c6ae99SBarry Smith /*@C
483aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
48447c6ae99SBarry Smith        representation.
48547c6ae99SBarry Smith 
48647c6ae99SBarry Smith     Collective on DMComposite
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith     Input Parameter:
48947c6ae99SBarry Smith +    dm - the packer object
49047c6ae99SBarry Smith .    gvec - the global vector
49147c6ae99SBarry Smith -    ... - the individual sequential or parallel objects (arrays or vectors)
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith     Level: advanced
49447c6ae99SBarry Smith 
4950c010503SBarry Smith .seealso  DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
4966eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
497aa219208SBarry Smith          DMCompositeRestoreAccess(), DMCompositeGetAccess()
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith @*/
50047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreAccess(DM dm,Vec gvec,...)
50147c6ae99SBarry Smith {
50247c6ae99SBarry Smith   va_list                Argp;
50347c6ae99SBarry Smith   PetscErrorCode         ierr;
50447c6ae99SBarry Smith   struct DMCompositeLink *next;
50547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
50647c6ae99SBarry Smith 
50747c6ae99SBarry Smith   PetscFunctionBegin;
50847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
50947c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
51047c6ae99SBarry Smith   next = com->next;
51147c6ae99SBarry Smith   if (!com->setup) {
512d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
51347c6ae99SBarry Smith   }
51447c6ae99SBarry Smith 
51547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
51647c6ae99SBarry Smith   va_start(Argp,gvec);
51747c6ae99SBarry Smith   while (next) {
51847c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
51947c6ae99SBarry Smith       PetscScalar **array;
52047c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
52147c6ae99SBarry Smith       ierr  = DMCompositeRestoreAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
52247c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
52347c6ae99SBarry Smith       Vec *vec;
52447c6ae99SBarry Smith       vec  = va_arg(Argp, Vec*);
52547c6ae99SBarry Smith       ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
52647c6ae99SBarry Smith     } else {
52747c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
52847c6ae99SBarry Smith     }
52947c6ae99SBarry Smith     next = next->next;
53047c6ae99SBarry Smith   }
53147c6ae99SBarry Smith   va_end(Argp);
53247c6ae99SBarry Smith   PetscFunctionReturn(0);
53347c6ae99SBarry Smith }
53447c6ae99SBarry Smith 
53547c6ae99SBarry Smith #undef __FUNCT__
53647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter"
53747c6ae99SBarry Smith /*@C
53847c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
53947c6ae99SBarry Smith 
54047c6ae99SBarry Smith     Collective on DMComposite
54147c6ae99SBarry Smith 
54247c6ae99SBarry Smith     Input Parameter:
54347c6ae99SBarry Smith +    dm - the packer object
54447c6ae99SBarry Smith .    gvec - the global vector
54547c6ae99SBarry Smith -    ... - the individual sequential objects (arrays or vectors)
54647c6ae99SBarry Smith 
54747c6ae99SBarry Smith     Level: advanced
54847c6ae99SBarry Smith 
5490c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
5506eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
55147c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
55247c6ae99SBarry Smith 
55347c6ae99SBarry Smith @*/
55447c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeScatter(DM dm,Vec gvec,...)
55547c6ae99SBarry Smith {
55647c6ae99SBarry Smith   va_list                Argp;
55747c6ae99SBarry Smith   PetscErrorCode         ierr;
55847c6ae99SBarry Smith   struct DMCompositeLink *next;
55947c6ae99SBarry Smith   PetscInt               cnt = 3;
56047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith   PetscFunctionBegin;
56347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
56447c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
56547c6ae99SBarry Smith   next = com->next;
56647c6ae99SBarry Smith   if (!com->setup) {
567d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
56847c6ae99SBarry Smith   }
56947c6ae99SBarry Smith 
57047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
57147c6ae99SBarry Smith   va_start(Argp,gvec);
57247c6ae99SBarry Smith   while (next) {
57347c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
57447c6ae99SBarry Smith       PetscScalar *array;
57547c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar*);
57647c6ae99SBarry Smith       ierr = DMCompositeScatter_Array(dm,next,gvec,array);CHKERRQ(ierr);
57747c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
57847c6ae99SBarry Smith       Vec vec;
57947c6ae99SBarry Smith       vec = va_arg(Argp, Vec);
58047c6ae99SBarry Smith       PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt);
58147c6ae99SBarry Smith       ierr = DMCompositeScatter_DM(dm,next,gvec,vec);CHKERRQ(ierr);
58247c6ae99SBarry Smith     } else {
58347c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
58447c6ae99SBarry Smith     }
58547c6ae99SBarry Smith     cnt++;
58647c6ae99SBarry Smith     next = next->next;
58747c6ae99SBarry Smith   }
58847c6ae99SBarry Smith   va_end(Argp);
58947c6ae99SBarry Smith   PetscFunctionReturn(0);
59047c6ae99SBarry Smith }
59147c6ae99SBarry Smith 
59247c6ae99SBarry Smith #undef __FUNCT__
59347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather"
59447c6ae99SBarry Smith /*@C
59547c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
59647c6ae99SBarry Smith 
59747c6ae99SBarry Smith     Collective on DMComposite
59847c6ae99SBarry Smith 
59947c6ae99SBarry Smith     Input Parameter:
60047c6ae99SBarry Smith +    dm - the packer object
60147c6ae99SBarry Smith .    gvec - the global vector
60247c6ae99SBarry Smith -    ... - the individual sequential objects (arrays or vectors)
60347c6ae99SBarry Smith 
60447c6ae99SBarry Smith     Level: advanced
60547c6ae99SBarry Smith 
6060c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
6076eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
60847c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
60947c6ae99SBarry Smith 
61047c6ae99SBarry Smith @*/
611df0c820aSJed Brown PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
61247c6ae99SBarry Smith {
61347c6ae99SBarry Smith   va_list                Argp;
61447c6ae99SBarry Smith   PetscErrorCode         ierr;
61547c6ae99SBarry Smith   struct DMCompositeLink *next;
61647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
61747c6ae99SBarry Smith 
61847c6ae99SBarry Smith   PetscFunctionBegin;
61947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
62047c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
62147c6ae99SBarry Smith   next = com->next;
62247c6ae99SBarry Smith   if (!com->setup) {
623d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
62447c6ae99SBarry Smith   }
62547c6ae99SBarry Smith 
62647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
627df0c820aSJed Brown   va_start(Argp,imode);
62847c6ae99SBarry Smith   while (next) {
62947c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
63047c6ae99SBarry Smith       PetscScalar *array;
63147c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar*);
632df0c820aSJed Brown       ierr  = DMCompositeGather_Array(dm,next,gvec,imode,array);CHKERRQ(ierr);
63347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
63447c6ae99SBarry Smith       Vec vec;
63547c6ae99SBarry Smith       vec = va_arg(Argp, Vec);
63647c6ae99SBarry Smith       PetscValidHeaderSpecific(vec,VEC_CLASSID,3);
637df0c820aSJed Brown       ierr = DMCompositeGather_DM(dm,next,gvec,imode,vec);CHKERRQ(ierr);
63847c6ae99SBarry Smith     } else {
63947c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
64047c6ae99SBarry Smith     }
64147c6ae99SBarry Smith     next = next->next;
64247c6ae99SBarry Smith   }
64347c6ae99SBarry Smith   va_end(Argp);
64447c6ae99SBarry Smith   PetscFunctionReturn(0);
64547c6ae99SBarry Smith }
64647c6ae99SBarry Smith 
64747c6ae99SBarry Smith #undef __FUNCT__
64847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddArray"
64947c6ae99SBarry Smith /*@C
65047c6ae99SBarry Smith     DMCompositeAddArray - adds an "redundant" array to a DMComposite. The array values will
65147c6ae99SBarry Smith        be stored in part of the array on process orank.
65247c6ae99SBarry Smith 
65347c6ae99SBarry Smith     Collective on DMComposite
65447c6ae99SBarry Smith 
65547c6ae99SBarry Smith     Input Parameter:
65647c6ae99SBarry Smith +    dm - the packer object
65747c6ae99SBarry Smith .    orank - the process on which the array entries officially live, this number must be
65847c6ae99SBarry Smith              the same on all processes.
65947c6ae99SBarry Smith -    n - the length of the array
66047c6ae99SBarry Smith 
66147c6ae99SBarry Smith     Level: advanced
66247c6ae99SBarry Smith 
6630c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
6646eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
66547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith @*/
66847c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeAddArray(DM dm,PetscMPIInt orank,PetscInt n)
66947c6ae99SBarry Smith {
67047c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
67147c6ae99SBarry Smith   PetscErrorCode         ierr;
67247c6ae99SBarry Smith   PetscMPIInt            rank;
67347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
67447c6ae99SBarry Smith 
67547c6ae99SBarry Smith   PetscFunctionBegin;
67647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
67747c6ae99SBarry Smith   next = com->next;
67847c6ae99SBarry Smith   if (com->setup) {
67947c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add an array once you have used the DMComposite");
68047c6ae99SBarry Smith   }
68147c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
68247c6ae99SBarry Smith   {
68347c6ae99SBarry Smith     PetscMPIInt orankmax;
68447c6ae99SBarry Smith     ierr = MPI_Allreduce(&orank,&orankmax,1,MPI_INT,MPI_MAX,((PetscObject)dm)->comm);CHKERRQ(ierr);
68547c6ae99SBarry 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);
68647c6ae99SBarry Smith   }
68747c6ae99SBarry Smith #endif
68847c6ae99SBarry Smith 
68947c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
69047c6ae99SBarry Smith   /* create new link */
69147c6ae99SBarry Smith   ierr                = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
69247c6ae99SBarry Smith   mine->n             = n;
69347c6ae99SBarry Smith   mine->rank          = orank;
69447c6ae99SBarry Smith   mine->dm            = PETSC_NULL;
69547c6ae99SBarry Smith   mine->type          = DMCOMPOSITE_ARRAY;
69647c6ae99SBarry Smith   mine->next          = PETSC_NULL;
69747c6ae99SBarry Smith   if (rank == mine->rank) {com->n += n;com->nmine++;}
69847c6ae99SBarry Smith 
69947c6ae99SBarry Smith   /* add to end of list */
70047c6ae99SBarry Smith   if (!next) {
70147c6ae99SBarry Smith     com->next = mine;
70247c6ae99SBarry Smith   } else {
70347c6ae99SBarry Smith     while (next->next) next = next->next;
70447c6ae99SBarry Smith     next->next = mine;
70547c6ae99SBarry Smith   }
70647c6ae99SBarry Smith   com->nredundant++;
70747c6ae99SBarry Smith   PetscFunctionReturn(0);
70847c6ae99SBarry Smith }
70947c6ae99SBarry Smith 
71047c6ae99SBarry Smith #undef __FUNCT__
71147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddDM"
71247c6ae99SBarry Smith /*@C
713aa219208SBarry Smith     DMCompositeAddDM - adds a DM  vector to a DMComposite
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith     Collective on DMComposite
71647c6ae99SBarry Smith 
71747c6ae99SBarry Smith     Input Parameter:
71847c6ae99SBarry Smith +    dm - the packer object
71947c6ae99SBarry Smith -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
72047c6ae99SBarry Smith 
72147c6ae99SBarry Smith     Level: advanced
72247c6ae99SBarry Smith 
7230c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
7246eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
72547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
72647c6ae99SBarry Smith 
72747c6ae99SBarry Smith @*/
72847c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeAddDM(DM dmc,DM dm)
72947c6ae99SBarry Smith {
73047c6ae99SBarry Smith   PetscErrorCode         ierr;
73147c6ae99SBarry Smith   PetscInt               n;
73247c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
73347c6ae99SBarry Smith   Vec                    global;
73447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmc->data;
73547c6ae99SBarry Smith 
73647c6ae99SBarry Smith   PetscFunctionBegin;
73747c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
73847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
73947c6ae99SBarry Smith   next = com->next;
740aa219208SBarry Smith   if (com->setup) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
74147c6ae99SBarry Smith 
74247c6ae99SBarry Smith   /* create new link */
74347c6ae99SBarry Smith   ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
74447c6ae99SBarry Smith   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
74547c6ae99SBarry Smith   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
74647c6ae99SBarry Smith   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
74747c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
74847c6ae99SBarry Smith   mine->n      = n;
74947c6ae99SBarry Smith   mine->dm     = dm;
75047c6ae99SBarry Smith   mine->type   = DMCOMPOSITE_DM;
75147c6ae99SBarry Smith   mine->next   = PETSC_NULL;
75247c6ae99SBarry Smith   com->n       += n;
75347c6ae99SBarry Smith 
75447c6ae99SBarry Smith   /* add to end of list */
75547c6ae99SBarry Smith   if (!next) {
75647c6ae99SBarry Smith     com->next = mine;
75747c6ae99SBarry Smith   } else {
75847c6ae99SBarry Smith     while (next->next) next = next->next;
75947c6ae99SBarry Smith     next->next = mine;
76047c6ae99SBarry Smith   }
76147c6ae99SBarry Smith   com->nDM++;
76247c6ae99SBarry Smith   com->nmine++;
76347c6ae99SBarry Smith   PetscFunctionReturn(0);
76447c6ae99SBarry Smith }
76547c6ae99SBarry Smith 
76647c6ae99SBarry Smith extern PetscErrorCode PETSCDM_DLLEXPORT VecView_MPI(Vec,PetscViewer);
76747c6ae99SBarry Smith EXTERN_C_BEGIN
76847c6ae99SBarry Smith #undef __FUNCT__
76947c6ae99SBarry Smith #define __FUNCT__ "VecView_DMComposite"
77047c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT VecView_DMComposite(Vec gvec,PetscViewer viewer)
77147c6ae99SBarry Smith {
77247c6ae99SBarry Smith   DM                     dm;
77347c6ae99SBarry Smith   PetscErrorCode         ierr;
77447c6ae99SBarry Smith   struct DMCompositeLink *next;
77547c6ae99SBarry Smith   PetscBool              isdraw;
776cef07954SSatish Balay   DM_Composite           *com;
77747c6ae99SBarry Smith 
77847c6ae99SBarry Smith   PetscFunctionBegin;
77947c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)gvec,"DMComposite",(PetscObject*)&dm);CHKERRQ(ierr);
78047c6ae99SBarry Smith   if (!dm) SETERRQ(((PetscObject)gvec)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
78147c6ae99SBarry Smith   com = (DM_Composite*)dm->data;
78247c6ae99SBarry Smith   next = com->next;
78347c6ae99SBarry Smith 
78447c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
78547c6ae99SBarry Smith   if (!isdraw) {
78647c6ae99SBarry Smith     /* do I really want to call this? */
78747c6ae99SBarry Smith     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
78847c6ae99SBarry Smith   } else {
78947c6ae99SBarry Smith     PetscInt cnt = 0;
79047c6ae99SBarry Smith 
79147c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
79247c6ae99SBarry Smith     while (next) {
79347c6ae99SBarry Smith       if (next->type == DMCOMPOSITE_ARRAY) {
79447c6ae99SBarry Smith 	PetscScalar *array;
79547c6ae99SBarry Smith 	ierr  = DMCompositeGetAccess_Array(dm,next,gvec,&array);CHKERRQ(ierr);
79647c6ae99SBarry Smith 
79747c6ae99SBarry Smith 	/*skip it for now */
79847c6ae99SBarry Smith       } else if (next->type == DMCOMPOSITE_DM) {
79947c6ae99SBarry Smith 	Vec      vec;
80047c6ae99SBarry Smith         PetscInt bs;
80147c6ae99SBarry Smith 
80247c6ae99SBarry Smith 	ierr = DMCompositeGetAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
80347c6ae99SBarry Smith 	ierr = VecView(vec,viewer);CHKERRQ(ierr);
80447c6ae99SBarry Smith         ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
80547c6ae99SBarry Smith 	ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
80647c6ae99SBarry Smith         ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
80747c6ae99SBarry Smith         cnt += bs;
80847c6ae99SBarry Smith       } else {
80947c6ae99SBarry Smith 	SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
81047c6ae99SBarry Smith       }
81147c6ae99SBarry Smith       next = next->next;
81247c6ae99SBarry Smith     }
81347c6ae99SBarry Smith     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
81447c6ae99SBarry Smith   }
81547c6ae99SBarry Smith   PetscFunctionReturn(0);
81647c6ae99SBarry Smith }
81747c6ae99SBarry Smith EXTERN_C_END
81847c6ae99SBarry Smith 
81947c6ae99SBarry Smith 
82047c6ae99SBarry Smith #undef __FUNCT__
8210c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Composite"
8220c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
82347c6ae99SBarry Smith {
82447c6ae99SBarry Smith   PetscErrorCode         ierr;
82547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
82647c6ae99SBarry Smith 
82747c6ae99SBarry Smith   PetscFunctionBegin;
82847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
82947c6ae99SBarry Smith   if (!com->setup) {
830d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
83147c6ae99SBarry Smith   }
83247c6ae99SBarry Smith   ierr = VecCreateMPI(((PetscObject)dm)->comm,com->n,com->N,gvec);CHKERRQ(ierr);
83347c6ae99SBarry Smith   ierr = PetscObjectCompose((PetscObject)*gvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr);
83447c6ae99SBarry Smith   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void(*)(void))VecView_DMComposite);CHKERRQ(ierr);
83547c6ae99SBarry Smith   PetscFunctionReturn(0);
83647c6ae99SBarry Smith }
83747c6ae99SBarry Smith 
83847c6ae99SBarry Smith #undef __FUNCT__
8390c010503SBarry Smith #define __FUNCT__ "DMCreateLocalVector_Composite"
8400c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateLocalVector_Composite(DM dm,Vec *lvec)
84147c6ae99SBarry Smith {
84247c6ae99SBarry Smith   PetscErrorCode         ierr;
84347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
84447c6ae99SBarry Smith 
84547c6ae99SBarry Smith   PetscFunctionBegin;
84647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
84747c6ae99SBarry Smith   if (!com->setup) {
848d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
84947c6ae99SBarry Smith   }
85047c6ae99SBarry Smith   ierr = VecCreateSeq(((PetscObject)dm)->comm,com->nghost,lvec);CHKERRQ(ierr);
85147c6ae99SBarry Smith   ierr = PetscObjectCompose((PetscObject)*lvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr);
85247c6ae99SBarry Smith   PetscFunctionReturn(0);
85347c6ae99SBarry Smith }
85447c6ae99SBarry Smith 
85547c6ae99SBarry Smith #undef __FUNCT__
8566eb61c8cSJed Brown #define __FUNCT__ "DMCompositeGetISLocalToGlobalMappings"
85747c6ae99SBarry Smith /*@C
8586eb61c8cSJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM/array in the DMComposite, maps to the composite global space
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith     Collective on DMComposite
86147c6ae99SBarry Smith 
86247c6ae99SBarry Smith     Input Parameter:
86347c6ae99SBarry Smith .    dm - the packer object
86447c6ae99SBarry Smith 
86547c6ae99SBarry Smith     Output Parameters:
8666eb61c8cSJed Brown .    ltogs - the individual mappings for each packed vector/array. Note that this includes
867aa219208SBarry Smith            all the ghost points that individual ghosted DMDA's may have. Also each process has an
8686eb61c8cSJed Brown            mapping for EACH redundant array (not just the local redundant arrays).
86947c6ae99SBarry Smith 
87047c6ae99SBarry Smith     Level: advanced
87147c6ae99SBarry Smith 
87247c6ae99SBarry Smith     Notes:
8736eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
87447c6ae99SBarry Smith 
8750c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
87647c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
87747c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
87847c6ae99SBarry Smith 
87947c6ae99SBarry Smith @*/
8806eb61c8cSJed Brown PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
88147c6ae99SBarry Smith {
88247c6ae99SBarry Smith   PetscErrorCode         ierr;
88347c6ae99SBarry Smith   PetscInt               i,*idx,n,cnt;
88447c6ae99SBarry Smith   struct DMCompositeLink *next;
88547c6ae99SBarry Smith   PetscMPIInt            rank;
88647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
88747c6ae99SBarry Smith 
88847c6ae99SBarry Smith   PetscFunctionBegin;
88947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
890acc1e270SJed Brown   ierr = PetscMalloc((com->nDM+com->nredundant)*sizeof(ISLocalToGlobalMapping),ltogs);CHKERRQ(ierr);
89147c6ae99SBarry Smith   next = com->next;
89247c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
89347c6ae99SBarry Smith 
89447c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
89547c6ae99SBarry Smith   cnt = 0;
89647c6ae99SBarry Smith   while (next) {
89747c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
89847c6ae99SBarry Smith       ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
89947c6ae99SBarry Smith       if (rank == next->rank) {
9006eb61c8cSJed Brown         for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
90147c6ae99SBarry Smith       }
90247c6ae99SBarry Smith       ierr = MPI_Bcast(idx,next->n,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
9036eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
90447c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
9056eb61c8cSJed Brown       ISLocalToGlobalMapping ltog;
9066eb61c8cSJed Brown       PetscMPIInt            size;
9076eb61c8cSJed Brown       const PetscInt         *suboff;
9086eb61c8cSJed Brown       Vec                    global;
90947c6ae99SBarry Smith 
9106eb61c8cSJed Brown       /* Get sub-DM global indices for each local dof */
9116eb61c8cSJed Brown       ierr = DMDAGetISLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr); /* This function should become generic to DM */
9126eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
91347c6ae99SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
9146eb61c8cSJed Brown       for (i=0; i<n; i++) idx[i] = i;
9156eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingApply(ltog,n,idx,idx);CHKERRQ(ierr); /* This would be nicer with an ISLocalToGlobalMappingGetIndices() */
91647c6ae99SBarry Smith 
9176eb61c8cSJed Brown       /* Get the offsets for the sub-DM global vector */
9186eb61c8cSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
9196eb61c8cSJed Brown       ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
9206eb61c8cSJed Brown       ierr = MPI_Comm_size(((PetscObject)global)->comm,&size);CHKERRQ(ierr);
9216eb61c8cSJed Brown 
9226eb61c8cSJed Brown       /* Shift the sub-DM definition of the global space to the composite global space */
9236eb61c8cSJed Brown       for (i=0; i<n; i++) {
9246eb61c8cSJed Brown         PetscInt subi = idx[i],lo = 0,hi = size,t;
9256eb61c8cSJed Brown         /* Binary search to find which rank owns subi */
9266eb61c8cSJed Brown         while (hi-lo > 1) {
9276eb61c8cSJed Brown           t = lo + (hi-lo)/2;
9286eb61c8cSJed Brown           if (suboff[t] > subi) hi = t;
9296eb61c8cSJed Brown           else                  lo = t;
9306eb61c8cSJed Brown         }
9316eb61c8cSJed Brown         idx[i] = subi - suboff[lo] + next->grstarts[lo];
9326eb61c8cSJed Brown       }
9336eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
9346eb61c8cSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
9356eb61c8cSJed Brown     } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
93647c6ae99SBarry Smith     next = next->next;
93747c6ae99SBarry Smith     cnt++;
93847c6ae99SBarry Smith   }
93947c6ae99SBarry Smith   PetscFunctionReturn(0);
94047c6ae99SBarry Smith }
94147c6ae99SBarry Smith 
94247c6ae99SBarry Smith #undef __FUNCT__
94387c85e80SJed Brown #define __FUNCT__ "DMCompositeGetLocalISs"
94487c85e80SJed Brown /*@C
94587c85e80SJed Brown    DMCompositeGetLocalISs - Gets index sets for each DM/array component of a composite local vector
94687c85e80SJed Brown 
94787c85e80SJed Brown    Not Collective
94887c85e80SJed Brown 
94987c85e80SJed Brown    Input Arguments:
95087c85e80SJed Brown . dm - composite DM
95187c85e80SJed Brown 
95287c85e80SJed Brown    Output Arguments:
95387c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
95487c85e80SJed Brown 
95587c85e80SJed Brown    Level: intermediate
95687c85e80SJed Brown 
95787c85e80SJed Brown    Notes:
95887c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
95987c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
96087c85e80SJed Brown    scatter to a composite local vector.
96187c85e80SJed Brown 
96287c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
96387c85e80SJed Brown 
96487c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
96587c85e80SJed Brown 
96687c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
96787c85e80SJed Brown 
96887c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
96987c85e80SJed Brown @*/
97087c85e80SJed Brown PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetLocalISs(DM dm,IS **is)
97187c85e80SJed Brown {
97287c85e80SJed Brown   PetscErrorCode         ierr;
97387c85e80SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
97487c85e80SJed Brown   struct DMCompositeLink *link;
97587c85e80SJed Brown   PetscInt cnt,start;
97687c85e80SJed Brown 
97787c85e80SJed Brown   PetscFunctionBegin;
97887c85e80SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
97987c85e80SJed Brown   PetscValidPointer(is,2);
98087c85e80SJed Brown   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
98187c85e80SJed Brown   for (cnt=0,start=0,link=com->next; link; start+=link->n,cnt++,link=link->next) {
98287c85e80SJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,link->n,start,1,&(*is)[cnt]);CHKERRQ(ierr);
983*520db06cSJed Brown     if (link->type == DMCOMPOSITE_DM) {
984*520db06cSJed Brown       Vec lvec;
985*520db06cSJed Brown       PetscInt bs;
986*520db06cSJed Brown       ierr = DMGetLocalVector(link->dm,&lvec);CHKERRQ(ierr);
987*520db06cSJed Brown       ierr = VecGetBlockSize(lvec,&bs);CHKERRQ(ierr);
988*520db06cSJed Brown       ierr = DMRestoreLocalVector(link->dm,&lvec);CHKERRQ(ierr);
989*520db06cSJed Brown       ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
990*520db06cSJed Brown     }
99187c85e80SJed Brown   }
99287c85e80SJed Brown   PetscFunctionReturn(0);
99387c85e80SJed Brown }
99487c85e80SJed Brown 
99587c85e80SJed Brown #undef __FUNCT__
99647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetGlobalISs"
99747c6ae99SBarry Smith /*@C
99847c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
99947c6ae99SBarry Smith 
100047c6ae99SBarry Smith     Collective on DMComposite
100147c6ae99SBarry Smith 
100247c6ae99SBarry Smith     Input Parameter:
100347c6ae99SBarry Smith .    dm - the packer object
100447c6ae99SBarry Smith 
100547c6ae99SBarry Smith     Output Parameters:
100647c6ae99SBarry Smith .    is - the array of index sets
100747c6ae99SBarry Smith 
100847c6ae99SBarry Smith     Level: advanced
100947c6ae99SBarry Smith 
101047c6ae99SBarry Smith     Notes:
101147c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
101247c6ae99SBarry Smith 
101347c6ae99SBarry Smith        The number of IS on each process will/may be different when redundant arrays are used
101447c6ae99SBarry Smith 
101547c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
101647c6ae99SBarry Smith 
10176eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
10186eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
10196eb61c8cSJed Brown        indices.
102047c6ae99SBarry Smith 
10210c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
102247c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
102347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
102447c6ae99SBarry Smith 
102547c6ae99SBarry Smith @*/
10266eb61c8cSJed Brown 
102747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetGlobalISs(DM dm,IS *is[])
102847c6ae99SBarry Smith {
102947c6ae99SBarry Smith   PetscErrorCode         ierr;
103047c6ae99SBarry Smith   PetscInt               cnt = 0,*idx,i;
103147c6ae99SBarry Smith   struct DMCompositeLink *next;
103247c6ae99SBarry Smith   PetscMPIInt            rank;
103347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
103447c6ae99SBarry Smith 
103547c6ae99SBarry Smith   PetscFunctionBegin;
103647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
103747c6ae99SBarry Smith   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
103847c6ae99SBarry Smith   next = com->next;
103947c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
104047c6ae99SBarry Smith 
104147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
104247c6ae99SBarry Smith   while (next) {
104347c6ae99SBarry Smith 
104447c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
104547c6ae99SBarry Smith 
104647c6ae99SBarry Smith       if (rank == next->rank) {
104747c6ae99SBarry Smith         ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
104847c6ae99SBarry Smith         for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
104947c6ae99SBarry Smith         ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
105047c6ae99SBarry Smith         cnt++;
105147c6ae99SBarry Smith       }
105247c6ae99SBarry Smith 
105347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
105447c6ae99SBarry Smith       ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
105547c6ae99SBarry Smith       for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
105647c6ae99SBarry Smith       ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
105747c6ae99SBarry Smith       cnt++;
105847c6ae99SBarry Smith     } else {
105947c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
106047c6ae99SBarry Smith     }
106147c6ae99SBarry Smith     next = next->next;
106247c6ae99SBarry Smith   }
106347c6ae99SBarry Smith   PetscFunctionReturn(0);
106447c6ae99SBarry Smith }
106547c6ae99SBarry Smith 
106647c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
106747c6ae99SBarry Smith #undef __FUNCT__
106847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_Array"
106947c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
107047c6ae99SBarry Smith {
107147c6ae99SBarry Smith   PetscErrorCode ierr;
107247c6ae99SBarry Smith   PetscFunctionBegin;
107347c6ae99SBarry Smith   if (array) {
107447c6ae99SBarry Smith     ierr = PetscMalloc(mine->n*sizeof(PetscScalar),array);CHKERRQ(ierr);
107547c6ae99SBarry Smith   }
107647c6ae99SBarry Smith   PetscFunctionReturn(0);
107747c6ae99SBarry Smith }
107847c6ae99SBarry Smith 
107947c6ae99SBarry Smith #undef __FUNCT__
108047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_DM"
108147c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
108247c6ae99SBarry Smith {
108347c6ae99SBarry Smith   PetscErrorCode ierr;
108447c6ae99SBarry Smith   PetscFunctionBegin;
108547c6ae99SBarry Smith   if (local) {
108647c6ae99SBarry Smith     ierr = DMGetLocalVector(mine->dm,local);CHKERRQ(ierr);
108747c6ae99SBarry Smith   }
108847c6ae99SBarry Smith   PetscFunctionReturn(0);
108947c6ae99SBarry Smith }
109047c6ae99SBarry Smith 
109147c6ae99SBarry Smith #undef __FUNCT__
109247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_Array"
109347c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
109447c6ae99SBarry Smith {
109547c6ae99SBarry Smith   PetscErrorCode ierr;
109647c6ae99SBarry Smith   PetscFunctionBegin;
109747c6ae99SBarry Smith   if (array) {
109847c6ae99SBarry Smith     ierr = PetscFree(*array);CHKERRQ(ierr);
109947c6ae99SBarry Smith   }
110047c6ae99SBarry Smith   PetscFunctionReturn(0);
110147c6ae99SBarry Smith }
110247c6ae99SBarry Smith 
110347c6ae99SBarry Smith #undef __FUNCT__
110447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_DM"
110547c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
110647c6ae99SBarry Smith {
110747c6ae99SBarry Smith   PetscErrorCode ierr;
110847c6ae99SBarry Smith   PetscFunctionBegin;
110947c6ae99SBarry Smith   if (local) {
111047c6ae99SBarry Smith     ierr = DMRestoreLocalVector(mine->dm,local);CHKERRQ(ierr);
111147c6ae99SBarry Smith   }
111247c6ae99SBarry Smith   PetscFunctionReturn(0);
111347c6ae99SBarry Smith }
111447c6ae99SBarry Smith 
111547c6ae99SBarry Smith #undef __FUNCT__
111647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors"
111747c6ae99SBarry Smith /*@C
111847c6ae99SBarry Smith     DMCompositeGetLocalVectors - Gets local vectors and arrays for each part of a DMComposite.
111947c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
112047c6ae99SBarry Smith 
112147c6ae99SBarry Smith     Not Collective
112247c6ae99SBarry Smith 
112347c6ae99SBarry Smith     Input Parameter:
112447c6ae99SBarry Smith .    dm - the packer object
112547c6ae99SBarry Smith 
112647c6ae99SBarry Smith     Output Parameter:
112747c6ae99SBarry Smith .    ... - the individual sequential objects (arrays or vectors)
112847c6ae99SBarry Smith 
112947c6ae99SBarry Smith     Level: advanced
113047c6ae99SBarry Smith 
11310c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
11326eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
113347c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
113447c6ae99SBarry Smith 
113547c6ae99SBarry Smith @*/
113647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetLocalVectors(DM dm,...)
113747c6ae99SBarry Smith {
113847c6ae99SBarry Smith   va_list                Argp;
113947c6ae99SBarry Smith   PetscErrorCode         ierr;
114047c6ae99SBarry Smith   struct DMCompositeLink *next;
114147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
114247c6ae99SBarry Smith 
114347c6ae99SBarry Smith   PetscFunctionBegin;
114447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
114547c6ae99SBarry Smith   next = com->next;
114647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
114747c6ae99SBarry Smith   va_start(Argp,dm);
114847c6ae99SBarry Smith   while (next) {
114947c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
115047c6ae99SBarry Smith       PetscScalar **array;
115147c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
115247c6ae99SBarry Smith       ierr = DMCompositeGetLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
115347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
115447c6ae99SBarry Smith       Vec *vec;
115547c6ae99SBarry Smith       vec = va_arg(Argp, Vec*);
115647c6ae99SBarry Smith       ierr = DMCompositeGetLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
115747c6ae99SBarry Smith     } else {
115847c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
115947c6ae99SBarry Smith     }
116047c6ae99SBarry Smith     next = next->next;
116147c6ae99SBarry Smith   }
116247c6ae99SBarry Smith   va_end(Argp);
116347c6ae99SBarry Smith   PetscFunctionReturn(0);
116447c6ae99SBarry Smith }
116547c6ae99SBarry Smith 
116647c6ae99SBarry Smith #undef __FUNCT__
116747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors"
116847c6ae99SBarry Smith /*@C
116947c6ae99SBarry Smith     DMCompositeRestoreLocalVectors - Restores local vectors and arrays for each part of a DMComposite.
117047c6ae99SBarry Smith 
117147c6ae99SBarry Smith     Not Collective
117247c6ae99SBarry Smith 
117347c6ae99SBarry Smith     Input Parameter:
117447c6ae99SBarry Smith .    dm - the packer object
117547c6ae99SBarry Smith 
117647c6ae99SBarry Smith     Output Parameter:
117747c6ae99SBarry Smith .    ... - the individual sequential objects (arrays or vectors)
117847c6ae99SBarry Smith 
117947c6ae99SBarry Smith     Level: advanced
118047c6ae99SBarry Smith 
11810c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
11826eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
118347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
118447c6ae99SBarry Smith 
118547c6ae99SBarry Smith @*/
118647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreLocalVectors(DM dm,...)
118747c6ae99SBarry Smith {
118847c6ae99SBarry Smith   va_list                Argp;
118947c6ae99SBarry Smith   PetscErrorCode         ierr;
119047c6ae99SBarry Smith   struct DMCompositeLink *next;
119147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
119247c6ae99SBarry Smith 
119347c6ae99SBarry Smith   PetscFunctionBegin;
119447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
119547c6ae99SBarry Smith   next = com->next;
119647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
119747c6ae99SBarry Smith   va_start(Argp,dm);
119847c6ae99SBarry Smith   while (next) {
119947c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
120047c6ae99SBarry Smith       PetscScalar **array;
120147c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
120247c6ae99SBarry Smith       ierr = DMCompositeRestoreLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
120347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
120447c6ae99SBarry Smith       Vec *vec;
120547c6ae99SBarry Smith       vec = va_arg(Argp, Vec*);
120647c6ae99SBarry Smith       ierr = DMCompositeRestoreLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
120747c6ae99SBarry Smith     } else {
120847c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
120947c6ae99SBarry Smith     }
121047c6ae99SBarry Smith     next = next->next;
121147c6ae99SBarry Smith   }
121247c6ae99SBarry Smith   va_end(Argp);
121347c6ae99SBarry Smith   PetscFunctionReturn(0);
121447c6ae99SBarry Smith }
121547c6ae99SBarry Smith 
121647c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
121747c6ae99SBarry Smith #undef __FUNCT__
121847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_Array"
121947c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_Array(DM dm,struct DMCompositeLink *mine,PetscInt *n)
122047c6ae99SBarry Smith {
122147c6ae99SBarry Smith   PetscFunctionBegin;
122247c6ae99SBarry Smith   if (n) *n = mine->n;
122347c6ae99SBarry Smith   PetscFunctionReturn(0);
122447c6ae99SBarry Smith }
122547c6ae99SBarry Smith 
122647c6ae99SBarry Smith #undef __FUNCT__
122747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_DM"
122847c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_DM(DM dmi,struct DMCompositeLink *mine,DM *dm)
122947c6ae99SBarry Smith {
123047c6ae99SBarry Smith   PetscFunctionBegin;
123147c6ae99SBarry Smith   if (dm) *dm = mine->dm;
123247c6ae99SBarry Smith   PetscFunctionReturn(0);
123347c6ae99SBarry Smith }
123447c6ae99SBarry Smith 
123547c6ae99SBarry Smith #undef __FUNCT__
123647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries"
123747c6ae99SBarry Smith /*@C
1238aa219208SBarry Smith     DMCompositeGetEntries - Gets the DM, redundant size, etc for each entry in a DMComposite.
123947c6ae99SBarry Smith 
124047c6ae99SBarry Smith     Not Collective
124147c6ae99SBarry Smith 
124247c6ae99SBarry Smith     Input Parameter:
124347c6ae99SBarry Smith .    dm - the packer object
124447c6ae99SBarry Smith 
124547c6ae99SBarry Smith     Output Parameter:
1246aa219208SBarry Smith .    ... - the individual entries, DMs or integer sizes)
124747c6ae99SBarry Smith 
124847c6ae99SBarry Smith     Level: advanced
124947c6ae99SBarry Smith 
12500c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
12516eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
125247c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
125347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
125447c6ae99SBarry Smith 
125547c6ae99SBarry Smith @*/
125647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetEntries(DM dm,...)
125747c6ae99SBarry Smith {
125847c6ae99SBarry Smith   va_list                Argp;
125947c6ae99SBarry Smith   PetscErrorCode         ierr;
126047c6ae99SBarry Smith   struct DMCompositeLink *next;
126147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
126247c6ae99SBarry Smith 
126347c6ae99SBarry Smith   PetscFunctionBegin;
126447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
126547c6ae99SBarry Smith   next = com->next;
126647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
126747c6ae99SBarry Smith   va_start(Argp,dm);
126847c6ae99SBarry Smith   while (next) {
126947c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
127047c6ae99SBarry Smith       PetscInt *n;
127147c6ae99SBarry Smith       n = va_arg(Argp, PetscInt*);
127247c6ae99SBarry Smith       ierr = DMCompositeGetEntries_Array(dm,next,n);CHKERRQ(ierr);
127347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
127447c6ae99SBarry Smith       DM *dmn;
127547c6ae99SBarry Smith       dmn = va_arg(Argp, DM*);
127647c6ae99SBarry Smith       ierr = DMCompositeGetEntries_DM(dm,next,dmn);CHKERRQ(ierr);
127747c6ae99SBarry Smith     } else {
127847c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
127947c6ae99SBarry Smith     }
128047c6ae99SBarry Smith     next = next->next;
128147c6ae99SBarry Smith   }
128247c6ae99SBarry Smith   va_end(Argp);
128347c6ae99SBarry Smith   PetscFunctionReturn(0);
128447c6ae99SBarry Smith }
128547c6ae99SBarry Smith 
128647c6ae99SBarry Smith #undef __FUNCT__
12870c010503SBarry Smith #define __FUNCT__ "DMRefine_Composite"
12880c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
128947c6ae99SBarry Smith {
129047c6ae99SBarry Smith   PetscErrorCode         ierr;
129147c6ae99SBarry Smith   struct DMCompositeLink *next;
129247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmi->data;
129347c6ae99SBarry Smith   DM                     dm;
129447c6ae99SBarry Smith 
129547c6ae99SBarry Smith   PetscFunctionBegin;
129647c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
129747c6ae99SBarry Smith   next = com->next;
129847c6ae99SBarry Smith   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
129947c6ae99SBarry Smith 
130047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
130147c6ae99SBarry Smith   while (next) {
130247c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
130347c6ae99SBarry Smith       ierr = DMCompositeAddArray(*fine,next->rank,next->n);CHKERRQ(ierr);
130447c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
130547c6ae99SBarry Smith       ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
130647c6ae99SBarry Smith       ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
130747c6ae99SBarry Smith       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
130847c6ae99SBarry Smith     } else {
130947c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
131047c6ae99SBarry Smith     }
131147c6ae99SBarry Smith     next = next->next;
131247c6ae99SBarry Smith   }
131347c6ae99SBarry Smith   PetscFunctionReturn(0);
131447c6ae99SBarry Smith }
131547c6ae99SBarry Smith 
131647c6ae99SBarry Smith #include "petscmat.h"
131747c6ae99SBarry Smith 
131847c6ae99SBarry Smith struct MatPackLink {
131947c6ae99SBarry Smith   Mat                A;
132047c6ae99SBarry Smith   struct MatPackLink *next;
132147c6ae99SBarry Smith };
132247c6ae99SBarry Smith 
132347c6ae99SBarry Smith struct MatPack {
132447c6ae99SBarry Smith   DM                 right,left;
132547c6ae99SBarry Smith   struct MatPackLink *next;
132647c6ae99SBarry Smith };
132747c6ae99SBarry Smith 
132847c6ae99SBarry Smith #undef __FUNCT__
132947c6ae99SBarry Smith #define __FUNCT__ "MatMultBoth_Shell_Pack"
133047c6ae99SBarry Smith PetscErrorCode MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscBool  add)
133147c6ae99SBarry Smith {
133247c6ae99SBarry Smith   struct MatPack         *mpack;
133347c6ae99SBarry Smith   struct DMCompositeLink *xnext,*ynext;
133447c6ae99SBarry Smith   struct MatPackLink     *anext;
133547c6ae99SBarry Smith   PetscScalar            *xarray,*yarray;
133647c6ae99SBarry Smith   PetscErrorCode         ierr;
133747c6ae99SBarry Smith   PetscInt               i;
133847c6ae99SBarry Smith   Vec                    xglobal,yglobal;
133947c6ae99SBarry Smith   PetscMPIInt            rank;
134047c6ae99SBarry Smith   DM_Composite           *comright;
134147c6ae99SBarry Smith   DM_Composite           *comleft;
134247c6ae99SBarry Smith 
134347c6ae99SBarry Smith   PetscFunctionBegin;
134447c6ae99SBarry Smith   ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
134547c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
134647c6ae99SBarry Smith   comright = (DM_Composite*)mpack->right->data;
134747c6ae99SBarry Smith   comleft = (DM_Composite*)mpack->left->data;
134847c6ae99SBarry Smith   xnext = comright->next;
134947c6ae99SBarry Smith   ynext = comleft->next;
135047c6ae99SBarry Smith   anext = mpack->next;
135147c6ae99SBarry Smith 
135247c6ae99SBarry Smith   while (xnext) {
135347c6ae99SBarry Smith     if (xnext->type == DMCOMPOSITE_ARRAY) {
135447c6ae99SBarry Smith       if (rank == xnext->rank) {
135547c6ae99SBarry Smith         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
135647c6ae99SBarry Smith         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
135747c6ae99SBarry Smith         if (add) {
135847c6ae99SBarry Smith           for (i=0; i<xnext->n; i++) {
135947c6ae99SBarry Smith             yarray[ynext->rstart+i] += xarray[xnext->rstart+i];
136047c6ae99SBarry Smith           }
136147c6ae99SBarry Smith         } else {
136247c6ae99SBarry Smith           ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
136347c6ae99SBarry Smith         }
136447c6ae99SBarry Smith         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
136547c6ae99SBarry Smith         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
136647c6ae99SBarry Smith       }
136747c6ae99SBarry Smith     } else if (xnext->type == DMCOMPOSITE_DM) {
136847c6ae99SBarry Smith       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
136947c6ae99SBarry Smith       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
137047c6ae99SBarry Smith       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
137147c6ae99SBarry Smith       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
137247c6ae99SBarry Smith       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
137347c6ae99SBarry Smith       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
137447c6ae99SBarry Smith       if (add) {
137547c6ae99SBarry Smith         ierr  = MatMultAdd(anext->A,xglobal,yglobal,yglobal);CHKERRQ(ierr);
137647c6ae99SBarry Smith       } else {
137747c6ae99SBarry Smith         ierr  = MatMult(anext->A,xglobal,yglobal);CHKERRQ(ierr);
137847c6ae99SBarry Smith       }
137947c6ae99SBarry Smith       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
138047c6ae99SBarry Smith       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
138147c6ae99SBarry Smith       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
138247c6ae99SBarry Smith       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
138347c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
138447c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
138547c6ae99SBarry Smith       anext = anext->next;
138647c6ae99SBarry Smith     } else {
138747c6ae99SBarry Smith       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
138847c6ae99SBarry Smith     }
138947c6ae99SBarry Smith     xnext = xnext->next;
139047c6ae99SBarry Smith     ynext = ynext->next;
139147c6ae99SBarry Smith   }
139247c6ae99SBarry Smith   PetscFunctionReturn(0);
139347c6ae99SBarry Smith }
139447c6ae99SBarry Smith 
139547c6ae99SBarry Smith #undef __FUNCT__
139647c6ae99SBarry Smith #define __FUNCT__ "MatMultAdd_Shell_Pack"
139747c6ae99SBarry Smith PetscErrorCode MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z)
139847c6ae99SBarry Smith {
139947c6ae99SBarry Smith   PetscErrorCode ierr;
140047c6ae99SBarry Smith   PetscFunctionBegin;
140147c6ae99SBarry Smith   if (z != y) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Handles y == z only");
140247c6ae99SBarry Smith   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
140347c6ae99SBarry Smith   PetscFunctionReturn(0);
140447c6ae99SBarry Smith }
140547c6ae99SBarry Smith 
140647c6ae99SBarry Smith #undef __FUNCT__
140747c6ae99SBarry Smith #define __FUNCT__ "MatMult_Shell_Pack"
140847c6ae99SBarry Smith PetscErrorCode MatMult_Shell_Pack(Mat A,Vec x,Vec y)
140947c6ae99SBarry Smith {
141047c6ae99SBarry Smith   PetscErrorCode ierr;
141147c6ae99SBarry Smith   PetscFunctionBegin;
141247c6ae99SBarry Smith   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
141347c6ae99SBarry Smith   PetscFunctionReturn(0);
141447c6ae99SBarry Smith }
141547c6ae99SBarry Smith 
141647c6ae99SBarry Smith #undef __FUNCT__
141747c6ae99SBarry Smith #define __FUNCT__ "MatMultTranspose_Shell_Pack"
141847c6ae99SBarry Smith PetscErrorCode MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y)
141947c6ae99SBarry Smith {
142047c6ae99SBarry Smith   struct MatPack         *mpack;
142147c6ae99SBarry Smith   struct DMCompositeLink *xnext,*ynext;
142247c6ae99SBarry Smith   struct MatPackLink     *anext;
142347c6ae99SBarry Smith   PetscScalar            *xarray,*yarray;
142447c6ae99SBarry Smith   PetscErrorCode         ierr;
142547c6ae99SBarry Smith   Vec                    xglobal,yglobal;
142647c6ae99SBarry Smith   PetscMPIInt            rank;
142747c6ae99SBarry Smith   DM_Composite           *comright;
142847c6ae99SBarry Smith   DM_Composite           *comleft;
142947c6ae99SBarry Smith 
143047c6ae99SBarry Smith   PetscFunctionBegin;
143147c6ae99SBarry Smith   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
143247c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
143347c6ae99SBarry Smith   comright = (DM_Composite*)mpack->right->data;
143447c6ae99SBarry Smith   comleft = (DM_Composite*)mpack->left->data;
143547c6ae99SBarry Smith   ynext = comright->next;
143647c6ae99SBarry Smith   xnext = comleft->next;
143747c6ae99SBarry Smith   anext = mpack->next;
143847c6ae99SBarry Smith 
143947c6ae99SBarry Smith   while (xnext) {
144047c6ae99SBarry Smith     if (xnext->type == DMCOMPOSITE_ARRAY) {
144147c6ae99SBarry Smith       if (rank == ynext->rank) {
144247c6ae99SBarry Smith         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
144347c6ae99SBarry Smith         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
144447c6ae99SBarry Smith         ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
144547c6ae99SBarry Smith         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
144647c6ae99SBarry Smith         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
144747c6ae99SBarry Smith       }
144847c6ae99SBarry Smith     } else if (xnext->type == DMCOMPOSITE_DM) {
144947c6ae99SBarry Smith       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
145047c6ae99SBarry Smith       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
145147c6ae99SBarry Smith       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
145247c6ae99SBarry Smith       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
145347c6ae99SBarry Smith       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
145447c6ae99SBarry Smith       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
145547c6ae99SBarry Smith       ierr  = MatMultTranspose(anext->A,xglobal,yglobal);CHKERRQ(ierr);
145647c6ae99SBarry Smith       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
145747c6ae99SBarry Smith       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
145847c6ae99SBarry Smith       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
145947c6ae99SBarry Smith       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
146047c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
146147c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
146247c6ae99SBarry Smith       anext = anext->next;
146347c6ae99SBarry Smith     } else {
146447c6ae99SBarry Smith       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
146547c6ae99SBarry Smith     }
146647c6ae99SBarry Smith     xnext = xnext->next;
146747c6ae99SBarry Smith     ynext = ynext->next;
146847c6ae99SBarry Smith   }
146947c6ae99SBarry Smith   PetscFunctionReturn(0);
147047c6ae99SBarry Smith }
147147c6ae99SBarry Smith 
147247c6ae99SBarry Smith #undef __FUNCT__
147347c6ae99SBarry Smith #define __FUNCT__ "MatDestroy_Shell_Pack"
147447c6ae99SBarry Smith PetscErrorCode MatDestroy_Shell_Pack(Mat A)
147547c6ae99SBarry Smith {
147647c6ae99SBarry Smith   struct MatPack     *mpack;
147747c6ae99SBarry Smith   struct MatPackLink *anext,*oldanext;
147847c6ae99SBarry Smith   PetscErrorCode     ierr;
147947c6ae99SBarry Smith 
148047c6ae99SBarry Smith   PetscFunctionBegin;
148147c6ae99SBarry Smith   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
148247c6ae99SBarry Smith   anext = mpack->next;
148347c6ae99SBarry Smith 
148447c6ae99SBarry Smith   while (anext) {
148547c6ae99SBarry Smith     ierr     = MatDestroy(anext->A);CHKERRQ(ierr);
148647c6ae99SBarry Smith     oldanext = anext;
148747c6ae99SBarry Smith     anext    = anext->next;
148847c6ae99SBarry Smith     ierr     = PetscFree(oldanext);CHKERRQ(ierr);
148947c6ae99SBarry Smith   }
149047c6ae99SBarry Smith   ierr = PetscFree(mpack);CHKERRQ(ierr);
149147c6ae99SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
149247c6ae99SBarry Smith   PetscFunctionReturn(0);
149347c6ae99SBarry Smith }
149447c6ae99SBarry Smith 
149547c6ae99SBarry Smith #undef __FUNCT__
14960c010503SBarry Smith #define __FUNCT__ "DMGetInterpolation_Composite"
14970c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
149847c6ae99SBarry Smith {
149947c6ae99SBarry Smith   PetscErrorCode         ierr;
150047c6ae99SBarry Smith   PetscInt               m,n,M,N;
150147c6ae99SBarry Smith   struct DMCompositeLink *nextc;
150247c6ae99SBarry Smith   struct DMCompositeLink *nextf;
150347c6ae99SBarry Smith   struct MatPackLink     *nextmat,*pnextmat = 0;
150447c6ae99SBarry Smith   struct MatPack         *mpack;
150547c6ae99SBarry Smith   Vec                    gcoarse,gfine;
150647c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
150747c6ae99SBarry Smith   DM_Composite           *comfine = (DM_Composite*)fine->data;
150847c6ae99SBarry Smith 
150947c6ae99SBarry Smith   PetscFunctionBegin;
151047c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
151147c6ae99SBarry Smith   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
151247c6ae99SBarry Smith   nextc = comcoarse->next;
151347c6ae99SBarry Smith   nextf = comfine->next;
151447c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
15150c010503SBarry Smith   ierr = DMCreateGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
15160c010503SBarry Smith   ierr = DMCreateGlobalVector(fine,&gfine);CHKERRQ(ierr);
151747c6ae99SBarry Smith   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
151847c6ae99SBarry Smith   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
151947c6ae99SBarry Smith   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
152047c6ae99SBarry Smith   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
152147c6ae99SBarry Smith   ierr = VecDestroy(gcoarse);CHKERRQ(ierr);
152247c6ae99SBarry Smith   ierr = VecDestroy(gfine);CHKERRQ(ierr);
152347c6ae99SBarry Smith 
152447c6ae99SBarry Smith   ierr         = PetscNew(struct MatPack,&mpack);CHKERRQ(ierr);
152547c6ae99SBarry Smith   mpack->right = coarse;
152647c6ae99SBarry Smith   mpack->left  = fine;
152747c6ae99SBarry Smith   ierr  = MatCreate(((PetscObject)fine)->comm,A);CHKERRQ(ierr);
152847c6ae99SBarry Smith   ierr  = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
152947c6ae99SBarry Smith   ierr  = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
153047c6ae99SBarry Smith   ierr  = MatShellSetContext(*A,mpack);CHKERRQ(ierr);
153147c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);CHKERRQ(ierr);
153247c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);CHKERRQ(ierr);
153347c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);CHKERRQ(ierr);
153447c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);CHKERRQ(ierr);
153547c6ae99SBarry Smith 
153647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
153747c6ae99SBarry Smith   while (nextc) {
153847c6ae99SBarry Smith     if (nextc->type != nextf->type) SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Two DM have different layout");
153947c6ae99SBarry Smith 
154047c6ae99SBarry Smith     if (nextc->type == DMCOMPOSITE_ARRAY) {
154147c6ae99SBarry Smith       ;
154247c6ae99SBarry Smith     } else if (nextc->type == DMCOMPOSITE_DM) {
154347c6ae99SBarry Smith       ierr          = PetscNew(struct MatPackLink,&nextmat);CHKERRQ(ierr);
154447c6ae99SBarry Smith       nextmat->next = 0;
154547c6ae99SBarry Smith       if (pnextmat) {
154647c6ae99SBarry Smith         pnextmat->next = nextmat;
154747c6ae99SBarry Smith         pnextmat       = nextmat;
154847c6ae99SBarry Smith       } else {
154947c6ae99SBarry Smith         pnextmat    = nextmat;
155047c6ae99SBarry Smith         mpack->next = nextmat;
155147c6ae99SBarry Smith       }
155247c6ae99SBarry Smith       ierr = DMGetInterpolation(nextc->dm,nextf->dm,&nextmat->A,PETSC_NULL);CHKERRQ(ierr);
155347c6ae99SBarry Smith     } else {
155447c6ae99SBarry Smith       SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
155547c6ae99SBarry Smith     }
155647c6ae99SBarry Smith     nextc = nextc->next;
155747c6ae99SBarry Smith     nextf = nextf->next;
155847c6ae99SBarry Smith   }
155947c6ae99SBarry Smith   PetscFunctionReturn(0);
156047c6ae99SBarry Smith }
156147c6ae99SBarry Smith 
156247c6ae99SBarry Smith #undef __FUNCT__
15630c010503SBarry Smith #define __FUNCT__ "DMGetMatrix_Composite"
15640c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix_Composite(DM dm, const MatType mtype,Mat *J)
156547c6ae99SBarry Smith {
156647c6ae99SBarry Smith   PetscErrorCode         ierr;
156747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
156847c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
156947c6ae99SBarry Smith   PetscInt               m,*dnz,*onz,i,j,mA;
157047c6ae99SBarry Smith   Mat                    Atmp;
157147c6ae99SBarry Smith   PetscMPIInt            rank;
157247c6ae99SBarry Smith   PetscScalar            zero = 0.0;
157347c6ae99SBarry Smith   PetscBool              dense = PETSC_FALSE;
157447c6ae99SBarry Smith 
157547c6ae99SBarry Smith   PetscFunctionBegin;
157647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
157747c6ae99SBarry Smith 
157847c6ae99SBarry Smith   /* use global vector to determine layout needed for matrix */
157947c6ae99SBarry Smith   m = com->n;
158047c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
158147c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr);
158247c6ae99SBarry Smith   ierr = MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
158347c6ae99SBarry Smith   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
158447c6ae99SBarry Smith 
158547c6ae99SBarry Smith   /*
158647c6ae99SBarry Smith      Extremely inefficient but will compute entire Jacobian for testing
158747c6ae99SBarry Smith   */
1588671f6225SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
158947c6ae99SBarry Smith   if (dense) {
159047c6ae99SBarry Smith     PetscInt    rstart,rend,*indices;
159147c6ae99SBarry Smith     PetscScalar *values;
159247c6ae99SBarry Smith 
159347c6ae99SBarry Smith     mA    = com->N;
159447c6ae99SBarry Smith     ierr = MatMPIAIJSetPreallocation(*J,mA,PETSC_NULL,mA-m,PETSC_NULL);CHKERRQ(ierr);
159547c6ae99SBarry Smith     ierr = MatSeqAIJSetPreallocation(*J,mA,PETSC_NULL);CHKERRQ(ierr);
159647c6ae99SBarry Smith 
159747c6ae99SBarry Smith     ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr);
159847c6ae99SBarry Smith     ierr = PetscMalloc2(mA,PetscScalar,&values,mA,PetscInt,&indices);CHKERRQ(ierr);
159947c6ae99SBarry Smith     ierr = PetscMemzero(values,mA*sizeof(PetscScalar));CHKERRQ(ierr);
160047c6ae99SBarry Smith     for (i=0; i<mA; i++) indices[i] = i;
160147c6ae99SBarry Smith     for (i=rstart; i<rend; i++) {
160247c6ae99SBarry Smith       ierr = MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES);CHKERRQ(ierr);
160347c6ae99SBarry Smith     }
160447c6ae99SBarry Smith     ierr = PetscFree2(values,indices);CHKERRQ(ierr);
160547c6ae99SBarry Smith     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160647c6ae99SBarry Smith     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160747c6ae99SBarry Smith     PetscFunctionReturn(0);
160847c6ae99SBarry Smith   }
160947c6ae99SBarry Smith 
161047c6ae99SBarry Smith   ierr = MatPreallocateInitialize(((PetscObject)dm)->comm,m,m,dnz,onz);CHKERRQ(ierr);
161147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
161247c6ae99SBarry Smith   next = com->next;
161347c6ae99SBarry Smith   while (next) {
161447c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
161547c6ae99SBarry Smith       if (rank == next->rank) {  /* zero the "little" block */
161647c6ae99SBarry Smith         for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) {
161747c6ae99SBarry Smith           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
161847c6ae99SBarry Smith             ierr = MatPreallocateSet(j,1,&i,dnz,onz);CHKERRQ(ierr);
161947c6ae99SBarry Smith           }
162047c6ae99SBarry Smith         }
162147c6ae99SBarry Smith       }
162247c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
162347c6ae99SBarry Smith       PetscInt       nc,rstart,*ccols,maxnc;
162447c6ae99SBarry Smith       const PetscInt *cols,*rstarts;
162547c6ae99SBarry Smith       PetscMPIInt    proc;
162647c6ae99SBarry Smith 
162747c6ae99SBarry Smith       ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr);
162847c6ae99SBarry Smith       ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr);
162947c6ae99SBarry Smith       ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
163047c6ae99SBarry Smith       ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr);
163147c6ae99SBarry Smith 
163247c6ae99SBarry Smith       maxnc = 0;
163347c6ae99SBarry Smith       for (i=0; i<mA; i++) {
163447c6ae99SBarry Smith         ierr  = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
163547c6ae99SBarry Smith         ierr  = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
163647c6ae99SBarry Smith         maxnc = PetscMax(nc,maxnc);
163747c6ae99SBarry Smith       }
163847c6ae99SBarry Smith       ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
163947c6ae99SBarry Smith       for (i=0; i<mA; i++) {
164047c6ae99SBarry Smith         ierr = MatGetRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr);
164147c6ae99SBarry Smith         /* remap the columns taking into how much they are shifted on each process */
164247c6ae99SBarry Smith         for (j=0; j<nc; j++) {
164347c6ae99SBarry Smith           proc = 0;
164447c6ae99SBarry Smith           while (cols[j] >= rstarts[proc+1]) proc++;
164547c6ae99SBarry Smith           ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
164647c6ae99SBarry Smith         }
164747c6ae99SBarry Smith         ierr = MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz);CHKERRQ(ierr);
164847c6ae99SBarry Smith         ierr = MatRestoreRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr);
164947c6ae99SBarry Smith       }
165047c6ae99SBarry Smith       ierr = PetscFree(ccols);CHKERRQ(ierr);
165147c6ae99SBarry Smith       ierr = MatDestroy(Atmp);CHKERRQ(ierr);
165247c6ae99SBarry Smith     } else {
165347c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
165447c6ae99SBarry Smith     }
165547c6ae99SBarry Smith     next = next->next;
165647c6ae99SBarry Smith   }
165747c6ae99SBarry Smith   if (com->FormCoupleLocations) {
165847c6ae99SBarry Smith     ierr = (*com->FormCoupleLocations)(dm,PETSC_NULL,dnz,onz,__rstart,__nrows,__start,__end);CHKERRQ(ierr);
165947c6ae99SBarry Smith   }
166047c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);CHKERRQ(ierr);
166147c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(*J,0,dnz);CHKERRQ(ierr);
166247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
166347c6ae99SBarry Smith 
166447c6ae99SBarry Smith   next = com->next;
166547c6ae99SBarry Smith   while (next) {
166647c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
166747c6ae99SBarry Smith       if (rank == next->rank) {
166847c6ae99SBarry Smith         for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) {
166947c6ae99SBarry Smith           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
167047c6ae99SBarry Smith             ierr = MatSetValues(*J,1,&j,1,&i,&zero,INSERT_VALUES);CHKERRQ(ierr);
167147c6ae99SBarry Smith           }
167247c6ae99SBarry Smith         }
167347c6ae99SBarry Smith       }
167447c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
167547c6ae99SBarry Smith       PetscInt          nc,rstart,row,maxnc,*ccols;
167647c6ae99SBarry Smith       const PetscInt    *cols,*rstarts;
167747c6ae99SBarry Smith       const PetscScalar *values;
167847c6ae99SBarry Smith       PetscMPIInt       proc;
167947c6ae99SBarry Smith 
168047c6ae99SBarry Smith       ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr);
168147c6ae99SBarry Smith       ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr);
168247c6ae99SBarry Smith       ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
168347c6ae99SBarry Smith       ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr);
168447c6ae99SBarry Smith       maxnc = 0;
168547c6ae99SBarry Smith       for (i=0; i<mA; i++) {
168647c6ae99SBarry Smith         ierr  = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
168747c6ae99SBarry Smith         ierr  = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
168847c6ae99SBarry Smith         maxnc = PetscMax(nc,maxnc);
168947c6ae99SBarry Smith       }
169047c6ae99SBarry Smith       ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
169147c6ae99SBarry Smith       for (i=0; i<mA; i++) {
169247c6ae99SBarry Smith         ierr = MatGetRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr);
169347c6ae99SBarry Smith         for (j=0; j<nc; j++) {
169447c6ae99SBarry Smith           proc = 0;
169547c6ae99SBarry Smith           while (cols[j] >= rstarts[proc+1]) proc++;
169647c6ae99SBarry Smith           ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
169747c6ae99SBarry Smith         }
169847c6ae99SBarry Smith         row  = com->rstart+next->rstart+i;
169947c6ae99SBarry Smith         ierr = MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES);CHKERRQ(ierr);
170047c6ae99SBarry Smith         ierr = MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr);
170147c6ae99SBarry Smith       }
170247c6ae99SBarry Smith       ierr = PetscFree(ccols);CHKERRQ(ierr);
170347c6ae99SBarry Smith       ierr = MatDestroy(Atmp);CHKERRQ(ierr);
170447c6ae99SBarry Smith     } else {
170547c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
170647c6ae99SBarry Smith     }
170747c6ae99SBarry Smith     next = next->next;
170847c6ae99SBarry Smith   }
170947c6ae99SBarry Smith   if (com->FormCoupleLocations) {
171047c6ae99SBarry Smith     PetscInt __rstart;
171147c6ae99SBarry Smith     ierr = MatGetOwnershipRange(*J,&__rstart,PETSC_NULL);CHKERRQ(ierr);
171247c6ae99SBarry Smith     ierr = (*com->FormCoupleLocations)(dm,*J,PETSC_NULL,PETSC_NULL,__rstart,0,0,0);CHKERRQ(ierr);
171347c6ae99SBarry Smith   }
171447c6ae99SBarry Smith   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
171547c6ae99SBarry Smith   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
171647c6ae99SBarry Smith   PetscFunctionReturn(0);
171747c6ae99SBarry Smith }
171847c6ae99SBarry Smith 
171947c6ae99SBarry Smith #undef __FUNCT__
17200c010503SBarry Smith #define __FUNCT__ "DMGetColoring_Composite"
17210c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
172247c6ae99SBarry Smith {
172347c6ae99SBarry Smith   PetscErrorCode         ierr;
172447c6ae99SBarry Smith   PetscInt               n,i,cnt;
172547c6ae99SBarry Smith   ISColoringValue        *colors;
172647c6ae99SBarry Smith   PetscBool              dense = PETSC_FALSE;
172747c6ae99SBarry Smith   ISColoringValue        maxcol = 0;
172847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
172947c6ae99SBarry Smith 
173047c6ae99SBarry Smith   PetscFunctionBegin;
173147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
173247c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED) {
173347c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" );
173447c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL) {
173547c6ae99SBarry Smith     n = com->n;
173647c6ae99SBarry Smith   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
173747c6ae99SBarry Smith   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
173847c6ae99SBarry Smith 
1739671f6225SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
174047c6ae99SBarry Smith   if (dense) {
174147c6ae99SBarry Smith     for (i=0; i<n; i++) {
174247c6ae99SBarry Smith       colors[i] = (ISColoringValue)(com->rstart + i);
174347c6ae99SBarry Smith     }
174447c6ae99SBarry Smith     maxcol = com->N;
174547c6ae99SBarry Smith   } else {
174647c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
174747c6ae99SBarry Smith     PetscMPIInt            rank;
174847c6ae99SBarry Smith 
174947c6ae99SBarry Smith     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
175047c6ae99SBarry Smith     cnt  = 0;
175147c6ae99SBarry Smith     while (next) {
175247c6ae99SBarry Smith       if (next->type == DMCOMPOSITE_ARRAY) {
175347c6ae99SBarry Smith         if (rank == next->rank) {  /* each column gets is own color */
175447c6ae99SBarry Smith           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
175547c6ae99SBarry Smith             colors[cnt++] = maxcol++;
175647c6ae99SBarry Smith           }
175747c6ae99SBarry Smith         }
175847c6ae99SBarry Smith         ierr = MPI_Bcast(&maxcol,1,MPIU_COLORING_VALUE,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
175947c6ae99SBarry Smith       } else if (next->type == DMCOMPOSITE_DM) {
176047c6ae99SBarry Smith         ISColoring     lcoloring;
176147c6ae99SBarry Smith 
176247c6ae99SBarry Smith         ierr = DMGetColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
176347c6ae99SBarry Smith         for (i=0; i<lcoloring->N; i++) {
176447c6ae99SBarry Smith           colors[cnt++] = maxcol + lcoloring->colors[i];
176547c6ae99SBarry Smith         }
176647c6ae99SBarry Smith         maxcol += lcoloring->n;
176747c6ae99SBarry Smith         ierr = ISColoringDestroy(lcoloring);CHKERRQ(ierr);
176847c6ae99SBarry Smith       } else {
176947c6ae99SBarry Smith         SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
177047c6ae99SBarry Smith       }
177147c6ae99SBarry Smith       next = next->next;
177247c6ae99SBarry Smith     }
177347c6ae99SBarry Smith   }
177447c6ae99SBarry Smith   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
177547c6ae99SBarry Smith   PetscFunctionReturn(0);
177647c6ae99SBarry Smith }
177747c6ae99SBarry Smith 
177847c6ae99SBarry Smith #undef __FUNCT__
17790c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
17800c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
178147c6ae99SBarry Smith {
178247c6ae99SBarry Smith   PetscErrorCode         ierr;
178347c6ae99SBarry Smith   struct DMCompositeLink *next;
178447c6ae99SBarry Smith   PetscInt               cnt = 3;
178547c6ae99SBarry Smith   PetscMPIInt            rank;
178647c6ae99SBarry Smith   PetscScalar            *garray,*larray;
178747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
178847c6ae99SBarry Smith 
178947c6ae99SBarry Smith   PetscFunctionBegin;
179047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
179147c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
179247c6ae99SBarry Smith   next = com->next;
179347c6ae99SBarry Smith   if (!com->setup) {
1794d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
179547c6ae99SBarry Smith   }
179647c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
179747c6ae99SBarry Smith   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
179847c6ae99SBarry Smith   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
179947c6ae99SBarry Smith 
180047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
180147c6ae99SBarry Smith   while (next) {
180247c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
180347c6ae99SBarry Smith       if (rank == next->rank) {
180447c6ae99SBarry Smith         ierr    = PetscMemcpy(larray,garray,next->n*sizeof(PetscScalar));CHKERRQ(ierr);
180547c6ae99SBarry Smith         garray += next->n;
180647c6ae99SBarry Smith       }
180747c6ae99SBarry Smith       /* does not handle ADD_VALUES */
180847c6ae99SBarry Smith       ierr = MPI_Bcast(larray,next->n,MPIU_SCALAR,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
180947c6ae99SBarry Smith       larray += next->n;
181047c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
181147c6ae99SBarry Smith       Vec      local,global;
181247c6ae99SBarry Smith       PetscInt N;
181347c6ae99SBarry Smith 
181447c6ae99SBarry Smith       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
181547c6ae99SBarry Smith       ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
181647c6ae99SBarry Smith       ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
181747c6ae99SBarry Smith       ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
181847c6ae99SBarry Smith       ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
181947c6ae99SBarry Smith       ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
182047c6ae99SBarry Smith       ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
182147c6ae99SBarry Smith       ierr = VecResetArray(global);CHKERRQ(ierr);
182247c6ae99SBarry Smith       ierr = VecResetArray(local);CHKERRQ(ierr);
182347c6ae99SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
182447c6ae99SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
182547c6ae99SBarry Smith       larray += next->n;
182647c6ae99SBarry Smith     } else {
182747c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
182847c6ae99SBarry Smith     }
182947c6ae99SBarry Smith     cnt++;
183047c6ae99SBarry Smith     next = next->next;
183147c6ae99SBarry Smith   }
183247c6ae99SBarry Smith 
183347c6ae99SBarry Smith   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
183447c6ae99SBarry Smith   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
183547c6ae99SBarry Smith   PetscFunctionReturn(0);
183647c6ae99SBarry Smith }
183747c6ae99SBarry Smith 
183847c6ae99SBarry Smith #undef __FUNCT__
18390c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
18400c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
18410c010503SBarry Smith {
18420c010503SBarry Smith   PetscFunctionBegin;
18430c010503SBarry Smith   PetscFunctionReturn(0);
18440c010503SBarry Smith }
184547c6ae99SBarry Smith 
1846a4121054SBarry Smith EXTERN_C_BEGIN
1847a4121054SBarry Smith #undef __FUNCT__
1848a4121054SBarry Smith #define __FUNCT__ "DMCreate_Composite"
1849a4121054SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreate_Composite(DM p)
1850a4121054SBarry Smith {
1851a4121054SBarry Smith   PetscErrorCode ierr;
1852a4121054SBarry Smith   DM_Composite   *com;
1853a4121054SBarry Smith 
1854a4121054SBarry Smith   PetscFunctionBegin;
1855a4121054SBarry Smith   ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1856a4121054SBarry Smith   p->data = com;
1857a4121054SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1858a4121054SBarry Smith   com->n            = 0;
1859a4121054SBarry Smith   com->next         = PETSC_NULL;
1860a4121054SBarry Smith   com->nredundant   = 0;
1861a4121054SBarry Smith   com->nDM          = 0;
1862a4121054SBarry Smith 
1863a4121054SBarry Smith   p->ops->createglobalvector = DMCreateGlobalVector_Composite;
1864a4121054SBarry Smith   p->ops->createlocalvector  = DMCreateLocalVector_Composite;
1865a4121054SBarry Smith   p->ops->refine             = DMRefine_Composite;
1866a4121054SBarry Smith   p->ops->getinterpolation   = DMGetInterpolation_Composite;
1867a4121054SBarry Smith   p->ops->getmatrix          = DMGetMatrix_Composite;
1868a4121054SBarry Smith   p->ops->getcoloring        = DMGetColoring_Composite;
1869a4121054SBarry Smith   p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Composite;
1870a4121054SBarry Smith   p->ops->globaltolocalend   = DMGlobalToLocalEnd_Composite;
1871a4121054SBarry Smith   p->ops->destroy            = DMDestroy_Composite;
1872a4121054SBarry Smith   p->ops->view               = DMView_Composite;
1873a4121054SBarry Smith   p->ops->setup              = DMSetUp_Composite;
1874a4121054SBarry Smith   PetscFunctionReturn(0);
1875a4121054SBarry Smith }
1876a4121054SBarry Smith EXTERN_C_END
1877a4121054SBarry Smith 
18780c010503SBarry Smith #undef __FUNCT__
18790c010503SBarry Smith #define __FUNCT__ "DMCompositeCreate"
18800c010503SBarry Smith /*@C
18810c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
18820c010503SBarry Smith       vectors made up of several subvectors.
18830c010503SBarry Smith 
18840c010503SBarry Smith     Collective on MPI_Comm
188547c6ae99SBarry Smith 
188647c6ae99SBarry Smith     Input Parameter:
18870c010503SBarry Smith .   comm - the processors that will share the global vector
18880c010503SBarry Smith 
18890c010503SBarry Smith     Output Parameters:
18900c010503SBarry Smith .   packer - the packer object
189147c6ae99SBarry Smith 
189247c6ae99SBarry Smith     Level: advanced
189347c6ae99SBarry Smith 
18940c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCompositeScatter(),
18956eb61c8cSJed Brown          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
189647c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
189747c6ae99SBarry Smith 
189847c6ae99SBarry Smith @*/
18990c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeCreate(MPI_Comm comm,DM *packer)
190047c6ae99SBarry Smith {
19010c010503SBarry Smith   PetscErrorCode ierr;
19020c010503SBarry Smith 
190347c6ae99SBarry Smith   PetscFunctionBegin;
19040c010503SBarry Smith   PetscValidPointer(packer,2);
1905a4121054SBarry Smith   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1906a4121054SBarry Smith   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
190747c6ae99SBarry Smith   PetscFunctionReturn(0);
190847c6ae99SBarry Smith }
1909