Actual source code: test15f.F90
slepc-3.22.2 2024-12-02
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: ! Program usage: mpiexec -n <np> ./test15f [-help] [-n <n>] [all SLEPc options]
11: !
12: ! Description: Tests custom monitors from Fortran.
13: !
14: ! The command line options are:
15: ! -n <n>, where <n> = number of grid points = matrix size
16: ! -my_eps_monitor, activates the custom monitor
17: !
18: ! ----------------------------------------------------------------------
19: !
20: program main
21: #include <slepc/finclude/slepceps.h>
22: use slepceps
23: implicit none
25: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26: ! Declarations
27: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
28: !
29: ! Variables:
30: ! A operator matrix
31: ! eps eigenproblem solver context
33: Mat A
34: EPS eps
35: EPSType tname
36: PetscInt n, i, Istart, Iend
37: PetscInt nev
38: PetscInt col(3)
39: PetscInt i1,i2,i3
40: PetscMPIInt rank
41: PetscErrorCode ierr
42: PetscBool flg
43: PetscScalar value(3)
45: ! Note: Any user-defined Fortran routines (such as MyEPSMonitor)
46: ! MUST be declared as external.
48: external MyEPSMonitor
50: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: ! Beginning of program
52: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER,ierr))
55: PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
56: n = 30
57: PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr))
59: if (rank .eq. 0) then
60: write(*,100) n
61: endif
62: 100 format (/'1-D Laplacian Eigenproblem, n =',I3,' (Fortran)')
64: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: ! Compute the operator matrix that defines the eigensystem, Ax=kx
66: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: PetscCallA(MatCreate(PETSC_COMM_WORLD,A,ierr))
69: PetscCallA(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr))
70: PetscCallA(MatSetFromOptions(A,ierr))
72: i1 = 1
73: i2 = 2
74: i3 = 3
75: PetscCallA(MatGetOwnershipRange(A,Istart,Iend,ierr))
76: if (Istart .eq. 0) then
77: i = 0
78: col(1) = 0
79: col(2) = 1
80: value(1) = 2.0
81: value(2) = -1.0
82: PetscCallA(MatSetValues(A,i1,[i],i2,col,value,INSERT_VALUES,ierr))
83: Istart = Istart+1
84: endif
85: if (Iend .eq. n) then
86: i = n-1
87: col(1) = n-2
88: col(2) = n-1
89: value(1) = -1.0
90: value(2) = 2.0
91: PetscCallA(MatSetValues(A,i1,[i],i2,col,value,INSERT_VALUES,ierr))
92: Iend = Iend-1
93: endif
94: value(1) = -1.0
95: value(2) = 2.0
96: value(3) = -1.0
97: do i=Istart,Iend-1
98: col(1) = i-1
99: col(2) = i
100: col(3) = i+1
101: PetscCallA(MatSetValues(A,i1,[i],i3,col,value,INSERT_VALUES,ierr))
102: enddo
104: PetscCallA(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr))
105: PetscCallA(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr))
107: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108: ! Create the eigensolver and display info
109: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: ! ** Create eigensolver context
112: PetscCallA(EPSCreate(PETSC_COMM_WORLD,eps,ierr))
114: ! ** Set operators. In this case, it is a standard eigenvalue problem
115: PetscCallA(EPSSetOperators(eps,A,PETSC_NULL_MAT,ierr))
116: PetscCallA(EPSSetProblemType(eps,EPS_HEP,ierr))
118: ! ** Set user-defined monitor
119: PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my_eps_monitor',flg,ierr))
120: if (flg) then
121: PetscCallA(EPSMonitorSet(eps,MyEPSMonitor,0,PETSC_NULL_FUNCTION,ierr))
122: endif
124: ! ** Set solver parameters at runtime
125: PetscCallA(EPSSetFromOptions(eps,ierr))
127: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: ! Solve the eigensystem
129: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131: PetscCallA(EPSSolve(eps,ierr))
133: ! ** Optional: Get some information from the solver and display it
134: PetscCallA(EPSGetType(eps,tname,ierr))
135: if (rank .eq. 0) then
136: write(*,120) tname
137: endif
138: 120 format (' Solution method: ',A)
139: PetscCallA(EPSGetDimensions(eps,nev,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr))
140: if (rank .eq. 0) then
141: write(*,130) nev
142: endif
143: 130 format (' Number of requested eigenvalues:',I2)
145: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: ! Display solution and clean up
147: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149: PetscCallA(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_NULL_VIEWER,ierr))
150: PetscCallA(EPSDestroy(eps,ierr))
151: PetscCallA(MatDestroy(A,ierr))
153: PetscCallA(SlepcFinalize(ierr))
154: end
156: ! --------------------------------------------------------------
157: !
158: ! MyEPSMonitor - This is a user-defined routine for monitoring
159: ! the EPS iterative solvers.
160: !
161: ! Input Parameters:
162: ! eps - eigensolver context
163: ! its - iteration number
164: ! nconv - number of converged eigenpairs
165: ! eigr - real part of the eigenvalues
166: ! eigi - imaginary part of the eigenvalues
167: ! errest- relative error estimates for each eigenpair
168: ! nest - number of error estimates
169: ! dummy - optional user-defined monitor context (unused here)
170: !
171: subroutine MyEPSMonitor(eps,its,nconv,eigr,eigi,errest,nest,dummy,ierr)
172: #include <slepc/finclude/slepceps.h>
173: use slepceps
174: implicit none
176: EPS eps
177: PetscErrorCode ierr
178: PetscInt its,nconv,nest,dummy
179: PetscScalar eigr(*),eigi(*)
180: PetscReal re,errest(*)
181: PetscMPIInt rank
183: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
184: if (its .gt. 0 .and. rank .eq. 0) then
185: re = PetscRealPart(eigr(nconv+1))
186: write(6,140) its,nconv,re,errest(nconv+1)
187: endif
189: 140 format(i3,' EPS nconv=',i2,' first unconverged value (error) ',f7.4,' (',g10.3,')')
190: ierr = 0
191: end
193: !/*TEST
194: !
195: ! test:
196: ! suffix: 1
197: ! args: -my_eps_monitor
198: ! requires: double
199: !
200: !TEST*/