Actual source code: ex1f90.F90

slepc-3.7.4 2017-05-17
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> ./ex1f90 [-help] [-n <n>] [all SLEPc options] 
 21: !
 22: !  Description: Simple example that solves an eigensystem with the EPS object.
 23: !  The standard symmetric eigenvalue problem to be solved corresponds to the 
 24: !  Laplacian operator in 1 dimension. 
 25: !
 26: !  The command line options are:
 27: !    -n <n>, where <n> = number of grid points = matrix size
 28: !
 29: ! ---------------------------------------------------------------------- 
 30: !
 31:       program main

 33: #include <slepc/finclude/slepcepsdef.h>
 34:       use slepceps

 36:       implicit none

 38: ! For usage without modules, uncomment the following lines and remove 
 39: ! the previous lines between 'program main' and 'implicit none'
 40: !
 41: !#include <petsc/finclude/petsc.h>
 42: !#include <slepc/finclude/slepc.h>

 44: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 45: !     Declarations
 46: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 47: !
 48: !  Variables:
 49: !     A      operator matrix
 50: !     eps    eigenproblem solver context

 52: #if defined(PETSC_USE_FORTRAN_DATATYPES)
 53:       type(Mat)      A
 54:       type(EPS)      eps
 55: #else
 56:       Mat            A
 57:       EPS            eps
 58: #endif
 59:       EPSType        tname
 60:       PetscInt       n, i, Istart, Iend, one, two, three
 61:       PetscInt       nev
 62:       PetscInt       row(1), col(3)
 63:       PetscMPIInt    rank
 64:       PetscErrorCode ierr
 65:       PetscBool      flg, terse
 66:       PetscScalar    value(3)

 68: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 69: !     Beginning of program
 70: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

 72:       one = 1
 73:       two = 2
 74:       three = 3
 75:       call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
 76:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 77:       n = 30
 78:       call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,   &
 79:      &                        '-n',n,flg,ierr)

 81:       if (rank .eq. 0) then
 82:         write(*,100) n
 83:       endif
 84:  100  format (/'1-D Laplacian Eigenproblem, n =',I4,' (Fortran)')

 86: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 87: !     Compute the operator matrix that defines the eigensystem, Ax=kx
 88: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

 90:       call MatCreate(PETSC_COMM_WORLD,A,ierr)
 91:       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr)
 92:       call MatSetFromOptions(A,ierr)
 93:       call MatSetUp(A,ierr)

 95:       call MatGetOwnershipRange(A,Istart,Iend,ierr)
 96:       if (Istart .eq. 0) then 
 97:         row(1) = 0
 98:         col(1) = 0
 99:         col(2) = 1
100:         value(1) =  2.0
101:         value(2) = -1.0
102:         call MatSetValues(A,one,row,two,col,value,INSERT_VALUES,ierr)
103:         Istart = Istart+1
104:       endif
105:       if (Iend .eq. n) then 
106:         row(1) = n-1
107:         col(1) = n-2
108:         col(2) = n-1
109:         value(1) = -1.0
110:         value(2) =  2.0
111:         call MatSetValues(A,one,row,two,col,value,INSERT_VALUES,ierr)
112:         Iend = Iend-1
113:       endif
114:       value(1) = -1.0
115:       value(2) =  2.0
116:       value(3) = -1.0
117:       do i=Istart,Iend-1
118:         row(1) = i
119:         col(1) = i-1
120:         col(2) = i
121:         col(3) = i+1
122:         call MatSetValues(A,one,row,three,col,value,INSERT_VALUES,ierr)
123:       enddo

125:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
126:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

128: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
129: !     Create the eigensolver and display info
130: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

132: !     ** Create eigensolver context
133:       call EPSCreate(PETSC_COMM_WORLD,eps,ierr)

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

139: !     ** Set solver parameters at runtime
140:       call EPSSetFromOptions(eps,ierr)

142: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
143: !     Solve the eigensystem
144: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

146:       call EPSSolve(eps,ierr) 

148: !     ** Optional: Get some information from the solver and display it
149:       call EPSGetType(eps,tname,ierr)
150:       if (rank .eq. 0) then
151:         write(*,120) tname
152:       endif
153:  120  format (' Solution method: ',A)
154:       call EPSGetDimensions(eps,nev,PETSC_NULL_INTEGER,                 &
155:      &                      PETSC_NULL_INTEGER,ierr)
156:       if (rank .eq. 0) then
157:         write(*,130) nev
158:       endif
159:  130  format (' Number of requested eigenvalues:',I4)

161: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
162: !     Display solution and clean up
163: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

165: !     ** show detailed info unless -terse option is given by user
166:       call PetscOptionsHasName(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,  &
167:      &                        '-terse',terse,ierr)
168:       if (terse) then
169:         call EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_NULL_OBJECT,ierr)
170:       else
171:         call PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,           &
172:      &                   PETSC_VIEWER_ASCII_INFO_DETAIL,ierr)
173:         call EPSReasonView(eps,PETSC_VIEWER_STDOUT_WORLD,ierr)
174:         call EPSErrorView(eps,EPS_ERROR_RELATIVE,                       &
175:      &                   PETSC_VIEWER_STDOUT_WORLD,ierr)
176:         call PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr)
177:       endif
178:       call EPSDestroy(eps,ierr)
179:       call MatDestroy(A,ierr)

181:       call SlepcFinalize(ierr)
182:       end