Actual source code: ex2f.F

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: !
  2: !
  3: !  Description: Builds a parallel vector with 1 component on the first
  4: !               processor, 2 on the second, etc.  Then each processor adds
  5: !               one to all elements except the last rank.
  6: !
  7: !/*T
  8: !   Concepts: vectors^assembling
  9: !   Processors: n
 10: !T*/
 11: ! -----------------------------------------------------------------------

 13:       program main
 14:       implicit none

 16: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 17: !                    Include files
 18: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 19: !
 20: !  The following include statements are required for Fortran programs
 21: !  that use PETSc vectors:
 22: !     petscsys.h       - base PETSc routines
 23: !     petscvec.h    - vectors
 24: !  Additional include statements may be needed if using additional
 25: !  PETSc routines in a Fortran program, e.g.,
 26: !     petscviewer.h - viewers
 27: !     petscis.h     - index sets
 28: !
 29: #include <petsc/finclude/petscsys.h>
 30: #include <petsc/finclude/petscvec.h>

 32: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 33: !                 Beginning of program
 34: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 36:       Vec     x
 37:       PetscInt N,i,ione
 38:       PetscErrorCode ierr
 39:       PetscMPIInt rank
 40:       PetscScalar  one

 42:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 43:       one   = 1.0
 44:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 46: !  Create a parallel vector.
 47: !   - In this case, we specify the size of the local portion on
 48: !     each processor, and PETSc computes the global size.  Alternatively,
 49: !     if we pass the global size and use PETSC_DECIDE for the
 50: !     local size PETSc will choose a reasonable partition trying
 51: !     to put nearly an equal number of elements on each processor.

 53:       N = rank + 1
 54:       call VecCreateMPI(PETSC_COMM_WORLD,N,PETSC_DECIDE,x,ierr)
 55:       call VecGetSize(x,N,ierr)
 56:       call VecSet(x,one,ierr)

 58: !  Set the vector elements.
 59: !   - Note that VecSetValues() uses 0-based row and column numbers
 60: !     in Fortran as well as in C.
 61: !   - Always specify global locations of vector entries.
 62: !   - Each processor can contribute any vector entries,
 63: !     regardless of which processor "owns" them; any nonlocal
 64: !     contributions will be transferred to the appropriate processor
 65: !     during the assembly process.
 66: !   - In this example, the flag ADD_VALUES indicates that all
 67: !     contributions will be added together.

 69:       ione = 1
 70:       do 100 i=0,N-rank-1
 71:          call VecSetValues(x,ione,i,one,ADD_VALUES,ierr)
 72:  100  continue

 74: !  Assemble vector, using the 2-step process:
 75: !    VecAssemblyBegin(), VecAssemblyEnd()
 76: !  Computations can be done while messages are in transition
 77: !  by placing code between these two statements.

 79:       call VecAssemblyBegin(x,ierr)
 80:       call VecAssemblyEnd(x,ierr)

 82: !  View the vector; then destroy it.

 84:       call VecView(x,PETSC_VIEWER_STDOUT_WORLD,ierr)
 85:       call VecDestroy(x,ierr)

 87:       call PetscFinalize(ierr)
 88:       end