Actual source code: ex20f90.F90
1: !
2: !
3: !
4: ! -----------------------------------------------------------------------
6: program main
8: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
10: !
11: !
12: ! This examples uses Fortran 90 MODULES instead of include files
13: !
14: #include <petsc/finclude/petscvec.h>
15: use petscvec
16: implicit none
18: !
19: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: ! Variable declarations
21: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: !
23: ! Variables:
24: ! x, y, w - vectors
25: ! z - array of vectors
26: !
27: type(tVec) x,y,w
28: type(tVec), pointer :: z(:)
30: PetscReal norm,v,v1,v2,tol
31: PetscInt n,ithree
32: PetscErrorCode ierr
33: PetscMPIInt rank
34: PetscBool flg
35: PetscScalar one,two,three
36: PetscScalar dots(3),dot
37: PetscReal nfloat
39: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40: ! Beginning of program
41: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43: PetscCallA(PetscInitialize(ierr))
44: tol = 1.e-10_PETSC_REAL_KIND
45: one = 1.0
46: two = 2.0
47: three = 3.0
48: n = 20
49: ithree = 3
51: PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr))
52: nfloat = n
53: PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
55: ! Create a vector, specifying only its global dimension.
56: ! When using VecCreate(), VecSetSizes() and VecSetFromOptions(),
57: ! the vector format (currently parallel
58: ! or sequential) is determined at runtime. Also, the parallel
59: ! partitioning of the vector is determined by PETSc at runtime.
60: !
61: ! Routines for creating particular vector types directly are:
62: ! VecCreateSeq() - uniprocessor vector
63: ! VecCreateMPI() - distributed vector, where the user can
64: ! determine the parallel partitioning
66: PetscCallA(VecCreate(PETSC_COMM_WORLD,x,ierr))
67: PetscCallA(VecSetSizes(x,PETSC_DECIDE,n,ierr))
68: PetscCallA(VecSetFromOptions(x,ierr))
70: ! Duplicate some work vectors (of the same format and
71: ! partitioning as the initial vector).
73: PetscCallA(VecDuplicate(x,y,ierr))
74: PetscCallA(VecDuplicate(x,w,ierr))
76: ! Duplicate more work vectors (of the same format and
77: ! partitioning as the initial vector). Here we duplicate
78: ! an array of vectors, which is often more convenient than
79: ! duplicating individual ones.
81: PetscCallA(VecDuplicateVecsF90(x,ithree,z,ierr))
83: ! Set the vectors to entries to a constant value.
85: PetscCallA(VecSet(x,one,ierr))
86: PetscCallA(VecSet(y,two,ierr))
87: PetscCallA(VecSet(z(1),one,ierr))
88: PetscCallA(VecSet(z(2),two,ierr))
89: PetscCallA(VecSet(z(3),three,ierr))
91: ! Demonstrate various basic vector routines.
93: PetscCallA(VecDot(x,x,dot,ierr))
94: PetscCallA(VecMDot(x,ithree,z,dots,ierr))
96: ! Note: If using a complex numbers version of PETSc, then
97: ! PETSC_USE_COMPLEX is defined in the makefiles; otherwise,
98: ! (when using real numbers) it is undefined.
100: if (rank .eq. 0) then
101: #if defined(PETSC_USE_COMPLEX)
102: write(6,100) int(PetscRealPart(dot))
103: write(6,110) int(PetscRealPart(dots(1))),int(PetscRealPart(dots(2))),int(PetscRealPart(dots(3)))
104: #else
105: write(6,100) int(dot)
106: write(6,110) int(dots(1)),int(dots(2)),int(dots(3))
107: #endif
108: write(6,120)
109: endif
110: 100 format ('Vector length ',i6)
111: 110 format ('Vector length ',3(i6))
112: 120 format ('All other values should be near zero')
114: PetscCallA(VecScale(x,two,ierr))
115: PetscCallA(VecNorm(x,NORM_2,norm,ierr))
116: v = abs(norm-2.0*sqrt(nfloat))
117: if (v .gt. -tol .and. v .lt. tol) v = 0.0
118: if (rank .eq. 0) write(6,130) v
119: 130 format ('VecScale ',1pe9.2)
121: PetscCallA(VecCopy(x,w,ierr))
122: PetscCallA(VecNorm(w,NORM_2,norm,ierr))
123: v = abs(norm-2.0*sqrt(nfloat))
124: if (v .gt. -tol .and. v .lt. tol) v = 0.0
125: if (rank .eq. 0) write(6,140) v
126: 140 format ('VecCopy ',1pe9.2)
128: PetscCallA(VecAXPY(y,three,x,ierr))
129: PetscCallA(VecNorm(y,NORM_2,norm,ierr))
130: v = abs(norm-8.0*sqrt(nfloat))
131: if (v .gt. -tol .and. v .lt. tol) v = 0.0
132: if (rank .eq. 0) write(6,150) v
133: 150 format ('VecAXPY ',1pe9.2)
135: PetscCallA(VecAYPX(y,two,x,ierr))
136: PetscCallA(VecNorm(y,NORM_2,norm,ierr))
137: v = abs(norm-18.0*sqrt(nfloat))
138: if (v .gt. -tol .and. v .lt. tol) v = 0.0
139: if (rank .eq. 0) write(6,160) v
140: 160 format ('VecAYXP ',1pe9.2)
142: PetscCallA(VecSwap(x,y,ierr))
143: PetscCallA(VecNorm(y,NORM_2,norm,ierr))
144: v = abs(norm-2.0*sqrt(nfloat))
145: if (v .gt. -tol .and. v .lt. tol) v = 0.0
146: if (rank .eq. 0) write(6,170) v
147: 170 format ('VecSwap ',1pe9.2)
149: PetscCallA(VecNorm(x,NORM_2,norm,ierr))
150: v = abs(norm-18.0*sqrt(nfloat))
151: if (v .gt. -tol .and. v .lt. tol) v = 0.0
152: if (rank .eq. 0) write(6,180) v
153: 180 format ('VecSwap ',1pe9.2)
155: PetscCallA(VecWAXPY(w,two,x,y,ierr))
156: PetscCallA(VecNorm(w,NORM_2,norm,ierr))
157: v = abs(norm-38.0*sqrt(nfloat))
158: if (v .gt. -tol .and. v .lt. tol) v = 0.0
159: if (rank .eq. 0) write(6,190) v
160: 190 format ('VecWAXPY ',1pe9.2)
162: PetscCallA(VecPointwiseMult(w,y,x,ierr))
163: PetscCallA(VecNorm(w,NORM_2,norm,ierr))
164: v = abs(norm-36.0*sqrt(nfloat))
165: if (v .gt. -tol .and. v .lt. tol) v = 0.0
166: if (rank .eq. 0) write(6,200) v
167: 200 format ('VecPointwiseMult ',1pe9.2)
169: PetscCallA(VecPointwiseDivide(w,x,y,ierr))
170: PetscCallA(VecNorm(w,NORM_2,norm,ierr))
171: v = abs(norm-9.0*sqrt(nfloat))
172: if (v .gt. -tol .and. v .lt. tol) v = 0.0
173: if (rank .eq. 0) write(6,210) v
174: 210 format ('VecPointwiseDivide ',1pe9.2)
176: dots(1) = one
177: dots(2) = three
178: dots(3) = two
179: PetscCallA(VecSet(x,one,ierr))
180: PetscCallA(VecMAXPY(x,ithree,dots,z,ierr))
181: PetscCallA(VecNorm(z(1),NORM_2,norm,ierr))
182: v = abs(norm-sqrt(nfloat))
183: if (v .gt. -tol .and. v .lt. tol) v = 0.0
184: PetscCallA(VecNorm(z(2),NORM_2,norm,ierr))
185: v1 = abs(norm-2.0*sqrt(nfloat))
186: if (v1 .gt. -tol .and. v1 .lt. tol) v1 = 0.0
187: PetscCallA(VecNorm(z(3),NORM_2,norm,ierr))
188: v2 = abs(norm-3.0*sqrt(nfloat))
189: if (v2 .gt. -tol .and. v2 .lt. tol) v2 = 0.0
190: if (rank .eq. 0) write(6,220) v,v1,v2
191: 220 format ('VecMAXPY ',3(1pe9.2))
193: ! Free work space. All PETSc objects should be destroyed when they
194: ! are no longer needed.
196: PetscCallA(VecDestroy(x,ierr))
197: PetscCallA(VecDestroy(y,ierr))
198: PetscCallA(VecDestroy(w,ierr))
199: PetscCallA(VecDestroyVecsF90(ithree,z,ierr))
200: PetscCallA(PetscFinalize(ierr))
202: end
204: !/*TEST
205: !
206: ! test:
207: !
208: !TEST*/