Actual source code: ex15f.F
slepc-3.7.3 2016-09-29
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> ./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_OBJECT,PETSC_NULL_CHARACTER, &
66: & '-n',n,flg,ierr)
67: mu = PETSC_SQRT_MACHINE_EPSILON
68: call PetscOptionsGetReal(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER, &
69: & '-mu',mu,flg,ierr)
71: if (rank .eq. 0) then
72: write(*,100) n, mu
73: endif
74: 100 format (/'Lauchli SVD, n =',I3,', mu=',E12.4,' (Fortran)')
76: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77: ! Build the Lauchli matrix
78: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: call MatCreate(PETSC_COMM_WORLD,A,ierr)
81: call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n+1,n,ierr)
82: call MatSetFromOptions(A,ierr)
83: call MatSetUp(A,ierr)
85: call MatGetOwnershipRange(A,Istart,Iend,ierr)
86: one = 1.0
87: do i=Istart,Iend-1
88: if (i .eq. 0) then
89: do j=0,n-1
90: call MatSetValue(A,i,j,one,INSERT_VALUES,ierr)
91: end do
92: else
93: call MatSetValue(A,i,i-1,mu,INSERT_VALUES,ierr)
94: end if
95: enddo
97: call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
98: call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
100: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: ! Create the singular value solver and display info
102: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: ! ** Create singular value solver context
105: call SVDCreate(PETSC_COMM_WORLD,svd,ierr)
107: ! ** Set operator
108: call SVDSetOperator(svd,A,ierr)
110: ! ** Use thick-restart Lanczos as default solver
111: call SVDSetType(svd,SVDTRLANCZOS,ierr)
113: ! ** Set solver parameters at runtime
114: call SVDSetFromOptions(svd,ierr)
116: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: ! Solve the singular value system
118: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: call SVDSolve(svd,ierr)
121: call SVDGetIterationNumber(svd,its,ierr)
122: if (rank .eq. 0) then
123: write(*,110) its
124: endif
125: 110 format (/' Number of iterations of the method:',I4)
127: ! ** Optional: Get some information from the solver and display it
128: call SVDGetType(svd,tname,ierr)
129: if (rank .eq. 0) then
130: write(*,120) tname
131: endif
132: 120 format (' Solution method: ',A)
133: call SVDGetDimensions(svd,nsv,PETSC_NULL_INTEGER, &
134: & PETSC_NULL_INTEGER,ierr)
135: if (rank .eq. 0) then
136: write(*,130) nsv
137: endif
138: 130 format (' Number of requested singular values:',I2)
139: call SVDGetTolerances(svd,tol,maxit,ierr)
140: if (rank .eq. 0) then
141: write(*,140) tol, maxit
142: endif
143: 140 format (' Stopping condition: tol=',1P,E10.4,', maxit=',I4)
145: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: ! Display solution and clean up
147: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149: ! ** Get number of converged singular triplets
150: call SVDGetConverged(svd,nconv,ierr)
151: if (rank .eq. 0) then
152: write(*,150) nconv
153: endif
154: 150 format (' Number of converged approximate singular triplets:',I2/)
156: ! ** Display singular values and relative errors
157: if (nconv.gt.0) then
158: if (rank .eq. 0) then
159: write(*,*) ' sigma relative error'
160: write(*,*) ' ----------------- ------------------'
161: endif
162: do i=0,nconv-1
163: ! ** Get converged singular triplet: i-th singular value is stored in sigma
164: call SVDGetSingularTriplet(svd,i,sigma,PETSC_NULL_OBJECT, &
165: & PETSC_NULL_OBJECT,ierr)
167: ! ** Compute the relative error associated to each eigenpair
168: call SVDComputeError(svd,i,SVD_ERROR_RELATIVE,error,ierr)
169: if (rank .eq. 0) then
170: write(*,160) sigma, error
171: endif
172: 160 format (1P,' ',E12.4,' ',E12.4)
174: enddo
175: if (rank .eq. 0) then
176: write(*,*)
177: endif
178: endif
180: ! ** Free work space
181: call SVDDestroy(svd,ierr)
182: call MatDestroy(A,ierr)
184: call SlepcFinalize(ierr)
185: end