1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96
|
/*
ARPACK++ v1.2 2/18/2000
c++ interface to ARPACK code.
MODULE DSVD.cc.
Example program that illustrates how to determine the largest
singular values of a matrix using arpack++.
1) Problem description:
In this example, Arpack++ is called to solve the symmetric problem:
(A'*A)*v = sigma*v
where A is an m by n real matrix.
This formulation is appropriate when m >= n.
The roles of A and A' must be reversed in the case that m < n.
2) Data structure used to represent the matrix:
A is stored columnwise in a vector called A.
3) Included header files:
File Contents
----------- --------------------------------------------
dnmatrxw.h MatrixW, a function that generates matrix A.
ardnsmat.h The ARdsNonSymMatrix class definition.
arssym.h The ARSymStdEig class definition.
4) ARPACK Authors:
Richard Lehoucq
Kristyn Maschhoff
Danny Sorensen
Chao Yang
Dept. of Computational & Applied Mathematics
Rice University
Houston, Texas
*/
#include <iostream>
#include <cmath>
#include "arssym.h"
#include "dnmatrxw.h"
#include "ardnsmat.h"
int main()
{
// Defining variables;
int i;
int m; // Number of rows in A.
int n; // Number of columns in A.
double* valA; // Pointer to an array that stores the elements of A.
double* svalue = new double[4];
// Creating a matrix.
m = 500;
n = 100;
MatrixW(m, n, valA);
// Using ARdsNonSymMatrix to store matrix information and to
// perform the product A'Ax (LU decomposition is not used).
ARdsNonSymMatrix<double, double> A(m, n, valA);
// Defining what we need: eigenvalues with largest magnitude.
ARSymStdEig<double, ARdsNonSymMatrix<double, double> >
dprob(n, 4L, &A, &ARdsNonSymMatrix<double, double>::MultMtMv);
// Finding eigenvalues.
dprob.Eigenvalues(svalue);
// Calculating singular values.
for (i=0; i<4; i++) {
svalue[i] = sqrt(svalue[i]);
}
// Printing some information about the problem.
std::cout << std::endl << "Testing ARPACK++ class ARSymStdEig" << std::endl;
std::cout << "Obtaining singular values by solving (A'*A)*v = sigma*v" << std::endl;
std::cout << std::endl << "greatest singular values: " << std::endl;
for (i=0; i<4; i++) {
std::cout << " sigma [" << i+1 << "]: " << svalue[i] << std::endl;
}
} // main.
|