Actual source code: ex16.c

petsc-3.3-p7 2013-05-11
  2: /* Program usage:  mpiexec ex1 [-help] [all PETSc options] */

  4: static char help[] = "Demonstrates VecStrideScatter() and VecStrideGather() with subvectors that are also strided.\n\n";

  6: /*T
  7:    Concepts: vectors^sub-vectors;
  8:    Processors: n
  9: T*/

 11: /* 
 12:   Include "petscvec.h" so that we can use vectors.  Note that this file
 13:   automatically includes:
 14:      petscsys.h       - base PETSc routines   petscis.h     - index sets
 15:      petscviewer.h - viewers
 16: */

 18: #include <petscvec.h>

 22: int main(int argc,char **argv)
 23: {
 24:   Vec            v,s,r,vecs[2];               /* vectors */
 25:   PetscInt       i,start,end,n = 20;
 27:   PetscScalar    value;

 29:   PetscInitialize(&argc,&argv,(char*)0,help);
 30:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

 32:   /* 
 33:       Create multi-component vector with 2 components
 34:   */
 35:   VecCreate(PETSC_COMM_WORLD,&v);
 36:   VecSetSizes(v,PETSC_DECIDE,n);
 37:   VecSetBlockSize(v,4);
 38:   VecSetFromOptions(v);

 40:   /* 
 41:       Create double-component vectors
 42:   */
 43:   VecCreate(PETSC_COMM_WORLD,&s);
 44:   VecSetSizes(s,PETSC_DECIDE,n/2);
 45:   VecSetBlockSize(s,2);
 46:   VecSetFromOptions(s);
 47:   VecDuplicate(s,&r);

 49:   vecs[0] = s;
 50:   vecs[1] = r;
 51:   /*
 52:      Set the vector values
 53:   */
 54:   VecGetOwnershipRange(v,&start,&end);
 55:   for (i=start; i<end; i++) {
 56:     value = i;
 57:     VecSetValues(v,1,&i,&value,INSERT_VALUES);
 58:   }

 60:   /*
 61:      Get the components from the multi-component vector to the other vectors
 62:   */
 63:   VecStrideGatherAll(v,vecs,INSERT_VALUES);

 65:   VecView(s,PETSC_VIEWER_STDOUT_WORLD);
 66:   VecView(r,PETSC_VIEWER_STDOUT_WORLD);

 68:   VecStrideScatterAll(vecs,v,ADD_VALUES);

 70:   VecView(v,PETSC_VIEWER_STDOUT_WORLD);

 72:   /* 
 73:      Free work space.  All PETSc objects should be destroyed when they
 74:      are no longer needed.
 75:   */
 76:   VecDestroy(&v);
 77:   VecDestroy(&s);
 78:   VecDestroy(&r);
 79:   PetscFinalize();
 80:   return 0;
 81: }
 82: