26#ifndef SCIMATH_MATRIXMATHLA_H 
   27#define SCIMATH_MATRIXMATHLA_H 
   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> 
  101#if !defined(NEED_FORTRAN_UNDERSCORES) 
  102#define NEED_FORTRAN_UNDERSCORES 1 
  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_ 
  118#define spotri spotri_ 
  119#define dpotri dpotri_ 
  120#define cpotri cpotri_ 
  121#define zpotri zpotri_ 
  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);
 
  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);
 
  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);
 
  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); }
 
 
  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); }  
 
 
  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); }  
 
 
  220#ifndef CASACORE_NO_AUTO_TEMPLATES 
  221#include <casacore/scimath/Mathematics/MatrixMathLA.tcc> 
this file contains all the compiler specific defines
 
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)
 
T determinate(const Matrix< T > &in)
 
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)
 
Matrix< T > invert(const Matrix< T > &in)