Actual source code: ex16f90.F90

  1: !  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  2: !  SLEPc - Scalable Library for Eigenvalue Problem Computations
  3: !  Copyright (c) 2002-2012, 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: mpirun -np n ex16f90 [-help] [-n <n>] [-m <m>] [SLEPc opts]
 21: !
 22: !  Description: Simple example that solves a quadratic eigensystem with the
 23: !  QEP object. This is the Fortran90 equivalent to ex16.c
 24: !
 25: !  The command line options are:
 26: !    -n <n>, where <n> = number of grid subdivisions in x dimension
 27: !    -m <m>, where <m> = number of grid subdivisions in y dimension
 28: !
 29: ! ----------------------------------------------------------------------
 30: !
 31:       program main

 33: #include <finclude/slepcqepdef.h>
 34:       use slepcqep

 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 <finclude/petsc.h>
 42: !#include <finclude/slepc.h>

 44: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 45: !     Declarations
 46: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 47: !
 48: !  Variables:
 49: !     M,C,K  problem matrices
 50: !     solver quadratic eigenproblem solver context

 52: #if defined(PETSC_USE_FORTRAN_DATATYPES)
 53:       type(Mat)      M, C, K
 54:       type(QEP)      solver
 55: #else
 56:       Mat            M, C, K
 57:       QEP            solver
 58: #endif
 59:       QEPType        tname
 60:       PetscReal      tol
 61:       PetscInt       N, nx, ny, i, j, Istart, Iend, II
 62:       PetscInt       nev, maxit
 63:       PetscMPIInt    rank
 64:       PetscErrorCode ierr
 65:       PetscBool      flg

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

 71:       call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
 72:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 73:       nx = 10
 74:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',nx,flg,ierr)
 75:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-m',ny,flg,ierr)
 76:       if (.not. flg) then
 77:         ny = nx
 78:       endif
 79:       N = nx*ny
 80:       if (rank .eq. 0) then
 81:         write(*,100) N, nx, ny
 82:       endif
 83:  100  format (/'Quadratic Eigenproblem, N=',I6,' (',I4,'x',I4,' grid)')

 85: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 86: !     Compute the matrices that define the eigensystem, (k^2*K+k*X+M)x=0
 87: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 89: !     ** K is the 2-D Laplacian
 90:       call MatCreate(PETSC_COMM_WORLD,K,ierr)
 91:       call MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N,ierr)
 92:       call MatSetFromOptions(K,ierr)
 93:       call MatSetUp(K,ierr)
 94:       call MatGetOwnershipRange(K,Istart,Iend,ierr)
 95:       do II=Istart,Iend-1
 96:         i = II/nx
 97:         j = II-i*nx
 98:         if (i .gt. 0) then
 99:           call MatSetValue(K,II,II-nx,-1.D0,INSERT_VALUES,ierr)
100:         endif
101:         if (i .lt. ny-1) then
102:           call MatSetValue(K,II,II+nx,-1.D0,INSERT_VALUES,ierr)
103:         endif
104:         if (j .gt. 0) then
105:           call MatSetValue(K,II,II-1,-1.D0,INSERT_VALUES,ierr)
106:         endif
107:         if (j .lt. nx-1) then
108:           call MatSetValue(K,II,II+1,-1.D0,INSERT_VALUES,ierr)
109:         endif
110:         call MatSetValue(K,II,II,4.D0,INSERT_VALUES,ierr)
111:       end do
112:       call MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY,ierr)
113:       call MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY,ierr)

115: !     ** C is the zero matrix
116:       call MatCreate(PETSC_COMM_WORLD,C,ierr)
117:       call MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N,ierr)
118:       call MatSetFromOptions(C,ierr)
119:       call MatSetUp(C,ierr)
120:       call MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY,ierr)
121:       call MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY,ierr)

123: !     ** M is the identity matrix
124:       call MatCreate(PETSC_COMM_WORLD,M,ierr)
125:       call MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N,ierr)
126:       call MatSetFromOptions(M,ierr)
127:       call MatSetUp(M,ierr)
128:       call MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY,ierr)
129:       call MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY,ierr)
130:       call MatShift(M,1.D0,ierr)

132: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133: !     Create the eigensolver and set various options
134: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

136: !     ** Create eigensolver context
137:       call QEPCreate(PETSC_COMM_WORLD,solver,ierr)

139: !     ** Set matrices and problem type
140:       call QEPSetOperators(solver,M,C,K,ierr)
141:       call QEPSetProblemType(solver,QEP_GENERAL,ierr)

143: !     ** Set solver parameters at runtime
144:       call QEPSetFromOptions(solver,ierr)

146: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: !     Solve the eigensystem
148: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

150:       call QEPSolve(solver,ierr)
151: 
152: !     ** Optional: Get some information from the solver and display it
153:       call QEPGetType(solver,tname,ierr)
154:       if (rank .eq. 0) then
155:         write(*,120) tname
156:       endif
157:  120  format (' Solution method: ',A)
158:       call QEPGetDimensions(solver,nev,PETSC_NULL_INTEGER,              &
159:      &                      PETSC_NULL_INTEGER,ierr)
160:       if (rank .eq. 0) then
161:         write(*,130) nev
162:       endif
163:  130  format (' Number of requested eigenvalues:',I4)
164:       call QEPGetTolerances(solver,tol,maxit,ierr)
165:       if (rank .eq. 0) then
166:         write(*,140) tol, maxit
167:       endif
168:  140  format (' Stopping condition: tol=',1P,E10.4,', maxit=',I4)

170: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171: !     Display solution and clean up
172: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

174:       call QEPPrintSolution(solver,PETSC_NULL_OBJECT,ierr)
175:       call QEPDestroy(solver,ierr)
176:       call MatDestroy(K,ierr)
177:       call MatDestroy(C,ierr)
178:       call MatDestroy(M,ierr)
179:       call SlepcFinalize(ierr)
180:       end