Actual source code: ex15f.F

  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 ex15f [-help] [-n <n>] [-mu <mu>] [all SLEPc options]
 21: !
 22: !  Description: Singular value decomposition of the Lauchli matrix.
 23: !
 24: !  The command line options are:
 25: !    -n <n>, where <n> = matrix dimension.
 26: !    -mu <mu>, where <mu> = subdiagonal value.
 27: !
 28: ! ----------------------------------------------------------------------
 29: !
 30:       program main
 31:       implicit none

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

 39: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 40: !     Declarations
 41: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 42: !
 43: !  Variables:
 44: !     A     operator matrix
 45: !     svd   singular value solver context

 47:       Mat            A
 48:       SVD            svd
 49:       SVDType        tname
 50:       PetscReal      tol, error, sigma, mu
 51:       PetscInt       n, i, j, Istart, Iend
 52:       PetscInt       nsv, maxit, its, nconv
 53:       PetscMPIInt    rank
 54:       PetscErrorCode ierr
 55:       PetscBool      flg

 57: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58: !     Beginning of program
 59: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 61:       call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
 62:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
 63:       n = 100
 64:       call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)
 65:       mu = PETSC_SQRT_MACHINE_EPSILON
 66:       call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-mu',mu,flg,ierr)

 68:       if (rank .eq. 0) then
 69:         write(*,100) n, mu
 70:       endif
 71:  100  format (/'Lauchli SVD, n =',I3,', mu=',E12.4,' (Fortran)')

 73: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 74: !     Build the Lauchli matrix
 75: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

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

 82:       call MatGetOwnershipRange(A,Istart,Iend,ierr)
 83:       do i=Istart,Iend-1
 84:         if (i .eq. 0) then
 85:           do j=0,n-1
 86:             call MatSetValue(A,i,j,1.d0,INSERT_VALUES,ierr)
 87:           end do
 88:         else
 89:           call MatSetValue(A,i,i-1,mu,INSERT_VALUES,ierr)
 90:         end if
 91:       enddo

 93:       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
 94:       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)

 96: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97: !     Create the singular value solver and display info
 98: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

100: !     ** Create singular value solver context
101:       call SVDCreate(PETSC_COMM_WORLD,svd,ierr)

103: !     ** Set operator
104:       call SVDSetOperator(svd,A,ierr)

106: !     ** Use thick-restart Lanczos as default solver
107:       call SVDSetType(svd,SVDTRLANCZOS,ierr)

109: !     ** Set solver parameters at runtime
110:       call SVDSetFromOptions(svd,ierr)

112: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: !     Solve the singular value system
114: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

116:       call SVDSolve(svd,ierr)
117:       call SVDGetIterationNumber(svd,its,ierr)
118:       if (rank .eq. 0) then
119:         write(*,110) its
120:       endif
121:  110  format (/' Number of iterations of the method:',I4)
122: 
123: !     ** Optional: Get some information from the solver and display it
124:       call SVDGetType(svd,tname,ierr)
125:       if (rank .eq. 0) then
126:         write(*,120) tname
127:       endif
128:  120  format (' Solution method: ',A)
129:       call SVDGetDimensions(svd,nsv,PETSC_NULL_INTEGER,                 &
130:      &                      PETSC_NULL_INTEGER,ierr)
131:       if (rank .eq. 0) then
132:         write(*,130) nsv
133:       endif
134:  130  format (' Number of requested singular values:',I2)
135:       call SVDGetTolerances(svd,tol,maxit,ierr)
136:       if (rank .eq. 0) then
137:         write(*,140) tol, maxit
138:       endif
139:  140  format (' Stopping condition: tol=',1P,E10.4,', maxit=',I4)

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

145: !     ** Get number of converged singular triplets
146:       call SVDGetConverged(svd,nconv,ierr)
147:       if (rank .eq. 0) then
148:         write(*,150) nconv
149:       endif
150:  150  format (' Number of converged approximate singular triplets:',I2/)

152: !     ** Display singular values and relative errors
153:       if (nconv.gt.0) then
154:         if (rank .eq. 0) then
155:           write(*,*) '       sigma          relative error'
156:           write(*,*) ' ----------------- ------------------'
157:         endif
158:         do i=0,nconv-1
159: !         ** Get converged singular triplet: i-th singular value is stored in sigma
160:           call SVDGetSingularTriplet(svd,i,sigma,PETSC_NULL_OBJECT,     &
161:      &         PETSC_NULL_OBJECT,ierr)

163: !         ** Compute the relative error associated to each eigenpair
164:           call SVDComputeRelativeError(svd,i,error,ierr)
165:           if (rank .eq. 0) then
166:             write(*,160) sigma, error
167:           endif
168:  160      format (1P,'   ',E12.4,'       ',E12.4)

170:         enddo
171:         if (rank .eq. 0) then
172:           write(*,*)
173:         endif
174:       endif

176: !     ** Free work space
177:       call SVDDestroy(svd,ierr)
178:       call MatDestroy(A,ierr)

180:       call SlepcFinalize(ierr)
181:       end