getrs.h
Go to the documentation of this file.
1 //=================================================================================================
33 //=================================================================================================
34 
35 #ifndef _BLAZE_MATH_LAPACK_CLAPACK_GETRS_H_
36 #define _BLAZE_MATH_LAPACK_CLAPACK_GETRS_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 sgetrs_( char* trans, int* n, int* nrhs, float* A, int* lda, int* ipiv, float* B, int* ldb, int* info );
58 void dgetrs_( char* trans, int* n, int* nrhs, double* A, int* lda, int* ipiv, double* B, int* ldb, int* info );
59 void cgetrs_( char* trans, int* n, int* nrhs, float* A, int* lda, int* ipiv, float* B, int* ldb, int* info );
60 void zgetrs_( char* trans, int* n, int* nrhs, double* A, int* lda, int* ipiv, double* B, int* ldb, int* info );
61 
62 }
64 //*************************************************************************************************
65 
66 
67 
68 
69 namespace blaze {
70 
71 //=================================================================================================
72 //
73 // LAPACK LU-BASED SUBSTITUTION FUNCTIONS (GETRS)
74 //
75 //=================================================================================================
76 
77 //*************************************************************************************************
80 inline void getrs( char trans, int n, int nrhs, const float* A, int lda, const int* ipiv,
81  float* B, int ldb, int* info );
82 
83 inline void getrs( char trans, int n, int nrhs, const double* A, int lda, const int* ipiv,
84  double* B, int ldb, int* info );
85 
86 inline void getrs( char trans, int n, int nrhs, const complex<float>* A, int lda,
87  const int* ipiv, complex<float>* B, int ldb, int* info );
88 
89 inline void getrs( char trans, int n, int nrhs, const complex<double>* A, int lda,
90  const int* ipiv, complex<double>* B, int ldb, int* info );
92 //*************************************************************************************************
93 
94 
95 //*************************************************************************************************
134 inline void getrs( char trans, int n, int nrhs, const float* A, int lda,
135  const int* ipiv, float* B, int ldb, int* info )
136 {
137  sgetrs_( &trans, &n, &nrhs, const_cast<float*>( A ), &lda,
138  const_cast<int*>( ipiv ), B, &ldb, info );
139 }
140 //*************************************************************************************************
141 
142 
143 //*************************************************************************************************
182 inline void getrs( char trans, int n, int nrhs, const double* A, int lda,
183  const int* ipiv, double* B, int ldb, int* info )
184 {
185  dgetrs_( &trans, &n, &nrhs, const_cast<double*>( A ), &lda,
186  const_cast<int*>( ipiv ), B, &ldb, info );
187 }
188 //*************************************************************************************************
189 
190 
191 //*************************************************************************************************
230 inline void getrs( char trans, int n, int nrhs, const complex<float>* A, int lda,
231  const int* ipiv, complex<float>* B, int ldb, int* info )
232 {
233  BLAZE_STATIC_ASSERT( sizeof( complex<float> ) == 2UL*sizeof( float ) );
234 
235  cgetrs_( &trans, &n, &nrhs, const_cast<float*>( reinterpret_cast<const float*>( A ) ),
236  &lda, const_cast<int*>( ipiv ), reinterpret_cast<float*>( B ), &ldb, info );
237 }
238 //*************************************************************************************************
239 
240 
241 //*************************************************************************************************
280 inline void getrs( char trans, int n, int nrhs, const complex<double>* A, int lda,
281  const int* ipiv, complex<double>* B, int ldb, int* info )
282 {
283  BLAZE_STATIC_ASSERT( sizeof( complex<double> ) == 2UL*sizeof( double ) );
284 
285  zgetrs_( &trans, &n, &nrhs, const_cast<double*>( reinterpret_cast<const double*>( A ) ),
286  &lda, const_cast<int*>( ipiv ), reinterpret_cast<double*>( B ), &ldb, info );
287 }
288 //*************************************************************************************************
289 
290 } // namespace blaze
291 
292 #endif
Namespace of the Blaze C++ math library.
Definition: Blaze.h:57
Compile time assertion.
void getrs(char trans, int n, int nrhs, const float *A, int lda, const int *ipiv, float *B, int ldb, int *info)
LAPACK kernel for the substitution step of solving a general single precision linear system of equati...
Definition: getrs.h:134
const DMatTransExpr< MT,!SO > trans(const DenseMatrix< MT, SO > &dm)
Calculation of the transpose of the given dense matrix.
Definition: DMatTransExpr.h:733
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