Actual source code: ex1f90.F
1: !
2: ! "$Id: ex1f90.F,v 1.26 2001/01/17 22:21:32 bsmith Exp $"
3: !
4: !/*T
5: ! Concepts: vectors^using basic vector routines;
6: ! Concepts: Fortran90;
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: Scalar 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() and VecSetFromOptions(), the vector format (currently parallel
62: ! or sequential) is determined at runtime. Also, the parallel
63: ! partitioning of the vector is determined by PETSc at runtime.
64: !
65: ! Routines for creating particular vector types directly are:
66: ! VecCreateSeq() - uniprocessor vector
67: ! VecCreateMPI() - distributed vector, where the user can
68: ! determine the parallel partitioning
70: call VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,n,x,ierr)
71: call VecSetFromOptions(x,ierr)
73: ! Duplicate some work vectors (of the same format and
74: ! partitioning as the initial vector).
76: call VecDuplicate(x,y,ierr)
77: call VecDuplicate(x,w,ierr)
79: ! Duplicate more work vectors (of the same format and
80: ! partitioning as the initial vector). Here we duplicate
81: ! an array of vectors, which is often more convenient than
82: ! duplicating individual ones.
84: call VecDuplicateVecsF90(x,3,z,ierr)
86: ! Set the vectors to entries to a constant value.
88: call VecSet(one,x,ierr)
89: call VecSet(two,y,ierr)
90: call VecSet(one,z(1),ierr)
91: call VecSet(two,z(2),ierr)
92: call VecSet(three,z(3),ierr)
94: ! Demonstrate various basic vector routines.
96: call VecDot(x,x,dot,ierr)
97: call VecMDot(3,x,z,dots,ierr)
99: ! Note: If using a complex numbers version of PETSc, then
100: ! PETSC_USE_COMPLEX is defined in the makefiles; otherwise,
101: ! (when using real numbers) it is undefined.
103: if (rank .eq. 0) then
104: #if defined(PETSC_USE_COMPLEX)
105: write(6,100) int(PetscRealPart(dot))
106: write(6,110) int(PetscRealPart(dots(1))), &
107: & int(PetscRealPart(dots(2))), &
108: & int(PetscRealPart(dots(3)))
109: #else
110: write(6,100) int(dot)
111: write(6,110) int(dots(1)),int(dots(2)),int(dots(3))
112: #endif
113: write(6,120)
114: endif
115: 100 format ("Vector length ",i6)
116: 110 format ("Vector length ",3(i6))
117: 120 format ("All other values should be near zero")
119: call VecScale(two,x,ierr)
120: call VecNorm(x,NORM_2,norm,ierr)
121: v = norm-2.0*sqrt(dble(n))
122: if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
123: if (rank .eq. 0) write(6,130) v
124: 130 format ("VecScale ",1pe8.2)
126: call VecCopy(x,w,ierr)
127: call VecNorm(w,NORM_2,norm,ierr)
128: v = norm-2.0*sqrt(dble(n))
129: if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
130: if (rank .eq. 0) write(6,140) v
131: 140 format ("VecCopy ",1pe8.2)
133: call VecAXPY(three,x,y,ierr)
134: call VecNorm(y,NORM_2,norm,ierr)
135: v = norm-8.0*sqrt(dble(n))
136: if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
137: if (rank .eq. 0) write(6,150) v
138: 150 format ("VecAXPY ",1pe8.2)
140: call VecAYPX(two,x,y,ierr)
141: call VecNorm(y,NORM_2,norm,ierr)
142: v = norm-18.0*sqrt(dble(n))
143: if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
144: if (rank .eq. 0) write(6,160) v
145: 160 format ("VecAYXP ",1pe8.2)
147: call VecSwap(x,y,ierr)
148: call VecNorm(y,NORM_2,norm,ierr)
149: v = norm-2.0*sqrt(dble(n))
150: if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
151: if (rank .eq. 0) write(6,170) v
152: 170 format ("VecSwap ",1pe8.2)
154: call VecNorm(x,NORM_2,norm,ierr)
155: v = norm-18.0*sqrt(dble(n))
156: if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
157: if (rank .eq. 0) write(6,180) v
158: 180 format ("VecSwap ",1pe8.2)
160: call VecWAXPY(two,x,y,w,ierr)
161: call VecNorm(w,NORM_2,norm,ierr)
162: v = norm-38.0*sqrt(dble(n))
163: if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
164: if (rank .eq. 0) write(6,190) v
165: 190 format ("VecWAXPY ",1pe8.2)
167: call VecPointwiseMult(y,x,w,ierr)
168: call VecNorm(w,NORM_2,norm,ierr)
169: v = norm-36.0*sqrt(dble(n))
170: if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
171: if (rank .eq. 0) write(6,200) v
172: 200 format ("VecPointwiseMult ",1pe8.2)
174: call VecPointwiseDivide(x,y,w,ierr)
175: call VecNorm(w,NORM_2,norm,ierr)
176: v = norm-9.0*sqrt(dble(n))
177: if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
178: if (rank .eq. 0) write(6,210) v
179: 210 format ("VecPointwiseDivide ",1pe8.2)
181:
182: dots(1) = one
183: dots(2) = three
184: dots(3) = two
185: call VecSet(one,x,ierr)
186: call VecMAXPY(3,dots,x,z,ierr)
187: call VecNorm(z(1),NORM_2,norm,ierr)
188: v = norm-sqrt(dble(n))
189: if (v .gt. -1.d-10 .and. v .lt. 1.d-10) v = 0.0
190: call VecNorm(z(2),NORM_2,norm,ierr)
191: v1 = norm-2.0*sqrt(dble(n))
192: if (v1 .gt. -1.d-10 .and. v1 .lt. 1.d-10) v1 = 0.0
193: call VecNorm(z(3),NORM_2,norm,ierr)
194: v2 = norm-3.0*sqrt(dble(n))
195: if (v2 .gt. -1.d-10 .and. v2 .lt. 1.d-10) v2 = 0.0
196: if (rank .eq. 0) write(6,220) v,v1,v2
197: 220 format ("VecMAXPY ",3(1pe8.2))
200: ! Test whether vector has been corrupted (just to demonstrate this
201: ! routine) not needed in most application codes.
203: call VecValid(x,flg,ierr)
204: if (flg .ne. PETSC_TRUE) then
205: if (rank .eq. 0) then
206: write(6,*) 'Corrupted vector!'
207: endif
208: SETERRQ(1,' ',ierr)
209: endif
211: ! Free work space. All PETSc objects should be destroyed when they
212: ! are no longer needed.
214: call VecDestroy(x,ierr)
215: call VecDestroy(y,ierr)
216: call VecDestroy(w,ierr)
217: call VecDestroyVecsF90(z,3,ierr)
218: call PetscFinalize(ierr)
220: end
221: