casacore
Loading...
Searching...
No Matches
MatrixMathLA.h
Go to the documentation of this file.
1//# MatrixMath.h: The Casacore linear algebra functions
2//# Copyright (C) 1994,1995,1996,1999,2000,2002
3//# Associated Universities, Inc. Washington DC, USA.
4//#
5//# This library is free software; you can redistribute it and/or modify it
6//# under the terms of the GNU Library General Public License as published by
7//# the Free Software Foundation; either version 2 of the License, or (at your
8//# option) any later version.
9//#
10//# This library is distributed in the hope that it will be useful, but WITHOUT
11//# ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
12//# FITNESS FOR A PARTICULAR PURPOSE. See the GNU Library General Public
13//# License for more details.
14//#
15//# You should have received a copy of the GNU Library General Public License
16//# along with this library; if not, write to the Free Software Foundation,
17//# Inc., 675 Massachusetts Ave, Cambridge, MA 02139, USA.
18//#
19//# Correspondence concerning AIPS++ should be addressed as follows:
20//# Internet email: casa-feedback@nrao.edu.
21//# Postal address: AIPS++ Project Office
22//# National Radio Astronomy Observatory
23//# 520 Edgemont Road
24//# Charlottesville, VA 22903-2475 USA
25
26#ifndef SCIMATH_MATRIXMATHLA_H
27#define SCIMATH_MATRIXMATHLA_H
28
29
30#include <casacore/casa/aips.h>
31#include <casacore/casa/Arrays/Vector.h>
32#include <casacore/casa/Arrays/Matrix.h>
33#include <casacore/casa/BasicSL/Complex.h>
34
35
36namespace casacore { //# NAMESPACE CASACORE - BEGIN
37
38//<summary>
39// Linear algebra functions on Vectors and Matrices.
40// </summary>
41//
42// <reviewed reviewer="UNKNOWN" date="before2004/08/25" tests="tLinAlgebra">
43//
44// <linkfrom anchor="Linear Algebra" classes="Vector Matrix">
45// <here>Linear Algebra</here> -- Linear algebra functions
46// on Vectors and Matrices.
47// </linkfrom>
48//
49//<group name="Linear Algebra">
51// Routines which calculate the inverse of a matrix. The inverse is very
52// often the worst way to do a calculation. Nevertheless it is often
53// convenient. The present implementation uses LU decomposition implemented
54// by LAPACK. The determinate can be calculated "for free" as it is the
55// product of the diagonal terms after decomposition. If the input matrix is
56// singular, a matrix with no rows or columns is returned. <src>in</src>
57// must be a square matrix.
58// <note role="warning">This function will only work for complex types if
59// Complex and DComplex map onto their FORTRAN counterparts.</note>
60//# We could special case small matrices for efficiency.
61//<group>
62template<class T> void invert(Matrix<T> & out, T& determinate,
63 const Matrix<T> &in);
64template<class T> Matrix<T> invert(const Matrix<T> &in);
65template<class T> T determinate(const Matrix<T> &in);
66//</group>
67
68// This function inverts a symmetric positive definite matrix. It is
69// written in C++, so it should work with any data type for which
70// operators +, -, *, /, =, and sqrt are defined. The function uses
71// the Cholesky decomposition method to invert the matrix. Cholesky
72// decomposition is about a factor of 2 better than LU decomposition
73// where symmetry is ignored.
74template<class T> void invertSymPosDef(Matrix<T> & out, T& determinate,
75 const Matrix<T> &in);
76template<class T> Matrix<T> invertSymPosDef(const Matrix<T> &in);
77
78//</group>
79
80//# These are actually used by invertSymPosDef. They will not
81//# normally be called by the end user.
82
83//# This function performs Cholesky decomposition.
84//# A is a positive-definite symmetric matrix. Only the upper triangle of
85//# A is needed on input. On output, the lower triangle of A contains the
86//# Cholesky factor L. The diagonal elements of L are returned in vector
87//# diag.
88template<class T> void CholeskyDecomp(Matrix<T> &A, Vector<T> &diag);
89
90//# Solve linear equation A*x = b, where A positive-definite symmetric.
91//# On input, A contains Cholesky factor L in its low triangle except the
92//# diagonal elements which are in vector diag. On return x contains the
93//# solution. b and x can be the same vector to save memory space.
94template<class T> void CholeskySolve(Matrix<T> &A, Vector<T> &diag,
95 Vector<T> &b, Vector<T> &x);
96
97
98//# These are the LAPACK routines actually used by invert. They will not
99//# normally be called by the end user.
100
101#if !defined(NEED_FORTRAN_UNDERSCORES)
102#define NEED_FORTRAN_UNDERSCORES 1
103#endif
104
105#if NEED_FORTRAN_UNDERSCORES
106#define sgetrf sgetrf_
107#define dgetrf dgetrf_
108#define cgetrf cgetrf_
109#define zgetrf zgetrf_
110#define sgetri sgetri_
111#define dgetri dgetri_
112#define cgetri cgetri_
113#define zgetri zgetri_
114#define sposv sposv_
115#define dposv dposv_
116#define cposv cposv_
117#define zposv zposv_
118#define spotri spotri_
119#define dpotri dpotri_
120#define cpotri cpotri_
121#define zpotri zpotri_
122#endif
123
124extern "C" {
125 void sgetrf(const int *m, const int *n, float *a, const int *lda,
126 int *ipiv, int *info);
127 void dgetrf(const int *m, const int *n, double *a, const int *lda,
128 int *ipiv, int *info);
129 void cgetrf(const int *m, const int *n, Complex *a, const int *lda,
130 int *ipiv, int *info);
131 void zgetrf(const int *m, const int *n, DComplex *a, const int *lda,
132 int *ipiv, int *info);
133 void sgetri(const int *m, float *a, const int *lda, const int *ipiv,
134 float *work, const int *lwork, int *info);
135 void dgetri(const int *m, double *a, const int *lda, const int *ipiv,
136 double *work, const int *lwork, int *info);
137 void cgetri(const int *m, Complex *a, const int *lda, const int *ipiv,
138 Complex *work, const int *lwork, int *info);
139 void zgetri(const int *m, DComplex *a, const int *lda, const int *ipiv,
140 DComplex *work, const int *lwork, int *info);
141
142
143 void sposv(const char *uplo, const int *n, const int* nrhs, float *a,
144 const int *lda, float *b, const int *ldb, int *info);
145 void dposv(const char *uplo, const int *n, const int* nrhs, double *a,
146 const int *lda, double *b, const int *ldb, int *info);
147 void cposv(const char *uplo, const int *n, const int* nrhs, Complex *a,
148 const int *lda, Complex *b, const int *ldb, int *info);
149 void zposv(const char *uplo, const int *n, const int* nrhs, DComplex *a,
150 const int *lda, DComplex *b, const int *ldb, int *info);
151
152
153 void spotri(const char *uplo, const int *n, float *a,
154 const int *lda, int *info);
155 void dpotri(const char *uplo, const int *n, double *a,
156 const int *lda, int *info);
157 void cpotri(const char *uplo, const int *n, Complex *a,
158 const int *lda, int *info);
159 void zpotri(const char *uplo, const int *n, DComplex *a,
160 const int *lda, int *info);
161
162}
163
164//# Overloaded versions of the above to make templating work more easily
165inline void getrf(const int *m, const int *n, float *a, const int *lda,
166 int *ipiv, int *info)
167 { sgetrf(m, n, a, lda, ipiv, info); }
168inline void getrf(const int *m, const int *n, double *a, const int *lda,
169 int *ipiv, int *info)
170 { dgetrf(m, n, a, lda, ipiv, info); }
171inline void getrf(const int *m, const int *n, Complex *a, const int *lda,
172 int *ipiv, int *info)
173 { cgetrf(m, n, a, lda, ipiv, info); }
174inline void getrf(const int *m, const int *n, DComplex *a, const int *lda,
175 int *ipiv, int *info)
176 { zgetrf(m, n, a, lda, ipiv, info); }
177inline void getri(const int *m, float *a, const int *lda, const int *ipiv,
178 float *work, const int *lwork, int *info)
179 { sgetri(m, a, lda, ipiv, work, lwork, info); }
180inline void getri(const int *m, double *a, const int *lda, const int *ipiv,
181 double *work, const int *lwork, int *info)
182 { dgetri(m, a, lda, ipiv, work, lwork, info); }
183inline void getri(const int *m, Complex *a, const int *lda, const int *ipiv,
184 Complex *work, const int *lwork, int *info)
185 { cgetri(m, a, lda, ipiv, work, lwork, info); }
186inline void getri(const int *m, DComplex *a, const int *lda, const int *ipiv,
187 DComplex *work, const int *lwork, int *info)
188 { zgetri(m, a, lda, ipiv, work, lwork, info); }
189
190inline void posv(const char *uplo, const int *n, const int* nrhs, float *a,
191 const int *lda, float *b, const int *ldb, int *info)
192 { sposv(uplo, n, nrhs, a, lda, b, ldb, info); }
193inline void posv(const char *uplo, const int *n, const int* nrhs, double *a,
194 const int *lda, double *b, const int *ldb, int *info)
195 { dposv(uplo, n, nrhs, a, lda, b, ldb, info); }
196inline void posv(const char *uplo, const int *n, const int* nrhs, Complex *a,
197 const int *lda, Complex *b, const int *ldb, int *info)
198 { cposv(uplo, n, nrhs, a, lda, b, ldb, info); }
199inline void posv(const char *uplo, const int *n, const int* nrhs, DComplex *a,
200 const int *lda, DComplex *b, const int *ldb, int *info)
201 { zposv(uplo, n, nrhs, a, lda, b, ldb, info); }
202
203inline void potri(const char *uplo, const int *n, float *a,
204 const int *lda, int *info)
205 { spotri(uplo, n, a, lda, info); }
206inline void potri(const char *uplo, const int *n, double *a,
207 const int *lda, int *info)
208 { dpotri(uplo, n, a, lda, info); }
209inline void potri(const char *uplo, const int *n, Complex *a,
210 const int *lda, int *info)
211 { cpotri(uplo, n, a, lda, info); }
212inline void potri(const char *uplo, const int *n, DComplex *a,
213 const int *lda, int *info)
214 { zpotri(uplo, n, a, lda, info); }
215
216
217
218} //# NAMESPACE CASACORE - END
219
220#ifndef CASACORE_NO_AUTO_TEMPLATES
221#include <casacore/scimath/Mathematics/MatrixMathLA.tcc>
222#endif //# CASACORE_NO_AUTO_TEMPLATES
223#endif
#define zpotri
#define sgetrf
#define zgetri
#define sgetri
#define dgetri
#define cposv
#define spotri
#define cgetrf
#define dposv
#define zgetrf
#define cpotri
#define sposv
#define dgetrf
#define zposv
#define cgetri
#define dpotri
this file contains all the compiler specific defines
Definition mainpage.dox:28
void CholeskySolve(Matrix< T > &A, Vector< T > &diag, Vector< T > &b, Vector< T > &x)
void CholeskyDecomp(Matrix< T > &A, Vector< T > &diag)
void getri(const int *m, float *a, const int *lda, const int *ipiv, float *work, const int *lwork, int *info)
void potri(const char *uplo, const int *n, float *a, const int *lda, int *info)
void posv(const char *uplo, const int *n, const int *nrhs, float *a, const int *lda, float *b, const int *ldb, int *info)
void getrf(const int *m, const int *n, float *a, const int *lda, int *ipiv, int *info)
void invert(Matrix< T > &out, T &determinate, const Matrix< T > &in)
Routines which calculate the inverse of a matrix.
void invertSymPosDef(Matrix< T > &out, T &determinate, const Matrix< T > &in)
This function inverts a symmetric positive definite matrix.
Matrix< T > invertSymPosDef(const Matrix< T > &in)