getri.h
Go to the documentation of this file.
1 //=================================================================================================
33 //=================================================================================================
34 
35 #ifndef _BLAZE_MATH_LAPACK_CLAPACK_GETRI_H_
36 #define _BLAZE_MATH_LAPACK_CLAPACK_GETRI_H_
37 
38 
39 //*************************************************************************************************
40 // Includes
41 //*************************************************************************************************
42 
43 #include <blaze/util/Complex.h>
45 
46 
47 //=================================================================================================
48 //
49 // LAPACK FORWARD DECLARATIONS
50 //
51 //=================================================================================================
52 
53 //*************************************************************************************************
55 extern "C" {
56 
57 void sgetri_( int* n, float* A, int* lda, int* ipiv, float* work, int* lwork, int* info );
58 void dgetri_( int* n, double* A, int* lda, int* ipiv, double* work, int* lwork, int* info );
59 void cgetri_( int* n, float* A, int* lda, int* ipiv, float* work, int* lwork, int* info );
60 void zgetri_( int* n, double* A, int* lda, int* ipiv, double* work, int* lwork, int* info );
61 
62 }
64 //*************************************************************************************************
65 
66 
67 
68 
69 namespace blaze {
70 
71 //=================================================================================================
72 //
73 // LAPACK LU-BASED INVERSION FUNCTIONS (GETRI)
74 //
75 //=================================================================================================
76 
77 //*************************************************************************************************
80 inline void getri( int n, float* A, int lda, const int* ipiv, float* work, int lwork, int* info );
81 
82 inline void getri( int n, double* A, int lda, const int* ipiv, double* work, int lwork, int* info );
83 
84 inline void getri( int n, complex<float>* A, int lda, const int* ipiv,
85  complex<float>* work, int lwork, int* info );
86 
87 inline void getri( int n, complex<double>* A, int lda, const int* ipiv,
88  complex<double>* work, int lwork, int* info );
90 //*************************************************************************************************
91 
92 
93 //*************************************************************************************************
128 inline void getri( int n, float* A, int lda, const int* ipiv, float* work, int lwork, int* info )
129 {
130  sgetri_( &n, A, &lda, const_cast<int*>( ipiv ), work, &lwork, info );
131 }
132 //*************************************************************************************************
133 
134 
135 //*************************************************************************************************
170 inline void getri( int n, double* A, int lda, const int* ipiv, double* work, int lwork, int* info )
171 {
172  dgetri_( &n, A, &lda, const_cast<int*>( ipiv ), work, &lwork, info );
173 }
174 //*************************************************************************************************
175 
176 
177 //*************************************************************************************************
212 inline void getri( int n, complex<float>* A, int lda, const int* ipiv,
213  complex<float>* work, int lwork, int* info )
214 {
215  BLAZE_STATIC_ASSERT( sizeof( complex<float> ) == 2UL*sizeof( float ) );
216 
217  cgetri_( &n, reinterpret_cast<float*>( A ), &lda, const_cast<int*>( ipiv ),
218  reinterpret_cast<float*>( work ), &lwork, info );
219 }
220 //*************************************************************************************************
221 
222 
223 //*************************************************************************************************
258 inline void getri( int n, complex<double>* A, int lda, const int* ipiv,
259  complex<double>* work, int lwork, int* info )
260 {
261  BLAZE_STATIC_ASSERT( sizeof( complex<double> ) == 2UL*sizeof( double ) );
262 
263  zgetri_( &n, reinterpret_cast<double*>( A ), &lda, const_cast<int*>( ipiv ),
264  reinterpret_cast<double*>( work ), &lwork, info );
265 }
266 //*************************************************************************************************
267 
268 } // namespace blaze
269 
270 #endif
Log level for high-level information.
Definition: LogLevel.h:80
Namespace of the Blaze C++ math library.
Definition: Blaze.h:57
Compile time assertion.
void getri(int n, float *A, int lda, const int *ipiv, float *work, int lwork, int *info)
LAPACK kernel for the inversion of the given dense general single precision column-major square matri...
Definition: getri.h:128
Header file for the complex data type.
#define BLAZE_STATIC_ASSERT(expr)
Compile time assertion macro.In case of an invalid compile time expression, a compilation error is cr...
Definition: StaticAssert.h:112