Actual source code: ex9.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: static char help[] = "Demonstrates use of VecCreateGhost().\n\n";

  4: /*T
  5:    Concepts: vectors^assembling vectors;
  6:    Concepts: vectors^ghost padding;
  7:    Processors: n

  9:    Description: Ghost padding is one way to handle local calculations that
 10:       involve values from other processors. VecCreateGhost() provides
 11:       a way to create vectors with extra room at the end of the vector
 12:       array to contain the needed ghost values from other processors,
 13:       vector computations are otherwise unaffected.
 14: T*/

 16: /*
 17:   Include "petscvec.h" so that we can use vectors.  Note that this file
 18:   automatically includes:
 19:      petscsys.h       - base PETSc routines   petscis.h     - index sets
 20:      petscviewer.h - viewers
 21: */
 22: #include <petscvec.h>

 26: int main(int argc,char **argv)
 27: {
 28:   PetscMPIInt    rank,size;
 29:   PetscInt       nlocal = 6,nghost = 2,ifrom[2],i,rstart,rend;
 31:   PetscBool      flg,flg2;
 32:   PetscScalar    value,*array,*tarray=0;
 33:   Vec            lx,gx,gxs;

 35:   PetscInitialize(&argc,&argv,(char*)0,help);
 36:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 37:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 38:   if (size != 2) SETERRQ(PETSC_COMM_SELF,1,"Must run example with two processors\n");

 40:   /*
 41:      Construct a two dimensional graph connecting nlocal degrees of
 42:      freedom per processor. From this we will generate the global
 43:      indices of needed ghost values

 45:      For simplicity we generate the entire graph on each processor:
 46:      in real application the graph would stored in parallel, but this
 47:      example is only to demonstrate the management of ghost padding
 48:      with VecCreateGhost().

 50:      In this example we consider the vector as representing
 51:      degrees of freedom in a one dimensional grid with periodic
 52:      boundary conditions.

 54:         ----Processor  1---------  ----Processor 2 --------
 55:          0    1   2   3   4    5    6    7   8   9   10   11
 56:                                |----|
 57:          |-------------------------------------------------|

 59:   */

 61:   if (!rank) {
 62:     ifrom[0] = 11; ifrom[1] = 6;
 63:   } else {
 64:     ifrom[0] = 0;  ifrom[1] = 5;
 65:   }

 67:   /*
 68:      Create the vector with two slots for ghost points. Note that both
 69:      the local vector (lx) and the global vector (gx) share the same
 70:      array for storing vector values.
 71:   */
 72:   PetscOptionsHasName(NULL,"-allocate",&flg);
 73:   PetscOptionsHasName(NULL,"-vecmpisetghost",&flg2);
 74:   if (flg) {
 75:     PetscMalloc1((nlocal+nghost),&tarray);
 76:     VecCreateGhostWithArray(PETSC_COMM_WORLD,nlocal,PETSC_DECIDE,nghost,ifrom,tarray,&gxs);
 77:   } else if (flg2) {
 78:     VecCreate(PETSC_COMM_WORLD,&gxs);
 79:     VecSetType(gxs,VECMPI);
 80:     VecSetSizes(gxs,nlocal,PETSC_DECIDE);
 81:     VecMPISetGhost(gxs,nghost,ifrom);
 82:   } else {
 83:     VecCreateGhost(PETSC_COMM_WORLD,nlocal,PETSC_DECIDE,nghost,ifrom,&gxs);
 84:   }

 86:   /*
 87:       Test VecDuplicate()
 88:   */
 89:   VecDuplicate(gxs,&gx);
 90:   VecDestroy(&gxs);

 92:   /*
 93:      Access the local representation
 94:   */
 95:   VecGhostGetLocalForm(gx,&lx);

 97:   /*
 98:      Set the values from 0 to 12 into the "global" vector
 99:   */
100:   VecGetOwnershipRange(gx,&rstart,&rend);
101:   for (i=rstart; i<rend; i++) {
102:     value = (PetscScalar) i;
103:     VecSetValues(gx,1,&i,&value,INSERT_VALUES);
104:   }
105:   VecAssemblyBegin(gx);
106:   VecAssemblyEnd(gx);

108:   VecGhostUpdateBegin(gx,INSERT_VALUES,SCATTER_FORWARD);
109:   VecGhostUpdateEnd(gx,INSERT_VALUES,SCATTER_FORWARD);

111:   /*
112:      Print out each vector, including the ghost padding region.
113:   */
114:   VecGetArray(lx,&array);
115:   for (i=0; i<nlocal+nghost; i++) {
116:     PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%D %g\n",i,(double)PetscRealPart(array[i]));
117:   }
118:   VecRestoreArray(lx,&array);
119:   PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);

121:   VecGhostRestoreLocalForm(gx,&lx);
122:   VecDestroy(&gx);
123:   if (flg) {PetscFree(tarray);}
124:   PetscFinalize();
125:   return 0;
126: }