1 //=================================================================================================
2 /*!
3 //  \file blaze/math/lapack/clapack/orgqr.h
4 //  \brief Header file for the CLAPACK orgqr 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_ORGQR_H_
36 #define _BLAZE_MATH_LAPACK_CLAPACK_ORGQR_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 sorgqr_( 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 dorgqr_( 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 QR DECOMPOSITION (ORGQR)
79 //
80 //=================================================================================================
81 
82 //*************************************************************************************************
83 /*!\name LAPACK functions to reconstruct Q from a QR decomposition (orgqr) */
84 //@{
85 void orgqr( blas_int_t m, blas_int_t n, blas_int_t k, float* A, blas_int_t lda,
86             const float* tau, float* work, blas_int_t lwork, blas_int_t* info );
87 
88 void orgqr( blas_int_t m, blas_int_t n, blas_int_t k, double* A, blas_int_t lda,
89             const double* tau, double* work, blas_int_t lwork, blas_int_t* info );
90 //@}
91 //*************************************************************************************************
92 
93 
94 //*************************************************************************************************
95 /*!\brief LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition.
96 // \ingroup lapack_decomposition
97 //
98 // \param m The number of rows of the given matrix \f$[0..\infty)\f$.
99 // \param n The number of columns of the given matrix \f$[0..m)\f$.
100 // \param k The number of elementary reflectors, whose product defines the matrix \f$[0..n)\f$.
101 // \param A Pointer to the first element of the single precision column-major matrix.
102 // \param lda The total number of elements between two columns of the matrix \f$[0..\infty)\f$.
103 // \param tau Array for the scalar factors of the elementary reflectors; size >= min( \a m, \a n ).
104 // \param work Auxiliary array; size >= max( 1, \a lwork ).
105 // \param lwork The dimension of the array \a work; size >= max( 1, \a n ).
106 // \param info Return code of the function call.
107 // \return void
108 //
109 // This function generates all or part of the orthogonal matrix Q from a QR decomposition based on
110 // the LAPACK sorgqr() function for single precision column-major matrices that have already been
111 // factorized by the sgeqrf() function. The \a info argument provides feedback on the success of
112 // the function call:
113 //
114 //   - = 0: The decomposition finished successfully.
115 //   - < 0: The i-th argument had an illegal value.
116 //
117 // For more information on the sorgqr() function, see the LAPACK online documentation browser:
118 //
119 //        http://www.netlib.org/lapack/explore-html/
120 //
121 // \note This function can only be used if a fitting LAPACK library, which supports this function,
122 // is available and linked to the executable. Otherwise a call to this function will result in a
123 // linker error.
124 */
orgqr(blas_int_t m,blas_int_t n,blas_int_t k,float * A,blas_int_t lda,const float * tau,float * work,blas_int_t lwork,blas_int_t * info)125 inline void orgqr( blas_int_t m, blas_int_t n, blas_int_t k, float* A, blas_int_t lda,
126                    const float* tau, float* work, blas_int_t lwork, blas_int_t* info )
127 {
128 #if defined(INTEL_MKL_VERSION)
129    BLAZE_STATIC_ASSERT( sizeof( MKL_INT ) == sizeof( blas_int_t ) );
130 #endif
131 
132    sorgqr_( &m, &n, &k, A, &lda, const_cast<float*>( tau ), work, &lwork, info );
133 }
134 //*************************************************************************************************
135 
136 
137 //*************************************************************************************************
138 /*!\brief LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition.
139 // \ingroup lapack_decomposition
140 //
141 // \param m The number of rows of the given matrix \f$[0..\infty)\f$.
142 // \param n The number of columns of the given matrix \f$[0..m)\f$.
143 // \param k The number of elementary reflectors, whose product defines the matrix \f$[0..n)\f$.
144 // \param A Pointer to the first element of the double precision column-major matrix.
145 // \param lda The total number of elements between two columns of the matrix \f$[0..\infty)\f$.
146 // \param tau Array for the scalar factors of the elementary reflectors; size >= min( \a m, \a n ).
147 // \param work Auxiliary array; size >= max( 1, \a lwork ).
148 // \param lwork The dimension of the array \a work; size >= max( 1, \a n ).
149 // \param info Return code of the function call.
150 // \return void
151 //
152 // This function generates all or part of the orthogonal matrix Q from a QR decomposition based on
153 // the LAPACK dorgqr() function for double precision column-major matrices that have already been
154 // factorized by the dgeqrf() function. The \a info argument provides feedback on the success of
155 // the function call:
156 //
157 //   - = 0: The decomposition finished successfully.
158 //   - < 0: The i-th argument had an illegal value.
159 //
160 // For more information on the dorgqr() function, see the LAPACK online documentation browser:
161 //
162 //        http://www.netlib.org/lapack/explore-html/
163 //
164 // \note This function can only be used if a fitting LAPACK library, which supports this function,
165 // is available and linked to the executable. Otherwise a call to this function will result in a
166 // linker error.
167 */
orgqr(blas_int_t m,blas_int_t n,blas_int_t k,double * A,blas_int_t lda,const double * tau,double * work,blas_int_t lwork,blas_int_t * info)168 inline void orgqr( blas_int_t m, blas_int_t n, blas_int_t k, double* A, blas_int_t lda,
169                    const double* tau, double* work, blas_int_t lwork, blas_int_t* info )
170 {
171 #if defined(INTEL_MKL_VERSION)
172    BLAZE_STATIC_ASSERT( sizeof( MKL_INT ) == sizeof( blas_int_t ) );
173 #endif
174 
175    dorgqr_( &m, &n, &k, A, &lda, const_cast<double*>( tau ), work, &lwork, info );
176 }
177 //*************************************************************************************************
178 
179 } // namespace blaze
180 
181 #endif
182