Actual source code: ex16f90.F90

slepc-3.7.1 2016-05-27
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> ./ex16f90 [-help] [-n <n>] [-m <m>] [SLEPc opts] 
 21: !
 22: !  Description: Simple example that solves a quadratic eigensystem with the
 23: !  PEP 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 <slepc/finclude/slepcpepdef.h>
 34:       use slepcpep

 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: !     M,C,K  problem matrices
 50: !     pep    polynomial eigenproblem solver context

 52: #if defined(PETSC_USE_FORTRAN_DATATYPES)
 53:       type(Mat)      M, C, K, A(3)
 54:       type(PEP)      pep
 55: #else
 56:       Mat            M, C, K, A(3)
 57:       PEP            pep
 58: #endif
 59:       PEPType        tname
 60:       PetscInt       N, nx, ny, i, j, Istart, Iend, II
 61:       PetscInt       nev, ithree
 62:       PetscMPIInt    rank
 63:       PetscErrorCode ierr
 64:       PetscBool      flg, terse
 65:       PetscScalar    mone, two, four, val

 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_OBJECT,PETSC_NULL_CHARACTER,   &
 75:      &                        '-n',nx,flg,ierr)
 76:       call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,   &
 77:      &                        '-m',ny,flg,ierr)
 78:       if (.not. flg) then
 79:         ny = nx
 80:       endif
 81:       N = nx*ny
 82:       if (rank .eq. 0) then
 83:         write(*,100) N, nx, ny
 84:       endif
 85:  100  format (/'Quadratic Eigenproblem, N=',I6,' (',I4,'x',I4,' grid)')

 87: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 88: !     Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 89: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

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

119: !     ** C is the 1-D Laplacian on horizontal lines
120:       call MatCreate(PETSC_COMM_WORLD,C,ierr)
121:       call MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N,ierr)
122:       call MatSetFromOptions(C,ierr)
123:       call MatSetUp(C,ierr)
124:       call MatGetOwnershipRange(C,Istart,Iend,ierr)
125:       two = 2.0
126:       do II=Istart,Iend-1
127:         i = II/nx
128:         j = II-i*nx
129:         if (j .gt. 0) then 
130:           call MatSetValue(C,II,II-1,mone,INSERT_VALUES,ierr)
131:         endif
132:         if (j .lt. nx-1) then 
133:           call MatSetValue(C,II,II+1,mone,INSERT_VALUES,ierr)
134:         endif
135:         call MatSetValue(C,II,II,two,INSERT_VALUES,ierr)
136:       end do
137:       call MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY,ierr)
138:       call MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY,ierr)

140: !     ** M is a diagonal matrix
141:       call MatCreate(PETSC_COMM_WORLD,M,ierr)
142:       call MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N,ierr)
143:       call MatSetFromOptions(M,ierr)
144:       call MatSetUp(M,ierr)
145:       call MatGetOwnershipRange(M,Istart,Iend,ierr)
146:       do II=Istart,Iend-1
147:         val = II+1
148:         call MatSetValue(M,II,II,val,INSERT_VALUES,ierr)
149:       end do
150:       call MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY,ierr)
151:       call MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY,ierr)

153: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
154: !     Create the eigensolver and set various options
155: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

157: !     ** Create eigensolver context
158:       call PEPCreate(PETSC_COMM_WORLD,pep,ierr)

160: !     ** Set matrices and problem type
161:       A(1) = K
162:       A(2) = C
163:       A(3) = M
164:       ithree = 3
165:       call PEPSetOperators(pep,ithree,A,ierr)
166:       call PEPSetProblemType(pep,PEP_GENERAL,ierr)

168: !     ** Set solver parameters at runtime
169:       call PEPSetFromOptions(pep,ierr)

171: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
172: !     Solve the eigensystem
173: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

175:       call PEPSolve(pep,ierr) 

177: !     ** Optional: Get some information from the solver and display it
178:       call PEPGetType(pep,tname,ierr)
179:       if (rank .eq. 0) then
180:         write(*,120) tname
181:       endif
182:  120  format (' Solution method: ',A)
183:       call PEPGetDimensions(pep,nev,PETSC_NULL_INTEGER,                 &
184:      &                      PETSC_NULL_INTEGER,ierr)
185:       if (rank .eq. 0) then
186:         write(*,130) nev
187:       endif
188:  130  format (' Number of requested eigenvalues:',I4)

190: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
191: !     Display solution and clean up
192: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

194: !     ** show detailed info unless -terse option is given by user
195:       call PetscOptionsHasName(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER,  &
196:      &                        '-terse',terse,ierr)
197:       if (terse) then
198:         call PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_NULL_OBJECT,ierr)
199:       else
200:         call PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,           &
201:      &                   PETSC_VIEWER_ASCII_INFO_DETAIL,ierr)
202:         call PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD,ierr)
203:         call PEPErrorView(pep,PEP_ERROR_BACKWARD,                       &
204:      &                   PETSC_VIEWER_STDOUT_WORLD,ierr)
205:         call PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr)
206:       endif
207:       call PEPDestroy(pep,ierr)
208:       call MatDestroy(K,ierr)
209:       call MatDestroy(C,ierr)
210:       call MatDestroy(M,ierr)
211:       call SlepcFinalize(ierr)
212:       end