Actual source code: test14f.F90

  1: !
  2: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3: !  SLEPc - Scalable Library for Eigenvalue Problem Computations
  4: !  Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
  5: !
  6: !  This file is part of SLEPc.
  7: !  SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: !
 10: !  Description: Simple example to test the EPS Fortran interface.
 11: !
 12: ! ----------------------------------------------------------------------
 13: !
 14: #include <slepc/finclude/slepceps.h>
 15: program test14f
 16:   use slepceps
 17:   implicit none

 19: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: ! Declarations
 21: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22:   Mat                  :: A, B
 23:   EPS                  :: eps
 24:   ST                   :: st
 25:   KSP                  :: ksp
 26:   DS                   :: ds
 27:   PetscReal            :: cut, tol, tolabs
 28:   PetscScalar          :: tget, val
 29:   PetscInt             :: n, i, its, Istart, Iend
 30:   PetscInt             :: nev, ncv, mpd
 31:   PetscBool            :: flg
 32:   EPSConvergedReason   :: reason
 33:   EPSType              :: tname
 34:   EPSExtraction        :: extr
 35:   EPSBalance           :: bal
 36:   EPSWhich             :: which
 37:   EPSConv              :: conv
 38:   EPSStop              :: stp
 39:   EPSProblemType       :: ptype
 40:   PetscMPIInt          :: rank
 41:   PetscErrorCode       :: ierr
 42:   PetscViewerAndFormat :: vf

 44: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 45: ! Beginning of program
 46: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 48:   PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER, ierr))
 49:   PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
 50:   n = 20
 51:   if (rank == 0) then
 52:     write (*, '(/a,i3,a)') 'Diagonal Eigenproblem, n =', n, ' (Fortran)'
 53:   end if

 55:   PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr))
 56:   PetscCallA(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n, ierr))
 57:   PetscCallA(MatSetFromOptions(A, ierr))
 58:   PetscCallA(MatGetOwnershipRange(A, Istart, Iend, ierr))
 59:   do i = Istart, Iend - 1
 60:     val = i + 1
 61:     PetscCallA(MatSetValue(A, i, i, val, INSERT_VALUES, ierr))
 62:   end do
 63:   PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr))
 64:   PetscCallA(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY, ierr))

 66: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67: ! Create eigensolver and test interface functions
 68: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 70:   PetscCallA(EPSCreate(PETSC_COMM_WORLD, eps, ierr))
 71:   PetscCallA(EPSSetOperators(eps, A, PETSC_NULL_MAT, ierr))
 72:   PetscCallA(EPSGetOperators(eps, B, PETSC_NULL_MAT, ierr))
 73:   PetscCallA(MatView(B, PETSC_NULL_VIEWER, ierr))

 75:   PetscCallA(EPSSetType(eps, EPSKRYLOVSCHUR, ierr))
 76:   PetscCallA(EPSGetType(eps, tname, ierr))
 77:   if (rank == 0) then
 78:     write (*, '(a,a)') ' Type set to ', tname
 79:   end if

 81:   PetscCallA(EPSGetProblemType(eps, ptype, ierr))
 82:   if (rank == 0) then
 83:     write (*, '(a,i2)') ' Problem type before changing = ', ptype
 84:   end if
 85:   PetscCallA(EPSSetProblemType(eps, EPS_HEP, ierr))
 86:   PetscCallA(EPSGetProblemType(eps, ptype, ierr))
 87:   if (rank == 0) then
 88:     write (*, '(a,i2)') ' ... changed to ', ptype
 89:   end if
 90:   PetscCallA(EPSIsGeneralized(eps, flg, ierr))
 91:   if (flg .and. rank == 0) then
 92:     write (*, *) 'generalized'
 93:   end if
 94:   PetscCallA(EPSIsHermitian(eps, flg, ierr))
 95:   if (flg .and. rank == 0) then
 96:     write (*, *) 'hermitian'
 97:   end if
 98:   PetscCallA(EPSIsPositive(eps, flg, ierr))
 99:   if (flg .and. rank == 0) then
100:     write (*, *) 'positive'
101:   end if

103:   PetscCallA(EPSGetExtraction(eps, extr, ierr))
104:   if (rank == 0) then
105:     write (*, '(a,i2)') ' Extraction before changing = ', extr
106:   end if
107:   PetscCallA(EPSSetExtraction(eps, EPS_HARMONIC, ierr))
108:   PetscCallA(EPSGetExtraction(eps, extr, ierr))
109:   if (rank == 0) then
110:     write (*, '(a,i2)') ' ... changed to ', extr
111:   end if

113:   its = 8
114:   cut = 2.0e-6
115:   bal = EPS_BALANCE_ONESIDE
116:   PetscCallA(EPSSetBalance(eps, bal, its, cut, ierr))
117:   PetscCallA(EPSGetBalance(eps, bal, its, cut, ierr))
118:   if (rank == 0) then
119:     write (*, '(a,i2,a,i2,a,f9.6)') ' Balance: ', bal, ', its=', its, ', cutoff=', cut
120:   end if

122:   tget = 4.8
123:   PetscCallA(EPSSetTarget(eps, tget, ierr))
124:   PetscCallA(EPSGetTarget(eps, tget, ierr))
125:   PetscCallA(EPSSetWhichEigenpairs(eps, EPS_TARGET_MAGNITUDE, ierr))
126:   PetscCallA(EPSGetWhichEigenpairs(eps, which, ierr))
127:   if (rank == 0) then
128:     write (*, '(a,i2,a,f4.1)') ' Which = ', which, ', target = ', PetscRealPart(tget)
129:   end if

131:   nev = 4
132:   PetscCallA(EPSSetDimensions(eps, nev, PETSC_DETERMINE_INTEGER, PETSC_DETERMINE_INTEGER, ierr))
133:   PetscCallA(EPSGetDimensions(eps, nev, ncv, mpd, ierr))
134:   if (rank == 0) then
135:     write (*, '(a,i2,a,i2,a,i2)') ' Dimensions: nev=', nev, ', ncv=', ncv, ', mpd=', mpd
136:   end if

138:   tol = 2.2e-4
139:   its = 200
140:   PetscCallA(EPSSetTolerances(eps, tol, its, ierr))
141:   PetscCallA(EPSGetTolerances(eps, tol, its, ierr))
142:   if (rank == 0) then
143:     write (*, '(a,f8.5,a,i4)') ' Tolerance =', tol, ', max_its =', its
144:   end if

146:   PetscCallA(EPSSetConvergenceTest(eps, EPS_CONV_ABS, ierr))
147:   PetscCallA(EPSGetConvergenceTest(eps, conv, ierr))
148:   PetscCallA(EPSSetStoppingTest(eps, EPS_STOP_BASIC, ierr))
149:   PetscCallA(EPSGetStoppingTest(eps, stp, ierr))
150:   if (rank == 0) then
151:     write (*, '(a,i2,a,i2)') ' Convergence test =', conv, ', stopping test =', stp
152:   end if

154:   PetscCallA(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, vf, ierr))
155:   PetscCallA(EPSMonitorSet(eps, EPSMONITORFIRST, vf, PetscViewerAndFormatDestroy, ierr))
156:   PetscCallA(EPSMonitorConvergedCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, PETSC_NULL, vf, ierr))
157:   PetscCallA(EPSMonitorSet(eps, EPSMONITORCONVERGED, vf, PetscViewerAndFormatDestroy, ierr))
158:   PetscCallA(EPSMonitorCancel(eps, ierr))

160:   PetscCallA(EPSGetST(eps, st, ierr))
161:   PetscCallA(STGetKSP(st, ksp, ierr))
162:   tol = 1.e-8
163:   tolabs = 1.e-35
164:   PetscCallA(KSPSetTolerances(ksp, tol, tolabs, PETSC_CURRENT_REAL, PETSC_CURRENT_INTEGER, ierr))
165:   PetscCallA(STView(st, PETSC_NULL_VIEWER, ierr))
166:   PetscCallA(EPSGetDS(eps, ds, ierr))
167:   PetscCallA(DSView(ds, PETSC_NULL_VIEWER, ierr))

169:   PetscCallA(EPSSetFromOptions(eps, ierr))
170:   PetscCallA(EPSSolve(eps, ierr))
171:   PetscCallA(EPSGetConvergedReason(eps, reason, ierr))
172:   PetscCallA(EPSGetIterationNumber(eps, its, ierr))
173:   if (rank == 0) then
174:     write (*, '(a,i2,a,i4)') ' Finished - converged reason =', reason, ', its = ', its
175:   end if

177: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178: ! Display solution and clean up
179: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180:   PetscCallA(EPSErrorView(eps, EPS_ERROR_RELATIVE, PETSC_NULL_VIEWER, ierr))
181:   PetscCallA(EPSDestroy(eps, ierr))
182:   PetscCallA(MatDestroy(A, ierr))

184:   PetscCallA(SlepcFinalize(ierr))
185: end program test14f

187: !/*TEST
188: !
189: !   test:
190: !      suffix: 1
191: !      args: -eps_ncv 14
192: !      filter: sed -e "s/00001/00000/" | sed -e "s/4.99999/5.00000/" | sed -e "s/5.99999/6.00000/"
193: !
194: !TEST*/