xref: /petsc/src/dm/impls/composite/pack.c (revision 3c0c59f39ee3bd561a1b4dcbd0e58366bfa60aba)
147c6ae99SBarry Smith 
23c48a1e8SJed Brown #include "packimpl.h" /*I   "petscdmcomposite.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #undef __FUNCT__
547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeSetCoupling"
647c6ae99SBarry Smith /*@C
747c6ae99SBarry Smith     DMCompositeSetCoupling - Sets user provided routines that compute the coupling between the
8aa219208SBarry Smith       seperate components (DMDA's and arrays) in a DMto build the correct matrix nonzero structure.
947c6ae99SBarry Smith 
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith     Logically Collective on MPI_Comm
1247c6ae99SBarry Smith 
1347c6ae99SBarry Smith     Input Parameter:
1447c6ae99SBarry Smith +   dm - the composite object
1547c6ae99SBarry Smith -   formcouplelocations - routine to set the nonzero locations in the matrix
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith     Level: advanced
1847c6ae99SBarry Smith 
1947c6ae99SBarry Smith     Notes: See DMCompositeSetContext() and DMCompositeGetContext() for how to get user information into
2047c6ae99SBarry Smith         this routine
2147c6ae99SBarry Smith 
2247c6ae99SBarry Smith @*/
237087cfbeSBarry Smith PetscErrorCode  DMCompositeSetCoupling(DM dm,PetscErrorCode (*FormCoupleLocations)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt))
2447c6ae99SBarry Smith {
2547c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   PetscFunctionBegin;
2847c6ae99SBarry Smith   com->FormCoupleLocations = FormCoupleLocations;
2947c6ae99SBarry Smith   PetscFunctionReturn(0);
3047c6ae99SBarry Smith }
3147c6ae99SBarry Smith 
3247c6ae99SBarry Smith #undef __FUNCT__
3347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeSetContext"
3447c6ae99SBarry Smith /*@
3547c6ae99SBarry Smith     DMCompositeSetContext - Allows user to stash data they may need within the form coupling routine they
3647c6ae99SBarry Smith       set with DMCompositeSetCoupling()
3747c6ae99SBarry Smith 
3847c6ae99SBarry Smith 
3947c6ae99SBarry Smith     Not Collective
4047c6ae99SBarry Smith 
4147c6ae99SBarry Smith     Input Parameter:
4247c6ae99SBarry Smith +   dm - the composite object
4347c6ae99SBarry Smith -   ctx - the user supplied context
4447c6ae99SBarry Smith 
4547c6ae99SBarry Smith     Level: advanced
4647c6ae99SBarry Smith 
4747c6ae99SBarry Smith     Notes: Use DMCompositeGetContext() to retrieve the context when needed.
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith @*/
507087cfbeSBarry Smith PetscErrorCode  DMCompositeSetContext(DM dm,void *ctx)
5147c6ae99SBarry Smith {
5247c6ae99SBarry Smith   PetscFunctionBegin;
5347c6ae99SBarry Smith   dm->ctx = ctx;
5447c6ae99SBarry Smith   PetscFunctionReturn(0);
5547c6ae99SBarry Smith }
5647c6ae99SBarry Smith 
5747c6ae99SBarry Smith #undef __FUNCT__
5847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetContext"
5947c6ae99SBarry Smith /*@
6047c6ae99SBarry Smith     DMCompositeGetContext - Access the context set with DMCompositeSetContext()
6147c6ae99SBarry Smith 
6247c6ae99SBarry Smith 
6347c6ae99SBarry Smith     Not Collective
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith     Input Parameter:
6647c6ae99SBarry Smith .   dm - the composite object
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith     Output Parameter:
6947c6ae99SBarry Smith .    ctx - the user supplied context
7047c6ae99SBarry Smith 
7147c6ae99SBarry Smith     Level: advanced
7247c6ae99SBarry Smith 
7347c6ae99SBarry Smith     Notes: Use DMCompositeGetContext() to retrieve the context when needed.
7447c6ae99SBarry Smith 
7547c6ae99SBarry Smith @*/
767087cfbeSBarry Smith PetscErrorCode  DMCompositeGetContext(DM dm,void **ctx)
7747c6ae99SBarry Smith {
7847c6ae99SBarry Smith   PetscFunctionBegin;
7947c6ae99SBarry Smith   *ctx = dm->ctx;
8047c6ae99SBarry Smith   PetscFunctionReturn(0);
8147c6ae99SBarry Smith }
8247c6ae99SBarry Smith 
8347c6ae99SBarry Smith #undef __FUNCT__
840c010503SBarry Smith #define __FUNCT__ "DMDestroy_Composite"
856bf464f9SBarry Smith PetscErrorCode  DMDestroy_Composite(DM dm)
8647c6ae99SBarry Smith {
8747c6ae99SBarry Smith   PetscErrorCode         ierr;
8847c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
8947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
9047c6ae99SBarry Smith 
9147c6ae99SBarry Smith   PetscFunctionBegin;
9247c6ae99SBarry Smith   next = com->next;
9347c6ae99SBarry Smith   while (next) {
9447c6ae99SBarry Smith     prev = next;
9547c6ae99SBarry Smith     next = next->next;
9647c6ae99SBarry Smith     if (prev->type == DMCOMPOSITE_DM) {
97fcfd50ebSBarry Smith       ierr = DMDestroy(&prev->dm);CHKERRQ(ierr);
9847c6ae99SBarry Smith     }
9947c6ae99SBarry Smith     ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
10047c6ae99SBarry Smith     ierr = PetscFree(prev);CHKERRQ(ierr);
10147c6ae99SBarry Smith   }
10247c6ae99SBarry Smith   PetscFunctionReturn(0);
10347c6ae99SBarry Smith }
10447c6ae99SBarry Smith 
10547c6ae99SBarry Smith #undef __FUNCT__
1060c010503SBarry Smith #define __FUNCT__ "DMView_Composite"
1077087cfbeSBarry Smith PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
10847c6ae99SBarry Smith {
10947c6ae99SBarry Smith   PetscErrorCode ierr;
11047c6ae99SBarry Smith   PetscBool      iascii;
11147c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite *)dm->data;
11247c6ae99SBarry Smith 
11347c6ae99SBarry Smith   PetscFunctionBegin;
11447c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
11547c6ae99SBarry Smith   if (iascii) {
11647c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
11747c6ae99SBarry Smith     PetscInt               i;
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix?((PetscObject)dm)->prefix:"no prefix");CHKERRQ(ierr);
12047c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"  contains %d DMs and %d redundant arrays\n",com->nDM,com->nredundant);CHKERRQ(ierr);
12147c6ae99SBarry Smith     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
12247c6ae99SBarry Smith     for (i=0; lnk; lnk=lnk->next,i++) {
12347c6ae99SBarry Smith       if (lnk->dm) {
12447c6ae99SBarry Smith         ierr = PetscViewerASCIIPrintf(v,"Link %d: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
12547c6ae99SBarry Smith         ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
12647c6ae99SBarry Smith         ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
12747c6ae99SBarry Smith         ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
12847c6ae99SBarry Smith       } else {
12906ebdd98SJed Brown         ierr = PetscViewerASCIIPrintf(v,"Link %d: Redundant array of size %d owned by rank %d\n",i,lnk->nlocal,lnk->rank);CHKERRQ(ierr);
13047c6ae99SBarry Smith       }
13147c6ae99SBarry Smith     }
13247c6ae99SBarry Smith     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
13347c6ae99SBarry Smith   }
13447c6ae99SBarry Smith   PetscFunctionReturn(0);
13547c6ae99SBarry Smith }
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
13847c6ae99SBarry Smith #undef __FUNCT__
139d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp_Composite"
1407087cfbeSBarry Smith PetscErrorCode  DMSetUp_Composite(DM dm)
14147c6ae99SBarry Smith {
14247c6ae99SBarry Smith   PetscErrorCode         ierr;
14347c6ae99SBarry Smith   PetscInt               nprev = 0;
14447c6ae99SBarry Smith   PetscMPIInt            rank,size;
14547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
14647c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
14747c6ae99SBarry Smith   PetscLayout            map;
14847c6ae99SBarry Smith 
14947c6ae99SBarry Smith   PetscFunctionBegin;
15047c6ae99SBarry Smith   if (com->setup) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
15147c6ae99SBarry Smith   ierr = PetscLayoutCreate(((PetscObject)dm)->comm,&map);CHKERRQ(ierr);
15247c6ae99SBarry Smith   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
15347c6ae99SBarry Smith   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
15447c6ae99SBarry Smith   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
15547c6ae99SBarry Smith   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
15647c6ae99SBarry Smith   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
15747c6ae99SBarry Smith   ierr = PetscLayoutGetRange(map,&com->rstart,PETSC_NULL);CHKERRQ(ierr);
158fcfd50ebSBarry Smith   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
15947c6ae99SBarry Smith 
16047c6ae99SBarry Smith   /* now set the rstart for each linked array/vector */
16147c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
16247c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)dm)->comm,&size);CHKERRQ(ierr);
16347c6ae99SBarry Smith   while (next) {
16447c6ae99SBarry Smith     next->rstart = nprev;
16506ebdd98SJed Brown     nprev += next->n;
16647c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
16747c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
16847c6ae99SBarry Smith       ierr = MPI_Bcast(&next->grstart,1,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
16947c6ae99SBarry Smith     } else {
17047c6ae99SBarry Smith       ierr = PetscMalloc(size*sizeof(PetscInt),&next->grstarts);CHKERRQ(ierr);
17147c6ae99SBarry Smith       ierr = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
17247c6ae99SBarry Smith     }
17347c6ae99SBarry Smith     next = next->next;
17447c6ae99SBarry Smith   }
17547c6ae99SBarry Smith   com->setup = PETSC_TRUE;
17647c6ae99SBarry Smith   PetscFunctionReturn(0);
17747c6ae99SBarry Smith }
17847c6ae99SBarry Smith 
17947c6ae99SBarry Smith 
18047c6ae99SBarry Smith #undef __FUNCT__
18147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess_Array"
18247c6ae99SBarry Smith PetscErrorCode DMCompositeGetAccess_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar **array)
18347c6ae99SBarry Smith {
18447c6ae99SBarry Smith   PetscErrorCode ierr;
18547c6ae99SBarry Smith   PetscScalar    *varray;
18647c6ae99SBarry Smith   PetscMPIInt    rank;
18747c6ae99SBarry Smith 
18847c6ae99SBarry Smith   PetscFunctionBegin;
18947c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
19047c6ae99SBarry Smith   if (array) {
19147c6ae99SBarry Smith     if (rank == mine->rank) {
19247c6ae99SBarry Smith       ierr    = VecGetArray(vec,&varray);CHKERRQ(ierr);
19347c6ae99SBarry Smith       *array  = varray + mine->rstart;
19447c6ae99SBarry Smith       ierr    = VecRestoreArray(vec,&varray);CHKERRQ(ierr);
19547c6ae99SBarry Smith     } else {
19647c6ae99SBarry Smith       *array = 0;
19747c6ae99SBarry Smith     }
19847c6ae99SBarry Smith   }
19947c6ae99SBarry Smith   PetscFunctionReturn(0);
20047c6ae99SBarry Smith }
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith #undef __FUNCT__
20347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess_DM"
20447c6ae99SBarry Smith PetscErrorCode DMCompositeGetAccess_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec *global)
20547c6ae99SBarry Smith {
20647c6ae99SBarry Smith   PetscErrorCode ierr;
20747c6ae99SBarry Smith   PetscScalar    *array;
20847c6ae99SBarry Smith 
20947c6ae99SBarry Smith   PetscFunctionBegin;
21047c6ae99SBarry Smith   if (global) {
21147c6ae99SBarry Smith     ierr    = DMGetGlobalVector(mine->dm,global);CHKERRQ(ierr);
21247c6ae99SBarry Smith     ierr    = VecGetArray(vec,&array);CHKERRQ(ierr);
21347c6ae99SBarry Smith     ierr    = VecPlaceArray(*global,array+mine->rstart);CHKERRQ(ierr);
21447c6ae99SBarry Smith     ierr    = VecRestoreArray(vec,&array);CHKERRQ(ierr);
21547c6ae99SBarry Smith   }
21647c6ae99SBarry Smith   PetscFunctionReturn(0);
21747c6ae99SBarry Smith }
21847c6ae99SBarry Smith 
21947c6ae99SBarry Smith #undef __FUNCT__
22047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess_Array"
22147c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreAccess_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar **array)
22247c6ae99SBarry Smith {
22347c6ae99SBarry Smith   PetscFunctionBegin;
22447c6ae99SBarry Smith   PetscFunctionReturn(0);
22547c6ae99SBarry Smith }
22647c6ae99SBarry Smith 
22747c6ae99SBarry Smith #undef __FUNCT__
22847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess_DM"
22947c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreAccess_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec *global)
23047c6ae99SBarry Smith {
23147c6ae99SBarry Smith   PetscErrorCode ierr;
23247c6ae99SBarry Smith 
23347c6ae99SBarry Smith   PetscFunctionBegin;
23447c6ae99SBarry Smith   if (global) {
23547c6ae99SBarry Smith     ierr = VecResetArray(*global);CHKERRQ(ierr);
23647c6ae99SBarry Smith     ierr = DMRestoreGlobalVector(mine->dm,global);CHKERRQ(ierr);
23747c6ae99SBarry Smith   }
23847c6ae99SBarry Smith   PetscFunctionReturn(0);
23947c6ae99SBarry Smith }
24047c6ae99SBarry Smith 
24147c6ae99SBarry Smith #undef __FUNCT__
24247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter_Array"
24347c6ae99SBarry Smith PetscErrorCode DMCompositeScatter_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar *array)
24447c6ae99SBarry Smith {
24547c6ae99SBarry Smith   PetscErrorCode ierr;
24647c6ae99SBarry Smith   PetscScalar    *varray;
24747c6ae99SBarry Smith   PetscMPIInt    rank;
24847c6ae99SBarry Smith 
24947c6ae99SBarry Smith   PetscFunctionBegin;
25047c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
25147c6ae99SBarry Smith   if (rank == mine->rank) {
25247c6ae99SBarry Smith     ierr    = VecGetArray(vec,&varray);CHKERRQ(ierr);
25347c6ae99SBarry Smith     ierr    = PetscMemcpy(array,varray+mine->rstart,mine->n*sizeof(PetscScalar));CHKERRQ(ierr);
25447c6ae99SBarry Smith     ierr    = VecRestoreArray(vec,&varray);CHKERRQ(ierr);
25547c6ae99SBarry Smith   }
25606ebdd98SJed Brown   ierr    = MPI_Bcast(array,mine->nlocal,MPIU_SCALAR,mine->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
25747c6ae99SBarry Smith   PetscFunctionReturn(0);
25847c6ae99SBarry Smith }
25947c6ae99SBarry Smith 
26047c6ae99SBarry Smith #undef __FUNCT__
26147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter_DM"
26247c6ae99SBarry Smith PetscErrorCode DMCompositeScatter_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec local)
26347c6ae99SBarry Smith {
26447c6ae99SBarry Smith   PetscErrorCode ierr;
26547c6ae99SBarry Smith   PetscScalar    *array;
26647c6ae99SBarry Smith   Vec            global;
26747c6ae99SBarry Smith 
26847c6ae99SBarry Smith   PetscFunctionBegin;
26947c6ae99SBarry Smith   ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr);
27047c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
27147c6ae99SBarry Smith   ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr);
27247c6ae99SBarry Smith   ierr = DMGlobalToLocalBegin(mine->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
27347c6ae99SBarry Smith   ierr = DMGlobalToLocalEnd(mine->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
27447c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
27547c6ae99SBarry Smith   ierr = VecResetArray(global);CHKERRQ(ierr);
27647c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr);
27747c6ae99SBarry Smith   PetscFunctionReturn(0);
27847c6ae99SBarry Smith }
27947c6ae99SBarry Smith 
28047c6ae99SBarry Smith #undef __FUNCT__
28147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_Array"
282acc1e270SJed Brown PetscErrorCode DMCompositeGather_Array(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,const PetscScalar *array)
28347c6ae99SBarry Smith {
28447c6ae99SBarry Smith   PetscErrorCode ierr;
28547c6ae99SBarry Smith   PetscScalar    *varray;
28647c6ae99SBarry Smith   PetscMPIInt    rank;
28747c6ae99SBarry Smith 
28847c6ae99SBarry Smith   PetscFunctionBegin;
28947c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
29047c6ae99SBarry Smith   if (rank == mine->rank) {
29147c6ae99SBarry Smith     ierr = VecGetArray(vec,&varray);CHKERRQ(ierr);
29247c6ae99SBarry Smith     if (varray+mine->rstart == array) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"You need not DMCompositeGather() into objects obtained via DMCompositeGetAccess()");
293acc1e270SJed Brown   }
294df0c820aSJed Brown   switch (imode) {
295df0c820aSJed Brown   case INSERT_VALUES:
296acc1e270SJed Brown     if (rank == mine->rank) {
29747c6ae99SBarry Smith       ierr = PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(PetscScalar));CHKERRQ(ierr);
298acc1e270SJed Brown     }
299df0c820aSJed Brown     break;
300df0c820aSJed Brown   case ADD_VALUES: {
301df0c820aSJed Brown     PetscInt          i;
302760fd489SMatthew G Knepley     void             *source;
303acc1e270SJed Brown     PetscScalar       *buffer,*dest;
304acc1e270SJed Brown     if (rank == mine->rank) {
305acc1e270SJed Brown       dest = &varray[mine->rstart];
306acc1e270SJed Brown #if defined(PETSC_HAVE_MPI_IN_PLACE)
307acc1e270SJed Brown       buffer = dest;
308acc1e270SJed Brown       source = MPI_IN_PLACE;
309acc1e270SJed Brown #else
31006ebdd98SJed Brown       ierr = PetscMalloc(mine->nlocal*sizeof(PetscScalar),&buffer);CHKERRQ(ierr);
311760fd489SMatthew G Knepley       source = (void *) buffer;
312acc1e270SJed Brown #endif
31306ebdd98SJed Brown       for (i=0; i<mine->nlocal; i++) buffer[i] = varray[mine->rstart+i] + array[i];
314acc1e270SJed Brown     } else {
315760fd489SMatthew G Knepley       source = (void *) array;
316acc1e270SJed Brown       dest   = PETSC_NULL;
317acc1e270SJed Brown     }
31806ebdd98SJed Brown     ierr = MPI_Reduce(source,dest,mine->nlocal,MPIU_SCALAR,MPI_SUM,mine->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
319acc1e270SJed Brown #if !defined(PETSC_HAVE_MPI_IN_PLACE)
320acc1e270SJed Brown     if (rank == mine->rank) {ierr = PetscFree(source);CHKERRQ(ierr);}
321acc1e270SJed Brown #endif
322df0c820aSJed Brown   } break;
323df0c820aSJed Brown   default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"imode");
324df0c820aSJed Brown   }
325acc1e270SJed Brown   if (rank == mine->rank) {ierr = VecRestoreArray(vec,&varray);CHKERRQ(ierr);}
32647c6ae99SBarry Smith   PetscFunctionReturn(0);
32747c6ae99SBarry Smith }
32847c6ae99SBarry Smith 
32947c6ae99SBarry Smith #undef __FUNCT__
33047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_DM"
331df0c820aSJed Brown PetscErrorCode DMCompositeGather_DM(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,Vec local)
33247c6ae99SBarry Smith {
33347c6ae99SBarry Smith   PetscErrorCode ierr;
33447c6ae99SBarry Smith   PetscScalar    *array;
33547c6ae99SBarry Smith   Vec            global;
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith   PetscFunctionBegin;
33847c6ae99SBarry Smith   ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr);
33947c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
34047c6ae99SBarry Smith   ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr);
341df0c820aSJed Brown   ierr = DMLocalToGlobalBegin(mine->dm,local,imode,global);CHKERRQ(ierr);
342df0c820aSJed Brown   ierr = DMLocalToGlobalEnd(mine->dm,local,imode,global);CHKERRQ(ierr);
34347c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
34447c6ae99SBarry Smith   ierr = VecResetArray(global);CHKERRQ(ierr);
34547c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr);
34647c6ae99SBarry Smith   PetscFunctionReturn(0);
34747c6ae99SBarry Smith }
34847c6ae99SBarry Smith 
34947c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
35047c6ae99SBarry Smith 
35147c6ae99SBarry Smith #include <stdarg.h>
35247c6ae99SBarry Smith 
35347c6ae99SBarry Smith #undef __FUNCT__
35447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetNumberDM"
35547c6ae99SBarry Smith /*@C
35647c6ae99SBarry Smith     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
35747c6ae99SBarry Smith        representation.
35847c6ae99SBarry Smith 
35947c6ae99SBarry Smith     Not Collective
36047c6ae99SBarry Smith 
36147c6ae99SBarry Smith     Input Parameter:
36247c6ae99SBarry Smith .    dm - the packer object
36347c6ae99SBarry Smith 
36447c6ae99SBarry Smith     Output Parameter:
36547c6ae99SBarry Smith .     nDM - the number of DMs
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith     Level: beginner
36847c6ae99SBarry Smith 
36947c6ae99SBarry Smith @*/
3707087cfbeSBarry Smith PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
37147c6ae99SBarry Smith {
37247c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
37347c6ae99SBarry Smith   PetscFunctionBegin;
37447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
37547c6ae99SBarry Smith   *nDM = com->nDM;
37647c6ae99SBarry Smith   PetscFunctionReturn(0);
37747c6ae99SBarry Smith }
37847c6ae99SBarry Smith 
37947c6ae99SBarry Smith 
38047c6ae99SBarry Smith #undef __FUNCT__
38147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess"
38247c6ae99SBarry Smith /*@C
38347c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
38447c6ae99SBarry Smith        representation.
38547c6ae99SBarry Smith 
38647c6ae99SBarry Smith     Collective on DMComposite
38747c6ae99SBarry Smith 
38847c6ae99SBarry Smith     Input Parameter:
38947c6ae99SBarry Smith +    dm - the packer object
39047c6ae99SBarry Smith .    gvec - the global vector
39147c6ae99SBarry Smith -    ... - the individual sequential or parallel objects (arrays or vectors)
39247c6ae99SBarry Smith 
39347c6ae99SBarry Smith     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
39447c6ae99SBarry Smith 
39547c6ae99SBarry Smith     Level: advanced
39647c6ae99SBarry Smith 
39747c6ae99SBarry Smith @*/
3987087cfbeSBarry Smith PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
39947c6ae99SBarry Smith {
40047c6ae99SBarry Smith   va_list                Argp;
40147c6ae99SBarry Smith   PetscErrorCode         ierr;
40247c6ae99SBarry Smith   struct DMCompositeLink *next;
40347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
40447c6ae99SBarry Smith 
40547c6ae99SBarry Smith   PetscFunctionBegin;
40647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
40747c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
40847c6ae99SBarry Smith   next = com->next;
40947c6ae99SBarry Smith   if (!com->setup) {
410d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
41147c6ae99SBarry Smith   }
41247c6ae99SBarry Smith 
41347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
41447c6ae99SBarry Smith   va_start(Argp,gvec);
41547c6ae99SBarry Smith   while (next) {
41647c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
41747c6ae99SBarry Smith       PetscScalar **array;
41847c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
41947c6ae99SBarry Smith       ierr  = DMCompositeGetAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
42047c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
42147c6ae99SBarry Smith       Vec *vec;
42247c6ae99SBarry Smith       vec  = va_arg(Argp, Vec*);
42347c6ae99SBarry Smith       ierr = DMCompositeGetAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
42447c6ae99SBarry Smith     } else {
42547c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
42647c6ae99SBarry Smith     }
42747c6ae99SBarry Smith     next = next->next;
42847c6ae99SBarry Smith   }
42947c6ae99SBarry Smith   va_end(Argp);
43047c6ae99SBarry Smith   PetscFunctionReturn(0);
43147c6ae99SBarry Smith }
43247c6ae99SBarry Smith 
43347c6ae99SBarry Smith #undef __FUNCT__
43447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess"
43547c6ae99SBarry Smith /*@C
436aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
43747c6ae99SBarry Smith        representation.
43847c6ae99SBarry Smith 
43947c6ae99SBarry Smith     Collective on DMComposite
44047c6ae99SBarry Smith 
44147c6ae99SBarry Smith     Input Parameter:
44247c6ae99SBarry Smith +    dm - the packer object
44347c6ae99SBarry Smith .    gvec - the global vector
44447c6ae99SBarry Smith -    ... - the individual sequential or parallel objects (arrays or vectors)
44547c6ae99SBarry Smith 
44647c6ae99SBarry Smith     Level: advanced
44747c6ae99SBarry Smith 
4480c010503SBarry Smith .seealso  DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
4496eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
450aa219208SBarry Smith          DMCompositeRestoreAccess(), DMCompositeGetAccess()
45147c6ae99SBarry Smith 
45247c6ae99SBarry Smith @*/
4537087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
45447c6ae99SBarry Smith {
45547c6ae99SBarry Smith   va_list                Argp;
45647c6ae99SBarry Smith   PetscErrorCode         ierr;
45747c6ae99SBarry Smith   struct DMCompositeLink *next;
45847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
45947c6ae99SBarry Smith 
46047c6ae99SBarry Smith   PetscFunctionBegin;
46147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
46247c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
46347c6ae99SBarry Smith   next = com->next;
46447c6ae99SBarry Smith   if (!com->setup) {
465d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
46647c6ae99SBarry Smith   }
46747c6ae99SBarry Smith 
46847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
46947c6ae99SBarry Smith   va_start(Argp,gvec);
47047c6ae99SBarry Smith   while (next) {
47147c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
47247c6ae99SBarry Smith       PetscScalar **array;
47347c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
47447c6ae99SBarry Smith       ierr  = DMCompositeRestoreAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
47547c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
47647c6ae99SBarry Smith       Vec *vec;
47747c6ae99SBarry Smith       vec  = va_arg(Argp, Vec*);
47847c6ae99SBarry Smith       ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
47947c6ae99SBarry Smith     } else {
48047c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
48147c6ae99SBarry Smith     }
48247c6ae99SBarry Smith     next = next->next;
48347c6ae99SBarry Smith   }
48447c6ae99SBarry Smith   va_end(Argp);
48547c6ae99SBarry Smith   PetscFunctionReturn(0);
48647c6ae99SBarry Smith }
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith #undef __FUNCT__
48947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter"
49047c6ae99SBarry Smith /*@C
49147c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith     Collective on DMComposite
49447c6ae99SBarry Smith 
49547c6ae99SBarry Smith     Input Parameter:
49647c6ae99SBarry Smith +    dm - the packer object
49747c6ae99SBarry Smith .    gvec - the global vector
4988fd8f222SJed Brown -    ... - the individual sequential objects (arrays or vectors), PETSC_NULL for those that are not needed
49947c6ae99SBarry Smith 
50047c6ae99SBarry Smith     Level: advanced
50147c6ae99SBarry Smith 
5020c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
5036eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
50447c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith @*/
5077087cfbeSBarry Smith PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
50847c6ae99SBarry Smith {
50947c6ae99SBarry Smith   va_list                Argp;
51047c6ae99SBarry Smith   PetscErrorCode         ierr;
51147c6ae99SBarry Smith   struct DMCompositeLink *next;
5128fd8f222SJed Brown   PetscInt               cnt;
51347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
51447c6ae99SBarry Smith 
51547c6ae99SBarry Smith   PetscFunctionBegin;
51647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
51747c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
51847c6ae99SBarry Smith   if (!com->setup) {
519d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
52047c6ae99SBarry Smith   }
52147c6ae99SBarry Smith 
52247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
52347c6ae99SBarry Smith   va_start(Argp,gvec);
5248fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
52547c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
52647c6ae99SBarry Smith       PetscScalar *array;
52747c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar*);
528c72eef85SJed Brown       if (array) PetscValidScalarPointer(array,cnt);
529c72eef85SJed Brown       PetscValidLogicalCollectiveBool(dm,!!array,cnt);
5308fd8f222SJed Brown       if (!array) continue;
53147c6ae99SBarry Smith       ierr = DMCompositeScatter_Array(dm,next,gvec,array);CHKERRQ(ierr);
53247c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
53347c6ae99SBarry Smith       Vec vec;
53447c6ae99SBarry Smith       vec = va_arg(Argp, Vec);
5358fd8f222SJed Brown       if (!vec) continue;
53647c6ae99SBarry Smith       PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt);
53747c6ae99SBarry Smith       ierr = DMCompositeScatter_DM(dm,next,gvec,vec);CHKERRQ(ierr);
53847c6ae99SBarry Smith     } else {
53947c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
54047c6ae99SBarry Smith     }
54147c6ae99SBarry Smith   }
54247c6ae99SBarry Smith   va_end(Argp);
54347c6ae99SBarry Smith   PetscFunctionReturn(0);
54447c6ae99SBarry Smith }
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith #undef __FUNCT__
54747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather"
54847c6ae99SBarry Smith /*@C
54947c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
55047c6ae99SBarry Smith 
55147c6ae99SBarry Smith     Collective on DMComposite
55247c6ae99SBarry Smith 
55347c6ae99SBarry Smith     Input Parameter:
55447c6ae99SBarry Smith +    dm - the packer object
55547c6ae99SBarry Smith .    gvec - the global vector
5568fd8f222SJed Brown -    ... - the individual sequential objects (arrays or vectors), PETSC_NULL for any that are not needed
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith     Level: advanced
55947c6ae99SBarry Smith 
5600c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
5616eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
56247c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
56347c6ae99SBarry Smith 
56447c6ae99SBarry Smith @*/
5657087cfbeSBarry Smith PetscErrorCode  DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
56647c6ae99SBarry Smith {
56747c6ae99SBarry Smith   va_list                Argp;
56847c6ae99SBarry Smith   PetscErrorCode         ierr;
56947c6ae99SBarry Smith   struct DMCompositeLink *next;
57047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
5718fd8f222SJed Brown   PetscInt               cnt;
57247c6ae99SBarry Smith 
57347c6ae99SBarry Smith   PetscFunctionBegin;
57447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
57547c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
57647c6ae99SBarry Smith   if (!com->setup) {
577d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
57847c6ae99SBarry Smith   }
57947c6ae99SBarry Smith 
58047c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
581df0c820aSJed Brown   va_start(Argp,imode);
5828fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
58347c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
58447c6ae99SBarry Smith       PetscScalar *array;
58547c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar*);
5868fd8f222SJed Brown       if (!array) continue;
5878fd8f222SJed Brown       PetscValidScalarPointer(array,cnt);
588df0c820aSJed Brown       ierr  = DMCompositeGather_Array(dm,next,gvec,imode,array);CHKERRQ(ierr);
58947c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
59047c6ae99SBarry Smith       Vec vec;
59147c6ae99SBarry Smith       vec = va_arg(Argp, Vec);
5928fd8f222SJed Brown       if (!vec) continue;
5938fd8f222SJed Brown       PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt);
594df0c820aSJed Brown       ierr = DMCompositeGather_DM(dm,next,gvec,imode,vec);CHKERRQ(ierr);
59547c6ae99SBarry Smith     } else {
59647c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
59747c6ae99SBarry Smith     }
59847c6ae99SBarry Smith   }
59947c6ae99SBarry Smith   va_end(Argp);
60047c6ae99SBarry Smith   PetscFunctionReturn(0);
60147c6ae99SBarry Smith }
60247c6ae99SBarry Smith 
60347c6ae99SBarry Smith #undef __FUNCT__
60447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddArray"
60547c6ae99SBarry Smith /*@C
60647c6ae99SBarry Smith     DMCompositeAddArray - adds an "redundant" array to a DMComposite. The array values will
60747c6ae99SBarry Smith        be stored in part of the array on process orank.
60847c6ae99SBarry Smith 
60947c6ae99SBarry Smith     Collective on DMComposite
61047c6ae99SBarry Smith 
61147c6ae99SBarry Smith     Input Parameter:
61247c6ae99SBarry Smith +    dm - the packer object
61347c6ae99SBarry Smith .    orank - the process on which the array entries officially live, this number must be
61447c6ae99SBarry Smith              the same on all processes.
61547c6ae99SBarry Smith -    n - the length of the array
61647c6ae99SBarry Smith 
61747c6ae99SBarry Smith     Level: advanced
61847c6ae99SBarry Smith 
6190c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
6206eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
62147c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith @*/
6247087cfbeSBarry Smith PetscErrorCode  DMCompositeAddArray(DM dm,PetscMPIInt orank,PetscInt n)
62547c6ae99SBarry Smith {
62647c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
62747c6ae99SBarry Smith   PetscErrorCode         ierr;
62847c6ae99SBarry Smith   PetscMPIInt            rank;
62947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
63047c6ae99SBarry Smith 
63147c6ae99SBarry Smith   PetscFunctionBegin;
63247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
63347c6ae99SBarry Smith   next = com->next;
63447c6ae99SBarry Smith   if (com->setup) {
63547c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add an array once you have used the DMComposite");
63647c6ae99SBarry Smith   }
63747c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
63847c6ae99SBarry Smith   {
63947c6ae99SBarry Smith     PetscMPIInt orankmax;
64047c6ae99SBarry Smith     ierr = MPI_Allreduce(&orank,&orankmax,1,MPI_INT,MPI_MAX,((PetscObject)dm)->comm);CHKERRQ(ierr);
64147c6ae99SBarry 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);
64247c6ae99SBarry Smith   }
64347c6ae99SBarry Smith #endif
64447c6ae99SBarry Smith 
64547c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
64647c6ae99SBarry Smith   /* create new link */
64747c6ae99SBarry Smith   ierr                = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
64806ebdd98SJed Brown   mine->nlocal        = n;
64906ebdd98SJed Brown   mine->n             = (rank == orank) ? n : 0;
65047c6ae99SBarry Smith   mine->rank          = orank;
65147c6ae99SBarry Smith   mine->dm            = PETSC_NULL;
65247c6ae99SBarry Smith   mine->type          = DMCOMPOSITE_ARRAY;
65347c6ae99SBarry Smith   mine->next          = PETSC_NULL;
65447c6ae99SBarry Smith   if (rank == mine->rank) {com->n += n;com->nmine++;}
65547c6ae99SBarry Smith 
65647c6ae99SBarry Smith   /* add to end of list */
65747c6ae99SBarry Smith   if (!next) {
65847c6ae99SBarry Smith     com->next = mine;
65947c6ae99SBarry Smith   } else {
66047c6ae99SBarry Smith     while (next->next) next = next->next;
66147c6ae99SBarry Smith     next->next = mine;
66247c6ae99SBarry Smith   }
66347c6ae99SBarry Smith   com->nredundant++;
66447c6ae99SBarry Smith   PetscFunctionReturn(0);
66547c6ae99SBarry Smith }
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith #undef __FUNCT__
66847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddDM"
66947c6ae99SBarry Smith /*@C
670aa219208SBarry Smith     DMCompositeAddDM - adds a DM  vector to a DMComposite
67147c6ae99SBarry Smith 
67247c6ae99SBarry Smith     Collective on DMComposite
67347c6ae99SBarry Smith 
67447c6ae99SBarry Smith     Input Parameter:
67547c6ae99SBarry Smith +    dm - the packer object
67647c6ae99SBarry Smith -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
67747c6ae99SBarry Smith 
67847c6ae99SBarry Smith     Level: advanced
67947c6ae99SBarry Smith 
6800c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
6816eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
68247c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
68347c6ae99SBarry Smith 
68447c6ae99SBarry Smith @*/
6857087cfbeSBarry Smith PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
68647c6ae99SBarry Smith {
68747c6ae99SBarry Smith   PetscErrorCode         ierr;
68806ebdd98SJed Brown   PetscInt               n,nlocal;
68947c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
69006ebdd98SJed Brown   Vec                    global,local;
69147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmc->data;
69247c6ae99SBarry Smith 
69347c6ae99SBarry Smith   PetscFunctionBegin;
69447c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
69547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
69647c6ae99SBarry Smith   next = com->next;
697aa219208SBarry Smith   if (com->setup) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
69847c6ae99SBarry Smith 
69947c6ae99SBarry Smith   /* create new link */
70047c6ae99SBarry Smith   ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
70147c6ae99SBarry Smith   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
70247c6ae99SBarry Smith   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
70347c6ae99SBarry Smith   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
70447c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
70506ebdd98SJed Brown   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
70606ebdd98SJed Brown   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
70706ebdd98SJed Brown   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
70847c6ae99SBarry Smith   mine->n      = n;
70906ebdd98SJed Brown   mine->nlocal = nlocal;
71047c6ae99SBarry Smith   mine->dm     = dm;
71147c6ae99SBarry Smith   mine->type   = DMCOMPOSITE_DM;
71247c6ae99SBarry Smith   mine->next   = PETSC_NULL;
71347c6ae99SBarry Smith   com->n       += n;
71447c6ae99SBarry Smith 
71547c6ae99SBarry Smith   /* add to end of list */
71647c6ae99SBarry Smith   if (!next) {
71747c6ae99SBarry Smith     com->next = mine;
71847c6ae99SBarry Smith   } else {
71947c6ae99SBarry Smith     while (next->next) next = next->next;
72047c6ae99SBarry Smith     next->next = mine;
72147c6ae99SBarry Smith   }
72247c6ae99SBarry Smith   com->nDM++;
72347c6ae99SBarry Smith   com->nmine++;
72447c6ae99SBarry Smith   PetscFunctionReturn(0);
72547c6ae99SBarry Smith }
72647c6ae99SBarry Smith 
7277087cfbeSBarry Smith extern PetscErrorCode  VecView_MPI(Vec,PetscViewer);
72847c6ae99SBarry Smith EXTERN_C_BEGIN
72947c6ae99SBarry Smith #undef __FUNCT__
73047c6ae99SBarry Smith #define __FUNCT__ "VecView_DMComposite"
7317087cfbeSBarry Smith PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
73247c6ae99SBarry Smith {
73347c6ae99SBarry Smith   DM                     dm;
73447c6ae99SBarry Smith   PetscErrorCode         ierr;
73547c6ae99SBarry Smith   struct DMCompositeLink *next;
73647c6ae99SBarry Smith   PetscBool              isdraw;
737cef07954SSatish Balay   DM_Composite           *com;
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith   PetscFunctionBegin;
740*3c0c59f3SBarry Smith   ierr = PetscObjectQuery((PetscObject)gvec,"DM",(PetscObject*)&dm);CHKERRQ(ierr);
74147c6ae99SBarry Smith   if (!dm) SETERRQ(((PetscObject)gvec)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
74247c6ae99SBarry Smith   com = (DM_Composite*)dm->data;
74347c6ae99SBarry Smith   next = com->next;
74447c6ae99SBarry Smith 
74547c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
74647c6ae99SBarry Smith   if (!isdraw) {
74747c6ae99SBarry Smith     /* do I really want to call this? */
74847c6ae99SBarry Smith     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
74947c6ae99SBarry Smith   } else {
75047c6ae99SBarry Smith     PetscInt cnt = 0;
75147c6ae99SBarry Smith 
75247c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
75347c6ae99SBarry Smith     while (next) {
75447c6ae99SBarry Smith       if (next->type == DMCOMPOSITE_ARRAY) {
75547c6ae99SBarry Smith 	PetscScalar *array;
75647c6ae99SBarry Smith 	ierr  = DMCompositeGetAccess_Array(dm,next,gvec,&array);CHKERRQ(ierr);
75747c6ae99SBarry Smith 
75847c6ae99SBarry Smith 	/*skip it for now */
75947c6ae99SBarry Smith       } else if (next->type == DMCOMPOSITE_DM) {
76047c6ae99SBarry Smith 	Vec      vec;
76147c6ae99SBarry Smith         PetscInt bs;
76247c6ae99SBarry Smith 
76347c6ae99SBarry Smith 	ierr = DMCompositeGetAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
76447c6ae99SBarry Smith 	ierr = VecView(vec,viewer);CHKERRQ(ierr);
76547c6ae99SBarry Smith         ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
76647c6ae99SBarry Smith 	ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
76747c6ae99SBarry Smith         ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
76847c6ae99SBarry Smith         cnt += bs;
76947c6ae99SBarry Smith       } else {
77047c6ae99SBarry Smith 	SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
77147c6ae99SBarry Smith       }
77247c6ae99SBarry Smith       next = next->next;
77347c6ae99SBarry Smith     }
77447c6ae99SBarry Smith     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
77547c6ae99SBarry Smith   }
77647c6ae99SBarry Smith   PetscFunctionReturn(0);
77747c6ae99SBarry Smith }
77847c6ae99SBarry Smith EXTERN_C_END
77947c6ae99SBarry Smith 
78047c6ae99SBarry Smith 
78147c6ae99SBarry Smith #undef __FUNCT__
7820c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Composite"
7837087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
78447c6ae99SBarry Smith {
78547c6ae99SBarry Smith   PetscErrorCode         ierr;
78647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
78747c6ae99SBarry Smith 
78847c6ae99SBarry Smith   PetscFunctionBegin;
78947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
79047c6ae99SBarry Smith   if (!com->setup) {
791d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
79247c6ae99SBarry Smith   }
79347c6ae99SBarry Smith   ierr = VecCreateMPI(((PetscObject)dm)->comm,com->n,com->N,gvec);CHKERRQ(ierr);
794*3c0c59f3SBarry Smith   ierr = PetscObjectCompose((PetscObject)*gvec,"DM",(PetscObject)dm);CHKERRQ(ierr);
79547c6ae99SBarry Smith   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void(*)(void))VecView_DMComposite);CHKERRQ(ierr);
79647c6ae99SBarry Smith   PetscFunctionReturn(0);
79747c6ae99SBarry Smith }
79847c6ae99SBarry Smith 
79947c6ae99SBarry Smith #undef __FUNCT__
8000c010503SBarry Smith #define __FUNCT__ "DMCreateLocalVector_Composite"
8017087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
80247c6ae99SBarry Smith {
80347c6ae99SBarry Smith   PetscErrorCode         ierr;
80447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
80547c6ae99SBarry Smith 
80647c6ae99SBarry Smith   PetscFunctionBegin;
80747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
80847c6ae99SBarry Smith   if (!com->setup) {
809d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
81047c6ae99SBarry Smith   }
81147c6ae99SBarry Smith   ierr = VecCreateSeq(((PetscObject)dm)->comm,com->nghost,lvec);CHKERRQ(ierr);
812*3c0c59f3SBarry Smith   ierr = PetscObjectCompose((PetscObject)*lvec,"DM",(PetscObject)dm);CHKERRQ(ierr);
81347c6ae99SBarry Smith   PetscFunctionReturn(0);
81447c6ae99SBarry Smith }
81547c6ae99SBarry Smith 
81647c6ae99SBarry Smith #undef __FUNCT__
8176eb61c8cSJed Brown #define __FUNCT__ "DMCompositeGetISLocalToGlobalMappings"
81847c6ae99SBarry Smith /*@C
8196eb61c8cSJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM/array in the DMComposite, maps to the composite global space
82047c6ae99SBarry Smith 
82106ebdd98SJed Brown     Collective on DM
82247c6ae99SBarry Smith 
82347c6ae99SBarry Smith     Input Parameter:
82447c6ae99SBarry Smith .    dm - the packer object
82547c6ae99SBarry Smith 
82647c6ae99SBarry Smith     Output Parameters:
8276eb61c8cSJed Brown .    ltogs - the individual mappings for each packed vector/array. Note that this includes
828aa219208SBarry Smith            all the ghost points that individual ghosted DMDA's may have. Also each process has an
8296eb61c8cSJed Brown            mapping for EACH redundant array (not just the local redundant arrays).
83047c6ae99SBarry Smith 
83147c6ae99SBarry Smith     Level: advanced
83247c6ae99SBarry Smith 
83347c6ae99SBarry Smith     Notes:
8346eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
83547c6ae99SBarry Smith 
8360c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
83747c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
83847c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
83947c6ae99SBarry Smith 
84047c6ae99SBarry Smith @*/
8417087cfbeSBarry Smith PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
84247c6ae99SBarry Smith {
84347c6ae99SBarry Smith   PetscErrorCode         ierr;
84447c6ae99SBarry Smith   PetscInt               i,*idx,n,cnt;
84547c6ae99SBarry Smith   struct DMCompositeLink *next;
84647c6ae99SBarry Smith   PetscMPIInt            rank;
84747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
84847c6ae99SBarry Smith 
84947c6ae99SBarry Smith   PetscFunctionBegin;
85047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
851acc1e270SJed Brown   ierr = PetscMalloc((com->nDM+com->nredundant)*sizeof(ISLocalToGlobalMapping),ltogs);CHKERRQ(ierr);
85247c6ae99SBarry Smith   next = com->next;
85347c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
85447c6ae99SBarry Smith 
85547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
85647c6ae99SBarry Smith   cnt = 0;
85747c6ae99SBarry Smith   while (next) {
85847c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
85906ebdd98SJed Brown       ierr = PetscMalloc(next->nlocal*sizeof(PetscInt),&idx);CHKERRQ(ierr);
86047c6ae99SBarry Smith       if (rank == next->rank) {
86106ebdd98SJed Brown         for (i=0; i<next->nlocal; i++) idx[i] = next->grstart + i;
86247c6ae99SBarry Smith       }
86306ebdd98SJed Brown       ierr = MPI_Bcast(idx,next->nlocal,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
86406ebdd98SJed Brown       ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,next->nlocal,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
86547c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
8666eb61c8cSJed Brown       ISLocalToGlobalMapping ltog;
8676eb61c8cSJed Brown       PetscMPIInt            size;
86886994e45SJed Brown       const PetscInt         *suboff,*indices;
8696eb61c8cSJed Brown       Vec                    global;
87047c6ae99SBarry Smith 
8716eb61c8cSJed Brown       /* Get sub-DM global indices for each local dof */
8721411c6eeSJed Brown       ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
8736eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
87486994e45SJed Brown       ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
87547c6ae99SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
87647c6ae99SBarry Smith 
8776eb61c8cSJed Brown       /* Get the offsets for the sub-DM global vector */
8786eb61c8cSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
8796eb61c8cSJed Brown       ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
8806eb61c8cSJed Brown       ierr = MPI_Comm_size(((PetscObject)global)->comm,&size);CHKERRQ(ierr);
8816eb61c8cSJed Brown 
8826eb61c8cSJed Brown       /* Shift the sub-DM definition of the global space to the composite global space */
8836eb61c8cSJed Brown       for (i=0; i<n; i++) {
88486994e45SJed Brown         PetscInt subi = indices[i],lo = 0,hi = size,t;
8856eb61c8cSJed Brown         /* Binary search to find which rank owns subi */
8866eb61c8cSJed Brown         while (hi-lo > 1) {
8876eb61c8cSJed Brown           t = lo + (hi-lo)/2;
8886eb61c8cSJed Brown           if (suboff[t] > subi) hi = t;
8896eb61c8cSJed Brown           else                  lo = t;
8906eb61c8cSJed Brown         }
8916eb61c8cSJed Brown         idx[i] = subi - suboff[lo] + next->grstarts[lo];
8926eb61c8cSJed Brown       }
89386994e45SJed Brown       ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
8946eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
8956eb61c8cSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
8966eb61c8cSJed Brown     } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
89747c6ae99SBarry Smith     next = next->next;
89847c6ae99SBarry Smith     cnt++;
89947c6ae99SBarry Smith   }
90047c6ae99SBarry Smith   PetscFunctionReturn(0);
90147c6ae99SBarry Smith }
90247c6ae99SBarry Smith 
90347c6ae99SBarry Smith #undef __FUNCT__
90487c85e80SJed Brown #define __FUNCT__ "DMCompositeGetLocalISs"
90587c85e80SJed Brown /*@C
90687c85e80SJed Brown    DMCompositeGetLocalISs - Gets index sets for each DM/array component of a composite local vector
90787c85e80SJed Brown 
90887c85e80SJed Brown    Not Collective
90987c85e80SJed Brown 
91087c85e80SJed Brown    Input Arguments:
91187c85e80SJed Brown . dm - composite DM
91287c85e80SJed Brown 
91387c85e80SJed Brown    Output Arguments:
91487c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
91587c85e80SJed Brown 
91687c85e80SJed Brown    Level: intermediate
91787c85e80SJed Brown 
91887c85e80SJed Brown    Notes:
91987c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
92087c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
92187c85e80SJed Brown    scatter to a composite local vector.
92287c85e80SJed Brown 
92387c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
92487c85e80SJed Brown 
92587c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
92687c85e80SJed Brown 
92787c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
92887c85e80SJed Brown 
92987c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
93087c85e80SJed Brown @*/
9317087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
93287c85e80SJed Brown {
93387c85e80SJed Brown   PetscErrorCode         ierr;
93487c85e80SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
93587c85e80SJed Brown   struct DMCompositeLink *link;
93687c85e80SJed Brown   PetscInt cnt,start;
93787c85e80SJed Brown 
93887c85e80SJed Brown   PetscFunctionBegin;
93987c85e80SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
94087c85e80SJed Brown   PetscValidPointer(is,2);
94187c85e80SJed Brown   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
94206ebdd98SJed Brown   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
94306ebdd98SJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
944520db06cSJed Brown     if (link->type == DMCOMPOSITE_DM) {
945520db06cSJed Brown       PetscInt bs;
9461411c6eeSJed Brown       ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
947520db06cSJed Brown       ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
948520db06cSJed Brown     }
94987c85e80SJed Brown   }
95087c85e80SJed Brown   PetscFunctionReturn(0);
95187c85e80SJed Brown }
95287c85e80SJed Brown 
95387c85e80SJed Brown #undef __FUNCT__
95447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetGlobalISs"
95547c6ae99SBarry Smith /*@C
95647c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
95747c6ae99SBarry Smith 
95847c6ae99SBarry Smith     Collective on DMComposite
95947c6ae99SBarry Smith 
96047c6ae99SBarry Smith     Input Parameter:
96147c6ae99SBarry Smith .    dm - the packer object
96247c6ae99SBarry Smith 
96347c6ae99SBarry Smith     Output Parameters:
96447c6ae99SBarry Smith .    is - the array of index sets
96547c6ae99SBarry Smith 
96647c6ae99SBarry Smith     Level: advanced
96747c6ae99SBarry Smith 
96847c6ae99SBarry Smith     Notes:
96947c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
97047c6ae99SBarry Smith 
97147c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
97247c6ae99SBarry Smith 
9736eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
9746eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
9756eb61c8cSJed Brown        indices.
97647c6ae99SBarry Smith 
9770c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
97847c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
97947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
98047c6ae99SBarry Smith 
98147c6ae99SBarry Smith @*/
9826eb61c8cSJed Brown 
9837087cfbeSBarry Smith PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
98447c6ae99SBarry Smith {
98547c6ae99SBarry Smith   PetscErrorCode         ierr;
98647c6ae99SBarry Smith   PetscInt               cnt = 0,*idx,i;
98747c6ae99SBarry Smith   struct DMCompositeLink *next;
98847c6ae99SBarry Smith   PetscMPIInt            rank;
98947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
99047c6ae99SBarry Smith 
99147c6ae99SBarry Smith   PetscFunctionBegin;
99247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
99306ebdd98SJed Brown   ierr = PetscMalloc((com->nDM+com->nredundant)*sizeof(IS),is);CHKERRQ(ierr);
99447c6ae99SBarry Smith   next = com->next;
99547c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
99647c6ae99SBarry Smith 
99747c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
99847c6ae99SBarry Smith   while (next) {
99947c6ae99SBarry Smith     ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
100047c6ae99SBarry Smith     for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
100147c6ae99SBarry Smith     ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
100247c6ae99SBarry Smith     cnt++;
100347c6ae99SBarry Smith     next = next->next;
100447c6ae99SBarry Smith   }
100547c6ae99SBarry Smith   PetscFunctionReturn(0);
100647c6ae99SBarry Smith }
100747c6ae99SBarry Smith 
100847c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
100947c6ae99SBarry Smith #undef __FUNCT__
101047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_Array"
101147c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
101247c6ae99SBarry Smith {
101347c6ae99SBarry Smith   PetscErrorCode ierr;
101447c6ae99SBarry Smith   PetscFunctionBegin;
101547c6ae99SBarry Smith   if (array) {
1016c72eef85SJed Brown     ierr = PetscMalloc(mine->nlocal*sizeof(PetscScalar),array);CHKERRQ(ierr);
101747c6ae99SBarry Smith   }
101847c6ae99SBarry Smith   PetscFunctionReturn(0);
101947c6ae99SBarry Smith }
102047c6ae99SBarry Smith 
102147c6ae99SBarry Smith #undef __FUNCT__
102247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_DM"
102347c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
102447c6ae99SBarry Smith {
102547c6ae99SBarry Smith   PetscErrorCode ierr;
102647c6ae99SBarry Smith   PetscFunctionBegin;
102747c6ae99SBarry Smith   if (local) {
102847c6ae99SBarry Smith     ierr = DMGetLocalVector(mine->dm,local);CHKERRQ(ierr);
102947c6ae99SBarry Smith   }
103047c6ae99SBarry Smith   PetscFunctionReturn(0);
103147c6ae99SBarry Smith }
103247c6ae99SBarry Smith 
103347c6ae99SBarry Smith #undef __FUNCT__
103447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_Array"
103547c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
103647c6ae99SBarry Smith {
103747c6ae99SBarry Smith   PetscErrorCode ierr;
103847c6ae99SBarry Smith   PetscFunctionBegin;
103947c6ae99SBarry Smith   if (array) {
104047c6ae99SBarry Smith     ierr = PetscFree(*array);CHKERRQ(ierr);
104147c6ae99SBarry Smith   }
104247c6ae99SBarry Smith   PetscFunctionReturn(0);
104347c6ae99SBarry Smith }
104447c6ae99SBarry Smith 
104547c6ae99SBarry Smith #undef __FUNCT__
104647c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_DM"
104747c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
104847c6ae99SBarry Smith {
104947c6ae99SBarry Smith   PetscErrorCode ierr;
105047c6ae99SBarry Smith   PetscFunctionBegin;
105147c6ae99SBarry Smith   if (local) {
105247c6ae99SBarry Smith     ierr = DMRestoreLocalVector(mine->dm,local);CHKERRQ(ierr);
105347c6ae99SBarry Smith   }
105447c6ae99SBarry Smith   PetscFunctionReturn(0);
105547c6ae99SBarry Smith }
105647c6ae99SBarry Smith 
105747c6ae99SBarry Smith #undef __FUNCT__
105847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors"
105947c6ae99SBarry Smith /*@C
106047c6ae99SBarry Smith     DMCompositeGetLocalVectors - Gets local vectors and arrays for each part of a DMComposite.
106147c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
106247c6ae99SBarry Smith 
106347c6ae99SBarry Smith     Not Collective
106447c6ae99SBarry Smith 
106547c6ae99SBarry Smith     Input Parameter:
106647c6ae99SBarry Smith .    dm - the packer object
106747c6ae99SBarry Smith 
106847c6ae99SBarry Smith     Output Parameter:
106947c6ae99SBarry Smith .    ... - the individual sequential objects (arrays or vectors)
107047c6ae99SBarry Smith 
107147c6ae99SBarry Smith     Level: advanced
107247c6ae99SBarry Smith 
10730c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
10746eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
107547c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
107647c6ae99SBarry Smith 
107747c6ae99SBarry Smith @*/
10787087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
107947c6ae99SBarry Smith {
108047c6ae99SBarry Smith   va_list                Argp;
108147c6ae99SBarry Smith   PetscErrorCode         ierr;
108247c6ae99SBarry Smith   struct DMCompositeLink *next;
108347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
108447c6ae99SBarry Smith 
108547c6ae99SBarry Smith   PetscFunctionBegin;
108647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
108747c6ae99SBarry Smith   next = com->next;
108847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
108947c6ae99SBarry Smith   va_start(Argp,dm);
109047c6ae99SBarry Smith   while (next) {
109147c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
109247c6ae99SBarry Smith       PetscScalar **array;
109347c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
109447c6ae99SBarry Smith       ierr = DMCompositeGetLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
109547c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
109647c6ae99SBarry Smith       Vec *vec;
109747c6ae99SBarry Smith       vec = va_arg(Argp, Vec*);
109847c6ae99SBarry Smith       ierr = DMCompositeGetLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
109947c6ae99SBarry Smith     } else {
110047c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
110147c6ae99SBarry Smith     }
110247c6ae99SBarry Smith     next = next->next;
110347c6ae99SBarry Smith   }
110447c6ae99SBarry Smith   va_end(Argp);
110547c6ae99SBarry Smith   PetscFunctionReturn(0);
110647c6ae99SBarry Smith }
110747c6ae99SBarry Smith 
110847c6ae99SBarry Smith #undef __FUNCT__
110947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors"
111047c6ae99SBarry Smith /*@C
111147c6ae99SBarry Smith     DMCompositeRestoreLocalVectors - Restores local vectors and arrays for each part of a DMComposite.
111247c6ae99SBarry Smith 
111347c6ae99SBarry Smith     Not Collective
111447c6ae99SBarry Smith 
111547c6ae99SBarry Smith     Input Parameter:
111647c6ae99SBarry Smith .    dm - the packer object
111747c6ae99SBarry Smith 
111847c6ae99SBarry Smith     Output Parameter:
111947c6ae99SBarry Smith .    ... - the individual sequential objects (arrays or vectors)
112047c6ae99SBarry Smith 
112147c6ae99SBarry Smith     Level: advanced
112247c6ae99SBarry Smith 
11230c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
11246eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
112547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
112647c6ae99SBarry Smith 
112747c6ae99SBarry Smith @*/
11287087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
112947c6ae99SBarry Smith {
113047c6ae99SBarry Smith   va_list                Argp;
113147c6ae99SBarry Smith   PetscErrorCode         ierr;
113247c6ae99SBarry Smith   struct DMCompositeLink *next;
113347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
113447c6ae99SBarry Smith 
113547c6ae99SBarry Smith   PetscFunctionBegin;
113647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
113747c6ae99SBarry Smith   next = com->next;
113847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
113947c6ae99SBarry Smith   va_start(Argp,dm);
114047c6ae99SBarry Smith   while (next) {
114147c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
114247c6ae99SBarry Smith       PetscScalar **array;
114347c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
114447c6ae99SBarry Smith       ierr = DMCompositeRestoreLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
114547c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
114647c6ae99SBarry Smith       Vec *vec;
114747c6ae99SBarry Smith       vec = va_arg(Argp, Vec*);
114847c6ae99SBarry Smith       ierr = DMCompositeRestoreLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
114947c6ae99SBarry Smith     } else {
115047c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
115147c6ae99SBarry Smith     }
115247c6ae99SBarry Smith     next = next->next;
115347c6ae99SBarry Smith   }
115447c6ae99SBarry Smith   va_end(Argp);
115547c6ae99SBarry Smith   PetscFunctionReturn(0);
115647c6ae99SBarry Smith }
115747c6ae99SBarry Smith 
115847c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
115947c6ae99SBarry Smith #undef __FUNCT__
116047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_Array"
116147c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_Array(DM dm,struct DMCompositeLink *mine,PetscInt *n)
116247c6ae99SBarry Smith {
116347c6ae99SBarry Smith   PetscFunctionBegin;
116406ebdd98SJed Brown   if (n) *n = mine->nlocal;
116547c6ae99SBarry Smith   PetscFunctionReturn(0);
116647c6ae99SBarry Smith }
116747c6ae99SBarry Smith 
116847c6ae99SBarry Smith #undef __FUNCT__
116947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_DM"
117047c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_DM(DM dmi,struct DMCompositeLink *mine,DM *dm)
117147c6ae99SBarry Smith {
117247c6ae99SBarry Smith   PetscFunctionBegin;
117347c6ae99SBarry Smith   if (dm) *dm = mine->dm;
117447c6ae99SBarry Smith   PetscFunctionReturn(0);
117547c6ae99SBarry Smith }
117647c6ae99SBarry Smith 
117747c6ae99SBarry Smith #undef __FUNCT__
117847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries"
117947c6ae99SBarry Smith /*@C
1180aa219208SBarry Smith     DMCompositeGetEntries - Gets the DM, redundant size, etc for each entry in a DMComposite.
118147c6ae99SBarry Smith 
118247c6ae99SBarry Smith     Not Collective
118347c6ae99SBarry Smith 
118447c6ae99SBarry Smith     Input Parameter:
118547c6ae99SBarry Smith .    dm - the packer object
118647c6ae99SBarry Smith 
118747c6ae99SBarry Smith     Output Parameter:
1188aa219208SBarry Smith .    ... - the individual entries, DMs or integer sizes)
118947c6ae99SBarry Smith 
119047c6ae99SBarry Smith     Level: advanced
119147c6ae99SBarry Smith 
11920c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
11936eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
119447c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
119547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
119647c6ae99SBarry Smith 
119747c6ae99SBarry Smith @*/
11987087cfbeSBarry Smith PetscErrorCode  DMCompositeGetEntries(DM dm,...)
119947c6ae99SBarry Smith {
120047c6ae99SBarry Smith   va_list                Argp;
120147c6ae99SBarry Smith   PetscErrorCode         ierr;
120247c6ae99SBarry Smith   struct DMCompositeLink *next;
120347c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
120447c6ae99SBarry Smith 
120547c6ae99SBarry Smith   PetscFunctionBegin;
120647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
120747c6ae99SBarry Smith   next = com->next;
120847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
120947c6ae99SBarry Smith   va_start(Argp,dm);
121047c6ae99SBarry Smith   while (next) {
121147c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
121247c6ae99SBarry Smith       PetscInt *n;
121347c6ae99SBarry Smith       n = va_arg(Argp, PetscInt*);
121447c6ae99SBarry Smith       ierr = DMCompositeGetEntries_Array(dm,next,n);CHKERRQ(ierr);
121547c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
121647c6ae99SBarry Smith       DM *dmn;
121747c6ae99SBarry Smith       dmn = va_arg(Argp, DM*);
121847c6ae99SBarry Smith       ierr = DMCompositeGetEntries_DM(dm,next,dmn);CHKERRQ(ierr);
121947c6ae99SBarry Smith     } else {
122047c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
122147c6ae99SBarry Smith     }
122247c6ae99SBarry Smith     next = next->next;
122347c6ae99SBarry Smith   }
122447c6ae99SBarry Smith   va_end(Argp);
122547c6ae99SBarry Smith   PetscFunctionReturn(0);
122647c6ae99SBarry Smith }
122747c6ae99SBarry Smith 
122847c6ae99SBarry Smith #undef __FUNCT__
12290c010503SBarry Smith #define __FUNCT__ "DMRefine_Composite"
12307087cfbeSBarry Smith PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
123147c6ae99SBarry Smith {
123247c6ae99SBarry Smith   PetscErrorCode         ierr;
123347c6ae99SBarry Smith   struct DMCompositeLink *next;
123447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmi->data;
123547c6ae99SBarry Smith   DM                     dm;
123647c6ae99SBarry Smith 
123747c6ae99SBarry Smith   PetscFunctionBegin;
123847c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
123947c6ae99SBarry Smith   next = com->next;
124047c6ae99SBarry Smith   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
124147c6ae99SBarry Smith 
124247c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
124347c6ae99SBarry Smith   while (next) {
124447c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
124506ebdd98SJed Brown       ierr = DMCompositeAddArray(*fine,next->rank,next->nlocal);CHKERRQ(ierr);
124647c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
124747c6ae99SBarry Smith       ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
124847c6ae99SBarry Smith       ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
124947c6ae99SBarry Smith       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
125047c6ae99SBarry Smith     } else {
125147c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
125247c6ae99SBarry Smith     }
125347c6ae99SBarry Smith     next = next->next;
125447c6ae99SBarry Smith   }
125547c6ae99SBarry Smith   PetscFunctionReturn(0);
125647c6ae99SBarry Smith }
125747c6ae99SBarry Smith 
1258c6db04a5SJed Brown #include <petscmat.h>
125947c6ae99SBarry Smith 
126047c6ae99SBarry Smith struct MatPackLink {
126147c6ae99SBarry Smith   Mat                A;
126247c6ae99SBarry Smith   struct MatPackLink *next;
126347c6ae99SBarry Smith };
126447c6ae99SBarry Smith 
126547c6ae99SBarry Smith struct MatPack {
126647c6ae99SBarry Smith   DM                 right,left;
126747c6ae99SBarry Smith   struct MatPackLink *next;
126847c6ae99SBarry Smith };
126947c6ae99SBarry Smith 
127047c6ae99SBarry Smith #undef __FUNCT__
127147c6ae99SBarry Smith #define __FUNCT__ "MatMultBoth_Shell_Pack"
127247c6ae99SBarry Smith PetscErrorCode MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscBool  add)
127347c6ae99SBarry Smith {
127447c6ae99SBarry Smith   struct MatPack         *mpack;
127547c6ae99SBarry Smith   struct DMCompositeLink *xnext,*ynext;
127647c6ae99SBarry Smith   struct MatPackLink     *anext;
127747c6ae99SBarry Smith   PetscScalar            *xarray,*yarray;
127847c6ae99SBarry Smith   PetscErrorCode         ierr;
127947c6ae99SBarry Smith   PetscInt               i;
128047c6ae99SBarry Smith   Vec                    xglobal,yglobal;
128147c6ae99SBarry Smith   PetscMPIInt            rank;
128247c6ae99SBarry Smith   DM_Composite           *comright;
128347c6ae99SBarry Smith   DM_Composite           *comleft;
128447c6ae99SBarry Smith 
128547c6ae99SBarry Smith   PetscFunctionBegin;
128647c6ae99SBarry Smith   ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
128747c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
128847c6ae99SBarry Smith   comright = (DM_Composite*)mpack->right->data;
128947c6ae99SBarry Smith   comleft = (DM_Composite*)mpack->left->data;
129047c6ae99SBarry Smith   xnext = comright->next;
129147c6ae99SBarry Smith   ynext = comleft->next;
129247c6ae99SBarry Smith   anext = mpack->next;
129347c6ae99SBarry Smith 
129447c6ae99SBarry Smith   while (xnext) {
129547c6ae99SBarry Smith     if (xnext->type == DMCOMPOSITE_ARRAY) {
129647c6ae99SBarry Smith       if (rank == xnext->rank) {
129747c6ae99SBarry Smith         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
129847c6ae99SBarry Smith         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
129947c6ae99SBarry Smith         if (add) {
130047c6ae99SBarry Smith           for (i=0; i<xnext->n; i++) {
130147c6ae99SBarry Smith             yarray[ynext->rstart+i] += xarray[xnext->rstart+i];
130247c6ae99SBarry Smith           }
130347c6ae99SBarry Smith         } else {
130447c6ae99SBarry Smith           ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
130547c6ae99SBarry Smith         }
130647c6ae99SBarry Smith         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
130747c6ae99SBarry Smith         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
130847c6ae99SBarry Smith       }
130947c6ae99SBarry Smith     } else if (xnext->type == DMCOMPOSITE_DM) {
131047c6ae99SBarry Smith       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
131147c6ae99SBarry Smith       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
131247c6ae99SBarry Smith       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
131347c6ae99SBarry Smith       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
131447c6ae99SBarry Smith       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
131547c6ae99SBarry Smith       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
131647c6ae99SBarry Smith       if (add) {
131747c6ae99SBarry Smith         ierr  = MatMultAdd(anext->A,xglobal,yglobal,yglobal);CHKERRQ(ierr);
131847c6ae99SBarry Smith       } else {
131947c6ae99SBarry Smith         ierr  = MatMult(anext->A,xglobal,yglobal);CHKERRQ(ierr);
132047c6ae99SBarry Smith       }
132147c6ae99SBarry Smith       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
132247c6ae99SBarry Smith       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
132347c6ae99SBarry Smith       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
132447c6ae99SBarry Smith       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
132547c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
132647c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
132747c6ae99SBarry Smith       anext = anext->next;
132847c6ae99SBarry Smith     } else {
132947c6ae99SBarry Smith       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
133047c6ae99SBarry Smith     }
133147c6ae99SBarry Smith     xnext = xnext->next;
133247c6ae99SBarry Smith     ynext = ynext->next;
133347c6ae99SBarry Smith   }
133447c6ae99SBarry Smith   PetscFunctionReturn(0);
133547c6ae99SBarry Smith }
133647c6ae99SBarry Smith 
133747c6ae99SBarry Smith #undef __FUNCT__
133847c6ae99SBarry Smith #define __FUNCT__ "MatMultAdd_Shell_Pack"
133947c6ae99SBarry Smith PetscErrorCode MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z)
134047c6ae99SBarry Smith {
134147c6ae99SBarry Smith   PetscErrorCode ierr;
134247c6ae99SBarry Smith   PetscFunctionBegin;
134347c6ae99SBarry Smith   if (z != y) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Handles y == z only");
134447c6ae99SBarry Smith   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
134547c6ae99SBarry Smith   PetscFunctionReturn(0);
134647c6ae99SBarry Smith }
134747c6ae99SBarry Smith 
134847c6ae99SBarry Smith #undef __FUNCT__
134947c6ae99SBarry Smith #define __FUNCT__ "MatMult_Shell_Pack"
135047c6ae99SBarry Smith PetscErrorCode MatMult_Shell_Pack(Mat A,Vec x,Vec y)
135147c6ae99SBarry Smith {
135247c6ae99SBarry Smith   PetscErrorCode ierr;
135347c6ae99SBarry Smith   PetscFunctionBegin;
135447c6ae99SBarry Smith   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
135547c6ae99SBarry Smith   PetscFunctionReturn(0);
135647c6ae99SBarry Smith }
135747c6ae99SBarry Smith 
135847c6ae99SBarry Smith #undef __FUNCT__
135947c6ae99SBarry Smith #define __FUNCT__ "MatMultTranspose_Shell_Pack"
136047c6ae99SBarry Smith PetscErrorCode MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y)
136147c6ae99SBarry Smith {
136247c6ae99SBarry Smith   struct MatPack         *mpack;
136347c6ae99SBarry Smith   struct DMCompositeLink *xnext,*ynext;
136447c6ae99SBarry Smith   struct MatPackLink     *anext;
136547c6ae99SBarry Smith   PetscScalar            *xarray,*yarray;
136647c6ae99SBarry Smith   PetscErrorCode         ierr;
136747c6ae99SBarry Smith   Vec                    xglobal,yglobal;
136847c6ae99SBarry Smith   PetscMPIInt            rank;
136947c6ae99SBarry Smith   DM_Composite           *comright;
137047c6ae99SBarry Smith   DM_Composite           *comleft;
137147c6ae99SBarry Smith 
137247c6ae99SBarry Smith   PetscFunctionBegin;
137347c6ae99SBarry Smith   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
137447c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
137547c6ae99SBarry Smith   comright = (DM_Composite*)mpack->right->data;
137647c6ae99SBarry Smith   comleft = (DM_Composite*)mpack->left->data;
137747c6ae99SBarry Smith   ynext = comright->next;
137847c6ae99SBarry Smith   xnext = comleft->next;
137947c6ae99SBarry Smith   anext = mpack->next;
138047c6ae99SBarry Smith 
138147c6ae99SBarry Smith   while (xnext) {
138247c6ae99SBarry Smith     if (xnext->type == DMCOMPOSITE_ARRAY) {
138347c6ae99SBarry Smith       if (rank == ynext->rank) {
138447c6ae99SBarry Smith         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
138547c6ae99SBarry Smith         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
138647c6ae99SBarry Smith         ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
138747c6ae99SBarry Smith         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
138847c6ae99SBarry Smith         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
138947c6ae99SBarry Smith       }
139047c6ae99SBarry Smith     } else if (xnext->type == DMCOMPOSITE_DM) {
139147c6ae99SBarry Smith       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
139247c6ae99SBarry Smith       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
139347c6ae99SBarry Smith       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
139447c6ae99SBarry Smith       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
139547c6ae99SBarry Smith       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
139647c6ae99SBarry Smith       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
139747c6ae99SBarry Smith       ierr  = MatMultTranspose(anext->A,xglobal,yglobal);CHKERRQ(ierr);
139847c6ae99SBarry Smith       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
139947c6ae99SBarry Smith       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
140047c6ae99SBarry Smith       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
140147c6ae99SBarry Smith       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
140247c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
140347c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
140447c6ae99SBarry Smith       anext = anext->next;
140547c6ae99SBarry Smith     } else {
140647c6ae99SBarry Smith       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
140747c6ae99SBarry Smith     }
140847c6ae99SBarry Smith     xnext = xnext->next;
140947c6ae99SBarry Smith     ynext = ynext->next;
141047c6ae99SBarry Smith   }
141147c6ae99SBarry Smith   PetscFunctionReturn(0);
141247c6ae99SBarry Smith }
141347c6ae99SBarry Smith 
141447c6ae99SBarry Smith #undef __FUNCT__
141547c6ae99SBarry Smith #define __FUNCT__ "MatDestroy_Shell_Pack"
141647c6ae99SBarry Smith PetscErrorCode MatDestroy_Shell_Pack(Mat A)
141747c6ae99SBarry Smith {
141847c6ae99SBarry Smith   struct MatPack     *mpack;
141947c6ae99SBarry Smith   struct MatPackLink *anext,*oldanext;
142047c6ae99SBarry Smith   PetscErrorCode     ierr;
142147c6ae99SBarry Smith 
142247c6ae99SBarry Smith   PetscFunctionBegin;
142347c6ae99SBarry Smith   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
142447c6ae99SBarry Smith   anext = mpack->next;
142547c6ae99SBarry Smith 
142647c6ae99SBarry Smith   while (anext) {
1427fcfd50ebSBarry Smith     ierr     = MatDestroy(&anext->A);CHKERRQ(ierr);
142847c6ae99SBarry Smith     oldanext = anext;
142947c6ae99SBarry Smith     anext    = anext->next;
143047c6ae99SBarry Smith     ierr     = PetscFree(oldanext);CHKERRQ(ierr);
143147c6ae99SBarry Smith   }
143247c6ae99SBarry Smith   ierr = PetscFree(mpack);CHKERRQ(ierr);
143347c6ae99SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
143447c6ae99SBarry Smith   PetscFunctionReturn(0);
143547c6ae99SBarry Smith }
143647c6ae99SBarry Smith 
143747c6ae99SBarry Smith #undef __FUNCT__
14380c010503SBarry Smith #define __FUNCT__ "DMGetInterpolation_Composite"
14397087cfbeSBarry Smith PetscErrorCode  DMGetInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
144047c6ae99SBarry Smith {
144147c6ae99SBarry Smith   PetscErrorCode         ierr;
144247c6ae99SBarry Smith   PetscInt               m,n,M,N;
144347c6ae99SBarry Smith   struct DMCompositeLink *nextc;
144447c6ae99SBarry Smith   struct DMCompositeLink *nextf;
144547c6ae99SBarry Smith   struct MatPackLink     *nextmat,*pnextmat = 0;
144647c6ae99SBarry Smith   struct MatPack         *mpack;
144747c6ae99SBarry Smith   Vec                    gcoarse,gfine;
144847c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
144947c6ae99SBarry Smith   DM_Composite           *comfine = (DM_Composite*)fine->data;
145047c6ae99SBarry Smith 
145147c6ae99SBarry Smith   PetscFunctionBegin;
145247c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
145347c6ae99SBarry Smith   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
145447c6ae99SBarry Smith   nextc = comcoarse->next;
145547c6ae99SBarry Smith   nextf = comfine->next;
145647c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14570c010503SBarry Smith   ierr = DMCreateGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
14580c010503SBarry Smith   ierr = DMCreateGlobalVector(fine,&gfine);CHKERRQ(ierr);
145947c6ae99SBarry Smith   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
146047c6ae99SBarry Smith   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
146147c6ae99SBarry Smith   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
146247c6ae99SBarry Smith   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
1463fcfd50ebSBarry Smith   ierr = VecDestroy(&gcoarse);CHKERRQ(ierr);
1464fcfd50ebSBarry Smith   ierr = VecDestroy(&gfine);CHKERRQ(ierr);
146547c6ae99SBarry Smith 
146647c6ae99SBarry Smith   ierr         = PetscNew(struct MatPack,&mpack);CHKERRQ(ierr);
146747c6ae99SBarry Smith   mpack->right = coarse;
146847c6ae99SBarry Smith   mpack->left  = fine;
146947c6ae99SBarry Smith   ierr  = MatCreate(((PetscObject)fine)->comm,A);CHKERRQ(ierr);
147047c6ae99SBarry Smith   ierr  = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
147147c6ae99SBarry Smith   ierr  = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
147247c6ae99SBarry Smith   ierr  = MatShellSetContext(*A,mpack);CHKERRQ(ierr);
147347c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);CHKERRQ(ierr);
147447c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);CHKERRQ(ierr);
147547c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);CHKERRQ(ierr);
147647c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);CHKERRQ(ierr);
147747c6ae99SBarry Smith 
147847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
147947c6ae99SBarry Smith   while (nextc) {
148047c6ae99SBarry Smith     if (nextc->type != nextf->type) SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Two DM have different layout");
148147c6ae99SBarry Smith 
148247c6ae99SBarry Smith     if (nextc->type == DMCOMPOSITE_ARRAY) {
148347c6ae99SBarry Smith       ;
148447c6ae99SBarry Smith     } else if (nextc->type == DMCOMPOSITE_DM) {
148547c6ae99SBarry Smith       ierr          = PetscNew(struct MatPackLink,&nextmat);CHKERRQ(ierr);
148647c6ae99SBarry Smith       nextmat->next = 0;
148747c6ae99SBarry Smith       if (pnextmat) {
148847c6ae99SBarry Smith         pnextmat->next = nextmat;
148947c6ae99SBarry Smith         pnextmat       = nextmat;
149047c6ae99SBarry Smith       } else {
149147c6ae99SBarry Smith         pnextmat    = nextmat;
149247c6ae99SBarry Smith         mpack->next = nextmat;
149347c6ae99SBarry Smith       }
149447c6ae99SBarry Smith       ierr = DMGetInterpolation(nextc->dm,nextf->dm,&nextmat->A,PETSC_NULL);CHKERRQ(ierr);
149547c6ae99SBarry Smith     } else {
149647c6ae99SBarry Smith       SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
149747c6ae99SBarry Smith     }
149847c6ae99SBarry Smith     nextc = nextc->next;
149947c6ae99SBarry Smith     nextf = nextf->next;
150047c6ae99SBarry Smith   }
150147c6ae99SBarry Smith   PetscFunctionReturn(0);
150247c6ae99SBarry Smith }
150347c6ae99SBarry Smith 
150447c6ae99SBarry Smith #undef __FUNCT__
15051411c6eeSJed Brown #define __FUNCT__ "DMCreateLocalToGlobalMapping_Composite"
15061411c6eeSJed Brown static PetscErrorCode DMCreateLocalToGlobalMapping_Composite(DM dm)
15071411c6eeSJed Brown {
15081411c6eeSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
15091411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1510a03f9cedSJed Brown   PetscInt               i,cnt,m,*idx;
15111411c6eeSJed Brown   PetscErrorCode         ierr;
15121411c6eeSJed Brown 
15131411c6eeSJed Brown   PetscFunctionBegin;
15141411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
15151411c6eeSJed Brown   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
15161411c6eeSJed Brown   for (cnt=0,i=0; i<(com->nDM+com->nredundant); i++) {
1517a03f9cedSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1518a03f9cedSJed Brown     cnt += m;
1519a03f9cedSJed Brown   }
1520a03f9cedSJed Brown   ierr = PetscMalloc(cnt*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1521a03f9cedSJed Brown   for (cnt=0,i=0; i<(com->nDM+com->nredundant); i++) {
15221411c6eeSJed Brown     const PetscInt *subidx;
15231411c6eeSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
15241411c6eeSJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
15251411c6eeSJed Brown     ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr);
15261411c6eeSJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
15271411c6eeSJed Brown     cnt += m;
15281411c6eeSJed Brown   }
15291411c6eeSJed Brown   ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,cnt,idx,PETSC_OWN_POINTER,&dm->ltogmap);CHKERRQ(ierr);
1530fcfd50ebSBarry Smith   for (i=0; i<com->nDM+com->nredundant; i++) {ierr = ISLocalToGlobalMappingDestroy(&ltogs[i]);CHKERRQ(ierr);}
15311411c6eeSJed Brown   ierr = PetscFree(ltogs);CHKERRQ(ierr);
15321411c6eeSJed Brown   PetscFunctionReturn(0);
15331411c6eeSJed Brown }
15341411c6eeSJed Brown 
15351411c6eeSJed Brown 
15361411c6eeSJed Brown #undef __FUNCT__
15370c010503SBarry Smith #define __FUNCT__ "DMGetColoring_Composite"
15387087cfbeSBarry Smith PetscErrorCode  DMGetColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
153947c6ae99SBarry Smith {
154047c6ae99SBarry Smith   PetscErrorCode         ierr;
154147c6ae99SBarry Smith   PetscInt               n,i,cnt;
154247c6ae99SBarry Smith   ISColoringValue        *colors;
154347c6ae99SBarry Smith   PetscBool              dense = PETSC_FALSE;
154447c6ae99SBarry Smith   ISColoringValue        maxcol = 0;
154547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
154647c6ae99SBarry Smith 
154747c6ae99SBarry Smith   PetscFunctionBegin;
154847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
154947c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED) {
155047c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" );
155147c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL) {
155247c6ae99SBarry Smith     n = com->n;
155347c6ae99SBarry Smith   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
155447c6ae99SBarry Smith   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
155547c6ae99SBarry Smith 
1556671f6225SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
155747c6ae99SBarry Smith   if (dense) {
155847c6ae99SBarry Smith     for (i=0; i<n; i++) {
155947c6ae99SBarry Smith       colors[i] = (ISColoringValue)(com->rstart + i);
156047c6ae99SBarry Smith     }
156147c6ae99SBarry Smith     maxcol = com->N;
156247c6ae99SBarry Smith   } else {
156347c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
156447c6ae99SBarry Smith     PetscMPIInt            rank;
156547c6ae99SBarry Smith 
156647c6ae99SBarry Smith     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
156747c6ae99SBarry Smith     cnt  = 0;
156847c6ae99SBarry Smith     while (next) {
156947c6ae99SBarry Smith       if (next->type == DMCOMPOSITE_ARRAY) {
157047c6ae99SBarry Smith         if (rank == next->rank) {  /* each column gets is own color */
157147c6ae99SBarry Smith           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
157247c6ae99SBarry Smith             colors[cnt++] = maxcol++;
157347c6ae99SBarry Smith           }
157447c6ae99SBarry Smith         }
157547c6ae99SBarry Smith         ierr = MPI_Bcast(&maxcol,1,MPIU_COLORING_VALUE,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
157647c6ae99SBarry Smith       } else if (next->type == DMCOMPOSITE_DM) {
157747c6ae99SBarry Smith         ISColoring     lcoloring;
157847c6ae99SBarry Smith 
157947c6ae99SBarry Smith         ierr = DMGetColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
158047c6ae99SBarry Smith         for (i=0; i<lcoloring->N; i++) {
158147c6ae99SBarry Smith           colors[cnt++] = maxcol + lcoloring->colors[i];
158247c6ae99SBarry Smith         }
158347c6ae99SBarry Smith         maxcol += lcoloring->n;
1584fcfd50ebSBarry Smith         ierr = ISColoringDestroy(&lcoloring);CHKERRQ(ierr);
158596e147daSBarry Smith       } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
158647c6ae99SBarry Smith       next = next->next;
158747c6ae99SBarry Smith     }
158847c6ae99SBarry Smith   }
158947c6ae99SBarry Smith   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
159047c6ae99SBarry Smith   PetscFunctionReturn(0);
159147c6ae99SBarry Smith }
159247c6ae99SBarry Smith 
159347c6ae99SBarry Smith #undef __FUNCT__
15940c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
15957087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
159647c6ae99SBarry Smith {
159747c6ae99SBarry Smith   PetscErrorCode         ierr;
159847c6ae99SBarry Smith   struct DMCompositeLink *next;
159947c6ae99SBarry Smith   PetscInt               cnt = 3;
160047c6ae99SBarry Smith   PetscMPIInt            rank;
160147c6ae99SBarry Smith   PetscScalar            *garray,*larray;
160247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
160347c6ae99SBarry Smith 
160447c6ae99SBarry Smith   PetscFunctionBegin;
160547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
160647c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
160747c6ae99SBarry Smith   next = com->next;
160847c6ae99SBarry Smith   if (!com->setup) {
1609d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
161047c6ae99SBarry Smith   }
161147c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
161247c6ae99SBarry Smith   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
161347c6ae99SBarry Smith   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
161447c6ae99SBarry Smith 
161547c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
161647c6ae99SBarry Smith   while (next) {
161747c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
161847c6ae99SBarry Smith       if (rank == next->rank) {
161947c6ae99SBarry Smith         ierr    = PetscMemcpy(larray,garray,next->n*sizeof(PetscScalar));CHKERRQ(ierr);
162047c6ae99SBarry Smith         garray += next->n;
162147c6ae99SBarry Smith       }
162247c6ae99SBarry Smith       /* does not handle ADD_VALUES */
162306ebdd98SJed Brown       ierr = MPI_Bcast(larray,next->nlocal,MPIU_SCALAR,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
162447c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
162547c6ae99SBarry Smith       Vec      local,global;
162647c6ae99SBarry Smith       PetscInt N;
162747c6ae99SBarry Smith 
162847c6ae99SBarry Smith       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
162947c6ae99SBarry Smith       ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
163047c6ae99SBarry Smith       ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
163147c6ae99SBarry Smith       ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
163247c6ae99SBarry Smith       ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
163347c6ae99SBarry Smith       ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
163447c6ae99SBarry Smith       ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
163547c6ae99SBarry Smith       ierr = VecResetArray(global);CHKERRQ(ierr);
163647c6ae99SBarry Smith       ierr = VecResetArray(local);CHKERRQ(ierr);
163747c6ae99SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
163847c6ae99SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
163947c6ae99SBarry Smith     } else {
164047c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
164147c6ae99SBarry Smith     }
164247c6ae99SBarry Smith     cnt++;
164306ebdd98SJed Brown     larray += next->nlocal;
164447c6ae99SBarry Smith     next    = next->next;
164547c6ae99SBarry Smith   }
164647c6ae99SBarry Smith 
164747c6ae99SBarry Smith   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
164847c6ae99SBarry Smith   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
164947c6ae99SBarry Smith   PetscFunctionReturn(0);
165047c6ae99SBarry Smith }
165147c6ae99SBarry Smith 
165247c6ae99SBarry Smith #undef __FUNCT__
16530c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
16547087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
16550c010503SBarry Smith {
16560c010503SBarry Smith   PetscFunctionBegin;
16570c010503SBarry Smith   PetscFunctionReturn(0);
16580c010503SBarry Smith }
165947c6ae99SBarry Smith 
1660a4121054SBarry Smith EXTERN_C_BEGIN
1661a4121054SBarry Smith #undef __FUNCT__
1662a4121054SBarry Smith #define __FUNCT__ "DMCreate_Composite"
16637087cfbeSBarry Smith PetscErrorCode  DMCreate_Composite(DM p)
1664a4121054SBarry Smith {
1665a4121054SBarry Smith   PetscErrorCode ierr;
1666a4121054SBarry Smith   DM_Composite   *com;
1667a4121054SBarry Smith 
1668a4121054SBarry Smith   PetscFunctionBegin;
1669a4121054SBarry Smith   ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1670a4121054SBarry Smith   p->data = com;
1671a4121054SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1672a4121054SBarry Smith   com->n            = 0;
1673a4121054SBarry Smith   com->next         = PETSC_NULL;
1674a4121054SBarry Smith   com->nredundant   = 0;
1675a4121054SBarry Smith   com->nDM          = 0;
1676a4121054SBarry Smith 
1677a4121054SBarry Smith   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1678a4121054SBarry Smith   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
16791411c6eeSJed Brown   p->ops->createlocaltoglobalmapping      = DMCreateLocalToGlobalMapping_Composite;
16801411c6eeSJed Brown   p->ops->createlocaltoglobalmappingblock = 0;
1681a4121054SBarry Smith   p->ops->refine                          = DMRefine_Composite;
1682a4121054SBarry Smith   p->ops->getinterpolation                = DMGetInterpolation_Composite;
1683a4121054SBarry Smith   p->ops->getmatrix                       = DMGetMatrix_Composite;
1684a4121054SBarry Smith   p->ops->getcoloring                     = DMGetColoring_Composite;
1685a4121054SBarry Smith   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1686a4121054SBarry Smith   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1687a4121054SBarry Smith   p->ops->destroy                         = DMDestroy_Composite;
1688a4121054SBarry Smith   p->ops->view                            = DMView_Composite;
1689a4121054SBarry Smith   p->ops->setup                           = DMSetUp_Composite;
1690a4121054SBarry Smith   PetscFunctionReturn(0);
1691a4121054SBarry Smith }
1692a4121054SBarry Smith EXTERN_C_END
1693a4121054SBarry Smith 
16940c010503SBarry Smith #undef __FUNCT__
16950c010503SBarry Smith #define __FUNCT__ "DMCompositeCreate"
16960c010503SBarry Smith /*@C
16970c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
16980c010503SBarry Smith       vectors made up of several subvectors.
16990c010503SBarry Smith 
17000c010503SBarry Smith     Collective on MPI_Comm
170147c6ae99SBarry Smith 
170247c6ae99SBarry Smith     Input Parameter:
17030c010503SBarry Smith .   comm - the processors that will share the global vector
17040c010503SBarry Smith 
17050c010503SBarry Smith     Output Parameters:
17060c010503SBarry Smith .   packer - the packer object
170747c6ae99SBarry Smith 
170847c6ae99SBarry Smith     Level: advanced
170947c6ae99SBarry Smith 
17100c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCompositeScatter(),
17116eb61c8cSJed Brown          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
171247c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
171347c6ae99SBarry Smith 
171447c6ae99SBarry Smith @*/
17157087cfbeSBarry Smith PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
171647c6ae99SBarry Smith {
17170c010503SBarry Smith   PetscErrorCode ierr;
17180c010503SBarry Smith 
171947c6ae99SBarry Smith   PetscFunctionBegin;
17200c010503SBarry Smith   PetscValidPointer(packer,2);
1721a4121054SBarry Smith   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1722a4121054SBarry Smith   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
172347c6ae99SBarry Smith   PetscFunctionReturn(0);
172447c6ae99SBarry Smith }
1725