casa  $Rev:20696$
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
MatrixMathLA.h
Go to the documentation of this file.
00001 //# MatrixMath.h: The AIPS++ linear algebra functions
00002 //# Copyright (C) 1994,1995,1996,1999,2000,2002
00003 //# Associated Universities, Inc. Washington DC, USA.
00004 //#
00005 //# This library is free software; you can redistribute it and/or modify it
00006 //# under the terms of the GNU Library General Public License as published by
00007 //# the Free Software Foundation; either version 2 of the License, or (at your
00008 //# option) any later version.
00009 //#
00010 //# This library is distributed in the hope that it will be useful, but WITHOUT
00011 //# ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00012 //# FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Library General Public
00013 //# License for more details.
00014 //#
00015 //# You should have received a copy of the GNU Library General Public License
00016 //# along with this library; if not, write to the Free Software Foundation,
00017 //# Inc., 675 Massachusetts Ave, Cambridge, MA 02139, USA.
00018 //#
00019 //# Correspondence concerning AIPS++ should be addressed as follows:
00020 //#        Internet email: aips2-request@nrao.edu.
00021 //#        Postal address: AIPS++ Project Office
00022 //#                        National Radio Astronomy Observatory
00023 //#                        520 Edgemont Road
00024 //#                        Charlottesville, VA 22903-2475 USA
00025 //#
00026 //# $Id: MatrixMathLA.h 21024 2011-03-01 11:46:18Z gervandiepen $
00027 
00028 #ifndef SCIMATH_MATRIXMATHLA_H
00029 #define SCIMATH_MATRIXMATHLA_H
00030 
00031 
00032 #include <casa/aips.h>
00033 #include <casa/Arrays/Vector.h>
00034 #include <casa/Arrays/Matrix.h>
00035 #include <casa/BasicSL/Complex.h>
00036 
00037 
00038 namespace casa { //# NAMESPACE CASA - BEGIN
00039 
00040 //<summary>
00041 //    Linear algebra functions on Vectors and Matrices.
00042 // </summary>
00043 //
00044 // <reviewed reviewer="UNKNOWN" date="before2004/08/25" tests="tLinAlgebra">
00045 //
00046 // <linkfrom anchor="Linear Algebra" classes="Vector Matrix">
00047 //    <here>Linear Algebra</here> -- Linear algebra functions
00048 //     on Vectors and Matrices.
00049 // </linkfrom>
00050 //
00051 //<group name="Linear Algebra">
00052 
00053 // Routines which calculate the inverse of a matrix.  The inverse is very
00054 // often the worst way to do a calculation.  Nevertheless it is often
00055 // convenient. The present implementation uses LU decomposition implemented
00056 // by LAPACK. The determinate can be calculated "for free" as it is the
00057 // product of the diagonal terms after decomposition.  If the input matrix is
00058 // singular, a matrix with no rows or columns is returned.  <src>in</src>
00059 // must be a square matrix.  
00060 // <note role="warning">This function will only work for complex types if
00061 // Complex and DComplex map onto their FORTRAN counterparts.</note>
00062 //# We could special case small matrices for efficiency. 
00063 //<group>
00064 template<class T> void invert(Matrix<T> & out, T& determinate, 
00065                               const Matrix<T> &in);
00066 template<class T> Matrix<T> invert(const Matrix<T> &in);
00067 template<class T> T determinate(const Matrix<T> &in);
00068 //</group>
00069 
00070 // This function inverts a symmetric positive definite matrix.  It is
00071 // written in C++, so it should work with any data type for which
00072 // operators +, -, *, /, =, and sqrt are defined.  The function uses
00073 // the Cholesky decomposition method to invert the matrix.  Cholesky
00074 // decomposition is about a factor of 2 better than LU decomposition
00075 // where symmetry is ignored. 
00076 template<class T> void invertSymPosDef(Matrix<T> & out, T& determinate, 
00077                                        const Matrix<T> &in);
00078 template<class T> Matrix<T> invertSymPosDef(const Matrix<T> &in);
00079 
00080 //</group>
00081 
00082 //# These are actually used by invertSymPosDef. They will not
00083 //# normally be called by the end user.
00084 
00085 //# This function performs Cholesky decomposition.
00086 //# A is a positive-definite symmetric matrix. Only the upper triangle of
00087 //# A is needed on input. On output, the lower triangle of A contains the
00088 //# Cholesky factor L.  The diagonal elements of L are returned in vector
00089 //# diag.
00090 template<class T> void CholeskyDecomp(Matrix<T> &A, Vector<T> &diag);
00091 
00092 //# Solve linear equation A*x = b, where A positive-definite symmetric.
00093 //# On input, A contains Cholesky factor L in its low triangle except the
00094 //# diagonal elements which are in vector diag.  On return x contains the
00095 //# solution.  b and x can be the same vector to save memory space.
00096 template<class T> void CholeskySolve(Matrix<T> &A, Vector<T> &diag, 
00097                                      Vector<T> &b, Vector<T> &x);
00098 
00099 
00100 //# These are the LAPACK routines actually used by invert. They will not
00101 //# normally be called by the end user.
00102 
00103 #if !defined(NEED_FORTRAN_UNDERSCORES)
00104 #define NEED_FORTRAN_UNDERSCORES 1
00105 #endif
00106 
00107 #if NEED_FORTRAN_UNDERSCORES
00108 #define sgetrf sgetrf_
00109 #define dgetrf dgetrf_
00110 #define cgetrf cgetrf_
00111 #define zgetrf zgetrf_
00112 #define sgetri sgetri_
00113 #define dgetri dgetri_
00114 #define cgetri cgetri_
00115 #define zgetri zgetri_
00116 #define sposv sposv_
00117 #define dposv dposv_
00118 #define cposv cposv_
00119 #define zposv zposv_
00120 #define spotri spotri_
00121 #define dpotri dpotri_
00122 #define cpotri cpotri_
00123 #define zpotri zpotri_
00124 #endif
00125 
00126 extern "C" {
00127       void sgetrf(const int *m, const int *n, float *a, const int *lda,
00128                   int *ipiv, int *info);
00129       void dgetrf(const int *m, const int *n, double *a, const int *lda,
00130                   int *ipiv, int *info);
00131       void cgetrf(const int *m, const int *n, Complex *a, const int *lda,
00132                   int *ipiv, int *info);
00133       void zgetrf(const int *m, const int *n, DComplex *a, const int *lda,
00134                   int *ipiv, int *info);
00135       void sgetri(const int *m, float *a, const int *lda, const int *ipiv,
00136                   float *work, const int *lwork, int *info);
00137       void dgetri(const int *m, double *a, const int *lda, const int *ipiv,
00138                   double *work, const int *lwork, int *info);
00139       void cgetri(const int *m, Complex *a, const int *lda, const int *ipiv,
00140                   Complex *work, const int *lwork, int *info);
00141       void zgetri(const int *m, DComplex *a, const int *lda, const int *ipiv,
00142                   DComplex *work, const int *lwork, int *info);
00143 
00144 
00145       void sposv(const char *uplo, const int *n, const int* nrhs, float *a, 
00146                  const int *lda, float *b, const int *ldb, int *info);
00147       void dposv(const char *uplo, const int *n, const int* nrhs, double *a, 
00148                  const int *lda, double *b, const int *ldb, int *info);
00149       void cposv(const char *uplo, const int *n, const int* nrhs, Complex *a, 
00150                  const int *lda, Complex *b, const int *ldb, int *info);
00151       void zposv(const char *uplo, const int *n, const int* nrhs, DComplex *a, 
00152                  const int *lda, DComplex *b, const int *ldb, int *info);
00153 
00154 
00155       void spotri(const char *uplo, const int *n, float *a, 
00156                   const int *lda, int *info);
00157       void dpotri(const char *uplo, const int *n, double *a, 
00158                   const int *lda, int *info);
00159       void cpotri(const char *uplo, const int *n, Complex *a, 
00160                   const int *lda, int *info);
00161       void zpotri(const char *uplo, const int *n, DComplex *a, 
00162                   const int *lda, int *info);
00163 
00164 }
00165 
00166 //# Overloaded versions of the above to make templating work more easily
00167 inline void getrf(const int *m, const int *n, float *a, const int *lda,
00168                   int *ipiv, int *info)
00169    { sgetrf(m, n, a, lda, ipiv, info); }
00170 inline void getrf(const int *m, const int *n, double *a, const int *lda,
00171                   int *ipiv, int *info)
00172    { dgetrf(m, n, a, lda, ipiv, info); }
00173 inline void getrf(const int *m, const int *n, Complex *a, const int *lda,
00174                   int *ipiv, int *info)
00175    { cgetrf(m, n, a, lda, ipiv, info); }
00176 inline void getrf(const int *m, const int *n, DComplex *a, const int *lda,
00177                   int *ipiv, int *info)
00178    { zgetrf(m, n, a, lda, ipiv, info); }
00179 inline void getri(const int *m, float *a, const int *lda, const int *ipiv,
00180                   float *work, const int *lwork, int *info)
00181    { sgetri(m, a, lda, ipiv, work, lwork, info); }
00182 inline void getri(const int *m, double *a, const int *lda, const int *ipiv,
00183                   double *work, const int *lwork, int *info)
00184    { dgetri(m, a, lda, ipiv, work, lwork, info); }
00185 inline void getri(const int *m, Complex *a, const int *lda, const int *ipiv,
00186                   Complex *work, const int *lwork, int *info)
00187    { cgetri(m, a, lda, ipiv, work, lwork, info); }
00188 inline void getri(const int *m, DComplex *a, const int *lda, const int *ipiv,
00189                   DComplex *work, const int *lwork, int *info)
00190    { zgetri(m, a, lda, ipiv, work, lwork, info); }
00191 
00192 inline void posv(const char *uplo, const int *n, const int* nrhs, float *a, 
00193                  const int *lda, float *b, const int *ldb, int *info)
00194    { sposv(uplo, n, nrhs, a, lda, b, ldb, info); }  
00195 inline void posv(const char *uplo, const int *n, const int* nrhs, double *a, 
00196                  const int *lda, double *b, const int *ldb, int *info)
00197    { dposv(uplo, n, nrhs, a, lda, b, ldb, info); }  
00198 inline void posv(const char *uplo, const int *n, const int* nrhs, Complex *a, 
00199                  const int *lda, Complex *b, const int *ldb, int *info)
00200    { cposv(uplo, n, nrhs, a, lda, b, ldb, info); }  
00201 inline void posv(const char *uplo, const int *n, const int* nrhs, DComplex *a, 
00202                  const int *lda, DComplex *b, const int *ldb, int *info)
00203    { zposv(uplo, n, nrhs, a, lda, b, ldb, info); }  
00204 
00205 inline void potri(const char *uplo, const int *n, float *a, 
00206                   const int *lda, int *info)
00207    { spotri(uplo, n, a, lda, info); }  
00208 inline void potri(const char *uplo, const int *n, double *a, 
00209                   const int *lda, int *info)
00210    { dpotri(uplo, n, a, lda, info); }  
00211 inline void potri(const char *uplo, const int *n, Complex *a, 
00212                   const int *lda, int *info)
00213    { cpotri(uplo, n, a, lda, info); }  
00214 inline void potri(const char *uplo, const int *n, DComplex *a, 
00215                   const int *lda, int *info)
00216    { zpotri(uplo, n, a, lda, info); }  
00217 
00218 
00219 
00220 } //# NAMESPACE CASA - END
00221 
00222 #ifndef CASACORE_NO_AUTO_TEMPLATES
00223 #include <scimath/Mathematics/MatrixMathLA.tcc>
00224 #endif //# CASACORE_NO_AUTO_TEMPLATES
00225 #endif