C++ main module for mmsd Package  1.0
Functions
dfullmatrix_functions.h File Reference
#include "lapack_functions.h"
Include dependency graph for dfullmatrix_functions.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

void DoubleFullMatrixVectorProduct (const tLVectorIndex &lx, const tLVectorIncrement &incx, const double *x, const tBoolean &isTransposed, const tLVectorIndex &nRows, const tLVectorIndex &nCols, const tLVectorIndex &ldaA, const double *A, const lapack_real &alpha, const lapack_real &beta, const tLVectorIndex &ly, const tLVectorIncrement &incy, double *y)
 compute Y=beta.Y + alpha op(A). X where op(A)=A or op(A)=tA More...
 
void DoubleFullUpperMatrixVectorProduct (const tBoolean &isTransA, const tLVectorIndex &nRowsA, const tLVectorIndex &ldA, const double *A, const tLVectorIndex &nX, const tLVectorIncrement &incX, double *x)
 compute X=op(A). X where op(A)=A or op(A)=tA where A is Upper matrix More...
 
void DoubleFullMatrixMatrixProduct (const tLVectorIndex &nARows, const tLVectorIndex &nACols, const tLVectorIndex &ldA, const double *A, const tLVectorIndex &nBRows, const tLVectorIndex &nBCols, const tLVectorIndex &ldB, const double *B, const tBoolean &isTrA, const tBoolean &isTrB, const lapack_real &alpha, const lapack_real &beta, const tLVectorIndex &nCRows, const tLVectorIndex &nCCols, const tLVectorIndex &ldC, double *C)
 
int DoubleFullMatrixEigenValues (const tLVectorIndex &nRows, const tLVectorIndex &nCols, const tLVectorIndex &ld, double *A, const tLVectorIndex &nU, double *U, double *C)
 
int DoubleFullMatrixLUFactorization (const tLVectorIndex &nRows, const tLVectorIndex &nCols, const tLVectorIndex &ld, double *A, lapack_int *LU_pivots)
 
bool DoubleFullMatrixInverse (const tLVectorIndex &nRows, const tLVectorIndex &nCols, const tLVectorIndex &ld, const lapack_int *ipiv, double *A)
 
bool DoubleFullMatrixSolveMatrixEquation (const tLVectorIndex &nARows, const tLVectorIndex &nACols, const tLVectorIndex &nAMaxRows, const tLVectorIndex &nAMaxCols, const tLVectorIndex &ldA, const tBoolean &isTransA, const lapack_int *ipiv, const double *A, const tLVectorIndex &nBRows, const tLVectorIndex &nBCols, const tLVectorIndex &ldB, double *B)
 
int DoubleFullMatrixQRFactorization (const tLVectorIndex &nRows, const tLVectorIndex &nCols, const tLVectorIndex &nMaxRows, const tLVectorIndex &nMaxCols, const tLVectorIndex &ld, double *A, const tBoolean &isOnlyQ)
 
int DoubleFullMatrixPQRFactorization (const tLVectorIndex &nRows, const tLVectorIndex &nCols, const tLVectorIndex &nMaxRows, const tLVectorIndex &nMaxCols, const tLVectorIndex &ld, double *A, const tBoolean &isOnlyQ)
 

Function Documentation

int DoubleFullMatrixEigenValues ( const tLVectorIndex nRows,
const tLVectorIndex nCols,
const tLVectorIndex ld,
double *  A,
const tLVectorIndex nU,
double *  U,
double *  C 
)

References dgeev(), F77NAME, lapack_char, lapack_int, lapack_real, and null.

Referenced by LAP_DoubleFullGeneralMatrix::computeEigenValueDecomposition(), and LAP_DoubleFullGeneralMatrix::computeEigenValues().

Here is the call graph for this function:

Here is the caller graph for this function:

bool DoubleFullMatrixInverse ( const tLVectorIndex nRows,
const tLVectorIndex nCols,
const tLVectorIndex ld,
const lapack_int ipiv,
double *  A 
)

References dgetri(), F77NAME, getBlockSizeEnvironment(), lapack_int, and lapack_real.

Referenced by LAP_DoubleFullGeneralMatrix::inverse().

Here is the call graph for this function:

Here is the caller graph for this function:

int DoubleFullMatrixLUFactorization ( const tLVectorIndex nRows,
const tLVectorIndex nCols,
const tLVectorIndex ld,
double *  A,
lapack_int LU_pivots 
)

References dgetrf(), F77NAME, and lapack_int.

Here is the call graph for this function:

void DoubleFullMatrixMatrixProduct ( const tLVectorIndex nARows,
const tLVectorIndex nACols,
const tLVectorIndex ldA,
const double *  A,
const tLVectorIndex nBRows,
const tLVectorIndex nBCols,
const tLVectorIndex ldB,
const double *  B,
const tBoolean isTrA,
const tBoolean isTrB,
const lapack_real alpha,
const lapack_real beta,
const tLVectorIndex nCRows,
const tLVectorIndex nCCols,
const tLVectorIndex ldC,
double *  C 
)

References dgemm(), F77NAME, lapack_char, lapack_int, and tLVectorIndex.

Referenced by LAP_DoubleFullGeneralMatrix::matrixProduct().

Here is the call graph for this function:

Here is the caller graph for this function:

int DoubleFullMatrixPQRFactorization ( const tLVectorIndex nRows,
const tLVectorIndex nCols,
const tLVectorIndex nMaxRows,
const tLVectorIndex nMaxCols,
const tLVectorIndex ld,
double *  A,
const tBoolean isOnlyQ 
)

References dgeqp3(), dorgqr(), F77NAME, getBlockSizeEnvironment(), lapack_int, and lapack_real.

Here is the call graph for this function:

int DoubleFullMatrixQRFactorization ( const tLVectorIndex nRows,
const tLVectorIndex nCols,
const tLVectorIndex nMaxRows,
const tLVectorIndex nMaxCols,
const tLVectorIndex ld,
double *  A,
const tBoolean isOnlyQ 
)

References dgeqrf(), dorgqr(), F77NAME, getBlockSizeEnvironment(), lapack_int, and lapack_real.

Here is the call graph for this function:

bool DoubleFullMatrixSolveMatrixEquation ( const tLVectorIndex nARows,
const tLVectorIndex nACols,
const tLVectorIndex nAMaxRows,
const tLVectorIndex nAMaxCols,
const tLVectorIndex ldA,
const tBoolean isTransA,
const lapack_int ipiv,
const double *  A,
const tLVectorIndex nBRows,
const tLVectorIndex nBCols,
const tLVectorIndex ldB,
double *  B 
)

References dgetrs(), F77NAME, lapack_char, and lapack_int.

Referenced by LAP_DoubleFullGeneralMatrix::solve().

Here is the call graph for this function:

Here is the caller graph for this function:

void DoubleFullMatrixVectorProduct ( const tLVectorIndex lx,
const tLVectorIncrement incx,
const double *  x,
const tBoolean isTransposed,
const tLVectorIndex nRows,
const tLVectorIndex nCols,
const tLVectorIndex ldaA,
const double *  A,
const lapack_real alpha,
const lapack_real beta,
const tLVectorIndex ly,
const tLVectorIncrement incy,
double *  y 
)

compute Y=beta.Y + alpha op(A). X where op(A)=A or op(A)=tA

Parameters
lxsize of x
incxincrement of x
x: values of x
isTransposedif true op(A)=tA otherwise op(A)=A
nRowsnumber of rows of A
nColsnumber of cols of A
ldAdistance in memory netween A[i][j] and A[i][j+1]
Avalues of A stored by column
alphaalpha value
betabeta value
lysize of y
incy:incrementof y
yvalues of y

if (op(A)=A)

  • lx=nCols
  • x array must contains lx*incx double values
  • A array must constains ldA*nCols*nRows double values
  • ly=nRows if (opA=tA)
  • lx=nRows
  • x array must contains lx*incx double values
  • A array must constains ldA*nCols*nRows double values
  • ly=nCols

References dgemv(), F77NAME, lapack_char, and lapack_int.

Referenced by LAP_DoubleFullGeneralMatrix::vectorProduct().

Here is the call graph for this function:

Here is the caller graph for this function:

void DoubleFullUpperMatrixVectorProduct ( const tBoolean isTransA,
const tLVectorIndex nRowsA,
const tLVectorIndex ldA,
const double *  A,
const tLVectorIndex nX,
const tLVectorIncrement incX,
double *  x 
)

compute X=op(A). X where op(A)=A or op(A)=tA where A is Upper matrix

Parameters
isTransAif true op(A)=tA otherwise op(A)=A
nRowsAnumber of rows of A
ldAdistance in memory netween A[i][j] and A[i][j+1]
Avalues of A stored by column
nXsize of x
incX:incrementof x
xvalues of x

Referenced by LAP_DoubleFullUpperMatrix::rankProduct(), and LAP_DoubleFullUpperMatrix::vectorProduct().

Here is the caller graph for this function: