1 //=================================================================================================
2 /*!
3 // \file blaze/math/lapack/clapack/ungrq.h
4 // \brief Header file for the CLAPACK ungrq wrapper functions
5 //
6 // Copyright (C) 2012-2020 Klaus Iglberger - All Rights Reserved
7 //
8 // This file is part of the Blaze library. You can redistribute it and/or modify it under
9 // the terms of the New (Revised) BSD License. Redistribution and use in source and binary
10 // forms, with or without modification, are permitted provided that the following conditions
11 // are met:
12 //
13 // 1. Redistributions of source code must retain the above copyright notice, this list of
14 // conditions and the following disclaimer.
15 // 2. Redistributions in binary form must reproduce the above copyright notice, this list
16 // of conditions and the following disclaimer in the documentation and/or other materials
17 // provided with the distribution.
18 // 3. Neither the names of the Blaze development group nor the names of its contributors
19 // may be used to endorse or promote products derived from this software without specific
20 // prior written permission.
21 //
22 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY
23 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES
24 // OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT
25 // SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
26 // INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
27 // TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR
28 // BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
29 // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 // ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH
31 // DAMAGE.
32 */
33 //=================================================================================================
34
35 #ifndef _BLAZE_MATH_LAPACK_CLAPACK_UNGRQ_H_
36 #define _BLAZE_MATH_LAPACK_CLAPACK_UNGRQ_H_
37
38
39 //*************************************************************************************************
40 // Includes
41 //*************************************************************************************************
42
43 #include <blaze/math/blas/Types.h>
44 #include <blaze/util/Complex.h>
45 #include <blaze/util/StaticAssert.h>
46
47
48 //=================================================================================================
49 //
50 // LAPACK FORWARD DECLARATIONS
51 //
52 //=================================================================================================
53
54 //*************************************************************************************************
55 /*! \cond BLAZE_INTERNAL */
56 #if !defined(INTEL_MKL_VERSION)
57 extern "C" {
58
59 void cungrq_( blaze::blas_int_t* m, blaze::blas_int_t* n, blaze::blas_int_t* k, float* A,
60 blaze::blas_int_t* lda, float* tau, float* work, blaze::blas_int_t* lwork,
61 blaze::blas_int_t* info );
62 void zungrq_( blaze::blas_int_t* m, blaze::blas_int_t* n, blaze::blas_int_t* k, double* A,
63 blaze::blas_int_t* lda, double* tau, double* work, blaze::blas_int_t* lwork,
64 blaze::blas_int_t* info );
65
66 }
67 #endif
68 /*! \endcond */
69 //*************************************************************************************************
70
71
72
73
74 namespace blaze {
75
76 //=================================================================================================
77 //
78 // LAPACK FUNCTIONS TO RECONSTRUCT Q FROM A RQ DECOMPOSITION (UNGRQ)
79 //
80 //=================================================================================================
81
82 //*************************************************************************************************
83 /*!\name LAPACK functions to reconstruct Q from a RQ decomposition (ungrq) */
84 //@{
85 void ungrq( blas_int_t m, blas_int_t n, blas_int_t k, complex<float>* A,
86 blas_int_t lda, const complex<float>* tau, complex<float>* work,
87 blas_int_t lwork, blas_int_t* info );
88
89 void ungrq( blas_int_t m, blas_int_t n, blas_int_t k, complex<double>* A,
90 blas_int_t lda, const complex<double>* tau, complex<double>* work,
91 blas_int_t lwork, blas_int_t* info );
92 //@}
93 //*************************************************************************************************
94
95
96 //*************************************************************************************************
97 /*!\brief LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition.
98 // \ingroup lapack_decomposition
99 //
100 // \param m The number of rows of the given matrix \f$[0..n)\f$.
101 // \param n The number of columns of the given matrix \f$[0..\infty)\f$.
102 // \param k The number of elementary reflectors, whose product defines the matrix \f$[0..m)\f$.
103 // \param A Pointer to the first element of the single precision column-major matrix.
104 // \param lda The total number of elements between two columns of the matrix \f$[0..\infty)\f$.
105 // \param tau Array for the scalar factors of the elementary reflectors; size >= min( \a m, \a n ).
106 // \param work Auxiliary array; size >= max( 1, \a lwork ).
107 // \param lwork The dimension of the array \a work; size >= max( 1, \a n ).
108 // \param info Return code of the function call.
109 // \return void
110 //
111 // This function generates all or part of the orthogonal matrix Q from a RQ decomposition based on
112 // the LAPACK cungrq() function for single precision column-major matrices that have already been
113 // factorized by the sgerqf() function. The \a info argument provides feedback on the success of
114 // the function call:
115 //
116 // - = 0: The decomposition finished successfully.
117 // - < 0: The i-th argument had an illegal value.
118 //
119 // For more information on the cungrq() function, see the LAPACK online documentation browser:
120 //
121 // http://www.netlib.org/lapack/explore-html/
122 //
123 // \note This function can only be used if a fitting LAPACK library, which supports this function,
124 // is available and linked to the executable. Otherwise a call to this function will result in a
125 // linker error.
126 */
ungrq(blas_int_t m,blas_int_t n,blas_int_t k,complex<float> * A,blas_int_t lda,const complex<float> * tau,complex<float> * work,blas_int_t lwork,blas_int_t * info)127 inline void ungrq( blas_int_t m, blas_int_t n, blas_int_t k, complex<float>* A,
128 blas_int_t lda, const complex<float>* tau, complex<float>* work,
129 blas_int_t lwork, blas_int_t* info )
130 {
131 BLAZE_STATIC_ASSERT( sizeof( complex<float> ) == 2UL*sizeof( float ) );
132
133 #if defined(INTEL_MKL_VERSION)
134 BLAZE_STATIC_ASSERT( sizeof( MKL_INT ) == sizeof( blas_int_t ) );
135 BLAZE_STATIC_ASSERT( sizeof( MKL_Complex8 ) == sizeof( complex<float> ) );
136 using ET = MKL_Complex8;
137 #else
138 using ET = float;
139 #endif
140
141 cungrq_( &m, &n, &k, reinterpret_cast<ET*>( A ), &lda,
142 const_cast<ET*>( reinterpret_cast<const ET*>( tau ) ),
143 reinterpret_cast<ET*>( work ), &lwork, info );
144 }
145 //*************************************************************************************************
146
147
148 //*************************************************************************************************
149 /*!\brief LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition.
150 // \ingroup lapack_decomposition
151 //
152 // \param m The number of rows of the given matrix \f$[0..n)\f$.
153 // \param n The number of columns of the given matrix \f$[0..\infty)\f$.
154 // \param k The number of elementary reflectors, whose product defines the matrix \f$[0..m)\f$.
155 // \param A Pointer to the first element of the double precision column-major matrix.
156 // \param lda The total number of elements between two columns of the matrix \f$[0..\infty)\f$.
157 // \param tau Array for the scalar factors of the elementary reflectors; size >= min( \a m, \a n ).
158 // \param work Auxiliary array; size >= max( 1, \a lwork ).
159 // \param lwork The dimension of the array \a work; size >= max( 1, \a n ).
160 // \param info Return code of the function call.
161 // \return void
162 //
163 // This function generates all or part of the orthogonal matrix Q from a RQ decomposition based on
164 // the LAPACK zungrq() function for double precision column-major matrices that have already been
165 // factorized by the dgerqf() function. The \a info argument provides feedback on the success of
166 // the function call:
167 //
168 // - = 0: The decomposition finished successfully.
169 // - < 0: The i-th argument had an illegal value.
170 //
171 // For more information on the zungrq() function, see the LAPACK online documentation browser:
172 //
173 // http://www.netlib.org/lapack/explore-html/
174 //
175 // \note This function can only be used if a fitting LAPACK library, which supports this function,
176 // is available and linked to the executable. Otherwise a call to this function will result in a
177 // linker error.
178 */
ungrq(blas_int_t m,blas_int_t n,blas_int_t k,complex<double> * A,blas_int_t lda,const complex<double> * tau,complex<double> * work,blas_int_t lwork,blas_int_t * info)179 inline void ungrq( blas_int_t m, blas_int_t n, blas_int_t k, complex<double>* A,
180 blas_int_t lda, const complex<double>* tau, complex<double>* work,
181 blas_int_t lwork, blas_int_t* info )
182 {
183 BLAZE_STATIC_ASSERT( sizeof( complex<double> ) == 2UL*sizeof( double ) );
184
185 #if defined(INTEL_MKL_VERSION)
186 BLAZE_STATIC_ASSERT( sizeof( MKL_INT ) == sizeof( blas_int_t ) );
187 BLAZE_STATIC_ASSERT( sizeof( MKL_Complex16 ) == sizeof( complex<double> ) );
188 using ET = MKL_Complex16;
189 #else
190 using ET = double;
191 #endif
192
193 zungrq_( &m, &n, &k, reinterpret_cast<ET*>( A ), &lda,
194 const_cast<ET*>( reinterpret_cast<const ET*>( tau ) ),
195 reinterpret_cast<ET*>( work ), &lwork, info );
196 }
197 //*************************************************************************************************
198
199 } // namespace blaze
200
201 #endif
202