tensor-0.1.0
 All Data Structures Namespaces Functions Variables Typedefs Enumerations Enumerator Groups Pages
linalg.h
1 // -*- mode: c++; fill-column: 80; c-basic-offset: 2; indent-tabs-mode: nil -*-
2 /*
3  Copyright (c) 2010 Juan Jose Garcia Ripoll
4 
5  Tensor is free software; you can redistribute it and/or modify it
6  under the terms of the GNU Library General Public License as published
7  by the Free Software Foundation; either version 2 of the License, or
8  (at your option) any later version.
9 
10  This program is distributed in the hope that it will be useful,
11  but WITHOUT ANY WARRANTY; without even the implied warranty of
12  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  GNU Library General Public License for more details.
14 
15  You should have received a copy of the GNU General Public License along
16  with this program; if not, write to the Free Software Foundation, Inc.,
17  51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
18 */
19 
20 #ifndef TENSOR_LINALG_H
21 #define TENSOR_LINALG_H
22 
23 #include <tensor/tensor.h>
24 #include <tensor/sparse.h>
25 
28 namespace linalg {
29 
30  using tensor::RTensor;
31  using tensor::CTensor;
32  using tensor::RSparse;
33  using tensor::CSparse;
34 
35  const RTensor solve(const RTensor &A, const RTensor &B);
36  const CTensor solve(const CTensor &A, const CTensor &B);
37 
38  const RTensor solve_with_svd(const RTensor &A, const RTensor &B, double tol = 0.0);
39  const CTensor solve_with_svd(const CTensor &A, const CTensor &B, double tol = 0.0);
40 
41  extern bool accurate_svd;
42 
43  #define SVD_ECONOMIC true
44  RTensor svd(RTensor A, RTensor *pU = 0, RTensor *pVT = 0, bool economic = 0);
45  RTensor svd(CTensor A, CTensor *pU = 0, CTensor *pVT = 0, bool economic = 0);
46 
48  const CTensor eig(const RTensor &A, CTensor *R = 0, CTensor *L = 0);
49 
51  const CTensor eig(const CTensor &A, CTensor *R = 0, CTensor *L = 0);
52 
53  double eig_power_right(const RTensor &A, RTensor *vector,
54  size_t iter = 0, double tol = 1e-11);
55  double eig_power_left(const RTensor &A, RTensor *vector,
56  size_t iter = 0, double tol = 1e-11);
57  tensor::cdouble eig_power_right(const CTensor &A, CTensor *vector,
58  size_t iter = 0, double tol = 1e-11);
59  tensor::cdouble eig_power_left(const CTensor &A, CTensor *vector,
60  size_t iter = 0, double tol = 1e-11);
61 
62  RTensor eig_sym(const RTensor &A, RTensor *pR = 0);
63  RTensor eig_sym(const CTensor &A, CTensor *pR = 0);
64 
65  const RTensor expm(const RTensor &A, unsigned int order = 7);
66  const CTensor expm(const CTensor &A, unsigned int order = 7);
67 
69  enum EigType {
78  };
79 
81  CTensor eigs(const CTensor &A, int eig_type, size_t neig,
82  CTensor *vectors = NULL,
83  const tensor::cdouble *initial_guess = NULL);
84 
86  CTensor eigs(const CSparse &A, int eig_type, size_t neig,
87  CTensor *vectors = NULL,
88  const tensor::cdouble *initial_guess = NULL);
89 
91  RTensor eigs(const RTensor &A, int eig_type, size_t neig,
92  RTensor *vectors = NULL,
93  const double *initial_guess = NULL);
94 
96  RTensor eigs(const RSparse &A, int eig_type, size_t neig,
97  RTensor *vectors = NULL,
98  const double *initial_guess = NULL);
99 
101  RTensor eigs_sym(const RTensor &A, int eig_type, size_t neig,
102  RTensor *vectors = NULL,
103  const double *initial_guess = NULL);
104 
106  RTensor eigs_sym(const CTensor &A, int eig_type, size_t neig,
107  CTensor *vectors = NULL,
108  const tensor::cdouble *initial_guess = NULL);
109 
111  RTensor eigs_sym(const RSparse &A, int eig_type, size_t neig,
112  RTensor *vectors = NULL,
113  const double *initial_guess = NULL);
114 
116  RTensor eigs_sym(const CSparse &A, int eig_type, size_t neig,
117  CTensor *vectors = NULL,
118  const tensor::cdouble *initial_guess = NULL);
119 
120 } // namespace linalg
121 
122 
123 #endif // !TENSOR_LINALG_H
A sparse matrix.
Definition: sparse.h:37
const CTensor eig(const RTensor &A, CTensor *R=0, CTensor *L=0)
Eigenvalue decomposition of a real matrix.
Definition: eig_d.cc:45
double eig_power_left(const RTensor &A, RTensor *vector, size_t iter=0, double tol=1e-11)
Left eigenvalue and eigenvector with the largest absolute value, computed using the power method...
Definition: eig_power_d.cc:41
const RTensor expm(const RTensor &A, unsigned int order=7)
Compute the exponential of a real matrix.
Definition: expm_d.cc:42
Real Tensor with elements of type "double".
Definition: tensor.h:349
CTensor eigs(const CTensor &A, int eig_type, size_t neig, CTensor *vectors=NULL, const tensor::cdouble *initial_guess=NULL)
Find out a few eigenvalues and eigenvectors of a complex nonsymmetric matrix.
RTensor eig_sym(const RTensor &A, RTensor *pR=0)
Eigenvalue decomposition of a real matrix.
Definition: eig_sym_d.cc:44
Complex Tensor with elements of type "cdouble".
Definition: tensor.h:435
RTensor eigs_sym(const RTensor &A, int eig_type, size_t neig, RTensor *vectors=NULL, const double *initial_guess=NULL)
Find out a few eigenvalues and eigenvectors of a symmetric real matrix.
const RTensor solve(const RTensor &A, const RTensor &B)
Solve a real linear system of equations by Gauss-Seidel method.
Definition: solve_d.cc:39
double eig_power_right(const RTensor &A, RTensor *vector, size_t iter=0, double tol=1e-11)
Right eigenvalue and eigenvector with the largest absolute value, computed using the power method...
Definition: eig_power_d.cc:30
EigType
Type of eigenvalues that eigs and Arpack compute.
Definition: linalg.h:69
const RTensor solve_with_svd(const RTensor &A, const RTensor &B, double tol=0.0)
Solution of a linear system of equations using Penrose's pseudoinvese.
RTensor svd(RTensor A, RTensor *pU=0, RTensor *pVT=0, bool economic=0)
Singular value decomposition of a real matrix.
Definition: svd_d.cc:50