Actual source code: ex1f90.F

  1: !
  2: !    "$Id: ex1f90.F,v 1.27 2001/08/07 03:02:34 balay Exp $"
  3: !
  4: !/*T
  5: !   Concepts: vectors^using basic vector routines;
  6: !   Concepts: Fortran90^using basic vector routines;
  7: !   Processors: n
  8: !T*/
  9: !
 10: ! -----------------------------------------------------------------------

 12:       program main
 13:       implicit none

 15: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 16: !                    Include files
 17: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 18: !
 19: !  The following include statements are required for Fortran programs
 20: !  that use PETSc vectors:
 21: !     petsc.h       - base PETSc routines
 22: !     petscvec.h    - vectors
 23: !     petscvec.h90  - to allow access to Fortran90 features of vectors
 24: !
 25: !  Additional include statements may be needed if using additional
 26: !  PETSc routines in a Fortran program, e.g.,
 27: !     petscviewer.h - viewers
 28: !     petscis.h     - index sets
 29: !
 30:  #include include/finclude/petsc.h
 31:  #include include/finclude/petscvec.h
 32: #include "include/finclude/petscvec.h90"
 33: !
 34: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 35: !                   Variable declarations
 36: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 37: !
 38: !  Variables:
 39: !     x, y, w - vectors
 40: !     z       - array of vectors
 41: !
 42:       Vec              x,y,w
 43:       Vec, pointer :: z(:)
 44:       double precision norm,v,v1,v2
 45:       integer          n,ierr,flg,rank
 46:       PetscScalar      one,two,three,dots(3),dot

 48: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49: !                 Beginning of program
 50: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 52:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 53:       one   = 1.0
 54:       two   = 2.0
 55:       three = 3.0
 56:       n     = 20
 57:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)
 58:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 60: !  Create a vector, specifying only its global dimension.
 61: !  When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
 62: !  the vector format (currently parallel
 63: !  or sequential) is determined at runtime.  Also, the parallel
 64: !  partitioning of the vector is determined by PETSc at runtime.
 65: !
 66: !  Routines for creating particular vector types directly are:
 67: !     VecCreateSeq() - uniprocessor vector
 68: !     VecCreateMPI() - distributed vector, where the user can
 69: !                      determine the parallel partitioning

 71:       call VecCreate(PETSC_COMM_WORLD,x,ierr)
 72:       call VecSetSizes(x,PETSC_DECIDE,n,ierr)
 73:       call VecSetFromOptions(x,ierr)

 75: !  Duplicate some work vectors (of the same format and
 76: !  partitioning as the initial vector).

 78:       call VecDuplicate(x,y,ierr)
 79:       call VecDuplicate(x,w,ierr)

 81: !  Duplicate more work vectors (of the same format and
 82: !  partitioning as the initial vector).  Here we duplicate
 83: !  an array of vectors, which is often more convenient than
 84: !  duplicating individual ones.

 86:       call VecDuplicateVecsF90(x,3,z,ierr)

 88: !  Set the vectors to entries to a constant value.

 90:       call VecSet(one,x,ierr)
 91:       call VecSet(two,y,ierr)
 92:       call VecSet(one,z(1),ierr)
 93:       call VecSet(two,z(2),ierr)
 94:       call VecSet(three,z(3),ierr)

 96: !  Demonstrate various basic vector routines.

 98:       call VecDot(x,x,dot,ierr)
 99:       call VecMDot(3,x,z,dots,ierr)

101: !  Note: If using a complex numbers version of PETSc, then
102: !  PETSC_USE_COMPLEX is defined in the makefiles; otherwise,
103: !  (when using real numbers) it is undefined.

105:       if (rank .eq. 0) then
106: #if defined(PETSC_USE_COMPLEX)
107:          write(6,100) int(PetscRealPart(dot))
108:          write(6,110) int(PetscRealPart(dots(1))),                               &
109:      &                int(PetscRealPart(dots(2))),                               &
110:      &                int(PetscRealPart(dots(3)))
111: #else
112:          write(6,100) int(dot)
113:          write(6,110) int(dots(1)),int(dots(2)),int(dots(3))
114: #endif
115:          write(6,120)
116:       endif
117:  100  format ("Vector length ",i6)
118:  110  format ("Vector length ",3(i6))
119:  120  format ("All other values should be near zero")

121:       call VecScale(two,x,ierr)
122:       call VecNorm(x,NORM_2,norm,ierr)
123:       v = norm-2.0*sqrt(dble(n))
124:       if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
125:       if (rank .eq. 0) write(6,130) v
126:  130  format ("VecScale ",1pe8.2)

128:       call VecCopy(x,w,ierr)
129:       call VecNorm(w,NORM_2,norm,ierr)
130:       v = norm-2.0*sqrt(dble(n))
131:       if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
132:       if (rank .eq. 0) write(6,140) v
133:  140  format ("VecCopy ",1pe8.2)

135:       call VecAXPY(three,x,y,ierr)
136:       call VecNorm(y,NORM_2,norm,ierr)
137:       v = norm-8.0*sqrt(dble(n))
138:       if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
139:       if (rank .eq. 0) write(6,150) v
140:  150  format ("VecAXPY ",1pe8.2)

142:       call VecAYPX(two,x,y,ierr)
143:       call VecNorm(y,NORM_2,norm,ierr)
144:       v = norm-18.0*sqrt(dble(n))
145:       if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
146:       if (rank .eq. 0) write(6,160) v
147:  160  format ("VecAYXP ",1pe8.2)

149:       call VecSwap(x,y,ierr)
150:       call VecNorm(y,NORM_2,norm,ierr)
151:       v = norm-2.0*sqrt(dble(n))
152:       if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
153:       if (rank .eq. 0) write(6,170) v
154:  170  format ("VecSwap ",1pe8.2)

156:       call VecNorm(x,NORM_2,norm,ierr)
157:       v = norm-18.0*sqrt(dble(n))
158:       if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
159:       if (rank .eq. 0) write(6,180) v
160:  180  format ("VecSwap ",1pe8.2)

162:       call VecWAXPY(two,x,y,w,ierr)
163:       call VecNorm(w,NORM_2,norm,ierr)
164:       v = norm-38.0*sqrt(dble(n))
165:       if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
166:       if (rank .eq. 0) write(6,190) v
167:  190  format ("VecWAXPY ",1pe8.2)

169:       call VecPointwiseMult(y,x,w,ierr)
170:       call VecNorm(w,NORM_2,norm,ierr)
171:       v = norm-36.0*sqrt(dble(n))
172:       if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
173:       if (rank .eq. 0) write(6,200) v
174:  200  format ("VecPointwiseMult ",1pe8.2)

176:       call VecPointwiseDivide(x,y,w,ierr)
177:       call VecNorm(w,NORM_2,norm,ierr)
178:       v = norm-9.0*sqrt(dble(n))
179:       if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
180:       if (rank .eq. 0) write(6,210) v
181:  210  format ("VecPointwiseDivide ",1pe8.2)

183: 
184:       dots(1) = one
185:       dots(2) = three
186:       dots(3) = two
187:       call VecSet(one,x,ierr)
188:       call VecMAXPY(3,dots,x,z,ierr)
189:       call VecNorm(z(1),NORM_2,norm,ierr)
190:       v = norm-sqrt(dble(n))
191:       if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
192:       call VecNorm(z(2),NORM_2,norm,ierr)
193:       v1 = norm-2.0*sqrt(dble(n))
194:       if (v1 .gt. -1.d-10 .and. v1 .lt. 1.d-10) v1 = 0.0
195:       call VecNorm(z(3),NORM_2,norm,ierr)
196:       v2 = norm-3.0*sqrt(dble(n))
197:       if (v2 .gt. -1.d-10 .and. v2 .lt. 1.d-10) v2 = 0.0
198:       if (rank .eq. 0) write(6,220) v,v1,v2
199:  220  format ("VecMAXPY ",3(1pe8.2))


202: !  Test whether vector has been corrupted (just to demonstrate this
203: !  routine) not needed in most application codes.

205:       call VecValid(x,flg,ierr)
206:       if (flg .ne. PETSC_TRUE) then
207:          if (rank .eq. 0) then
208:             write(6,*) 'Corrupted vector!'
209:          endif
210:          SETERRQ(1,' ',ierr)
211:       endif

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

216:       call VecDestroy(x,ierr)
217:       call VecDestroy(y,ierr)
218:       call VecDestroy(w,ierr)
219:       call VecDestroyVecsF90(z,3,ierr)
220:       call PetscFinalize(ierr)

222:       end
223: