Actual source code: ex16.c

petsc-3.9.4 2018-09-11
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>

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

 25:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 26:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);

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

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

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

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

 61:   VecView(s,PETSC_VIEWER_STDOUT_WORLD);
 62:   VecView(r,PETSC_VIEWER_STDOUT_WORLD);

 64:   VecStrideScatterAll(vecs,v,ADD_VALUES);

 66:   VecView(v,PETSC_VIEWER_STDOUT_WORLD);

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

 79: /*TEST

 81:      test:
 82:        nsize: 2

 84: TEST*/