accelerInt  v0.1
Functions
complexInverse.c File Reference

Implementation of LU factorization of complex (variable-sized) matricies. More...

#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <float.h>
#include <complex.h>
#include <string.h>
#include "solver_props.h"
#include "lapack_dfns.h"
Include dependency graph for complexInverse.c:

Go to the source code of this file.

Functions

static int getHessenbergLU (const int n, double complex *__restrict__ A, int *__restrict__ indPivot)
 Computes the LU factorization of a (n x STRIDE) Hessenberg Matrix using partial pivoting with row interchanges. More...
 
void scaleComplex (const int n, const double complex val, double complex *__restrict__ arrX)
 scaleComplex scales a vector (with increment equal to one) by a constant val. More...
 
void multiplyComplexUpperMV (const int n, double complex *x, const int lda, const double complex *A)
 Performs the matrix-vector operation \(x_v:= A*x_v\). More...
 
void complexGEMV (const int m, const int n, const int lda, const double alpha, const double complex *__restrict__ A, const double complex *arrX, double complex *arrY)
 Computes the matrix-vector operation \(alpha*A*x + y\) where alpha is a scalar, x and y are vectors and A is an m by n matrix. More...
 
int getComplexInverseHessenbergLU (const int n, double complex *__restrict__ A, const int *__restrict__ indPivot)
 getComplexInverseHessenbergLU computes the inverse of a matrix using the LU factorization computed by getHessenbergLU. More...
 
void getComplexInverseHessenberg (const int n, double complex *__restrict__ A, int *__restrict__ ipiv, int *__restrict__ info)
 getComplexInverseHessenberg computes the inverse of an upper Hessenberg matrix A using a LU factorization method More...
 

Detailed Description

Implementation of LU factorization of complex (variable-sized) matricies.

Adapted from Lapack LU factorization and inversion routines

Definition in file complexInverse.c.

Function Documentation

◆ complexGEMV()

void complexGEMV ( const int  m,
const int  n,
const int  lda,
const double  alpha,
const double complex *__restrict__  A,
const double complex *  arrX,
double complex *  arrY 
)

Computes the matrix-vector operation \(alpha*A*x + y\) where alpha is a scalar, x and y are vectors and A is an m by n matrix.

Parameters
[in]mOn entry, M specifies the number of rows of the matrix A. Must be >= 0
[out]nOn entry, N specifies the number of columns of the matrix A. Must be >= 0
[in]ldaThe stride of the matrix
See also
STRIDE
Parameters
[in]alphaThe scalar value
[in]AA is an array of dimension (lda, n). Before entry, the leading m by n part of the array A must contain the matrix of coefficients.
[in]arrXarrX is an array of dimension at least (n) Before entry, the incremented array arrX must contain the vector x.
[in,out]arrYarrY is an array of dimension at least (m). Before entry, the incremented array arrY must contain the vector y. On exit, arrY is overwritten by the updated vector y.

Note: These pointers cannot use the __restrict__ modifier, as they may overlap

Definition at line 146 of file complexInverse.c.

◆ getComplexInverseHessenberg()

void getComplexInverseHessenberg ( const int  n,
double complex *__restrict__  A,
int *__restrict__  ipiv,
int *__restrict__  info 
)

getComplexInverseHessenberg computes the inverse of an upper Hessenberg matrix A using a LU factorization method

Parameters
[in]nThe order of the matrix A. n >= 0.
[in,out]AThe array, dimension (STRIDE, n)
See also
STRIDE
Parameters
[out]ipivipiv is an array of dimension (n). The pivot indices from getHessenbergLU; for 0<=i<=n-1, row i of the matrix was interchanged with row indPiv[i].
[out]infoIf not zero, an error occured during factorization

Definition at line 238 of file complexInverse.c.

◆ getComplexInverseHessenbergLU()

int getComplexInverseHessenbergLU ( const int  n,
double complex *__restrict__  A,
const int *__restrict__  indPivot 
)

getComplexInverseHessenbergLU computes the inverse of a matrix using the LU factorization computed by getHessenbergLU.

This method inverts U and then computes inv(A) by solving the system inv(A)*L = inv(U) for inv(A).

Parameters
[in]nThe order of the matrix A. n >= 0.
[in,out]AThe array, dimension (STRIDE, n)
See also
STRIDE
Parameters
[in]indPivotindPivot is an array of dimension (n). The pivot indices from getHessenbergLU; for 0<=i<=n-1, row i of the matrix was interchanged with row indPiv[i].

Definition at line 178 of file complexInverse.c.

◆ getHessenbergLU()

static int getHessenbergLU ( const int  n,
double complex *__restrict__  A,
int *__restrict__  indPivot 
)
inlinestatic

Computes the LU factorization of a (n x STRIDE) Hessenberg Matrix using partial pivoting with row interchanges.

See also
STRIDE
Parameters
[in]nThe matrix size
[in,out]AThe matrix to factorize (nxn) with stride defined in solver_props.h
See also
STRIDE
Parameters
[out]indPivotindPivot is an array of dimension (n). The pivot indices from getHessenbergLU; for 0<=i<=n-1, row i of the matrix was interchanged with row indPiv[i].

The factorization has the form: \(A = P * L * U\) where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n). For full reference see: G. W. Stewart, Matrix Algorithms: Volume 1: Basic Decompositions, SIAM, Philadelphia, 1998. doi:10.1137/1.9781611971408.

Definition at line 35 of file complexInverse.c.

◆ multiplyComplexUpperMV()

void multiplyComplexUpperMV ( const int  n,
double complex *  x,
const int  lda,
const double complex *  A 
)

Performs the matrix-vector operation \(x_v:= A*x_v\).

Parameters
[in]nOn entry, n specifies the order of the matrix A. n must be at least zero.
[out]xx is an array of dimension at least (n). Before entry, the incremented array X must contain the n element vector \(x_v\). On exit, X is overwritten with the transformed vector \(x_v\).
[in]ldaThe stride of the matrix
See also
STRIDE
Parameters
[in]AA is an array of dimension (lda, n). Before entry the leading n by n upper triangular part of the array A must contain the upper triangular matrix and the strictly lower triangular part of A is not referenced.

Note: These pointers can't use the __restrict__ attribute, as they may overlap

Definition at line 108 of file complexInverse.c.

◆ scaleComplex()

void scaleComplex ( const int  n,
const double complex  val,
double complex *__restrict__  arrX 
)

scaleComplex scales a vector (with increment equal to one) by a constant val.

Parameters
[in]nThe vector size
[out]valThe value to scale by
[out]arrXThe vector to scale

Definition at line 83 of file complexInverse.c.