Actual source code: ex1f90.F90

slepc-3.20.1 2023-11-27
Report Typos and Errors
  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> ./ex1f90 [-help] [-n <n>] [all SLEPc options]
 11: !
 12: !  Description: Simple example that solves an eigensystem with the EPS object.
 13: !  The standard symmetric eigenvalue problem to be solved corresponds to the
 14: !  Laplacian operator in 1 dimension.
 15: !
 16: !  The command line options are:
 17: !    -n <n>, where <n> = number of grid points = matrix size
 18: !
 19: ! ----------------------------------------------------------------------
 20: !
 21:       program main
 22: #include <slepc/finclude/slepceps.h>
 23:       use slepceps
 24:       use,intrinsic :: iso_c_binding
 25:       implicit none

 27: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 28: !     Declarations
 29: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 30: !
 31: !  Variables:
 32: !     A      operator matrix
 33: !     eps    eigenproblem solver context

 35:       Mat            A
 36:       EPS            eps
 37:       EPSType        tname
 38:       PetscInt       n, i, Istart, Iend, one, two, three
 39:       PetscInt       nev
 40:       PetscInt       row(1), col(3)
 41:       PetscMPIInt    rank
 42:       PetscErrorCode ierr
 43:       PetscBool      flg, terse
 44:       PetscScalar    val(3)

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

 50:       one = 1
 51:       two = 2
 52:       three = 3
 53:       PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER,"ex1f90 test"//c_new_line,ierr))
 54:       if (ierr .ne. 0) then
 55:         print*,'SlepcInitialize failed'
 56:         stop
 57:       endif
 58:       PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
 59:       n = 30
 60:       PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr))

 62:       if (rank .eq. 0) then
 63:         write(*,100) n
 64:       endif
 65:  100  format (/'1-D Laplacian Eigenproblem, n =',I4,' (Fortran)')

 67: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 68: !     Compute the operator matrix that defines the eigensystem, Ax=kx
 69: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 71:       PetscCallA(MatCreate(PETSC_COMM_WORLD,A,ierr))
 72:       PetscCallA(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr))
 73:       PetscCallA(MatSetFromOptions(A,ierr))
 74:       PetscCallA(MatSetUp(A,ierr))

 76:       PetscCallA(MatGetOwnershipRange(A,Istart,Iend,ierr))
 77:       if (Istart .eq. 0) then
 78:         row(1) = 0
 79:         col(1) = 0
 80:         col(2) = 1
 81:         val(1) =  2.0
 82:         val(2) = -1.0
 83:         PetscCallA(MatSetValues(A,one,row,two,col,val,INSERT_VALUES,ierr))
 84:         Istart = Istart+1
 85:       endif
 86:       if (Iend .eq. n) then
 87:         row(1) = n-1
 88:         col(1) = n-2
 89:         col(2) = n-1
 90:         val(1) = -1.0
 91:         val(2) =  2.0
 92:         PetscCallA(MatSetValues(A,one,row,two,col,val,INSERT_VALUES,ierr))
 93:         Iend = Iend-1
 94:       endif
 95:       val(1) = -1.0
 96:       val(2) =  2.0
 97:       val(3) = -1.0
 98:       do i=Istart,Iend-1
 99:         row(1) = i
100:         col(1) = i-1
101:         col(2) = i
102:         col(3) = i+1
103:         PetscCallA(MatSetValues(A,one,row,three,col,val,INSERT_VALUES,ierr))
104:       enddo

106:       PetscCallA(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr))
107:       PetscCallA(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr))

109: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: !     Create the eigensolver and display info
111: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

113: !     ** Create eigensolver context
114:       PetscCallA(EPSCreate(PETSC_COMM_WORLD,eps,ierr))

116: !     ** Set operators. In this case, it is a standard eigenvalue problem
117:       PetscCallA(EPSSetOperators(eps,A,PETSC_NULL_MAT,ierr))
118:       PetscCallA(EPSSetProblemType(eps,EPS_HEP,ierr))

120: !     ** Set solver parameters at runtime
121:       PetscCallA(EPSSetFromOptions(eps,ierr))

123: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124: !     Solve the eigensystem
125: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

127:       PetscCallA(EPSSolve(eps,ierr))

129: !     ** Optional: Get some information from the solver and display it
130:       PetscCallA(EPSGetType(eps,tname,ierr))
131:       if (rank .eq. 0) then
132:         write(*,120) tname
133:       endif
134:  120  format (' Solution method: ',A)
135:       PetscCallA(EPSGetDimensions(eps,nev,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr))
136:       if (rank .eq. 0) then
137:         write(*,130) nev
138:       endif
139:  130  format (' Number of requested eigenvalues:',I4)

141: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: !     Display solution and clean up
143: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

145: !     ** show detailed info unless -terse option is given by user
146:       PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-terse',terse,ierr))
147:       if (terse) then
148:         PetscCallA(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_NULL_VIEWER,ierr))
149:       else
150:         PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL,ierr))
151:         PetscCallA(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD,ierr))
152:         PetscCallA(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD,ierr))
153:         PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr))
154:       endif
155:       PetscCallA(EPSDestroy(eps,ierr))
156:       PetscCallA(MatDestroy(A,ierr))

158:       PetscCallA(SlepcFinalize(ierr))
159:       end

161: !/*TEST
162: !
163: !   build:
164: !      requires: defined(PETSC_USING_F2003) defined(PETSC_USING_F90FREEFORM)
165: !
166: !   test:
167: !      args: -eps_nev 4 -terse
168: !      filter: sed -e "s/3.83791/3.83792/"
169: !
170: !TEST*/