Actual source code: test3f.F

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  2: !  SLEPc - Scalable Library for Eigenvalue Problem Computations
  3: !  Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
  4: !
  5: !  This file is part of SLEPc.
  6: !     
  7: !  SLEPc is free software: you can redistribute it and/or modify it under  the
  8: !  terms of version 3 of the GNU Lesser General Public License as published by
  9: !  the Free Software Foundation.
 10: !
 11: !  SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 12: !  WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 13: !  FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 14: !  more details.
 15: !
 16: !  You  should have received a copy of the GNU Lesser General  Public  License
 17: !  along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 18: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 19: !
 20: !  Description: Simple example to test the PEP Fortran interface.
 21: !
 22: ! ---------------------------------------------------------------------- 
 23: !
 24:       program main
 25:       implicit none

 27: #include <petsc/finclude/petscsys.h>
 28: #include <petsc/finclude/petscvec.h>
 29: #include <petsc/finclude/petscmat.h>
 30: #include <petsc/finclude/petscviewer.h>
 31: #include <slepc/finclude/slepcsys.h>
 32: #include <slepc/finclude/slepcpep.h>

 34: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 35: !     Declarations
 36: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 37:       Mat                A(3),B
 38:       PEP                pep
 39:       ST                 st
 40:       KSP                ksp
 41:       DS                 ds
 42:       PetscReal          tol,tolabs,alpha,lambda
 43:       PetscScalar        tget,val
 44:       PetscInt           n,i,its,Istart,Iend
 45:       PetscInt           nev,ncv,mpd,nmat,np
 46:       PEPWhich           which
 47:       PEPConvergedReason reason
 48:       PEPType            tname
 49:       PEPExtract         extr
 50:       PEPBasis           basis
 51:       PEPScale           scal
 52:       PEPRefine          refine
 53:       PEPRefineScheme    rscheme
 54:       PEPConv            conv
 55:       PEPStop            stp
 56:       PEPProblemType     ptype
 57:       PetscMPIInt        rank
 58:       PetscErrorCode     ierr
 59:       SlepcConvMonitor   ctx
 60:       PetscViewerAndFormat vf

 62: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63: !     Beginning of program
 64: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

 66:       call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
 67:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 68:       n = 20
 69:       if (rank .eq. 0) then
 70:         write(*,100) n
 71:       endif
 72:  100  format (/'Diagonal Quadratic Eigenproblem, n =',I3,' (Fortran)')

 74:       call MatCreate(PETSC_COMM_WORLD,A(1),ierr)
 75:       call MatSetSizes(A(1),PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
 76:       call MatSetFromOptions(A(1),ierr)
 77:       call MatSetUp(A(1),ierr)
 78:       call MatGetOwnershipRange(A(1),Istart,Iend,ierr)
 79:       do i=Istart,Iend-1
 80:         val = i+1.
 81:         call MatSetValue(A(1),i,i,val,INSERT_VALUES,ierr)
 82:       enddo
 83:       call MatAssemblyBegin(A(1),MAT_FINAL_ASSEMBLY,ierr)
 84:       call MatAssemblyEnd(A(1),MAT_FINAL_ASSEMBLY,ierr)

 86:       call MatCreate(PETSC_COMM_WORLD,A(2),ierr)
 87:       call MatSetSizes(A(2),PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
 88:       call MatSetFromOptions(A(2),ierr)
 89:       call MatSetUp(A(2),ierr)
 90:       call MatGetOwnershipRange(A(2),Istart,Iend,ierr)
 91:       do i=Istart,Iend-1
 92:         val = 1
 93:         call MatSetValue(A(2),i,i,val,INSERT_VALUES,ierr)
 94:       enddo
 95:       call MatAssemblyBegin(A(2),MAT_FINAL_ASSEMBLY,ierr)
 96:       call MatAssemblyEnd(A(2),MAT_FINAL_ASSEMBLY,ierr)

 98:       call MatCreate(PETSC_COMM_WORLD,A(3),ierr)
 99:       call MatSetSizes(A(3),PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
100:       call MatSetFromOptions(A(3),ierr)
101:       call MatSetUp(A(3),ierr)
102:       call MatGetOwnershipRange(A(3),Istart,Iend,ierr)
103:       do i=Istart,Iend-1
104:         val = n/(i+1.)
105:         call MatSetValue(A(3),i,i,val,INSERT_VALUES,ierr)
106:       enddo
107:       call MatAssemblyBegin(A(3),MAT_FINAL_ASSEMBLY,ierr)
108:       call MatAssemblyEnd(A(3),MAT_FINAL_ASSEMBLY,ierr)

110: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
111: !     Create eigensolver and test interface functions
112: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

114:       call PEPCreate(PETSC_COMM_WORLD,pep,ierr)
115:       nmat = 3
116:       call PEPSetOperators(pep,nmat,A,ierr)
117:       call PEPGetNumMatrices(pep,nmat,ierr)
118:       if (rank .eq. 0) then
119:         write(*,110) nmat-1
120:       endif
121:  110  format (' Polynomial of degree ',I2)
122:       i = 0
123:       call PEPGetOperators(pep,i,B,ierr)
124:       call MatView(B,PETSC_NULL_OBJECT,ierr)

126:       call PEPSetType(pep,PEPTOAR,ierr)
127:       call PEPGetType(pep,tname,ierr)
128:       if (rank .eq. 0) then
129:         write(*,120) tname
130:       endif
131:  120  format (' Type set to ',A)

133:       call PEPGetProblemType(pep,ptype,ierr)
134:       if (rank .eq. 0) then
135:         write(*,130) ptype
136:       endif
137:  130  format (' Problem type before changing = ',I2)
138:       call PEPSetProblemType(pep,PEP_HERMITIAN,ierr)
139:       call PEPGetProblemType(pep,ptype,ierr)
140:       if (rank .eq. 0) then
141:         write(*,140) ptype
142:       endif
143:  140  format (' ... changed to ',I2)

145:       call PEPGetExtract(pep,extr,ierr)
146:       if (rank .eq. 0) then
147:         write(*,150) extr
148:       endif
149:  150  format (' Extraction before changing = ',I2)
150:       call PEPSetExtract(pep,PEP_EXTRACT_STRUCTURED,ierr)
151:       call PEPGetExtract(pep,extr,ierr)
152:       if (rank .eq. 0) then
153:         write(*,160) extr
154:       endif
155:  160  format (' ... changed to ',I2)

157:       alpha = .1
158:       its = 5
159:       lambda = 1.
160:       scal = PEP_SCALE_SCALAR
161:       call PEPSetScale(pep,scal,alpha,PETSC_NULL_OBJECT,                &
162:      &                 PETSC_NULL_OBJECT,its,lambda,ierr)
163:       call PEPGetScale(pep,scal,alpha,PETSC_NULL_OBJECT,                &
164:      &                 PETSC_NULL_OBJECT,its,lambda,ierr)
165:       if (rank .eq. 0) then
166:         write(*,170) scal,alpha,its
167:       endif
168:  170  format (' Scaling: ',I2,', alpha=',F6.4,', its=',I2)

170:       basis = PEP_BASIS_CHEBYSHEV1
171:       call PEPSetBasis(pep,basis,ierr)
172:       call PEPGetBasis(pep,basis,ierr)
173:       if (rank .eq. 0) then
174:         write(*,180) basis
175:       endif
176:  180  format (' Polynomial basis: ',I2)
177:  
178:       np = 1
179:       tol = 1e-9
180:       its = 2
181:       refine = PEP_REFINE_SIMPLE
182:       rscheme = PEP_REFINE_SCHEME_SCHUR
183:       call PEPSetRefine(pep,refine,np,tol,its,rscheme,ierr)
184:       call PEPGetRefine(pep,refine,np,tol,its,rscheme,ierr)
185:       if (rank .eq. 0) then
186:         write(*,190) refine,tol,its,rscheme
187:       endif
188:  190  format (' Refinement: ',I2,', tol=',F6.4,', its=',I2', schem=',I2)

190:       tget = 4.8
191:       call PEPSetTarget(pep,tget,ierr)
192:       call PEPGetTarget(pep,tget,ierr)
193:       call PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE,ierr)
194:       call PEPGetWhichEigenpairs(pep,which,ierr)
195:       if (rank .eq. 0) then
196:         write(*,200) which,PetscRealPart(tget)
197:       endif
198:  200  format (' Which = ',I2,', target = ',F3.1)

200:       nev = 4
201:       call PEPSetDimensions(pep,nev,PETSC_DEFAULT_INTEGER,              &
202:      &                      PETSC_DEFAULT_INTEGER,ierr)
203:       call PEPGetDimensions(pep,nev,ncv,mpd,ierr)
204:       if (rank .eq. 0) then
205:         write(*,210) nev,ncv,mpd
206:       endif
207:  210  format (' Dimensions: nev=',I2,', ncv=',I2,', mpd=',I2)

209:       tol = 2.2e-4
210:       its = 200
211:       call PEPSetTolerances(pep,tol,its,ierr)
212:       call PEPGetTolerances(pep,tol,its,ierr)
213:       if (rank .eq. 0) then
214:         write(*,220) tol,its
215:       endif
216:  220  format (' Tolerance =',F7.5,', max_its =',I4)

218:       call PEPSetConvergenceTest(pep,PEP_CONV_ABS,ierr)
219:       call PEPGetConvergenceTest(pep,conv,ierr)
220:       call PEPSetStoppingTest(pep,PEP_STOP_BASIC,ierr)
221:       call PEPGetStoppingTest(pep,stp,ierr)
222:       if (rank .eq. 0) then
223:         write(*,230) conv,stp
224:       endif
225:  230  format (' Convergence test =',I2,', stopping test =',I2)

227:       call PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,        &
228:      &                   PETSC_VIEWER_DEFAULT,vf,ierr)
229:       call PEPMonitorSet(pep,PEPMONITORFIRST,vf,                        &
230:      &                   PetscViewerAndFormatDestroy,ierr)
231:       call SlepcConvMonitorCreate(PETSC_VIEWER_STDOUT_WORLD,            &
232:      &                   PETSC_VIEWER_DEFAULT,ctx,ierr)
233:       call PEPMonitorSet(pep,PEPMONITORCONVERGED,ctx,                   &
234:      &                   SlepcConvMonitorDestroy,ierr)
235:       call PEPMonitorCancel(pep,ierr)

237:       call PEPGetST(pep,st,ierr) 
238:       call STGetKSP(st,ksp,ierr) 
239:       tol = 1.e-8
240:       tolabs = 1.e-35
241:       call KSPSetTolerances(ksp,tol,tolabs,PETSC_DEFAULT_REAL,          &
242:      &                      PETSC_DEFAULT_INTEGER,ierr)
243:       call STView(st,PETSC_NULL_OBJECT,ierr) 
244:       call PEPGetDS(pep,ds,ierr) 
245:       call DSView(ds,PETSC_NULL_OBJECT,ierr) 

247:       call PEPSetFromOptions(pep,ierr)
248:       call PEPSolve(pep,ierr) 
249:       call PEPGetConvergedReason(pep,reason,ierr)
250:       call PEPGetIterationNumber(pep,its,ierr)
251:       if (rank .eq. 0) then
252:         write(*,240) reason,its
253:       endif
254:  240  format (' Finished - converged reason =',I2,', its=',I4)

256: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
257: !     Display solution and clean up
258: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
259:       call PEPErrorView(pep,PEP_ERROR_RELATIVE,PETSC_NULL_OBJECT,ierr)
260:       call PEPDestroy(pep,ierr)
261:       call MatDestroy(A(1),ierr)
262:       call MatDestroy(A(2),ierr)
263:       call MatDestroy(A(3),ierr)

265:       call SlepcFinalize(ierr)
266:       end