xref: /petsc/src/dm/impls/composite/pack.c (revision 7087cfbefd1a42b179f217f9994fb6cb0d0c1824)
147c6ae99SBarry Smith #define PETSCDM_DLL
247c6ae99SBarry Smith 
3ec9811eeSJed Brown #include "packimpl.h" /*I   "petscdm.h"   I*/
447c6ae99SBarry Smith 
547c6ae99SBarry Smith #undef __FUNCT__
647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeSetCoupling"
747c6ae99SBarry Smith /*@C
847c6ae99SBarry Smith     DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
9aa219208SBarry Smith       seperate components (DMDA's and arrays) in a DMto build the correct matrix nonzero structure.
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith 
1247c6ae99SBarry Smith     Logically Collective on MPI_Comm
1347c6ae99SBarry Smith 
1447c6ae99SBarry Smith     Input Parameter:
1547c6ae99SBarry Smith +   dm - the composite object
1647c6ae99SBarry Smith -   formcouplelocations - routine to set the nonzero locations in the matrix
1747c6ae99SBarry Smith 
1847c6ae99SBarry Smith     Level: advanced
1947c6ae99SBarry Smith 
2047c6ae99SBarry Smith     Notes: See DMCompositeSetContext() and DMCompositeGetContext() for how to get user information into
2147c6ae99SBarry Smith         this routine
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith @*/
24*7087cfbeSBarry Smith PetscErrorCode  DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt))
2547c6ae99SBarry Smith {
2647c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith   PetscFunctionBegin;
2947c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
3047c6ae99SBarry Smith   PetscFunctionReturn(0);
3147c6ae99SBarry Smith }
3247c6ae99SBarry Smith 
3347c6ae99SBarry Smith #undef __FUNCT__
3447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeSetContext"
3547c6ae99SBarry Smith /*@
3647c6ae99SBarry Smith     DMCompositeSetContext - Allows user to stash data they may need within the form coupling routine they
3747c6ae99SBarry Smith       set with DMCompositeSetCoupling()
3847c6ae99SBarry Smith 
3947c6ae99SBarry Smith 
4047c6ae99SBarry Smith     Not Collective
4147c6ae99SBarry Smith 
4247c6ae99SBarry Smith     Input Parameter:
4347c6ae99SBarry Smith +   dm - the composite object
4447c6ae99SBarry Smith -   ctx - the user supplied context
4547c6ae99SBarry Smith 
4647c6ae99SBarry Smith     Level: advanced
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith     Notes: Use DMCompositeGetContext() to retrieve the context when needed.
4947c6ae99SBarry Smith 
5047c6ae99SBarry Smith @*/
51*7087cfbeSBarry Smith PetscErrorCode  DMCompositeSetContext(DM dm,void *ctx)
5247c6ae99SBarry Smith {
5347c6ae99SBarry Smith   PetscFunctionBegin;
5447c6ae99SBarry Smith   dm->ctx = ctx;
5547c6ae99SBarry Smith   PetscFunctionReturn(0);
5647c6ae99SBarry Smith }
5747c6ae99SBarry Smith 
5847c6ae99SBarry Smith #undef __FUNCT__
5947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetContext"
6047c6ae99SBarry Smith /*@
6147c6ae99SBarry Smith     DMCompositeGetContext - Access the context set with DMCompositeSetContext()
6247c6ae99SBarry Smith 
6347c6ae99SBarry Smith 
6447c6ae99SBarry Smith     Not Collective
6547c6ae99SBarry Smith 
6647c6ae99SBarry Smith     Input Parameter:
6747c6ae99SBarry Smith .   dm - the composite object
6847c6ae99SBarry Smith 
6947c6ae99SBarry Smith     Output Parameter:
7047c6ae99SBarry Smith .    ctx - the user supplied context
7147c6ae99SBarry Smith 
7247c6ae99SBarry Smith     Level: advanced
7347c6ae99SBarry Smith 
7447c6ae99SBarry Smith     Notes: Use DMCompositeGetContext() to retrieve the context when needed.
7547c6ae99SBarry Smith 
7647c6ae99SBarry Smith @*/
77*7087cfbeSBarry Smith PetscErrorCode  DMCompositeGetContext(DM dm,void **ctx)
7847c6ae99SBarry Smith {
7947c6ae99SBarry Smith   PetscFunctionBegin;
8047c6ae99SBarry Smith   *ctx = dm->ctx;
8147c6ae99SBarry Smith   PetscFunctionReturn(0);
8247c6ae99SBarry Smith }
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith 
8647c6ae99SBarry Smith extern PetscErrorCode DMDestroy_Private(DM,PetscBool *);
8747c6ae99SBarry Smith 
8847c6ae99SBarry Smith #undef __FUNCT__
890c010503SBarry Smith #define __FUNCT__ "DMDestroy_Composite"
90*7087cfbeSBarry Smith PetscErrorCode  DMDestroy_Composite(DM dm)
9147c6ae99SBarry Smith {
9247c6ae99SBarry Smith   PetscErrorCode         ierr;
9347c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
9447c6ae99SBarry Smith   PetscBool              done;
9547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
9647c6ae99SBarry Smith 
9747c6ae99SBarry Smith   PetscFunctionBegin;
9847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
9947c6ae99SBarry Smith   ierr = DMDestroy_Private((DM)dm,&done);CHKERRQ(ierr);
10047c6ae99SBarry Smith   if (!done) PetscFunctionReturn(0);
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith   next = com->next;
10347c6ae99SBarry Smith   while (next) {
10447c6ae99SBarry Smith     prev = next;
10547c6ae99SBarry Smith     next = next->next;
10647c6ae99SBarry Smith     if (prev->type == DMCOMPOSITE_DM) {
10747c6ae99SBarry Smith       ierr = DMDestroy(prev->dm);CHKERRQ(ierr);
10847c6ae99SBarry Smith     }
10947c6ae99SBarry Smith     if (prev->grstarts) {
11047c6ae99SBarry Smith       ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
11147c6ae99SBarry Smith     }
11247c6ae99SBarry Smith     ierr = PetscFree(prev);CHKERRQ(ierr);
11347c6ae99SBarry Smith   }
11447c6ae99SBarry Smith   ierr = PetscFree(com);CHKERRQ(ierr);
11547c6ae99SBarry Smith   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
11647c6ae99SBarry Smith   PetscFunctionReturn(0);
11747c6ae99SBarry Smith }
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith #undef __FUNCT__
1200c010503SBarry Smith #define __FUNCT__ "DMView_Composite"
121*7087cfbeSBarry Smith PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
12247c6ae99SBarry Smith {
12347c6ae99SBarry Smith   PetscErrorCode ierr;
12447c6ae99SBarry Smith   PetscBool      iascii;
12547c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite *)dm->data;
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith   PetscFunctionBegin;
12847c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
12947c6ae99SBarry Smith   if (iascii) {
13047c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
13147c6ae99SBarry Smith     PetscInt               i;
13247c6ae99SBarry Smith 
13347c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix?((PetscObject)dm)->prefix:"no prefix");CHKERRQ(ierr);
13447c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"  contains %d DMs and %d redundant arrays\n",com->nDM,com->nredundant);CHKERRQ(ierr);
13547c6ae99SBarry Smith     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
13647c6ae99SBarry Smith     for (i=0; lnk; lnk=lnk->next,i++) {
13747c6ae99SBarry Smith       if (lnk->dm) {
13847c6ae99SBarry Smith         ierr = PetscViewerASCIIPrintf(v,"Link %d: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
13947c6ae99SBarry Smith         ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
14047c6ae99SBarry Smith         ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
14147c6ae99SBarry Smith         ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
14247c6ae99SBarry Smith       } else {
14306ebdd98SJed Brown         ierr = PetscViewerASCIIPrintf(v,"Link %d: Redundant array of size %d owned by rank %d\n",i,lnk->nlocal,lnk->rank);CHKERRQ(ierr);
14447c6ae99SBarry Smith       }
14547c6ae99SBarry Smith     }
14647c6ae99SBarry Smith     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
14747c6ae99SBarry Smith   }
14847c6ae99SBarry Smith   PetscFunctionReturn(0);
14947c6ae99SBarry Smith }
15047c6ae99SBarry Smith 
15147c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
15247c6ae99SBarry Smith #undef __FUNCT__
153d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp_Composite"
154*7087cfbeSBarry Smith PetscErrorCode  DMSetUp_Composite(DM dm)
15547c6ae99SBarry Smith {
15647c6ae99SBarry Smith   PetscErrorCode         ierr;
15747c6ae99SBarry Smith   PetscInt               nprev = 0;
15847c6ae99SBarry Smith   PetscMPIInt            rank,size;
15947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
16047c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
16147c6ae99SBarry Smith   PetscLayout            map;
16247c6ae99SBarry Smith 
16347c6ae99SBarry Smith   PetscFunctionBegin;
16447c6ae99SBarry Smith   if (com->setup) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
16547c6ae99SBarry Smith   ierr = PetscLayoutCreate(((PetscObject)dm)->comm,&map);CHKERRQ(ierr);
16647c6ae99SBarry Smith   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
16747c6ae99SBarry Smith   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
16847c6ae99SBarry Smith   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
16947c6ae99SBarry Smith   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
17047c6ae99SBarry Smith   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
17147c6ae99SBarry Smith   ierr = PetscLayoutGetRange(map,&com->rstart,PETSC_NULL);CHKERRQ(ierr);
17247c6ae99SBarry Smith   ierr = PetscLayoutDestroy(map);CHKERRQ(ierr);
17347c6ae99SBarry Smith 
17447c6ae99SBarry Smith   /* now set the rstart for each linked array/vector */
17547c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
17647c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)dm)->comm,&size);CHKERRQ(ierr);
17747c6ae99SBarry Smith   while (next) {
17847c6ae99SBarry Smith     next->rstart = nprev;
17906ebdd98SJed Brown     nprev += next->n;
18047c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
18147c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
18247c6ae99SBarry Smith       ierr = MPI_Bcast(&next->grstart,1,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
18347c6ae99SBarry Smith     } else {
18447c6ae99SBarry Smith       ierr = PetscMalloc(size*sizeof(PetscInt),&next->grstarts);CHKERRQ(ierr);
18547c6ae99SBarry Smith       ierr = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
18647c6ae99SBarry Smith     }
18747c6ae99SBarry Smith     next = next->next;
18847c6ae99SBarry Smith   }
18947c6ae99SBarry Smith   com->setup = PETSC_TRUE;
19047c6ae99SBarry Smith   PetscFunctionReturn(0);
19147c6ae99SBarry Smith }
19247c6ae99SBarry Smith 
19347c6ae99SBarry Smith 
19447c6ae99SBarry Smith #undef __FUNCT__
19547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess_Array"
19647c6ae99SBarry Smith PetscErrorCode DMCompositeGetAccess_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar **array)
19747c6ae99SBarry Smith {
19847c6ae99SBarry Smith   PetscErrorCode ierr;
19947c6ae99SBarry Smith   PetscScalar    *varray;
20047c6ae99SBarry Smith   PetscMPIInt    rank;
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith   PetscFunctionBegin;
20347c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
20447c6ae99SBarry Smith   if (array) {
20547c6ae99SBarry Smith     if (rank == mine->rank) {
20647c6ae99SBarry Smith       ierr    = VecGetArray(vec,&varray);CHKERRQ(ierr);
20747c6ae99SBarry Smith       *array  = varray + mine->rstart;
20847c6ae99SBarry Smith       ierr    = VecRestoreArray(vec,&varray);CHKERRQ(ierr);
20947c6ae99SBarry Smith     } else {
21047c6ae99SBarry Smith       *array = 0;
21147c6ae99SBarry Smith     }
21247c6ae99SBarry Smith   }
21347c6ae99SBarry Smith   PetscFunctionReturn(0);
21447c6ae99SBarry Smith }
21547c6ae99SBarry Smith 
21647c6ae99SBarry Smith #undef __FUNCT__
21747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess_DM"
21847c6ae99SBarry Smith PetscErrorCode DMCompositeGetAccess_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec *global)
21947c6ae99SBarry Smith {
22047c6ae99SBarry Smith   PetscErrorCode ierr;
22147c6ae99SBarry Smith   PetscScalar    *array;
22247c6ae99SBarry Smith 
22347c6ae99SBarry Smith   PetscFunctionBegin;
22447c6ae99SBarry Smith   if (global) {
22547c6ae99SBarry Smith     ierr    = DMGetGlobalVector(mine->dm,global);CHKERRQ(ierr);
22647c6ae99SBarry Smith     ierr    = VecGetArray(vec,&array);CHKERRQ(ierr);
22747c6ae99SBarry Smith     ierr    = VecPlaceArray(*global,array+mine->rstart);CHKERRQ(ierr);
22847c6ae99SBarry Smith     ierr    = VecRestoreArray(vec,&array);CHKERRQ(ierr);
22947c6ae99SBarry Smith   }
23047c6ae99SBarry Smith   PetscFunctionReturn(0);
23147c6ae99SBarry Smith }
23247c6ae99SBarry Smith 
23347c6ae99SBarry Smith #undef __FUNCT__
23447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess_Array"
23547c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreAccess_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar **array)
23647c6ae99SBarry Smith {
23747c6ae99SBarry Smith   PetscFunctionBegin;
23847c6ae99SBarry Smith   PetscFunctionReturn(0);
23947c6ae99SBarry Smith }
24047c6ae99SBarry Smith 
24147c6ae99SBarry Smith #undef __FUNCT__
24247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess_DM"
24347c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreAccess_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec *global)
24447c6ae99SBarry Smith {
24547c6ae99SBarry Smith   PetscErrorCode ierr;
24647c6ae99SBarry Smith 
24747c6ae99SBarry Smith   PetscFunctionBegin;
24847c6ae99SBarry Smith   if (global) {
24947c6ae99SBarry Smith     ierr = VecResetArray(*global);CHKERRQ(ierr);
25047c6ae99SBarry Smith     ierr = DMRestoreGlobalVector(mine->dm,global);CHKERRQ(ierr);
25147c6ae99SBarry Smith   }
25247c6ae99SBarry Smith   PetscFunctionReturn(0);
25347c6ae99SBarry Smith }
25447c6ae99SBarry Smith 
25547c6ae99SBarry Smith #undef __FUNCT__
25647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter_Array"
25747c6ae99SBarry Smith PetscErrorCode DMCompositeScatter_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar *array)
25847c6ae99SBarry Smith {
25947c6ae99SBarry Smith   PetscErrorCode ierr;
26047c6ae99SBarry Smith   PetscScalar    *varray;
26147c6ae99SBarry Smith   PetscMPIInt    rank;
26247c6ae99SBarry Smith 
26347c6ae99SBarry Smith   PetscFunctionBegin;
26447c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
26547c6ae99SBarry Smith   if (rank == mine->rank) {
26647c6ae99SBarry Smith     ierr    = VecGetArray(vec,&varray);CHKERRQ(ierr);
26747c6ae99SBarry Smith     ierr    = PetscMemcpy(array,varray+mine->rstart,mine->n*sizeof(PetscScalar));CHKERRQ(ierr);
26847c6ae99SBarry Smith     ierr    = VecRestoreArray(vec,&varray);CHKERRQ(ierr);
26947c6ae99SBarry Smith   }
27006ebdd98SJed Brown   ierr    = MPI_Bcast(array,mine->nlocal,MPIU_SCALAR,mine->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
27147c6ae99SBarry Smith   PetscFunctionReturn(0);
27247c6ae99SBarry Smith }
27347c6ae99SBarry Smith 
27447c6ae99SBarry Smith #undef __FUNCT__
27547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter_DM"
27647c6ae99SBarry Smith PetscErrorCode DMCompositeScatter_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec local)
27747c6ae99SBarry Smith {
27847c6ae99SBarry Smith   PetscErrorCode ierr;
27947c6ae99SBarry Smith   PetscScalar    *array;
28047c6ae99SBarry Smith   Vec            global;
28147c6ae99SBarry Smith 
28247c6ae99SBarry Smith   PetscFunctionBegin;
28347c6ae99SBarry Smith   ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr);
28447c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
28547c6ae99SBarry Smith   ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr);
28647c6ae99SBarry Smith   ierr = DMGlobalToLocalBegin(mine->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
28747c6ae99SBarry Smith   ierr = DMGlobalToLocalEnd(mine->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
28847c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
28947c6ae99SBarry Smith   ierr = VecResetArray(global);CHKERRQ(ierr);
29047c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr);
29147c6ae99SBarry Smith   PetscFunctionReturn(0);
29247c6ae99SBarry Smith }
29347c6ae99SBarry Smith 
29447c6ae99SBarry Smith #undef __FUNCT__
29547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_Array"
296acc1e270SJed Brown PetscErrorCode DMCompositeGather_Array(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,const PetscScalar *array)
29747c6ae99SBarry Smith {
29847c6ae99SBarry Smith   PetscErrorCode ierr;
29947c6ae99SBarry Smith   PetscScalar    *varray;
30047c6ae99SBarry Smith   PetscMPIInt    rank;
30147c6ae99SBarry Smith 
30247c6ae99SBarry Smith   PetscFunctionBegin;
30347c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
30447c6ae99SBarry Smith   if (rank == mine->rank) {
30547c6ae99SBarry Smith     ierr = VecGetArray(vec,&varray);CHKERRQ(ierr);
30647c6ae99SBarry Smith     if (varray+mine->rstart == array) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"You need not DMCompositeGather() into objects obtained via DMCompositeGetAccess()");
307acc1e270SJed Brown   }
308df0c820aSJed Brown   switch (imode) {
309df0c820aSJed Brown   case INSERT_VALUES:
310acc1e270SJed Brown     if (rank == mine->rank) {
31147c6ae99SBarry Smith       ierr = PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(PetscScalar));CHKERRQ(ierr);
312acc1e270SJed Brown     }
313df0c820aSJed Brown     break;
314df0c820aSJed Brown   case ADD_VALUES: {
315df0c820aSJed Brown     PetscInt          i;
316760fd489SMatthew G Knepley     void             *source;
317acc1e270SJed Brown     PetscScalar       *buffer,*dest;
318acc1e270SJed Brown     if (rank == mine->rank) {
319acc1e270SJed Brown       dest = &varray[mine->rstart];
320acc1e270SJed Brown #if defined(PETSC_HAVE_MPI_IN_PLACE)
321acc1e270SJed Brown       buffer = dest;
322acc1e270SJed Brown       source = MPI_IN_PLACE;
323acc1e270SJed Brown #else
32406ebdd98SJed Brown       ierr = PetscMalloc(mine->nlocal*sizeof(PetscScalar),&buffer);CHKERRQ(ierr);
325760fd489SMatthew G Knepley       source = (void *) buffer;
326acc1e270SJed Brown #endif
32706ebdd98SJed Brown       for (i=0; i<mine->nlocal; i++) buffer[i] = varray[mine->rstart+i] + array[i];
328acc1e270SJed Brown     } else {
329760fd489SMatthew G Knepley       source = (void *) array;
330acc1e270SJed Brown       dest   = PETSC_NULL;
331acc1e270SJed Brown     }
33206ebdd98SJed Brown     ierr = MPI_Reduce(source,dest,mine->nlocal,MPIU_SCALAR,MPI_SUM,mine->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
333acc1e270SJed Brown #if !defined(PETSC_HAVE_MPI_IN_PLACE)
334acc1e270SJed Brown     if (rank == mine->rank) {ierr = PetscFree(source);CHKERRQ(ierr);}
335acc1e270SJed Brown #endif
336df0c820aSJed Brown   } break;
337df0c820aSJed Brown   default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"imode");
338df0c820aSJed Brown   }
339acc1e270SJed Brown   if (rank == mine->rank) {ierr = VecRestoreArray(vec,&varray);CHKERRQ(ierr);}
34047c6ae99SBarry Smith   PetscFunctionReturn(0);
34147c6ae99SBarry Smith }
34247c6ae99SBarry Smith 
34347c6ae99SBarry Smith #undef __FUNCT__
34447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_DM"
345df0c820aSJed Brown PetscErrorCode DMCompositeGather_DM(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,Vec local)
34647c6ae99SBarry Smith {
34747c6ae99SBarry Smith   PetscErrorCode ierr;
34847c6ae99SBarry Smith   PetscScalar    *array;
34947c6ae99SBarry Smith   Vec            global;
35047c6ae99SBarry Smith 
35147c6ae99SBarry Smith   PetscFunctionBegin;
35247c6ae99SBarry Smith   ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr);
35347c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
35447c6ae99SBarry Smith   ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr);
355df0c820aSJed Brown   ierr = DMLocalToGlobalBegin(mine->dm,local,imode,global);CHKERRQ(ierr);
356df0c820aSJed Brown   ierr = DMLocalToGlobalEnd(mine->dm,local,imode,global);CHKERRQ(ierr);
35747c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
35847c6ae99SBarry Smith   ierr = VecResetArray(global);CHKERRQ(ierr);
35947c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr);
36047c6ae99SBarry Smith   PetscFunctionReturn(0);
36147c6ae99SBarry Smith }
36247c6ae99SBarry Smith 
36347c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
36447c6ae99SBarry Smith 
36547c6ae99SBarry Smith #include <stdarg.h>
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith #undef __FUNCT__
36847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetNumberDM"
36947c6ae99SBarry Smith /*@C
37047c6ae99SBarry Smith     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
37147c6ae99SBarry Smith        representation.
37247c6ae99SBarry Smith 
37347c6ae99SBarry Smith     Not Collective
37447c6ae99SBarry Smith 
37547c6ae99SBarry Smith     Input Parameter:
37647c6ae99SBarry Smith .    dm - the packer object
37747c6ae99SBarry Smith 
37847c6ae99SBarry Smith     Output Parameter:
37947c6ae99SBarry Smith .     nDM - the number of DMs
38047c6ae99SBarry Smith 
38147c6ae99SBarry Smith     Level: beginner
38247c6ae99SBarry Smith 
38347c6ae99SBarry Smith @*/
384*7087cfbeSBarry Smith PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
38547c6ae99SBarry Smith {
38647c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
38747c6ae99SBarry Smith   PetscFunctionBegin;
38847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
38947c6ae99SBarry Smith   *nDM = com->nDM;
39047c6ae99SBarry Smith   PetscFunctionReturn(0);
39147c6ae99SBarry Smith }
39247c6ae99SBarry Smith 
39347c6ae99SBarry Smith 
39447c6ae99SBarry Smith #undef __FUNCT__
39547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess"
39647c6ae99SBarry Smith /*@C
39747c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
39847c6ae99SBarry Smith        representation.
39947c6ae99SBarry Smith 
40047c6ae99SBarry Smith     Collective on DMComposite
40147c6ae99SBarry Smith 
40247c6ae99SBarry Smith     Input Parameter:
40347c6ae99SBarry Smith +    dm - the packer object
40447c6ae99SBarry Smith .    gvec - the global vector
40547c6ae99SBarry Smith -    ... - the individual sequential or parallel objects (arrays or vectors)
40647c6ae99SBarry Smith 
40747c6ae99SBarry Smith     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
40847c6ae99SBarry Smith 
40947c6ae99SBarry Smith     Level: advanced
41047c6ae99SBarry Smith 
41147c6ae99SBarry Smith @*/
412*7087cfbeSBarry Smith PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
41347c6ae99SBarry Smith {
41447c6ae99SBarry Smith   va_list                Argp;
41547c6ae99SBarry Smith   PetscErrorCode         ierr;
41647c6ae99SBarry Smith   struct DMCompositeLink *next;
41747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
41847c6ae99SBarry Smith 
41947c6ae99SBarry Smith   PetscFunctionBegin;
42047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
42147c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
42247c6ae99SBarry Smith   next = com->next;
42347c6ae99SBarry Smith   if (!com->setup) {
424d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
42547c6ae99SBarry Smith   }
42647c6ae99SBarry Smith 
42747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
42847c6ae99SBarry Smith   va_start(Argp,gvec);
42947c6ae99SBarry Smith   while (next) {
43047c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
43147c6ae99SBarry Smith       PetscScalar **array;
43247c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
43347c6ae99SBarry Smith       ierr  = DMCompositeGetAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
43447c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
43547c6ae99SBarry Smith       Vec *vec;
43647c6ae99SBarry Smith       vec  = va_arg(Argp, Vec*);
43747c6ae99SBarry Smith       ierr = DMCompositeGetAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
43847c6ae99SBarry Smith     } else {
43947c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
44047c6ae99SBarry Smith     }
44147c6ae99SBarry Smith     next = next->next;
44247c6ae99SBarry Smith   }
44347c6ae99SBarry Smith   va_end(Argp);
44447c6ae99SBarry Smith   PetscFunctionReturn(0);
44547c6ae99SBarry Smith }
44647c6ae99SBarry Smith 
44747c6ae99SBarry Smith #undef __FUNCT__
44847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess"
44947c6ae99SBarry Smith /*@C
450aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
45147c6ae99SBarry Smith        representation.
45247c6ae99SBarry Smith 
45347c6ae99SBarry Smith     Collective on DMComposite
45447c6ae99SBarry Smith 
45547c6ae99SBarry Smith     Input Parameter:
45647c6ae99SBarry Smith +    dm - the packer object
45747c6ae99SBarry Smith .    gvec - the global vector
45847c6ae99SBarry Smith -    ... - the individual sequential or parallel objects (arrays or vectors)
45947c6ae99SBarry Smith 
46047c6ae99SBarry Smith     Level: advanced
46147c6ae99SBarry Smith 
4620c010503SBarry Smith .seealso  DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
4636eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
464aa219208SBarry Smith          DMCompositeRestoreAccess(), DMCompositeGetAccess()
46547c6ae99SBarry Smith 
46647c6ae99SBarry Smith @*/
467*7087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
46847c6ae99SBarry Smith {
46947c6ae99SBarry Smith   va_list                Argp;
47047c6ae99SBarry Smith   PetscErrorCode         ierr;
47147c6ae99SBarry Smith   struct DMCompositeLink *next;
47247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
47347c6ae99SBarry Smith 
47447c6ae99SBarry Smith   PetscFunctionBegin;
47547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
47647c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
47747c6ae99SBarry Smith   next = com->next;
47847c6ae99SBarry Smith   if (!com->setup) {
479d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
48047c6ae99SBarry Smith   }
48147c6ae99SBarry Smith 
48247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
48347c6ae99SBarry Smith   va_start(Argp,gvec);
48447c6ae99SBarry Smith   while (next) {
48547c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
48647c6ae99SBarry Smith       PetscScalar **array;
48747c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
48847c6ae99SBarry Smith       ierr  = DMCompositeRestoreAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
48947c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
49047c6ae99SBarry Smith       Vec *vec;
49147c6ae99SBarry Smith       vec  = va_arg(Argp, Vec*);
49247c6ae99SBarry Smith       ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
49347c6ae99SBarry Smith     } else {
49447c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
49547c6ae99SBarry Smith     }
49647c6ae99SBarry Smith     next = next->next;
49747c6ae99SBarry Smith   }
49847c6ae99SBarry Smith   va_end(Argp);
49947c6ae99SBarry Smith   PetscFunctionReturn(0);
50047c6ae99SBarry Smith }
50147c6ae99SBarry Smith 
50247c6ae99SBarry Smith #undef __FUNCT__
50347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter"
50447c6ae99SBarry Smith /*@C
50547c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
50647c6ae99SBarry Smith 
50747c6ae99SBarry Smith     Collective on DMComposite
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith     Input Parameter:
51047c6ae99SBarry Smith +    dm - the packer object
51147c6ae99SBarry Smith .    gvec - the global vector
5128fd8f222SJed Brown -    ... - the individual sequential objects (arrays or vectors), PETSC_NULL for those that are not needed
51347c6ae99SBarry Smith 
51447c6ae99SBarry Smith     Level: advanced
51547c6ae99SBarry Smith 
5160c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
5176eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
51847c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith @*/
521*7087cfbeSBarry Smith PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
52247c6ae99SBarry Smith {
52347c6ae99SBarry Smith   va_list                Argp;
52447c6ae99SBarry Smith   PetscErrorCode         ierr;
52547c6ae99SBarry Smith   struct DMCompositeLink *next;
5268fd8f222SJed Brown   PetscInt               cnt;
52747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
52847c6ae99SBarry Smith 
52947c6ae99SBarry Smith   PetscFunctionBegin;
53047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
53147c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
53247c6ae99SBarry Smith   if (!com->setup) {
533d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
53447c6ae99SBarry Smith   }
53547c6ae99SBarry Smith 
53647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
53747c6ae99SBarry Smith   va_start(Argp,gvec);
5388fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
53947c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
54047c6ae99SBarry Smith       PetscScalar *array;
54147c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar*);
542c72eef85SJed Brown       if (array) PetscValidScalarPointer(array,cnt);
543c72eef85SJed Brown       PetscValidLogicalCollectiveBool(dm,!!array,cnt);
5448fd8f222SJed Brown       if (!array) continue;
54547c6ae99SBarry Smith       ierr = DMCompositeScatter_Array(dm,next,gvec,array);CHKERRQ(ierr);
54647c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
54747c6ae99SBarry Smith       Vec vec;
54847c6ae99SBarry Smith       vec = va_arg(Argp, Vec);
5498fd8f222SJed Brown       if (!vec) continue;
55047c6ae99SBarry Smith       PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt);
55147c6ae99SBarry Smith       ierr = DMCompositeScatter_DM(dm,next,gvec,vec);CHKERRQ(ierr);
55247c6ae99SBarry Smith     } else {
55347c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
55447c6ae99SBarry Smith     }
55547c6ae99SBarry Smith   }
55647c6ae99SBarry Smith   va_end(Argp);
55747c6ae99SBarry Smith   PetscFunctionReturn(0);
55847c6ae99SBarry Smith }
55947c6ae99SBarry Smith 
56047c6ae99SBarry Smith #undef __FUNCT__
56147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather"
56247c6ae99SBarry Smith /*@C
56347c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
56447c6ae99SBarry Smith 
56547c6ae99SBarry Smith     Collective on DMComposite
56647c6ae99SBarry Smith 
56747c6ae99SBarry Smith     Input Parameter:
56847c6ae99SBarry Smith +    dm - the packer object
56947c6ae99SBarry Smith .    gvec - the global vector
5708fd8f222SJed Brown -    ... - the individual sequential objects (arrays or vectors), PETSC_NULL for any that are not needed
57147c6ae99SBarry Smith 
57247c6ae99SBarry Smith     Level: advanced
57347c6ae99SBarry Smith 
5740c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
5756eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
57647c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
57747c6ae99SBarry Smith 
57847c6ae99SBarry Smith @*/
579*7087cfbeSBarry Smith PetscErrorCode  DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
58047c6ae99SBarry Smith {
58147c6ae99SBarry Smith   va_list                Argp;
58247c6ae99SBarry Smith   PetscErrorCode         ierr;
58347c6ae99SBarry Smith   struct DMCompositeLink *next;
58447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
5858fd8f222SJed Brown   PetscInt               cnt;
58647c6ae99SBarry Smith 
58747c6ae99SBarry Smith   PetscFunctionBegin;
58847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
58947c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
59047c6ae99SBarry Smith   if (!com->setup) {
591d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
59247c6ae99SBarry Smith   }
59347c6ae99SBarry Smith 
59447c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
595df0c820aSJed Brown   va_start(Argp,imode);
5968fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
59747c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
59847c6ae99SBarry Smith       PetscScalar *array;
59947c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar*);
6008fd8f222SJed Brown       if (!array) continue;
6018fd8f222SJed Brown       PetscValidScalarPointer(array,cnt);
602df0c820aSJed Brown       ierr  = DMCompositeGather_Array(dm,next,gvec,imode,array);CHKERRQ(ierr);
60347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
60447c6ae99SBarry Smith       Vec vec;
60547c6ae99SBarry Smith       vec = va_arg(Argp, Vec);
6068fd8f222SJed Brown       if (!vec) continue;
6078fd8f222SJed Brown       PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt);
608df0c820aSJed Brown       ierr = DMCompositeGather_DM(dm,next,gvec,imode,vec);CHKERRQ(ierr);
60947c6ae99SBarry Smith     } else {
61047c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
61147c6ae99SBarry Smith     }
61247c6ae99SBarry Smith   }
61347c6ae99SBarry Smith   va_end(Argp);
61447c6ae99SBarry Smith   PetscFunctionReturn(0);
61547c6ae99SBarry Smith }
61647c6ae99SBarry Smith 
61747c6ae99SBarry Smith #undef __FUNCT__
61847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddArray"
61947c6ae99SBarry Smith /*@C
62047c6ae99SBarry Smith     DMCompositeAddArray - adds an "redundant" array to a DMComposite. The array values will
62147c6ae99SBarry Smith        be stored in part of the array on process orank.
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith     Collective on DMComposite
62447c6ae99SBarry Smith 
62547c6ae99SBarry Smith     Input Parameter:
62647c6ae99SBarry Smith +    dm - the packer object
62747c6ae99SBarry Smith .    orank - the process on which the array entries officially live, this number must be
62847c6ae99SBarry Smith              the same on all processes.
62947c6ae99SBarry Smith -    n - the length of the array
63047c6ae99SBarry Smith 
63147c6ae99SBarry Smith     Level: advanced
63247c6ae99SBarry Smith 
6330c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
6346eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
63547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
63647c6ae99SBarry Smith 
63747c6ae99SBarry Smith @*/
638*7087cfbeSBarry Smith PetscErrorCode  DMCompositeAddArray(DM dm,PetscMPIInt orank,PetscInt n)
63947c6ae99SBarry Smith {
64047c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
64147c6ae99SBarry Smith   PetscErrorCode         ierr;
64247c6ae99SBarry Smith   PetscMPIInt            rank;
64347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
64447c6ae99SBarry Smith 
64547c6ae99SBarry Smith   PetscFunctionBegin;
64647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
64747c6ae99SBarry Smith   next = com->next;
64847c6ae99SBarry Smith   if (com->setup) {
64947c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add an array once you have used the DMComposite");
65047c6ae99SBarry Smith   }
65147c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
65247c6ae99SBarry Smith   {
65347c6ae99SBarry Smith     PetscMPIInt orankmax;
65447c6ae99SBarry Smith     ierr = MPI_Allreduce(&orank,&orankmax,1,MPI_INT,MPI_MAX,((PetscObject)dm)->comm);CHKERRQ(ierr);
65547c6ae99SBarry 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);
65647c6ae99SBarry Smith   }
65747c6ae99SBarry Smith #endif
65847c6ae99SBarry Smith 
65947c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
66047c6ae99SBarry Smith   /* create new link */
66147c6ae99SBarry Smith   ierr                = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
66206ebdd98SJed Brown   mine->nlocal        = n;
66306ebdd98SJed Brown   mine->n             = (rank == orank) ? n : 0;
66447c6ae99SBarry Smith   mine->rank          = orank;
66547c6ae99SBarry Smith   mine->dm            = PETSC_NULL;
66647c6ae99SBarry Smith   mine->type          = DMCOMPOSITE_ARRAY;
66747c6ae99SBarry Smith   mine->next          = PETSC_NULL;
66847c6ae99SBarry Smith   if (rank == mine->rank) {com->n += n;com->nmine++;}
66947c6ae99SBarry Smith 
67047c6ae99SBarry Smith   /* add to end of list */
67147c6ae99SBarry Smith   if (!next) {
67247c6ae99SBarry Smith     com->next = mine;
67347c6ae99SBarry Smith   } else {
67447c6ae99SBarry Smith     while (next->next) next = next->next;
67547c6ae99SBarry Smith     next->next = mine;
67647c6ae99SBarry Smith   }
67747c6ae99SBarry Smith   com->nredundant++;
67847c6ae99SBarry Smith   PetscFunctionReturn(0);
67947c6ae99SBarry Smith }
68047c6ae99SBarry Smith 
68147c6ae99SBarry Smith #undef __FUNCT__
68247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddDM"
68347c6ae99SBarry Smith /*@C
684aa219208SBarry Smith     DMCompositeAddDM - adds a DM  vector to a DMComposite
68547c6ae99SBarry Smith 
68647c6ae99SBarry Smith     Collective on DMComposite
68747c6ae99SBarry Smith 
68847c6ae99SBarry Smith     Input Parameter:
68947c6ae99SBarry Smith +    dm - the packer object
69047c6ae99SBarry Smith -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
69147c6ae99SBarry Smith 
69247c6ae99SBarry Smith     Level: advanced
69347c6ae99SBarry Smith 
6940c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
6956eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
69647c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
69747c6ae99SBarry Smith 
69847c6ae99SBarry Smith @*/
699*7087cfbeSBarry Smith PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
70047c6ae99SBarry Smith {
70147c6ae99SBarry Smith   PetscErrorCode         ierr;
70206ebdd98SJed Brown   PetscInt               n,nlocal;
70347c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
70406ebdd98SJed Brown   Vec                    global,local;
70547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmc->data;
70647c6ae99SBarry Smith 
70747c6ae99SBarry Smith   PetscFunctionBegin;
70847c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
70947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
71047c6ae99SBarry Smith   next = com->next;
711aa219208SBarry Smith   if (com->setup) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
71247c6ae99SBarry Smith 
71347c6ae99SBarry Smith   /* create new link */
71447c6ae99SBarry Smith   ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
71547c6ae99SBarry Smith   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
71647c6ae99SBarry Smith   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
71747c6ae99SBarry Smith   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
71847c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
71906ebdd98SJed Brown   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
72006ebdd98SJed Brown   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
72106ebdd98SJed Brown   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
72247c6ae99SBarry Smith   mine->n      = n;
72306ebdd98SJed Brown   mine->nlocal = nlocal;
72447c6ae99SBarry Smith   mine->dm     = dm;
72547c6ae99SBarry Smith   mine->type   = DMCOMPOSITE_DM;
72647c6ae99SBarry Smith   mine->next   = PETSC_NULL;
72747c6ae99SBarry Smith   com->n       += n;
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith   /* add to end of list */
73047c6ae99SBarry Smith   if (!next) {
73147c6ae99SBarry Smith     com->next = mine;
73247c6ae99SBarry Smith   } else {
73347c6ae99SBarry Smith     while (next->next) next = next->next;
73447c6ae99SBarry Smith     next->next = mine;
73547c6ae99SBarry Smith   }
73647c6ae99SBarry Smith   com->nDM++;
73747c6ae99SBarry Smith   com->nmine++;
73847c6ae99SBarry Smith   PetscFunctionReturn(0);
73947c6ae99SBarry Smith }
74047c6ae99SBarry Smith 
741*7087cfbeSBarry Smith extern PetscErrorCode  VecView_MPI(Vec,PetscViewer);
74247c6ae99SBarry Smith EXTERN_C_BEGIN
74347c6ae99SBarry Smith #undef __FUNCT__
74447c6ae99SBarry Smith #define __FUNCT__ "VecView_DMComposite"
745*7087cfbeSBarry Smith PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
74647c6ae99SBarry Smith {
74747c6ae99SBarry Smith   DM                     dm;
74847c6ae99SBarry Smith   PetscErrorCode         ierr;
74947c6ae99SBarry Smith   struct DMCompositeLink *next;
75047c6ae99SBarry Smith   PetscBool              isdraw;
751cef07954SSatish Balay   DM_Composite           *com;
75247c6ae99SBarry Smith 
75347c6ae99SBarry Smith   PetscFunctionBegin;
75447c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)gvec,"DMComposite",(PetscObject*)&dm);CHKERRQ(ierr);
75547c6ae99SBarry Smith   if (!dm) SETERRQ(((PetscObject)gvec)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
75647c6ae99SBarry Smith   com = (DM_Composite*)dm->data;
75747c6ae99SBarry Smith   next = com->next;
75847c6ae99SBarry Smith 
75947c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
76047c6ae99SBarry Smith   if (!isdraw) {
76147c6ae99SBarry Smith     /* do I really want to call this? */
76247c6ae99SBarry Smith     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
76347c6ae99SBarry Smith   } else {
76447c6ae99SBarry Smith     PetscInt cnt = 0;
76547c6ae99SBarry Smith 
76647c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
76747c6ae99SBarry Smith     while (next) {
76847c6ae99SBarry Smith       if (next->type == DMCOMPOSITE_ARRAY) {
76947c6ae99SBarry Smith 	PetscScalar *array;
77047c6ae99SBarry Smith 	ierr  = DMCompositeGetAccess_Array(dm,next,gvec,&array);CHKERRQ(ierr);
77147c6ae99SBarry Smith 
77247c6ae99SBarry Smith 	/*skip it for now */
77347c6ae99SBarry Smith       } else if (next->type == DMCOMPOSITE_DM) {
77447c6ae99SBarry Smith 	Vec      vec;
77547c6ae99SBarry Smith         PetscInt bs;
77647c6ae99SBarry Smith 
77747c6ae99SBarry Smith 	ierr = DMCompositeGetAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
77847c6ae99SBarry Smith 	ierr = VecView(vec,viewer);CHKERRQ(ierr);
77947c6ae99SBarry Smith         ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
78047c6ae99SBarry Smith 	ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
78147c6ae99SBarry Smith         ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
78247c6ae99SBarry Smith         cnt += bs;
78347c6ae99SBarry Smith       } else {
78447c6ae99SBarry Smith 	SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
78547c6ae99SBarry Smith       }
78647c6ae99SBarry Smith       next = next->next;
78747c6ae99SBarry Smith     }
78847c6ae99SBarry Smith     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
78947c6ae99SBarry Smith   }
79047c6ae99SBarry Smith   PetscFunctionReturn(0);
79147c6ae99SBarry Smith }
79247c6ae99SBarry Smith EXTERN_C_END
79347c6ae99SBarry Smith 
79447c6ae99SBarry Smith 
79547c6ae99SBarry Smith #undef __FUNCT__
7960c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Composite"
797*7087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
79847c6ae99SBarry Smith {
79947c6ae99SBarry Smith   PetscErrorCode         ierr;
80047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
80147c6ae99SBarry Smith 
80247c6ae99SBarry Smith   PetscFunctionBegin;
80347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
80447c6ae99SBarry Smith   if (!com->setup) {
805d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
80647c6ae99SBarry Smith   }
80747c6ae99SBarry Smith   ierr = VecCreateMPI(((PetscObject)dm)->comm,com->n,com->N,gvec);CHKERRQ(ierr);
80847c6ae99SBarry Smith   ierr = PetscObjectCompose((PetscObject)*gvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr);
80947c6ae99SBarry Smith   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void(*)(void))VecView_DMComposite);CHKERRQ(ierr);
81047c6ae99SBarry Smith   PetscFunctionReturn(0);
81147c6ae99SBarry Smith }
81247c6ae99SBarry Smith 
81347c6ae99SBarry Smith #undef __FUNCT__
8140c010503SBarry Smith #define __FUNCT__ "DMCreateLocalVector_Composite"
815*7087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
81647c6ae99SBarry Smith {
81747c6ae99SBarry Smith   PetscErrorCode         ierr;
81847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
81947c6ae99SBarry Smith 
82047c6ae99SBarry Smith   PetscFunctionBegin;
82147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
82247c6ae99SBarry Smith   if (!com->setup) {
823d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
82447c6ae99SBarry Smith   }
82547c6ae99SBarry Smith   ierr = VecCreateSeq(((PetscObject)dm)->comm,com->nghost,lvec);CHKERRQ(ierr);
82647c6ae99SBarry Smith   ierr = PetscObjectCompose((PetscObject)*lvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr);
82747c6ae99SBarry Smith   PetscFunctionReturn(0);
82847c6ae99SBarry Smith }
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith #undef __FUNCT__
8316eb61c8cSJed Brown #define __FUNCT__ "DMCompositeGetISLocalToGlobalMappings"
83247c6ae99SBarry Smith /*@C
8336eb61c8cSJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM/array in the DMComposite, maps to the composite global space
83447c6ae99SBarry Smith 
83506ebdd98SJed Brown     Collective on DM
83647c6ae99SBarry Smith 
83747c6ae99SBarry Smith     Input Parameter:
83847c6ae99SBarry Smith .    dm - the packer object
83947c6ae99SBarry Smith 
84047c6ae99SBarry Smith     Output Parameters:
8416eb61c8cSJed Brown .    ltogs - the individual mappings for each packed vector/array. Note that this includes
842aa219208SBarry Smith            all the ghost points that individual ghosted DMDA's may have. Also each process has an
8436eb61c8cSJed Brown            mapping for EACH redundant array (not just the local redundant arrays).
84447c6ae99SBarry Smith 
84547c6ae99SBarry Smith     Level: advanced
84647c6ae99SBarry Smith 
84747c6ae99SBarry Smith     Notes:
8486eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
84947c6ae99SBarry Smith 
8500c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
85147c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
85247c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith @*/
855*7087cfbeSBarry Smith PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
85647c6ae99SBarry Smith {
85747c6ae99SBarry Smith   PetscErrorCode         ierr;
85847c6ae99SBarry Smith   PetscInt               i,*idx,n,cnt;
85947c6ae99SBarry Smith   struct DMCompositeLink *next;
86047c6ae99SBarry Smith   PetscMPIInt            rank;
86147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
86247c6ae99SBarry Smith 
86347c6ae99SBarry Smith   PetscFunctionBegin;
86447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
865acc1e270SJed Brown   ierr = PetscMalloc((com->nDM+com->nredundant)*sizeof(ISLocalToGlobalMapping),ltogs);CHKERRQ(ierr);
86647c6ae99SBarry Smith   next = com->next;
86747c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
86847c6ae99SBarry Smith 
86947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
87047c6ae99SBarry Smith   cnt = 0;
87147c6ae99SBarry Smith   while (next) {
87247c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
87306ebdd98SJed Brown       ierr = PetscMalloc(next->nlocal*sizeof(PetscInt),&idx);CHKERRQ(ierr);
87447c6ae99SBarry Smith       if (rank == next->rank) {
87506ebdd98SJed Brown         for (i=0; i<next->nlocal; i++) idx[i] = next->grstart + i;
87647c6ae99SBarry Smith       }
87706ebdd98SJed Brown       ierr = MPI_Bcast(idx,next->nlocal,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
87806ebdd98SJed Brown       ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,next->nlocal,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
87947c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
8806eb61c8cSJed Brown       ISLocalToGlobalMapping ltog;
8816eb61c8cSJed Brown       PetscMPIInt            size;
88286994e45SJed Brown       const PetscInt         *suboff,*indices;
8836eb61c8cSJed Brown       Vec                    global;
88447c6ae99SBarry Smith 
8856eb61c8cSJed Brown       /* Get sub-DM global indices for each local dof */
8861411c6eeSJed Brown       ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
8876eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
88886994e45SJed Brown       ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
88947c6ae99SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
89047c6ae99SBarry Smith 
8916eb61c8cSJed Brown       /* Get the offsets for the sub-DM global vector */
8926eb61c8cSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
8936eb61c8cSJed Brown       ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
8946eb61c8cSJed Brown       ierr = MPI_Comm_size(((PetscObject)global)->comm,&size);CHKERRQ(ierr);
8956eb61c8cSJed Brown 
8966eb61c8cSJed Brown       /* Shift the sub-DM definition of the global space to the composite global space */
8976eb61c8cSJed Brown       for (i=0; i<n; i++) {
89886994e45SJed Brown         PetscInt subi = indices[i],lo = 0,hi = size,t;
8996eb61c8cSJed Brown         /* Binary search to find which rank owns subi */
9006eb61c8cSJed Brown         while (hi-lo > 1) {
9016eb61c8cSJed Brown           t = lo + (hi-lo)/2;
9026eb61c8cSJed Brown           if (suboff[t] > subi) hi = t;
9036eb61c8cSJed Brown           else                  lo = t;
9046eb61c8cSJed Brown         }
9056eb61c8cSJed Brown         idx[i] = subi - suboff[lo] + next->grstarts[lo];
9066eb61c8cSJed Brown       }
90786994e45SJed Brown       ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
9086eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
9096eb61c8cSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
9106eb61c8cSJed Brown     } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
91147c6ae99SBarry Smith     next = next->next;
91247c6ae99SBarry Smith     cnt++;
91347c6ae99SBarry Smith   }
91447c6ae99SBarry Smith   PetscFunctionReturn(0);
91547c6ae99SBarry Smith }
91647c6ae99SBarry Smith 
91747c6ae99SBarry Smith #undef __FUNCT__
91887c85e80SJed Brown #define __FUNCT__ "DMCompositeGetLocalISs"
91987c85e80SJed Brown /*@C
92087c85e80SJed Brown    DMCompositeGetLocalISs - Gets index sets for each DM/array component of a composite local vector
92187c85e80SJed Brown 
92287c85e80SJed Brown    Not Collective
92387c85e80SJed Brown 
92487c85e80SJed Brown    Input Arguments:
92587c85e80SJed Brown . dm - composite DM
92687c85e80SJed Brown 
92787c85e80SJed Brown    Output Arguments:
92887c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
92987c85e80SJed Brown 
93087c85e80SJed Brown    Level: intermediate
93187c85e80SJed Brown 
93287c85e80SJed Brown    Notes:
93387c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
93487c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
93587c85e80SJed Brown    scatter to a composite local vector.
93687c85e80SJed Brown 
93787c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
93887c85e80SJed Brown 
93987c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
94087c85e80SJed Brown 
94187c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
94287c85e80SJed Brown 
94387c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
94487c85e80SJed Brown @*/
945*7087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
94687c85e80SJed Brown {
94787c85e80SJed Brown   PetscErrorCode         ierr;
94887c85e80SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
94987c85e80SJed Brown   struct DMCompositeLink *link;
95087c85e80SJed Brown   PetscInt cnt,start;
95187c85e80SJed Brown 
95287c85e80SJed Brown   PetscFunctionBegin;
95387c85e80SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
95487c85e80SJed Brown   PetscValidPointer(is,2);
95587c85e80SJed Brown   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
95606ebdd98SJed Brown   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
95706ebdd98SJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
958520db06cSJed Brown     if (link->type == DMCOMPOSITE_DM) {
959520db06cSJed Brown       PetscInt bs;
9601411c6eeSJed Brown       ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
961520db06cSJed Brown       ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
962520db06cSJed Brown     }
96387c85e80SJed Brown   }
96487c85e80SJed Brown   PetscFunctionReturn(0);
96587c85e80SJed Brown }
96687c85e80SJed Brown 
96787c85e80SJed Brown #undef __FUNCT__
96847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetGlobalISs"
96947c6ae99SBarry Smith /*@C
97047c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
97147c6ae99SBarry Smith 
97247c6ae99SBarry Smith     Collective on DMComposite
97347c6ae99SBarry Smith 
97447c6ae99SBarry Smith     Input Parameter:
97547c6ae99SBarry Smith .    dm - the packer object
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith     Output Parameters:
97847c6ae99SBarry Smith .    is - the array of index sets
97947c6ae99SBarry Smith 
98047c6ae99SBarry Smith     Level: advanced
98147c6ae99SBarry Smith 
98247c6ae99SBarry Smith     Notes:
98347c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
98447c6ae99SBarry Smith 
98547c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
98647c6ae99SBarry Smith 
9876eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
9886eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
9896eb61c8cSJed Brown        indices.
99047c6ae99SBarry Smith 
9910c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
99247c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
99347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
99447c6ae99SBarry Smith 
99547c6ae99SBarry Smith @*/
9966eb61c8cSJed Brown 
997*7087cfbeSBarry Smith PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
99847c6ae99SBarry Smith {
99947c6ae99SBarry Smith   PetscErrorCode         ierr;
100047c6ae99SBarry Smith   PetscInt               cnt = 0,*idx,i;
100147c6ae99SBarry Smith   struct DMCompositeLink *next;
100247c6ae99SBarry Smith   PetscMPIInt            rank;
100347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
100447c6ae99SBarry Smith 
100547c6ae99SBarry Smith   PetscFunctionBegin;
100647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
100706ebdd98SJed Brown   ierr = PetscMalloc((com->nDM+com->nredundant)*sizeof(IS),is);CHKERRQ(ierr);
100847c6ae99SBarry Smith   next = com->next;
100947c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
101047c6ae99SBarry Smith 
101147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
101247c6ae99SBarry Smith   while (next) {
101347c6ae99SBarry Smith     ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
101447c6ae99SBarry Smith     for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
101547c6ae99SBarry Smith     ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
101647c6ae99SBarry Smith     cnt++;
101747c6ae99SBarry Smith     next = next->next;
101847c6ae99SBarry Smith   }
101947c6ae99SBarry Smith   PetscFunctionReturn(0);
102047c6ae99SBarry Smith }
102147c6ae99SBarry Smith 
102247c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
102347c6ae99SBarry Smith #undef __FUNCT__
102447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_Array"
102547c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
102647c6ae99SBarry Smith {
102747c6ae99SBarry Smith   PetscErrorCode ierr;
102847c6ae99SBarry Smith   PetscFunctionBegin;
102947c6ae99SBarry Smith   if (array) {
1030c72eef85SJed Brown     ierr = PetscMalloc(mine->nlocal*sizeof(PetscScalar),array);CHKERRQ(ierr);
103147c6ae99SBarry Smith   }
103247c6ae99SBarry Smith   PetscFunctionReturn(0);
103347c6ae99SBarry Smith }
103447c6ae99SBarry Smith 
103547c6ae99SBarry Smith #undef __FUNCT__
103647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_DM"
103747c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
103847c6ae99SBarry Smith {
103947c6ae99SBarry Smith   PetscErrorCode ierr;
104047c6ae99SBarry Smith   PetscFunctionBegin;
104147c6ae99SBarry Smith   if (local) {
104247c6ae99SBarry Smith     ierr = DMGetLocalVector(mine->dm,local);CHKERRQ(ierr);
104347c6ae99SBarry Smith   }
104447c6ae99SBarry Smith   PetscFunctionReturn(0);
104547c6ae99SBarry Smith }
104647c6ae99SBarry Smith 
104747c6ae99SBarry Smith #undef __FUNCT__
104847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_Array"
104947c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
105047c6ae99SBarry Smith {
105147c6ae99SBarry Smith   PetscErrorCode ierr;
105247c6ae99SBarry Smith   PetscFunctionBegin;
105347c6ae99SBarry Smith   if (array) {
105447c6ae99SBarry Smith     ierr = PetscFree(*array);CHKERRQ(ierr);
105547c6ae99SBarry Smith   }
105647c6ae99SBarry Smith   PetscFunctionReturn(0);
105747c6ae99SBarry Smith }
105847c6ae99SBarry Smith 
105947c6ae99SBarry Smith #undef __FUNCT__
106047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_DM"
106147c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
106247c6ae99SBarry Smith {
106347c6ae99SBarry Smith   PetscErrorCode ierr;
106447c6ae99SBarry Smith   PetscFunctionBegin;
106547c6ae99SBarry Smith   if (local) {
106647c6ae99SBarry Smith     ierr = DMRestoreLocalVector(mine->dm,local);CHKERRQ(ierr);
106747c6ae99SBarry Smith   }
106847c6ae99SBarry Smith   PetscFunctionReturn(0);
106947c6ae99SBarry Smith }
107047c6ae99SBarry Smith 
107147c6ae99SBarry Smith #undef __FUNCT__
107247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors"
107347c6ae99SBarry Smith /*@C
107447c6ae99SBarry Smith     DMCompositeGetLocalVectors - Gets local vectors and arrays for each part of a DMComposite.
107547c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
107647c6ae99SBarry Smith 
107747c6ae99SBarry Smith     Not Collective
107847c6ae99SBarry Smith 
107947c6ae99SBarry Smith     Input Parameter:
108047c6ae99SBarry Smith .    dm - the packer object
108147c6ae99SBarry Smith 
108247c6ae99SBarry Smith     Output Parameter:
108347c6ae99SBarry Smith .    ... - the individual sequential objects (arrays or vectors)
108447c6ae99SBarry Smith 
108547c6ae99SBarry Smith     Level: advanced
108647c6ae99SBarry Smith 
10870c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
10886eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
108947c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
109047c6ae99SBarry Smith 
109147c6ae99SBarry Smith @*/
1092*7087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
109347c6ae99SBarry Smith {
109447c6ae99SBarry Smith   va_list                Argp;
109547c6ae99SBarry Smith   PetscErrorCode         ierr;
109647c6ae99SBarry Smith   struct DMCompositeLink *next;
109747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
109847c6ae99SBarry Smith 
109947c6ae99SBarry Smith   PetscFunctionBegin;
110047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
110147c6ae99SBarry Smith   next = com->next;
110247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
110347c6ae99SBarry Smith   va_start(Argp,dm);
110447c6ae99SBarry Smith   while (next) {
110547c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
110647c6ae99SBarry Smith       PetscScalar **array;
110747c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
110847c6ae99SBarry Smith       ierr = DMCompositeGetLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
110947c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
111047c6ae99SBarry Smith       Vec *vec;
111147c6ae99SBarry Smith       vec = va_arg(Argp, Vec*);
111247c6ae99SBarry Smith       ierr = DMCompositeGetLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
111347c6ae99SBarry Smith     } else {
111447c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
111547c6ae99SBarry Smith     }
111647c6ae99SBarry Smith     next = next->next;
111747c6ae99SBarry Smith   }
111847c6ae99SBarry Smith   va_end(Argp);
111947c6ae99SBarry Smith   PetscFunctionReturn(0);
112047c6ae99SBarry Smith }
112147c6ae99SBarry Smith 
112247c6ae99SBarry Smith #undef __FUNCT__
112347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors"
112447c6ae99SBarry Smith /*@C
112547c6ae99SBarry Smith     DMCompositeRestoreLocalVectors - Restores local vectors and arrays for each part of a DMComposite.
112647c6ae99SBarry Smith 
112747c6ae99SBarry Smith     Not Collective
112847c6ae99SBarry Smith 
112947c6ae99SBarry Smith     Input Parameter:
113047c6ae99SBarry Smith .    dm - the packer object
113147c6ae99SBarry Smith 
113247c6ae99SBarry Smith     Output Parameter:
113347c6ae99SBarry Smith .    ... - the individual sequential objects (arrays or vectors)
113447c6ae99SBarry Smith 
113547c6ae99SBarry Smith     Level: advanced
113647c6ae99SBarry Smith 
11370c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
11386eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
113947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
114047c6ae99SBarry Smith 
114147c6ae99SBarry Smith @*/
1142*7087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
114347c6ae99SBarry Smith {
114447c6ae99SBarry Smith   va_list                Argp;
114547c6ae99SBarry Smith   PetscErrorCode         ierr;
114647c6ae99SBarry Smith   struct DMCompositeLink *next;
114747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
114847c6ae99SBarry Smith 
114947c6ae99SBarry Smith   PetscFunctionBegin;
115047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
115147c6ae99SBarry Smith   next = com->next;
115247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
115347c6ae99SBarry Smith   va_start(Argp,dm);
115447c6ae99SBarry Smith   while (next) {
115547c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
115647c6ae99SBarry Smith       PetscScalar **array;
115747c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
115847c6ae99SBarry Smith       ierr = DMCompositeRestoreLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
115947c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
116047c6ae99SBarry Smith       Vec *vec;
116147c6ae99SBarry Smith       vec = va_arg(Argp, Vec*);
116247c6ae99SBarry Smith       ierr = DMCompositeRestoreLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
116347c6ae99SBarry Smith     } else {
116447c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
116547c6ae99SBarry Smith     }
116647c6ae99SBarry Smith     next = next->next;
116747c6ae99SBarry Smith   }
116847c6ae99SBarry Smith   va_end(Argp);
116947c6ae99SBarry Smith   PetscFunctionReturn(0);
117047c6ae99SBarry Smith }
117147c6ae99SBarry Smith 
117247c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
117347c6ae99SBarry Smith #undef __FUNCT__
117447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_Array"
117547c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_Array(DM dm,struct DMCompositeLink *mine,PetscInt *n)
117647c6ae99SBarry Smith {
117747c6ae99SBarry Smith   PetscFunctionBegin;
117806ebdd98SJed Brown   if (n) *n = mine->nlocal;
117947c6ae99SBarry Smith   PetscFunctionReturn(0);
118047c6ae99SBarry Smith }
118147c6ae99SBarry Smith 
118247c6ae99SBarry Smith #undef __FUNCT__
118347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_DM"
118447c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_DM(DM dmi,struct DMCompositeLink *mine,DM *dm)
118547c6ae99SBarry Smith {
118647c6ae99SBarry Smith   PetscFunctionBegin;
118747c6ae99SBarry Smith   if (dm) *dm = mine->dm;
118847c6ae99SBarry Smith   PetscFunctionReturn(0);
118947c6ae99SBarry Smith }
119047c6ae99SBarry Smith 
119147c6ae99SBarry Smith #undef __FUNCT__
119247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries"
119347c6ae99SBarry Smith /*@C
1194aa219208SBarry Smith     DMCompositeGetEntries - Gets the DM, redundant size, etc for each entry in a DMComposite.
119547c6ae99SBarry Smith 
119647c6ae99SBarry Smith     Not Collective
119747c6ae99SBarry Smith 
119847c6ae99SBarry Smith     Input Parameter:
119947c6ae99SBarry Smith .    dm - the packer object
120047c6ae99SBarry Smith 
120147c6ae99SBarry Smith     Output Parameter:
1202aa219208SBarry Smith .    ... - the individual entries, DMs or integer sizes)
120347c6ae99SBarry Smith 
120447c6ae99SBarry Smith     Level: advanced
120547c6ae99SBarry Smith 
12060c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
12076eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
120847c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
120947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
121047c6ae99SBarry Smith 
121147c6ae99SBarry Smith @*/
1212*7087cfbeSBarry Smith PetscErrorCode  DMCompositeGetEntries(DM dm,...)
121347c6ae99SBarry Smith {
121447c6ae99SBarry Smith   va_list                Argp;
121547c6ae99SBarry Smith   PetscErrorCode         ierr;
121647c6ae99SBarry Smith   struct DMCompositeLink *next;
121747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
121847c6ae99SBarry Smith 
121947c6ae99SBarry Smith   PetscFunctionBegin;
122047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
122147c6ae99SBarry Smith   next = com->next;
122247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
122347c6ae99SBarry Smith   va_start(Argp,dm);
122447c6ae99SBarry Smith   while (next) {
122547c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
122647c6ae99SBarry Smith       PetscInt *n;
122747c6ae99SBarry Smith       n = va_arg(Argp, PetscInt*);
122847c6ae99SBarry Smith       ierr = DMCompositeGetEntries_Array(dm,next,n);CHKERRQ(ierr);
122947c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
123047c6ae99SBarry Smith       DM *dmn;
123147c6ae99SBarry Smith       dmn = va_arg(Argp, DM*);
123247c6ae99SBarry Smith       ierr = DMCompositeGetEntries_DM(dm,next,dmn);CHKERRQ(ierr);
123347c6ae99SBarry Smith     } else {
123447c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
123547c6ae99SBarry Smith     }
123647c6ae99SBarry Smith     next = next->next;
123747c6ae99SBarry Smith   }
123847c6ae99SBarry Smith   va_end(Argp);
123947c6ae99SBarry Smith   PetscFunctionReturn(0);
124047c6ae99SBarry Smith }
124147c6ae99SBarry Smith 
124247c6ae99SBarry Smith #undef __FUNCT__
12430c010503SBarry Smith #define __FUNCT__ "DMRefine_Composite"
1244*7087cfbeSBarry Smith PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
124547c6ae99SBarry Smith {
124647c6ae99SBarry Smith   PetscErrorCode         ierr;
124747c6ae99SBarry Smith   struct DMCompositeLink *next;
124847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmi->data;
124947c6ae99SBarry Smith   DM                     dm;
125047c6ae99SBarry Smith 
125147c6ae99SBarry Smith   PetscFunctionBegin;
125247c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
125347c6ae99SBarry Smith   next = com->next;
125447c6ae99SBarry Smith   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
125547c6ae99SBarry Smith 
125647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
125747c6ae99SBarry Smith   while (next) {
125847c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
125906ebdd98SJed Brown       ierr = DMCompositeAddArray(*fine,next->rank,next->nlocal);CHKERRQ(ierr);
126047c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
126147c6ae99SBarry Smith       ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
126247c6ae99SBarry Smith       ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
126347c6ae99SBarry Smith       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
126447c6ae99SBarry Smith     } else {
126547c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
126647c6ae99SBarry Smith     }
126747c6ae99SBarry Smith     next = next->next;
126847c6ae99SBarry Smith   }
126947c6ae99SBarry Smith   PetscFunctionReturn(0);
127047c6ae99SBarry Smith }
127147c6ae99SBarry Smith 
127247c6ae99SBarry Smith #include "petscmat.h"
127347c6ae99SBarry Smith 
127447c6ae99SBarry Smith struct MatPackLink {
127547c6ae99SBarry Smith   Mat                A;
127647c6ae99SBarry Smith   struct MatPackLink *next;
127747c6ae99SBarry Smith };
127847c6ae99SBarry Smith 
127947c6ae99SBarry Smith struct MatPack {
128047c6ae99SBarry Smith   DM                 right,left;
128147c6ae99SBarry Smith   struct MatPackLink *next;
128247c6ae99SBarry Smith };
128347c6ae99SBarry Smith 
128447c6ae99SBarry Smith #undef __FUNCT__
128547c6ae99SBarry Smith #define __FUNCT__ "MatMultBoth_Shell_Pack"
128647c6ae99SBarry Smith PetscErrorCode MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscBool  add)
128747c6ae99SBarry Smith {
128847c6ae99SBarry Smith   struct MatPack         *mpack;
128947c6ae99SBarry Smith   struct DMCompositeLink *xnext,*ynext;
129047c6ae99SBarry Smith   struct MatPackLink     *anext;
129147c6ae99SBarry Smith   PetscScalar            *xarray,*yarray;
129247c6ae99SBarry Smith   PetscErrorCode         ierr;
129347c6ae99SBarry Smith   PetscInt               i;
129447c6ae99SBarry Smith   Vec                    xglobal,yglobal;
129547c6ae99SBarry Smith   PetscMPIInt            rank;
129647c6ae99SBarry Smith   DM_Composite           *comright;
129747c6ae99SBarry Smith   DM_Composite           *comleft;
129847c6ae99SBarry Smith 
129947c6ae99SBarry Smith   PetscFunctionBegin;
130047c6ae99SBarry Smith   ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
130147c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
130247c6ae99SBarry Smith   comright = (DM_Composite*)mpack->right->data;
130347c6ae99SBarry Smith   comleft = (DM_Composite*)mpack->left->data;
130447c6ae99SBarry Smith   xnext = comright->next;
130547c6ae99SBarry Smith   ynext = comleft->next;
130647c6ae99SBarry Smith   anext = mpack->next;
130747c6ae99SBarry Smith 
130847c6ae99SBarry Smith   while (xnext) {
130947c6ae99SBarry Smith     if (xnext->type == DMCOMPOSITE_ARRAY) {
131047c6ae99SBarry Smith       if (rank == xnext->rank) {
131147c6ae99SBarry Smith         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
131247c6ae99SBarry Smith         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
131347c6ae99SBarry Smith         if (add) {
131447c6ae99SBarry Smith           for (i=0; i<xnext->n; i++) {
131547c6ae99SBarry Smith             yarray[ynext->rstart+i] += xarray[xnext->rstart+i];
131647c6ae99SBarry Smith           }
131747c6ae99SBarry Smith         } else {
131847c6ae99SBarry Smith           ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
131947c6ae99SBarry Smith         }
132047c6ae99SBarry Smith         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
132147c6ae99SBarry Smith         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
132247c6ae99SBarry Smith       }
132347c6ae99SBarry Smith     } else if (xnext->type == DMCOMPOSITE_DM) {
132447c6ae99SBarry Smith       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
132547c6ae99SBarry Smith       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
132647c6ae99SBarry Smith       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
132747c6ae99SBarry Smith       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
132847c6ae99SBarry Smith       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
132947c6ae99SBarry Smith       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
133047c6ae99SBarry Smith       if (add) {
133147c6ae99SBarry Smith         ierr  = MatMultAdd(anext->A,xglobal,yglobal,yglobal);CHKERRQ(ierr);
133247c6ae99SBarry Smith       } else {
133347c6ae99SBarry Smith         ierr  = MatMult(anext->A,xglobal,yglobal);CHKERRQ(ierr);
133447c6ae99SBarry Smith       }
133547c6ae99SBarry Smith       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
133647c6ae99SBarry Smith       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
133747c6ae99SBarry Smith       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
133847c6ae99SBarry Smith       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
133947c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
134047c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
134147c6ae99SBarry Smith       anext = anext->next;
134247c6ae99SBarry Smith     } else {
134347c6ae99SBarry Smith       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
134447c6ae99SBarry Smith     }
134547c6ae99SBarry Smith     xnext = xnext->next;
134647c6ae99SBarry Smith     ynext = ynext->next;
134747c6ae99SBarry Smith   }
134847c6ae99SBarry Smith   PetscFunctionReturn(0);
134947c6ae99SBarry Smith }
135047c6ae99SBarry Smith 
135147c6ae99SBarry Smith #undef __FUNCT__
135247c6ae99SBarry Smith #define __FUNCT__ "MatMultAdd_Shell_Pack"
135347c6ae99SBarry Smith PetscErrorCode MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z)
135447c6ae99SBarry Smith {
135547c6ae99SBarry Smith   PetscErrorCode ierr;
135647c6ae99SBarry Smith   PetscFunctionBegin;
135747c6ae99SBarry Smith   if (z != y) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Handles y == z only");
135847c6ae99SBarry Smith   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
135947c6ae99SBarry Smith   PetscFunctionReturn(0);
136047c6ae99SBarry Smith }
136147c6ae99SBarry Smith 
136247c6ae99SBarry Smith #undef __FUNCT__
136347c6ae99SBarry Smith #define __FUNCT__ "MatMult_Shell_Pack"
136447c6ae99SBarry Smith PetscErrorCode MatMult_Shell_Pack(Mat A,Vec x,Vec y)
136547c6ae99SBarry Smith {
136647c6ae99SBarry Smith   PetscErrorCode ierr;
136747c6ae99SBarry Smith   PetscFunctionBegin;
136847c6ae99SBarry Smith   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
136947c6ae99SBarry Smith   PetscFunctionReturn(0);
137047c6ae99SBarry Smith }
137147c6ae99SBarry Smith 
137247c6ae99SBarry Smith #undef __FUNCT__
137347c6ae99SBarry Smith #define __FUNCT__ "MatMultTranspose_Shell_Pack"
137447c6ae99SBarry Smith PetscErrorCode MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y)
137547c6ae99SBarry Smith {
137647c6ae99SBarry Smith   struct MatPack         *mpack;
137747c6ae99SBarry Smith   struct DMCompositeLink *xnext,*ynext;
137847c6ae99SBarry Smith   struct MatPackLink     *anext;
137947c6ae99SBarry Smith   PetscScalar            *xarray,*yarray;
138047c6ae99SBarry Smith   PetscErrorCode         ierr;
138147c6ae99SBarry Smith   Vec                    xglobal,yglobal;
138247c6ae99SBarry Smith   PetscMPIInt            rank;
138347c6ae99SBarry Smith   DM_Composite           *comright;
138447c6ae99SBarry Smith   DM_Composite           *comleft;
138547c6ae99SBarry Smith 
138647c6ae99SBarry Smith   PetscFunctionBegin;
138747c6ae99SBarry Smith   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
138847c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
138947c6ae99SBarry Smith   comright = (DM_Composite*)mpack->right->data;
139047c6ae99SBarry Smith   comleft = (DM_Composite*)mpack->left->data;
139147c6ae99SBarry Smith   ynext = comright->next;
139247c6ae99SBarry Smith   xnext = comleft->next;
139347c6ae99SBarry Smith   anext = mpack->next;
139447c6ae99SBarry Smith 
139547c6ae99SBarry Smith   while (xnext) {
139647c6ae99SBarry Smith     if (xnext->type == DMCOMPOSITE_ARRAY) {
139747c6ae99SBarry Smith       if (rank == ynext->rank) {
139847c6ae99SBarry Smith         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
139947c6ae99SBarry Smith         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
140047c6ae99SBarry Smith         ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
140147c6ae99SBarry Smith         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
140247c6ae99SBarry Smith         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
140347c6ae99SBarry Smith       }
140447c6ae99SBarry Smith     } else if (xnext->type == DMCOMPOSITE_DM) {
140547c6ae99SBarry Smith       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
140647c6ae99SBarry Smith       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
140747c6ae99SBarry Smith       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
140847c6ae99SBarry Smith       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
140947c6ae99SBarry Smith       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
141047c6ae99SBarry Smith       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
141147c6ae99SBarry Smith       ierr  = MatMultTranspose(anext->A,xglobal,yglobal);CHKERRQ(ierr);
141247c6ae99SBarry Smith       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
141347c6ae99SBarry Smith       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
141447c6ae99SBarry Smith       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
141547c6ae99SBarry Smith       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
141647c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
141747c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
141847c6ae99SBarry Smith       anext = anext->next;
141947c6ae99SBarry Smith     } else {
142047c6ae99SBarry Smith       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
142147c6ae99SBarry Smith     }
142247c6ae99SBarry Smith     xnext = xnext->next;
142347c6ae99SBarry Smith     ynext = ynext->next;
142447c6ae99SBarry Smith   }
142547c6ae99SBarry Smith   PetscFunctionReturn(0);
142647c6ae99SBarry Smith }
142747c6ae99SBarry Smith 
142847c6ae99SBarry Smith #undef __FUNCT__
142947c6ae99SBarry Smith #define __FUNCT__ "MatDestroy_Shell_Pack"
143047c6ae99SBarry Smith PetscErrorCode MatDestroy_Shell_Pack(Mat A)
143147c6ae99SBarry Smith {
143247c6ae99SBarry Smith   struct MatPack     *mpack;
143347c6ae99SBarry Smith   struct MatPackLink *anext,*oldanext;
143447c6ae99SBarry Smith   PetscErrorCode     ierr;
143547c6ae99SBarry Smith 
143647c6ae99SBarry Smith   PetscFunctionBegin;
143747c6ae99SBarry Smith   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
143847c6ae99SBarry Smith   anext = mpack->next;
143947c6ae99SBarry Smith 
144047c6ae99SBarry Smith   while (anext) {
144147c6ae99SBarry Smith     ierr     = MatDestroy(anext->A);CHKERRQ(ierr);
144247c6ae99SBarry Smith     oldanext = anext;
144347c6ae99SBarry Smith     anext    = anext->next;
144447c6ae99SBarry Smith     ierr     = PetscFree(oldanext);CHKERRQ(ierr);
144547c6ae99SBarry Smith   }
144647c6ae99SBarry Smith   ierr = PetscFree(mpack);CHKERRQ(ierr);
144747c6ae99SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
144847c6ae99SBarry Smith   PetscFunctionReturn(0);
144947c6ae99SBarry Smith }
145047c6ae99SBarry Smith 
145147c6ae99SBarry Smith #undef __FUNCT__
14520c010503SBarry Smith #define __FUNCT__ "DMGetInterpolation_Composite"
1453*7087cfbeSBarry Smith PetscErrorCode  DMGetInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
145447c6ae99SBarry Smith {
145547c6ae99SBarry Smith   PetscErrorCode         ierr;
145647c6ae99SBarry Smith   PetscInt               m,n,M,N;
145747c6ae99SBarry Smith   struct DMCompositeLink *nextc;
145847c6ae99SBarry Smith   struct DMCompositeLink *nextf;
145947c6ae99SBarry Smith   struct MatPackLink     *nextmat,*pnextmat = 0;
146047c6ae99SBarry Smith   struct MatPack         *mpack;
146147c6ae99SBarry Smith   Vec                    gcoarse,gfine;
146247c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
146347c6ae99SBarry Smith   DM_Composite           *comfine = (DM_Composite*)fine->data;
146447c6ae99SBarry Smith 
146547c6ae99SBarry Smith   PetscFunctionBegin;
146647c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
146747c6ae99SBarry Smith   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
146847c6ae99SBarry Smith   nextc = comcoarse->next;
146947c6ae99SBarry Smith   nextf = comfine->next;
147047c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14710c010503SBarry Smith   ierr = DMCreateGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
14720c010503SBarry Smith   ierr = DMCreateGlobalVector(fine,&gfine);CHKERRQ(ierr);
147347c6ae99SBarry Smith   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
147447c6ae99SBarry Smith   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
147547c6ae99SBarry Smith   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
147647c6ae99SBarry Smith   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
147747c6ae99SBarry Smith   ierr = VecDestroy(gcoarse);CHKERRQ(ierr);
147847c6ae99SBarry Smith   ierr = VecDestroy(gfine);CHKERRQ(ierr);
147947c6ae99SBarry Smith 
148047c6ae99SBarry Smith   ierr         = PetscNew(struct MatPack,&mpack);CHKERRQ(ierr);
148147c6ae99SBarry Smith   mpack->right = coarse;
148247c6ae99SBarry Smith   mpack->left  = fine;
148347c6ae99SBarry Smith   ierr  = MatCreate(((PetscObject)fine)->comm,A);CHKERRQ(ierr);
148447c6ae99SBarry Smith   ierr  = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
148547c6ae99SBarry Smith   ierr  = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
148647c6ae99SBarry Smith   ierr  = MatShellSetContext(*A,mpack);CHKERRQ(ierr);
148747c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);CHKERRQ(ierr);
148847c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);CHKERRQ(ierr);
148947c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);CHKERRQ(ierr);
149047c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);CHKERRQ(ierr);
149147c6ae99SBarry Smith 
149247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
149347c6ae99SBarry Smith   while (nextc) {
149447c6ae99SBarry Smith     if (nextc->type != nextf->type) SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Two DM have different layout");
149547c6ae99SBarry Smith 
149647c6ae99SBarry Smith     if (nextc->type == DMCOMPOSITE_ARRAY) {
149747c6ae99SBarry Smith       ;
149847c6ae99SBarry Smith     } else if (nextc->type == DMCOMPOSITE_DM) {
149947c6ae99SBarry Smith       ierr          = PetscNew(struct MatPackLink,&nextmat);CHKERRQ(ierr);
150047c6ae99SBarry Smith       nextmat->next = 0;
150147c6ae99SBarry Smith       if (pnextmat) {
150247c6ae99SBarry Smith         pnextmat->next = nextmat;
150347c6ae99SBarry Smith         pnextmat       = nextmat;
150447c6ae99SBarry Smith       } else {
150547c6ae99SBarry Smith         pnextmat    = nextmat;
150647c6ae99SBarry Smith         mpack->next = nextmat;
150747c6ae99SBarry Smith       }
150847c6ae99SBarry Smith       ierr = DMGetInterpolation(nextc->dm,nextf->dm,&nextmat->A,PETSC_NULL);CHKERRQ(ierr);
150947c6ae99SBarry Smith     } else {
151047c6ae99SBarry Smith       SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
151147c6ae99SBarry Smith     }
151247c6ae99SBarry Smith     nextc = nextc->next;
151347c6ae99SBarry Smith     nextf = nextf->next;
151447c6ae99SBarry Smith   }
151547c6ae99SBarry Smith   PetscFunctionReturn(0);
151647c6ae99SBarry Smith }
151747c6ae99SBarry Smith 
151847c6ae99SBarry Smith #undef __FUNCT__
15191411c6eeSJed Brown #define __FUNCT__ "DMCreateLocalToGlobalMapping_Composite"
15201411c6eeSJed Brown static PetscErrorCode DMCreateLocalToGlobalMapping_Composite(DM dm)
15211411c6eeSJed Brown {
15221411c6eeSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
15231411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1524a03f9cedSJed Brown   PetscInt               i,cnt,m,*idx;
15251411c6eeSJed Brown   PetscErrorCode         ierr;
15261411c6eeSJed Brown 
15271411c6eeSJed Brown   PetscFunctionBegin;
15281411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
15291411c6eeSJed Brown   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
15301411c6eeSJed Brown   for (cnt=0,i=0; i<(com->nDM+com->nredundant); i++) {
1531a03f9cedSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1532a03f9cedSJed Brown     cnt += m;
1533a03f9cedSJed Brown   }
1534a03f9cedSJed Brown   ierr = PetscMalloc(cnt*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1535a03f9cedSJed Brown   for (cnt=0,i=0; i<(com->nDM+com->nredundant); i++) {
15361411c6eeSJed Brown     const PetscInt *subidx;
15371411c6eeSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
15381411c6eeSJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
15391411c6eeSJed Brown     ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr);
15401411c6eeSJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
15411411c6eeSJed Brown     cnt += m;
15421411c6eeSJed Brown   }
15431411c6eeSJed Brown   ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,cnt,idx,PETSC_OWN_POINTER,&dm->ltogmap);CHKERRQ(ierr);
15441411c6eeSJed Brown   for (i=0; i<com->nDM+com->nredundant; i++) {ierr = ISLocalToGlobalMappingDestroy(ltogs[i]);CHKERRQ(ierr);}
15451411c6eeSJed Brown   ierr = PetscFree(ltogs);CHKERRQ(ierr);
15461411c6eeSJed Brown   PetscFunctionReturn(0);
15471411c6eeSJed Brown }
15481411c6eeSJed Brown 
15491411c6eeSJed Brown 
15501411c6eeSJed Brown #undef __FUNCT__
15510c010503SBarry Smith #define __FUNCT__ "DMGetColoring_Composite"
1552*7087cfbeSBarry Smith PetscErrorCode  DMGetColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
155347c6ae99SBarry Smith {
155447c6ae99SBarry Smith   PetscErrorCode         ierr;
155547c6ae99SBarry Smith   PetscInt               n,i,cnt;
155647c6ae99SBarry Smith   ISColoringValue        *colors;
155747c6ae99SBarry Smith   PetscBool              dense = PETSC_FALSE;
155847c6ae99SBarry Smith   ISColoringValue        maxcol = 0;
155947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
156047c6ae99SBarry Smith 
156147c6ae99SBarry Smith   PetscFunctionBegin;
156247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
156347c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED) {
156447c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" );
156547c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL) {
156647c6ae99SBarry Smith     n = com->n;
156747c6ae99SBarry Smith   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
156847c6ae99SBarry Smith   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
156947c6ae99SBarry Smith 
1570671f6225SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
157147c6ae99SBarry Smith   if (dense) {
157247c6ae99SBarry Smith     for (i=0; i<n; i++) {
157347c6ae99SBarry Smith       colors[i] = (ISColoringValue)(com->rstart + i);
157447c6ae99SBarry Smith     }
157547c6ae99SBarry Smith     maxcol = com->N;
157647c6ae99SBarry Smith   } else {
157747c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
157847c6ae99SBarry Smith     PetscMPIInt            rank;
157947c6ae99SBarry Smith 
158047c6ae99SBarry Smith     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
158147c6ae99SBarry Smith     cnt  = 0;
158247c6ae99SBarry Smith     while (next) {
158347c6ae99SBarry Smith       if (next->type == DMCOMPOSITE_ARRAY) {
158447c6ae99SBarry Smith         if (rank == next->rank) {  /* each column gets is own color */
158547c6ae99SBarry Smith           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
158647c6ae99SBarry Smith             colors[cnt++] = maxcol++;
158747c6ae99SBarry Smith           }
158847c6ae99SBarry Smith         }
158947c6ae99SBarry Smith         ierr = MPI_Bcast(&maxcol,1,MPIU_COLORING_VALUE,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
159047c6ae99SBarry Smith       } else if (next->type == DMCOMPOSITE_DM) {
159147c6ae99SBarry Smith         ISColoring     lcoloring;
159247c6ae99SBarry Smith 
159347c6ae99SBarry Smith         ierr = DMGetColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
159447c6ae99SBarry Smith         for (i=0; i<lcoloring->N; i++) {
159547c6ae99SBarry Smith           colors[cnt++] = maxcol + lcoloring->colors[i];
159647c6ae99SBarry Smith         }
159747c6ae99SBarry Smith         maxcol += lcoloring->n;
159847c6ae99SBarry Smith         ierr = ISColoringDestroy(lcoloring);CHKERRQ(ierr);
159947c6ae99SBarry Smith       } else {
160047c6ae99SBarry Smith         SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
160147c6ae99SBarry Smith       }
160247c6ae99SBarry Smith       next = next->next;
160347c6ae99SBarry Smith     }
160447c6ae99SBarry Smith   }
160547c6ae99SBarry Smith   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
160647c6ae99SBarry Smith   PetscFunctionReturn(0);
160747c6ae99SBarry Smith }
160847c6ae99SBarry Smith 
160947c6ae99SBarry Smith #undef __FUNCT__
16100c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
1611*7087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
161247c6ae99SBarry Smith {
161347c6ae99SBarry Smith   PetscErrorCode         ierr;
161447c6ae99SBarry Smith   struct DMCompositeLink *next;
161547c6ae99SBarry Smith   PetscInt               cnt = 3;
161647c6ae99SBarry Smith   PetscMPIInt            rank;
161747c6ae99SBarry Smith   PetscScalar            *garray,*larray;
161847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
161947c6ae99SBarry Smith 
162047c6ae99SBarry Smith   PetscFunctionBegin;
162147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
162247c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
162347c6ae99SBarry Smith   next = com->next;
162447c6ae99SBarry Smith   if (!com->setup) {
1625d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
162647c6ae99SBarry Smith   }
162747c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
162847c6ae99SBarry Smith   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
162947c6ae99SBarry Smith   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
163047c6ae99SBarry Smith 
163147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
163247c6ae99SBarry Smith   while (next) {
163347c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
163447c6ae99SBarry Smith       if (rank == next->rank) {
163547c6ae99SBarry Smith         ierr    = PetscMemcpy(larray,garray,next->n*sizeof(PetscScalar));CHKERRQ(ierr);
163647c6ae99SBarry Smith         garray += next->n;
163747c6ae99SBarry Smith       }
163847c6ae99SBarry Smith       /* does not handle ADD_VALUES */
163906ebdd98SJed Brown       ierr = MPI_Bcast(larray,next->nlocal,MPIU_SCALAR,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
164047c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
164147c6ae99SBarry Smith       Vec      local,global;
164247c6ae99SBarry Smith       PetscInt N;
164347c6ae99SBarry Smith 
164447c6ae99SBarry Smith       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
164547c6ae99SBarry Smith       ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
164647c6ae99SBarry Smith       ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
164747c6ae99SBarry Smith       ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
164847c6ae99SBarry Smith       ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
164947c6ae99SBarry Smith       ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
165047c6ae99SBarry Smith       ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
165147c6ae99SBarry Smith       ierr = VecResetArray(global);CHKERRQ(ierr);
165247c6ae99SBarry Smith       ierr = VecResetArray(local);CHKERRQ(ierr);
165347c6ae99SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
165447c6ae99SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
165547c6ae99SBarry Smith     } else {
165647c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
165747c6ae99SBarry Smith     }
165847c6ae99SBarry Smith     cnt++;
165906ebdd98SJed Brown     larray += next->nlocal;
166047c6ae99SBarry Smith     next    = next->next;
166147c6ae99SBarry Smith   }
166247c6ae99SBarry Smith 
166347c6ae99SBarry Smith   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
166447c6ae99SBarry Smith   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
166547c6ae99SBarry Smith   PetscFunctionReturn(0);
166647c6ae99SBarry Smith }
166747c6ae99SBarry Smith 
166847c6ae99SBarry Smith #undef __FUNCT__
16690c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
1670*7087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
16710c010503SBarry Smith {
16720c010503SBarry Smith   PetscFunctionBegin;
16730c010503SBarry Smith   PetscFunctionReturn(0);
16740c010503SBarry Smith }
167547c6ae99SBarry Smith 
1676a4121054SBarry Smith EXTERN_C_BEGIN
1677a4121054SBarry Smith #undef __FUNCT__
1678a4121054SBarry Smith #define __FUNCT__ "DMCreate_Composite"
1679*7087cfbeSBarry Smith PetscErrorCode  DMCreate_Composite(DM p)
1680a4121054SBarry Smith {
1681a4121054SBarry Smith   PetscErrorCode ierr;
1682a4121054SBarry Smith   DM_Composite   *com;
1683a4121054SBarry Smith 
1684a4121054SBarry Smith   PetscFunctionBegin;
1685a4121054SBarry Smith   ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1686a4121054SBarry Smith   p->data = com;
1687a4121054SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1688a4121054SBarry Smith   com->n            = 0;
1689a4121054SBarry Smith   com->next         = PETSC_NULL;
1690a4121054SBarry Smith   com->nredundant   = 0;
1691a4121054SBarry Smith   com->nDM          = 0;
1692a4121054SBarry Smith 
1693a4121054SBarry Smith   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1694a4121054SBarry Smith   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
16951411c6eeSJed Brown   p->ops->createlocaltoglobalmapping      = DMCreateLocalToGlobalMapping_Composite;
16961411c6eeSJed Brown   p->ops->createlocaltoglobalmappingblock = 0;
1697a4121054SBarry Smith   p->ops->refine                          = DMRefine_Composite;
1698a4121054SBarry Smith   p->ops->getinterpolation                = DMGetInterpolation_Composite;
1699a4121054SBarry Smith   p->ops->getmatrix                       = DMGetMatrix_Composite;
1700a4121054SBarry Smith   p->ops->getcoloring                     = DMGetColoring_Composite;
1701a4121054SBarry Smith   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1702a4121054SBarry Smith   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1703a4121054SBarry Smith   p->ops->destroy                         = DMDestroy_Composite;
1704a4121054SBarry Smith   p->ops->view                            = DMView_Composite;
1705a4121054SBarry Smith   p->ops->setup                           = DMSetUp_Composite;
1706a4121054SBarry Smith   PetscFunctionReturn(0);
1707a4121054SBarry Smith }
1708a4121054SBarry Smith EXTERN_C_END
1709a4121054SBarry Smith 
17100c010503SBarry Smith #undef __FUNCT__
17110c010503SBarry Smith #define __FUNCT__ "DMCompositeCreate"
17120c010503SBarry Smith /*@C
17130c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
17140c010503SBarry Smith       vectors made up of several subvectors.
17150c010503SBarry Smith 
17160c010503SBarry Smith     Collective on MPI_Comm
171747c6ae99SBarry Smith 
171847c6ae99SBarry Smith     Input Parameter:
17190c010503SBarry Smith .   comm - the processors that will share the global vector
17200c010503SBarry Smith 
17210c010503SBarry Smith     Output Parameters:
17220c010503SBarry Smith .   packer - the packer object
172347c6ae99SBarry Smith 
172447c6ae99SBarry Smith     Level: advanced
172547c6ae99SBarry Smith 
17260c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCompositeScatter(),
17276eb61c8cSJed Brown          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
172847c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
172947c6ae99SBarry Smith 
173047c6ae99SBarry Smith @*/
1731*7087cfbeSBarry Smith PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
173247c6ae99SBarry Smith {
17330c010503SBarry Smith   PetscErrorCode ierr;
17340c010503SBarry Smith 
173547c6ae99SBarry Smith   PetscFunctionBegin;
17360c010503SBarry Smith   PetscValidPointer(packer,2);
1737a4121054SBarry Smith   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1738a4121054SBarry Smith   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
173947c6ae99SBarry Smith   PetscFunctionReturn(0);
174047c6ae99SBarry Smith }
1741