1 //=================================================================================================
2 /*!
3 //  \file blaze/math/lapack/clapack/unmlq.h
4 //  \brief Header file for the CLAPACK unmlq 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_UNMLQ_H_
36 #define _BLAZE_MATH_LAPACK_CLAPACK_UNMLQ_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 #include <blaze/util/Types.h>
47 
48 
49 //=================================================================================================
50 //
51 //  LAPACK FORWARD DECLARATIONS
52 //
53 //=================================================================================================
54 
55 //*************************************************************************************************
56 /*! \cond BLAZE_INTERNAL */
57 #if !defined(INTEL_MKL_VERSION)
58 extern "C" {
59 
60 void cunmlq_( char* side, char* trans, blaze::blas_int_t* m, blaze::blas_int_t* n,
61               blaze::blas_int_t* k, float* A, blaze::blas_int_t* lda, float* tau, float* C,
62               blaze::blas_int_t* ldc, float* work, blaze::blas_int_t* lwork, blaze::blas_int_t* info,
63               blaze::fortran_charlen_t nside, blaze::fortran_charlen_t ntrans );
64 void zunmlq_( char* side, char* trans, blaze::blas_int_t* m, blaze::blas_int_t* n,
65               blaze::blas_int_t* k, double* A, blaze::blas_int_t* lda, double* tau, double* C,
66               blaze::blas_int_t* ldc, double* work, blaze::blas_int_t* lwork, blaze::blas_int_t* info,
67               blaze::fortran_charlen_t nside, blaze::fortran_charlen_t ntrans );
68 
69 }
70 #endif
71 /*! \endcond */
72 //*************************************************************************************************
73 
74 
75 
76 
77 namespace blaze {
78 
79 //=================================================================================================
80 //
81 //  LAPACK FUNCTIONS TO MULTIPLY Q FROM A LQ DECOMPOSITION WITH A MATRIX (UNMLQ)
82 //
83 //=================================================================================================
84 
85 //*************************************************************************************************
86 /*!\name LAPACK functions to multiply Q from a LQ decomposition with a matrix (unmlq) */
87 //@{
88 void unmlq( char side, char trans, blas_int_t m, blas_int_t n,
89             blas_int_t k, const complex<float>* A, blas_int_t lda,
90             const complex<float>* tau, complex<float>* C, blas_int_t ldc,
91             complex<float>* work, blas_int_t lwork, blas_int_t* info );
92 
93 void unmlq( char side, char trans, blas_int_t m, blas_int_t n,
94             blas_int_t k, const complex<double>* A, blas_int_t lda,
95             const complex<double>* tau, complex<double>* C, blas_int_t ldc,
96             complex<double>* work, blas_int_t lwork, blas_int_t* info );
97 //@}
98 //*************************************************************************************************
99 
100 
101 //*************************************************************************************************
102 /*!\brief LAPACK kernel for the multiplication of the single precision complex Q from a LQ
103 //        decomposition with another matrix.
104 // \ingroup lapack_decomposition
105 //
106 // \param side \c 'L' to apply \f$ Q \f$ or \f$ Q^H \f$ from the left, \c 'R' to apply from the right.
107 // \param trans \c 'N' for \f$ Q \f$, \c 'C' for \f$ Q^H \f$.
108 // \param m The number of rows of the matrix \a C \f$[0..\infty)\f$.
109 // \param n The number of columns of the matrix \a C \f$[0..\infty)\f$.
110 // \param k The number of elementary reflectors, whose product defines the matrix \a Q.
111 // \param A Pointer to the first element of the LQ decomposed single precision complex column-major matrix.
112 // \param lda The total number of elements between two columns of the matrix \a A \f$[0..\infty)\f$.
113 // \param tau Array for the scalar factors of the elementary reflectors; size >= min( \a m, \a n ).
114 // \param C Pointer to the first element of the single precision complex column-major matrix multiplicator.
115 // \param ldc The total number of elements between two columns of the matrix \a C \f$[0..\infty)\f$.
116 // \param work Auxiliary array; size >= max( 1, \a lwork ).
117 // \param lwork The dimension of the array \a work.
118 // \param info Return code of the function call.
119 // \return void
120 //
121 // This function multiplies the \a Q matrix resulting from the LQ decomposition of the sgelqf()
122 // function with the given general single precision complex \a m-by-\a n matrix \a C. Depending
123 // on the settings of \a side and \a trans it overwrites \a C with
124 
125    \code
126                 | side = 'L'    | side = 'R'
127    -------------|---------------|---------------
128    trans = 'N': | Q * C         | C * Q
129    trans = 'C': | ctrans(Q) * C | C * ctrans(Q)
130    \endcode
131 
132 // In case \a side is specified as \c 'L', the \a Q matrix is expected to be of size \a m-by-\a m,
133 // in case \a side is set to \c 'R', \a Q is expected to be of size \a n-by-\a n.
134 //
135 // The \a info argument provides feedback on the success of the function call:
136 //
137 //   - = 0: The decomposition finished successfully.
138 //   - < 0: The i-th argument had an illegal value.
139 //
140 // For more information on the cunmlq() function, see the LAPACK online documentation browser:
141 //
142 //        http://www.netlib.org/lapack/explore-html/
143 //
144 // \note This function can only be used if a fitting LAPACK library, which supports this function,
145 // is available and linked to the executable. Otherwise a call to this function will result in a
146 // linker error.
147 */
unmlq(char side,char trans,blas_int_t m,blas_int_t n,blas_int_t k,const complex<float> * A,blas_int_t lda,const complex<float> * tau,complex<float> * C,blas_int_t ldc,complex<float> * work,blas_int_t lwork,blas_int_t * info)148 inline void unmlq( char side, char trans, blas_int_t m, blas_int_t n,
149                    blas_int_t k, const complex<float>* A, blas_int_t lda,
150                    const complex<float>* tau, complex<float>* C, blas_int_t ldc,
151                    complex<float>* work, blas_int_t lwork, blas_int_t* info )
152 {
153    BLAZE_STATIC_ASSERT( sizeof( complex<float> ) == 2UL*sizeof( float ) );
154 
155 #if defined(INTEL_MKL_VERSION)
156    BLAZE_STATIC_ASSERT( sizeof( MKL_INT ) == sizeof( blas_int_t ) );
157    BLAZE_STATIC_ASSERT( sizeof( MKL_Complex8 ) == sizeof( complex<float> ) );
158    using ET = MKL_Complex8;
159 #else
160    using ET = float;
161 #endif
162 
163    cunmlq_( &side, &trans, &m, &n, &k,
164             const_cast<ET*>( reinterpret_cast<const ET*>( A ) ), &lda,
165             const_cast<ET*>( reinterpret_cast<const ET*>( tau ) ),
166             reinterpret_cast<ET*>( C ), &ldc, reinterpret_cast<ET*>( work ),
167             &lwork, info
168 #if !defined(INTEL_MKL_VERSION)
169           , blaze::fortran_charlen_t(1), blaze::fortran_charlen_t(1)
170 #endif
171           );
172 }
173 //*************************************************************************************************
174 
175 
176 //*************************************************************************************************
177 /*!\brief LAPACK kernel for the multiplication of the double precision complex Q from a LQ
178 //        decomposition with another matrix.
179 // \ingroup lapack_decomposition
180 //
181 // \param side \c 'L' to apply \f$ Q \f$ or \f$ Q^H \f$ from the left, \c 'R' to apply from the right.
182 // \param trans \c 'N' for \f$ Q \f$, \c 'C' for \f$ Q^H \f$.
183 // \param m The number of rows of the matrix \a C \f$[0..\infty)\f$.
184 // \param n The number of columns of the matrix \a C \f$[0..\infty)\f$.
185 // \param k The number of elementary reflectors, whose product defines the matrix \a Q.
186 // \param A Pointer to the first element of the LQ decomposed double precision complex column-major matrix.
187 // \param lda The total number of elements between two columns of the matrix \a A \f$[0..\infty)\f$.
188 // \param tau Array for the scalar factors of the elementary reflectors; size >= min( \a m, \a n ).
189 // \param C Pointer to the first element of the double precision complex column-major matrix multiplicator.
190 // \param ldc The total number of elements between two columns of the matrix \a C \f$[0..\infty)\f$.
191 // \param work Auxiliary array; size >= max( 1, \a lwork ).
192 // \param lwork The dimension of the array \a work.
193 // \param info Return code of the function call.
194 // \return void
195 //
196 // This function multiplies the \a Q matrix resulting from the LQ decomposition of the sgelqf()
197 // function with the given general double precision complex \a m-by-\a n matrix \a C. Depending
198 // on the settings of \a side and \a trans it overwrites \a C with
199 
200    \code
201                 | side = 'L'    | side = 'R'
202    -------------|---------------|---------------
203    trans = 'N': | Q * C         | C * Q
204    trans = 'C': | ctrans(Q) * C | C * ctrans(Q)
205    \endcode
206 
207 // In case \a side is specified as \c 'L', the \a Q matrix is expected to be of size \a m-by-\a m,
208 // in case \a side is set to \c 'R', \a Q is expected to be of size \a n-by-\a n.
209 //
210 // The \a info argument provides feedback on the success of the function call:
211 //
212 //   - = 0: The decomposition finished successfully.
213 //   - < 0: The i-th argument had an illegal value.
214 //
215 // For more information on the zunmlq() function, see the LAPACK online documentation browser:
216 //
217 //        http://www.netlib.org/lapack/explore-html/
218 //
219 // \note This function can only be used if a fitting LAPACK library, which supports this function,
220 // is available and linked to the executable. Otherwise a call to this function will result in a
221 // linker error.
222 */
unmlq(char side,char trans,blas_int_t m,blas_int_t n,blas_int_t k,const complex<double> * A,blas_int_t lda,const complex<double> * tau,complex<double> * C,blas_int_t ldc,complex<double> * work,blas_int_t lwork,blas_int_t * info)223 inline void unmlq( char side, char trans, blas_int_t m, blas_int_t n,
224                    blas_int_t k, const complex<double>* A, blas_int_t lda,
225                    const complex<double>* tau, complex<double>* C, blas_int_t ldc,
226                    complex<double>* work, blas_int_t lwork, blas_int_t* info )
227 {
228    BLAZE_STATIC_ASSERT( sizeof( complex<double> ) == 2UL*sizeof( double ) );
229 
230 #if defined(INTEL_MKL_VERSION)
231    BLAZE_STATIC_ASSERT( sizeof( MKL_INT ) == sizeof( blas_int_t ) );
232    BLAZE_STATIC_ASSERT( sizeof( MKL_Complex16 ) == sizeof( complex<double> ) );
233    using ET = MKL_Complex16;
234 #else
235    using ET = double;
236 #endif
237 
238    zunmlq_( &side, &trans, &m, &n, &k,
239             const_cast<ET*>( reinterpret_cast<const ET*>( A ) ), &lda,
240             const_cast<ET*>( reinterpret_cast<const ET*>( tau ) ),
241             reinterpret_cast<ET*>( C ), &ldc, reinterpret_cast<ET*>( work ),
242             &lwork, info
243 #if !defined(INTEL_MKL_VERSION)
244           , blaze::fortran_charlen_t(1), blaze::fortran_charlen_t(1)
245 #endif
246           );
247 }
248 //*************************************************************************************************
249 
250 } // namespace blaze
251 
252 #endif
253