Actual source code: test7f.F

slepc-3.7.3 2016-09-29
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: !  Program usage: mpiexec -n <np> ./test7f [-help] [-n <n>] [all SLEPc options] 
 21: !
 22: !  Description: Simple example that solves an eigensystem with the EPS object.
 23: !  Same problem as ex1f but with simplified output.
 24: !
 25: !  The command line options are:
 26: !    -n <n>, where <n> = number of grid points = matrix size
 27: !
 28: ! ---------------------------------------------------------------------- 
 29: !
 30:       program main
 31:       implicit none

 33: #include <petsc/finclude/petscsys.h>
 34: #include <petsc/finclude/petscvec.h>
 35: #include <petsc/finclude/petscmat.h>
 36: #include <slepc/finclude/slepcsys.h>
 37: #include <slepc/finclude/slepceps.h>

 39: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 40: !     Declarations
 41: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 42: !
 43: !  Variables:
 44: !     A     operator matrix
 45: !     eps   eigenproblem solver context

 47:       Mat            A
 48:       EPS            eps
 49:       EPSType        tname
 50:       PetscInt       n, i, Istart, Iend
 51:       PetscInt       nev
 52:       PetscInt       col(3)
 53:       PetscInt       i1,i2,i3
 54:       PetscMPIInt    rank
 55:       PetscErrorCode ierr
 56:       PetscBool      flg
 57:       PetscScalar    value(3)

 59: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60: !     Beginning of program
 61: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

 63:       call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
 64:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 65:       n = 30
 66:       call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,   &
 67:      &                        '-n',n,flg,ierr)

 69:       if (rank .eq. 0) then
 70:         write(*,100) n
 71:       endif
 72:  100  format (/'1-D Laplacian Eigenproblem, n =',I3,' (Fortran)')

 74: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 75: !     Compute the operator matrix that defines the eigensystem, Ax=kx
 76: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

 78:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
 79:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
 80:       call MatSetFromOptions(A,ierr)
 81:       call MatSetUp(A,ierr)

 83:       i1 = 1
 84:       i2 = 2
 85:       i3 = 3
 86:       call MatGetOwnershipRange(A,Istart,Iend,ierr)
 87:       if (Istart .eq. 0) then 
 88:         i = 0
 89:         col(1) = 0
 90:         col(2) = 1
 91:         value(1) =  2.0
 92:         value(2) = -1.0
 93:         call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
 94:         Istart = Istart+1
 95:       endif
 96:       if (Iend .eq. n) then 
 97:         i = n-1
 98:         col(1) = n-2
 99:         col(2) = n-1
100:         value(1) = -1.0
101:         value(2) =  2.0
102:         call MatSetValues(A,i1,i,i2,col,value,INSERT_VALUES,ierr)
103:         Iend = Iend-1
104:       endif
105:       value(1) = -1.0
106:       value(2) =  2.0
107:       value(3) = -1.0
108:       do i=Istart,Iend-1
109:         col(1) = i-1
110:         col(2) = i
111:         col(3) = i+1
112:         call MatSetValues(A,i1,i,i3,col,value,INSERT_VALUES,ierr)
113:       enddo

115:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
116:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

118: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
119: !     Create the eigensolver and display info
120: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

122: !     ** Create eigensolver context
123:       call EPSCreate(PETSC_COMM_WORLD,eps,ierr)

125: !     ** Set operators. In this case, it is a standard eigenvalue problem
126:       call EPSSetOperators(eps,A,PETSC_NULL_OBJECT,ierr)
127:       call EPSSetProblemType(eps,EPS_HEP,ierr)

129: !     ** Set solver parameters at runtime
130:       call EPSSetFromOptions(eps,ierr)

132: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
133: !     Solve the eigensystem
134: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

136:       call EPSSolve(eps,ierr) 

138: !     ** Optional: Get some information from the solver and display it
139:       call EPSGetType(eps,tname,ierr)
140:       if (rank .eq. 0) then
141:         write(*,120) tname
142:       endif
143:  120  format (' Solution method: ',A)
144:       call EPSGetDimensions(eps,nev,PETSC_NULL_INTEGER,                 &
145:      &                      PETSC_NULL_INTEGER,ierr)
146:       if (rank .eq. 0) then
147:         write(*,130) nev
148:       endif
149:  130  format (' Number of requested eigenvalues:',I2)

151: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
152: !     Display solution and clean up
153: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

155:       call EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_NULL_OBJECT,ierr)
156:       call EPSDestroy(eps,ierr)
157:       call MatDestroy(A,ierr)

159:       call SlepcFinalize(ierr)
160:       end