Actual source code: ex15f.F
slepc-3.6.3 2016-03-29
1: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
2: ! SLEPc - Scalable Library for Eigenvalue Problem Computations
3: ! Copyright (c) 2002-2015, 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 <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/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
56: PetscScalar one
58: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: ! Beginning of program
60: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: call SlepcInitialize(PETSC_NULL_CHARACTER,ierr)
63: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
64: n = 100
65: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)
66: mu = PETSC_SQRT_MACHINE_EPSILON
67: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-mu',mu,flg,ierr)
69: if (rank .eq. 0) then
70: write(*,100) n, mu
71: endif
72: 100 format (/'Lauchli SVD, n =',I3,', mu=',E12.4,' (Fortran)')
74: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75: ! Build the Lauchli matrix
76: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: call MatCreate(PETSC_COMM_WORLD,A,ierr)
79: call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n+1,n,ierr)
80: call MatSetFromOptions(A,ierr)
81: call MatSetUp(A,ierr)
83: call MatGetOwnershipRange(A,Istart,Iend,ierr)
84: one = 1.0
85: do i=Istart,Iend-1
86: if (i .eq. 0) then
87: do j=0,n-1
88: call MatSetValue(A,i,j,one,INSERT_VALUES,ierr)
89: end do
90: else
91: call MatSetValue(A,i,i-1,mu,INSERT_VALUES,ierr)
92: end if
93: enddo
95: call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
96: call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
98: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99: ! Create the singular value solver and display info
100: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: ! ** Create singular value solver context
103: call SVDCreate(PETSC_COMM_WORLD,svd,ierr)
105: ! ** Set operator
106: call SVDSetOperator(svd,A,ierr)
108: ! ** Use thick-restart Lanczos as default solver
109: call SVDSetType(svd,SVDTRLANCZOS,ierr)
111: ! ** Set solver parameters at runtime
112: call SVDSetFromOptions(svd,ierr)
114: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: ! Solve the singular value system
116: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118: call SVDSolve(svd,ierr)
119: call SVDGetIterationNumber(svd,its,ierr)
120: if (rank .eq. 0) then
121: write(*,110) its
122: endif
123: 110 format (/' Number of iterations of the method:',I4)
125: ! ** Optional: Get some information from the solver and display it
126: call SVDGetType(svd,tname,ierr)
127: if (rank .eq. 0) then
128: write(*,120) tname
129: endif
130: 120 format (' Solution method: ',A)
131: call SVDGetDimensions(svd,nsv,PETSC_NULL_INTEGER, &
132: & PETSC_NULL_INTEGER,ierr)
133: if (rank .eq. 0) then
134: write(*,130) nsv
135: endif
136: 130 format (' Number of requested singular values:',I2)
137: call SVDGetTolerances(svd,tol,maxit,ierr)
138: if (rank .eq. 0) then
139: write(*,140) tol, maxit
140: endif
141: 140 format (' Stopping condition: tol=',1P,E10.4,', maxit=',I4)
143: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144: ! Display solution and clean up
145: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: ! ** Get number of converged singular triplets
148: call SVDGetConverged(svd,nconv,ierr)
149: if (rank .eq. 0) then
150: write(*,150) nconv
151: endif
152: 150 format (' Number of converged approximate singular triplets:',I2/)
154: ! ** Display singular values and relative errors
155: if (nconv.gt.0) then
156: if (rank .eq. 0) then
157: write(*,*) ' sigma relative error'
158: write(*,*) ' ----------------- ------------------'
159: endif
160: do i=0,nconv-1
161: ! ** Get converged singular triplet: i-th singular value is stored in sigma
162: call SVDGetSingularTriplet(svd,i,sigma,PETSC_NULL_OBJECT, &
163: & PETSC_NULL_OBJECT,ierr)
165: ! ** Compute the relative error associated to each eigenpair
166: call SVDComputeError(svd,i,SVD_ERROR_RELATIVE,error,ierr)
167: if (rank .eq. 0) then
168: write(*,160) sigma, error
169: endif
170: 160 format (1P,' ',E12.4,' ',E12.4)
172: enddo
173: if (rank .eq. 0) then
174: write(*,*)
175: endif
176: endif
178: ! ** Free work space
179: call SVDDestroy(svd,ierr)
180: call MatDestroy(A,ierr)
182: call SlepcFinalize(ierr)
183: end