Actual source code: ex9.c

  1: /*$Id: ex9.c,v 1.32 2001/08/07 03:02:34 balay Exp $*/

  3: static char help[] = "Demonstrates use of VecCreateGhost().nn";

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

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

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

 25: int main(int argc,char **argv)
 26: {
 27:   int          rank,nlocal = 6,nghost = 2,ifrom[2],size,ierr,i,rstart,rend;
 28:   PetscTruth   flg;
 29:   PetscScalar  value,*array,*tarray=0;
 30:   Vec          lx,gx,gxs;

 32:   PetscInitialize(&argc,&argv,(char *)0,help);
 33:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 34:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 35:   if (size != 2) SETERRQ(1,"Must run example with two processorsn");

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

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

 47:      In this example we consider the vector as representing 
 48:      degrees of freedom in a one dimensional grid with periodic 
 49:      boundary conditions.

 51:         ----Processor  1---------  ----Processor 2 --------
 52:          0    1   2   3   4    5    6    7   8   9   10   11
 53:                                |----| 
 54:          |-------------------------------------------------|

 56:   */

 58:   if (!rank) {
 59:     ifrom[0] = 11; ifrom[1] = 6;
 60:   } else {
 61:     ifrom[0] = 0;  ifrom[1] = 5;
 62:   }

 64:   /*
 65:      Create the vector with two slots for ghost points. Note that both 
 66:      the local vector (lx) and the global vector (gx) share the same 
 67:      array for storing vector values.
 68:   */
 69:   PetscOptionsHasName(PETSC_NULL,"-allocate",&flg);
 70:   if (flg) {
 71:     PetscMalloc((nlocal+nghost)*sizeof(PetscScalar),&tarray);
 72:     VecCreateGhostWithArray(PETSC_COMM_WORLD,nlocal,PETSC_DECIDE,nghost,ifrom,tarray,&gxs);
 73:   } else {
 74:     VecCreateGhost(PETSC_COMM_WORLD,nlocal,PETSC_DECIDE,nghost,ifrom,&gxs);
 75:   }

 77:   /*
 78:       Test VecDuplicate()
 79:   */
 80:   VecDuplicate(gxs,&gx);
 81:   VecDestroy(gxs);

 83:   /*
 84:      Access the local representation
 85:   */
 86:   VecGhostGetLocalForm(gx,&lx);

 88:   /*
 89:      Set the values from 0 to 12 into the "global" vector 
 90:   */
 91:   VecGetOwnershipRange(gx,&rstart,&rend);
 92:   for (i=rstart; i<rend; i++) {
 93:     value = (PetscScalar) i;
 94:     ierr  = VecSetValues(gx,1,&i,&value,INSERT_VALUES);
 95:   }
 96:   VecAssemblyBegin(gx);
 97:   VecAssemblyEnd(gx);

 99:   VecGhostUpdateBegin(gx,INSERT_VALUES,SCATTER_FORWARD);
100:   VecGhostUpdateEnd(gx,INSERT_VALUES,SCATTER_FORWARD);

102:   /*
103:      Print out each vector, including the ghost padding region. 
104:   */
105:   VecGetArray(lx,&array);
106:   for (i=0; i<nlocal+nghost; i++) {
107:     PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%d %gn",i,PetscRealPart(array[i]));
108:   }
109:   VecRestoreArray(lx,&array);
110:   PetscSynchronizedFlush(PETSC_COMM_WORLD);

112:   VecGhostRestoreLocalForm(gx,&lx);
113:   VecDestroy(gx);
114:   if (flg) {PetscFree(tarray);}
115:   PetscFinalize();
116:   return 0;
117: }
118: