Actual source code: test2f.F

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: !
  2: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3: !  SLEPc - Scalable Library for Eigenvalue Problem Computations
  4: !  Copyright (c) 2002-2021, 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 NEP Fortran interface.
 11: !
 12: ! ----------------------------------------------------------------------
 13: !
 14:       program main
 15: #include <slepc/finclude/slepcnep.h>
 16:       use slepcnep
 17:       implicit none

 19: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: !     Declarations
 21: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22:       Mat                A(3),B
 23:       FN                 f(3),g
 24:       NEP                nep
 25:       DS                 ds
 26:       RG                 rg
 27:       PetscReal          tol
 28:       PetscScalar        coeffs(2),tget,val
 29:       PetscInt           n,i,its,Istart,Iend
 30:       PetscInt           nev,ncv,mpd,nterm
 31:       PetscInt           nc,np
 32:       NEPWhich           which
 33:       NEPConvergedReason reason
 34:       NEPType            tname
 35:       NEPRefine          refine
 36:       NEPRefineScheme    rscheme
 37:       NEPConv            conv
 38:       NEPStop            stp
 39:       NEPProblemType     ptype
 40:       MatStructure       mstr
 41:       PetscMPIInt        rank
 42:       PetscErrorCode     ierr
 43:       PetscViewerAndFormat vf

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

 49:       call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
 50:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 51:       n = 20
 52:       if (rank .eq. 0) then
 53:         write(*,100) n
 54:       endif
 55:  100  format (/'Diagonal Nonlinear Eigenproblem, n =',I3,' (Fortran)')

 57: !     Matrices
 58:       call MatCreate(PETSC_COMM_WORLD,A(1),ierr)
 59:       call MatSetSizes(A(1),PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
 60:       call MatSetFromOptions(A(1),ierr)
 61:       call MatSetUp(A(1),ierr)
 62:       call MatGetOwnershipRange(A(1),Istart,Iend,ierr)
 63:       do i=Istart,Iend-1
 64:         val = i+1
 65:         call MatSetValue(A(1),i,i,val,INSERT_VALUES,ierr)
 66:       enddo
 67:       call MatAssemblyBegin(A(1),MAT_FINAL_ASSEMBLY,ierr)
 68:       call MatAssemblyEnd(A(1),MAT_FINAL_ASSEMBLY,ierr)

 70:       call MatCreate(PETSC_COMM_WORLD,A(2),ierr)
 71:       call MatSetSizes(A(2),PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
 72:       call MatSetFromOptions(A(2),ierr)
 73:       call MatSetUp(A(2),ierr)
 74:       call MatGetOwnershipRange(A(2),Istart,Iend,ierr)
 75:       do i=Istart,Iend-1
 76:         val = 1
 77:         call MatSetValue(A(2),i,i,val,INSERT_VALUES,ierr)
 78:       enddo
 79:       call MatAssemblyBegin(A(2),MAT_FINAL_ASSEMBLY,ierr)
 80:       call MatAssemblyEnd(A(2),MAT_FINAL_ASSEMBLY,ierr)

 82:       call MatCreate(PETSC_COMM_WORLD,A(3),ierr)
 83:       call MatSetSizes(A(3),PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
 84:       call MatSetFromOptions(A(3),ierr)
 85:       call MatSetUp(A(3),ierr)
 86:       call MatGetOwnershipRange(A(3),Istart,Iend,ierr)
 87:       do i=Istart,Iend-1
 88:         val = real(n)/real(i+1)
 89:         call MatSetValue(A(3),i,i,val,INSERT_VALUES,ierr)
 90:       enddo
 91:       call MatAssemblyBegin(A(3),MAT_FINAL_ASSEMBLY,ierr)
 92:       call MatAssemblyEnd(A(3),MAT_FINAL_ASSEMBLY,ierr)

 94: !     Functions: f0=-lambda, f1=1.0, f2=sqrt(lambda)
 95:       call FNCreate(PETSC_COMM_WORLD,f(1),ierr)
 96:       call FNSetType(f(1),FNRATIONAL,ierr)
 97:       nc = 2
 98:       coeffs(1) = -1.0
 99:       coeffs(2) = 0.0
100:       call FNRationalSetNumerator(f(1),nc,coeffs,ierr)

102:       call FNCreate(PETSC_COMM_WORLD,f(2),ierr)
103:       call FNSetType(f(2),FNRATIONAL,ierr)
104:       nc = 1
105:       coeffs(1) = 1.0
106:       call FNRationalSetNumerator(f(2),nc,coeffs,ierr)

108:       call FNCreate(PETSC_COMM_WORLD,f(3),ierr)
109:       call FNSetType(f(3),FNSQRT,ierr)

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

115:       call NEPCreate(PETSC_COMM_WORLD,nep,ierr)
116:       nterm = 3
117:       mstr = SAME_NONZERO_PATTERN
118:       call NEPSetSplitOperator(nep,nterm,A,f,mstr,ierr)
119:       call NEPGetSplitOperatorInfo(nep,nterm,mstr,ierr)
120:       if (rank .eq. 0) then
121:         write(*,110) nterm
122:       endif
123:  110  format (' Nonlinear function with ',I2,' terms')
124:       i = 0
125:       call NEPGetSplitOperatorTerm(nep,i,B,g,ierr)
126:       call MatView(B,PETSC_NULL_VIEWER,ierr)
127:       call FNView(g,PETSC_NULL_VIEWER,ierr)

129:       call NEPSetType(nep,NEPRII,ierr)
130:       call NEPGetType(nep,tname,ierr)
131:       if (rank .eq. 0) then
132:         write(*,120) tname
133:       endif
134:  120  format (' Type set to ',A)

136:       call NEPGetProblemType(nep,ptype,ierr)
137:       if (rank .eq. 0) then
138:         write(*,130) ptype
139:       endif
140:  130  format (' Problem type before changing = ',I2)
141:       call NEPSetProblemType(nep,NEP_RATIONAL,ierr)
142:       call NEPGetProblemType(nep,ptype,ierr)
143:       if (rank .eq. 0) then
144:         write(*,140) ptype
145:       endif
146:  140  format (' ... changed to ',I2)

148:       np = 1
149:       tol = 1e-9
150:       its = 2
151:       refine = NEP_REFINE_SIMPLE
152:       rscheme = NEP_REFINE_SCHEME_EXPLICIT
153:       call NEPSetRefine(nep,refine,np,tol,its,rscheme,ierr)
154:       call NEPGetRefine(nep,refine,np,tol,its,rscheme,ierr)
155:       if (rank .eq. 0) then
156:         write(*,190) refine,tol,its,rscheme
157:       endif
158:  190  format (' Refinement: ',I2,', tol=',F12.9,', its=',I2,', scheme=',&
159:      &   I2)

161:       tget = 1.1
162:       call NEPSetTarget(nep,tget,ierr)
163:       call NEPGetTarget(nep,tget,ierr)
164:       call NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE,ierr)
165:       call NEPGetWhichEigenpairs(nep,which,ierr)
166:       if (rank .eq. 0) then
167:         write(*,200) which,PetscRealPart(tget)
168:       endif
169:  200  format (' Which = ',I2,', target = ',F4.1)

171:       nev = 1
172:       ncv = 12
173:       call NEPSetDimensions(nep,nev,ncv,PETSC_DEFAULT_INTEGER,ierr)
174:       call NEPGetDimensions(nep,nev,ncv,mpd,ierr)
175:       if (rank .eq. 0) then
176:         write(*,210) nev,ncv,mpd
177:       endif
178:  210  format (' Dimensions: nev=',I2,', ncv=',I2,', mpd=',I2)

180:       tol = 1.0e-6
181:       its = 200
182:       call NEPSetTolerances(nep,tol,its,ierr)
183:       call NEPGetTolerances(nep,tol,its,ierr)
184:       if (rank .eq. 0) then
185:         write(*,220) tol,its
186:       endif
187:  220  format (' Tolerance =',F9.6,', max_its =',I4)

189:       call NEPSetConvergenceTest(nep,NEP_CONV_ABS,ierr)
190:       call NEPGetConvergenceTest(nep,conv,ierr)
191:       call NEPSetStoppingTest(nep,NEP_STOP_BASIC,ierr)
192:       call NEPGetStoppingTest(nep,stp,ierr)
193:       if (rank .eq. 0) then
194:         write(*,230) conv,stp
195:       endif
196:  230  format (' Convergence test =',I2,', stopping test =',I2)

198:       call PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,        &
199:      &                   PETSC_VIEWER_DEFAULT,vf,ierr)
200:       call NEPMonitorSet(nep,NEPMONITORFIRST,vf,                        &
201:      &                   PetscViewerAndFormatDestroy,ierr)
202:       call NEPMonitorConvergedCreate(PETSC_VIEWER_STDOUT_WORLD,         &
203:      &                   PETSC_VIEWER_DEFAULT,PETSC_NULL_VEC,vf,ierr)
204:       call NEPMonitorSet(nep,NEPMONITORCONVERGED,vf,                    &
205:      &                   NEPMonitorConvergedDestroy,ierr)
206:       call NEPMonitorCancel(nep,ierr)

208:       call NEPGetDS(nep,ds,ierr)
209:       call DSView(ds,PETSC_NULL_VIEWER,ierr)
210:       call NEPSetFromOptions(nep,ierr)

212:       call NEPGetRG(nep,rg,ierr)
213:       call RGView(rg,PETSC_NULL_VIEWER,ierr)

215:       call NEPSolve(nep,ierr)
216:       call NEPGetConvergedReason(nep,reason,ierr)
217:       if (rank .eq. 0) then
218:         write(*,240) reason
219:       endif
220:  240  format (' Finished - converged reason =',I2)

222: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223: !     Display solution and clean up
224: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225:       call NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_NULL_VIEWER,ierr)
226:       call NEPDestroy(nep,ierr)
227:       call MatDestroy(A(1),ierr)
228:       call MatDestroy(A(2),ierr)
229:       call MatDestroy(A(3),ierr)
230:       call FNDestroy(f(1),ierr)
231:       call FNDestroy(f(2),ierr)
232:       call FNDestroy(f(3),ierr)

234:       call SlepcFinalize(ierr)
235:       end

237: !/*TEST
238: !
239: !   test:
240: !      suffix: 1
241: !      requires: !single
242: !
243: !TEST*/