1 //=================================================================================================
2 /*!
3 //  \file blaze/math/lapack/clapack/potri.h
4 //  \brief Header file for the CLAPACK potri 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; ORx
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_POTRI_H_
36 #define _BLAZE_MATH_LAPACK_CLAPACK_POTRI_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) && !defined(BLAS_H)
58 extern "C" {
59 
60 void spotri_( char* uplo, blaze::blas_int_t* n, float* A, blaze::blas_int_t* lda,
61               blaze::blas_int_t* info, blaze::fortran_charlen_t nuplo );
62 void dpotri_( char* uplo, blaze::blas_int_t* n, double* A, blaze::blas_int_t* lda,
63               blaze::blas_int_t* info, blaze::fortran_charlen_t nuplo );
64 void cpotri_( char* uplo, blaze::blas_int_t* n, float* A, blaze::blas_int_t* lda,
65               blaze::blas_int_t* info, blaze::fortran_charlen_t nuplo );
66 void zpotri_( char* uplo, blaze::blas_int_t* n, double* A, blaze::blas_int_t* lda,
67               blaze::blas_int_t* info, blaze::fortran_charlen_t nuplo );
68 
69 }
70 #endif
71 /*! \endcond */
72 //*************************************************************************************************
73 
74 
75 
76 
77 namespace blaze {
78 
79 //=================================================================================================
80 //
81 //  LAPACK LLH-BASED INVERSION FUNCTIONS (POTRI)
82 //
83 //=================================================================================================
84 
85 //*************************************************************************************************
86 /*!\name LAPACK LLH-based inversion functions (potri) */
87 //@{
88 void potri( char uplo, blas_int_t n, float* A, blas_int_t lda, blas_int_t* info );
89 
90 void potri( char uplo, blas_int_t n, double* A, blas_int_t lda, blas_int_t* info );
91 
92 void potri( char uplo, blas_int_t n, complex<float>* A, blas_int_t lda, blas_int_t* info );
93 
94 void potri( char uplo, blas_int_t n, complex<double>* A, blas_int_t lda, blas_int_t* info );
95 //@}
96 //*************************************************************************************************
97 
98 
99 //*************************************************************************************************
100 /*!\brief LAPACK kernel for the inversion of the given dense positive definite single precision
101 //        column-major square matrix.
102 // \ingroup lapack_inversion
103 //
104 // \param uplo \c 'L' to use the lower part of the matrix, \c 'U' to use the upper part.
105 // \param n The number of rows/columns of the matrix \f$[0..\infty)\f$.
106 // \param A Pointer to the first element of the single precision column-major matrix.
107 // \param lda The total number of elements between two columns of the matrix \f$[0..\infty)\f$.
108 // \param info Return code of the function call.
109 // \return void
110 //
111 // This function performs the dense matrix inversion based on the LAPACK spotri() function for
112 // positive-definite single precision column-major matrices that have already been factorized
113 // by the spotrf() function. The resulting symmetric inverse of \a A is stored either in the
114 // lower part of \a A (\a uplo = \c 'L') or in the upper part (\a uplo = \c 'U').
115 //
116 // The \a info argument provides feedback on the success of the function call:
117 //
118 //   - = 0: The inversion finished successfully.
119 //   - < 0: If \a info = -i, the i-th argument had an illegal value.
120 //   - > 0: If \a info = i, element (i,i) of U or L is zero and the inverse could not be computed.
121 //
122 // For more information on the spotri() function, see the LAPACK online documentation browser:
123 //
124 //        http://www.netlib.org/lapack/explore-html/
125 //
126 // \note This function can only be used if a fitting LAPACK library, which supports this function,
127 // is available and linked to the executable. Otherwise a call to this function will result in a
128 // linker error.
129 */
potri(char uplo,blas_int_t n,float * A,blas_int_t lda,blas_int_t * info)130 inline void potri( char uplo, blas_int_t n, float* A, blas_int_t lda, blas_int_t* info )
131 {
132 #if defined(INTEL_MKL_VERSION)
133    BLAZE_STATIC_ASSERT( sizeof( MKL_INT ) == sizeof( blas_int_t ) );
134 #endif
135 
136    spotri_( &uplo, &n, A, &lda, info
137 #if !defined(INTEL_MKL_VERSION) && !defined(BLAS_H)
138           , blaze::fortran_charlen_t(1)
139 #endif
140           );
141 }
142 //*************************************************************************************************
143 
144 
145 //*************************************************************************************************
146 /*!\brief LAPACK kernel for the inversion of the given dense positive definite double precision
147 //        column-major square matrix.
148 // \ingroup lapack_inversion
149 //
150 // \param uplo \c 'L' to use the lower part of the matrix, \c 'U' to use the upper part.
151 // \param n The number of rows/columns of the matrix \f$[0..\infty)\f$.
152 // \param A Pointer to the first element of the double precision column-major matrix.
153 // \param lda The total number of elements between two columns of the matrix \f$[0..\infty)\f$.
154 // \param info Return code of the function call.
155 // \return void
156 //
157 // This function performs the dense matrix inversion based on the LAPACK dpotri() function for
158 // positive-definite double precision column-major matrices that have already been factorized
159 // by the dpotrf() function. The resulting symmetric inverse of \a A is stored either in the
160 // lower part of \a A (\a uplo = \c 'L') or in the upper part (\a uplo = \c 'U').
161 //
162 // The \a info argument provides feedback on the success of the function call:
163 //
164 //   - = 0: The inversion finished successfully.
165 //   - < 0: If \a info = -i, the i-th argument had an illegal value.
166 //   - > 0: If \a info = i, element (i,i) of U or L is zero and the inverse could not be computed.
167 //
168 // For more information on the spotri() function, see the LAPACK online documentation browser:
169 //
170 //        http://www.netlib.org/lapack/explore-html/
171 //
172 // \note This function can only be used if a fitting LAPACK library, which supports this function,
173 // is available and linked to the executable. Otherwise a call to this function will result in a
174 // linker error.
175 */
potri(char uplo,blas_int_t n,double * A,blas_int_t lda,blas_int_t * info)176 inline void potri( char uplo, blas_int_t n, double* A, blas_int_t lda, blas_int_t* info )
177 {
178 #if defined(INTEL_MKL_VERSION)
179    BLAZE_STATIC_ASSERT( sizeof( MKL_INT ) == sizeof( blas_int_t ) );
180 #endif
181 
182    dpotri_( &uplo, &n, A, &lda, info
183 #if !defined(INTEL_MKL_VERSION) && !defined(BLAS_H)
184           , blaze::fortran_charlen_t(1)
185 #endif
186           );
187 }
188 //*************************************************************************************************
189 
190 
191 //*************************************************************************************************
192 /*!\brief LAPACK kernel for the inversion of the given dense positive definite single precision
193 //        complex column-major square matrix.
194 // \ingroup lapack_inversion
195 //
196 // \param uplo \c 'L' to use the lower part of the matrix, \c 'U' to use the upper part.
197 // \param n The number of rows/columns of the matrix \f$[0..\infty)\f$.
198 // \param A Pointer to the first element of the single precision complex column-major matrix.
199 // \param lda The total number of elements between two columns of the matrix \f$[0..\infty)\f$.
200 // \param info Return code of the function call.
201 // \return void
202 //
203 // This function performs the dense matrix inversion based on the LAPACK cpotri() function for
204 // positive-definite single precision complex column-major matrices that have already been
205 // factorized by the cpotrf() function. The resulting symmetric inverse of \a A is stored either
206 // in the lower part of \a A (\a uplo = \c 'L') or in the upper part (\a uplo = \c 'U').
207 //
208 // The \a info argument provides feedback on the success of the function call:
209 //
210 //   - = 0: The inversion finished successfully.
211 //   - < 0: If \a info = -i, the i-th argument had an illegal value.
212 //   - > 0: If \a info = i, element (i,i) of U or L is zero and the inverse could not be computed.
213 //
214 // For more information on the cpotri() function, see the LAPACK online documentation browser:
215 //
216 //        http://www.netlib.org/lapack/explore-html/
217 //
218 // \note This function can only be used if a fitting LAPACK library, which supports this function,
219 // is available and linked to the executable. Otherwise a call to this function will result in a
220 // linker error.
221 */
potri(char uplo,blas_int_t n,complex<float> * A,blas_int_t lda,blas_int_t * info)222 inline void potri( char uplo, blas_int_t n, complex<float>* A, blas_int_t lda, blas_int_t* info )
223 {
224    BLAZE_STATIC_ASSERT( sizeof( complex<float> ) == 2UL*sizeof( float ) );
225 
226 #if defined(INTEL_MKL_VERSION)
227    BLAZE_STATIC_ASSERT( sizeof( MKL_INT ) == sizeof( blas_int_t ) );
228    BLAZE_STATIC_ASSERT( sizeof( MKL_Complex8 ) == sizeof( complex<float> ) );
229    using ET = MKL_Complex8;
230 #else
231    using ET = float;
232 #endif
233 
234    cpotri_( &uplo, &n, reinterpret_cast<ET*>( A ), &lda, info
235 #if !defined(INTEL_MKL_VERSION) && !defined(BLAS_H)
236           , blaze::fortran_charlen_t(1)
237 #endif
238           );
239 }
240 //*************************************************************************************************
241 
242 
243 //*************************************************************************************************
244 /*!\brief LAPACK kernel for the inversion of the given dense positive definite double precision
245 //        complex column-major square matrix.
246 // \ingroup lapack_inversion
247 //
248 // \param uplo \c 'L' to use the lower part of the matrix, \c 'U' to use the upper part.
249 // \param n The number of rows/columns of the matrix \f$[0..\infty)\f$.
250 // \param A Pointer to the first element of the double precision complex column-major matrix.
251 // \param lda The total number of elements between two columns of the matrix \f$[0..\infty)\f$.
252 // \param info Return code of the function call.
253 // \return void
254 //
255 // This function performs the dense matrix inversion based on the LAPACK zpotri() function for
256 // positive-definite double precision complex column-major matrices that have already been
257 // factorized by the zpotrf() function. The resulting symmetric inverse of \a A is stored either
258 // in the lower part of \a A (\a uplo = \c 'L') or in the upper part (\a uplo = \c 'U').
259 //
260 // The \a info argument provides feedback on the success of the function call:
261 //
262 //   - = 0: The inversion finished successfully.
263 //   - < 0: If \a info = -i, the i-th argument had an illegal value.
264 //   - > 0: If \a info = i, element (i,i) of U or L is zero and the inverse could not be computed.
265 //
266 // For more information on the zpotri() function, see the LAPACK online documentation browser:
267 //
268 //        http://www.netlib.org/lapack/explore-html/
269 //
270 // \note This function can only be used if a fitting LAPACK library, which supports this function,
271 // is available and linked to the executable. Otherwise a call to this function will result in a
272 // linker error.
273 */
potri(char uplo,blas_int_t n,complex<double> * A,blas_int_t lda,blas_int_t * info)274 inline void potri( char uplo, blas_int_t n, complex<double>* A, blas_int_t lda, blas_int_t* info )
275 {
276    BLAZE_STATIC_ASSERT( sizeof( complex<double> ) == 2UL*sizeof( double ) );
277 
278 #if defined(INTEL_MKL_VERSION)
279    BLAZE_STATIC_ASSERT( sizeof( MKL_INT ) == sizeof( blas_int_t ) );
280    BLAZE_STATIC_ASSERT( sizeof( MKL_Complex16 ) == sizeof( complex<double> ) );
281    using ET = MKL_Complex16;
282 #else
283    using ET = double;
284 #endif
285 
286    zpotri_( &uplo, &n, reinterpret_cast<ET*>( A ), &lda, info
287 #if !defined(INTEL_MKL_VERSION) && !defined(BLAS_H)
288           , blaze::fortran_charlen_t(1)
289 #endif
290           );
291 }
292 //*************************************************************************************************
293 
294 } // namespace blaze
295 
296 #endif
297