Actual source code: ex3f.F

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: !
  2: !
  3: !  Description: Displays a vector visually.
  4: !
  5: !/*T
  6: !   Concepts: vectors^drawing vectors;
  7: !   Processors: n
  8: !T*/
  9: ! -----------------------------------------------------------------------

 11:       program main
 12:       implicit none

 14: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 15: !                    Include files
 16: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 17: !
 18: !  The following include statements are required for Fortran programs
 19: !  that use PETSc vectors:
 20: !     petscsys.h       - base PETSc routines
 21: !     petscvec.h    - vectors
 22: !  Include petscviewer.h so that we can use the PETSc viewers.
 23: !
 24: #include <petsc/finclude/petscsys.h>
 25: #include <petsc/finclude/petscvec.h>
 26: #include <petsc/finclude/petscviewer.h>

 28: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 29: !                 Beginning of program
 30: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 32:       Vec     x
 33:       PetscViewer  viewer
 34:       PetscScalar  v
 35:       PetscInt i,istart,iend,n,ione
 36:       PetscErrorCode ierr
 37:       PetscBool  flg

 39:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 40:       n = 50
 41:       call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,      &
 42:      &                        '-n',n,flg,ierr)

 44: !  Create a vector, specifying only its global dimension.
 45: !  When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
 46: !  the vector format (currently parallel
 47: !  or sequential) is determined at runtime.  Also, the parallel
 48: !  partitioning of the vector is determined by PETSc at runtime.
 49:       call VecCreate(PETSC_COMM_WORLD,x,ierr)
 50:       call VecSetSizes(x,PETSC_DECIDE,n,ierr)
 51:       call VecSetFromOptions(x,ierr)

 53: !  Currently, all PETSc parallel vectors are partitioned by
 54: !  contiguous chunks of rows across the processors.  Determine
 55: !  which vector are locally owned.
 56:       call VecGetOwnershipRange(x,istart,iend,ierr)

 58: !  Set the vector elements.
 59: !   - Always specify global locations of vector entries.
 60: !   - Each processor needs to insert only elements that it owns locally.
 61:       ione = 1
 62:       do 100 i=istart,iend-1
 63:          v = 1.0*real(i)
 64:          call VecSetValues(x,ione,i,v,INSERT_VALUES,ierr)
 65:  100  continue

 67: !  Assemble vector, using the 2-step process:
 68: !    VecAssemblyBegin(), VecAssemblyEnd()
 69: !  Computations can be done while messages are in transition
 70: !  by placing code between these two statements.
 71:       call VecAssemblyBegin(x,ierr)
 72:       call VecAssemblyEnd(x,ierr)

 74: !  Open an X-window viewer.  Note that we specify the same communicator
 75: !  for the viewer as we used for the distributed vector (PETSC_COMM_WORLD).
 76: !    - Helpful runtime option:
 77: !         -draw_pause <pause> : sets time (in seconds) that the
 78: !               program pauses after PetscDrawPause() has been called
 79: !              (0 is default, -1 implies until user input).

 81:       call PetscViewerDrawOpen(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,        &
 82:      &                   PETSC_NULL_CHARACTER,0,0,300,300,viewer,ierr)

 84: !  View the vector
 85:       call VecView(x,viewer,ierr)

 87: !  Free work space.  All PETSc objects should be destroyed when they
 88: !  are no longer needed.

 90:       call PetscViewerDestroy(viewer,ierr)
 91:       call VecDestroy(x,ierr)

 93:       call PetscFinalize(ierr)
 94:       end