Xmipp  v3.23.11-Nereus
Classes | Functions
matrix2d.h File Reference
#include <complex>
#include "xmipp_random_mode.h"
#include "xmipp_macros.h"
#include "xmipp_error.h"
#include "xmipp_strings.h"
#include <sys/mman.h>
Include dependency graph for matrix2d.h:

Go to the source code of this file.

Classes

class  Matrix1D< T >
 
class  Matrix2D< T >
 
class  Matrix2D< T >
 

Functions

template<typename T >
void ludcmp (const Matrix2D< T > &A, Matrix2D< T > &LU, Matrix1D< int > &indx, T &d)
 
template<typename T >
void lubksb (const Matrix2D< T > &LU, Matrix1D< int > &indx, Matrix1D< T > &b)
 
template<typename T >
void svdcmp (const Matrix2D< T > &a, Matrix2D< double > &u, Matrix1D< double > &w, Matrix2D< double > &v)
 
void svbksb (Matrix2D< double > &u, Matrix1D< double > &w, Matrix2D< double > &v, Matrix1D< double > &b, Matrix1D< double > &x)
 
void cholesky (const Matrix2D< double > &M, Matrix2D< double > &L)
 
void schur (const Matrix2D< double > &M, Matrix2D< double > &O, Matrix2D< double > &T)
 

Matrices speed up macros

#define MATRIX2D_ARRAY(m)   ((m).mdata)
 
#define FOR_ALL_ELEMENTS_IN_MATRIX2D(m)
 
#define MAT_ELEM(m, i, j)   ((m).mdata[((size_t)i)*(m).mdimx+((size_t)j)])
 
#define MAT_XSIZE(m)   ((m).mdimx)
 
#define MAT_YSIZE(m)   ((m).mdimy)
 
#define MAT_SIZE(m)   ((m).mdim)
 
#define dMij(m, i, j)   MAT_ELEM(m, i, j)
 
#define dMn(m, n)   ((m).mdata[(n)])
 
#define M3x3_BY_V3x1(a, M, b)
 
#define M3x3_BY_M3x3(A, B, C)
 
#define M2x2_BY_V2x1(a, M, b)
 
#define M2x2_BY_CT(M2, M1, k)
 
#define M3x3_BY_CT(M2, M1, k)
 
#define M4x4_BY_CT(M2, M1, k)
 
#define M2x2_INV(Ainv, A)
 
#define M3x3_INV(Ainv, A)
 
#define M4x4_INV(Ainv, A)
 
typedef Matrix2D< double > DMatrix
 
typedef Matrix2D< int > IMatrix
 
template<typename T >
bool operator== (const Matrix2D< T > &op1, const Matrix2D< T > &op2)
 

Matrix Related functions

These functions are not methods of Matrix2D

#define VIA_BILIB
 
void generalizedEigs (const Matrix2D< double > &A, const Matrix2D< double > &B, Matrix1D< double > &D, Matrix2D< double > &P)
 
void firstEigs (const Matrix2D< double > &A, size_t M, Matrix1D< double > &D, Matrix2D< double > &P, bool Pneeded=true)
 
void lastEigs (const Matrix2D< double > &A, size_t M, Matrix1D< double > &D, Matrix2D< double > &P)
 
void eigsBetween (const Matrix2D< double > &A, size_t I1, size_t I2, Matrix1D< double > &D, Matrix2D< double > &P)
 
void allEigs (const Matrix2D< double > &A, std::vector< std::complex< double > > &eigs)
 
void connectedComponentsOfUndirectedGraph (const Matrix2D< double > &G, Matrix1D< int > &component)
 
template<typename T1 , typename T2 >
void typeCast (const Matrix2D< T1 > &v1, Matrix2D< T2 > &v2)
 
template<typename T1 >
void typeCast (const Matrix2D< T1 > &v1, Matrix2D< T1 > &v2)
 
void orthogonalizeColumnsGramSchmidt (Matrix2D< double > &M)
 
void normalizeColumns (Matrix2D< double > &A)
 
void normalizeColumnsBetween0and1 (Matrix2D< double > &A)
 
void subtractColumnMeans (Matrix2D< double > &A)
 
void matrixOperation_AtA (const Matrix2D< double > &A, Matrix2D< double > &B)
 
void matrixOperation_AAt (const Matrix2D< double > &A, Matrix2D< double > &C)
 
void matrixOperation_AB (const Matrix2D< double > &A, const Matrix2D< double > &B, Matrix2D< double > &C)
 
void matrixOperation_Ax (const Matrix2D< double > &A, const Matrix1D< double > &x, Matrix1D< double > &y)
 
void matrixOperation_ABt (const Matrix2D< double > &A, const Matrix2D< double > &B, Matrix2D< double > &C)
 
void matrixOperation_AtB (const Matrix2D< double > &A, const Matrix2D< double > &B, Matrix2D< double > &C)
 
void matrixOperation_Atx (const Matrix2D< double > &A, const Matrix1D< double > &x, Matrix1D< double > &y)
 
void matrixOperation_AtBt (const Matrix2D< double > &A, const Matrix2D< double > &B, Matrix2D< double > &C)
 
void matrixOperation_XtAX_symmetric (const Matrix2D< double > &X, const Matrix2D< double > &A, Matrix2D< double > &B)
 
void matrixOperation_IplusA (Matrix2D< double > &A)
 
void matrixOperation_IminusA (Matrix2D< double > &A)
 
void eraseFirstColumn (Matrix2D< double > &A)
 
void keepColumns (Matrix2D< double > &A, int j0, int jF)
 

Function Documentation

◆ cholesky()

void cholesky ( const Matrix2D< double > &  M,
Matrix2D< double > &  L 
)

Cholesky decomposition. Given M, this function decomposes M as M=L*L^t where L is a lower triangular matrix. M must be positive semi-definite.

Definition at line 160 of file matrix2d.cpp.

161 {
162  L=M;
164  p.initZeros(MAT_XSIZE(M));
167  if (i==j)
168  MAT_ELEM(L,i,j)=VEC_ELEM(p,i);
169  else if (i<j)
170  MAT_ELEM(L,i,j)=0.0;
171 }
#define FOR_ALL_ELEMENTS_IN_MATRIX2D(m)
Definition: matrix2d.h:104
#define VEC_ELEM(v, i)
Definition: matrix1d.h:245
void choldc(double *a, int n, double *p)
T * adaptForNumericalRecipes2() const
Definition: matrix2d.h:875
#define i
#define MAT_ELEM(m, i, j)
Definition: matrix2d.h:116
void initZeros()
Definition: matrix1d.h:592
#define j
#define MAT_XSIZE(m)
Definition: matrix2d.h:120
T * adaptForNumericalRecipes() const
Definition: matrix1d.h:844

◆ schur()

void schur ( const Matrix2D< double > &  M,
Matrix2D< double > &  O,
Matrix2D< double > &  T 
)

Schur decomposition. Given M, this function decomposes M as M = O*T*O' where O is an orthogonal matrix.

Definition at line 251 of file matrix2d.cpp.

252 {
255  bool ok=rmatrixschur(a, MAT_YSIZE(M), s);
256  if (!ok)
257  REPORT_ERROR(ERR_NUMERICAL,"Could not perform Schur decomposition");
258  O.resizeNoCopy(M);
259  T.resizeNoCopy(M);
261  {
262  MAT_ELEM(O,i,j)=s(i,j);
263  MAT_ELEM(T,i,j)=a(i,j);
264  }
265 }
#define FOR_ALL_ELEMENTS_IN_MATRIX2D(m)
Definition: matrix2d.h:104
#define MAT_YSIZE(m)
Definition: matrix2d.h:124
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
#define i
void resizeNoCopy(int Ydim, int Xdim)
Definition: matrix2d.h:534
#define MAT_ELEM(m, i, j)
Definition: matrix2d.h:116
bool rmatrixschur(real_2d_array &a, const ae_int_t n, real_2d_array &s)
Definition: linalg.cpp:7050
Error related to numerical calculation.
Definition: xmipp_error.h:179
#define j
#define MAT_XSIZE(m)
Definition: matrix2d.h:120
void setcontent(ae_int_t irows, ae_int_t icols, const double *pContent)
Definition: ap.cpp:6660
#define MATRIX2D_ARRAY(m)
Definition: matrix2d.h:89
doublereal * a