xref: /petsc/src/dm/impls/composite/pack.c (revision 3c48a1e8da19189ff2402a4e41a2fc082d52c349)
147c6ae99SBarry Smith 
2*3c48a1e8SJed 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 
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith extern PetscErrorCode DMDestroy_Private(DM,PetscBool *);
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith #undef __FUNCT__
880c010503SBarry Smith #define __FUNCT__ "DMDestroy_Composite"
897087cfbeSBarry Smith PetscErrorCode  DMDestroy_Composite(DM dm)
9047c6ae99SBarry Smith {
9147c6ae99SBarry Smith   PetscErrorCode         ierr;
9247c6ae99SBarry Smith   struct DMCompositeLink *next, *prev;
9347c6ae99SBarry Smith   PetscBool              done;
9447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite *)dm->data;
9547c6ae99SBarry Smith 
9647c6ae99SBarry Smith   PetscFunctionBegin;
9747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
9847c6ae99SBarry Smith   ierr = DMDestroy_Private((DM)dm,&done);CHKERRQ(ierr);
9947c6ae99SBarry Smith   if (!done) PetscFunctionReturn(0);
10047c6ae99SBarry Smith 
10147c6ae99SBarry Smith   next = com->next;
10247c6ae99SBarry Smith   while (next) {
10347c6ae99SBarry Smith     prev = next;
10447c6ae99SBarry Smith     next = next->next;
10547c6ae99SBarry Smith     if (prev->type == DMCOMPOSITE_DM) {
10647c6ae99SBarry Smith       ierr = DMDestroy(prev->dm);CHKERRQ(ierr);
10747c6ae99SBarry Smith     }
10847c6ae99SBarry Smith     ierr = PetscFree(prev->grstarts);CHKERRQ(ierr);
10947c6ae99SBarry Smith     ierr = PetscFree(prev);CHKERRQ(ierr);
11047c6ae99SBarry Smith   }
11147c6ae99SBarry Smith   ierr = PetscFree(com);CHKERRQ(ierr);
11247c6ae99SBarry Smith   ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr);
11347c6ae99SBarry Smith   PetscFunctionReturn(0);
11447c6ae99SBarry Smith }
11547c6ae99SBarry Smith 
11647c6ae99SBarry Smith #undef __FUNCT__
1170c010503SBarry Smith #define __FUNCT__ "DMView_Composite"
1187087cfbeSBarry Smith PetscErrorCode  DMView_Composite(DM dm,PetscViewer v)
11947c6ae99SBarry Smith {
12047c6ae99SBarry Smith   PetscErrorCode ierr;
12147c6ae99SBarry Smith   PetscBool      iascii;
12247c6ae99SBarry Smith   DM_Composite   *com = (DM_Composite *)dm->data;
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith   PetscFunctionBegin;
12547c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)v,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
12647c6ae99SBarry Smith   if (iascii) {
12747c6ae99SBarry Smith     struct DMCompositeLink *lnk = com->next;
12847c6ae99SBarry Smith     PetscInt               i;
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"DM (%s)\n",((PetscObject)dm)->prefix?((PetscObject)dm)->prefix:"no prefix");CHKERRQ(ierr);
13147c6ae99SBarry Smith     ierr = PetscViewerASCIIPrintf(v,"  contains %d DMs and %d redundant arrays\n",com->nDM,com->nredundant);CHKERRQ(ierr);
13247c6ae99SBarry Smith     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
13347c6ae99SBarry Smith     for (i=0; lnk; lnk=lnk->next,i++) {
13447c6ae99SBarry Smith       if (lnk->dm) {
13547c6ae99SBarry Smith         ierr = PetscViewerASCIIPrintf(v,"Link %d: DM of type %s\n",i,((PetscObject)lnk->dm)->type_name);CHKERRQ(ierr);
13647c6ae99SBarry Smith         ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
13747c6ae99SBarry Smith         ierr = DMView(lnk->dm,v);CHKERRQ(ierr);
13847c6ae99SBarry Smith         ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
13947c6ae99SBarry Smith       } else {
14006ebdd98SJed Brown         ierr = PetscViewerASCIIPrintf(v,"Link %d: Redundant array of size %d owned by rank %d\n",i,lnk->nlocal,lnk->rank);CHKERRQ(ierr);
14147c6ae99SBarry Smith       }
14247c6ae99SBarry Smith     }
14347c6ae99SBarry Smith     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
14447c6ae99SBarry Smith   }
14547c6ae99SBarry Smith   PetscFunctionReturn(0);
14647c6ae99SBarry Smith }
14747c6ae99SBarry Smith 
14847c6ae99SBarry Smith /* --------------------------------------------------------------------------------------*/
14947c6ae99SBarry Smith #undef __FUNCT__
150d7bf68aeSBarry Smith #define __FUNCT__ "DMSetUp_Composite"
1517087cfbeSBarry Smith PetscErrorCode  DMSetUp_Composite(DM dm)
15247c6ae99SBarry Smith {
15347c6ae99SBarry Smith   PetscErrorCode         ierr;
15447c6ae99SBarry Smith   PetscInt               nprev = 0;
15547c6ae99SBarry Smith   PetscMPIInt            rank,size;
15647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
15747c6ae99SBarry Smith   struct DMCompositeLink *next = com->next;
15847c6ae99SBarry Smith   PetscLayout            map;
15947c6ae99SBarry Smith 
16047c6ae99SBarry Smith   PetscFunctionBegin;
16147c6ae99SBarry Smith   if (com->setup) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Packer has already been setup");
16247c6ae99SBarry Smith   ierr = PetscLayoutCreate(((PetscObject)dm)->comm,&map);CHKERRQ(ierr);
16347c6ae99SBarry Smith   ierr = PetscLayoutSetLocalSize(map,com->n);CHKERRQ(ierr);
16447c6ae99SBarry Smith   ierr = PetscLayoutSetSize(map,PETSC_DETERMINE);CHKERRQ(ierr);
16547c6ae99SBarry Smith   ierr = PetscLayoutSetBlockSize(map,1);CHKERRQ(ierr);
16647c6ae99SBarry Smith   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
16747c6ae99SBarry Smith   ierr = PetscLayoutGetSize(map,&com->N);CHKERRQ(ierr);
16847c6ae99SBarry Smith   ierr = PetscLayoutGetRange(map,&com->rstart,PETSC_NULL);CHKERRQ(ierr);
16947c6ae99SBarry Smith   ierr = PetscLayoutDestroy(map);CHKERRQ(ierr);
17047c6ae99SBarry Smith 
17147c6ae99SBarry Smith   /* now set the rstart for each linked array/vector */
17247c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
17347c6ae99SBarry Smith   ierr = MPI_Comm_size(((PetscObject)dm)->comm,&size);CHKERRQ(ierr);
17447c6ae99SBarry Smith   while (next) {
17547c6ae99SBarry Smith     next->rstart = nprev;
17606ebdd98SJed Brown     nprev += next->n;
17747c6ae99SBarry Smith     next->grstart = com->rstart + next->rstart;
17847c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
17947c6ae99SBarry Smith       ierr = MPI_Bcast(&next->grstart,1,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
18047c6ae99SBarry Smith     } else {
18147c6ae99SBarry Smith       ierr = PetscMalloc(size*sizeof(PetscInt),&next->grstarts);CHKERRQ(ierr);
18247c6ae99SBarry Smith       ierr = MPI_Allgather(&next->grstart,1,MPIU_INT,next->grstarts,1,MPIU_INT,((PetscObject)dm)->comm);CHKERRQ(ierr);
18347c6ae99SBarry Smith     }
18447c6ae99SBarry Smith     next = next->next;
18547c6ae99SBarry Smith   }
18647c6ae99SBarry Smith   com->setup = PETSC_TRUE;
18747c6ae99SBarry Smith   PetscFunctionReturn(0);
18847c6ae99SBarry Smith }
18947c6ae99SBarry Smith 
19047c6ae99SBarry Smith 
19147c6ae99SBarry Smith #undef __FUNCT__
19247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess_Array"
19347c6ae99SBarry Smith PetscErrorCode DMCompositeGetAccess_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar **array)
19447c6ae99SBarry Smith {
19547c6ae99SBarry Smith   PetscErrorCode ierr;
19647c6ae99SBarry Smith   PetscScalar    *varray;
19747c6ae99SBarry Smith   PetscMPIInt    rank;
19847c6ae99SBarry Smith 
19947c6ae99SBarry Smith   PetscFunctionBegin;
20047c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
20147c6ae99SBarry Smith   if (array) {
20247c6ae99SBarry Smith     if (rank == mine->rank) {
20347c6ae99SBarry Smith       ierr    = VecGetArray(vec,&varray);CHKERRQ(ierr);
20447c6ae99SBarry Smith       *array  = varray + mine->rstart;
20547c6ae99SBarry Smith       ierr    = VecRestoreArray(vec,&varray);CHKERRQ(ierr);
20647c6ae99SBarry Smith     } else {
20747c6ae99SBarry Smith       *array = 0;
20847c6ae99SBarry Smith     }
20947c6ae99SBarry Smith   }
21047c6ae99SBarry Smith   PetscFunctionReturn(0);
21147c6ae99SBarry Smith }
21247c6ae99SBarry Smith 
21347c6ae99SBarry Smith #undef __FUNCT__
21447c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess_DM"
21547c6ae99SBarry Smith PetscErrorCode DMCompositeGetAccess_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec *global)
21647c6ae99SBarry Smith {
21747c6ae99SBarry Smith   PetscErrorCode ierr;
21847c6ae99SBarry Smith   PetscScalar    *array;
21947c6ae99SBarry Smith 
22047c6ae99SBarry Smith   PetscFunctionBegin;
22147c6ae99SBarry Smith   if (global) {
22247c6ae99SBarry Smith     ierr    = DMGetGlobalVector(mine->dm,global);CHKERRQ(ierr);
22347c6ae99SBarry Smith     ierr    = VecGetArray(vec,&array);CHKERRQ(ierr);
22447c6ae99SBarry Smith     ierr    = VecPlaceArray(*global,array+mine->rstart);CHKERRQ(ierr);
22547c6ae99SBarry Smith     ierr    = VecRestoreArray(vec,&array);CHKERRQ(ierr);
22647c6ae99SBarry Smith   }
22747c6ae99SBarry Smith   PetscFunctionReturn(0);
22847c6ae99SBarry Smith }
22947c6ae99SBarry Smith 
23047c6ae99SBarry Smith #undef __FUNCT__
23147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess_Array"
23247c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreAccess_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar **array)
23347c6ae99SBarry Smith {
23447c6ae99SBarry Smith   PetscFunctionBegin;
23547c6ae99SBarry Smith   PetscFunctionReturn(0);
23647c6ae99SBarry Smith }
23747c6ae99SBarry Smith 
23847c6ae99SBarry Smith #undef __FUNCT__
23947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess_DM"
24047c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreAccess_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec *global)
24147c6ae99SBarry Smith {
24247c6ae99SBarry Smith   PetscErrorCode ierr;
24347c6ae99SBarry Smith 
24447c6ae99SBarry Smith   PetscFunctionBegin;
24547c6ae99SBarry Smith   if (global) {
24647c6ae99SBarry Smith     ierr = VecResetArray(*global);CHKERRQ(ierr);
24747c6ae99SBarry Smith     ierr = DMRestoreGlobalVector(mine->dm,global);CHKERRQ(ierr);
24847c6ae99SBarry Smith   }
24947c6ae99SBarry Smith   PetscFunctionReturn(0);
25047c6ae99SBarry Smith }
25147c6ae99SBarry Smith 
25247c6ae99SBarry Smith #undef __FUNCT__
25347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter_Array"
25447c6ae99SBarry Smith PetscErrorCode DMCompositeScatter_Array(DM dm,struct DMCompositeLink *mine,Vec vec,PetscScalar *array)
25547c6ae99SBarry Smith {
25647c6ae99SBarry Smith   PetscErrorCode ierr;
25747c6ae99SBarry Smith   PetscScalar    *varray;
25847c6ae99SBarry Smith   PetscMPIInt    rank;
25947c6ae99SBarry Smith 
26047c6ae99SBarry Smith   PetscFunctionBegin;
26147c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
26247c6ae99SBarry Smith   if (rank == mine->rank) {
26347c6ae99SBarry Smith     ierr    = VecGetArray(vec,&varray);CHKERRQ(ierr);
26447c6ae99SBarry Smith     ierr    = PetscMemcpy(array,varray+mine->rstart,mine->n*sizeof(PetscScalar));CHKERRQ(ierr);
26547c6ae99SBarry Smith     ierr    = VecRestoreArray(vec,&varray);CHKERRQ(ierr);
26647c6ae99SBarry Smith   }
26706ebdd98SJed Brown   ierr    = MPI_Bcast(array,mine->nlocal,MPIU_SCALAR,mine->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
26847c6ae99SBarry Smith   PetscFunctionReturn(0);
26947c6ae99SBarry Smith }
27047c6ae99SBarry Smith 
27147c6ae99SBarry Smith #undef __FUNCT__
27247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter_DM"
27347c6ae99SBarry Smith PetscErrorCode DMCompositeScatter_DM(DM dm,struct DMCompositeLink *mine,Vec vec,Vec local)
27447c6ae99SBarry Smith {
27547c6ae99SBarry Smith   PetscErrorCode ierr;
27647c6ae99SBarry Smith   PetscScalar    *array;
27747c6ae99SBarry Smith   Vec            global;
27847c6ae99SBarry Smith 
27947c6ae99SBarry Smith   PetscFunctionBegin;
28047c6ae99SBarry Smith   ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr);
28147c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
28247c6ae99SBarry Smith   ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr);
28347c6ae99SBarry Smith   ierr = DMGlobalToLocalBegin(mine->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
28447c6ae99SBarry Smith   ierr = DMGlobalToLocalEnd(mine->dm,global,INSERT_VALUES,local);CHKERRQ(ierr);
28547c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
28647c6ae99SBarry Smith   ierr = VecResetArray(global);CHKERRQ(ierr);
28747c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr);
28847c6ae99SBarry Smith   PetscFunctionReturn(0);
28947c6ae99SBarry Smith }
29047c6ae99SBarry Smith 
29147c6ae99SBarry Smith #undef __FUNCT__
29247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_Array"
293acc1e270SJed Brown PetscErrorCode DMCompositeGather_Array(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,const PetscScalar *array)
29447c6ae99SBarry Smith {
29547c6ae99SBarry Smith   PetscErrorCode ierr;
29647c6ae99SBarry Smith   PetscScalar    *varray;
29747c6ae99SBarry Smith   PetscMPIInt    rank;
29847c6ae99SBarry Smith 
29947c6ae99SBarry Smith   PetscFunctionBegin;
30047c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
30147c6ae99SBarry Smith   if (rank == mine->rank) {
30247c6ae99SBarry Smith     ierr = VecGetArray(vec,&varray);CHKERRQ(ierr);
30347c6ae99SBarry Smith     if (varray+mine->rstart == array) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"You need not DMCompositeGather() into objects obtained via DMCompositeGetAccess()");
304acc1e270SJed Brown   }
305df0c820aSJed Brown   switch (imode) {
306df0c820aSJed Brown   case INSERT_VALUES:
307acc1e270SJed Brown     if (rank == mine->rank) {
30847c6ae99SBarry Smith       ierr = PetscMemcpy(varray+mine->rstart,array,mine->n*sizeof(PetscScalar));CHKERRQ(ierr);
309acc1e270SJed Brown     }
310df0c820aSJed Brown     break;
311df0c820aSJed Brown   case ADD_VALUES: {
312df0c820aSJed Brown     PetscInt          i;
313760fd489SMatthew G Knepley     void             *source;
314acc1e270SJed Brown     PetscScalar       *buffer,*dest;
315acc1e270SJed Brown     if (rank == mine->rank) {
316acc1e270SJed Brown       dest = &varray[mine->rstart];
317acc1e270SJed Brown #if defined(PETSC_HAVE_MPI_IN_PLACE)
318acc1e270SJed Brown       buffer = dest;
319acc1e270SJed Brown       source = MPI_IN_PLACE;
320acc1e270SJed Brown #else
32106ebdd98SJed Brown       ierr = PetscMalloc(mine->nlocal*sizeof(PetscScalar),&buffer);CHKERRQ(ierr);
322760fd489SMatthew G Knepley       source = (void *) buffer;
323acc1e270SJed Brown #endif
32406ebdd98SJed Brown       for (i=0; i<mine->nlocal; i++) buffer[i] = varray[mine->rstart+i] + array[i];
325acc1e270SJed Brown     } else {
326760fd489SMatthew G Knepley       source = (void *) array;
327acc1e270SJed Brown       dest   = PETSC_NULL;
328acc1e270SJed Brown     }
32906ebdd98SJed Brown     ierr = MPI_Reduce(source,dest,mine->nlocal,MPIU_SCALAR,MPI_SUM,mine->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
330acc1e270SJed Brown #if !defined(PETSC_HAVE_MPI_IN_PLACE)
331acc1e270SJed Brown     if (rank == mine->rank) {ierr = PetscFree(source);CHKERRQ(ierr);}
332acc1e270SJed Brown #endif
333df0c820aSJed Brown   } break;
334df0c820aSJed Brown   default: SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"imode");
335df0c820aSJed Brown   }
336acc1e270SJed Brown   if (rank == mine->rank) {ierr = VecRestoreArray(vec,&varray);CHKERRQ(ierr);}
33747c6ae99SBarry Smith   PetscFunctionReturn(0);
33847c6ae99SBarry Smith }
33947c6ae99SBarry Smith 
34047c6ae99SBarry Smith #undef __FUNCT__
34147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather_DM"
342df0c820aSJed Brown PetscErrorCode DMCompositeGather_DM(DM dm,struct DMCompositeLink *mine,Vec vec,InsertMode imode,Vec local)
34347c6ae99SBarry Smith {
34447c6ae99SBarry Smith   PetscErrorCode ierr;
34547c6ae99SBarry Smith   PetscScalar    *array;
34647c6ae99SBarry Smith   Vec            global;
34747c6ae99SBarry Smith 
34847c6ae99SBarry Smith   PetscFunctionBegin;
34947c6ae99SBarry Smith   ierr = DMGetGlobalVector(mine->dm,&global);CHKERRQ(ierr);
35047c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
35147c6ae99SBarry Smith   ierr = VecPlaceArray(global,array+mine->rstart);CHKERRQ(ierr);
352df0c820aSJed Brown   ierr = DMLocalToGlobalBegin(mine->dm,local,imode,global);CHKERRQ(ierr);
353df0c820aSJed Brown   ierr = DMLocalToGlobalEnd(mine->dm,local,imode,global);CHKERRQ(ierr);
35447c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
35547c6ae99SBarry Smith   ierr = VecResetArray(global);CHKERRQ(ierr);
35647c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(mine->dm,&global);CHKERRQ(ierr);
35747c6ae99SBarry Smith   PetscFunctionReturn(0);
35847c6ae99SBarry Smith }
35947c6ae99SBarry Smith 
36047c6ae99SBarry Smith /* ----------------------------------------------------------------------------------*/
36147c6ae99SBarry Smith 
36247c6ae99SBarry Smith #include <stdarg.h>
36347c6ae99SBarry Smith 
36447c6ae99SBarry Smith #undef __FUNCT__
36547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetNumberDM"
36647c6ae99SBarry Smith /*@C
36747c6ae99SBarry Smith     DMCompositeGetNumberDM - Get's the number of DM objects in the DMComposite
36847c6ae99SBarry Smith        representation.
36947c6ae99SBarry Smith 
37047c6ae99SBarry Smith     Not Collective
37147c6ae99SBarry Smith 
37247c6ae99SBarry Smith     Input Parameter:
37347c6ae99SBarry Smith .    dm - the packer object
37447c6ae99SBarry Smith 
37547c6ae99SBarry Smith     Output Parameter:
37647c6ae99SBarry Smith .     nDM - the number of DMs
37747c6ae99SBarry Smith 
37847c6ae99SBarry Smith     Level: beginner
37947c6ae99SBarry Smith 
38047c6ae99SBarry Smith @*/
3817087cfbeSBarry Smith PetscErrorCode  DMCompositeGetNumberDM(DM dm,PetscInt *nDM)
38247c6ae99SBarry Smith {
38347c6ae99SBarry Smith   DM_Composite *com = (DM_Composite*)dm->data;
38447c6ae99SBarry Smith   PetscFunctionBegin;
38547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
38647c6ae99SBarry Smith   *nDM = com->nDM;
38747c6ae99SBarry Smith   PetscFunctionReturn(0);
38847c6ae99SBarry Smith }
38947c6ae99SBarry Smith 
39047c6ae99SBarry Smith 
39147c6ae99SBarry Smith #undef __FUNCT__
39247c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetAccess"
39347c6ae99SBarry Smith /*@C
39447c6ae99SBarry Smith     DMCompositeGetAccess - Allows one to access the individual packed vectors in their global
39547c6ae99SBarry Smith        representation.
39647c6ae99SBarry Smith 
39747c6ae99SBarry Smith     Collective on DMComposite
39847c6ae99SBarry Smith 
39947c6ae99SBarry Smith     Input Parameter:
40047c6ae99SBarry Smith +    dm - the packer object
40147c6ae99SBarry Smith .    gvec - the global vector
40247c6ae99SBarry Smith -    ... - the individual sequential or parallel objects (arrays or vectors)
40347c6ae99SBarry Smith 
40447c6ae99SBarry Smith     Notes: Use DMCompositeRestoreAccess() to return the vectors when you no longer need them
40547c6ae99SBarry Smith 
40647c6ae99SBarry Smith     Level: advanced
40747c6ae99SBarry Smith 
40847c6ae99SBarry Smith @*/
4097087cfbeSBarry Smith PetscErrorCode  DMCompositeGetAccess(DM dm,Vec gvec,...)
41047c6ae99SBarry Smith {
41147c6ae99SBarry Smith   va_list                Argp;
41247c6ae99SBarry Smith   PetscErrorCode         ierr;
41347c6ae99SBarry Smith   struct DMCompositeLink *next;
41447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
41547c6ae99SBarry Smith 
41647c6ae99SBarry Smith   PetscFunctionBegin;
41747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
41847c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
41947c6ae99SBarry Smith   next = com->next;
42047c6ae99SBarry Smith   if (!com->setup) {
421d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
42247c6ae99SBarry Smith   }
42347c6ae99SBarry Smith 
42447c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
42547c6ae99SBarry Smith   va_start(Argp,gvec);
42647c6ae99SBarry Smith   while (next) {
42747c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
42847c6ae99SBarry Smith       PetscScalar **array;
42947c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
43047c6ae99SBarry Smith       ierr  = DMCompositeGetAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
43147c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
43247c6ae99SBarry Smith       Vec *vec;
43347c6ae99SBarry Smith       vec  = va_arg(Argp, Vec*);
43447c6ae99SBarry Smith       ierr = DMCompositeGetAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
43547c6ae99SBarry Smith     } else {
43647c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
43747c6ae99SBarry Smith     }
43847c6ae99SBarry Smith     next = next->next;
43947c6ae99SBarry Smith   }
44047c6ae99SBarry Smith   va_end(Argp);
44147c6ae99SBarry Smith   PetscFunctionReturn(0);
44247c6ae99SBarry Smith }
44347c6ae99SBarry Smith 
44447c6ae99SBarry Smith #undef __FUNCT__
44547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreAccess"
44647c6ae99SBarry Smith /*@C
447aa219208SBarry Smith     DMCompositeRestoreAccess - Returns the vectors obtained with DMCompositeGetAccess()
44847c6ae99SBarry Smith        representation.
44947c6ae99SBarry Smith 
45047c6ae99SBarry Smith     Collective on DMComposite
45147c6ae99SBarry Smith 
45247c6ae99SBarry Smith     Input Parameter:
45347c6ae99SBarry Smith +    dm - the packer object
45447c6ae99SBarry Smith .    gvec - the global vector
45547c6ae99SBarry Smith -    ... - the individual sequential or parallel objects (arrays or vectors)
45647c6ae99SBarry Smith 
45747c6ae99SBarry Smith     Level: advanced
45847c6ae99SBarry Smith 
4590c010503SBarry Smith .seealso  DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
4606eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeScatter(),
461aa219208SBarry Smith          DMCompositeRestoreAccess(), DMCompositeGetAccess()
46247c6ae99SBarry Smith 
46347c6ae99SBarry Smith @*/
4647087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreAccess(DM dm,Vec gvec,...)
46547c6ae99SBarry Smith {
46647c6ae99SBarry Smith   va_list                Argp;
46747c6ae99SBarry Smith   PetscErrorCode         ierr;
46847c6ae99SBarry Smith   struct DMCompositeLink *next;
46947c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
47047c6ae99SBarry Smith 
47147c6ae99SBarry Smith   PetscFunctionBegin;
47247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
47347c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
47447c6ae99SBarry Smith   next = com->next;
47547c6ae99SBarry Smith   if (!com->setup) {
476d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
47747c6ae99SBarry Smith   }
47847c6ae99SBarry Smith 
47947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
48047c6ae99SBarry Smith   va_start(Argp,gvec);
48147c6ae99SBarry Smith   while (next) {
48247c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
48347c6ae99SBarry Smith       PetscScalar **array;
48447c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
48547c6ae99SBarry Smith       ierr  = DMCompositeRestoreAccess_Array(dm,next,gvec,array);CHKERRQ(ierr);
48647c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
48747c6ae99SBarry Smith       Vec *vec;
48847c6ae99SBarry Smith       vec  = va_arg(Argp, Vec*);
48947c6ae99SBarry Smith       ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,vec);CHKERRQ(ierr);
49047c6ae99SBarry Smith     } else {
49147c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
49247c6ae99SBarry Smith     }
49347c6ae99SBarry Smith     next = next->next;
49447c6ae99SBarry Smith   }
49547c6ae99SBarry Smith   va_end(Argp);
49647c6ae99SBarry Smith   PetscFunctionReturn(0);
49747c6ae99SBarry Smith }
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith #undef __FUNCT__
50047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeScatter"
50147c6ae99SBarry Smith /*@C
50247c6ae99SBarry Smith     DMCompositeScatter - Scatters from a global packed vector into its individual local vectors
50347c6ae99SBarry Smith 
50447c6ae99SBarry Smith     Collective on DMComposite
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith     Input Parameter:
50747c6ae99SBarry Smith +    dm - the packer object
50847c6ae99SBarry Smith .    gvec - the global vector
5098fd8f222SJed Brown -    ... - the individual sequential objects (arrays or vectors), PETSC_NULL for those that are not needed
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith     Level: advanced
51247c6ae99SBarry Smith 
5130c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
5146eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
51547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
51647c6ae99SBarry Smith 
51747c6ae99SBarry Smith @*/
5187087cfbeSBarry Smith PetscErrorCode  DMCompositeScatter(DM dm,Vec gvec,...)
51947c6ae99SBarry Smith {
52047c6ae99SBarry Smith   va_list                Argp;
52147c6ae99SBarry Smith   PetscErrorCode         ierr;
52247c6ae99SBarry Smith   struct DMCompositeLink *next;
5238fd8f222SJed Brown   PetscInt               cnt;
52447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
52547c6ae99SBarry Smith 
52647c6ae99SBarry Smith   PetscFunctionBegin;
52747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
52847c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
52947c6ae99SBarry Smith   if (!com->setup) {
530d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
53147c6ae99SBarry Smith   }
53247c6ae99SBarry Smith 
53347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
53447c6ae99SBarry Smith   va_start(Argp,gvec);
5358fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
53647c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
53747c6ae99SBarry Smith       PetscScalar *array;
53847c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar*);
539c72eef85SJed Brown       if (array) PetscValidScalarPointer(array,cnt);
540c72eef85SJed Brown       PetscValidLogicalCollectiveBool(dm,!!array,cnt);
5418fd8f222SJed Brown       if (!array) continue;
54247c6ae99SBarry Smith       ierr = DMCompositeScatter_Array(dm,next,gvec,array);CHKERRQ(ierr);
54347c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
54447c6ae99SBarry Smith       Vec vec;
54547c6ae99SBarry Smith       vec = va_arg(Argp, Vec);
5468fd8f222SJed Brown       if (!vec) continue;
54747c6ae99SBarry Smith       PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt);
54847c6ae99SBarry Smith       ierr = DMCompositeScatter_DM(dm,next,gvec,vec);CHKERRQ(ierr);
54947c6ae99SBarry Smith     } else {
55047c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
55147c6ae99SBarry Smith     }
55247c6ae99SBarry Smith   }
55347c6ae99SBarry Smith   va_end(Argp);
55447c6ae99SBarry Smith   PetscFunctionReturn(0);
55547c6ae99SBarry Smith }
55647c6ae99SBarry Smith 
55747c6ae99SBarry Smith #undef __FUNCT__
55847c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGather"
55947c6ae99SBarry Smith /*@C
56047c6ae99SBarry Smith     DMCompositeGather - Gathers into a global packed vector from its individual local vectors
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith     Collective on DMComposite
56347c6ae99SBarry Smith 
56447c6ae99SBarry Smith     Input Parameter:
56547c6ae99SBarry Smith +    dm - the packer object
56647c6ae99SBarry Smith .    gvec - the global vector
5678fd8f222SJed Brown -    ... - the individual sequential objects (arrays or vectors), PETSC_NULL for any that are not needed
56847c6ae99SBarry Smith 
56947c6ae99SBarry Smith     Level: advanced
57047c6ae99SBarry Smith 
5710c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
5726eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
57347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
57447c6ae99SBarry Smith 
57547c6ae99SBarry Smith @*/
5767087cfbeSBarry Smith PetscErrorCode  DMCompositeGather(DM dm,Vec gvec,InsertMode imode,...)
57747c6ae99SBarry Smith {
57847c6ae99SBarry Smith   va_list                Argp;
57947c6ae99SBarry Smith   PetscErrorCode         ierr;
58047c6ae99SBarry Smith   struct DMCompositeLink *next;
58147c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
5828fd8f222SJed Brown   PetscInt               cnt;
58347c6ae99SBarry Smith 
58447c6ae99SBarry Smith   PetscFunctionBegin;
58547c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
58647c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
58747c6ae99SBarry Smith   if (!com->setup) {
588d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
58947c6ae99SBarry Smith   }
59047c6ae99SBarry Smith 
59147c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
592df0c820aSJed Brown   va_start(Argp,imode);
5938fd8f222SJed Brown   for (cnt=3,next=com->next; next; cnt++,next=next->next) {
59447c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
59547c6ae99SBarry Smith       PetscScalar *array;
59647c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar*);
5978fd8f222SJed Brown       if (!array) continue;
5988fd8f222SJed Brown       PetscValidScalarPointer(array,cnt);
599df0c820aSJed Brown       ierr  = DMCompositeGather_Array(dm,next,gvec,imode,array);CHKERRQ(ierr);
60047c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
60147c6ae99SBarry Smith       Vec vec;
60247c6ae99SBarry Smith       vec = va_arg(Argp, Vec);
6038fd8f222SJed Brown       if (!vec) continue;
6048fd8f222SJed Brown       PetscValidHeaderSpecific(vec,VEC_CLASSID,cnt);
605df0c820aSJed Brown       ierr = DMCompositeGather_DM(dm,next,gvec,imode,vec);CHKERRQ(ierr);
60647c6ae99SBarry Smith     } else {
60747c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
60847c6ae99SBarry Smith     }
60947c6ae99SBarry Smith   }
61047c6ae99SBarry Smith   va_end(Argp);
61147c6ae99SBarry Smith   PetscFunctionReturn(0);
61247c6ae99SBarry Smith }
61347c6ae99SBarry Smith 
61447c6ae99SBarry Smith #undef __FUNCT__
61547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddArray"
61647c6ae99SBarry Smith /*@C
61747c6ae99SBarry Smith     DMCompositeAddArray - adds an "redundant" array to a DMComposite. The array values will
61847c6ae99SBarry Smith        be stored in part of the array on process orank.
61947c6ae99SBarry Smith 
62047c6ae99SBarry Smith     Collective on DMComposite
62147c6ae99SBarry Smith 
62247c6ae99SBarry Smith     Input Parameter:
62347c6ae99SBarry Smith +    dm - the packer object
62447c6ae99SBarry Smith .    orank - the process on which the array entries officially live, this number must be
62547c6ae99SBarry Smith              the same on all processes.
62647c6ae99SBarry Smith -    n - the length of the array
62747c6ae99SBarry Smith 
62847c6ae99SBarry Smith     Level: advanced
62947c6ae99SBarry Smith 
6300c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
6316eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
63247c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
63347c6ae99SBarry Smith 
63447c6ae99SBarry Smith @*/
6357087cfbeSBarry Smith PetscErrorCode  DMCompositeAddArray(DM dm,PetscMPIInt orank,PetscInt n)
63647c6ae99SBarry Smith {
63747c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
63847c6ae99SBarry Smith   PetscErrorCode         ierr;
63947c6ae99SBarry Smith   PetscMPIInt            rank;
64047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
64147c6ae99SBarry Smith 
64247c6ae99SBarry Smith   PetscFunctionBegin;
64347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
64447c6ae99SBarry Smith   next = com->next;
64547c6ae99SBarry Smith   if (com->setup) {
64647c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add an array once you have used the DMComposite");
64747c6ae99SBarry Smith   }
64847c6ae99SBarry Smith #if defined(PETSC_USE_DEBUG)
64947c6ae99SBarry Smith   {
65047c6ae99SBarry Smith     PetscMPIInt orankmax;
65147c6ae99SBarry Smith     ierr = MPI_Allreduce(&orank,&orankmax,1,MPI_INT,MPI_MAX,((PetscObject)dm)->comm);CHKERRQ(ierr);
65247c6ae99SBarry 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);
65347c6ae99SBarry Smith   }
65447c6ae99SBarry Smith #endif
65547c6ae99SBarry Smith 
65647c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
65747c6ae99SBarry Smith   /* create new link */
65847c6ae99SBarry Smith   ierr                = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
65906ebdd98SJed Brown   mine->nlocal        = n;
66006ebdd98SJed Brown   mine->n             = (rank == orank) ? n : 0;
66147c6ae99SBarry Smith   mine->rank          = orank;
66247c6ae99SBarry Smith   mine->dm            = PETSC_NULL;
66347c6ae99SBarry Smith   mine->type          = DMCOMPOSITE_ARRAY;
66447c6ae99SBarry Smith   mine->next          = PETSC_NULL;
66547c6ae99SBarry Smith   if (rank == mine->rank) {com->n += n;com->nmine++;}
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith   /* add to end of list */
66847c6ae99SBarry Smith   if (!next) {
66947c6ae99SBarry Smith     com->next = mine;
67047c6ae99SBarry Smith   } else {
67147c6ae99SBarry Smith     while (next->next) next = next->next;
67247c6ae99SBarry Smith     next->next = mine;
67347c6ae99SBarry Smith   }
67447c6ae99SBarry Smith   com->nredundant++;
67547c6ae99SBarry Smith   PetscFunctionReturn(0);
67647c6ae99SBarry Smith }
67747c6ae99SBarry Smith 
67847c6ae99SBarry Smith #undef __FUNCT__
67947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeAddDM"
68047c6ae99SBarry Smith /*@C
681aa219208SBarry Smith     DMCompositeAddDM - adds a DM  vector to a DMComposite
68247c6ae99SBarry Smith 
68347c6ae99SBarry Smith     Collective on DMComposite
68447c6ae99SBarry Smith 
68547c6ae99SBarry Smith     Input Parameter:
68647c6ae99SBarry Smith +    dm - the packer object
68747c6ae99SBarry Smith -    dm - the DM object, if the DM is a da you will need to caste it with a (DM)
68847c6ae99SBarry Smith 
68947c6ae99SBarry Smith     Level: advanced
69047c6ae99SBarry Smith 
6910c010503SBarry Smith .seealso DMDestroy(), DMCompositeGather(), DMCompositeAddDM(), DMCreateGlobalVector(),
6926eb61c8cSJed Brown          DMCompositeScatter(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
69347c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
69447c6ae99SBarry Smith 
69547c6ae99SBarry Smith @*/
6967087cfbeSBarry Smith PetscErrorCode  DMCompositeAddDM(DM dmc,DM dm)
69747c6ae99SBarry Smith {
69847c6ae99SBarry Smith   PetscErrorCode         ierr;
69906ebdd98SJed Brown   PetscInt               n,nlocal;
70047c6ae99SBarry Smith   struct DMCompositeLink *mine,*next;
70106ebdd98SJed Brown   Vec                    global,local;
70247c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmc->data;
70347c6ae99SBarry Smith 
70447c6ae99SBarry Smith   PetscFunctionBegin;
70547c6ae99SBarry Smith   PetscValidHeaderSpecific(dmc,DM_CLASSID,1);
70647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,2);
70747c6ae99SBarry Smith   next = com->next;
708aa219208SBarry Smith   if (com->setup) SETERRQ(((PetscObject)dmc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot add a DM once you have used the DMComposite");
70947c6ae99SBarry Smith 
71047c6ae99SBarry Smith   /* create new link */
71147c6ae99SBarry Smith   ierr = PetscNew(struct DMCompositeLink,&mine);CHKERRQ(ierr);
71247c6ae99SBarry Smith   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
71347c6ae99SBarry Smith   ierr = DMGetGlobalVector(dm,&global);CHKERRQ(ierr);
71447c6ae99SBarry Smith   ierr = VecGetLocalSize(global,&n);CHKERRQ(ierr);
71547c6ae99SBarry Smith   ierr = DMRestoreGlobalVector(dm,&global);CHKERRQ(ierr);
71606ebdd98SJed Brown   ierr = DMGetLocalVector(dm,&local);CHKERRQ(ierr);
71706ebdd98SJed Brown   ierr = VecGetSize(local,&nlocal);CHKERRQ(ierr);
71806ebdd98SJed Brown   ierr = DMRestoreLocalVector(dm,&local);CHKERRQ(ierr);
71947c6ae99SBarry Smith   mine->n      = n;
72006ebdd98SJed Brown   mine->nlocal = nlocal;
72147c6ae99SBarry Smith   mine->dm     = dm;
72247c6ae99SBarry Smith   mine->type   = DMCOMPOSITE_DM;
72347c6ae99SBarry Smith   mine->next   = PETSC_NULL;
72447c6ae99SBarry Smith   com->n       += n;
72547c6ae99SBarry Smith 
72647c6ae99SBarry Smith   /* add to end of list */
72747c6ae99SBarry Smith   if (!next) {
72847c6ae99SBarry Smith     com->next = mine;
72947c6ae99SBarry Smith   } else {
73047c6ae99SBarry Smith     while (next->next) next = next->next;
73147c6ae99SBarry Smith     next->next = mine;
73247c6ae99SBarry Smith   }
73347c6ae99SBarry Smith   com->nDM++;
73447c6ae99SBarry Smith   com->nmine++;
73547c6ae99SBarry Smith   PetscFunctionReturn(0);
73647c6ae99SBarry Smith }
73747c6ae99SBarry Smith 
7387087cfbeSBarry Smith extern PetscErrorCode  VecView_MPI(Vec,PetscViewer);
73947c6ae99SBarry Smith EXTERN_C_BEGIN
74047c6ae99SBarry Smith #undef __FUNCT__
74147c6ae99SBarry Smith #define __FUNCT__ "VecView_DMComposite"
7427087cfbeSBarry Smith PetscErrorCode  VecView_DMComposite(Vec gvec,PetscViewer viewer)
74347c6ae99SBarry Smith {
74447c6ae99SBarry Smith   DM                     dm;
74547c6ae99SBarry Smith   PetscErrorCode         ierr;
74647c6ae99SBarry Smith   struct DMCompositeLink *next;
74747c6ae99SBarry Smith   PetscBool              isdraw;
748cef07954SSatish Balay   DM_Composite           *com;
74947c6ae99SBarry Smith 
75047c6ae99SBarry Smith   PetscFunctionBegin;
75147c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)gvec,"DMComposite",(PetscObject*)&dm);CHKERRQ(ierr);
75247c6ae99SBarry Smith   if (!dm) SETERRQ(((PetscObject)gvec)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMComposite");
75347c6ae99SBarry Smith   com = (DM_Composite*)dm->data;
75447c6ae99SBarry Smith   next = com->next;
75547c6ae99SBarry Smith 
75647c6ae99SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
75747c6ae99SBarry Smith   if (!isdraw) {
75847c6ae99SBarry Smith     /* do I really want to call this? */
75947c6ae99SBarry Smith     ierr = VecView_MPI(gvec,viewer);CHKERRQ(ierr);
76047c6ae99SBarry Smith   } else {
76147c6ae99SBarry Smith     PetscInt cnt = 0;
76247c6ae99SBarry Smith 
76347c6ae99SBarry Smith     /* loop over packed objects, handling one at at time */
76447c6ae99SBarry Smith     while (next) {
76547c6ae99SBarry Smith       if (next->type == DMCOMPOSITE_ARRAY) {
76647c6ae99SBarry Smith 	PetscScalar *array;
76747c6ae99SBarry Smith 	ierr  = DMCompositeGetAccess_Array(dm,next,gvec,&array);CHKERRQ(ierr);
76847c6ae99SBarry Smith 
76947c6ae99SBarry Smith 	/*skip it for now */
77047c6ae99SBarry Smith       } else if (next->type == DMCOMPOSITE_DM) {
77147c6ae99SBarry Smith 	Vec      vec;
77247c6ae99SBarry Smith         PetscInt bs;
77347c6ae99SBarry Smith 
77447c6ae99SBarry Smith 	ierr = DMCompositeGetAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
77547c6ae99SBarry Smith 	ierr = VecView(vec,viewer);CHKERRQ(ierr);
77647c6ae99SBarry Smith         ierr = VecGetBlockSize(vec,&bs);CHKERRQ(ierr);
77747c6ae99SBarry Smith 	ierr = DMCompositeRestoreAccess_DM(dm,next,gvec,&vec);CHKERRQ(ierr);
77847c6ae99SBarry Smith         ierr = PetscViewerDrawBaseAdd(viewer,bs);CHKERRQ(ierr);
77947c6ae99SBarry Smith         cnt += bs;
78047c6ae99SBarry Smith       } else {
78147c6ae99SBarry Smith 	SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
78247c6ae99SBarry Smith       }
78347c6ae99SBarry Smith       next = next->next;
78447c6ae99SBarry Smith     }
78547c6ae99SBarry Smith     ierr = PetscViewerDrawBaseAdd(viewer,-cnt);CHKERRQ(ierr);
78647c6ae99SBarry Smith   }
78747c6ae99SBarry Smith   PetscFunctionReturn(0);
78847c6ae99SBarry Smith }
78947c6ae99SBarry Smith EXTERN_C_END
79047c6ae99SBarry Smith 
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith #undef __FUNCT__
7930c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Composite"
7947087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector_Composite(DM dm,Vec *gvec)
79547c6ae99SBarry Smith {
79647c6ae99SBarry Smith   PetscErrorCode         ierr;
79747c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
79847c6ae99SBarry Smith 
79947c6ae99SBarry Smith   PetscFunctionBegin;
80047c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
80147c6ae99SBarry Smith   if (!com->setup) {
802d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
80347c6ae99SBarry Smith   }
80447c6ae99SBarry Smith   ierr = VecCreateMPI(((PetscObject)dm)->comm,com->n,com->N,gvec);CHKERRQ(ierr);
80547c6ae99SBarry Smith   ierr = PetscObjectCompose((PetscObject)*gvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr);
80647c6ae99SBarry Smith   ierr = VecSetOperation(*gvec,VECOP_VIEW,(void(*)(void))VecView_DMComposite);CHKERRQ(ierr);
80747c6ae99SBarry Smith   PetscFunctionReturn(0);
80847c6ae99SBarry Smith }
80947c6ae99SBarry Smith 
81047c6ae99SBarry Smith #undef __FUNCT__
8110c010503SBarry Smith #define __FUNCT__ "DMCreateLocalVector_Composite"
8127087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_Composite(DM dm,Vec *lvec)
81347c6ae99SBarry Smith {
81447c6ae99SBarry Smith   PetscErrorCode         ierr;
81547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
81647c6ae99SBarry Smith 
81747c6ae99SBarry Smith   PetscFunctionBegin;
81847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
81947c6ae99SBarry Smith   if (!com->setup) {
820d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
82147c6ae99SBarry Smith   }
82247c6ae99SBarry Smith   ierr = VecCreateSeq(((PetscObject)dm)->comm,com->nghost,lvec);CHKERRQ(ierr);
82347c6ae99SBarry Smith   ierr = PetscObjectCompose((PetscObject)*lvec,"DMComposite",(PetscObject)dm);CHKERRQ(ierr);
82447c6ae99SBarry Smith   PetscFunctionReturn(0);
82547c6ae99SBarry Smith }
82647c6ae99SBarry Smith 
82747c6ae99SBarry Smith #undef __FUNCT__
8286eb61c8cSJed Brown #define __FUNCT__ "DMCompositeGetISLocalToGlobalMappings"
82947c6ae99SBarry Smith /*@C
8306eb61c8cSJed Brown     DMCompositeGetISLocalToGlobalMappings - gets an ISLocalToGlobalMapping for each DM/array in the DMComposite, maps to the composite global space
83147c6ae99SBarry Smith 
83206ebdd98SJed Brown     Collective on DM
83347c6ae99SBarry Smith 
83447c6ae99SBarry Smith     Input Parameter:
83547c6ae99SBarry Smith .    dm - the packer object
83647c6ae99SBarry Smith 
83747c6ae99SBarry Smith     Output Parameters:
8386eb61c8cSJed Brown .    ltogs - the individual mappings for each packed vector/array. Note that this includes
839aa219208SBarry Smith            all the ghost points that individual ghosted DMDA's may have. Also each process has an
8406eb61c8cSJed Brown            mapping for EACH redundant array (not just the local redundant arrays).
84147c6ae99SBarry Smith 
84247c6ae99SBarry Smith     Level: advanced
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith     Notes:
8456eb61c8cSJed Brown        Each entry of ltogs should be destroyed with ISLocalToGlobalMappingDestroy(), the ltogs array should be freed with PetscFree().
84647c6ae99SBarry Smith 
8470c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
84847c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
84947c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
85047c6ae99SBarry Smith 
85147c6ae99SBarry Smith @*/
8527087cfbeSBarry Smith PetscErrorCode  DMCompositeGetISLocalToGlobalMappings(DM dm,ISLocalToGlobalMapping **ltogs)
85347c6ae99SBarry Smith {
85447c6ae99SBarry Smith   PetscErrorCode         ierr;
85547c6ae99SBarry Smith   PetscInt               i,*idx,n,cnt;
85647c6ae99SBarry Smith   struct DMCompositeLink *next;
85747c6ae99SBarry Smith   PetscMPIInt            rank;
85847c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith   PetscFunctionBegin;
86147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
862acc1e270SJed Brown   ierr = PetscMalloc((com->nDM+com->nredundant)*sizeof(ISLocalToGlobalMapping),ltogs);CHKERRQ(ierr);
86347c6ae99SBarry Smith   next = com->next;
86447c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
86547c6ae99SBarry Smith 
86647c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
86747c6ae99SBarry Smith   cnt = 0;
86847c6ae99SBarry Smith   while (next) {
86947c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
87006ebdd98SJed Brown       ierr = PetscMalloc(next->nlocal*sizeof(PetscInt),&idx);CHKERRQ(ierr);
87147c6ae99SBarry Smith       if (rank == next->rank) {
87206ebdd98SJed Brown         for (i=0; i<next->nlocal; i++) idx[i] = next->grstart + i;
87347c6ae99SBarry Smith       }
87406ebdd98SJed Brown       ierr = MPI_Bcast(idx,next->nlocal,MPIU_INT,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
87506ebdd98SJed Brown       ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,next->nlocal,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
87647c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
8776eb61c8cSJed Brown       ISLocalToGlobalMapping ltog;
8786eb61c8cSJed Brown       PetscMPIInt            size;
87986994e45SJed Brown       const PetscInt         *suboff,*indices;
8806eb61c8cSJed Brown       Vec                    global;
88147c6ae99SBarry Smith 
8826eb61c8cSJed Brown       /* Get sub-DM global indices for each local dof */
8831411c6eeSJed Brown       ierr = DMGetLocalToGlobalMapping(next->dm,&ltog);CHKERRQ(ierr);
8846eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingGetSize(ltog,&n);CHKERRQ(ierr);
88586994e45SJed Brown       ierr = ISLocalToGlobalMappingGetIndices(ltog,&indices);CHKERRQ(ierr);
88647c6ae99SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
88747c6ae99SBarry Smith 
8886eb61c8cSJed Brown       /* Get the offsets for the sub-DM global vector */
8896eb61c8cSJed Brown       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
8906eb61c8cSJed Brown       ierr = VecGetOwnershipRanges(global,&suboff);CHKERRQ(ierr);
8916eb61c8cSJed Brown       ierr = MPI_Comm_size(((PetscObject)global)->comm,&size);CHKERRQ(ierr);
8926eb61c8cSJed Brown 
8936eb61c8cSJed Brown       /* Shift the sub-DM definition of the global space to the composite global space */
8946eb61c8cSJed Brown       for (i=0; i<n; i++) {
89586994e45SJed Brown         PetscInt subi = indices[i],lo = 0,hi = size,t;
8966eb61c8cSJed Brown         /* Binary search to find which rank owns subi */
8976eb61c8cSJed Brown         while (hi-lo > 1) {
8986eb61c8cSJed Brown           t = lo + (hi-lo)/2;
8996eb61c8cSJed Brown           if (suboff[t] > subi) hi = t;
9006eb61c8cSJed Brown           else                  lo = t;
9016eb61c8cSJed Brown         }
9026eb61c8cSJed Brown         idx[i] = subi - suboff[lo] + next->grstarts[lo];
9036eb61c8cSJed Brown       }
90486994e45SJed Brown       ierr = ISLocalToGlobalMappingRestoreIndices(ltog,&indices);CHKERRQ(ierr);
9056eb61c8cSJed Brown       ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,n,idx,PETSC_OWN_POINTER,&(*ltogs)[cnt]);CHKERRQ(ierr);
9066eb61c8cSJed Brown       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
9076eb61c8cSJed Brown     } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
90847c6ae99SBarry Smith     next = next->next;
90947c6ae99SBarry Smith     cnt++;
91047c6ae99SBarry Smith   }
91147c6ae99SBarry Smith   PetscFunctionReturn(0);
91247c6ae99SBarry Smith }
91347c6ae99SBarry Smith 
91447c6ae99SBarry Smith #undef __FUNCT__
91587c85e80SJed Brown #define __FUNCT__ "DMCompositeGetLocalISs"
91687c85e80SJed Brown /*@C
91787c85e80SJed Brown    DMCompositeGetLocalISs - Gets index sets for each DM/array component of a composite local vector
91887c85e80SJed Brown 
91987c85e80SJed Brown    Not Collective
92087c85e80SJed Brown 
92187c85e80SJed Brown    Input Arguments:
92287c85e80SJed Brown . dm - composite DM
92387c85e80SJed Brown 
92487c85e80SJed Brown    Output Arguments:
92587c85e80SJed Brown . is - array of serial index sets for each each component of the DMComposite
92687c85e80SJed Brown 
92787c85e80SJed Brown    Level: intermediate
92887c85e80SJed Brown 
92987c85e80SJed Brown    Notes:
93087c85e80SJed Brown    At present, a composite local vector does not normally exist.  This function is used to provide index sets for
93187c85e80SJed Brown    MatGetLocalSubMatrix().  In the future, the scatters for each entry in the DMComposite may be be merged into a single
93287c85e80SJed Brown    scatter to a composite local vector.
93387c85e80SJed Brown 
93487c85e80SJed Brown    To get the composite global indices at all local points (including ghosts), use DMCompositeGetISLocalToGlobalMappings().
93587c85e80SJed Brown 
93687c85e80SJed Brown    To get index sets for pieces of the composite global vector, use DMCompositeGetGlobalISs().
93787c85e80SJed Brown 
93887c85e80SJed Brown    Each returned IS should be destroyed with ISDestroy(), the array should be freed with PetscFree().
93987c85e80SJed Brown 
94087c85e80SJed Brown .seealso: DMCompositeGetGlobalISs(), DMCompositeGetISLocalToGlobalMappings(), MatGetLocalSubMatrix(), MatCreateLocalRef()
94187c85e80SJed Brown @*/
9427087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalISs(DM dm,IS **is)
94387c85e80SJed Brown {
94487c85e80SJed Brown   PetscErrorCode         ierr;
94587c85e80SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
94687c85e80SJed Brown   struct DMCompositeLink *link;
94787c85e80SJed Brown   PetscInt cnt,start;
94887c85e80SJed Brown 
94987c85e80SJed Brown   PetscFunctionBegin;
95087c85e80SJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
95187c85e80SJed Brown   PetscValidPointer(is,2);
95287c85e80SJed Brown   ierr = PetscMalloc(com->nmine*sizeof(IS),is);CHKERRQ(ierr);
95306ebdd98SJed Brown   for (cnt=0,start=0,link=com->next; link; start+=link->nlocal,cnt++,link=link->next) {
95406ebdd98SJed Brown     ierr = ISCreateStride(PETSC_COMM_SELF,link->nlocal,start,1,&(*is)[cnt]);CHKERRQ(ierr);
955520db06cSJed Brown     if (link->type == DMCOMPOSITE_DM) {
956520db06cSJed Brown       PetscInt bs;
9571411c6eeSJed Brown       ierr = DMGetBlockSize(link->dm,&bs);CHKERRQ(ierr);
958520db06cSJed Brown       ierr = ISSetBlockSize((*is)[cnt],bs);CHKERRQ(ierr);
959520db06cSJed Brown     }
96087c85e80SJed Brown   }
96187c85e80SJed Brown   PetscFunctionReturn(0);
96287c85e80SJed Brown }
96387c85e80SJed Brown 
96487c85e80SJed Brown #undef __FUNCT__
96547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetGlobalISs"
96647c6ae99SBarry Smith /*@C
96747c6ae99SBarry Smith     DMCompositeGetGlobalISs - Gets the index sets for each composed object
96847c6ae99SBarry Smith 
96947c6ae99SBarry Smith     Collective on DMComposite
97047c6ae99SBarry Smith 
97147c6ae99SBarry Smith     Input Parameter:
97247c6ae99SBarry Smith .    dm - the packer object
97347c6ae99SBarry Smith 
97447c6ae99SBarry Smith     Output Parameters:
97547c6ae99SBarry Smith .    is - the array of index sets
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith     Level: advanced
97847c6ae99SBarry Smith 
97947c6ae99SBarry Smith     Notes:
98047c6ae99SBarry Smith        The is entries should be destroyed with ISDestroy(), the is array should be freed with PetscFree()
98147c6ae99SBarry Smith 
98247c6ae99SBarry Smith        These could be used to extract a subset of vector entries for a "multi-physics" preconditioner
98347c6ae99SBarry Smith 
9846eb61c8cSJed Brown        Use DMCompositeGetLocalISs() for index sets in the packed local numbering, and
9856eb61c8cSJed Brown        DMCompositeGetISLocalToGlobalMappings() for to map local sub-DM (including ghost) indices to packed global
9866eb61c8cSJed Brown        indices.
98747c6ae99SBarry Smith 
9880c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
98947c6ae99SBarry Smith          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetAccess(), DMCompositeScatter(),
99047c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(),DMCompositeGetEntries()
99147c6ae99SBarry Smith 
99247c6ae99SBarry Smith @*/
9936eb61c8cSJed Brown 
9947087cfbeSBarry Smith PetscErrorCode  DMCompositeGetGlobalISs(DM dm,IS *is[])
99547c6ae99SBarry Smith {
99647c6ae99SBarry Smith   PetscErrorCode         ierr;
99747c6ae99SBarry Smith   PetscInt               cnt = 0,*idx,i;
99847c6ae99SBarry Smith   struct DMCompositeLink *next;
99947c6ae99SBarry Smith   PetscMPIInt            rank;
100047c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
100147c6ae99SBarry Smith 
100247c6ae99SBarry Smith   PetscFunctionBegin;
100347c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
100406ebdd98SJed Brown   ierr = PetscMalloc((com->nDM+com->nredundant)*sizeof(IS),is);CHKERRQ(ierr);
100547c6ae99SBarry Smith   next = com->next;
100647c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
100747c6ae99SBarry Smith 
100847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
100947c6ae99SBarry Smith   while (next) {
101047c6ae99SBarry Smith     ierr = PetscMalloc(next->n*sizeof(PetscInt),&idx);CHKERRQ(ierr);
101147c6ae99SBarry Smith     for (i=0; i<next->n; i++) idx[i] = next->grstart + i;
101247c6ae99SBarry Smith     ierr = ISCreateGeneral(((PetscObject)dm)->comm,next->n,idx,PETSC_OWN_POINTER,&(*is)[cnt]);CHKERRQ(ierr);
101347c6ae99SBarry Smith     cnt++;
101447c6ae99SBarry Smith     next = next->next;
101547c6ae99SBarry Smith   }
101647c6ae99SBarry Smith   PetscFunctionReturn(0);
101747c6ae99SBarry Smith }
101847c6ae99SBarry Smith 
101947c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
102047c6ae99SBarry Smith #undef __FUNCT__
102147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_Array"
102247c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
102347c6ae99SBarry Smith {
102447c6ae99SBarry Smith   PetscErrorCode ierr;
102547c6ae99SBarry Smith   PetscFunctionBegin;
102647c6ae99SBarry Smith   if (array) {
1027c72eef85SJed Brown     ierr = PetscMalloc(mine->nlocal*sizeof(PetscScalar),array);CHKERRQ(ierr);
102847c6ae99SBarry Smith   }
102947c6ae99SBarry Smith   PetscFunctionReturn(0);
103047c6ae99SBarry Smith }
103147c6ae99SBarry Smith 
103247c6ae99SBarry Smith #undef __FUNCT__
103347c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors_DM"
103447c6ae99SBarry Smith PetscErrorCode DMCompositeGetLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
103547c6ae99SBarry Smith {
103647c6ae99SBarry Smith   PetscErrorCode ierr;
103747c6ae99SBarry Smith   PetscFunctionBegin;
103847c6ae99SBarry Smith   if (local) {
103947c6ae99SBarry Smith     ierr = DMGetLocalVector(mine->dm,local);CHKERRQ(ierr);
104047c6ae99SBarry Smith   }
104147c6ae99SBarry Smith   PetscFunctionReturn(0);
104247c6ae99SBarry Smith }
104347c6ae99SBarry Smith 
104447c6ae99SBarry Smith #undef __FUNCT__
104547c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_Array"
104647c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_Array(DM dm,struct DMCompositeLink *mine,PetscScalar **array)
104747c6ae99SBarry Smith {
104847c6ae99SBarry Smith   PetscErrorCode ierr;
104947c6ae99SBarry Smith   PetscFunctionBegin;
105047c6ae99SBarry Smith   if (array) {
105147c6ae99SBarry Smith     ierr = PetscFree(*array);CHKERRQ(ierr);
105247c6ae99SBarry Smith   }
105347c6ae99SBarry Smith   PetscFunctionReturn(0);
105447c6ae99SBarry Smith }
105547c6ae99SBarry Smith 
105647c6ae99SBarry Smith #undef __FUNCT__
105747c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors_DM"
105847c6ae99SBarry Smith PetscErrorCode DMCompositeRestoreLocalVectors_DM(DM dm,struct DMCompositeLink *mine,Vec *local)
105947c6ae99SBarry Smith {
106047c6ae99SBarry Smith   PetscErrorCode ierr;
106147c6ae99SBarry Smith   PetscFunctionBegin;
106247c6ae99SBarry Smith   if (local) {
106347c6ae99SBarry Smith     ierr = DMRestoreLocalVector(mine->dm,local);CHKERRQ(ierr);
106447c6ae99SBarry Smith   }
106547c6ae99SBarry Smith   PetscFunctionReturn(0);
106647c6ae99SBarry Smith }
106747c6ae99SBarry Smith 
106847c6ae99SBarry Smith #undef __FUNCT__
106947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetLocalVectors"
107047c6ae99SBarry Smith /*@C
107147c6ae99SBarry Smith     DMCompositeGetLocalVectors - Gets local vectors and arrays for each part of a DMComposite.
107247c6ae99SBarry Smith        Use DMCompositeRestoreLocalVectors() to return them.
107347c6ae99SBarry Smith 
107447c6ae99SBarry Smith     Not Collective
107547c6ae99SBarry Smith 
107647c6ae99SBarry Smith     Input Parameter:
107747c6ae99SBarry Smith .    dm - the packer object
107847c6ae99SBarry Smith 
107947c6ae99SBarry Smith     Output Parameter:
108047c6ae99SBarry Smith .    ... - the individual sequential objects (arrays or vectors)
108147c6ae99SBarry Smith 
108247c6ae99SBarry Smith     Level: advanced
108347c6ae99SBarry Smith 
10840c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
10856eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
108647c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
108747c6ae99SBarry Smith 
108847c6ae99SBarry Smith @*/
10897087cfbeSBarry Smith PetscErrorCode  DMCompositeGetLocalVectors(DM dm,...)
109047c6ae99SBarry Smith {
109147c6ae99SBarry Smith   va_list                Argp;
109247c6ae99SBarry Smith   PetscErrorCode         ierr;
109347c6ae99SBarry Smith   struct DMCompositeLink *next;
109447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
109547c6ae99SBarry Smith 
109647c6ae99SBarry Smith   PetscFunctionBegin;
109747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
109847c6ae99SBarry Smith   next = com->next;
109947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
110047c6ae99SBarry Smith   va_start(Argp,dm);
110147c6ae99SBarry Smith   while (next) {
110247c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
110347c6ae99SBarry Smith       PetscScalar **array;
110447c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
110547c6ae99SBarry Smith       ierr = DMCompositeGetLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
110647c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
110747c6ae99SBarry Smith       Vec *vec;
110847c6ae99SBarry Smith       vec = va_arg(Argp, Vec*);
110947c6ae99SBarry Smith       ierr = DMCompositeGetLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
111047c6ae99SBarry Smith     } else {
111147c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
111247c6ae99SBarry Smith     }
111347c6ae99SBarry Smith     next = next->next;
111447c6ae99SBarry Smith   }
111547c6ae99SBarry Smith   va_end(Argp);
111647c6ae99SBarry Smith   PetscFunctionReturn(0);
111747c6ae99SBarry Smith }
111847c6ae99SBarry Smith 
111947c6ae99SBarry Smith #undef __FUNCT__
112047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeRestoreLocalVectors"
112147c6ae99SBarry Smith /*@C
112247c6ae99SBarry Smith     DMCompositeRestoreLocalVectors - Restores local vectors and arrays for each part of a DMComposite.
112347c6ae99SBarry Smith 
112447c6ae99SBarry Smith     Not Collective
112547c6ae99SBarry Smith 
112647c6ae99SBarry Smith     Input Parameter:
112747c6ae99SBarry Smith .    dm - the packer object
112847c6ae99SBarry Smith 
112947c6ae99SBarry Smith     Output Parameter:
113047c6ae99SBarry Smith .    ... - the individual sequential objects (arrays or vectors)
113147c6ae99SBarry Smith 
113247c6ae99SBarry Smith     Level: advanced
113347c6ae99SBarry Smith 
11340c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
11356eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
113647c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeScatter(), DMCompositeGetEntries()
113747c6ae99SBarry Smith 
113847c6ae99SBarry Smith @*/
11397087cfbeSBarry Smith PetscErrorCode  DMCompositeRestoreLocalVectors(DM dm,...)
114047c6ae99SBarry Smith {
114147c6ae99SBarry Smith   va_list                Argp;
114247c6ae99SBarry Smith   PetscErrorCode         ierr;
114347c6ae99SBarry Smith   struct DMCompositeLink *next;
114447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
114547c6ae99SBarry Smith 
114647c6ae99SBarry Smith   PetscFunctionBegin;
114747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
114847c6ae99SBarry Smith   next = com->next;
114947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
115047c6ae99SBarry Smith   va_start(Argp,dm);
115147c6ae99SBarry Smith   while (next) {
115247c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
115347c6ae99SBarry Smith       PetscScalar **array;
115447c6ae99SBarry Smith       array = va_arg(Argp, PetscScalar**);
115547c6ae99SBarry Smith       ierr = DMCompositeRestoreLocalVectors_Array(dm,next,array);CHKERRQ(ierr);
115647c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
115747c6ae99SBarry Smith       Vec *vec;
115847c6ae99SBarry Smith       vec = va_arg(Argp, Vec*);
115947c6ae99SBarry Smith       ierr = DMCompositeRestoreLocalVectors_DM(dm,next,vec);CHKERRQ(ierr);
116047c6ae99SBarry Smith     } else {
116147c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
116247c6ae99SBarry Smith     }
116347c6ae99SBarry Smith     next = next->next;
116447c6ae99SBarry Smith   }
116547c6ae99SBarry Smith   va_end(Argp);
116647c6ae99SBarry Smith   PetscFunctionReturn(0);
116747c6ae99SBarry Smith }
116847c6ae99SBarry Smith 
116947c6ae99SBarry Smith /* -------------------------------------------------------------------------------------*/
117047c6ae99SBarry Smith #undef __FUNCT__
117147c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_Array"
117247c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_Array(DM dm,struct DMCompositeLink *mine,PetscInt *n)
117347c6ae99SBarry Smith {
117447c6ae99SBarry Smith   PetscFunctionBegin;
117506ebdd98SJed Brown   if (n) *n = mine->nlocal;
117647c6ae99SBarry Smith   PetscFunctionReturn(0);
117747c6ae99SBarry Smith }
117847c6ae99SBarry Smith 
117947c6ae99SBarry Smith #undef __FUNCT__
118047c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries_DM"
118147c6ae99SBarry Smith PetscErrorCode DMCompositeGetEntries_DM(DM dmi,struct DMCompositeLink *mine,DM *dm)
118247c6ae99SBarry Smith {
118347c6ae99SBarry Smith   PetscFunctionBegin;
118447c6ae99SBarry Smith   if (dm) *dm = mine->dm;
118547c6ae99SBarry Smith   PetscFunctionReturn(0);
118647c6ae99SBarry Smith }
118747c6ae99SBarry Smith 
118847c6ae99SBarry Smith #undef __FUNCT__
118947c6ae99SBarry Smith #define __FUNCT__ "DMCompositeGetEntries"
119047c6ae99SBarry Smith /*@C
1191aa219208SBarry Smith     DMCompositeGetEntries - Gets the DM, redundant size, etc for each entry in a DMComposite.
119247c6ae99SBarry Smith 
119347c6ae99SBarry Smith     Not Collective
119447c6ae99SBarry Smith 
119547c6ae99SBarry Smith     Input Parameter:
119647c6ae99SBarry Smith .    dm - the packer object
119747c6ae99SBarry Smith 
119847c6ae99SBarry Smith     Output Parameter:
1199aa219208SBarry Smith .    ... - the individual entries, DMs or integer sizes)
120047c6ae99SBarry Smith 
120147c6ae99SBarry Smith     Level: advanced
120247c6ae99SBarry Smith 
12030c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCreateGlobalVector(),
12046eb61c8cSJed Brown          DMCompositeGather(), DMCompositeCreate(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess(),
120547c6ae99SBarry Smith          DMCompositeRestoreLocalVectors(), DMCompositeGetLocalVectors(),  DMCompositeScatter(),
120647c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors()
120747c6ae99SBarry Smith 
120847c6ae99SBarry Smith @*/
12097087cfbeSBarry Smith PetscErrorCode  DMCompositeGetEntries(DM dm,...)
121047c6ae99SBarry Smith {
121147c6ae99SBarry Smith   va_list                Argp;
121247c6ae99SBarry Smith   PetscErrorCode         ierr;
121347c6ae99SBarry Smith   struct DMCompositeLink *next;
121447c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
121547c6ae99SBarry Smith 
121647c6ae99SBarry Smith   PetscFunctionBegin;
121747c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
121847c6ae99SBarry Smith   next = com->next;
121947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
122047c6ae99SBarry Smith   va_start(Argp,dm);
122147c6ae99SBarry Smith   while (next) {
122247c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
122347c6ae99SBarry Smith       PetscInt *n;
122447c6ae99SBarry Smith       n = va_arg(Argp, PetscInt*);
122547c6ae99SBarry Smith       ierr = DMCompositeGetEntries_Array(dm,next,n);CHKERRQ(ierr);
122647c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
122747c6ae99SBarry Smith       DM *dmn;
122847c6ae99SBarry Smith       dmn = va_arg(Argp, DM*);
122947c6ae99SBarry Smith       ierr = DMCompositeGetEntries_DM(dm,next,dmn);CHKERRQ(ierr);
123047c6ae99SBarry Smith     } else {
123147c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
123247c6ae99SBarry Smith     }
123347c6ae99SBarry Smith     next = next->next;
123447c6ae99SBarry Smith   }
123547c6ae99SBarry Smith   va_end(Argp);
123647c6ae99SBarry Smith   PetscFunctionReturn(0);
123747c6ae99SBarry Smith }
123847c6ae99SBarry Smith 
123947c6ae99SBarry Smith #undef __FUNCT__
12400c010503SBarry Smith #define __FUNCT__ "DMRefine_Composite"
12417087cfbeSBarry Smith PetscErrorCode  DMRefine_Composite(DM dmi,MPI_Comm comm,DM *fine)
124247c6ae99SBarry Smith {
124347c6ae99SBarry Smith   PetscErrorCode         ierr;
124447c6ae99SBarry Smith   struct DMCompositeLink *next;
124547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dmi->data;
124647c6ae99SBarry Smith   DM                     dm;
124747c6ae99SBarry Smith 
124847c6ae99SBarry Smith   PetscFunctionBegin;
124947c6ae99SBarry Smith   PetscValidHeaderSpecific(dmi,DM_CLASSID,1);
125047c6ae99SBarry Smith   next = com->next;
125147c6ae99SBarry Smith   ierr = DMCompositeCreate(comm,fine);CHKERRQ(ierr);
125247c6ae99SBarry Smith 
125347c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
125447c6ae99SBarry Smith   while (next) {
125547c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
125606ebdd98SJed Brown       ierr = DMCompositeAddArray(*fine,next->rank,next->nlocal);CHKERRQ(ierr);
125747c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
125847c6ae99SBarry Smith       ierr = DMRefine(next->dm,comm,&dm);CHKERRQ(ierr);
125947c6ae99SBarry Smith       ierr = DMCompositeAddDM(*fine,dm);CHKERRQ(ierr);
126047c6ae99SBarry Smith       ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr);
126147c6ae99SBarry Smith     } else {
126247c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
126347c6ae99SBarry Smith     }
126447c6ae99SBarry Smith     next = next->next;
126547c6ae99SBarry Smith   }
126647c6ae99SBarry Smith   PetscFunctionReturn(0);
126747c6ae99SBarry Smith }
126847c6ae99SBarry Smith 
126947c6ae99SBarry Smith #include "petscmat.h"
127047c6ae99SBarry Smith 
127147c6ae99SBarry Smith struct MatPackLink {
127247c6ae99SBarry Smith   Mat                A;
127347c6ae99SBarry Smith   struct MatPackLink *next;
127447c6ae99SBarry Smith };
127547c6ae99SBarry Smith 
127647c6ae99SBarry Smith struct MatPack {
127747c6ae99SBarry Smith   DM                 right,left;
127847c6ae99SBarry Smith   struct MatPackLink *next;
127947c6ae99SBarry Smith };
128047c6ae99SBarry Smith 
128147c6ae99SBarry Smith #undef __FUNCT__
128247c6ae99SBarry Smith #define __FUNCT__ "MatMultBoth_Shell_Pack"
128347c6ae99SBarry Smith PetscErrorCode MatMultBoth_Shell_Pack(Mat A,Vec x,Vec y,PetscBool  add)
128447c6ae99SBarry Smith {
128547c6ae99SBarry Smith   struct MatPack         *mpack;
128647c6ae99SBarry Smith   struct DMCompositeLink *xnext,*ynext;
128747c6ae99SBarry Smith   struct MatPackLink     *anext;
128847c6ae99SBarry Smith   PetscScalar            *xarray,*yarray;
128947c6ae99SBarry Smith   PetscErrorCode         ierr;
129047c6ae99SBarry Smith   PetscInt               i;
129147c6ae99SBarry Smith   Vec                    xglobal,yglobal;
129247c6ae99SBarry Smith   PetscMPIInt            rank;
129347c6ae99SBarry Smith   DM_Composite           *comright;
129447c6ae99SBarry Smith   DM_Composite           *comleft;
129547c6ae99SBarry Smith 
129647c6ae99SBarry Smith   PetscFunctionBegin;
129747c6ae99SBarry Smith   ierr = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
129847c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
129947c6ae99SBarry Smith   comright = (DM_Composite*)mpack->right->data;
130047c6ae99SBarry Smith   comleft = (DM_Composite*)mpack->left->data;
130147c6ae99SBarry Smith   xnext = comright->next;
130247c6ae99SBarry Smith   ynext = comleft->next;
130347c6ae99SBarry Smith   anext = mpack->next;
130447c6ae99SBarry Smith 
130547c6ae99SBarry Smith   while (xnext) {
130647c6ae99SBarry Smith     if (xnext->type == DMCOMPOSITE_ARRAY) {
130747c6ae99SBarry Smith       if (rank == xnext->rank) {
130847c6ae99SBarry Smith         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
130947c6ae99SBarry Smith         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
131047c6ae99SBarry Smith         if (add) {
131147c6ae99SBarry Smith           for (i=0; i<xnext->n; i++) {
131247c6ae99SBarry Smith             yarray[ynext->rstart+i] += xarray[xnext->rstart+i];
131347c6ae99SBarry Smith           }
131447c6ae99SBarry Smith         } else {
131547c6ae99SBarry Smith           ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
131647c6ae99SBarry Smith         }
131747c6ae99SBarry Smith         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
131847c6ae99SBarry Smith         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
131947c6ae99SBarry Smith       }
132047c6ae99SBarry Smith     } else if (xnext->type == DMCOMPOSITE_DM) {
132147c6ae99SBarry Smith       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
132247c6ae99SBarry Smith       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
132347c6ae99SBarry Smith       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
132447c6ae99SBarry Smith       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
132547c6ae99SBarry Smith       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
132647c6ae99SBarry Smith       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
132747c6ae99SBarry Smith       if (add) {
132847c6ae99SBarry Smith         ierr  = MatMultAdd(anext->A,xglobal,yglobal,yglobal);CHKERRQ(ierr);
132947c6ae99SBarry Smith       } else {
133047c6ae99SBarry Smith         ierr  = MatMult(anext->A,xglobal,yglobal);CHKERRQ(ierr);
133147c6ae99SBarry Smith       }
133247c6ae99SBarry Smith       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
133347c6ae99SBarry Smith       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
133447c6ae99SBarry Smith       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
133547c6ae99SBarry Smith       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
133647c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
133747c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
133847c6ae99SBarry Smith       anext = anext->next;
133947c6ae99SBarry Smith     } else {
134047c6ae99SBarry Smith       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
134147c6ae99SBarry Smith     }
134247c6ae99SBarry Smith     xnext = xnext->next;
134347c6ae99SBarry Smith     ynext = ynext->next;
134447c6ae99SBarry Smith   }
134547c6ae99SBarry Smith   PetscFunctionReturn(0);
134647c6ae99SBarry Smith }
134747c6ae99SBarry Smith 
134847c6ae99SBarry Smith #undef __FUNCT__
134947c6ae99SBarry Smith #define __FUNCT__ "MatMultAdd_Shell_Pack"
135047c6ae99SBarry Smith PetscErrorCode MatMultAdd_Shell_Pack(Mat A,Vec x,Vec y,Vec z)
135147c6ae99SBarry Smith {
135247c6ae99SBarry Smith   PetscErrorCode ierr;
135347c6ae99SBarry Smith   PetscFunctionBegin;
135447c6ae99SBarry Smith   if (z != y) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Handles y == z only");
135547c6ae99SBarry Smith   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
135647c6ae99SBarry Smith   PetscFunctionReturn(0);
135747c6ae99SBarry Smith }
135847c6ae99SBarry Smith 
135947c6ae99SBarry Smith #undef __FUNCT__
136047c6ae99SBarry Smith #define __FUNCT__ "MatMult_Shell_Pack"
136147c6ae99SBarry Smith PetscErrorCode MatMult_Shell_Pack(Mat A,Vec x,Vec y)
136247c6ae99SBarry Smith {
136347c6ae99SBarry Smith   PetscErrorCode ierr;
136447c6ae99SBarry Smith   PetscFunctionBegin;
136547c6ae99SBarry Smith   ierr = MatMultBoth_Shell_Pack(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
136647c6ae99SBarry Smith   PetscFunctionReturn(0);
136747c6ae99SBarry Smith }
136847c6ae99SBarry Smith 
136947c6ae99SBarry Smith #undef __FUNCT__
137047c6ae99SBarry Smith #define __FUNCT__ "MatMultTranspose_Shell_Pack"
137147c6ae99SBarry Smith PetscErrorCode MatMultTranspose_Shell_Pack(Mat A,Vec x,Vec y)
137247c6ae99SBarry Smith {
137347c6ae99SBarry Smith   struct MatPack         *mpack;
137447c6ae99SBarry Smith   struct DMCompositeLink *xnext,*ynext;
137547c6ae99SBarry Smith   struct MatPackLink     *anext;
137647c6ae99SBarry Smith   PetscScalar            *xarray,*yarray;
137747c6ae99SBarry Smith   PetscErrorCode         ierr;
137847c6ae99SBarry Smith   Vec                    xglobal,yglobal;
137947c6ae99SBarry Smith   PetscMPIInt            rank;
138047c6ae99SBarry Smith   DM_Composite           *comright;
138147c6ae99SBarry Smith   DM_Composite           *comleft;
138247c6ae99SBarry Smith 
138347c6ae99SBarry Smith   PetscFunctionBegin;
138447c6ae99SBarry Smith   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
138547c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)mpack->right)->comm,&rank);CHKERRQ(ierr);
138647c6ae99SBarry Smith   comright = (DM_Composite*)mpack->right->data;
138747c6ae99SBarry Smith   comleft = (DM_Composite*)mpack->left->data;
138847c6ae99SBarry Smith   ynext = comright->next;
138947c6ae99SBarry Smith   xnext = comleft->next;
139047c6ae99SBarry Smith   anext = mpack->next;
139147c6ae99SBarry Smith 
139247c6ae99SBarry Smith   while (xnext) {
139347c6ae99SBarry Smith     if (xnext->type == DMCOMPOSITE_ARRAY) {
139447c6ae99SBarry Smith       if (rank == ynext->rank) {
139547c6ae99SBarry Smith         ierr    = VecGetArray(x,&xarray);CHKERRQ(ierr);
139647c6ae99SBarry Smith         ierr    = VecGetArray(y,&yarray);CHKERRQ(ierr);
139747c6ae99SBarry Smith         ierr    = PetscMemcpy(yarray+ynext->rstart,xarray+xnext->rstart,xnext->n*sizeof(PetscScalar));CHKERRQ(ierr);
139847c6ae99SBarry Smith         ierr    = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
139947c6ae99SBarry Smith         ierr    = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
140047c6ae99SBarry Smith       }
140147c6ae99SBarry Smith     } else if (xnext->type == DMCOMPOSITE_DM) {
140247c6ae99SBarry Smith       ierr  = VecGetArray(x,&xarray);CHKERRQ(ierr);
140347c6ae99SBarry Smith       ierr  = VecGetArray(y,&yarray);CHKERRQ(ierr);
140447c6ae99SBarry Smith       ierr  = DMGetGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
140547c6ae99SBarry Smith       ierr  = DMGetGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
140647c6ae99SBarry Smith       ierr  = VecPlaceArray(xglobal,xarray+xnext->rstart);CHKERRQ(ierr);
140747c6ae99SBarry Smith       ierr  = VecPlaceArray(yglobal,yarray+ynext->rstart);CHKERRQ(ierr);
140847c6ae99SBarry Smith       ierr  = MatMultTranspose(anext->A,xglobal,yglobal);CHKERRQ(ierr);
140947c6ae99SBarry Smith       ierr  = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
141047c6ae99SBarry Smith       ierr  = VecRestoreArray(y,&yarray);CHKERRQ(ierr);
141147c6ae99SBarry Smith       ierr  = VecResetArray(xglobal);CHKERRQ(ierr);
141247c6ae99SBarry Smith       ierr  = VecResetArray(yglobal);CHKERRQ(ierr);
141347c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(xnext->dm,&xglobal);CHKERRQ(ierr);
141447c6ae99SBarry Smith       ierr  = DMRestoreGlobalVector(ynext->dm,&yglobal);CHKERRQ(ierr);
141547c6ae99SBarry Smith       anext = anext->next;
141647c6ae99SBarry Smith     } else {
141747c6ae99SBarry Smith       SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
141847c6ae99SBarry Smith     }
141947c6ae99SBarry Smith     xnext = xnext->next;
142047c6ae99SBarry Smith     ynext = ynext->next;
142147c6ae99SBarry Smith   }
142247c6ae99SBarry Smith   PetscFunctionReturn(0);
142347c6ae99SBarry Smith }
142447c6ae99SBarry Smith 
142547c6ae99SBarry Smith #undef __FUNCT__
142647c6ae99SBarry Smith #define __FUNCT__ "MatDestroy_Shell_Pack"
142747c6ae99SBarry Smith PetscErrorCode MatDestroy_Shell_Pack(Mat A)
142847c6ae99SBarry Smith {
142947c6ae99SBarry Smith   struct MatPack     *mpack;
143047c6ae99SBarry Smith   struct MatPackLink *anext,*oldanext;
143147c6ae99SBarry Smith   PetscErrorCode     ierr;
143247c6ae99SBarry Smith 
143347c6ae99SBarry Smith   PetscFunctionBegin;
143447c6ae99SBarry Smith   ierr  = MatShellGetContext(A,(void**)&mpack);CHKERRQ(ierr);
143547c6ae99SBarry Smith   anext = mpack->next;
143647c6ae99SBarry Smith 
143747c6ae99SBarry Smith   while (anext) {
143847c6ae99SBarry Smith     ierr     = MatDestroy(anext->A);CHKERRQ(ierr);
143947c6ae99SBarry Smith     oldanext = anext;
144047c6ae99SBarry Smith     anext    = anext->next;
144147c6ae99SBarry Smith     ierr     = PetscFree(oldanext);CHKERRQ(ierr);
144247c6ae99SBarry Smith   }
144347c6ae99SBarry Smith   ierr = PetscFree(mpack);CHKERRQ(ierr);
144447c6ae99SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
144547c6ae99SBarry Smith   PetscFunctionReturn(0);
144647c6ae99SBarry Smith }
144747c6ae99SBarry Smith 
144847c6ae99SBarry Smith #undef __FUNCT__
14490c010503SBarry Smith #define __FUNCT__ "DMGetInterpolation_Composite"
14507087cfbeSBarry Smith PetscErrorCode  DMGetInterpolation_Composite(DM coarse,DM fine,Mat *A,Vec *v)
145147c6ae99SBarry Smith {
145247c6ae99SBarry Smith   PetscErrorCode         ierr;
145347c6ae99SBarry Smith   PetscInt               m,n,M,N;
145447c6ae99SBarry Smith   struct DMCompositeLink *nextc;
145547c6ae99SBarry Smith   struct DMCompositeLink *nextf;
145647c6ae99SBarry Smith   struct MatPackLink     *nextmat,*pnextmat = 0;
145747c6ae99SBarry Smith   struct MatPack         *mpack;
145847c6ae99SBarry Smith   Vec                    gcoarse,gfine;
145947c6ae99SBarry Smith   DM_Composite           *comcoarse = (DM_Composite*)coarse->data;
146047c6ae99SBarry Smith   DM_Composite           *comfine = (DM_Composite*)fine->data;
146147c6ae99SBarry Smith 
146247c6ae99SBarry Smith   PetscFunctionBegin;
146347c6ae99SBarry Smith   PetscValidHeaderSpecific(coarse,DM_CLASSID,1);
146447c6ae99SBarry Smith   PetscValidHeaderSpecific(fine,DM_CLASSID,2);
146547c6ae99SBarry Smith   nextc = comcoarse->next;
146647c6ae99SBarry Smith   nextf = comfine->next;
146747c6ae99SBarry Smith   /* use global vectors only for determining matrix layout */
14680c010503SBarry Smith   ierr = DMCreateGlobalVector(coarse,&gcoarse);CHKERRQ(ierr);
14690c010503SBarry Smith   ierr = DMCreateGlobalVector(fine,&gfine);CHKERRQ(ierr);
147047c6ae99SBarry Smith   ierr = VecGetLocalSize(gcoarse,&n);CHKERRQ(ierr);
147147c6ae99SBarry Smith   ierr = VecGetLocalSize(gfine,&m);CHKERRQ(ierr);
147247c6ae99SBarry Smith   ierr = VecGetSize(gcoarse,&N);CHKERRQ(ierr);
147347c6ae99SBarry Smith   ierr = VecGetSize(gfine,&M);CHKERRQ(ierr);
147447c6ae99SBarry Smith   ierr = VecDestroy(gcoarse);CHKERRQ(ierr);
147547c6ae99SBarry Smith   ierr = VecDestroy(gfine);CHKERRQ(ierr);
147647c6ae99SBarry Smith 
147747c6ae99SBarry Smith   ierr         = PetscNew(struct MatPack,&mpack);CHKERRQ(ierr);
147847c6ae99SBarry Smith   mpack->right = coarse;
147947c6ae99SBarry Smith   mpack->left  = fine;
148047c6ae99SBarry Smith   ierr  = MatCreate(((PetscObject)fine)->comm,A);CHKERRQ(ierr);
148147c6ae99SBarry Smith   ierr  = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
148247c6ae99SBarry Smith   ierr  = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
148347c6ae99SBarry Smith   ierr  = MatShellSetContext(*A,mpack);CHKERRQ(ierr);
148447c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT,(void(*)(void))MatMult_Shell_Pack);CHKERRQ(ierr);
148547c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Shell_Pack);CHKERRQ(ierr);
148647c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_Shell_Pack);CHKERRQ(ierr);
148747c6ae99SBarry Smith   ierr  = MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_Shell_Pack);CHKERRQ(ierr);
148847c6ae99SBarry Smith 
148947c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
149047c6ae99SBarry Smith   while (nextc) {
149147c6ae99SBarry Smith     if (nextc->type != nextf->type) SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_ARG_INCOMP,"Two DM have different layout");
149247c6ae99SBarry Smith 
149347c6ae99SBarry Smith     if (nextc->type == DMCOMPOSITE_ARRAY) {
149447c6ae99SBarry Smith       ;
149547c6ae99SBarry Smith     } else if (nextc->type == DMCOMPOSITE_DM) {
149647c6ae99SBarry Smith       ierr          = PetscNew(struct MatPackLink,&nextmat);CHKERRQ(ierr);
149747c6ae99SBarry Smith       nextmat->next = 0;
149847c6ae99SBarry Smith       if (pnextmat) {
149947c6ae99SBarry Smith         pnextmat->next = nextmat;
150047c6ae99SBarry Smith         pnextmat       = nextmat;
150147c6ae99SBarry Smith       } else {
150247c6ae99SBarry Smith         pnextmat    = nextmat;
150347c6ae99SBarry Smith         mpack->next = nextmat;
150447c6ae99SBarry Smith       }
150547c6ae99SBarry Smith       ierr = DMGetInterpolation(nextc->dm,nextf->dm,&nextmat->A,PETSC_NULL);CHKERRQ(ierr);
150647c6ae99SBarry Smith     } else {
150747c6ae99SBarry Smith       SETERRQ(((PetscObject)fine)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
150847c6ae99SBarry Smith     }
150947c6ae99SBarry Smith     nextc = nextc->next;
151047c6ae99SBarry Smith     nextf = nextf->next;
151147c6ae99SBarry Smith   }
151247c6ae99SBarry Smith   PetscFunctionReturn(0);
151347c6ae99SBarry Smith }
151447c6ae99SBarry Smith 
151547c6ae99SBarry Smith #undef __FUNCT__
15161411c6eeSJed Brown #define __FUNCT__ "DMCreateLocalToGlobalMapping_Composite"
15171411c6eeSJed Brown static PetscErrorCode DMCreateLocalToGlobalMapping_Composite(DM dm)
15181411c6eeSJed Brown {
15191411c6eeSJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
15201411c6eeSJed Brown   ISLocalToGlobalMapping *ltogs;
1521a03f9cedSJed Brown   PetscInt               i,cnt,m,*idx;
15221411c6eeSJed Brown   PetscErrorCode         ierr;
15231411c6eeSJed Brown 
15241411c6eeSJed Brown   PetscFunctionBegin;
15251411c6eeSJed Brown   /* Set the ISLocalToGlobalMapping on the new matrix */
15261411c6eeSJed Brown   ierr = DMCompositeGetISLocalToGlobalMappings(dm,&ltogs);CHKERRQ(ierr);
15271411c6eeSJed Brown   for (cnt=0,i=0; i<(com->nDM+com->nredundant); i++) {
1528a03f9cedSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1529a03f9cedSJed Brown     cnt += m;
1530a03f9cedSJed Brown   }
1531a03f9cedSJed Brown   ierr = PetscMalloc(cnt*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1532a03f9cedSJed Brown   for (cnt=0,i=0; i<(com->nDM+com->nredundant); i++) {
15331411c6eeSJed Brown     const PetscInt *subidx;
15341411c6eeSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
15351411c6eeSJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
15361411c6eeSJed Brown     ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr);
15371411c6eeSJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
15381411c6eeSJed Brown     cnt += m;
15391411c6eeSJed Brown   }
15401411c6eeSJed Brown   ierr = ISLocalToGlobalMappingCreate(((PetscObject)dm)->comm,cnt,idx,PETSC_OWN_POINTER,&dm->ltogmap);CHKERRQ(ierr);
15411411c6eeSJed Brown   for (i=0; i<com->nDM+com->nredundant; i++) {ierr = ISLocalToGlobalMappingDestroy(ltogs[i]);CHKERRQ(ierr);}
15421411c6eeSJed Brown   ierr = PetscFree(ltogs);CHKERRQ(ierr);
15431411c6eeSJed Brown   PetscFunctionReturn(0);
15441411c6eeSJed Brown }
15451411c6eeSJed Brown 
15461411c6eeSJed Brown 
15471411c6eeSJed Brown #undef __FUNCT__
15480c010503SBarry Smith #define __FUNCT__ "DMGetColoring_Composite"
15497087cfbeSBarry Smith PetscErrorCode  DMGetColoring_Composite(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
155047c6ae99SBarry Smith {
155147c6ae99SBarry Smith   PetscErrorCode         ierr;
155247c6ae99SBarry Smith   PetscInt               n,i,cnt;
155347c6ae99SBarry Smith   ISColoringValue        *colors;
155447c6ae99SBarry Smith   PetscBool              dense = PETSC_FALSE;
155547c6ae99SBarry Smith   ISColoringValue        maxcol = 0;
155647c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
155747c6ae99SBarry Smith 
155847c6ae99SBarry Smith   PetscFunctionBegin;
155947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
156047c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED) {
156147c6ae99SBarry Smith     SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Currently you must use -dmmg_iscoloring_type global" );
156247c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL) {
156347c6ae99SBarry Smith     n = com->n;
156447c6ae99SBarry Smith   } else SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Unknown ISColoringType");
156547c6ae99SBarry Smith   ierr = PetscMalloc(n*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); /* freed in ISColoringDestroy() */
156647c6ae99SBarry Smith 
1567671f6225SBarry Smith   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
156847c6ae99SBarry Smith   if (dense) {
156947c6ae99SBarry Smith     for (i=0; i<n; i++) {
157047c6ae99SBarry Smith       colors[i] = (ISColoringValue)(com->rstart + i);
157147c6ae99SBarry Smith     }
157247c6ae99SBarry Smith     maxcol = com->N;
157347c6ae99SBarry Smith   } else {
157447c6ae99SBarry Smith     struct DMCompositeLink *next = com->next;
157547c6ae99SBarry Smith     PetscMPIInt            rank;
157647c6ae99SBarry Smith 
157747c6ae99SBarry Smith     ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
157847c6ae99SBarry Smith     cnt  = 0;
157947c6ae99SBarry Smith     while (next) {
158047c6ae99SBarry Smith       if (next->type == DMCOMPOSITE_ARRAY) {
158147c6ae99SBarry Smith         if (rank == next->rank) {  /* each column gets is own color */
158247c6ae99SBarry Smith           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
158347c6ae99SBarry Smith             colors[cnt++] = maxcol++;
158447c6ae99SBarry Smith           }
158547c6ae99SBarry Smith         }
158647c6ae99SBarry Smith         ierr = MPI_Bcast(&maxcol,1,MPIU_COLORING_VALUE,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
158747c6ae99SBarry Smith       } else if (next->type == DMCOMPOSITE_DM) {
158847c6ae99SBarry Smith         ISColoring     lcoloring;
158947c6ae99SBarry Smith 
159047c6ae99SBarry Smith         ierr = DMGetColoring(next->dm,IS_COLORING_GLOBAL,mtype,&lcoloring);CHKERRQ(ierr);
159147c6ae99SBarry Smith         for (i=0; i<lcoloring->N; i++) {
159247c6ae99SBarry Smith           colors[cnt++] = maxcol + lcoloring->colors[i];
159347c6ae99SBarry Smith         }
159447c6ae99SBarry Smith         maxcol += lcoloring->n;
159547c6ae99SBarry Smith         ierr = ISColoringDestroy(lcoloring);CHKERRQ(ierr);
159647c6ae99SBarry Smith       } else {
159747c6ae99SBarry Smith         SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
159847c6ae99SBarry Smith       }
159947c6ae99SBarry Smith       next = next->next;
160047c6ae99SBarry Smith     }
160147c6ae99SBarry Smith   }
160247c6ae99SBarry Smith   ierr = ISColoringCreate(((PetscObject)dm)->comm,maxcol,n,colors,coloring);CHKERRQ(ierr);
160347c6ae99SBarry Smith   PetscFunctionReturn(0);
160447c6ae99SBarry Smith }
160547c6ae99SBarry Smith 
160647c6ae99SBarry Smith #undef __FUNCT__
16070c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin_Composite"
16087087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalBegin_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
160947c6ae99SBarry Smith {
161047c6ae99SBarry Smith   PetscErrorCode         ierr;
161147c6ae99SBarry Smith   struct DMCompositeLink *next;
161247c6ae99SBarry Smith   PetscInt               cnt = 3;
161347c6ae99SBarry Smith   PetscMPIInt            rank;
161447c6ae99SBarry Smith   PetscScalar            *garray,*larray;
161547c6ae99SBarry Smith   DM_Composite           *com = (DM_Composite*)dm->data;
161647c6ae99SBarry Smith 
161747c6ae99SBarry Smith   PetscFunctionBegin;
161847c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
161947c6ae99SBarry Smith   PetscValidHeaderSpecific(gvec,VEC_CLASSID,2);
162047c6ae99SBarry Smith   next = com->next;
162147c6ae99SBarry Smith   if (!com->setup) {
1622d7bf68aeSBarry Smith     ierr = DMSetUp(dm);CHKERRQ(ierr);
162347c6ae99SBarry Smith   }
162447c6ae99SBarry Smith   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
162547c6ae99SBarry Smith   ierr = VecGetArray(gvec,&garray);CHKERRQ(ierr);
162647c6ae99SBarry Smith   ierr = VecGetArray(lvec,&larray);CHKERRQ(ierr);
162747c6ae99SBarry Smith 
162847c6ae99SBarry Smith   /* loop over packed objects, handling one at at time */
162947c6ae99SBarry Smith   while (next) {
163047c6ae99SBarry Smith     if (next->type == DMCOMPOSITE_ARRAY) {
163147c6ae99SBarry Smith       if (rank == next->rank) {
163247c6ae99SBarry Smith         ierr    = PetscMemcpy(larray,garray,next->n*sizeof(PetscScalar));CHKERRQ(ierr);
163347c6ae99SBarry Smith         garray += next->n;
163447c6ae99SBarry Smith       }
163547c6ae99SBarry Smith       /* does not handle ADD_VALUES */
163606ebdd98SJed Brown       ierr = MPI_Bcast(larray,next->nlocal,MPIU_SCALAR,next->rank,((PetscObject)dm)->comm);CHKERRQ(ierr);
163747c6ae99SBarry Smith     } else if (next->type == DMCOMPOSITE_DM) {
163847c6ae99SBarry Smith       Vec      local,global;
163947c6ae99SBarry Smith       PetscInt N;
164047c6ae99SBarry Smith 
164147c6ae99SBarry Smith       ierr = DMGetGlobalVector(next->dm,&global);CHKERRQ(ierr);
164247c6ae99SBarry Smith       ierr = VecGetLocalSize(global,&N);CHKERRQ(ierr);
164347c6ae99SBarry Smith       ierr = VecPlaceArray(global,garray);CHKERRQ(ierr);
164447c6ae99SBarry Smith       ierr = DMGetLocalVector(next->dm,&local);CHKERRQ(ierr);
164547c6ae99SBarry Smith       ierr = VecPlaceArray(local,larray);CHKERRQ(ierr);
164647c6ae99SBarry Smith       ierr = DMGlobalToLocalBegin(next->dm,global,mode,local);CHKERRQ(ierr);
164747c6ae99SBarry Smith       ierr = DMGlobalToLocalEnd(next->dm,global,mode,local);CHKERRQ(ierr);
164847c6ae99SBarry Smith       ierr = VecResetArray(global);CHKERRQ(ierr);
164947c6ae99SBarry Smith       ierr = VecResetArray(local);CHKERRQ(ierr);
165047c6ae99SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,&global);CHKERRQ(ierr);
165147c6ae99SBarry Smith       ierr = DMRestoreGlobalVector(next->dm,&local);CHKERRQ(ierr);
165247c6ae99SBarry Smith     } else {
165347c6ae99SBarry Smith       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
165447c6ae99SBarry Smith     }
165547c6ae99SBarry Smith     cnt++;
165606ebdd98SJed Brown     larray += next->nlocal;
165747c6ae99SBarry Smith     next    = next->next;
165847c6ae99SBarry Smith   }
165947c6ae99SBarry Smith 
166047c6ae99SBarry Smith   ierr = VecRestoreArray(gvec,PETSC_NULL);CHKERRQ(ierr);
166147c6ae99SBarry Smith   ierr = VecRestoreArray(lvec,PETSC_NULL);CHKERRQ(ierr);
166247c6ae99SBarry Smith   PetscFunctionReturn(0);
166347c6ae99SBarry Smith }
166447c6ae99SBarry Smith 
166547c6ae99SBarry Smith #undef __FUNCT__
16660c010503SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd_Composite"
16677087cfbeSBarry Smith PetscErrorCode  DMGlobalToLocalEnd_Composite(DM dm,Vec gvec,InsertMode mode,Vec lvec)
16680c010503SBarry Smith {
16690c010503SBarry Smith   PetscFunctionBegin;
16700c010503SBarry Smith   PetscFunctionReturn(0);
16710c010503SBarry Smith }
167247c6ae99SBarry Smith 
1673a4121054SBarry Smith EXTERN_C_BEGIN
1674a4121054SBarry Smith #undef __FUNCT__
1675a4121054SBarry Smith #define __FUNCT__ "DMCreate_Composite"
16767087cfbeSBarry Smith PetscErrorCode  DMCreate_Composite(DM p)
1677a4121054SBarry Smith {
1678a4121054SBarry Smith   PetscErrorCode ierr;
1679a4121054SBarry Smith   DM_Composite   *com;
1680a4121054SBarry Smith 
1681a4121054SBarry Smith   PetscFunctionBegin;
1682a4121054SBarry Smith   ierr = PetscNewLog(p,DM_Composite,&com);CHKERRQ(ierr);
1683a4121054SBarry Smith   p->data = com;
1684a4121054SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)p,"DMComposite");CHKERRQ(ierr);
1685a4121054SBarry Smith   com->n            = 0;
1686a4121054SBarry Smith   com->next         = PETSC_NULL;
1687a4121054SBarry Smith   com->nredundant   = 0;
1688a4121054SBarry Smith   com->nDM          = 0;
1689a4121054SBarry Smith 
1690a4121054SBarry Smith   p->ops->createglobalvector              = DMCreateGlobalVector_Composite;
1691a4121054SBarry Smith   p->ops->createlocalvector               = DMCreateLocalVector_Composite;
16921411c6eeSJed Brown   p->ops->createlocaltoglobalmapping      = DMCreateLocalToGlobalMapping_Composite;
16931411c6eeSJed Brown   p->ops->createlocaltoglobalmappingblock = 0;
1694a4121054SBarry Smith   p->ops->refine                          = DMRefine_Composite;
1695a4121054SBarry Smith   p->ops->getinterpolation                = DMGetInterpolation_Composite;
1696a4121054SBarry Smith   p->ops->getmatrix                       = DMGetMatrix_Composite;
1697a4121054SBarry Smith   p->ops->getcoloring                     = DMGetColoring_Composite;
1698a4121054SBarry Smith   p->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Composite;
1699a4121054SBarry Smith   p->ops->globaltolocalend                = DMGlobalToLocalEnd_Composite;
1700a4121054SBarry Smith   p->ops->destroy                         = DMDestroy_Composite;
1701a4121054SBarry Smith   p->ops->view                            = DMView_Composite;
1702a4121054SBarry Smith   p->ops->setup                           = DMSetUp_Composite;
1703a4121054SBarry Smith   PetscFunctionReturn(0);
1704a4121054SBarry Smith }
1705a4121054SBarry Smith EXTERN_C_END
1706a4121054SBarry Smith 
17070c010503SBarry Smith #undef __FUNCT__
17080c010503SBarry Smith #define __FUNCT__ "DMCompositeCreate"
17090c010503SBarry Smith /*@C
17100c010503SBarry Smith     DMCompositeCreate - Creates a vector packer, used to generate "composite"
17110c010503SBarry Smith       vectors made up of several subvectors.
17120c010503SBarry Smith 
17130c010503SBarry Smith     Collective on MPI_Comm
171447c6ae99SBarry Smith 
171547c6ae99SBarry Smith     Input Parameter:
17160c010503SBarry Smith .   comm - the processors that will share the global vector
17170c010503SBarry Smith 
17180c010503SBarry Smith     Output Parameters:
17190c010503SBarry Smith .   packer - the packer object
172047c6ae99SBarry Smith 
172147c6ae99SBarry Smith     Level: advanced
172247c6ae99SBarry Smith 
17230c010503SBarry Smith .seealso DMDestroy(), DMCompositeAddArray(), DMCompositeAddDM(), DMCompositeScatter(),
17246eb61c8cSJed Brown          DMCompositeGather(), DMCreateGlobalVector(), DMCompositeGetISLocalToGlobalMappings(), DMCompositeGetAccess()
172547c6ae99SBarry Smith          DMCompositeGetLocalVectors(), DMCompositeRestoreLocalVectors(), DMCompositeGetEntries()
172647c6ae99SBarry Smith 
172747c6ae99SBarry Smith @*/
17287087cfbeSBarry Smith PetscErrorCode  DMCompositeCreate(MPI_Comm comm,DM *packer)
172947c6ae99SBarry Smith {
17300c010503SBarry Smith   PetscErrorCode ierr;
17310c010503SBarry Smith 
173247c6ae99SBarry Smith   PetscFunctionBegin;
17330c010503SBarry Smith   PetscValidPointer(packer,2);
1734a4121054SBarry Smith   ierr = DMCreate(comm,packer);CHKERRQ(ierr);
1735a4121054SBarry Smith   ierr = DMSetType(*packer,DMCOMPOSITE);CHKERRQ(ierr);
173647c6ae99SBarry Smith   PetscFunctionReturn(0);
173747c6ae99SBarry Smith }
1738