Actual source code: ex16.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: static char help[] = "Demonstrates VecStrideScatter() and VecStrideGather() with subvectors that are also strided.\n\n";

  4: /*T
  5:    Concepts: vectors^sub-vectors;
  6:    Processors: n
  7: T*/

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

 16: #include <petscvec.h>

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

 27:   PetscInitialize(&argc,&argv,(char*)0,help);
 28:   PetscOptionsGetInt(NULL,"-n",&n,NULL);

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

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

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

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

 63:   VecView(s,PETSC_VIEWER_STDOUT_WORLD);
 64:   VecView(r,PETSC_VIEWER_STDOUT_WORLD);

 66:   VecStrideScatterAll(vecs,v,ADD_VALUES);

 68:   VecView(v,PETSC_VIEWER_STDOUT_WORLD);

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