xref: /petsc/src/dm/impls/composite/pack.c (revision 06ebdd98f26409ef314cec3e093dd7660efa9fcc)
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;
17*06ebdd98SJed Brown   PetscInt               n;             /* number of owned */
18*06ebdd98SJed Brown   PetscInt               rstart;        /* rstart is relative to this process */
1947c6ae99SBarry Smith   PetscInt               grstart;       /* grstart is relative to all processes */
20*06ebdd98SJed Brown   PetscInt               nlocal;
2147c6ae99SBarry Smith 
2247c6ae99SBarry Smith   /* only used for DMCOMPOSITE_DM */
2347c6ae99SBarry Smith   PetscInt               *grstarts;     /* global row for first unknown of this DM on each process */
2447c6ae99SBarry Smith   DM                     dm;
2547c6ae99SBarry Smith 
2647c6ae99SBarry Smith   /* only used for DMCOMPOSITE_ARRAY */
2747c6ae99SBarry Smith   PetscMPIInt            rank;          /* process where array unknowns live */
2847c6ae99SBarry Smith };
2947c6ae99SBarry Smith 
3047c6ae99SBarry Smith typedef struct {
3147c6ae99SBarry Smith   PetscInt               n,N,rstart;           /* rstart is relative to all processors, n unknowns owned by this process, N is total unknowns */
32aa219208SBarry Smith   PetscInt               nghost;               /* number of all local entries include DMDA ghost points and any shared redundant arrays */
3347c6ae99SBarry Smith   PetscInt               nDM,nredundant,nmine; /* how many DM's and seperate redundant arrays used to build DM(nmine is ones on this process) */
3447c6ae99SBarry Smith   PetscBool              setup;                /* after this is set, cannot add new links to the DM*/
3547c6ae99SBarry Smith   struct DMCompositeLink *next;
3647c6ae99SBarry Smith 
3747c6ae99SBarry Smith   PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt);
3847c6ae99SBarry Smith } DM_Composite;
3947c6ae99SBarry Smith 
4047c6ae99SBarry Smith #undef __FUNCT__
4147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeSetCoupling"
4247c6ae99SBarry Smith /*@C
4347c6ae99SBarry Smith     DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
44aa219208SBarry Smith       seperate components (DMDA's and arrays) in a DMto build the correct matrix nonzero structure.
4547c6ae99SBarry Smith 
4647c6ae99SBarry Smith 
4747c6ae99SBarry Smith     Logically Collective on MPI_Comm
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith     Input Parameter:
5047c6ae99SBarry Smith +   dm - the composite object
5147c6ae99SBarry Smith -   formcouplelocations - routine to set the nonzero locations in the matrix
5247c6ae99SBarry Smith 
5347c6ae99SBarry Smith     Level: advanced
5447c6ae99SBarry Smith 
5547c6ae99SBarry Smith     Notes: See DMCompositeSetContext() and DMCompositeGetContext() for how to get user information into
5647c6ae99SBarry Smith         this routine
5747c6ae99SBarry Smith 
5847c6ae99SBarry Smith @*/
5947c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt))
6047c6ae99SBarry Smith {
6147c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
6247c6ae99SBarry Smith 
6347c6ae99SBarry Smith   PetscFunctionBegin;
6447c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
6547c6ae99SBarry Smith   PetscFunctionReturn(0);
6647c6ae99SBarry Smith }
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith #undef __FUNCT__
6947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeSetContext"
7047c6ae99SBarry Smith /*@
7147c6ae99SBarry Smith     DMCompositeSetContext - Allows user to stash data they may need within the form coupling routine they
7247c6ae99SBarry Smith       set with DMCompositeSetCoupling()
7347c6ae99SBarry Smith 
7447c6ae99SBarry Smith 
7547c6ae99SBarry Smith     Not Collective
7647c6ae99SBarry Smith 
7747c6ae99SBarry Smith     Input Parameter:
7847c6ae99SBarry Smith +   dm - the composite object
7947c6ae99SBarry Smith -   ctx - the user supplied context
8047c6ae99SBarry Smith 
8147c6ae99SBarry Smith     Level: advanced
8247c6ae99SBarry Smith 
8347c6ae99SBarry Smith     Notes: Use DMCompositeGetContext() to retrieve the context when needed.
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith @*/
8647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeSetContext(DM dm,void *ctx)
8747c6ae99SBarry Smith {
8847c6ae99SBarry Smith   PetscFunctionBegin;
8947c6ae99SBarry Smith   dm->ctx = ctx;
9047c6ae99SBarry Smith   PetscFunctionReturn(0);
9147c6ae99SBarry Smith }
9247c6ae99SBarry Smith 
9347c6ae99SBarry Smith #undef __FUNCT__
9447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetContext"
9547c6ae99SBarry Smith /*@
9647c6ae99SBarry Smith     DMCompositeGetContext - Access the context set with DMCompositeSetContext()
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith 
9947c6ae99SBarry Smith     Not Collective
10047c6ae99SBarry Smith 
10147c6ae99SBarry Smith     Input Parameter:
10247c6ae99SBarry Smith .   dm - the composite object
10347c6ae99SBarry Smith 
10447c6ae99SBarry Smith     Output Parameter:
10547c6ae99SBarry Smith .    ctx - the user supplied context
10647c6ae99SBarry Smith 
10747c6ae99SBarry Smith     Level: advanced
10847c6ae99SBarry Smith 
10947c6ae99SBarry Smith     Notes: Use DMCompositeGetContext() to retrieve the context when needed.
11047c6ae99SBarry Smith 
11147c6ae99SBarry Smith @*/
11247c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetContext(DM dm,void **ctx)
11347c6ae99SBarry Smith {
11447c6ae99SBarry Smith   PetscFunctionBegin;
11547c6ae99SBarry Smith   *ctx = dm->ctx;
11647c6ae99SBarry Smith   PetscFunctionReturn(0);
11747c6ae99SBarry Smith }
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith 
12047c6ae99SBarry Smith 
12147c6ae99SBarry Smith extern PetscErrorCode DMDestroy_Private(DM,PetscBool *);
12247c6ae99SBarry Smith 
12347c6ae99SBarry Smith #undef __FUNCT__
1240c010503SBarry Smith #define __FUNCT__ "DMDestroy_Composite"
1250c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDestroy_Composite(DM dm)
12647c6ae99SBarry Smith {
12747c6ae99SBarry Smith   PetscErrorCode         ierr;
12847c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
12947c6ae99SBarry Smith   PetscBool              done;
13047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   PetscFunctionBegin;
13347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
13447c6ae99SBarry Smith   ierr = DMDestroy_Private((DM)dm,&done);CHKERRQ(ierr);
13547c6ae99SBarry Smith   if (!done) PetscFunctionReturn(0);
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith   next = com->next;
13847c6ae99SBarry Smith   while (next) {
13947c6ae99SBarry Smith     prev = next;
14047c6ae99SBarry Smith     next = next->next;
14147c6ae99SBarry Smith     if (prev->type == DMCOMPOSITE_DM) {
14247c6ae99SBarry Smith       ierr = DMDestroy(prev->dm);CHKERRQ(ierr);
14347c6ae99SBarry Smith     }
14447c6ae99SBarry Smith     if (prev->grstarts) {
14547c6ae99SBarry Smith       ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
14647c6ae99SBarry Smith     }
14747c6ae99SBarry Smith     ierr = PetscFree(prev);CHKERRQ(ierr);
14847c6ae99SBarry Smith   }
14947c6ae99SBarry Smith   ierr = PetscFree(com);CHKERRQ(ierr);
15047c6ae99SBarry Smith   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
15147c6ae99SBarry Smith   PetscFunctionReturn(0);
15247c6ae99SBarry Smith }
15347c6ae99SBarry Smith 
15447c6ae99SBarry Smith #undef __FUNCT__
1550c010503SBarry Smith #define __FUNCT__ "DMView_Composite"
1560c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMView_Composite(DM dm,PetscViewer v)
15747c6ae99SBarry Smith {
15847c6ae99SBarry Smith   PetscErrorCode ierr;
15947c6ae99SBarry Smith   PetscBool      iascii;
16047c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite *)dm->data;
16147c6ae99SBarry Smith 
16247c6ae99SBarry Smith   PetscFunctionBegin;
16347c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
16447c6ae99SBarry Smith   if (iascii) {
16547c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
16647c6ae99SBarry Smith     PetscInt               i;
16747c6ae99SBarry Smith 
16847c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix?((PetscObject)dm)->prefix:"no prefix");CHKERRQ(ierr);
16947c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"  contains %d DMs and %d redundant arrays\n",com->nDM,com->nredundant);CHKERRQ(ierr);
17047c6ae99SBarry Smith     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
17147c6ae99SBarry Smith     for (i=0; lnk; lnk=lnk->next,i++) {
17247c6ae99SBarry Smith       if (lnk->dm) {
17347c6ae99SBarry Smith         ierr = PetscViewerASCIIPrintf(v,"Link %d: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
17447c6ae99SBarry Smith         ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
17547c6ae99SBarry Smith         ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
17647c6ae99SBarry Smith         ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
17747c6ae99SBarry Smith       } else {
178*06ebdd98SJed Brown         ierr = PetscViewerASCIIPrintf(v,"Link %d: Redundant array of size %d owned by rank %d\n",i,lnk->nlocal,lnk->rank);CHKERRQ(ierr);
17947c6ae99SBarry Smith       }
18047c6ae99SBarry Smith     }
18147c6ae99SBarry Smith     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
18247c6ae99SBarry Smith   }
18347c6ae99SBarry Smith   PetscFunctionReturn(0);
18447c6ae99SBarry Smith }
18547c6ae99SBarry Smith 
18647c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
18747c6ae99SBarry Smith #undef __FUNCT__
188d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp_Composite"
189d7bf68aeSBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSetUp_Composite(DM dm)
19047c6ae99SBarry Smith {
19147c6ae99SBarry Smith   PetscErrorCode         ierr;
19247c6ae99SBarry Smith   PetscInt               nprev = 0;
19347c6ae99SBarry Smith   PetscMPIInt            rank,size;
19447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
19547c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
19647c6ae99SBarry Smith   PetscLayout            map;
19747c6ae99SBarry Smith 
19847c6ae99SBarry Smith   PetscFunctionBegin;
19947c6ae99SBarry Smith   if (com->setup) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
20047c6ae99SBarry Smith   ierr = PetscLayoutCreate(((PetscObject)dm)->comm,&map);CHKERRQ(ierr);
20147c6ae99SBarry Smith   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
20247c6ae99SBarry Smith   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
20347c6ae99SBarry Smith   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
20447c6ae99SBarry Smith   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
20547c6ae99SBarry Smith   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
20647c6ae99SBarry Smith   ierr = PetscLayoutGetRange(map,&com->rstart,PETSC_NULL);CHKERRQ(ierr);
20747c6ae99SBarry Smith   ierr = PetscLayoutDestroy(map);CHKERRQ(ierr);
20847c6ae99SBarry Smith 
20947c6ae99SBarry Smith   /* now set the rstart for each linked array/vector */
21047c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
21147c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)dm)->comm,&size);CHKERRQ(ierr);
21247c6ae99SBarry Smith   while (next) {
21347c6ae99SBarry Smith     next->rstart = nprev;
214*06ebdd98SJed Brown     nprev += next->n;
21547c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
21647c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
21747c6ae99SBarry Smith       ierr = MPI_Bcast(&next->grstart,1,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
21847c6ae99SBarry Smith     } else {
21947c6ae99SBarry Smith       ierr = PetscMalloc(size*sizeof(PetscInt),&next->grstarts);CHKERRQ(ierr);
22047c6ae99SBarry Smith       ierr = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
22147c6ae99SBarry Smith     }
22247c6ae99SBarry Smith     next = next->next;
22347c6ae99SBarry Smith   }
22447c6ae99SBarry Smith   com->setup = PETSC_TRUE;
22547c6ae99SBarry Smith   PetscFunctionReturn(0);
22647c6ae99SBarry Smith }
22747c6ae99SBarry Smith 
22847c6ae99SBarry Smith 
22947c6ae99SBarry Smith #undef __FUNCT__
23047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess_Array"
23147c6ae99SBarry Smith PetscErrorCode DMCompositeGetAccess_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar **array)
23247c6ae99SBarry Smith {
23347c6ae99SBarry Smith   PetscErrorCode ierr;
23447c6ae99SBarry Smith   PetscScalar    *varray;
23547c6ae99SBarry Smith   PetscMPIInt    rank;
23647c6ae99SBarry Smith 
23747c6ae99SBarry Smith   PetscFunctionBegin;
23847c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
23947c6ae99SBarry Smith   if (array) {
24047c6ae99SBarry Smith     if (rank == mine->rank) {
24147c6ae99SBarry Smith       ierr    = VecGetArray(vec,&varray);CHKERRQ(ierr);
24247c6ae99SBarry Smith       *array  = varray + mine->rstart;
24347c6ae99SBarry Smith       ierr    = VecRestoreArray(vec,&varray);CHKERRQ(ierr);
24447c6ae99SBarry Smith     } else {
24547c6ae99SBarry Smith       *array = 0;
24647c6ae99SBarry Smith     }
24747c6ae99SBarry Smith   }
24847c6ae99SBarry Smith   PetscFunctionReturn(0);
24947c6ae99SBarry Smith }
25047c6ae99SBarry Smith 
25147c6ae99SBarry Smith #undef __FUNCT__
25247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess_DM"
25347c6ae99SBarry Smith PetscErrorCode DMCompositeGetAccess_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec *global)
25447c6ae99SBarry Smith {
25547c6ae99SBarry Smith   PetscErrorCode ierr;
25647c6ae99SBarry Smith   PetscScalar    *array;
25747c6ae99SBarry Smith 
25847c6ae99SBarry Smith   PetscFunctionBegin;
25947c6ae99SBarry Smith   if (global) {
26047c6ae99SBarry Smith     ierr    = DMGetGlobalVector(mine->dm,global);CHKERRQ(ierr);
26147c6ae99SBarry Smith     ierr    = VecGetArray(vec,&array);CHKERRQ(ierr);
26247c6ae99SBarry Smith     ierr    = VecPlaceArray(*global,array+mine->rstart);CHKERRQ(ierr);
26347c6ae99SBarry Smith     ierr    = VecRestoreArray(vec,&array);CHKERRQ(ierr);
26447c6ae99SBarry Smith   }
26547c6ae99SBarry Smith   PetscFunctionReturn(0);
26647c6ae99SBarry Smith }
26747c6ae99SBarry Smith 
26847c6ae99SBarry Smith #undef __FUNCT__
26947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess_Array"
27047c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreAccess_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar **array)
27147c6ae99SBarry Smith {
27247c6ae99SBarry Smith   PetscFunctionBegin;
27347c6ae99SBarry Smith   PetscFunctionReturn(0);
27447c6ae99SBarry Smith }
27547c6ae99SBarry Smith 
27647c6ae99SBarry Smith #undef __FUNCT__
27747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess_DM"
27847c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreAccess_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec *global)
27947c6ae99SBarry Smith {
28047c6ae99SBarry Smith   PetscErrorCode ierr;
28147c6ae99SBarry Smith 
28247c6ae99SBarry Smith   PetscFunctionBegin;
28347c6ae99SBarry Smith   if (global) {
28447c6ae99SBarry Smith     ierr = VecResetArray(*global);CHKERRQ(ierr);
28547c6ae99SBarry Smith     ierr = DMRestoreGlobalVector(mine->dm,global);CHKERRQ(ierr);
28647c6ae99SBarry Smith   }
28747c6ae99SBarry Smith   PetscFunctionReturn(0);
28847c6ae99SBarry Smith }
28947c6ae99SBarry Smith 
29047c6ae99SBarry Smith #undef __FUNCT__
29147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter_Array"
29247c6ae99SBarry Smith PetscErrorCode DMCompositeScatter_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar *array)
29347c6ae99SBarry Smith {
29447c6ae99SBarry Smith   PetscErrorCode ierr;
29547c6ae99SBarry Smith   PetscScalar    *varray;
29647c6ae99SBarry Smith   PetscMPIInt    rank;
29747c6ae99SBarry Smith 
29847c6ae99SBarry Smith   PetscFunctionBegin;
29947c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
30047c6ae99SBarry Smith   if (rank == mine->rank) {
30147c6ae99SBarry Smith     ierr    = VecGetArray(vec,&varray);CHKERRQ(ierr);
30247c6ae99SBarry Smith     ierr    = PetscMemcpy(array,varray+mine->rstart,mine->n*sizeof(PetscScalar));CHKERRQ(ierr);
30347c6ae99SBarry Smith     ierr    = VecRestoreArray(vec,&varray);CHKERRQ(ierr);
30447c6ae99SBarry Smith   }
305*06ebdd98SJed Brown   ierr    = MPI_Bcast(array,mine->nlocal,MPIU_SCALAR,mine->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
30647c6ae99SBarry Smith   PetscFunctionReturn(0);
30747c6ae99SBarry Smith }
30847c6ae99SBarry Smith 
30947c6ae99SBarry Smith #undef __FUNCT__
31047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter_DM"
31147c6ae99SBarry Smith PetscErrorCode DMCompositeScatter_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec local)
31247c6ae99SBarry Smith {
31347c6ae99SBarry Smith   PetscErrorCode ierr;
31447c6ae99SBarry Smith   PetscScalar    *array;
31547c6ae99SBarry Smith   Vec            global;
31647c6ae99SBarry Smith 
31747c6ae99SBarry Smith   PetscFunctionBegin;
31847c6ae99SBarry Smith   ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr);
31947c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
32047c6ae99SBarry Smith   ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr);
32147c6ae99SBarry Smith   ierr = DMGlobalToLocalBegin(mine->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
32247c6ae99SBarry Smith   ierr = DMGlobalToLocalEnd(mine->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
32347c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
32447c6ae99SBarry Smith   ierr = VecResetArray(global);CHKERRQ(ierr);
32547c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr);
32647c6ae99SBarry Smith   PetscFunctionReturn(0);
32747c6ae99SBarry Smith }
32847c6ae99SBarry Smith 
32947c6ae99SBarry Smith #undef __FUNCT__
33047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_Array"
331acc1e270SJed Brown PetscErrorCode DMCompositeGather_Array(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,const PetscScalar *array)
33247c6ae99SBarry Smith {
33347c6ae99SBarry Smith   PetscErrorCode ierr;
33447c6ae99SBarry Smith   PetscScalar    *varray;
33547c6ae99SBarry Smith   PetscMPIInt    rank;
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith   PetscFunctionBegin;
33847c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
33947c6ae99SBarry Smith   if (rank == mine->rank) {
34047c6ae99SBarry Smith     ierr = VecGetArray(vec,&varray);CHKERRQ(ierr);
34147c6ae99SBarry Smith     if (varray+mine->rstart == array) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"You need not DMCompositeGather() into objects obtained via DMCompositeGetAccess()");
342acc1e270SJed Brown   }
343df0c820aSJed Brown   switch (imode) {
344df0c820aSJed Brown   case INSERT_VALUES:
345acc1e270SJed Brown     if (rank == mine->rank) {
34647c6ae99SBarry Smith       ierr = PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(PetscScalar));CHKERRQ(ierr);
347acc1e270SJed Brown     }
348df0c820aSJed Brown     break;
349df0c820aSJed Brown   case ADD_VALUES: {
350df0c820aSJed Brown     PetscInt          i;
351760fd489SMatthew G Knepley     void             *source;
352acc1e270SJed Brown     PetscScalar       *buffer,*dest;
353acc1e270SJed Brown     if (rank == mine->rank) {
354acc1e270SJed Brown       dest = &varray[mine->rstart];
355acc1e270SJed Brown #if defined(PETSC_HAVE_MPI_IN_PLACE)
356acc1e270SJed Brown       buffer = dest;
357acc1e270SJed Brown       source = MPI_IN_PLACE;
358acc1e270SJed Brown #else
359*06ebdd98SJed Brown       ierr = PetscMalloc(mine->nlocal*sizeof(PetscScalar),&buffer);CHKERRQ(ierr);
360760fd489SMatthew G Knepley       source = (void *) buffer;
361acc1e270SJed Brown #endif
362*06ebdd98SJed Brown       for (i=0; i<mine->nlocal; i++) buffer[i] = varray[mine->rstart+i] + array[i];
363acc1e270SJed Brown     } else {
364760fd489SMatthew G Knepley       source = (void *) array;
365acc1e270SJed Brown       dest   = PETSC_NULL;
366acc1e270SJed Brown     }
367*06ebdd98SJed Brown     ierr = MPI_Reduce(source,dest,mine->nlocal,MPIU_SCALAR,MPI_SUM,mine->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
368acc1e270SJed Brown #if !defined(PETSC_HAVE_MPI_IN_PLACE)
369acc1e270SJed Brown     if (rank == mine->rank) {ierr = PetscFree(source);CHKERRQ(ierr);}
370acc1e270SJed Brown #endif
371df0c820aSJed Brown   } break;
372df0c820aSJed Brown   default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"imode");
373df0c820aSJed Brown   }
374acc1e270SJed Brown   if (rank == mine->rank) {ierr = VecRestoreArray(vec,&varray);CHKERRQ(ierr);}
37547c6ae99SBarry Smith   PetscFunctionReturn(0);
37647c6ae99SBarry Smith }
37747c6ae99SBarry Smith 
37847c6ae99SBarry Smith #undef __FUNCT__
37947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_DM"
380df0c820aSJed Brown PetscErrorCode DMCompositeGather_DM(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,Vec local)
38147c6ae99SBarry Smith {
38247c6ae99SBarry Smith   PetscErrorCode ierr;
38347c6ae99SBarry Smith   PetscScalar    *array;
38447c6ae99SBarry Smith   Vec            global;
38547c6ae99SBarry Smith 
38647c6ae99SBarry Smith   PetscFunctionBegin;
38747c6ae99SBarry Smith   ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr);
38847c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
38947c6ae99SBarry Smith   ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr);
390df0c820aSJed Brown   ierr = DMLocalToGlobalBegin(mine->dm,local,imode,global);CHKERRQ(ierr);
391df0c820aSJed Brown   ierr = DMLocalToGlobalEnd(mine->dm,local,imode,global);CHKERRQ(ierr);
39247c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
39347c6ae99SBarry Smith   ierr = VecResetArray(global);CHKERRQ(ierr);
39447c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr);
39547c6ae99SBarry Smith   PetscFunctionReturn(0);
39647c6ae99SBarry Smith }
39747c6ae99SBarry Smith 
39847c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
39947c6ae99SBarry Smith 
40047c6ae99SBarry Smith #include <stdarg.h>
40147c6ae99SBarry Smith 
40247c6ae99SBarry Smith #undef __FUNCT__
40347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetNumberDM"
40447c6ae99SBarry Smith /*@C
40547c6ae99SBarry Smith     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
40647c6ae99SBarry Smith        representation.
40747c6ae99SBarry Smith 
40847c6ae99SBarry Smith     Not Collective
40947c6ae99SBarry Smith 
41047c6ae99SBarry Smith     Input Parameter:
41147c6ae99SBarry Smith .    dm - the packer object
41247c6ae99SBarry Smith 
41347c6ae99SBarry Smith     Output Parameter:
41447c6ae99SBarry Smith .     nDM - the number of DMs
41547c6ae99SBarry Smith 
41647c6ae99SBarry Smith     Level: beginner
41747c6ae99SBarry Smith 
41847c6ae99SBarry Smith @*/
41947c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
42047c6ae99SBarry Smith {
42147c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
42247c6ae99SBarry Smith   PetscFunctionBegin;
42347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
42447c6ae99SBarry Smith   *nDM = com->nDM;
42547c6ae99SBarry Smith   PetscFunctionReturn(0);
42647c6ae99SBarry Smith }
42747c6ae99SBarry Smith 
42847c6ae99SBarry Smith 
42947c6ae99SBarry Smith #undef __FUNCT__
43047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess"
43147c6ae99SBarry Smith /*@C
43247c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
43347c6ae99SBarry Smith        representation.
43447c6ae99SBarry Smith 
43547c6ae99SBarry Smith     Collective on DMComposite
43647c6ae99SBarry Smith 
43747c6ae99SBarry Smith     Input Parameter:
43847c6ae99SBarry Smith +    dm - the packer object
43947c6ae99SBarry Smith .    gvec - the global vector
44047c6ae99SBarry Smith -    ... - the individual sequential or parallel objects (arrays or vectors)
44147c6ae99SBarry Smith 
44247c6ae99SBarry Smith     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
44347c6ae99SBarry Smith 
44447c6ae99SBarry Smith     Level: advanced
44547c6ae99SBarry Smith 
44647c6ae99SBarry Smith @*/
44747c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetAccess(DM dm,Vec gvec,...)
44847c6ae99SBarry Smith {
44947c6ae99SBarry Smith   va_list                Argp;
45047c6ae99SBarry Smith   PetscErrorCode         ierr;
45147c6ae99SBarry Smith   struct DMCompositeLink *next;
45247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
45347c6ae99SBarry Smith 
45447c6ae99SBarry Smith   PetscFunctionBegin;
45547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
45647c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
45747c6ae99SBarry Smith   next = com->next;
45847c6ae99SBarry Smith   if (!com->setup) {
459d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
46047c6ae99SBarry Smith   }
46147c6ae99SBarry Smith 
46247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
46347c6ae99SBarry Smith   va_start(Argp,gvec);
46447c6ae99SBarry Smith   while (next) {
46547c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
46647c6ae99SBarry Smith       PetscScalar **array;
46747c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
46847c6ae99SBarry Smith       ierr  = DMCompositeGetAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
46947c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
47047c6ae99SBarry Smith       Vec *vec;
47147c6ae99SBarry Smith       vec  = va_arg(Argp, Vec*);
47247c6ae99SBarry Smith       ierr = DMCompositeGetAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
47347c6ae99SBarry Smith     } else {
47447c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
47547c6ae99SBarry Smith     }
47647c6ae99SBarry Smith     next = next->next;
47747c6ae99SBarry Smith   }
47847c6ae99SBarry Smith   va_end(Argp);
47947c6ae99SBarry Smith   PetscFunctionReturn(0);
48047c6ae99SBarry Smith }
48147c6ae99SBarry Smith 
48247c6ae99SBarry Smith #undef __FUNCT__
48347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess"
48447c6ae99SBarry Smith /*@C
485aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
48647c6ae99SBarry Smith        representation.
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith     Collective on DMComposite
48947c6ae99SBarry Smith 
49047c6ae99SBarry Smith     Input Parameter:
49147c6ae99SBarry Smith +    dm - the packer object
49247c6ae99SBarry Smith .    gvec - the global vector
49347c6ae99SBarry Smith -    ... - the individual sequential or parallel objects (arrays or vectors)
49447c6ae99SBarry Smith 
49547c6ae99SBarry Smith     Level: advanced
49647c6ae99SBarry Smith 
4970c010503SBarry Smith .seealso  DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
4986eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
499aa219208SBarry Smith          DMCompositeRestoreAccess(), DMCompositeGetAccess()
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith @*/
50247c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreAccess(DM dm,Vec gvec,...)
50347c6ae99SBarry Smith {
50447c6ae99SBarry Smith   va_list                Argp;
50547c6ae99SBarry Smith   PetscErrorCode         ierr;
50647c6ae99SBarry Smith   struct DMCompositeLink *next;
50747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith   PetscFunctionBegin;
51047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
51147c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
51247c6ae99SBarry Smith   next = com->next;
51347c6ae99SBarry Smith   if (!com->setup) {
514d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
51547c6ae99SBarry Smith   }
51647c6ae99SBarry Smith 
51747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
51847c6ae99SBarry Smith   va_start(Argp,gvec);
51947c6ae99SBarry Smith   while (next) {
52047c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
52147c6ae99SBarry Smith       PetscScalar **array;
52247c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
52347c6ae99SBarry Smith       ierr  = DMCompositeRestoreAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
52447c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
52547c6ae99SBarry Smith       Vec *vec;
52647c6ae99SBarry Smith       vec  = va_arg(Argp, Vec*);
52747c6ae99SBarry Smith       ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
52847c6ae99SBarry Smith     } else {
52947c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
53047c6ae99SBarry Smith     }
53147c6ae99SBarry Smith     next = next->next;
53247c6ae99SBarry Smith   }
53347c6ae99SBarry Smith   va_end(Argp);
53447c6ae99SBarry Smith   PetscFunctionReturn(0);
53547c6ae99SBarry Smith }
53647c6ae99SBarry Smith 
53747c6ae99SBarry Smith #undef __FUNCT__
53847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter"
53947c6ae99SBarry Smith /*@C
54047c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
54147c6ae99SBarry Smith 
54247c6ae99SBarry Smith     Collective on DMComposite
54347c6ae99SBarry Smith 
54447c6ae99SBarry Smith     Input Parameter:
54547c6ae99SBarry Smith +    dm - the packer object
54647c6ae99SBarry Smith .    gvec - the global vector
5478fd8f222SJed Brown -    ... - the individual sequential objects (arrays or vectors), PETSC_NULL for those that are not needed
54847c6ae99SBarry Smith 
54947c6ae99SBarry Smith     Level: advanced
55047c6ae99SBarry Smith 
5510c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
5526eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
55347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
55447c6ae99SBarry Smith 
55547c6ae99SBarry Smith @*/
55647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeScatter(DM dm,Vec gvec,...)
55747c6ae99SBarry Smith {
55847c6ae99SBarry Smith   va_list                Argp;
55947c6ae99SBarry Smith   PetscErrorCode         ierr;
56047c6ae99SBarry Smith   struct DMCompositeLink *next;
5618fd8f222SJed Brown   PetscInt               cnt;
56247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
56347c6ae99SBarry Smith 
56447c6ae99SBarry Smith   PetscFunctionBegin;
56547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
56647c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
56747c6ae99SBarry Smith   if (!com->setup) {
568d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
56947c6ae99SBarry Smith   }
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
57247c6ae99SBarry Smith   va_start(Argp,gvec);
5738fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
57447c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
57547c6ae99SBarry Smith       PetscScalar *array;
57647c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar*);
5778fd8f222SJed Brown       if (!array) continue;
5788fd8f222SJed Brown       PetscValidScalarPointer(array,cnt);
57947c6ae99SBarry Smith       ierr = DMCompositeScatter_Array(dm,next,gvec,array);CHKERRQ(ierr);
58047c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
58147c6ae99SBarry Smith       Vec vec;
58247c6ae99SBarry Smith       vec = va_arg(Argp, Vec);
5838fd8f222SJed Brown       if (!vec) continue;
58447c6ae99SBarry Smith       PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt);
58547c6ae99SBarry Smith       ierr = DMCompositeScatter_DM(dm,next,gvec,vec);CHKERRQ(ierr);
58647c6ae99SBarry Smith     } else {
58747c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
58847c6ae99SBarry Smith     }
58947c6ae99SBarry Smith   }
59047c6ae99SBarry Smith   va_end(Argp);
59147c6ae99SBarry Smith   PetscFunctionReturn(0);
59247c6ae99SBarry Smith }
59347c6ae99SBarry Smith 
59447c6ae99SBarry Smith #undef __FUNCT__
59547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather"
59647c6ae99SBarry Smith /*@C
59747c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
59847c6ae99SBarry Smith 
59947c6ae99SBarry Smith     Collective on DMComposite
60047c6ae99SBarry Smith 
60147c6ae99SBarry Smith     Input Parameter:
60247c6ae99SBarry Smith +    dm - the packer object
60347c6ae99SBarry Smith .    gvec - the global vector
6048fd8f222SJed Brown -    ... - the individual sequential objects (arrays or vectors), PETSC_NULL for any that are not needed
60547c6ae99SBarry Smith 
60647c6ae99SBarry Smith     Level: advanced
60747c6ae99SBarry Smith 
6080c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
6096eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
61047c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
61147c6ae99SBarry Smith 
61247c6ae99SBarry Smith @*/
613df0c820aSJed Brown PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
61447c6ae99SBarry Smith {
61547c6ae99SBarry Smith   va_list                Argp;
61647c6ae99SBarry Smith   PetscErrorCode         ierr;
61747c6ae99SBarry Smith   struct DMCompositeLink *next;
61847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
6198fd8f222SJed Brown   PetscInt               cnt;
62047c6ae99SBarry Smith 
62147c6ae99SBarry Smith   PetscFunctionBegin;
62247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
62347c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
62447c6ae99SBarry Smith   if (!com->setup) {
625d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
62647c6ae99SBarry Smith   }
62747c6ae99SBarry Smith 
62847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
629df0c820aSJed Brown   va_start(Argp,imode);
6308fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
63147c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
63247c6ae99SBarry Smith       PetscScalar *array;
63347c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar*);
6348fd8f222SJed Brown       if (!array) continue;
6358fd8f222SJed Brown       PetscValidScalarPointer(array,cnt);
636df0c820aSJed Brown       ierr  = DMCompositeGather_Array(dm,next,gvec,imode,array);CHKERRQ(ierr);
63747c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
63847c6ae99SBarry Smith       Vec vec;
63947c6ae99SBarry Smith       vec = va_arg(Argp, Vec);
6408fd8f222SJed Brown       if (!vec) continue;
6418fd8f222SJed Brown       PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt);
642df0c820aSJed Brown       ierr = DMCompositeGather_DM(dm,next,gvec,imode,vec);CHKERRQ(ierr);
64347c6ae99SBarry Smith     } else {
64447c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
64547c6ae99SBarry Smith     }
64647c6ae99SBarry Smith   }
64747c6ae99SBarry Smith   va_end(Argp);
64847c6ae99SBarry Smith   PetscFunctionReturn(0);
64947c6ae99SBarry Smith }
65047c6ae99SBarry Smith 
65147c6ae99SBarry Smith #undef __FUNCT__
65247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddArray"
65347c6ae99SBarry Smith /*@C
65447c6ae99SBarry Smith     DMCompositeAddArray - adds an "redundant" array to a DMComposite. The array values will
65547c6ae99SBarry Smith        be stored in part of the array on process orank.
65647c6ae99SBarry Smith 
65747c6ae99SBarry Smith     Collective on DMComposite
65847c6ae99SBarry Smith 
65947c6ae99SBarry Smith     Input Parameter:
66047c6ae99SBarry Smith +    dm - the packer object
66147c6ae99SBarry Smith .    orank - the process on which the array entries officially live, this number must be
66247c6ae99SBarry Smith              the same on all processes.
66347c6ae99SBarry Smith -    n - the length of the array
66447c6ae99SBarry Smith 
66547c6ae99SBarry Smith     Level: advanced
66647c6ae99SBarry Smith 
6670c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
6686eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
66947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
67047c6ae99SBarry Smith 
67147c6ae99SBarry Smith @*/
67247c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeAddArray(DM dm,PetscMPIInt orank,PetscInt n)
67347c6ae99SBarry Smith {
67447c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
67547c6ae99SBarry Smith   PetscErrorCode         ierr;
67647c6ae99SBarry Smith   PetscMPIInt            rank;
67747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
67847c6ae99SBarry Smith 
67947c6ae99SBarry Smith   PetscFunctionBegin;
68047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
68147c6ae99SBarry Smith   next = com->next;
68247c6ae99SBarry Smith   if (com->setup) {
68347c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add an array once you have used the DMComposite");
68447c6ae99SBarry Smith   }
68547c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
68647c6ae99SBarry Smith   {
68747c6ae99SBarry Smith     PetscMPIInt orankmax;
68847c6ae99SBarry Smith     ierr = MPI_Allreduce(&orank,&orankmax,1,MPI_INT,MPI_MAX,((PetscObject)dm)->comm);CHKERRQ(ierr);
68947c6ae99SBarry 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);
69047c6ae99SBarry Smith   }
69147c6ae99SBarry Smith #endif
69247c6ae99SBarry Smith 
69347c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
69447c6ae99SBarry Smith   /* create new link */
69547c6ae99SBarry Smith   ierr                = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
696*06ebdd98SJed Brown   mine->nlocal        = n;
697*06ebdd98SJed Brown   mine->n             = (rank == orank) ? n : 0;
69847c6ae99SBarry Smith   mine->rank          = orank;
69947c6ae99SBarry Smith   mine->dm            = PETSC_NULL;
70047c6ae99SBarry Smith   mine->type          = DMCOMPOSITE_ARRAY;
70147c6ae99SBarry Smith   mine->next          = PETSC_NULL;
70247c6ae99SBarry Smith   if (rank == mine->rank) {com->n += n;com->nmine++;}
70347c6ae99SBarry Smith 
70447c6ae99SBarry Smith   /* add to end of list */
70547c6ae99SBarry Smith   if (!next) {
70647c6ae99SBarry Smith     com->next = mine;
70747c6ae99SBarry Smith   } else {
70847c6ae99SBarry Smith     while (next->next) next = next->next;
70947c6ae99SBarry Smith     next->next = mine;
71047c6ae99SBarry Smith   }
71147c6ae99SBarry Smith   com->nredundant++;
71247c6ae99SBarry Smith   PetscFunctionReturn(0);
71347c6ae99SBarry Smith }
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith #undef __FUNCT__
71647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddDM"
71747c6ae99SBarry Smith /*@C
718aa219208SBarry Smith     DMCompositeAddDM - adds a DM  vector to a DMComposite
71947c6ae99SBarry Smith 
72047c6ae99SBarry Smith     Collective on DMComposite
72147c6ae99SBarry Smith 
72247c6ae99SBarry Smith     Input Parameter:
72347c6ae99SBarry Smith +    dm - the packer object
72447c6ae99SBarry Smith -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
72547c6ae99SBarry Smith 
72647c6ae99SBarry Smith     Level: advanced
72747c6ae99SBarry Smith 
7280c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
7296eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
73047c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
73147c6ae99SBarry Smith 
73247c6ae99SBarry Smith @*/
73347c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeAddDM(DM dmc,DM dm)
73447c6ae99SBarry Smith {
73547c6ae99SBarry Smith   PetscErrorCode         ierr;
736*06ebdd98SJed Brown   PetscInt               n,nlocal;
73747c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
738*06ebdd98SJed Brown   Vec                    global,local;
73947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmc->data;
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith   PetscFunctionBegin;
74247c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
74347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
74447c6ae99SBarry Smith   next = com->next;
745aa219208SBarry Smith   if (com->setup) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
74647c6ae99SBarry Smith 
74747c6ae99SBarry Smith   /* create new link */
74847c6ae99SBarry Smith   ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
74947c6ae99SBarry Smith   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
75047c6ae99SBarry Smith   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
75147c6ae99SBarry Smith   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
75247c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
753*06ebdd98SJed Brown   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
754*06ebdd98SJed Brown   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
755*06ebdd98SJed Brown   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
75647c6ae99SBarry Smith   mine->n      = n;
757*06ebdd98SJed Brown   mine->nlocal = nlocal;
75847c6ae99SBarry Smith   mine->dm     = dm;
75947c6ae99SBarry Smith   mine->type   = DMCOMPOSITE_DM;
76047c6ae99SBarry Smith   mine->next   = PETSC_NULL;
76147c6ae99SBarry Smith   com->n       += n;
76247c6ae99SBarry Smith 
76347c6ae99SBarry Smith   /* add to end of list */
76447c6ae99SBarry Smith   if (!next) {
76547c6ae99SBarry Smith     com->next = mine;
76647c6ae99SBarry Smith   } else {
76747c6ae99SBarry Smith     while (next->next) next = next->next;
76847c6ae99SBarry Smith     next->next = mine;
76947c6ae99SBarry Smith   }
77047c6ae99SBarry Smith   com->nDM++;
77147c6ae99SBarry Smith   com->nmine++;
77247c6ae99SBarry Smith   PetscFunctionReturn(0);
77347c6ae99SBarry Smith }
77447c6ae99SBarry Smith 
77547c6ae99SBarry Smith extern PetscErrorCode PETSCDM_DLLEXPORT VecView_MPI(Vec,PetscViewer);
77647c6ae99SBarry Smith EXTERN_C_BEGIN
77747c6ae99SBarry Smith #undef __FUNCT__
77847c6ae99SBarry Smith #define __FUNCT__ "VecView_DMComposite"
77947c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT VecView_DMComposite(Vec gvec,PetscViewer viewer)
78047c6ae99SBarry Smith {
78147c6ae99SBarry Smith   DM                     dm;
78247c6ae99SBarry Smith   PetscErrorCode         ierr;
78347c6ae99SBarry Smith   struct DMCompositeLink *next;
78447c6ae99SBarry Smith   PetscBool              isdraw;
785cef07954SSatish Balay   DM_Composite           *com;
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith   PetscFunctionBegin;
78847c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)gvec,"DMComposite",(PetscObject*)&dm);CHKERRQ(ierr);
78947c6ae99SBarry Smith   if (!dm) SETERRQ(((PetscObject)gvec)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
79047c6ae99SBarry Smith   com = (DM_Composite*)dm->data;
79147c6ae99SBarry Smith   next = com->next;
79247c6ae99SBarry Smith 
79347c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
79447c6ae99SBarry Smith   if (!isdraw) {
79547c6ae99SBarry Smith     /* do I really want to call this? */
79647c6ae99SBarry Smith     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
79747c6ae99SBarry Smith   } else {
79847c6ae99SBarry Smith     PetscInt cnt = 0;
79947c6ae99SBarry Smith 
80047c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
80147c6ae99SBarry Smith     while (next) {
80247c6ae99SBarry Smith       if (next->type == DMCOMPOSITE_ARRAY) {
80347c6ae99SBarry Smith 	PetscScalar *array;
80447c6ae99SBarry Smith 	ierr  = DMCompositeGetAccess_Array(dm,next,gvec,&array);CHKERRQ(ierr);
80547c6ae99SBarry Smith 
80647c6ae99SBarry Smith 	/*skip it for now */
80747c6ae99SBarry Smith       } else if (next->type == DMCOMPOSITE_DM) {
80847c6ae99SBarry Smith 	Vec      vec;
80947c6ae99SBarry Smith         PetscInt bs;
81047c6ae99SBarry Smith 
81147c6ae99SBarry Smith 	ierr = DMCompositeGetAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
81247c6ae99SBarry Smith 	ierr = VecView(vec,viewer);CHKERRQ(ierr);
81347c6ae99SBarry Smith         ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
81447c6ae99SBarry Smith 	ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
81547c6ae99SBarry Smith         ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
81647c6ae99SBarry Smith         cnt += bs;
81747c6ae99SBarry Smith       } else {
81847c6ae99SBarry Smith 	SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
81947c6ae99SBarry Smith       }
82047c6ae99SBarry Smith       next = next->next;
82147c6ae99SBarry Smith     }
82247c6ae99SBarry Smith     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
82347c6ae99SBarry Smith   }
82447c6ae99SBarry Smith   PetscFunctionReturn(0);
82547c6ae99SBarry Smith }
82647c6ae99SBarry Smith EXTERN_C_END
82747c6ae99SBarry Smith 
82847c6ae99SBarry Smith 
82947c6ae99SBarry Smith #undef __FUNCT__
8300c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Composite"
8310c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
83247c6ae99SBarry Smith {
83347c6ae99SBarry Smith   PetscErrorCode         ierr;
83447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
83547c6ae99SBarry Smith 
83647c6ae99SBarry Smith   PetscFunctionBegin;
83747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
83847c6ae99SBarry Smith   if (!com->setup) {
839d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
84047c6ae99SBarry Smith   }
84147c6ae99SBarry Smith   ierr = VecCreateMPI(((PetscObject)dm)->comm,com->n,com->N,gvec);CHKERRQ(ierr);
84247c6ae99SBarry Smith   ierr = PetscObjectCompose((PetscObject)*gvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr);
84347c6ae99SBarry Smith   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void(*)(void))VecView_DMComposite);CHKERRQ(ierr);
84447c6ae99SBarry Smith   PetscFunctionReturn(0);
84547c6ae99SBarry Smith }
84647c6ae99SBarry Smith 
84747c6ae99SBarry Smith #undef __FUNCT__
8480c010503SBarry Smith #define __FUNCT__ "DMCreateLocalVector_Composite"
8490c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateLocalVector_Composite(DM dm,Vec *lvec)
85047c6ae99SBarry Smith {
85147c6ae99SBarry Smith   PetscErrorCode         ierr;
85247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith   PetscFunctionBegin;
85547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
85647c6ae99SBarry Smith   if (!com->setup) {
857d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
85847c6ae99SBarry Smith   }
85947c6ae99SBarry Smith   ierr = VecCreateSeq(((PetscObject)dm)->comm,com->nghost,lvec);CHKERRQ(ierr);
86047c6ae99SBarry Smith   ierr = PetscObjectCompose((PetscObject)*lvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr);
86147c6ae99SBarry Smith   PetscFunctionReturn(0);
86247c6ae99SBarry Smith }
86347c6ae99SBarry Smith 
86447c6ae99SBarry Smith #undef __FUNCT__
8656eb61c8cSJed Brown #define __FUNCT__ "DMCompositeGetISLocalToGlobalMappings"
86647c6ae99SBarry Smith /*@C
8676eb61c8cSJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM/array in the DMComposite, maps to the composite global space
86847c6ae99SBarry Smith 
869*06ebdd98SJed Brown     Collective on DM
87047c6ae99SBarry Smith 
87147c6ae99SBarry Smith     Input Parameter:
87247c6ae99SBarry Smith .    dm - the packer object
87347c6ae99SBarry Smith 
87447c6ae99SBarry Smith     Output Parameters:
8756eb61c8cSJed Brown .    ltogs - the individual mappings for each packed vector/array. Note that this includes
876aa219208SBarry Smith            all the ghost points that individual ghosted DMDA's may have. Also each process has an
8776eb61c8cSJed Brown            mapping for EACH redundant array (not just the local redundant arrays).
87847c6ae99SBarry Smith 
87947c6ae99SBarry Smith     Level: advanced
88047c6ae99SBarry Smith 
88147c6ae99SBarry Smith     Notes:
8826eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
88347c6ae99SBarry Smith 
8840c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
88547c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
88647c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
88747c6ae99SBarry Smith 
88847c6ae99SBarry Smith @*/
8896eb61c8cSJed Brown PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
89047c6ae99SBarry Smith {
89147c6ae99SBarry Smith   PetscErrorCode         ierr;
89247c6ae99SBarry Smith   PetscInt               i,*idx,n,cnt;
89347c6ae99SBarry Smith   struct DMCompositeLink *next;
89447c6ae99SBarry Smith   PetscMPIInt            rank;
89547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
89647c6ae99SBarry Smith 
89747c6ae99SBarry Smith   PetscFunctionBegin;
89847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
899acc1e270SJed Brown   ierr = PetscMalloc((com->nDM+com->nredundant)*sizeof(ISLocalToGlobalMapping),ltogs);CHKERRQ(ierr);
90047c6ae99SBarry Smith   next = com->next;
90147c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
90247c6ae99SBarry Smith 
90347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
90447c6ae99SBarry Smith   cnt = 0;
90547c6ae99SBarry Smith   while (next) {
90647c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
907*06ebdd98SJed Brown       ierr = PetscMalloc(next->nlocal*sizeof(PetscInt),&idx);CHKERRQ(ierr);
90847c6ae99SBarry Smith       if (rank == next->rank) {
909*06ebdd98SJed Brown         for (i=0; i<next->nlocal; i++) idx[i] = next->grstart + i;
91047c6ae99SBarry Smith       }
911*06ebdd98SJed Brown       ierr = MPI_Bcast(idx,next->nlocal,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
912*06ebdd98SJed Brown       ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,next->nlocal,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
91347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
9146eb61c8cSJed Brown       ISLocalToGlobalMapping ltog;
9156eb61c8cSJed Brown       PetscMPIInt            size;
91686994e45SJed Brown       const PetscInt         *suboff,*indices;
9176eb61c8cSJed Brown       Vec                    global;
91847c6ae99SBarry Smith 
9196eb61c8cSJed Brown       /* Get sub-DM global indices for each local dof */
9201411c6eeSJed Brown       ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
9216eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
92286994e45SJed Brown       ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
92347c6ae99SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
92447c6ae99SBarry Smith 
9256eb61c8cSJed Brown       /* Get the offsets for the sub-DM global vector */
9266eb61c8cSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
9276eb61c8cSJed Brown       ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
9286eb61c8cSJed Brown       ierr = MPI_Comm_size(((PetscObject)global)->comm,&size);CHKERRQ(ierr);
9296eb61c8cSJed Brown 
9306eb61c8cSJed Brown       /* Shift the sub-DM definition of the global space to the composite global space */
9316eb61c8cSJed Brown       for (i=0; i<n; i++) {
93286994e45SJed Brown         PetscInt subi = indices[i],lo = 0,hi = size,t;
9336eb61c8cSJed Brown         /* Binary search to find which rank owns subi */
9346eb61c8cSJed Brown         while (hi-lo > 1) {
9356eb61c8cSJed Brown           t = lo + (hi-lo)/2;
9366eb61c8cSJed Brown           if (suboff[t] > subi) hi = t;
9376eb61c8cSJed Brown           else                  lo = t;
9386eb61c8cSJed Brown         }
9396eb61c8cSJed Brown         idx[i] = subi - suboff[lo] + next->grstarts[lo];
9406eb61c8cSJed Brown       }
94186994e45SJed Brown       ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
9426eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
9436eb61c8cSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
9446eb61c8cSJed Brown     } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
94547c6ae99SBarry Smith     next = next->next;
94647c6ae99SBarry Smith     cnt++;
94747c6ae99SBarry Smith   }
94847c6ae99SBarry Smith   PetscFunctionReturn(0);
94947c6ae99SBarry Smith }
95047c6ae99SBarry Smith 
95147c6ae99SBarry Smith #undef __FUNCT__
95287c85e80SJed Brown #define __FUNCT__ "DMCompositeGetLocalISs"
95387c85e80SJed Brown /*@C
95487c85e80SJed Brown    DMCompositeGetLocalISs - Gets index sets for each DM/array component of a composite local vector
95587c85e80SJed Brown 
95687c85e80SJed Brown    Not Collective
95787c85e80SJed Brown 
95887c85e80SJed Brown    Input Arguments:
95987c85e80SJed Brown . dm - composite DM
96087c85e80SJed Brown 
96187c85e80SJed Brown    Output Arguments:
96287c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
96387c85e80SJed Brown 
96487c85e80SJed Brown    Level: intermediate
96587c85e80SJed Brown 
96687c85e80SJed Brown    Notes:
96787c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
96887c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
96987c85e80SJed Brown    scatter to a composite local vector.
97087c85e80SJed Brown 
97187c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
97287c85e80SJed Brown 
97387c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
97487c85e80SJed Brown 
97587c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
97687c85e80SJed Brown 
97787c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
97887c85e80SJed Brown @*/
97987c85e80SJed Brown PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetLocalISs(DM dm,IS **is)
98087c85e80SJed Brown {
98187c85e80SJed Brown   PetscErrorCode         ierr;
98287c85e80SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
98387c85e80SJed Brown   struct DMCompositeLink *link;
98487c85e80SJed Brown   PetscInt cnt,start;
98587c85e80SJed Brown 
98687c85e80SJed Brown   PetscFunctionBegin;
98787c85e80SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
98887c85e80SJed Brown   PetscValidPointer(is,2);
98987c85e80SJed Brown   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
990*06ebdd98SJed Brown   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
991*06ebdd98SJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
992520db06cSJed Brown     if (link->type == DMCOMPOSITE_DM) {
993520db06cSJed Brown       PetscInt bs;
9941411c6eeSJed Brown       ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
995520db06cSJed Brown       ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
996520db06cSJed Brown     }
99787c85e80SJed Brown   }
99887c85e80SJed Brown   PetscFunctionReturn(0);
99987c85e80SJed Brown }
100087c85e80SJed Brown 
100187c85e80SJed Brown #undef __FUNCT__
100247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetGlobalISs"
100347c6ae99SBarry Smith /*@C
100447c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
100547c6ae99SBarry Smith 
100647c6ae99SBarry Smith     Collective on DMComposite
100747c6ae99SBarry Smith 
100847c6ae99SBarry Smith     Input Parameter:
100947c6ae99SBarry Smith .    dm - the packer object
101047c6ae99SBarry Smith 
101147c6ae99SBarry Smith     Output Parameters:
101247c6ae99SBarry Smith .    is - the array of index sets
101347c6ae99SBarry Smith 
101447c6ae99SBarry Smith     Level: advanced
101547c6ae99SBarry Smith 
101647c6ae99SBarry Smith     Notes:
101747c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
101847c6ae99SBarry Smith 
101947c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
102047c6ae99SBarry Smith 
10216eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
10226eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
10236eb61c8cSJed Brown        indices.
102447c6ae99SBarry Smith 
10250c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
102647c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
102747c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
102847c6ae99SBarry Smith 
102947c6ae99SBarry Smith @*/
10306eb61c8cSJed Brown 
103147c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetGlobalISs(DM dm,IS *is[])
103247c6ae99SBarry Smith {
103347c6ae99SBarry Smith   PetscErrorCode         ierr;
103447c6ae99SBarry Smith   PetscInt               cnt = 0,*idx,i;
103547c6ae99SBarry Smith   struct DMCompositeLink *next;
103647c6ae99SBarry Smith   PetscMPIInt            rank;
103747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
103847c6ae99SBarry Smith 
103947c6ae99SBarry Smith   PetscFunctionBegin;
104047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1041*06ebdd98SJed Brown   ierr = PetscMalloc((com->nDM+com->nredundant)*sizeof(IS),is);CHKERRQ(ierr);
104247c6ae99SBarry Smith   next = com->next;
104347c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
104447c6ae99SBarry Smith 
104547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
104647c6ae99SBarry Smith   while (next) {
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     next = next->next;
105247c6ae99SBarry Smith   }
105347c6ae99SBarry Smith   PetscFunctionReturn(0);
105447c6ae99SBarry Smith }
105547c6ae99SBarry Smith 
105647c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
105747c6ae99SBarry Smith #undef __FUNCT__
105847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_Array"
105947c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
106047c6ae99SBarry Smith {
106147c6ae99SBarry Smith   PetscErrorCode ierr;
106247c6ae99SBarry Smith   PetscFunctionBegin;
106347c6ae99SBarry Smith   if (array) {
106447c6ae99SBarry Smith     ierr = PetscMalloc(mine->n*sizeof(PetscScalar),array);CHKERRQ(ierr);
106547c6ae99SBarry Smith   }
106647c6ae99SBarry Smith   PetscFunctionReturn(0);
106747c6ae99SBarry Smith }
106847c6ae99SBarry Smith 
106947c6ae99SBarry Smith #undef __FUNCT__
107047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_DM"
107147c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
107247c6ae99SBarry Smith {
107347c6ae99SBarry Smith   PetscErrorCode ierr;
107447c6ae99SBarry Smith   PetscFunctionBegin;
107547c6ae99SBarry Smith   if (local) {
107647c6ae99SBarry Smith     ierr = DMGetLocalVector(mine->dm,local);CHKERRQ(ierr);
107747c6ae99SBarry Smith   }
107847c6ae99SBarry Smith   PetscFunctionReturn(0);
107947c6ae99SBarry Smith }
108047c6ae99SBarry Smith 
108147c6ae99SBarry Smith #undef __FUNCT__
108247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_Array"
108347c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
108447c6ae99SBarry Smith {
108547c6ae99SBarry Smith   PetscErrorCode ierr;
108647c6ae99SBarry Smith   PetscFunctionBegin;
108747c6ae99SBarry Smith   if (array) {
108847c6ae99SBarry Smith     ierr = PetscFree(*array);CHKERRQ(ierr);
108947c6ae99SBarry Smith   }
109047c6ae99SBarry Smith   PetscFunctionReturn(0);
109147c6ae99SBarry Smith }
109247c6ae99SBarry Smith 
109347c6ae99SBarry Smith #undef __FUNCT__
109447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_DM"
109547c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
109647c6ae99SBarry Smith {
109747c6ae99SBarry Smith   PetscErrorCode ierr;
109847c6ae99SBarry Smith   PetscFunctionBegin;
109947c6ae99SBarry Smith   if (local) {
110047c6ae99SBarry Smith     ierr = DMRestoreLocalVector(mine->dm,local);CHKERRQ(ierr);
110147c6ae99SBarry Smith   }
110247c6ae99SBarry Smith   PetscFunctionReturn(0);
110347c6ae99SBarry Smith }
110447c6ae99SBarry Smith 
110547c6ae99SBarry Smith #undef __FUNCT__
110647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors"
110747c6ae99SBarry Smith /*@C
110847c6ae99SBarry Smith     DMCompositeGetLocalVectors - Gets local vectors and arrays for each part of a DMComposite.
110947c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
111047c6ae99SBarry Smith 
111147c6ae99SBarry Smith     Not Collective
111247c6ae99SBarry Smith 
111347c6ae99SBarry Smith     Input Parameter:
111447c6ae99SBarry Smith .    dm - the packer object
111547c6ae99SBarry Smith 
111647c6ae99SBarry Smith     Output Parameter:
111747c6ae99SBarry Smith .    ... - the individual sequential objects (arrays or vectors)
111847c6ae99SBarry Smith 
111947c6ae99SBarry Smith     Level: advanced
112047c6ae99SBarry Smith 
11210c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
11226eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
112347c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
112447c6ae99SBarry Smith 
112547c6ae99SBarry Smith @*/
112647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetLocalVectors(DM dm,...)
112747c6ae99SBarry Smith {
112847c6ae99SBarry Smith   va_list                Argp;
112947c6ae99SBarry Smith   PetscErrorCode         ierr;
113047c6ae99SBarry Smith   struct DMCompositeLink *next;
113147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
113247c6ae99SBarry Smith 
113347c6ae99SBarry Smith   PetscFunctionBegin;
113447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
113547c6ae99SBarry Smith   next = com->next;
113647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
113747c6ae99SBarry Smith   va_start(Argp,dm);
113847c6ae99SBarry Smith   while (next) {
113947c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
114047c6ae99SBarry Smith       PetscScalar **array;
114147c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
114247c6ae99SBarry Smith       ierr = DMCompositeGetLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
114347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
114447c6ae99SBarry Smith       Vec *vec;
114547c6ae99SBarry Smith       vec = va_arg(Argp, Vec*);
114647c6ae99SBarry Smith       ierr = DMCompositeGetLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
114747c6ae99SBarry Smith     } else {
114847c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
114947c6ae99SBarry Smith     }
115047c6ae99SBarry Smith     next = next->next;
115147c6ae99SBarry Smith   }
115247c6ae99SBarry Smith   va_end(Argp);
115347c6ae99SBarry Smith   PetscFunctionReturn(0);
115447c6ae99SBarry Smith }
115547c6ae99SBarry Smith 
115647c6ae99SBarry Smith #undef __FUNCT__
115747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors"
115847c6ae99SBarry Smith /*@C
115947c6ae99SBarry Smith     DMCompositeRestoreLocalVectors - Restores local vectors and arrays for each part of a DMComposite.
116047c6ae99SBarry Smith 
116147c6ae99SBarry Smith     Not Collective
116247c6ae99SBarry Smith 
116347c6ae99SBarry Smith     Input Parameter:
116447c6ae99SBarry Smith .    dm - the packer object
116547c6ae99SBarry Smith 
116647c6ae99SBarry Smith     Output Parameter:
116747c6ae99SBarry Smith .    ... - the individual sequential objects (arrays or vectors)
116847c6ae99SBarry Smith 
116947c6ae99SBarry Smith     Level: advanced
117047c6ae99SBarry Smith 
11710c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
11726eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
117347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
117447c6ae99SBarry Smith 
117547c6ae99SBarry Smith @*/
117647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeRestoreLocalVectors(DM dm,...)
117747c6ae99SBarry Smith {
117847c6ae99SBarry Smith   va_list                Argp;
117947c6ae99SBarry Smith   PetscErrorCode         ierr;
118047c6ae99SBarry Smith   struct DMCompositeLink *next;
118147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
118247c6ae99SBarry Smith 
118347c6ae99SBarry Smith   PetscFunctionBegin;
118447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
118547c6ae99SBarry Smith   next = com->next;
118647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
118747c6ae99SBarry Smith   va_start(Argp,dm);
118847c6ae99SBarry Smith   while (next) {
118947c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
119047c6ae99SBarry Smith       PetscScalar **array;
119147c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
119247c6ae99SBarry Smith       ierr = DMCompositeRestoreLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
119347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
119447c6ae99SBarry Smith       Vec *vec;
119547c6ae99SBarry Smith       vec = va_arg(Argp, Vec*);
119647c6ae99SBarry Smith       ierr = DMCompositeRestoreLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
119747c6ae99SBarry Smith     } else {
119847c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
119947c6ae99SBarry Smith     }
120047c6ae99SBarry Smith     next = next->next;
120147c6ae99SBarry Smith   }
120247c6ae99SBarry Smith   va_end(Argp);
120347c6ae99SBarry Smith   PetscFunctionReturn(0);
120447c6ae99SBarry Smith }
120547c6ae99SBarry Smith 
120647c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
120747c6ae99SBarry Smith #undef __FUNCT__
120847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_Array"
120947c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_Array(DM dm,struct DMCompositeLink *mine,PetscInt *n)
121047c6ae99SBarry Smith {
121147c6ae99SBarry Smith   PetscFunctionBegin;
1212*06ebdd98SJed Brown   if (n) *n = mine->nlocal;
121347c6ae99SBarry Smith   PetscFunctionReturn(0);
121447c6ae99SBarry Smith }
121547c6ae99SBarry Smith 
121647c6ae99SBarry Smith #undef __FUNCT__
121747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_DM"
121847c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_DM(DM dmi,struct DMCompositeLink *mine,DM *dm)
121947c6ae99SBarry Smith {
122047c6ae99SBarry Smith   PetscFunctionBegin;
122147c6ae99SBarry Smith   if (dm) *dm = mine->dm;
122247c6ae99SBarry Smith   PetscFunctionReturn(0);
122347c6ae99SBarry Smith }
122447c6ae99SBarry Smith 
122547c6ae99SBarry Smith #undef __FUNCT__
122647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries"
122747c6ae99SBarry Smith /*@C
1228aa219208SBarry Smith     DMCompositeGetEntries - Gets the DM, redundant size, etc for each entry in a DMComposite.
122947c6ae99SBarry Smith 
123047c6ae99SBarry Smith     Not Collective
123147c6ae99SBarry Smith 
123247c6ae99SBarry Smith     Input Parameter:
123347c6ae99SBarry Smith .    dm - the packer object
123447c6ae99SBarry Smith 
123547c6ae99SBarry Smith     Output Parameter:
1236aa219208SBarry Smith .    ... - the individual entries, DMs or integer sizes)
123747c6ae99SBarry Smith 
123847c6ae99SBarry Smith     Level: advanced
123947c6ae99SBarry Smith 
12400c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
12416eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
124247c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
124347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
124447c6ae99SBarry Smith 
124547c6ae99SBarry Smith @*/
124647c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeGetEntries(DM dm,...)
124747c6ae99SBarry Smith {
124847c6ae99SBarry Smith   va_list                Argp;
124947c6ae99SBarry Smith   PetscErrorCode         ierr;
125047c6ae99SBarry Smith   struct DMCompositeLink *next;
125147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
125247c6ae99SBarry Smith 
125347c6ae99SBarry Smith   PetscFunctionBegin;
125447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
125547c6ae99SBarry Smith   next = com->next;
125647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
125747c6ae99SBarry Smith   va_start(Argp,dm);
125847c6ae99SBarry Smith   while (next) {
125947c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
126047c6ae99SBarry Smith       PetscInt *n;
126147c6ae99SBarry Smith       n = va_arg(Argp, PetscInt*);
126247c6ae99SBarry Smith       ierr = DMCompositeGetEntries_Array(dm,next,n);CHKERRQ(ierr);
126347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
126447c6ae99SBarry Smith       DM *dmn;
126547c6ae99SBarry Smith       dmn = va_arg(Argp, DM*);
126647c6ae99SBarry Smith       ierr = DMCompositeGetEntries_DM(dm,next,dmn);CHKERRQ(ierr);
126747c6ae99SBarry Smith     } else {
126847c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
126947c6ae99SBarry Smith     }
127047c6ae99SBarry Smith     next = next->next;
127147c6ae99SBarry Smith   }
127247c6ae99SBarry Smith   va_end(Argp);
127347c6ae99SBarry Smith   PetscFunctionReturn(0);
127447c6ae99SBarry Smith }
127547c6ae99SBarry Smith 
127647c6ae99SBarry Smith #undef __FUNCT__
12770c010503SBarry Smith #define __FUNCT__ "DMRefine_Composite"
12780c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
127947c6ae99SBarry Smith {
128047c6ae99SBarry Smith   PetscErrorCode         ierr;
128147c6ae99SBarry Smith   struct DMCompositeLink *next;
128247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmi->data;
128347c6ae99SBarry Smith   DM                     dm;
128447c6ae99SBarry Smith 
128547c6ae99SBarry Smith   PetscFunctionBegin;
128647c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
128747c6ae99SBarry Smith   next = com->next;
128847c6ae99SBarry Smith   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
128947c6ae99SBarry Smith 
129047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
129147c6ae99SBarry Smith   while (next) {
129247c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
1293*06ebdd98SJed Brown       ierr = DMCompositeAddArray(*fine,next->rank,next->nlocal);CHKERRQ(ierr);
129447c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
129547c6ae99SBarry Smith       ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
129647c6ae99SBarry Smith       ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
129747c6ae99SBarry Smith       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
129847c6ae99SBarry Smith     } else {
129947c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
130047c6ae99SBarry Smith     }
130147c6ae99SBarry Smith     next = next->next;
130247c6ae99SBarry Smith   }
130347c6ae99SBarry Smith   PetscFunctionReturn(0);
130447c6ae99SBarry Smith }
130547c6ae99SBarry Smith 
130647c6ae99SBarry Smith #include "petscmat.h"
130747c6ae99SBarry Smith 
130847c6ae99SBarry Smith struct MatPackLink {
130947c6ae99SBarry Smith   Mat                A;
131047c6ae99SBarry Smith   struct MatPackLink *next;
131147c6ae99SBarry Smith };
131247c6ae99SBarry Smith 
131347c6ae99SBarry Smith struct MatPack {
131447c6ae99SBarry Smith   DM                 right,left;
131547c6ae99SBarry Smith   struct MatPackLink *next;
131647c6ae99SBarry Smith };
131747c6ae99SBarry Smith 
131847c6ae99SBarry Smith #undef __FUNCT__
131947c6ae99SBarry Smith #define __FUNCT__ "MatMultBoth_Shell_Pack"
132047c6ae99SBarry Smith PetscErrorCode MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscBool  add)
132147c6ae99SBarry Smith {
132247c6ae99SBarry Smith   struct MatPack         *mpack;
132347c6ae99SBarry Smith   struct DMCompositeLink *xnext,*ynext;
132447c6ae99SBarry Smith   struct MatPackLink     *anext;
132547c6ae99SBarry Smith   PetscScalar            *xarray,*yarray;
132647c6ae99SBarry Smith   PetscErrorCode         ierr;
132747c6ae99SBarry Smith   PetscInt               i;
132847c6ae99SBarry Smith   Vec                    xglobal,yglobal;
132947c6ae99SBarry Smith   PetscMPIInt            rank;
133047c6ae99SBarry Smith   DM_Composite           *comright;
133147c6ae99SBarry Smith   DM_Composite           *comleft;
133247c6ae99SBarry Smith 
133347c6ae99SBarry Smith   PetscFunctionBegin;
133447c6ae99SBarry Smith   ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
133547c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
133647c6ae99SBarry Smith   comright = (DM_Composite*)mpack->right->data;
133747c6ae99SBarry Smith   comleft = (DM_Composite*)mpack->left->data;
133847c6ae99SBarry Smith   xnext = comright->next;
133947c6ae99SBarry Smith   ynext = comleft->next;
134047c6ae99SBarry Smith   anext = mpack->next;
134147c6ae99SBarry Smith 
134247c6ae99SBarry Smith   while (xnext) {
134347c6ae99SBarry Smith     if (xnext->type == DMCOMPOSITE_ARRAY) {
134447c6ae99SBarry Smith       if (rank == xnext->rank) {
134547c6ae99SBarry Smith         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
134647c6ae99SBarry Smith         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
134747c6ae99SBarry Smith         if (add) {
134847c6ae99SBarry Smith           for (i=0; i<xnext->n; i++) {
134947c6ae99SBarry Smith             yarray[ynext->rstart+i] += xarray[xnext->rstart+i];
135047c6ae99SBarry Smith           }
135147c6ae99SBarry Smith         } else {
135247c6ae99SBarry Smith           ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
135347c6ae99SBarry Smith         }
135447c6ae99SBarry Smith         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
135547c6ae99SBarry Smith         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
135647c6ae99SBarry Smith       }
135747c6ae99SBarry Smith     } else if (xnext->type == DMCOMPOSITE_DM) {
135847c6ae99SBarry Smith       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
135947c6ae99SBarry Smith       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
136047c6ae99SBarry Smith       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
136147c6ae99SBarry Smith       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
136247c6ae99SBarry Smith       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
136347c6ae99SBarry Smith       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
136447c6ae99SBarry Smith       if (add) {
136547c6ae99SBarry Smith         ierr  = MatMultAdd(anext->A,xglobal,yglobal,yglobal);CHKERRQ(ierr);
136647c6ae99SBarry Smith       } else {
136747c6ae99SBarry Smith         ierr  = MatMult(anext->A,xglobal,yglobal);CHKERRQ(ierr);
136847c6ae99SBarry Smith       }
136947c6ae99SBarry Smith       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
137047c6ae99SBarry Smith       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
137147c6ae99SBarry Smith       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
137247c6ae99SBarry Smith       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
137347c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
137447c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
137547c6ae99SBarry Smith       anext = anext->next;
137647c6ae99SBarry Smith     } else {
137747c6ae99SBarry Smith       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
137847c6ae99SBarry Smith     }
137947c6ae99SBarry Smith     xnext = xnext->next;
138047c6ae99SBarry Smith     ynext = ynext->next;
138147c6ae99SBarry Smith   }
138247c6ae99SBarry Smith   PetscFunctionReturn(0);
138347c6ae99SBarry Smith }
138447c6ae99SBarry Smith 
138547c6ae99SBarry Smith #undef __FUNCT__
138647c6ae99SBarry Smith #define __FUNCT__ "MatMultAdd_Shell_Pack"
138747c6ae99SBarry Smith PetscErrorCode MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z)
138847c6ae99SBarry Smith {
138947c6ae99SBarry Smith   PetscErrorCode ierr;
139047c6ae99SBarry Smith   PetscFunctionBegin;
139147c6ae99SBarry Smith   if (z != y) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Handles y == z only");
139247c6ae99SBarry Smith   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
139347c6ae99SBarry Smith   PetscFunctionReturn(0);
139447c6ae99SBarry Smith }
139547c6ae99SBarry Smith 
139647c6ae99SBarry Smith #undef __FUNCT__
139747c6ae99SBarry Smith #define __FUNCT__ "MatMult_Shell_Pack"
139847c6ae99SBarry Smith PetscErrorCode MatMult_Shell_Pack(Mat A,Vec x,Vec y)
139947c6ae99SBarry Smith {
140047c6ae99SBarry Smith   PetscErrorCode ierr;
140147c6ae99SBarry Smith   PetscFunctionBegin;
140247c6ae99SBarry Smith   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
140347c6ae99SBarry Smith   PetscFunctionReturn(0);
140447c6ae99SBarry Smith }
140547c6ae99SBarry Smith 
140647c6ae99SBarry Smith #undef __FUNCT__
140747c6ae99SBarry Smith #define __FUNCT__ "MatMultTranspose_Shell_Pack"
140847c6ae99SBarry Smith PetscErrorCode MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y)
140947c6ae99SBarry Smith {
141047c6ae99SBarry Smith   struct MatPack         *mpack;
141147c6ae99SBarry Smith   struct DMCompositeLink *xnext,*ynext;
141247c6ae99SBarry Smith   struct MatPackLink     *anext;
141347c6ae99SBarry Smith   PetscScalar            *xarray,*yarray;
141447c6ae99SBarry Smith   PetscErrorCode         ierr;
141547c6ae99SBarry Smith   Vec                    xglobal,yglobal;
141647c6ae99SBarry Smith   PetscMPIInt            rank;
141747c6ae99SBarry Smith   DM_Composite           *comright;
141847c6ae99SBarry Smith   DM_Composite           *comleft;
141947c6ae99SBarry Smith 
142047c6ae99SBarry Smith   PetscFunctionBegin;
142147c6ae99SBarry Smith   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
142247c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
142347c6ae99SBarry Smith   comright = (DM_Composite*)mpack->right->data;
142447c6ae99SBarry Smith   comleft = (DM_Composite*)mpack->left->data;
142547c6ae99SBarry Smith   ynext = comright->next;
142647c6ae99SBarry Smith   xnext = comleft->next;
142747c6ae99SBarry Smith   anext = mpack->next;
142847c6ae99SBarry Smith 
142947c6ae99SBarry Smith   while (xnext) {
143047c6ae99SBarry Smith     if (xnext->type == DMCOMPOSITE_ARRAY) {
143147c6ae99SBarry Smith       if (rank == ynext->rank) {
143247c6ae99SBarry Smith         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
143347c6ae99SBarry Smith         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
143447c6ae99SBarry Smith         ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
143547c6ae99SBarry Smith         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
143647c6ae99SBarry Smith         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
143747c6ae99SBarry Smith       }
143847c6ae99SBarry Smith     } else if (xnext->type == DMCOMPOSITE_DM) {
143947c6ae99SBarry Smith       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
144047c6ae99SBarry Smith       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
144147c6ae99SBarry Smith       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
144247c6ae99SBarry Smith       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
144347c6ae99SBarry Smith       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
144447c6ae99SBarry Smith       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
144547c6ae99SBarry Smith       ierr  = MatMultTranspose(anext->A,xglobal,yglobal);CHKERRQ(ierr);
144647c6ae99SBarry Smith       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
144747c6ae99SBarry Smith       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
144847c6ae99SBarry Smith       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
144947c6ae99SBarry Smith       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
145047c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
145147c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
145247c6ae99SBarry Smith       anext = anext->next;
145347c6ae99SBarry Smith     } else {
145447c6ae99SBarry Smith       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
145547c6ae99SBarry Smith     }
145647c6ae99SBarry Smith     xnext = xnext->next;
145747c6ae99SBarry Smith     ynext = ynext->next;
145847c6ae99SBarry Smith   }
145947c6ae99SBarry Smith   PetscFunctionReturn(0);
146047c6ae99SBarry Smith }
146147c6ae99SBarry Smith 
146247c6ae99SBarry Smith #undef __FUNCT__
146347c6ae99SBarry Smith #define __FUNCT__ "MatDestroy_Shell_Pack"
146447c6ae99SBarry Smith PetscErrorCode MatDestroy_Shell_Pack(Mat A)
146547c6ae99SBarry Smith {
146647c6ae99SBarry Smith   struct MatPack     *mpack;
146747c6ae99SBarry Smith   struct MatPackLink *anext,*oldanext;
146847c6ae99SBarry Smith   PetscErrorCode     ierr;
146947c6ae99SBarry Smith 
147047c6ae99SBarry Smith   PetscFunctionBegin;
147147c6ae99SBarry Smith   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
147247c6ae99SBarry Smith   anext = mpack->next;
147347c6ae99SBarry Smith 
147447c6ae99SBarry Smith   while (anext) {
147547c6ae99SBarry Smith     ierr     = MatDestroy(anext->A);CHKERRQ(ierr);
147647c6ae99SBarry Smith     oldanext = anext;
147747c6ae99SBarry Smith     anext    = anext->next;
147847c6ae99SBarry Smith     ierr     = PetscFree(oldanext);CHKERRQ(ierr);
147947c6ae99SBarry Smith   }
148047c6ae99SBarry Smith   ierr = PetscFree(mpack);CHKERRQ(ierr);
148147c6ae99SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
148247c6ae99SBarry Smith   PetscFunctionReturn(0);
148347c6ae99SBarry Smith }
148447c6ae99SBarry Smith 
148547c6ae99SBarry Smith #undef __FUNCT__
14860c010503SBarry Smith #define __FUNCT__ "DMGetInterpolation_Composite"
14870c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
148847c6ae99SBarry Smith {
148947c6ae99SBarry Smith   PetscErrorCode         ierr;
149047c6ae99SBarry Smith   PetscInt               m,n,M,N;
149147c6ae99SBarry Smith   struct DMCompositeLink *nextc;
149247c6ae99SBarry Smith   struct DMCompositeLink *nextf;
149347c6ae99SBarry Smith   struct MatPackLink     *nextmat,*pnextmat = 0;
149447c6ae99SBarry Smith   struct MatPack         *mpack;
149547c6ae99SBarry Smith   Vec                    gcoarse,gfine;
149647c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
149747c6ae99SBarry Smith   DM_Composite           *comfine = (DM_Composite*)fine->data;
149847c6ae99SBarry Smith 
149947c6ae99SBarry Smith   PetscFunctionBegin;
150047c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
150147c6ae99SBarry Smith   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
150247c6ae99SBarry Smith   nextc = comcoarse->next;
150347c6ae99SBarry Smith   nextf = comfine->next;
150447c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
15050c010503SBarry Smith   ierr = DMCreateGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
15060c010503SBarry Smith   ierr = DMCreateGlobalVector(fine,&gfine);CHKERRQ(ierr);
150747c6ae99SBarry Smith   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
150847c6ae99SBarry Smith   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
150947c6ae99SBarry Smith   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
151047c6ae99SBarry Smith   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
151147c6ae99SBarry Smith   ierr = VecDestroy(gcoarse);CHKERRQ(ierr);
151247c6ae99SBarry Smith   ierr = VecDestroy(gfine);CHKERRQ(ierr);
151347c6ae99SBarry Smith 
151447c6ae99SBarry Smith   ierr         = PetscNew(struct MatPack,&mpack);CHKERRQ(ierr);
151547c6ae99SBarry Smith   mpack->right = coarse;
151647c6ae99SBarry Smith   mpack->left  = fine;
151747c6ae99SBarry Smith   ierr  = MatCreate(((PetscObject)fine)->comm,A);CHKERRQ(ierr);
151847c6ae99SBarry Smith   ierr  = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
151947c6ae99SBarry Smith   ierr  = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
152047c6ae99SBarry Smith   ierr  = MatShellSetContext(*A,mpack);CHKERRQ(ierr);
152147c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);CHKERRQ(ierr);
152247c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);CHKERRQ(ierr);
152347c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);CHKERRQ(ierr);
152447c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);CHKERRQ(ierr);
152547c6ae99SBarry Smith 
152647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
152747c6ae99SBarry Smith   while (nextc) {
152847c6ae99SBarry Smith     if (nextc->type != nextf->type) SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Two DM have different layout");
152947c6ae99SBarry Smith 
153047c6ae99SBarry Smith     if (nextc->type == DMCOMPOSITE_ARRAY) {
153147c6ae99SBarry Smith       ;
153247c6ae99SBarry Smith     } else if (nextc->type == DMCOMPOSITE_DM) {
153347c6ae99SBarry Smith       ierr          = PetscNew(struct MatPackLink,&nextmat);CHKERRQ(ierr);
153447c6ae99SBarry Smith       nextmat->next = 0;
153547c6ae99SBarry Smith       if (pnextmat) {
153647c6ae99SBarry Smith         pnextmat->next = nextmat;
153747c6ae99SBarry Smith         pnextmat       = nextmat;
153847c6ae99SBarry Smith       } else {
153947c6ae99SBarry Smith         pnextmat    = nextmat;
154047c6ae99SBarry Smith         mpack->next = nextmat;
154147c6ae99SBarry Smith       }
154247c6ae99SBarry Smith       ierr = DMGetInterpolation(nextc->dm,nextf->dm,&nextmat->A,PETSC_NULL);CHKERRQ(ierr);
154347c6ae99SBarry Smith     } else {
154447c6ae99SBarry Smith       SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
154547c6ae99SBarry Smith     }
154647c6ae99SBarry Smith     nextc = nextc->next;
154747c6ae99SBarry Smith     nextf = nextf->next;
154847c6ae99SBarry Smith   }
154947c6ae99SBarry Smith   PetscFunctionReturn(0);
155047c6ae99SBarry Smith }
155147c6ae99SBarry Smith 
155247c6ae99SBarry Smith #undef __FUNCT__
15530c010503SBarry Smith #define __FUNCT__ "DMGetMatrix_Composite"
15540c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix_Composite(DM dm, const MatType mtype,Mat *J)
155547c6ae99SBarry Smith {
155647c6ae99SBarry Smith   PetscErrorCode         ierr;
155747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
155847c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
155947c6ae99SBarry Smith   PetscInt               m,*dnz,*onz,i,j,mA;
156047c6ae99SBarry Smith   Mat                    Atmp;
156147c6ae99SBarry Smith   PetscMPIInt            rank;
156247c6ae99SBarry Smith   PetscScalar            zero = 0.0;
156347c6ae99SBarry Smith   PetscBool              dense = PETSC_FALSE;
15641411c6eeSJed Brown   ISLocalToGlobalMapping ltogmap,ltogmapb;
156547c6ae99SBarry Smith 
156647c6ae99SBarry Smith   PetscFunctionBegin;
156747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
156847c6ae99SBarry Smith 
156947c6ae99SBarry Smith   /* use global vector to determine layout needed for matrix */
157047c6ae99SBarry Smith   m = com->n;
157147c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
157247c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr);
157347c6ae99SBarry Smith   ierr = MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
157447c6ae99SBarry Smith   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
157547c6ae99SBarry Smith 
157647c6ae99SBarry Smith   /*
157747c6ae99SBarry Smith      Extremely inefficient but will compute entire Jacobian for testing
157847c6ae99SBarry Smith   */
1579671f6225SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
158047c6ae99SBarry Smith   if (dense) {
158147c6ae99SBarry Smith     PetscInt    rstart,rend,*indices;
158247c6ae99SBarry Smith     PetscScalar *values;
158347c6ae99SBarry Smith 
158447c6ae99SBarry Smith     mA    = com->N;
158547c6ae99SBarry Smith     ierr = MatMPIAIJSetPreallocation(*J,mA,PETSC_NULL,mA-m,PETSC_NULL);CHKERRQ(ierr);
158647c6ae99SBarry Smith     ierr = MatSeqAIJSetPreallocation(*J,mA,PETSC_NULL);CHKERRQ(ierr);
158747c6ae99SBarry Smith 
158847c6ae99SBarry Smith     ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr);
158947c6ae99SBarry Smith     ierr = PetscMalloc2(mA,PetscScalar,&values,mA,PetscInt,&indices);CHKERRQ(ierr);
159047c6ae99SBarry Smith     ierr = PetscMemzero(values,mA*sizeof(PetscScalar));CHKERRQ(ierr);
159147c6ae99SBarry Smith     for (i=0; i<mA; i++) indices[i] = i;
159247c6ae99SBarry Smith     for (i=rstart; i<rend; i++) {
159347c6ae99SBarry Smith       ierr = MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES);CHKERRQ(ierr);
159447c6ae99SBarry Smith     }
159547c6ae99SBarry Smith     ierr = PetscFree2(values,indices);CHKERRQ(ierr);
159647c6ae99SBarry Smith     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
159747c6ae99SBarry Smith     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
159847c6ae99SBarry Smith     PetscFunctionReturn(0);
159947c6ae99SBarry Smith   }
160047c6ae99SBarry Smith 
160147c6ae99SBarry Smith   ierr = MatPreallocateInitialize(((PetscObject)dm)->comm,m,m,dnz,onz);CHKERRQ(ierr);
160247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
160347c6ae99SBarry Smith   next = com->next;
160447c6ae99SBarry Smith   while (next) {
160547c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
160647c6ae99SBarry Smith       if (rank == next->rank) {  /* zero the "little" block */
160747c6ae99SBarry Smith         for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) {
160847c6ae99SBarry Smith           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
160947c6ae99SBarry Smith             ierr = MatPreallocateSet(j,1,&i,dnz,onz);CHKERRQ(ierr);
161047c6ae99SBarry Smith           }
161147c6ae99SBarry Smith         }
161247c6ae99SBarry Smith       }
161347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
161447c6ae99SBarry Smith       PetscInt       nc,rstart,*ccols,maxnc;
161547c6ae99SBarry Smith       const PetscInt *cols,*rstarts;
161647c6ae99SBarry Smith       PetscMPIInt    proc;
161747c6ae99SBarry Smith 
161847c6ae99SBarry Smith       ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr);
161947c6ae99SBarry Smith       ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr);
162047c6ae99SBarry Smith       ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
162147c6ae99SBarry Smith       ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr);
162247c6ae99SBarry Smith 
162347c6ae99SBarry Smith       maxnc = 0;
162447c6ae99SBarry Smith       for (i=0; i<mA; i++) {
162547c6ae99SBarry Smith         ierr  = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
162647c6ae99SBarry Smith         ierr  = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
162747c6ae99SBarry Smith         maxnc = PetscMax(nc,maxnc);
162847c6ae99SBarry Smith       }
162947c6ae99SBarry Smith       ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
163047c6ae99SBarry Smith       for (i=0; i<mA; i++) {
163147c6ae99SBarry Smith         ierr = MatGetRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr);
163247c6ae99SBarry Smith         /* remap the columns taking into how much they are shifted on each process */
163347c6ae99SBarry Smith         for (j=0; j<nc; j++) {
163447c6ae99SBarry Smith           proc = 0;
163547c6ae99SBarry Smith           while (cols[j] >= rstarts[proc+1]) proc++;
163647c6ae99SBarry Smith           ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
163747c6ae99SBarry Smith         }
163847c6ae99SBarry Smith         ierr = MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz);CHKERRQ(ierr);
163947c6ae99SBarry Smith         ierr = MatRestoreRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr);
164047c6ae99SBarry Smith       }
164147c6ae99SBarry Smith       ierr = PetscFree(ccols);CHKERRQ(ierr);
164247c6ae99SBarry Smith       ierr = MatDestroy(Atmp);CHKERRQ(ierr);
164347c6ae99SBarry Smith     } else {
164447c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
164547c6ae99SBarry Smith     }
164647c6ae99SBarry Smith     next = next->next;
164747c6ae99SBarry Smith   }
164847c6ae99SBarry Smith   if (com->FormCoupleLocations) {
164947c6ae99SBarry Smith     ierr = (*com->FormCoupleLocations)(dm,PETSC_NULL,dnz,onz,__rstart,__nrows,__start,__end);CHKERRQ(ierr);
165047c6ae99SBarry Smith   }
165147c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);CHKERRQ(ierr);
165247c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(*J,0,dnz);CHKERRQ(ierr);
165347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
165447c6ae99SBarry Smith 
165547c6ae99SBarry Smith   next = com->next;
165647c6ae99SBarry Smith   while (next) {
165747c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
165847c6ae99SBarry Smith       if (rank == next->rank) {
165947c6ae99SBarry Smith         for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) {
166047c6ae99SBarry Smith           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
166147c6ae99SBarry Smith             ierr = MatSetValues(*J,1,&j,1,&i,&zero,INSERT_VALUES);CHKERRQ(ierr);
166247c6ae99SBarry Smith           }
166347c6ae99SBarry Smith         }
166447c6ae99SBarry Smith       }
166547c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
166647c6ae99SBarry Smith       PetscInt          nc,rstart,row,maxnc,*ccols;
166747c6ae99SBarry Smith       const PetscInt    *cols,*rstarts;
166847c6ae99SBarry Smith       const PetscScalar *values;
166947c6ae99SBarry Smith       PetscMPIInt       proc;
167047c6ae99SBarry Smith 
167147c6ae99SBarry Smith       ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr);
167247c6ae99SBarry Smith       ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr);
167347c6ae99SBarry Smith       ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
167447c6ae99SBarry Smith       ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr);
167547c6ae99SBarry Smith       maxnc = 0;
167647c6ae99SBarry Smith       for (i=0; i<mA; i++) {
167747c6ae99SBarry Smith         ierr  = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
167847c6ae99SBarry Smith         ierr  = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
167947c6ae99SBarry Smith         maxnc = PetscMax(nc,maxnc);
168047c6ae99SBarry Smith       }
168147c6ae99SBarry Smith       ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
168247c6ae99SBarry Smith       for (i=0; i<mA; i++) {
168347c6ae99SBarry Smith         ierr = MatGetRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr);
168447c6ae99SBarry Smith         for (j=0; j<nc; j++) {
168547c6ae99SBarry Smith           proc = 0;
168647c6ae99SBarry Smith           while (cols[j] >= rstarts[proc+1]) proc++;
168747c6ae99SBarry Smith           ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
168847c6ae99SBarry Smith         }
168947c6ae99SBarry Smith         row  = com->rstart+next->rstart+i;
169047c6ae99SBarry Smith         ierr = MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES);CHKERRQ(ierr);
169147c6ae99SBarry Smith         ierr = MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr);
169247c6ae99SBarry Smith       }
169347c6ae99SBarry Smith       ierr = PetscFree(ccols);CHKERRQ(ierr);
169447c6ae99SBarry Smith       ierr = MatDestroy(Atmp);CHKERRQ(ierr);
169547c6ae99SBarry Smith     } else {
169647c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
169747c6ae99SBarry Smith     }
169847c6ae99SBarry Smith     next = next->next;
169947c6ae99SBarry Smith   }
170047c6ae99SBarry Smith   if (com->FormCoupleLocations) {
170147c6ae99SBarry Smith     PetscInt __rstart;
170247c6ae99SBarry Smith     ierr = MatGetOwnershipRange(*J,&__rstart,PETSC_NULL);CHKERRQ(ierr);
170347c6ae99SBarry Smith     ierr = (*com->FormCoupleLocations)(dm,*J,PETSC_NULL,PETSC_NULL,__rstart,0,0,0);CHKERRQ(ierr);
170447c6ae99SBarry Smith   }
170547c6ae99SBarry Smith   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
170647c6ae99SBarry Smith   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17071411c6eeSJed Brown 
17081411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(dm,&ltogmap);CHKERRQ(ierr);
17091411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(dm,&ltogmapb);CHKERRQ(ierr);
17101411c6eeSJed Brown   ierr = MatSetLocalToGlobalMapping(*J,ltogmap,ltogmap);CHKERRQ(ierr);
17111411c6eeSJed Brown   ierr = MatSetLocalToGlobalMappingBlock(*J,ltogmapb,ltogmapb);CHKERRQ(ierr);
171247c6ae99SBarry Smith   PetscFunctionReturn(0);
171347c6ae99SBarry Smith }
171447c6ae99SBarry Smith 
171547c6ae99SBarry Smith #undef __FUNCT__
17161411c6eeSJed Brown #define __FUNCT__ "DMCreateLocalToGlobalMapping_Composite"
17171411c6eeSJed Brown static PetscErrorCode DMCreateLocalToGlobalMapping_Composite(DM dm)
17181411c6eeSJed Brown {
17191411c6eeSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
17201411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1721a03f9cedSJed Brown   PetscInt               i,cnt,m,*idx;
17221411c6eeSJed Brown   PetscErrorCode         ierr;
17231411c6eeSJed Brown 
17241411c6eeSJed Brown   PetscFunctionBegin;
17251411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
17261411c6eeSJed Brown   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
17271411c6eeSJed Brown   for (cnt=0,i=0; i<(com->nDM+com->nredundant); i++) {
1728a03f9cedSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1729a03f9cedSJed Brown     cnt += m;
1730a03f9cedSJed Brown   }
1731a03f9cedSJed Brown   ierr = PetscMalloc(cnt*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1732a03f9cedSJed Brown   for (cnt=0,i=0; i<(com->nDM+com->nredundant); i++) {
17331411c6eeSJed Brown     const PetscInt *subidx;
17341411c6eeSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
17351411c6eeSJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
17361411c6eeSJed Brown     ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr);
17371411c6eeSJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
17381411c6eeSJed Brown     cnt += m;
17391411c6eeSJed Brown   }
17401411c6eeSJed Brown   ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,cnt,idx,PETSC_OWN_POINTER,&dm->ltogmap);CHKERRQ(ierr);
17411411c6eeSJed Brown   for (i=0; i<com->nDM+com->nredundant; i++) {ierr = ISLocalToGlobalMappingDestroy(ltogs[i]);CHKERRQ(ierr);}
17421411c6eeSJed Brown   ierr = PetscFree(ltogs);CHKERRQ(ierr);
17431411c6eeSJed Brown   PetscFunctionReturn(0);
17441411c6eeSJed Brown }
17451411c6eeSJed Brown 
17461411c6eeSJed Brown 
17471411c6eeSJed Brown #undef __FUNCT__
17480c010503SBarry Smith #define __FUNCT__ "DMGetColoring_Composite"
17490c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
175047c6ae99SBarry Smith {
175147c6ae99SBarry Smith   PetscErrorCode         ierr;
175247c6ae99SBarry Smith   PetscInt               n,i,cnt;
175347c6ae99SBarry Smith   ISColoringValue        *colors;
175447c6ae99SBarry Smith   PetscBool              dense = PETSC_FALSE;
175547c6ae99SBarry Smith   ISColoringValue        maxcol = 0;
175647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
175747c6ae99SBarry Smith 
175847c6ae99SBarry Smith   PetscFunctionBegin;
175947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
176047c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED) {
176147c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" );
176247c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL) {
176347c6ae99SBarry Smith     n = com->n;
176447c6ae99SBarry Smith   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
176547c6ae99SBarry Smith   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
176647c6ae99SBarry Smith 
1767671f6225SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
176847c6ae99SBarry Smith   if (dense) {
176947c6ae99SBarry Smith     for (i=0; i<n; i++) {
177047c6ae99SBarry Smith       colors[i] = (ISColoringValue)(com->rstart + i);
177147c6ae99SBarry Smith     }
177247c6ae99SBarry Smith     maxcol = com->N;
177347c6ae99SBarry Smith   } else {
177447c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
177547c6ae99SBarry Smith     PetscMPIInt            rank;
177647c6ae99SBarry Smith 
177747c6ae99SBarry Smith     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
177847c6ae99SBarry Smith     cnt  = 0;
177947c6ae99SBarry Smith     while (next) {
178047c6ae99SBarry Smith       if (next->type == DMCOMPOSITE_ARRAY) {
178147c6ae99SBarry Smith         if (rank == next->rank) {  /* each column gets is own color */
178247c6ae99SBarry Smith           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
178347c6ae99SBarry Smith             colors[cnt++] = maxcol++;
178447c6ae99SBarry Smith           }
178547c6ae99SBarry Smith         }
178647c6ae99SBarry Smith         ierr = MPI_Bcast(&maxcol,1,MPIU_COLORING_VALUE,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
178747c6ae99SBarry Smith       } else if (next->type == DMCOMPOSITE_DM) {
178847c6ae99SBarry Smith         ISColoring     lcoloring;
178947c6ae99SBarry Smith 
179047c6ae99SBarry Smith         ierr = DMGetColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
179147c6ae99SBarry Smith         for (i=0; i<lcoloring->N; i++) {
179247c6ae99SBarry Smith           colors[cnt++] = maxcol + lcoloring->colors[i];
179347c6ae99SBarry Smith         }
179447c6ae99SBarry Smith         maxcol += lcoloring->n;
179547c6ae99SBarry Smith         ierr = ISColoringDestroy(lcoloring);CHKERRQ(ierr);
179647c6ae99SBarry Smith       } else {
179747c6ae99SBarry Smith         SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
179847c6ae99SBarry Smith       }
179947c6ae99SBarry Smith       next = next->next;
180047c6ae99SBarry Smith     }
180147c6ae99SBarry Smith   }
180247c6ae99SBarry Smith   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
180347c6ae99SBarry Smith   PetscFunctionReturn(0);
180447c6ae99SBarry Smith }
180547c6ae99SBarry Smith 
180647c6ae99SBarry Smith #undef __FUNCT__
18070c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
18080c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
180947c6ae99SBarry Smith {
181047c6ae99SBarry Smith   PetscErrorCode         ierr;
181147c6ae99SBarry Smith   struct DMCompositeLink *next;
181247c6ae99SBarry Smith   PetscInt               cnt = 3;
181347c6ae99SBarry Smith   PetscMPIInt            rank;
181447c6ae99SBarry Smith   PetscScalar            *garray,*larray;
181547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
181647c6ae99SBarry Smith 
181747c6ae99SBarry Smith   PetscFunctionBegin;
181847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
181947c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
182047c6ae99SBarry Smith   next = com->next;
182147c6ae99SBarry Smith   if (!com->setup) {
1822d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
182347c6ae99SBarry Smith   }
182447c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
182547c6ae99SBarry Smith   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
182647c6ae99SBarry Smith   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
182747c6ae99SBarry Smith 
182847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
182947c6ae99SBarry Smith   while (next) {
183047c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
183147c6ae99SBarry Smith       if (rank == next->rank) {
183247c6ae99SBarry Smith         ierr    = PetscMemcpy(larray,garray,next->n*sizeof(PetscScalar));CHKERRQ(ierr);
183347c6ae99SBarry Smith         garray += next->n;
183447c6ae99SBarry Smith       }
183547c6ae99SBarry Smith       /* does not handle ADD_VALUES */
1836*06ebdd98SJed Brown       ierr = MPI_Bcast(larray,next->nlocal,MPIU_SCALAR,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
183747c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
183847c6ae99SBarry Smith       Vec      local,global;
183947c6ae99SBarry Smith       PetscInt N;
184047c6ae99SBarry Smith 
184147c6ae99SBarry Smith       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
184247c6ae99SBarry Smith       ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
184347c6ae99SBarry Smith       ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
184447c6ae99SBarry Smith       ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
184547c6ae99SBarry Smith       ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
184647c6ae99SBarry Smith       ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
184747c6ae99SBarry Smith       ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
184847c6ae99SBarry Smith       ierr = VecResetArray(global);CHKERRQ(ierr);
184947c6ae99SBarry Smith       ierr = VecResetArray(local);CHKERRQ(ierr);
185047c6ae99SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
185147c6ae99SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
185247c6ae99SBarry Smith     } else {
185347c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
185447c6ae99SBarry Smith     }
185547c6ae99SBarry Smith     cnt++;
1856*06ebdd98SJed Brown     larray += next->nlocal;
185747c6ae99SBarry Smith     next    = next->next;
185847c6ae99SBarry Smith   }
185947c6ae99SBarry Smith 
186047c6ae99SBarry Smith   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
186147c6ae99SBarry Smith   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
186247c6ae99SBarry Smith   PetscFunctionReturn(0);
186347c6ae99SBarry Smith }
186447c6ae99SBarry Smith 
186547c6ae99SBarry Smith #undef __FUNCT__
18660c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
18670c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
18680c010503SBarry Smith {
18690c010503SBarry Smith   PetscFunctionBegin;
18700c010503SBarry Smith   PetscFunctionReturn(0);
18710c010503SBarry Smith }
187247c6ae99SBarry Smith 
1873a4121054SBarry Smith EXTERN_C_BEGIN
1874a4121054SBarry Smith #undef __FUNCT__
1875a4121054SBarry Smith #define __FUNCT__ "DMCreate_Composite"
1876a4121054SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreate_Composite(DM p)
1877a4121054SBarry Smith {
1878a4121054SBarry Smith   PetscErrorCode ierr;
1879a4121054SBarry Smith   DM_Composite   *com;
1880a4121054SBarry Smith 
1881a4121054SBarry Smith   PetscFunctionBegin;
1882a4121054SBarry Smith   ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1883a4121054SBarry Smith   p->data = com;
1884a4121054SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1885a4121054SBarry Smith   com->n            = 0;
1886a4121054SBarry Smith   com->next         = PETSC_NULL;
1887a4121054SBarry Smith   com->nredundant   = 0;
1888a4121054SBarry Smith   com->nDM          = 0;
1889a4121054SBarry Smith 
1890a4121054SBarry Smith   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1891a4121054SBarry Smith   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
18921411c6eeSJed Brown   p->ops->createlocaltoglobalmapping      = DMCreateLocalToGlobalMapping_Composite;
18931411c6eeSJed Brown   p->ops->createlocaltoglobalmappingblock = 0;
1894a4121054SBarry Smith   p->ops->refine                          = DMRefine_Composite;
1895a4121054SBarry Smith   p->ops->getinterpolation                = DMGetInterpolation_Composite;
1896a4121054SBarry Smith   p->ops->getmatrix                       = DMGetMatrix_Composite;
1897a4121054SBarry Smith   p->ops->getcoloring                     = DMGetColoring_Composite;
1898a4121054SBarry Smith   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1899a4121054SBarry Smith   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1900a4121054SBarry Smith   p->ops->destroy                         = DMDestroy_Composite;
1901a4121054SBarry Smith   p->ops->view                            = DMView_Composite;
1902a4121054SBarry Smith   p->ops->setup                           = DMSetUp_Composite;
1903a4121054SBarry Smith   PetscFunctionReturn(0);
1904a4121054SBarry Smith }
1905a4121054SBarry Smith EXTERN_C_END
1906a4121054SBarry Smith 
19070c010503SBarry Smith #undef __FUNCT__
19080c010503SBarry Smith #define __FUNCT__ "DMCompositeCreate"
19090c010503SBarry Smith /*@C
19100c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
19110c010503SBarry Smith       vectors made up of several subvectors.
19120c010503SBarry Smith 
19130c010503SBarry Smith     Collective on MPI_Comm
191447c6ae99SBarry Smith 
191547c6ae99SBarry Smith     Input Parameter:
19160c010503SBarry Smith .   comm - the processors that will share the global vector
19170c010503SBarry Smith 
19180c010503SBarry Smith     Output Parameters:
19190c010503SBarry Smith .   packer - the packer object
192047c6ae99SBarry Smith 
192147c6ae99SBarry Smith     Level: advanced
192247c6ae99SBarry Smith 
19230c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCompositeScatter(),
19246eb61c8cSJed Brown          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
192547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
192647c6ae99SBarry Smith 
192747c6ae99SBarry Smith @*/
19280c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCompositeCreate(MPI_Comm comm,DM *packer)
192947c6ae99SBarry Smith {
19300c010503SBarry Smith   PetscErrorCode ierr;
19310c010503SBarry Smith 
193247c6ae99SBarry Smith   PetscFunctionBegin;
19330c010503SBarry Smith   PetscValidPointer(packer,2);
1934a4121054SBarry Smith   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1935a4121054SBarry Smith   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
193647c6ae99SBarry Smith   PetscFunctionReturn(0);
193747c6ae99SBarry Smith }
1938