Teuchos - Trilinos Tools Package Version of the Day
Loading...
Searching...
No Matches
Teuchos_LAPACK.cpp
1// @HEADER
2// ***********************************************************************
3//
4// Teuchos: Common Tools Package
5// Copyright (2004) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38//
39// ***********************************************************************
40// @HEADER
41
42#include <vector>
45#include "Teuchos_LAPACK.hpp"
47#ifdef HAVE_TEUCHOSCORE_QUADMATH
48# include "Teuchos_Details_Lapack128.hpp" // impl for __float128
49#endif // HAVE_TEUCHOSCORE_QUADMATH
50#ifdef HAVE_TEUCHOS_LONG_DOUBLE
51# include "Teuchos_Details_LapackLongDouble.hpp" // impl for long double
52#endif // HAVE_TEUCHOS_LONG_DOUBLE
54
55/* for INTEL_CXML, the second arg may need to be changed to 'one'. If so
56the appropriate declaration of one will need to be added back into
57functions that include the macro:
58*/
59
60#ifdef CHAR_MACRO
61#undef CHAR_MACRO
62#endif
63#if defined (INTEL_CXML)
64#define CHAR_MACRO(char_var) &char_var, one
65#else
66#define CHAR_MACRO(char_var) &char_var
67#endif
68
69#ifdef CHARPTR_MACRO
70#undef CHARPR_MACRO
71#endif
72#if defined (INTEL_CXML)
73#define CHARPTR_MACRO(charptr_var) charptr_var, one
74#else
75#define CHARPTR_MACRO(charptr_var) charptr_var
76#endif
77
78namespace {
79
80#if defined (INTEL_CXML)
81 unsigned int one=1;
82#endif
83
84// Use a wrapper function to handle calling ILAENV(). This removes
85// duplicaiton and avoid name lookup problems with member functions called
86// ILAENV() trying to call nonmember functions called ILAENV() (which does not
87// work on Intel compiler on Windows, see Trilinos bug 5762).
88inline
89int ilaenv_wrapper(
90 const int* ispec, const char* name, const unsigned int& name_length,
91 const char* opts, const unsigned int& opts_length,
92 const int* N1, const int* N2, const int* N3, const int* N4 )
93{
94#if defined (INTEL_CXML)
95 return ILAENV_F77(ispec, name, name_length, opts, opts_length, N1, N2, N3, N4 );
96#else
97 return ILAENV_F77(ispec, name, opts, N1, N2, N3, N4, name_length, opts_length );
98#endif
99}
100
101} // namespace
102
103
104
105extern "C" {
106
107
108typedef int (*gees_nullfptr_t)(double*,double*);
109
110
111} // extern "C"
112
113
114namespace Teuchos
115{
116 // BEGIN INT, FLOAT SPECIALIZATION IMPLEMENTATION //
117
118 void LAPACK<int, float>::PTTRF(const int& n, float* d, float* e, int* info) const
119 { SPTTRF_F77(&n,d,e,info); }
120
121
122 void LAPACK<int, float>::PTTRS(const int& n, const int& nrhs, const float* d, const float* e, float* B, const int& ldb, int* info) const
123 { SPTTRS_F77(&n,&nrhs,d,e,B,&ldb,info); }
124
125
126 void LAPACK<int, float>::POTRF(const char& UPLO, const int& n, float* A, const int& lda, int* info) const
127 { SPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info); }
128
129
130 void LAPACK<int, float>::POTRS(const char& UPLO, const int& n, const int& nrhs, const float* A, const int& lda, float* B, const int& ldb, int* info) const
131 { SPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info); }
132
133
134 void LAPACK<int, float>::POTRI(const char& UPLO, const int& n, float* A, const int& lda, int* info) const
135 { SPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info); }
136
137
138 void LAPACK<int, float>::POCON(const char& UPLO, const int& n, const float* A, const int& lda, const float& anorm, float* rcond, float* WORK, int* IWORK, int* info) const
139 { SPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, IWORK, info); }
140
141
142 void LAPACK<int, float>::POSV(const char& UPLO, const int& n, const int& nrhs, float* A, const int& lda, float* B, const int& ldb, int* info) const
143 { SPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info); }
144
145
146 void LAPACK<int, float>::POEQU(const int& n, const float* A, const int& lda, float* S, float* scond, float* amax, int* info) const
147 { SPOEQU_F77(&n, A, &lda, S, scond, amax, info); }
148
149
150 void LAPACK<int, float>::PORFS(const char& UPLO, const int& n, const int& nrhs, float* A, const int& lda, const float* AF, const int& ldaf, const float* B, const int& ldb, float* X, const int& ldx, float* FERR, float* BERR, float* WORK, int* IWORK, int* info) const
151 { SPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info); }
152
153 void LAPACK<int, float>::POSVX(const char& FACT, const char& UPLO, const int& n, const int& nrhs, float* A, const int& lda, float* AF, const int& ldaf, char* EQUED, float* S, float* B, const int& ldb, float* X, const int& ldx, float* rcond, float* FERR, float* BERR, float* WORK, int* IWORK, int* info) const
154 { SPOSVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, CHARPTR_MACRO(EQUED), S, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, IWORK, info); }
155
156
157 void LAPACK<int,float>::GELS(const char& TRANS, const int& m, const int& n, const int& nrhs, float* A, const int& lda, float* B, const int& ldb, float* WORK, const int& lwork, int* info) const
158 { SGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, info); }
159
160 void LAPACK<int,float>::GELSS (const int& m, const int& n, const int& nrhs, float* A, const int& lda, float* B, const int& ldb, float* S, const float& rcond, int* rank, float* WORK, const int& lwork, float* rwork, int* info) const
161 {
162 (void) rwork;
163 SGELSS_F77(&m, &n, &nrhs, A, &lda, B, &ldb, S, &rcond, rank, WORK, &lwork, info);
164 }
165
166 void LAPACK<int,float>::GELSS(const int& m, const int& n, const int& nrhs, float* A, const int& lda, float* B, const int& ldb, float* S, const float& rcond, int* rank, float* WORK, const int& lwork, int* info) const
167 { SGELSS_F77(&m, &n, &nrhs, A, &lda, B, &ldb, S, &rcond, rank, WORK, &lwork, info); }
168
169
170 void LAPACK<int,float>::GGLSE(const int& m, const int& n, const int& p, float* A, const int& lda, float* B, const int& ldb, float* C, float* D, float* X, float* WORK, const int& lwork, int* info) const
171 { SGGLSE_F77(&m, &n, &p, A, &lda, B, &ldb, C, D, X, WORK, &lwork, info); }
172
173
174 void LAPACK<int,float>::GEQRF( const int& m, const int& n, float* A, const int& lda, float* TAU, float* WORK, const int& lwork, int* info) const
175 { SGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info); }
176
177 void LAPACK<int,float>::GEQR2 (const int& m, const int& n, float A[], const int& lda, float TAU[], float WORK[], int* const info) const
178 {
179 SGEQR2_F77(&m, &n, A, &lda, TAU, WORK, info);
180 }
181
182 void LAPACK<int,float>::GETRF(const int& m, const int& n, float* A, const int& lda, int* IPIV, int* info) const
183 { SGETRF_F77(&m, &n, A, &lda, IPIV, info); }
184
185
186 void LAPACK<int,float>::GETRS(const char& TRANS, const int& n, const int& nrhs, const float* A, const int& lda, const int* IPIV, float* B, const int& ldb, int* info) const
187 { SGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info); }
188
189
190 void LAPACK<int,float>::LASCL(const char& TYPE, const int& kl, const int& ku, const float& cfrom, const float& cto, const int& m, const int& n, float* A, const int& lda, int* info) const
191 { SLASCL_F77(CHAR_MACRO(TYPE), &kl, &ku, &cfrom, &cto, &m, &n, A, &lda, info); }
192
193 void
195 GEQP3 (const int& m,
196 const int& n,
197 float* A,
198 const int& lda,
199 int* jpvt,
200 float* TAU,
201 float* WORK,
202 const int& lwork,
203 float* RWORK,
204 int* info) const
205 {
206 (void) RWORK;
207 SGEQP3_F77(&m, &n, A, &lda, jpvt, TAU, WORK, &lwork, info);
208 }
209
211 LASWP (const int& N,
212 float A[],
213 const int& LDA,
214 const int& K1,
215 const int& K2,
216 const int IPIV[],
217 const int& INCX) const
218 {
219 SLASWP_F77(&N, A, &LDA, &K1, &K2, IPIV, &INCX);
220 }
221
222 void LAPACK<int,float>::GBTRF(const int& m, const int& n, const int& kl, const int& ku, float* A, const int& lda, int* IPIV, int* info) const
223 { SGBTRF_F77(&m, &n, &kl, &ku, A, &lda, IPIV, info); }
224
225
226 void LAPACK<int,float>::GBTRS(const char& TRANS, const int& n, const int& kl, const int& ku, const int& nrhs, const float* A, const int& lda, const int* IPIV, float* B, const int& ldb, int* info) const
227 { SGBTRS_F77(CHAR_MACRO(TRANS), &n, &kl, &ku, &nrhs, A, &lda, IPIV, B, &ldb, info); }
228
229
230 void LAPACK<int,float>::GTTRF(const int& n, float* dl, float* d, float* du, float* du2, int* IPIV, int* info) const
231 { SGTTRF_F77(&n, dl, d, du, du2, IPIV, info); }
232
233
234 void LAPACK<int,float>::GTTRS(const char& TRANS, const int& n, const int& nrhs, const float* dl, const float* d, const float* du, const float* du2, const int* IPIV, float* B, const int& ldb, int* info) const
235 { SGTTRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, dl, d, du, du2, IPIV, B, &ldb, info); }
236
237
238 void LAPACK<int,float>::GETRI(const int& n, float* A, const int& lda, const int* IPIV, float* WORK, const int& lwork, int* info) const
239 { SGETRI_F77(&n, A, &lda, IPIV, WORK, &lwork, info); }
240
241 void
242 LAPACK<int, float>::LATRS (const char& UPLO,
243 const char& TRANS,
244 const char& DIAG,
245 const char& NORMIN,
246 const int& N,
247 float* A,
248 const int& LDA,
249 float* X,
250 float* SCALE,
251 float* CNORM,
252 int* INFO) const
253 {
254 SLATRS_F77(CHAR_MACRO(UPLO),
255 CHAR_MACRO(TRANS),
256 CHAR_MACRO(DIAG),
257 CHAR_MACRO(NORMIN),
258 &N,
259 A,
260 &LDA,
261 X,
262 SCALE,
263 CNORM,
264 INFO);
265 }
266
267
268 void LAPACK<int,float>::GECON(const char& NORM, const int& n, const float* A, const int& lda, const float& anorm, float* rcond, float* WORK, int* IWORK, int* info) const
269 { SGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, IWORK, info); }
270
271
272 void LAPACK<int,float>::GBCON(const char& NORM, const int& n, const int& kl, const int& ku, const float* A, const int& lda, int* IPIV, const float& anorm, float* rcond, float* WORK, int* IWORK, int* info) const
273 { SGBCON_F77(CHAR_MACRO(NORM), &n, &kl, &ku, A, &lda, IPIV, &anorm, rcond, WORK, IWORK, info); }
274
275
276 float LAPACK<int,float>::LANGB(const char& NORM, const int& n, const int& kl, const int& ku, const float* A, const int& lda, float* WORK) const
277 { return( SLANGB_F77(CHAR_MACRO(NORM), &n, &kl, &ku, A, &lda, WORK) ); }
278
279
280 void LAPACK<int,float>::GESV(const int& n, const int& nrhs, float* A, const int& lda, int* IPIV, float* B, const int& ldb, int* info) const
281 { SGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info); }
282
283
284 void LAPACK<int,float>::GEEQU(const int& m, const int& n, const float* A, const int& lda, float* R, float* C, float* rowcond, float* colcond, float* amax, int* info) const
285 { SGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info); }
286
287
288 void LAPACK<int,float>::GERFS(const char& TRANS, const int& n, const int& nrhs, const float* A, const int& lda, const float* AF, const int& ldaf, const int* IPIV, const float* B, const int& ldb, float* X, const int& ldx, float* FERR, float* BERR, float* WORK, int* IWORK, int* info) const
289 { SGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info); }
290
291
292 void LAPACK<int,float>::GBEQU(const int& m, const int& n, const int& kl, const int& ku, const float* A, const int& lda, float* R, float* C, float* rowcond, float* colcond, float* amax, int* info) const
293 { SGBEQU_F77(&m, &n, &kl, &ku, A, &lda, R, C, rowcond, colcond, amax, info); }
294
295
296 void LAPACK<int,float>::GBRFS(const char& TRANS, const int& n, const int& kl, const int& ku, const int& nrhs, const float* A, const int& lda, const float* AF, const int& ldaf, const int* IPIV, const float* B, const int& ldb, float* X, const int& ldx, float* FERR, float* BERR, float* WORK, int* IWORK, int* info) const
297 { SGBRFS_F77(CHAR_MACRO(TRANS), &n, &kl, &ku, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info); }
298
299 void LAPACK<int,float>::GESVX(const char& FACT, const char& TRANS, const int& n, const int& nrhs, float* A, const int& lda, float* AF, const int& ldaf, int* IPIV, char* EQUED, float* R, float* C, float* B, const int& ldb, float* X, const int& ldx, float* rcond, float* FERR, float* BERR, float* WORK, int* IWORK, int* info) const
300 { SGESVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, CHARPTR_MACRO(EQUED), R, C, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, IWORK, info); }
301
302
303 void LAPACK<int,float>::SYTRD(const char& UPLO, const int& n, float* A, const int& lda, float* D, float* E, float* TAU, float* WORK, const int& lwork, int* info) const
304 { SSYTRD_F77(CHAR_MACRO(UPLO), &n, A, &lda, D, E, TAU, WORK, &lwork, info); }
305
306
307 void LAPACK<int,float>::GEHRD(const int& n, const int& ilo, const int& ihi, float* A, const int& lda, float* TAU, float* WORK, const int& lwork, int* info) const
308 { SGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info); }
309
310
311 void LAPACK<int,float>::TRTRS(const char& UPLO, const char& TRANS, const char& DIAG, const int& n, const int& nrhs, const float* A, const int& lda, float* B, const int& ldb, int* info) const
312 { STRTRS_F77(CHAR_MACRO(UPLO), CHAR_MACRO(TRANS), CHAR_MACRO(DIAG), &n, &nrhs, A, &lda, B, &ldb, info); }
313
314
315 void LAPACK<int,float>::TRTRI(const char& UPLO, const char& DIAG, const int& n, const float* A, const int& lda, int* info) const
316 { STRTRI_F77(CHAR_MACRO(UPLO), CHAR_MACRO(DIAG), &n, A, &lda, info); }
317
318
319 void LAPACK<int,float>::SPEV(const char& JOBZ, const char& UPLO, const int& n, float* AP, float* W, float* Z, const int& ldz, float* WORK, int* info) const
320 { SSPEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, AP, W, Z, &ldz, WORK, info); }
321
322
323 void LAPACK<int,float>::SYEV(const char& JOBZ, const char& UPLO, const int& n, float* A, const int& lda, float* W, float* WORK, const int& lwork, int* info) const
324 { SSYEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, W, WORK, &lwork, info); }
325
326
327 void LAPACK<int,float>::SYGV(const int& itype, const char& JOBZ, const char& UPLO, const int& n, float* A, const int& lda, float* B, const int& ldb, float* W, float* WORK, const int& lwork, int* info) const
328 { SSYGV_F77(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, B, &ldb, W, WORK, &lwork, info); }
329
330
331 void LAPACK<int,float>::HEEV(const char& JOBZ, const char& UPLO, const int& n, float* A, const int& lda, float* W, float* WORK, const int& lwork, float* /* RWORK */, int* info) const
332 { SSYEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, W, WORK, &lwork, info); }
333
334
335 void LAPACK<int,float>::HEGV(const int& itype, const char& JOBZ, const char& UPLO, const int& n, float* A, const int& lda, float* B, const int& ldb, float* W, float* WORK, const int& lwork, float* /* RWORK */, int* info) const
336 { SSYGV_F77(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, B, &ldb, W, WORK, &lwork, info); }
337
338
339 void LAPACK<int,float>::STEQR(const char& COMPZ, const int& n, float* D, float* E, float* Z, const int& ldz, float* WORK, int* info) const
340 { SSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info); }
341
342
343 void LAPACK<int,float>::PTEQR(const char& COMPZ, const int& n, float* D, float* E, float* Z, const int& ldz, float* WORK, int* info) const
344 { SPTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info); }
345
346
347 void LAPACK<int, float>::HSEQR(const char& JOB, const char& COMPZ, const int& n, const int& ilo, const int& ihi, float* H, const int& ldh, float* WR, float* WI, float* Z, const int& ldz, float* WORK, const int& lwork, int* info) const
348 { SHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, WR, WI, Z, &ldz, WORK, &lwork, info); }
349
350
351 void LAPACK<int, float>::GEES(const char& JOBVS, const char& SORT, int (*ptr2func)(float*, float*), const int& n, float* A, const int& lda, int* sdim, float* WR, float* WI, float* VS, const int& ldvs, float* WORK, const int& lwork, int* BWORK, int* info) const
352 { SGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), ptr2func, &n, A, &lda, sdim, WR, WI, VS, &ldvs, WORK, &lwork, BWORK, info); }
353
354
355 void LAPACK<int, float>::GEES(const char& JOBVS, const int& n, float* A, const int& lda, int* sdim, float* WR, float* WI, float* VS, const int& ldvs, float* WORK, const int& lwork, float* /* RWORK */, int* BWORK, int* info) const
356 {
357 int (*nullfptr)(float*,float*) = NULL;
358 const char sort = 'N';
359 SGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(sort), nullfptr, &n, A, &lda, sdim, WR, WI, VS, &ldvs, WORK, &lwork, BWORK, info);
360 }
361
362
363 void LAPACK<int, float>::GEEV(const char& JOBVL, const char& JOBVR, const int& n, float* A, const int& lda, float* WR, float* WI, float* VL, const int& ldvl, float* VR, const int& ldvr, float* WORK, const int& lwork, int* info) const
364 { SGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, WR, WI, VL, &ldvl, VR, &ldvr, WORK, &lwork, info); }
365
366 void LAPACK<int, float>::GEEV(const char& JOBVL, const char& JOBVR, const int& n, float* A, const int& lda, float* WR, float* WI, float* VL, const int& ldvl, float* VR, const int& ldvr, float* WORK, const int& lwork, float* /* RWORK */, int* info) const
367 {
368 GEEV (JOBVL, JOBVR, n, A, lda, WR, WI, VL, ldvl, VR, ldvr, WORK, lwork, info);
369 }
370
371
372 void LAPACK<int, float>::GESVD(const char& JOBU, const char& JOBVT, const int& m, const int& n, float* A, const int& lda, float* S, float* U, const int& ldu, float* V, const int& ldv, float* WORK, const int& lwork, float* /* RWORK */, int* info) const
373 { SGESVD_F77(CHAR_MACRO(JOBU), CHAR_MACRO(JOBVT), &m, &n, A, &lda, S, U, &ldu, V, &ldv, WORK, &lwork, info); }
374
375
376 void LAPACK<int,float>::GEEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, float* A, const int& lda, float* WR, float* WI, float* VL, const int& ldvl, float* VR, const int& ldvr, int* ilo, int* ihi, float* SCALE, float* abnrm, float* RCONDE, float* RCONDV, float* WORK, const int& lwork, int* IWORK, int* info) const
377 { SGEEVX_F77(CHAR_MACRO(BALANC), CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), CHAR_MACRO(SENSE), &n, A, &lda, WR, WI, VL, &ldvl, VR, &ldvr, ilo, ihi, SCALE, abnrm, RCONDE, RCONDV, WORK, &lwork, IWORK, info); }
378
379
380 void LAPACK<int,float>::GGEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, float* A, const int& lda, float* B, const int& ldb, float* ALPHAR, float* ALPHAI, float* BETA, float* VL, const int& ldvl, float* VR, const int& ldvr, int* ilo, int* ihi, float* lscale, float* rscale, float* abnrm, float* bbnrm, float* RCONDE, float* RCONDV, float* WORK, const int& lwork, int* IWORK, int* BWORK, int* info) const
381 { SGGEVX_F77(CHAR_MACRO(BALANC), CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), CHAR_MACRO(SENSE), &n, A, &lda, B, &ldb, ALPHAR, ALPHAI, BETA, VL, &ldvl, VR, &ldvr, ilo, ihi, lscale, rscale, abnrm, bbnrm, RCONDE, RCONDV, WORK, &lwork, IWORK, BWORK, info); }
382
383 void LAPACK<int,float>::GGEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, float* A, const int& lda, float* B, const int& ldb, float* ALPHAR, float* ALPHAI, float* BETA, float* VL, const int& ldvl, float* VR, const int& ldvr, int* ilo, int* ihi, float* lscale, float* rscale, float* abnrm, float* bbnrm, float* RCONDE, float* RCONDV, float* WORK, const int& lwork, float* /* RWORK */, int* IWORK, int* BWORK, int* info) const
384 {
386 }
387
388 void LAPACK<int, float>::GGEV(const char& JOBVL, const char& JOBVR, const int& n, float* A, const int& lda, float* B, const int& ldb, float* ALPHAR, float* ALPHAI, float* BETA, float* VL, const int& ldvl, float* VR, const int& ldvr, float* WORK, const int& lwork, int* info) const
389 { SGGEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, B, &ldb, ALPHAR, ALPHAI, BETA, VL, &ldvl, VR, &ldvr, WORK, &lwork, info); }
390
391
392 void LAPACK<int, float>::TRSEN(const char& JOB, const char& COMPQ, const int* SELECT, const int& n, float* T, const int& ldt, float* Q, const int& ldq, float* WR, float* WI, int* M, float* S, float* SEP, float* WORK, const int& lwork, int* IWORK, const int& liwork, int* info ) const
393 { STRSEN_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPQ), SELECT, &n, T, &ldt, Q, &ldq, WR, WI, M, S, SEP, WORK, &lwork, IWORK, &liwork, info); }
394
395
396 void LAPACK<int, float>::TGSEN(const int& ijob, const int& wantq, const int& wantz, const int* SELECT, const int& n, float* A, const int& lda, float* B, const int& ldb, float* ALPHAR, float* ALPHAI, float* BETA, float* Q, const int& ldq, float* Z, const int& ldz, int* M, float* PL, float* PR, float* DIF, float* WORK, const int& lwork, int* IWORK, const int& liwork, int* info ) const
397 { STGSEN_F77(&ijob, &wantq, &wantz, SELECT, &n, A, &lda, B, &ldb, ALPHAR, ALPHAI, BETA, Q, &ldq, Z, &ldz, M, PL, PR, DIF, WORK, &lwork, IWORK, &liwork, info); }
398
399
400 void LAPACK<int, float>::GGES(const char& JOBVL, const char& JOBVR, const char& SORT, int (*ptr2func)(float* , float* , float* ), const int& n, float* A, const int& lda, float* B, const int& ldb, int* sdim, float* ALPHAR, float* ALPHAI, float* BETA, float* VL, const int& ldvl, float* VR, const int& ldvr, float* WORK, const int& lwork, int* BWORK, int* info ) const
401 { SGGES_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), CHAR_MACRO(SORT), ptr2func, &n, A, &lda, B, &ldb, sdim, ALPHAR, ALPHAI, BETA, VL, &ldvl, VR, &ldvr, WORK, &lwork, BWORK, info); }
402
403
404 void LAPACK<int, float>::ORMQR(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& k, float* A, const int& lda, const float* TAU, float* C, const int& ldc, float* WORK, const int& lwork, int* info) const
405 { SORMQR_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &k, A, &lda, TAU, C, &ldc, WORK, &lwork, info); }
406
407
408 void LAPACK<int, float>::ORM2R(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& k, const float A[], const int& lda, const float TAU[], float C[], const int& ldc, float WORK[], int* const info) const
409 { SORM2R_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &k, A, &lda, TAU, C, &ldc, WORK, info); }
410
411
412 void LAPACK<int, float>::UNMQR(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& k, float* A, const int& lda, const float* TAU, float* C, const int& ldc, float* WORK, const int& lwork, int* info) const
413 {
414 // LAPACK only defines UNMQR for Z and C (complex*8
415 // resp. complex*16), but logically, UNMQR means the same thing as
416 // ORMQR for real arithmetic.
417 ORMQR (SIDE, TRANS, m, n, k, A, lda, TAU, C, ldc, WORK, lwork, info);
418 }
419
420 void LAPACK<int, float>::UNM2R (const char& SIDE, const char& TRANS, const int& M, const int& N, const int& K, const float A[], const int& LDA, const float TAU[], float C[], const int& LDC, float WORK[], int* const INFO) const
421 {
422 // LAPACK only defines UNM2R for Z and C (complex*8
423 // resp. complex*16), but logically, UNM2R means the same thing as
424 // ORM2R for real arithmetic.
425 ORM2R (SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, INFO);
426 }
427
428
429 void LAPACK<int, float>::ORGQR(const int& m, const int& n, const int& k, float* A, const int& lda, const float* TAU, float* WORK, const int& lwork, int* info) const
430 { SORGQR_F77( &m, &n, &k, A, &lda, TAU, WORK, &lwork, info); }
431
432
433 void LAPACK<int, float>::UNGQR(const int& m, const int& n, const int& k, float* A, const int& lda, const float* TAU, float* WORK, const int& lwork, int* info) const
434 { SORGQR_F77( &m, &n, &k, A, &lda, TAU, WORK, &lwork, info); }
435
436
437 void LAPACK<int, float>::ORGHR(const int& n, const int& ilo, const int& ihi, float* A, const int& lda, const float* TAU, float* WORK, const int& lwork, int* info) const
438 { SORGHR_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info); }
439
440
441 void LAPACK<int, float>::ORMHR(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& ilo, const int& ihi, const float* A, const int& lda, const float* TAU, float* C, const int& ldc, float* WORK, const int& lwork, int* info) const
442 { SORMHR_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &ilo, &ihi, A, &lda, TAU, C, &ldc, WORK, &lwork, info); }
443
444
445 void LAPACK<int, float>::TREVC(const char& SIDE, const char& HOWMNY, int* select, const int& n, const float* T, const int& ldt, float* VL, const int& ldvl, float* VR, const int& ldvr, const int& mm, int* m, float* WORK, int* info) const
446 { STREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), select, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, info); }
447
448
449 void LAPACK<int, float>::TREVC(const char& SIDE, const int& n, const float* T, const int& ldt, float* VL, const int& ldvl, float* VR, const int& ldvr, const int& mm, int* m, float* WORK, float* /* RWORK */, int* info) const
450 {
451 std::vector<int> select(1);
452 const char whch = 'A';
453 STREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(whch), &select[0], &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, info);
454 }
455
456 void LAPACK<int, float>::TREXC(const char& COMPQ, const int& n, float* T, const int& ldt, float* Q, const int& ldq, int* ifst, int* ilst, float* WORK, int* info) const
457 { STREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, ifst, ilst, WORK, info); }
458
459
460 void LAPACK<int, float>::TGEVC(const char& SIDE, const char& HOWMNY, const int* SELECT, const int& n, float* S, const int& lds, float* P, const int& ldp, float* VL, const int& ldvl, float* VR, const int& ldvr, const int& mm, int* M, float* WORK, int* info) const
461 { STGEVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), SELECT, &n, S, &lds, P, &ldp, VL, &ldvl, VR, &ldvr, &mm, M, WORK, info); }
462
463
464 void LAPACK<int, float>::LARTG( const float& f, const float& g, float* c, float* s, float* r ) const
465 { SLARTG_F77(&f, &g, c, s, r); }
466
467
468 void LAPACK<int, float>::LARFG( const int& n, float* alpha, float* x, const int& incx, float* tau ) const
469 { SLARFG_F77(&n, alpha, x, &incx, tau); }
470
471 void LAPACK<int, float>::GEBAL(const char& JOBZ, const int& n, float* A, const int& lda, int* ilo, int* ihi, float* scale, int* info) const
472 { SGEBAL_F77(CHAR_MACRO(JOBZ),&n, A, &lda, ilo, ihi, scale, info); }
473
474
475 void LAPACK<int, float>::GEBAK(const char& JOBZ, const char& SIDE, const int& n, const int& ilo, const int& ihi, const float* scale, const int& m, float* V, const int& ldv, int* info) const
476 { SGEBAK_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(SIDE), &n, &ilo, &ihi, scale, &m, V, &ldv, info); }
477
478#ifdef HAVE_TEUCHOS_LAPACKLARND
479 float LAPACK<int, float>::LARND( const int& idist, int* seed ) const
480 { return(SLARND_F77(&idist, seed)); }
481#endif
482
483 void LAPACK<int, float>::LARNV( const int& idist, int* seed, const int& n, float* v ) const
484 { SLARNV_F77(&idist, seed, &n, v); }
485
486
487 float LAPACK<int, float>::LAMCH(const char& CMACH) const
488 { return(SLAMCH_F77(CHAR_MACRO(CMACH))); }
489
490
491 int LAPACK<int, float>::ILAENV( const int& ispec, const std::string& NAME, const std::string& OPTS, const int& N1, const int& N2, const int& N3, const int& N4 ) const
492 {
493 unsigned int opts_length = OPTS.length();
494 // if user queries a Hermitian routine, change it to a symmetric routine
495 std::string temp_NAME = "s" + NAME;
496 if (temp_NAME.substr(1,2) == "he") {
497 temp_NAME.replace(1,2,"sy");
498 }
499 unsigned int name_length = temp_NAME.length();
500 return ilaenv_wrapper(&ispec, &temp_NAME[0], name_length, &OPTS[0], opts_length, &N1, &N2, &N3, &N4);
501 }
502
503
504 float LAPACK<int, float>::LAPY2(const float& x, const float& y) const
505 {
506#if defined(HAVE_TEUCHOS_BLASFLOAT)
507 return SLAPY2_F77(&x, &y);
508#else
509 typedef ScalarTraits<float> ST;
510 const float xabs = ST::magnitude(x);
511 const float yabs = ST::magnitude(y);
512 const float w = TEUCHOS_MAX(xabs, yabs);
513 const float z = TEUCHOS_MIN(xabs, yabs);
514 if ( z == 0.0 ) {
515 return w;
516 }
517 const float z_over_w = z/w;
518 return w*ST::squareroot( 1.0+(z_over_w*z_over_w));
519#endif
520 }
521
522 // END INT, FLOAT SPECIALIZATION IMPLEMENTATION //
523
524 // BEGIN INT, DOUBLE SPECIALIZATION IMPLEMENTATION //
525
526
527
528 void LAPACK<int, double>::PTTRF(const int& n, double* d, double* e, int* info) const
529 { DPTTRF_F77(&n,d,e,info); }
530
531
532 void LAPACK<int, double>::PTTRS(const int& n, const int& nrhs, const double* d, const double* e, double* B, const int& ldb, int* info) const
533 { DPTTRS_F77(&n,&nrhs,d,e,B,&ldb,info); }
534
535
536 void LAPACK<int, double>::POTRF(const char& UPLO, const int& n, double* A, const int& lda, int* info) const
537 { DPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info); }
538
539
540 void LAPACK<int, double>::POTRS(const char& UPLO, const int& n, const int& nrhs, const double* A, const int& lda, double* B, const int& ldb, int* info) const
541 { DPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info); }
542
543
544 void LAPACK<int, double>::POTRI(const char& UPLO, const int& n, double* A, const int& lda, int* info) const
545 { DPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info); }
546
547
548 void LAPACK<int, double>::POCON(const char& UPLO, const int& n, const double* A, const int& lda, const double& anorm, double* rcond, double* WORK, int* IWORK, int* info) const
549 { DPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, IWORK, info); }
550
551
552 void LAPACK<int, double>::POSV(const char& UPLO, const int& n, const int& nrhs, double* A, const int& lda, double* B, const int& ldb, int* info) const
553 { DPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info); }
554
555
556 void LAPACK<int, double>::POEQU(const int& n, const double* A, const int& lda, double* S, double* scond, double* amax, int* info) const
557 { DPOEQU_F77(&n, A, &lda, S, scond, amax, info); }
558
559
560 void LAPACK<int, double>::PORFS(const char& UPLO, const int& n, const int& nrhs, double* A, const int& lda, const double* AF, const int& ldaf, const double* B, const int& ldb, double* X, const int& ldx, double* FERR, double* BERR, double* WORK, int* IWORK, int* info) const
561 { DPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info); }
562
563 void LAPACK<int, double>::POSVX(const char& FACT, const char& UPLO, const int& n, const int& nrhs, double* A, const int& lda, double* AF, const int& ldaf, char* EQUED, double* S, double* B, const int& ldb, double* X, const int& ldx, double* rcond, double* FERR, double* BERR, double* WORK, int* IWORK, int* info) const
564 { DPOSVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, CHARPTR_MACRO(EQUED), S, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, IWORK, info); }
565
566
567 void LAPACK<int,double>::GELS(const char& TRANS, const int& m, const int& n, const int& nrhs, double* A, const int& lda, double* B, const int& ldb, double* WORK, const int& lwork, int* info) const
568 { DGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, info); }
569
570
571 void LAPACK<int,double>::GELSS(const int& m, const int& n, const int& nrhs, double* A, const int& lda, double* B, const int& ldb, double* S, const double& rcond, int* rank, double* WORK, const int& lwork, double* rwork, int* info) const
572 {
573 (void) rwork;
574 DGELSS_F77(&m, &n, &nrhs, A, &lda, B, &ldb, S, &rcond, rank, WORK, &lwork, info);
575 }
576
577
578 void LAPACK<int,double>::GELSS(const int& m, const int& n, const int& nrhs, double* A, const int& lda, double* B, const int& ldb, double* S, const double& rcond, int* rank, double* WORK, const int& lwork, int* info) const
579 { DGELSS_F77(&m, &n, &nrhs, A, &lda, B, &ldb, S, &rcond, rank, WORK, &lwork, info); }
580
581
582 void LAPACK<int,double>::GGLSE(const int& m, const int& n, const int& p, double* A, const int& lda, double* B, const int& ldb, double* C, double* D, double* X, double* WORK, const int& lwork, int* info) const
583 { DGGLSE_F77(&m, &n, &p, A, &lda, B, &ldb, C, D, X, WORK, &lwork, info); }
584
585
586 void LAPACK<int,double>::GEQRF( const int& m, const int& n, double* A, const int& lda, double* TAU, double* WORK, const int& lwork, int* info) const
587 { DGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info); }
588
589 void LAPACK<int,double>::GEQR2 (const int& m, const int& n, double A[], const int& lda, double TAU[], double WORK[], int* const info) const
590 {
591 DGEQR2_F77(&m, &n, A, &lda, TAU, WORK, info);
592 }
593
594 void LAPACK<int,double>::GETRF(const int& m, const int& n, double* A, const int& lda, int* IPIV, int* info) const
595 { DGETRF_F77(&m, &n, A, &lda, IPIV, info); }
596
597
598 void LAPACK<int,double>::GETRS(const char& TRANS, const int& n, const int& nrhs, const double* A, const int& lda, const int* IPIV, double* B, const int& ldb, int* info) const
599 { DGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info); }
600
601
602 void LAPACK<int,double>::LASCL(const char& TYPE, const int& kl, const int& ku, const double& cfrom, const double& cto, const int& m, const int& n, double* A, const int& lda, int* info) const
603 { DLASCL_F77(CHAR_MACRO(TYPE), &kl, &ku, &cfrom, &cto, &m, &n, A, &lda, info); }
604
605 void LAPACK<int,double>::GEQP3(const int& m, const int& n, double* A, const int& lda, int* jpvt, double* TAU, double* WORK, const int& lwork, double* RWORK, int* info ) const
606 {
607 (void) RWORK;
608 DGEQP3_F77(&m, &n, A, &lda, jpvt, TAU, WORK, &lwork, info);
609 }
610
612 LASWP (const int& N,
613 double A[],
614 const int& LDA,
615 const int& K1,
616 const int& K2,
617 const int IPIV[],
618 const int& INCX) const
619 {
620 DLASWP_F77(&N, A, &LDA, &K1, &K2, IPIV, &INCX);
621 }
622
623 void LAPACK<int,double>::GBTRF(const int& m, const int& n, const int& kl, const int& ku, double* A, const int& lda, int* IPIV, int* info) const
624 { DGBTRF_F77(&m, &n, &kl, &ku, A, &lda, IPIV, info); }
625
626
627 void LAPACK<int,double>::GBTRS(const char& TRANS, const int& n, const int& kl, const int& ku, const int& nrhs, const double* A, const int& lda, const int* IPIV, double* B, const int& ldb, int* info) const
628 { DGBTRS_F77(CHAR_MACRO(TRANS), &n, &kl, &ku, &nrhs, A, &lda, IPIV, B, &ldb, info); }
629
630
631 void LAPACK<int,double>::GTTRF(const int& n, double* dl, double* d, double* du, double* du2, int* IPIV, int* info) const
632 { DGTTRF_F77(&n, dl, d, du, du2, IPIV, info); }
633
634
635 void LAPACK<int,double>::GTTRS(const char& TRANS, const int& n, const int& nrhs, const double* dl, const double* d, const double* du, const double* du2, const int* IPIV, double* B, const int& ldb, int* info) const
636 { DGTTRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, dl, d, du, du2, IPIV, B, &ldb, info); }
637
638
639 void LAPACK<int,double>::GETRI(const int& n, double* A, const int& lda, const int* IPIV, double* WORK, const int& lwork, int* info) const
640 { DGETRI_F77(&n, A, &lda, IPIV, WORK, &lwork, info); }
641
642 void
643 LAPACK<int, double>::LATRS (const char& UPLO,
644 const char& TRANS,
645 const char& DIAG,
646 const char& NORMIN,
647 const int& N,
648 double* A,
649 const int& LDA,
650 double* X,
651 double* SCALE,
652 double* CNORM,
653 int* INFO) const
654 {
655 DLATRS_F77(CHAR_MACRO(UPLO),
656 CHAR_MACRO(TRANS),
657 CHAR_MACRO(DIAG),
658 CHAR_MACRO(NORMIN),
659 &N,
660 A,
661 &LDA,
662 X,
663 SCALE,
664 CNORM,
665 INFO);
666 }
667
668 void LAPACK<int,double>::GECON(const char& NORM, const int& n, const double* A, const int& lda, const double& anorm, double* rcond, double* WORK, int* IWORK, int* info) const
669 { DGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, IWORK, info); }
670
671
672 void LAPACK<int,double>::GBCON(const char& NORM, const int& n, const int& kl, const int& ku, const double* A, const int& lda, int* IPIV, const double& anorm, double* rcond, double* WORK, int* IWORK, int* info) const
673 { DGBCON_F77(CHAR_MACRO(NORM), &n, &kl, &ku, A, &lda, IPIV, &anorm, rcond, WORK, IWORK, info); }
674
675
676 double LAPACK<int,double>::LANGB(const char& NORM, const int& n, const int& kl, const int& ku, const double* A, const int& lda, double* WORK) const
677 { return( DLANGB_F77(CHAR_MACRO(NORM), &n, &kl, &ku, A, &lda, WORK) ); }
678
679
680 void LAPACK<int,double>::GESV(const int& n, const int& nrhs, double* A, const int& lda, int* IPIV, double* B, const int& ldb, int* info) const
681 { DGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info); }
682
683
684 void LAPACK<int,double>::GEEQU(const int& m, const int& n, const double* A, const int& lda, double* R, double* C, double* rowcond, double* colcond, double* amax, int* info) const
685 { DGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info); }
686
687
688 void LAPACK<int,double>::GERFS(const char& TRANS, const int& n, const int& nrhs, const double* A, const int& lda, const double* AF, const int& ldaf, const int* IPIV, const double* B, const int& ldb, double* X, const int& ldx, double* FERR, double* BERR, double* WORK, int* IWORK, int* info) const
689 { DGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info); }
690
691
692 void LAPACK<int,double>::GBEQU(const int& m, const int& n, const int& kl, const int& ku, const double* A, const int& lda, double* R, double* C, double* rowcond, double* colcond, double* amax, int* info) const
693 { DGBEQU_F77(&m, &n, &kl, &ku, A, &lda, R, C, rowcond, colcond, amax, info); }
694
695
696 void LAPACK<int,double>::GBRFS(const char& TRANS, const int& n, const int& kl, const int& ku, const int& nrhs, const double* A, const int& lda, const double* AF, const int& ldaf, const int* IPIV, const double* B, const int& ldb, double* X, const int& ldx, double* FERR, double* BERR, double* WORK, int* IWORK, int* info) const
697 { DGBRFS_F77(CHAR_MACRO(TRANS), &n, &kl, &ku, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, IWORK, info); }
698
699 void LAPACK<int,double>::GESVX(const char& FACT, const char& TRANS, const int& n, const int& nrhs, double* A, const int& lda, double* AF, const int& ldaf, int* IPIV, char* EQUED, double* R, double* C, double* B, const int& ldb, double* X, const int& ldx, double* rcond, double* FERR, double* BERR, double* WORK, int* IWORK, int* info) const
700 { DGESVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, CHARPTR_MACRO(EQUED), R, C, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, IWORK, info); }
701
702
703 void LAPACK<int,double>::SYTRD(const char& UPLO, const int& n, double* A, const int& lda, double* D, double* E, double* TAU, double* WORK, const int& lwork, int* info) const
704 { DSYTRD_F77(CHAR_MACRO(UPLO), &n, A, &lda, D, E, TAU, WORK, &lwork, info); }
705
706
707 void LAPACK<int, double>::GEHRD(const int& n, const int& ilo, const int& ihi, double* A, const int& lda, double* TAU, double* WORK, const int& lwork, int* info) const
708 { DGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info); }
709
710
711 void LAPACK<int,double>::TRTRS(const char& UPLO, const char& TRANS, const char& DIAG, const int& n, const int& nrhs, const double* A, const int& lda, double* B, const int& ldb, int* info) const
712 { DTRTRS_F77(CHAR_MACRO(UPLO), CHAR_MACRO(TRANS), CHAR_MACRO(DIAG), &n, &nrhs, A, &lda, B, &ldb, info); }
713
714
715 void LAPACK<int,double>::TRTRI(const char& UPLO, const char& DIAG, const int& n, const double* A, const int& lda, int* info) const
716 { DTRTRI_F77(CHAR_MACRO(UPLO), CHAR_MACRO(DIAG), &n, A, &lda, info); }
717
718
719 void LAPACK<int,double>::SPEV(const char& JOBZ, const char& UPLO, const int& n, double* AP, double* W, double* Z, const int& ldz, double* WORK, int* info) const
720 { DSPEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, AP, W, Z, &ldz, WORK, info); }
721
722
723 void LAPACK<int,double>::SYEV(const char& JOBZ, const char& UPLO, const int& n, double* A, const int& lda, double* W, double* WORK, const int& lwork, int* info) const
724 {
725 DSYEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, W, WORK, &lwork, info);
726 }
727
728
729 void LAPACK<int,double>::SYGV(const int& itype, const char& JOBZ, const char& UPLO, const int& n, double* A, const int& lda, double* B, const int& ldb, double* W, double* WORK, const int& lwork, int* info) const
730 {
731 DSYGV_F77(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, B, &ldb, W, WORK, &lwork, info);
732 }
733
734
735 void LAPACK<int,double>::HEEV(const char& JOBZ, const char& UPLO, const int& n, double* A, const int& lda, double* W, double* WORK, const int& lwork, double* /* RWORK */, int* info) const
736 {
737 DSYEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, W, WORK, &lwork, info);
738 }
739
740
741 void LAPACK<int,double>::HEGV(const int& itype, const char& JOBZ, const char& UPLO, const int& n, double* A, const int& lda, double* B, const int& ldb, double* W, double* WORK, const int& lwork, double* /* RWORK */, int* info) const
742 {
743 DSYGV_F77(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, B, &ldb, W, WORK, &lwork, info);
744 }
745
746
747 void LAPACK<int,double>::STEQR(const char& COMPZ, const int& n, double* D, double* E, double* Z, const int& ldz, double* WORK, int* info) const
748 { DSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info); }
749
750
751 void LAPACK<int,double>::PTEQR(const char& COMPZ, const int& n, double* D, double* E, double* Z, const int& ldz, double* WORK, int* info) const
752 { DPTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info); }
753
754
755 void LAPACK<int, double>::HSEQR(const char& JOB, const char& COMPZ, const int& n, const int& ilo, const int& ihi, double* H, const int& ldh, double* WR, double* WI, double* Z, const int& ldz, double* WORK, const int& lwork, int* info) const
756 {
757 DHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, WR, WI, Z, &ldz, WORK, &lwork, info);
758 }
759
760
761 void LAPACK<int, double>::GEES(const char& JOBVS, const char& SORT, int (*ptr2func)(double*, double*), const int& n, double* A, const int& lda, int* sdim, double* WR, double* WI, double* VS, const int& ldvs, double* WORK, const int& lwork, int* BWORK, int* info) const
762 {
763 DGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), ptr2func, &n, A, &lda, sdim, WR, WI, VS, &ldvs, WORK, &lwork, BWORK, info);
764 }
765
766
767 void LAPACK<int, double>::GEES(const char& JOBVS, const int& n, double* A, const int& lda, int* sdim, double* WR, double* WI, double* VS, const int& ldvs, double* WORK, const int& lwork, double* /* RWORK */, int* BWORK, int* info) const
768 {
769 //int (*nullfptr)(double*,double*) = NULL;
770 gees_nullfptr_t nullfptr = 0;
771 const char sort = 'N';
772 DGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(sort), nullfptr, &n, A, &lda, sdim, WR, WI, VS, &ldvs, WORK, &lwork, BWORK, info);
773 }
774
775
776 void LAPACK<int, double>::GEEV(const char& JOBVL, const char& JOBVR, const int& n, double* A, const int& lda, double* WR, double* WI, double* VL, const int& ldvl, double* VR, const int& ldvr, double* WORK, const int& lwork, int* info) const
777 {
778 DGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, WR, WI, VL, &ldvl, VR, &ldvr, WORK, &lwork, info);
779 }
780
781 void LAPACK<int, double>::GEEV(const char& JOBVL, const char& JOBVR, const int& n, double* A, const int& lda, double* WR, double* WI, double* VL, const int& ldvl, double* VR, const int& ldvr, double* WORK, const int& lwork, double* /* RWORK */, int* info) const
782 {
783 GEEV (JOBVL, JOBVR, n, A, lda, WR, WI, VL, ldvl, VR, ldvr, WORK, lwork, info);
784 }
785
786
787 void LAPACK<int, double>::GESVD(const char& JOBU, const char& JOBVT, const int& m, const int& n, double* A, const int& lda, double* S, double* U, const int& ldu, double* V, const int& ldv, double* WORK, const int& lwork, double* /* RWORK */, int* info) const {
788 DGESVD_F77(CHAR_MACRO(JOBU), CHAR_MACRO(JOBVT), &m, &n, A, &lda, S, U, &ldu, V, &ldv, WORK, &lwork, info);
789 }
790
791
792 void LAPACK<int,double>::GEEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, double* A, const int& lda, double* WR, double* WI, double* VL, const int& ldvl, double* VR, const int& ldvr, int* ilo, int* ihi, double* SCALE, double* abnrm, double* RCONDE, double* RCONDV, double* WORK, const int& lwork, int* IWORK, int* info) const
793 {
794 DGEEVX_F77(CHAR_MACRO(BALANC), CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), CHAR_MACRO(SENSE), &n, A, &lda, WR, WI, VL, &ldvl, VR, &ldvr, ilo, ihi, SCALE, abnrm, RCONDE, RCONDV, WORK, &lwork, IWORK, info);
795 }
796
797
798 void LAPACK<int, double>::GGEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, double* A, const int& lda, double* B, const int& ldb, double* ALPHAR, double* ALPHAI, double* BETA, double* VL, const int& ldvl, double* VR, const int& ldvr, int* ilo, int* ihi, double* lscale, double* rscale, double* abnrm, double* bbnrm, double* RCONDE, double* RCONDV, double* WORK, const int& lwork, int* IWORK, int* BWORK, int* info) const
799 {
800 DGGEVX_F77(CHAR_MACRO(BALANC), CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), CHAR_MACRO(SENSE), &n, A, &lda, B, &ldb, ALPHAR, ALPHAI, BETA, VL, &ldvl, VR, &ldvr, ilo, ihi, lscale, rscale, abnrm, bbnrm, RCONDE, RCONDV, WORK, &lwork, IWORK, BWORK, info);
801 }
802
803 void LAPACK<int, double>::GGEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, double* A, const int& lda, double* B, const int& ldb, double* ALPHAR, double* ALPHAI, double* BETA, double* VL, const int& ldvl, double* VR, const int& ldvr, int* ilo, int* ihi, double* lscale, double* rscale, double* abnrm, double* bbnrm, double* RCONDE, double* RCONDV, double* WORK, const int& lwork, double* /* RWORK */, int* IWORK, int* BWORK, int* info) const
804 {
805 GGEVX(BALANC, JOBVL, JOBVR, SENSE, n, A, lda, B, ldb, ALPHAR, ALPHAI, BETA, VL, ldvl, VR, ldvr, ilo, ihi, lscale, rscale, abnrm, bbnrm, RCONDE, RCONDV, WORK, lwork, IWORK, BWORK, info);
806 }
807
808 void LAPACK<int, double>::GGEV(const char& JOBVL, const char& JOBVR, const int& n, double* A, const int& lda, double* B, const int& ldb, double* ALPHAR, double* ALPHAI, double* BETA, double* VL, const int& ldvl, double* VR, const int& ldvr, double* WORK, const int& lwork, int* info) const
809 {
810 DGGEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, B, &ldb, ALPHAR, ALPHAI, BETA, VL, &ldvl, VR, &ldvr, WORK, &lwork, info);
811 }
812
813 void LAPACK<int, double>::TRSEN(const char& JOB, const char& COMPQ, const int* SELECT, const int& n, double* T, const int& ldt, double* Q, const int& ldq, double* WR, double* WI, int* M, double* S, double* SEP, double* WORK, const int& lwork, int* IWORK, const int& liwork, int* info ) const
814 { DTRSEN_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPQ), SELECT, &n, T, &ldt, Q, &ldq, WR, WI, M, S, SEP, WORK, &lwork, IWORK, &liwork, info); }
815
816
817 void LAPACK<int, double>::TGSEN(const int& ijob, const int& wantq, const int& wantz, const int* SELECT, const int& n, double* A, const int& lda, double* B, const int& ldb, double* ALPHAR, double* ALPHAI, double* BETA, double* Q, const int& ldq, double* Z, const int& ldz, int* M, double* PL, double* PR, double* DIF, double* WORK, const int& lwork, int* IWORK, const int& liwork, int* info ) const
818 { DTGSEN_F77(&ijob, &wantq, &wantz, SELECT, &n, A, &lda, B, &ldb, ALPHAR, ALPHAI, BETA, Q, &ldq, Z, &ldz, M, PL, PR, DIF, WORK, &lwork, IWORK, &liwork, info); }
819
820
821 void LAPACK<int, double>::GGES(const char& JOBVL, const char& JOBVR, const char& SORT, int (*ptr2func)(double* , double* , double* ), const int& n, double* A, const int& lda, double* B, const int& ldb, int* sdim, double* ALPHAR, double* ALPHAI, double* BETA, double* VL, const int& ldvl, double* VR, const int& ldvr, double* WORK, const int& lwork, int* BWORK, int* info ) const
822 { DGGES_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), CHAR_MACRO(SORT), ptr2func, &n, A, &lda, B, &ldb, sdim, ALPHAR, ALPHAI, BETA, VL, &ldvl, VR, &ldvr, WORK, &lwork, BWORK, info); }
823
824
825 void LAPACK<int, double>::ORMQR(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& k, double* A, const int& lda, const double* TAU, double* C, const int& ldc, double* WORK, const int& lwork, int* info) const
826 {
827 DORMQR_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &k, A, &lda, TAU, C, &ldc, WORK, &lwork, info);
828 }
829
830 void LAPACK<int, double>::ORM2R(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& k, const double A[], const int& lda, const double TAU[], double C[], const int& ldc, double WORK[], int* const info) const
831 {
832 DORM2R_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &k, A, &lda, TAU, C, &ldc, WORK, info);
833 }
834
835 void LAPACK<int, double>::UNMQR(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& k, double* A, const int& lda, const double* TAU, double* C, const int& ldc, double* WORK, const int& lwork, int* info) const
836 {
837 // LAPACK only defines UNMQR for Z and C (complex*8
838 // resp. complex*16), but logically, UNMQR means the same thing as
839 // ORMQR for real arithmetic.
840 ORMQR (SIDE, TRANS, m, n, k, A, lda, TAU, C, ldc, WORK, lwork, info);
841 }
842
843 void LAPACK<int, double>::UNM2R (const char& SIDE, const char& TRANS, const int& M, const int& N, const int& K, const double A[], const int& LDA, const double TAU[], double C[], const int& LDC, double WORK[], int* const INFO) const
844 {
845 // LAPACK only defines UNM2R for Z and C (complex*8
846 // resp. complex*16), but logically, UNM2R means the same thing as
847 // ORM2R for real arithmetic.
848 ORM2R (SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, INFO);
849 }
850
851 void LAPACK<int, double>::ORGQR(const int& m, const int& n, const int& k, double* A, const int& lda, const double* TAU, double* WORK, const int& lwork, int* info) const
852 {
853 DORGQR_F77( &m, &n, &k, A, &lda, TAU, WORK, &lwork, info);
854 }
855
856
857 void LAPACK<int, double>::UNGQR(const int& m, const int& n, const int& k, double* A, const int& lda, const double* TAU, double* WORK, const int& lwork, int* info) const
858 {
859 DORGQR_F77( &m, &n, &k, A, &lda, TAU, WORK, &lwork, info);
860 }
861
862
863 void LAPACK<int, double>::ORGHR(const int& n, const int& ilo, const int& ihi, double* A, const int& lda, const double* TAU, double* WORK, const int& lwork, int* info) const
864 {
865 DORGHR_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
866 }
867
868
869 void LAPACK<int, double>::ORMHR(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& ilo, const int& ihi, const double* A, const int& lda, const double* TAU, double* C, const int& ldc, double* WORK, const int& lwork, int* info) const
870 {
871 DORMHR_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &ilo, &ihi, A, &lda, TAU, C, &ldc, WORK, &lwork, info);
872 }
873
874
875 void LAPACK<int, double>::TREVC(const char& SIDE, const char& HOWMNY, int* select, const int& n, const double* T, const int& ldt, double* VL, const int& ldvl, double* VR, const int& ldvr, const int& mm, int* m, double* WORK, int* info) const
876 {
877 DTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), select, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, info);
878 }
879
880
881 void LAPACK<int, double>::TREVC(const char& SIDE, const int& n, const double* T, const int& ldt, double* VL, const int& ldvl, double* VR, const int& ldvr, const int& mm, int* m, double* WORK, double* /* RWORK */, int* info) const
882 {
883 std::vector<int> select(1);
884 const char whch = 'A';
885 DTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(whch), &select[0], &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, info);
886 }
887
888 void LAPACK<int, double>::TREXC(const char& COMPQ, const int& n, double* T, const int& ldt, double* Q, const int& ldq, int* ifst, int* ilst, double* WORK, int* info) const
889 {
890 DTREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, ifst, ilst, WORK, info);
891 }
892
893
894 void LAPACK<int, double>::TGEVC(const char& SIDE, const char& HOWMNY, const int* SELECT, const int& n, double* S, const int& lds, double* P, const int& ldp, double* VL, const int& ldvl, double* VR, const int& ldvr, const int& mm, int* M, double* WORK, int* info) const
895 { DTGEVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), SELECT, &n, S, &lds, P, &ldp, VL, &ldvl, VR, &ldvr, &mm, M, WORK, info); }
896
897
898 void LAPACK<int, double>::LARTG( const double& f, const double& g, double* c, double* s, double* r ) const
899 {
900 DLARTG_F77(&f, &g, c, s, r);
901 }
902
903
904 void LAPACK<int, double>::LARFG( const int& n, double* alpha, double* x, const int& incx, double* tau ) const
905 {
906 DLARFG_F77(&n, alpha, x, &incx, tau);
907 }
908
909 void LAPACK<int, double>::GEBAL(const char& JOBZ, const int& n, double* A, const int& lda, int* ilo, int* ihi, double* scale, int* info) const
910 {
911 DGEBAL_F77(CHAR_MACRO(JOBZ),&n, A, &lda, ilo, ihi, scale, info);
912 }
913
914
915 void LAPACK<int, double>::GEBAK(const char& JOBZ, const char& SIDE, const int& n, const int& ilo, const int& ihi, const double* scale, const int& m, double* V, const int& ldv, int* info) const
916 {
917 DGEBAK_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(SIDE), &n, &ilo, &ihi, scale, &m, V, &ldv, info);
918 }
919
920
921#ifdef HAVE_TEUCHOS_LAPACKLARND
922 double LAPACK<int, double>::LARND( const int& idist, int* seed ) const
923 {
924 return(DLARND_F77(&idist, seed));
925 }
926#endif
927
928 void LAPACK<int, double>::LARNV( const int& idist, int* seed, const int& n, double* v ) const
929 {
930 DLARNV_F77(&idist, seed, &n, v);
931 }
932
933
934 double LAPACK<int, double>::LAMCH(const char& CMACH) const
935 {
936 return(DLAMCH_F77(CHAR_MACRO(CMACH)));
937 }
938
939
940 int LAPACK<int, double>::ILAENV( const int& ispec, const std::string& NAME, const std::string& OPTS, const int& N1, const int& N2, const int& N3, const int& N4 ) const
941 {
942 unsigned int opts_length = OPTS.length();
943 // if user queries a Hermitian routine, change it to a symmetric routine
944 std::string temp_NAME = "d" + NAME;
945 if (temp_NAME.substr(1,2) == "he") {
946 temp_NAME.replace(1,2,"sy");
947 }
948 unsigned int name_length = temp_NAME.length();
949 return ilaenv_wrapper(&ispec, &temp_NAME[0], name_length, &OPTS[0], opts_length, &N1, &N2, &N3, &N4);
950 }
951
952
953 double LAPACK<int, double>::LAPY2(const double& x, const double& y) const
954 {
955 return DLAPY2_F77(&x, &y);
956 }
957
958 // END INT, DOUBLE SPECIALIZATION IMPLEMENTATION //
959
960#ifdef HAVE_TEUCHOS_COMPLEX
961
962 // BEGIN INT, COMPLEX<FLOAT> SPECIALIZATION IMPLEMENTATION //
963
964
965 void LAPACK<int, std::complex<float> >::PTTRF(const int& n, std::complex<float>* d, std::complex<float>* e, int* info) const
966 {
967 CPTTRF_F77(&n,d,e,info);
968 }
969
970
971 void LAPACK<int, std::complex<float> >::PTTRS(const int& n, const int& nrhs, const std::complex<float>* d, const std::complex<float>* e, std::complex<float>* B, const int& ldb, int* info) const
972 {
973 CPTTRS_F77(&n,&nrhs,d,e,B,&ldb,info);
974 }
975
976
977 void LAPACK<int, std::complex<float> >::POTRF(const char& UPLO, const int& n, std::complex<float>* A, const int& lda, int* info) const
978 {
979 CPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
980 }
981
982
983 void LAPACK<int, std::complex<float> >::POTRS(const char& UPLO, const int& n, const int& nrhs, const std::complex<float>* A, const int& lda, std::complex<float>* B, const int& ldb, int* info) const
984 {
985 CPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
986 }
987
988
989 void LAPACK<int, std::complex<float> >::POTRI(const char& UPLO, const int& n, std::complex<float>* A, const int& lda, int* info) const
990 {
991 CPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
992 }
993
994
995 void LAPACK<int, std::complex<float> >::POCON(const char& UPLO, const int& n, const std::complex<float>* A, const int& lda, const float& anorm, float* rcond, std::complex<float>* WORK, float* RWORK, int* info) const
996 {
997 CPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
998 }
999
1000
1001 void LAPACK<int, std::complex<float> >::POSV(const char& UPLO, const int& n, const int& nrhs, std::complex<float>* A, const int& lda, std::complex<float>* B, const int& ldb, int* info) const
1002 {
1003 CPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
1004 }
1005
1006
1007 void LAPACK<int, std::complex<float> >::POEQU(const int& n, const std::complex<float>* A, const int& lda, float* S, float* scond, float* amax, int* info) const
1008 {
1009 CPOEQU_F77(&n, A, &lda, S, scond, amax, info);
1010 }
1011
1012
1013 void LAPACK<int, std::complex<float> >::PORFS(const char& UPLO, const int& n, const int& nrhs, std::complex<float>* A, const int& lda, const std::complex<float>* AF, const int& ldaf, const std::complex<float>* B, const int& ldb, std::complex<float>* X, const int& ldx, float* FERR, float* BERR, std::complex<float>* WORK, float* RWORK, int* info) const
1014 {
1015 CPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
1016 }
1017
1018 void LAPACK<int, std::complex<float> >::POSVX(const char& FACT, const char& UPLO, const int& n, const int& nrhs, std::complex<float>* A, const int& lda, std::complex<float>* AF, const int& ldaf, char* EQUED, float* S, std::complex<float>* B, const int& ldb, std::complex<float>* X, const int& ldx, float* rcond, float* FERR, float* BERR, std::complex<float>* WORK, float* RWORK, int* info) const
1019 {
1020 CPOSVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, CHARPTR_MACRO(EQUED), S, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, RWORK, info);
1021 }
1022
1023
1024 void LAPACK<int,std::complex<float> >::GELS(const char& TRANS, const int& m, const int& n, const int& nrhs, std::complex<float>* A, const int& lda, std::complex<float>* B, const int& ldb, std::complex<float>* WORK, const int& lwork, int* info) const
1025 {
1026 CGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, info);
1027 }
1028
1029 void LAPACK<int, std::complex<float> >::GELSS(const int& m, const int& n, const int& nrhs, std::complex<float>* A, const int& lda, std::complex<float>* B, const int& ldb, float* S, const float& rcond, int* rank, std::complex<float>* WORK, const int& lwork, float* rwork, int* info) const
1030 {
1031 CGELSS_F77(&m, &n, &nrhs, A, &lda, B, &ldb, S, &rcond, rank, WORK, &lwork, rwork, info);
1032 }
1033
1034 void LAPACK<int,std::complex<float> >::GEQRF( const int& m, const int& n, std::complex<float>* A, const int& lda, std::complex<float>* TAU, std::complex<float>* WORK, const int& lwork, int* info) const
1035 {
1036 CGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info);
1037 }
1038
1039 void LAPACK<int,std::complex<float> >::GEQR2 (const int& m, const int& n, std::complex<float> A[], const int& lda, std::complex<float> TAU[], std::complex<float> WORK[], int* const info) const
1040 {
1041 CGEQR2_F77(&m, &n, A, &lda, TAU, WORK, info);
1042 }
1043
1044 void LAPACK<int,std::complex<float> >::UNGQR(const int& m, const int& n, const int& k, std::complex<float>* A, const int& lda, const std::complex<float>* TAU, std::complex<float>* WORK, const int& lwork, int* info) const
1045 {
1046 CUNGQR_F77( &m, &n, &k, A, &lda, TAU, WORK, &lwork, info);
1047 }
1048
1049 void LAPACK<int,std::complex<float> >::UNMQR(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& k, std::complex<float>* A, const int& lda, const std::complex<float>* TAU, std::complex<float>* C, const int& ldc, std::complex<float>* WORK, const int& lwork, int* info) const
1050 {
1051 CUNMQR_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &k, A, &lda, TAU, C, &ldc, WORK, &lwork, info);
1052 }
1053
1054 void LAPACK<int,std::complex<float> >::UNM2R (const char& SIDE, const char& TRANS, const int& M, const int& N, const int& K, const std::complex<float> A[], const int& LDA, const std::complex<float> TAU[], std::complex<float> C[], const int& LDC, std::complex<float> WORK[], int* const INFO) const
1055 {
1056 CUNM2R_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &M, &N, &K, A, &LDA, TAU, C, &LDC, WORK, INFO);
1057 }
1058
1059 void LAPACK<int,std::complex<float> >::GETRF(const int& m, const int& n, std::complex<float>* A, const int& lda, int* IPIV, int* info) const
1060 {
1061 CGETRF_F77(&m, &n, A, &lda, IPIV, info);
1062 }
1063
1064 void LAPACK<int,std::complex<float> >::GETRS(const char& TRANS, const int& n, const int& nrhs, const std::complex<float>* A, const int& lda, const int* IPIV, std::complex<float>* B , const int& ldb, int* info) const
1065 {
1066 CGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info);
1067 }
1068
1069 void LAPACK<int,std::complex<float> >::LASCL(const char& TYPE, const int& kl, const int& ku, const float& cfrom, const float& cto, const int& m, const int& n, std::complex<float>* A, const int& lda, int* info) const
1070 { CLASCL_F77(CHAR_MACRO(TYPE), &kl, &ku, &cfrom, &cto, &m, &n, A, &lda, info); }
1071
1072 void LAPACK<int,std::complex<float> >::GEQP3(const int& m, const int& n, std::complex<float>* A, const int& lda, int* jpvt, std::complex<float>* TAU, std::complex<float>* WORK, const int& lwork, float* RWORK, int* info ) const
1073 {
1074 CGEQP3_F77(&m, &n, A, &lda, jpvt, TAU, WORK, &lwork, RWORK, info);
1075 }
1076
1077 void LAPACK<int, std::complex<float> >::
1078 LASWP (const int& N,
1079 std::complex<float> A[],
1080 const int& LDA,
1081 const int& K1,
1082 const int& K2,
1083 const int IPIV[],
1084 const int& INCX) const
1085 {
1086 CLASWP_F77(&N, A, &LDA, &K1, &K2, IPIV, &INCX);
1087 }
1088
1089 void LAPACK<int,std::complex<float> >::GBTRF(const int& m, const int& n, const int& kl, const int& ku, std::complex<float>* A, const int& lda, int* IPIV, int* info) const
1090 {
1091 CGBTRF_F77(&m, &kl, &ku, &n, A, &lda, IPIV, info);
1092 }
1093
1094
1095 void LAPACK<int,std::complex<float> >::GBTRS(const char& TRANS, const int& n, const int& kl, const int& ku, const int& nrhs, const std::complex<float>* A, const int& lda, const int* IPIV, std::complex<float>* B , const int& ldb, int* info) const
1096 {
1097 CGBTRS_F77(CHAR_MACRO(TRANS), &n, &kl, &ku, &nrhs, A, &lda, IPIV, B, &ldb, info);
1098 }
1099
1100
1101 void LAPACK<int,std::complex<float> >::GTTRF(const int& n, std::complex<float>* dl, std::complex<float>* d, std::complex<float>* du, std::complex<float>* du2, int* IPIV, int* info) const
1102 {
1103 CGTTRF_F77(&n, dl, d, du, du2, IPIV, info);
1104 }
1105
1106
1107 void LAPACK<int,std::complex<float> >::GTTRS(const char& TRANS, const int& n, const int& nrhs, const std::complex<float>* dl, const std::complex<float>* d, const std::complex<float>* du, const std::complex<float>* du2, const int* IPIV, std::complex<float>* B, const int& ldb, int* info) const
1108 {
1109 CGTTRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, dl, d, du, du2, IPIV, B, &ldb, info);
1110 }
1111
1112
1113 void LAPACK<int,std::complex<float> >::GETRI(const int& n, std::complex<float>* A, const int& lda, const int* IPIV, std::complex<float>* WORK, const int& lwork, int* info) const
1114 {
1115 CGETRI_F77(&n, A, &lda, IPIV, WORK, &lwork, info);
1116 }
1117
1118
1119 void
1120 LAPACK<int, std::complex<float> >::LATRS (const char& UPLO,
1121 const char& TRANS,
1122 const char& DIAG,
1123 const char& NORMIN,
1124 const int& N,
1125 std::complex<float>* A,
1126 const int& LDA,
1127 std::complex<float>* X,
1128 float* SCALE,
1129 float* CNORM,
1130 int* INFO) const
1131 {
1132 CLATRS_F77(CHAR_MACRO(UPLO),
1133 CHAR_MACRO(TRANS),
1134 CHAR_MACRO(DIAG),
1135 CHAR_MACRO(NORMIN),
1136 &N,
1137 A,
1138 &LDA,
1139 X,
1140 SCALE,
1141 CNORM,
1142 INFO);
1143 }
1144
1145
1146 void LAPACK<int,std::complex<float> >::GECON(const char& NORM, const int& n, const std::complex<float>* A, const int& lda, const float& anorm, float* rcond, std::complex<float>* WORK, float* RWORK, int* info) const
1147 {
1148 CGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
1149 }
1150
1151
1152 void LAPACK<int,std::complex<float> >::GBCON(const char& NORM, const int& n, const int& kl, const int& ku, const std::complex<float>* A, const int& lda, int* IPIV, const float& anorm, float* rcond, std::complex<float>* WORK, float* RWORK, int* info) const
1153 {
1154 CGBCON_F77(CHAR_MACRO(NORM), &n, &kl, &ku, A, &lda, IPIV, &anorm, rcond, WORK, RWORK, info);
1155 }
1156
1157
1158 float LAPACK<int,std::complex<float> >::LANGB(const char& NORM, const int& n, const int& kl, const int& ku, const std::complex<float>* A, const int& lda, float* WORK) const
1159 {
1160 return( CLANGB_F77(CHAR_MACRO(NORM), &n, &kl, &ku, A, &lda, WORK) );
1161 }
1162
1163
1164 void LAPACK<int,std::complex<float> >::GESV(const int& n, const int& nrhs, std::complex<float>* A, const int& lda, int* IPIV, std::complex<float>* B, const int& ldb, int* info) const
1165 {
1166 CGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info);
1167 }
1168
1169
1170 void LAPACK<int,std::complex<float> >::GEEQU(const int& m, const int& n, const std::complex<float>* A, const int& lda, float* R, float* C, float* rowcond, float* colcond, float* amax, int* info) const
1171 {
1172 CGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info);
1173 }
1174
1175
1176 void LAPACK<int,std::complex<float> >::GERFS(const char& TRANS, const int& n, const int& nrhs, const std::complex<float>* A, const int& lda, const std::complex<float>* AF, const int& ldaf, const int* IPIV, const std::complex<float>* B, const int& ldb, std::complex<float>* X, const int& ldx, float* FERR, float* BERR, std::complex<float>* WORK, float* RWORK, int* info) const
1177 {
1178 CGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
1179 }
1180
1181
1182 void LAPACK<int,std::complex<float> >::GBEQU(const int& m, const int& n, const int& kl, const int& ku, const std::complex<float>* A, const int& lda, float* R, float* C, float* rowcond, float* colcond, float* amax, int* info) const
1183 {
1184 CGBEQU_F77(&m, &n, &kl, &ku, A, &lda, R, C, rowcond, colcond, amax, info);
1185 }
1186
1187
1188 void LAPACK<int,std::complex<float> >::GBRFS(const char& TRANS, const int& n, const int& kl, const int& ku, const int& nrhs, const std::complex<float>* A, const int& lda, const std::complex<float>* AF, const int& ldaf, const int* IPIV, const std::complex<float>* B, const int& ldb, std::complex<float>* X, const int& ldx, float* FERR, float* BERR, std::complex<float>* WORK, float* RWORK, int* info) const
1189 {
1190 CGBRFS_F77(CHAR_MACRO(TRANS), &n, &kl, &ku, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
1191 }
1192
1193 void LAPACK<int,std::complex<float> >::GESVX(const char& FACT, const char& TRANS, const int& n, const int& nrhs, std::complex<float>* A, const int& lda, std::complex<float>* AF, const int& ldaf, int* IPIV, char* EQUED, float* R, float* C, std::complex<float>* B, const int& ldb, std::complex<float>* X, const int& ldx, float* rcond, float* FERR, float* BERR, std::complex<float>* WORK, float* RWORK, int* info) const
1194 {
1195 CGESVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, CHARPTR_MACRO(EQUED), R, C, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, RWORK, info);
1196 }
1197
1198
1199 void LAPACK<int,std::complex<float> >::GEHRD(const int& n, const int& ilo, const int& ihi, std::complex<float>* A, const int& lda, std::complex<float>* TAU, std::complex<float>* WORK, const int& lwork, int* info) const
1200 {
1201 CGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
1202 }
1203
1204
1205 void LAPACK<int,std::complex<float> >::TRTRS(const char& UPLO, const char& TRANS, const char& DIAG, const int& n, const int& nrhs, const std::complex<float>* A, const int& lda, std::complex<float>* B, const int& ldb, int* info) const
1206 {
1207 CTRTRS_F77(CHAR_MACRO(UPLO), CHAR_MACRO(TRANS), CHAR_MACRO(DIAG), &n, &nrhs, A, &lda, B, &ldb, info);
1208 }
1209
1210
1211 void LAPACK<int,std::complex<float> >::TRTRI(const char& UPLO, const char& DIAG, const int& n, const std::complex<float>* A, const int& lda, int* info) const
1212 {
1213 CTRTRI_F77(CHAR_MACRO(UPLO), CHAR_MACRO(DIAG), &n, A, &lda, info);
1214 }
1215
1216
1217 void LAPACK<int,std::complex<float> >::STEQR(const char& COMPZ, const int& n, float* D, float* E, std::complex<float>* Z, const int& ldz, float* WORK, int* info) const
1218 {
1219 CSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
1220 }
1221
1222
1223 void LAPACK<int,std::complex<float> >::PTEQR(const char& COMPZ, const int& n, float* D, float* E, std::complex<float>* Z, const int& ldz, float* WORK, int* info) const
1224 {
1225 CPTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
1226 }
1227
1228
1229 void LAPACK<int,std::complex<float> >::HEEV(const char& JOBZ, const char& UPLO, const int& n, std::complex<float> * A, const int& lda, float* W, std::complex<float> * WORK, const int& lwork, float* RWORK, int* info) const
1230 {
1231 CHEEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, W, WORK, &lwork, RWORK, info);
1232 }
1233
1234
1235 void LAPACK<int,std::complex<float> >::HEGV(const int& itype, const char& JOBZ, const char& UPLO, const int& n, std::complex<float> * A, const int& lda, std::complex<float> * B, const int& ldb, float* W, std::complex<float> * WORK, const int& lwork, float* RWORK, int* info) const
1236 {
1237 CHEGV_F77(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, B, &ldb, W, WORK, &lwork, RWORK, info);
1238 }
1239
1240
1241 void LAPACK<int, std::complex<float> >::HSEQR(const char& JOB, const char& COMPZ, const int& n, const int& ilo, const int& ihi, std::complex<float>* H, const int& ldh, std::complex<float>* W, std::complex<float>* Z, const int& ldz, std::complex<float>* WORK, const int& lwork, int* info) const
1242 {
1243 CHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, W, Z, &ldz, WORK, &lwork, info);
1244 }
1245
1246
1247 void LAPACK<int, std::complex<float> >::GEES(const char& JOBVS, const char& SORT, int (*ptr2func)(std::complex<float>*), const int& n, std::complex<float>* A, const int& lda, int* sdim, std::complex<float>* W, std::complex<float>* VS, const int& ldvs, std::complex<float>* WORK, const int& lwork, float* RWORK, int* BWORK, int* info) const
1248 {
1249 CGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), ptr2func, &n, A, &lda, sdim, W, VS, &ldvs, WORK, &lwork, RWORK, BWORK, info);
1250 }
1251
1252
1253 void LAPACK<int, std::complex<float> >::GEES(const char& JOBVS, const int& n, std::complex<float>* A, const int& lda, int* sdim, float* WR, float* WI, std::complex<float>* VS, const int& ldvs, std::complex<float>* WORK, const int& lwork, float* RWORK, int* BWORK, int* info) const
1254 {
1255 int (*nullfptr)(std::complex<float>*) = NULL;
1256 std::vector< std::complex<float> > W(n);
1257 const char sort = 'N';
1258 CGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(sort), nullfptr, &n, A, &lda, sdim, &W[0], VS, &ldvs, WORK, &lwork, RWORK, BWORK, info);
1259 for (int i=0; i<n; i++) {
1260 WR[i] = W[i].real();
1261 WI[i] = W[i].imag();
1262 }
1263 }
1264
1265
1266 void LAPACK<int, std::complex<float> >::GEEV(const char& JOBVL, const char& JOBVR, const int& n, std::complex<float>* A, const int& lda, std::complex<float>* W, std::complex<float>* VL, const int& ldvl, std::complex<float>* VR, const int& ldvr, std::complex<float>* WORK, const int& lwork, float* RWORK, int* info) const
1267 {
1268 CGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, W, VL, &ldvl, VR, &ldvr, WORK, &lwork, RWORK, info);
1269 }
1270
1271 void LAPACK<int, std::complex<float> >::GEEV(const char& JOBVL, const char& JOBVR, const int& n, std::complex<float>* A, const int& lda, float* WR, float* WI, std::complex<float>* VL, const int& ldvl, std::complex<float>* VR, const int& ldvr, std::complex<float>* WORK, const int& lwork, float* RWORK, int* info) const
1272 {
1273 std::vector<std::complex<float> > w (n);
1274 std::complex<float>* w_rawPtr = (n == 0) ? NULL : &w[0];
1275 GEEV (JOBVL, JOBVR, n, A, lda, w_rawPtr, VL, ldvl, VR, ldvr, WORK, lwork, RWORK, info);
1276 if (*info == 0) {
1277 // The eigenvalues are only valid on output if INFO is zero.
1278 // Otherwise, we shouldn't even write to WR or WI.
1279 for (int k = 0; k < n; ++k) {
1280 WR[k] = w[k].real ();
1281 WI[k] = w[k].imag ();
1282 }
1283 }
1284 }
1285
1286 void LAPACK<int, std::complex<float> >::GESVD(const char& JOBU, const char& JOBVT, const int& m, const int& n, std::complex<float> * A, const int& lda, float* S, std::complex<float> * U, const int& ldu, std::complex<float> * V, const int& ldv, std::complex<float> * WORK, const int& lwork, float* RWORK, int* info) const {
1287 CGESVD_F77(CHAR_MACRO(JOBU), CHAR_MACRO(JOBVT), &m, &n, A, &lda, S, U, &ldu, V, &ldv, WORK, &lwork, RWORK, info);
1288 }
1289
1290
1291 void LAPACK<int, std::complex<float> >::GEEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, std::complex<float>* A, const int& lda, std::complex<float>* W, std::complex<float>* VL, const int& ldvl, std::complex<float>* VR, const int& ldvr, int* ilo, int* ihi, float* SCALE, float* abnrm, float* RCONDE, float* RCONDV, std::complex<float>* WORK, const int& lwork, float* RWORK, int* info) const
1292 {
1293 CGEEVX_F77(CHAR_MACRO(BALANC), CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), CHAR_MACRO(SENSE), &n, A, &lda, W, VL, &ldvl, VR, &ldvr, ilo, ihi, SCALE, abnrm, RCONDE, RCONDV, WORK, &lwork, RWORK, info);
1294 }
1295
1296
1297 void LAPACK<int, std::complex<float> >::GGEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, std::complex<float>* A, const int& lda, std::complex<float>* B, const int& ldb, std::complex<float>* ALPHA, std::complex<float>* BETA, std::complex<float>* VL, const int& ldvl, std::complex<float>* VR, const int& ldvr, int* ilo, int* ihi, float* lscale, float* rscale, float* abnrm, float* bbnrm, float* RCONDE, float* RCONDV, std::complex<float>* WORK, const int& lwork, float* RWORK, int* IWORK, int* BWORK, int* info) const
1298 {
1299 CGGEVX_F77(CHAR_MACRO(BALANC), CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), CHAR_MACRO(SENSE), &n, A, &lda, B, &ldb, ALPHA, BETA, VL, &ldvl, VR, &ldvr, ilo, ihi, lscale, rscale, abnrm, bbnrm, RCONDE, RCONDV, WORK, &lwork, RWORK, IWORK, BWORK, info);
1300 }
1301
1302 void LAPACK<int, std::complex<float> >::GGEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, std::complex<float>* A, const int& lda, std::complex<float>* B, const int& ldb, float* ALPHAR, float* ALPHAI, std::complex<float>* BETA, std::complex<float>* VL, const int& ldvl, std::complex<float>* VR, const int& ldvr, int* ilo, int* ihi, float* lscale, float* rscale, float* abnrm, float* bbnrm, float* RCONDE, float* RCONDV, std::complex<float>* WORK, const int& lwork, float* RWORK, int* IWORK, int* BWORK, int* info) const
1303 {
1304 std::vector<std::complex<float> > w (n);
1305 std::complex<float>* w_rawPtr = (n == 0) ? NULL : &w[0];
1306 GGEVX(BALANC, JOBVL, JOBVR, SENSE, n, A, lda, B, ldb, w_rawPtr, BETA, VL, ldvl, VR, ldvr, ilo, ihi, lscale, rscale, abnrm, bbnrm, RCONDE, RCONDV, WORK, lwork, RWORK, IWORK, BWORK, info);
1307 if (*info == 0) {
1308 // The eigenvalues are only valid on output if INFO is zero.
1309 // Otherwise, we shouldn't even write to WR or WI.
1310 for (int k = 0; k < n; ++k) {
1311 ALPHAR[k] = w[k].real ();
1312 ALPHAI[k] = w[k].imag ();
1313 }
1314 }
1315 }
1316
1317
1318 void LAPACK<int, std::complex<float> >::TREVC(const char& SIDE, const char& HOWMNY, int* select, const int& n, const std::complex<float>* T, const int& ldt, std::complex<float>* VL, const int& ldvl, std::complex<float>* VR, const int& ldvr, const int& mm, int* m, std::complex<float>* WORK, float* RWORK, int* info) const
1319 {
1320 CTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), select, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, RWORK, info);
1321 }
1322
1323
1324 void LAPACK<int, std::complex<float> >::TREVC(const char& SIDE, const int& n, const std::complex<float>* T, const int& ldt, std::complex<float>* VL, const int& ldvl, std::complex<float>* VR, const int& ldvr, const int& mm, int* m, std::complex<float>* WORK, float* RWORK, int* info) const
1325 {
1326 std::vector<int> select(1);
1327 const char& whch = 'A';
1328 CTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(whch), &select[0], &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, RWORK, info);
1329 }
1330
1331 void LAPACK<int, std::complex<float> >::TREXC(const char& COMPQ, const int& n, std::complex<float>* T, const int& ldt, std::complex<float>* Q, const int& ldq, int* ifst, int* ilst, std::complex<float>* WORK, int* info) const
1332 {
1333 CTREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, ifst, ilst, info);
1334 }
1335
1336
1337 void LAPACK<int, std::complex<float> >::LARTG( const std::complex<float> f, const std::complex<float> g, float* c, std::complex<float>* s, std::complex<float>* r ) const
1338 {
1339 CLARTG_F77(&f, &g, c, s, r);
1340 }
1341
1342
1343 void LAPACK<int, std::complex<float> >::LARFG( const int& n, std::complex<float>* alpha, std::complex<float>* x, const int& incx, std::complex<float>* tau ) const
1344 {
1345 CLARFG_F77(&n, alpha, x, &incx, tau);
1346 }
1347
1348 void LAPACK<int, std::complex<float> >::GEBAL(const char& JOBZ, const int& n, std::complex<float>* A, const int& lda, int* ilo, int* ihi, float* scale, int* info) const
1349 {
1350 CGEBAL_F77(CHAR_MACRO(JOBZ),&n, A, &lda, ilo, ihi, scale, info);
1351 }
1352
1353
1354 void LAPACK<int, std::complex<float> >::GEBAK(const char& JOBZ, const char& SIDE, const int& n, const int& ilo, const int& ihi, const float* scale, const int& m, std::complex<float>* V, const int& ldv, int* info) const
1355 {
1356 CGEBAK_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(SIDE), &n, &ilo, &ihi, scale, &m, V, &ldv, info);
1357 }
1358
1359
1360#ifdef HAVE_TEUCHOS_LAPACKLARND
1361 std::complex<float> LAPACK<int, std::complex<float> >::LARND( const int& idist, int* seed ) const
1362 {
1363 return(CLARND_F77(&idist, seed));
1364 }
1365#endif
1366
1367 void LAPACK<int, std::complex<float> >::LARNV( const int& idist, int* seed, const int& n, std::complex<float>* v ) const
1368 {
1369 CLARNV_F77(&idist, seed, &n, v);
1370 }
1371
1372
1373 int LAPACK<int, std::complex<float> >::ILAENV( const int& ispec, const std::string& NAME, const std::string& OPTS, const int& N1, const int& N2, const int& N3, const int& N4 ) const
1374 {
1375 unsigned int opts_length = OPTS.length();
1376 std::string temp_NAME = "c" + NAME;
1377 unsigned int name_length = temp_NAME.length();
1378 return ilaenv_wrapper(&ispec, &temp_NAME[0], name_length, &OPTS[0], opts_length, &N1, &N2, &N3, &N4);
1379 }
1380
1381 // END INT, COMPLEX<FLOAT> SPECIALIZATION IMPLEMENTATION //
1382
1383 // BEGIN INT, COMPLEX<DOUBLE> SPECIALIZATION IMPLEMENTATION //
1384
1385
1386 void LAPACK<int, std::complex<double> >::PTTRF(const int& n, std::complex<double>* d, std::complex<double>* e, int* info) const
1387 {
1388 ZPTTRF_F77(&n,d,e,info);
1389 }
1390
1391
1392 void LAPACK<int, std::complex<double> >::PTTRS(const int& n, const int& nrhs, const std::complex<double>* d, const std::complex<double>* e, std::complex<double>* B, const int& ldb, int* info) const
1393 {
1394 ZPTTRS_F77(&n,&nrhs,d,e,B,&ldb,info);
1395 }
1396
1397
1398 void LAPACK<int, std::complex<double> >::POTRF(const char& UPLO, const int& n, std::complex<double>* A, const int& lda, int* info) const
1399 {
1400 ZPOTRF_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
1401 }
1402
1403
1404 void LAPACK<int, std::complex<double> >::POTRS(const char& UPLO, const int& n, const int& nrhs, const std::complex<double>* A, const int& lda, std::complex<double>* B, const int& ldb, int* info) const
1405 {
1406 ZPOTRS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
1407 }
1408
1409
1410 void LAPACK<int, std::complex<double> >::POTRI(const char& UPLO, const int& n, std::complex<double>* A, const int& lda, int* info) const
1411 {
1412 ZPOTRI_F77(CHAR_MACRO(UPLO), &n, A, &lda, info);
1413 }
1414
1415
1416 void LAPACK<int, std::complex<double> >::POCON(const char& UPLO, const int& n, const std::complex<double>* A, const int& lda, const double& anorm, double* rcond, std::complex<double>* WORK, double* RWORK, int* info) const
1417 {
1418 ZPOCON_F77(CHAR_MACRO(UPLO), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
1419 }
1420
1421
1422 void LAPACK<int, std::complex<double> >::POSV(const char& UPLO, const int& n, const int& nrhs, std::complex<double>* A, const int& lda, std::complex<double>* B, const int& ldb, int* info) const
1423 {
1424 ZPOSV_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, B, &ldb, info);
1425 }
1426
1427
1428 void LAPACK<int, std::complex<double> >::POEQU(const int& n, const std::complex<double>* A, const int& lda, double* S, double* scond, double* amax, int* info) const
1429 {
1430 ZPOEQU_F77(&n, A, &lda, S, scond, amax, info);
1431 }
1432
1433
1434 void LAPACK<int, std::complex<double> >::PORFS(const char& UPLO, const int& n, const int& nrhs, std::complex<double>* A, const int& lda, const std::complex<double>* AF, const int& ldaf, const std::complex<double>* B, const int& ldb, std::complex<double>* X, const int& ldx, double* FERR, double* BERR, std::complex<double>* WORK, double* RWORK, int* info) const
1435 {
1436 ZPORFS_F77(CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
1437 }
1438
1439 void LAPACK<int, std::complex<double> >::POSVX(const char& FACT, const char& UPLO, const int& n, const int& nrhs, std::complex<double>* A, const int& lda, std::complex<double>* AF, const int& ldaf, char* EQUED, double* S, std::complex<double>* B, const int& ldb, std::complex<double>* X, const int& ldx, double* rcond, double* FERR, double* BERR, std::complex<double>* WORK, double* RWORK, int* info) const
1440 {
1441 ZPOSVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(UPLO), &n, &nrhs, A, &lda, AF, &ldaf, CHARPTR_MACRO(EQUED), S, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, RWORK, info);
1442 }
1443
1444
1445 void LAPACK<int,std::complex<double> >::GELS(const char& TRANS, const int& m, const int& n, const int& nrhs, std::complex<double>* A, const int& lda, std::complex<double>* B, const int& ldb, std::complex<double>* WORK, const int& lwork, int* info) const
1446 {
1447 ZGELS_F77(CHAR_MACRO(TRANS), &m, &n, &nrhs, A, &lda, B, &ldb, WORK, &lwork, info);
1448 }
1449
1450
1451 void LAPACK<int, std::complex<double> >::GELSS(const int& m, const int& n, const int& nrhs, std::complex<double>* A, const int& lda, std::complex<double>* B, const int& ldb, double* S, const double& rcond, int* rank, std::complex<double>* WORK, const int& lwork, double* rwork, int* info) const
1452 {
1453 ZGELSS_F77(&m, &n, &nrhs, A, &lda, B, &ldb, S, &rcond, rank, WORK, &lwork, rwork, info);
1454 }
1455
1456
1457 void LAPACK<int,std::complex<double> >::GEQRF( const int& m, const int& n, std::complex<double>* A, const int& lda, std::complex<double>* TAU, std::complex<double>* WORK, const int& lwork, int* info) const
1458 {
1459 ZGEQRF_F77(&m, &n, A, &lda, TAU, WORK, &lwork, info);
1460 }
1461
1462 void LAPACK<int,std::complex<double> >::GEQR2 (const int& m, const int& n, std::complex<double> A[], const int& lda, std::complex<double> TAU[], std::complex<double> WORK[], int* const info) const
1463 {
1464 ZGEQR2_F77(&m, &n, A, &lda, TAU, WORK, info);
1465 }
1466
1467 void LAPACK<int,std::complex<double> >::UNGQR(const int& m, const int& n, const int& k, std::complex<double>* A, const int& lda, const std::complex<double>* TAU, std::complex<double>* WORK, const int& lwork, int* info) const
1468 {
1469 ZUNGQR_F77( &m, &n, &k, A, &lda, TAU, WORK, &lwork, info);
1470 }
1471
1472
1473 void LAPACK<int,std::complex<double> >::UNMQR(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& k, std::complex<double>* A, const int& lda, const std::complex<double>* TAU, std::complex<double>* C, const int& ldc, std::complex<double>* WORK, const int& lwork, int* info) const
1474 {
1475 ZUNMQR_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &m, &n, &k, A, &lda, TAU, C, &ldc, WORK, &lwork, info);
1476 }
1477
1478 void LAPACK<int,std::complex<double> >::UNM2R (const char& SIDE, const char& TRANS, const int& M, const int& N, const int& K, const std::complex<double> A[], const int& LDA, const std::complex<double> TAU[], std::complex<double> C[], const int& LDC, std::complex<double> WORK[], int* const INFO) const
1479 {
1480 ZUNM2R_F77(CHAR_MACRO(SIDE), CHAR_MACRO(TRANS), &M, &N, &K, A, &LDA, TAU, C, &LDC, WORK, INFO);
1481 }
1482
1483 void LAPACK<int,std::complex<double> >::GETRF(const int& m, const int& n, std::complex<double>* A, const int& lda, int* IPIV, int* info) const
1484 {
1485 ZGETRF_F77(&m, &n, A, &lda, IPIV, info);
1486 }
1487
1488
1489 void LAPACK<int,std::complex<double> >::GETRS(const char& TRANS, const int& n, const int& nrhs, const std::complex<double>* A, const int& lda, const int* IPIV, std::complex<double>* B, const int& ldb, int* info) const
1490 {
1491 ZGETRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, IPIV, B, &ldb, info);
1492 }
1493
1494
1495 void LAPACK<int,std::complex<double> >::LASCL(const char& TYPE, const int& kl, const int& ku, const double& cfrom, const double& cto, const int& m, const int& n, std::complex<double>* A, const int& lda, int* info) const
1496 { ZLASCL_F77(CHAR_MACRO(TYPE), &kl, &ku, &cfrom, &cto, &m, &n, A, &lda, info); }
1497
1498 void LAPACK<int,std::complex<double> >::GEQP3(const int& m, const int& n, std::complex<double>* A, const int& lda, int* jpvt, std::complex<double>* TAU, std::complex<double>* WORK, const int& lwork, double* RWORK, int* info ) const
1499 {
1500 ZGEQP3_F77(&m, &n, A, &lda, jpvt, TAU, WORK, &lwork, RWORK, info);
1501 }
1502
1503 void LAPACK<int, std::complex<double> >::
1504 LASWP (const int& N,
1505 std::complex<double> A[],
1506 const int& LDA,
1507 const int& K1,
1508 const int& K2,
1509 const int IPIV[],
1510 const int& INCX) const
1511 {
1512 ZLASWP_F77(&N, A, &LDA, &K1, &K2, IPIV, &INCX);
1513 }
1514
1515 void LAPACK<int,std::complex<double> >::GBTRF(const int& m, const int& n, const int& kl, const int& ku, std::complex<double>* A, const int& lda, int* IPIV, int* info) const
1516 {
1517 ZGBTRF_F77(&m, &n, &kl, &ku, A, &lda, IPIV, info);
1518 }
1519
1520
1521 void LAPACK<int,std::complex<double> >::GBTRS(const char& TRANS, const int& n, const int& kl, const int& ku, const int& nrhs, const std::complex<double>* A, const int& lda, const int* IPIV, std::complex<double>* B, const int& ldb, int* info) const
1522 {
1523 ZGBTRS_F77(CHAR_MACRO(TRANS), &n, &kl, &ku, &nrhs, A, &lda, IPIV, B, &ldb, info);
1524 }
1525
1526
1527 void LAPACK<int,std::complex<double> >::GTTRF(const int& n, std::complex<double>* dl, std::complex<double>* d, std::complex<double>* du, std::complex<double>* du2, int* IPIV, int* info) const
1528 {
1529 ZGTTRF_F77(&n, dl, d, du, du2, IPIV, info);
1530 }
1531
1532
1533 void LAPACK<int,std::complex<double> >::GTTRS(const char& TRANS, const int& n, const int& nrhs, const std::complex<double>* dl, const std::complex<double>* d, const std::complex<double>* du, const std::complex<double>* du2, const int* IPIV, std::complex<double>* B, const int& ldb, int* info) const
1534 {
1535 ZGTTRS_F77(CHAR_MACRO(TRANS), &n, &nrhs, dl, d, du, du2, IPIV, B, &ldb, info);
1536 }
1537
1538
1539 void LAPACK<int,std::complex<double> >::GETRI(const int& n, std::complex<double>* A, const int& lda, const int* IPIV, std::complex<double>* WORK, const int& lwork, int* info) const
1540 {
1541 ZGETRI_F77(&n, A, &lda, IPIV, WORK, &lwork, info);
1542 }
1543
1544 void
1545 LAPACK<int, std::complex<double> >::LATRS (const char& UPLO,
1546 const char& TRANS,
1547 const char& DIAG,
1548 const char& NORMIN,
1549 const int& N,
1550 std::complex<double>* A,
1551 const int& LDA,
1552 std::complex<double>* X,
1553 double* SCALE,
1554 double* CNORM,
1555 int* INFO) const
1556 {
1557 ZLATRS_F77(CHAR_MACRO(UPLO),
1558 CHAR_MACRO(TRANS),
1559 CHAR_MACRO(DIAG),
1560 CHAR_MACRO(NORMIN),
1561 &N,
1562 A,
1563 &LDA,
1564 X,
1565 SCALE,
1566 CNORM,
1567 INFO);
1568 }
1569
1570 void LAPACK<int,std::complex<double> >::GECON(const char& NORM, const int& n, const std::complex<double>* A, const int& lda, const double& anorm, double* rcond, std::complex<double>* WORK, double* RWORK, int* info) const
1571 {
1572 ZGECON_F77(CHAR_MACRO(NORM), &n, A, &lda, &anorm, rcond, WORK, RWORK, info);
1573 }
1574
1575
1576 void LAPACK<int,std::complex<double> >::GBCON(const char& NORM, const int& n, const int& kl, const int& ku, const std::complex<double>* A, const int& lda, int* IPIV, const double& anorm, double* rcond, std::complex<double>* WORK, double* RWORK, int* info) const
1577 {
1578 ZGBCON_F77(CHAR_MACRO(NORM), &n, &kl, &ku, A, &lda, IPIV, &anorm, rcond, WORK, RWORK, info);
1579 }
1580
1581
1582 double LAPACK<int,std::complex<double> >::LANGB(const char& NORM, const int& n, const int& kl, const int& ku, const std::complex<double>* A, const int& lda, double* WORK) const
1583 {
1584 return( ZLANGB_F77(CHAR_MACRO(NORM), &n, &kl, &ku, A, &lda, WORK) );
1585 }
1586
1587
1588 void LAPACK<int,std::complex<double> >::GESV(const int& n, const int& nrhs, std::complex<double>* A, const int& lda, int* IPIV, std::complex<double>* B, const int& ldb, int* info) const
1589 {
1590 ZGESV_F77(&n, &nrhs, A, &lda, IPIV, B, &ldb, info);
1591 }
1592
1593
1594 void LAPACK<int,std::complex<double> >::GEEQU(const int& m, const int& n, const std::complex<double>* A, const int& lda, double* R, double* C, double* rowcond, double* colcond, double* amax, int* info) const
1595 {
1596 ZGEEQU_F77(&m, &n, A, &lda, R, C, rowcond, colcond, amax, info);
1597 }
1598
1599
1600 void LAPACK<int,std::complex<double> >::GERFS(const char& TRANS, const int& n, const int& nrhs, const std::complex<double>* A, const int& lda, const std::complex<double>* AF, const int& ldaf, const int* IPIV, const std::complex<double>* B, const int& ldb, std::complex<double>* X, const int& ldx, double* FERR, double* BERR, std::complex<double>* WORK, double* RWORK, int* info) const
1601 {
1602 ZGERFS_F77(CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
1603 }
1604
1605
1606 void LAPACK<int,std::complex<double> >::GBEQU(const int& m, const int& n, const int& kl, const int& ku, const std::complex<double>* A, const int& lda, double* R, double* C, double* rowcond, double* colcond, double* amax, int* info) const
1607 {
1608 ZGBEQU_F77(&m, &n, &kl, &ku, A, &lda, R, C, rowcond, colcond, amax, info);
1609 }
1610
1611
1612 void LAPACK<int,std::complex<double> >::GBRFS(const char& TRANS, const int& n, const int& kl, const int& ku, const int& nrhs, const std::complex<double>* A, const int& lda, const std::complex<double>* AF, const int& ldaf, const int* IPIV, const std::complex<double>* B, const int& ldb, std::complex<double>* X, const int& ldx, double* FERR, double* BERR, std::complex<double>* WORK, double* RWORK, int* info) const
1613 {
1614 ZGBRFS_F77(CHAR_MACRO(TRANS), &n, &kl, &ku, &nrhs, A, &lda, AF, &ldaf, IPIV, B, &ldb, X, &ldx, FERR, BERR, WORK, RWORK, info);
1615 }
1616
1617 void LAPACK<int,std::complex<double> >::GESVX(const char& FACT, const char& TRANS, const int& n, const int& nrhs, std::complex<double>* A, const int& lda, std::complex<double>* AF, const int& ldaf, int* IPIV, char* EQUED, double* R, double* C, std::complex<double>* B, const int& ldb, std::complex<double>* X, const int& ldx, double* rcond, double* FERR, double* BERR, std::complex<double>* WORK, double* RWORK, int* info) const
1618 {
1619 ZGESVX_F77(CHAR_MACRO(FACT), CHAR_MACRO(TRANS), &n, &nrhs, A, &lda, AF, &ldaf, IPIV, CHARPTR_MACRO(EQUED), R, C, B, &ldb, X, &ldx, rcond, FERR, BERR, WORK, RWORK, info);
1620 }
1621
1622
1623 void LAPACK<int,std::complex<double> >::GEHRD(const int& n, const int& ilo, const int& ihi, std::complex<double>* A, const int& lda, std::complex<double>* TAU, std::complex<double>* WORK, const int& lwork, int* info) const
1624 {
1625 ZGEHRD_F77(&n, &ilo, &ihi, A, &lda, TAU, WORK, &lwork, info);
1626 }
1627
1628
1629 void LAPACK<int,std::complex<double> >::TRTRS(const char& UPLO, const char& TRANS, const char& DIAG, const int& n, const int& nrhs, const std::complex<double>* A, const int& lda, std::complex<double>* B, const int& ldb, int* info) const
1630 {
1631 ZTRTRS_F77(CHAR_MACRO(UPLO), CHAR_MACRO(TRANS), CHAR_MACRO(DIAG), &n, &nrhs, A, &lda, B, &ldb, info);
1632 }
1633
1634
1635 void LAPACK<int,std::complex<double> >::TRTRI(const char& UPLO, const char& DIAG, const int& n, const std::complex<double>* A, const int& lda, int* info) const
1636 {
1637 ZTRTRI_F77(CHAR_MACRO(UPLO), CHAR_MACRO(DIAG), &n, A, &lda, info);
1638 }
1639
1640
1641 void LAPACK<int,std::complex<double> >::STEQR(const char& COMPZ, const int& n, double* D, double* E, std::complex<double>* Z, const int& ldz, double* WORK, int* info) const
1642 {
1643 ZSTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
1644 }
1645
1646
1647 void LAPACK<int,std::complex<double> >::PTEQR(const char& COMPZ, const int& n, double* D, double* E, std::complex<double>* Z, const int& ldz, double* WORK, int* info) const
1648 {
1649 ZPTEQR_F77(CHAR_MACRO(COMPZ), &n, D, E, Z, &ldz, WORK, info);
1650 }
1651
1652
1653 void LAPACK<int,std::complex<double> >::HEEV(const char& JOBZ, const char& UPLO, const int& n, std::complex<double> * A, const int& lda, double* W, std::complex<double> * WORK, const int& lwork, double* RWORK, int* info) const
1654 {
1655 ZHEEV_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, W, WORK, &lwork, RWORK, info);
1656 }
1657
1658
1659 void LAPACK<int,std::complex<double> >::HEGV(const int& itype, const char& JOBZ, const char& UPLO, const int& n, std::complex<double> * A, const int& lda, std::complex<double> * B, const int& ldb, double* W, std::complex<double> * WORK, const int& lwork, double* RWORK, int* info) const
1660 {
1661 ZHEGV_F77(&itype, CHAR_MACRO(JOBZ), CHAR_MACRO(UPLO), &n, A, &lda, B, &ldb, W, WORK, &lwork, RWORK, info);
1662 }
1663
1664
1665 void LAPACK<int, std::complex<double> >::HSEQR(const char& JOB, const char& COMPZ, const int& n, const int& ilo, const int& ihi, std::complex<double>* H, const int& ldh, std::complex<double>* W, std::complex<double>* Z, const int& ldz, std::complex<double>* WORK, const int& lwork, int* info) const
1666 {
1667 ZHSEQR_F77(CHAR_MACRO(JOB), CHAR_MACRO(COMPZ), &n, &ilo, &ihi, H, &ldh, W, Z, &ldz, WORK, &lwork, info);
1668 }
1669
1670
1671 void LAPACK<int, std::complex<double> >::GEES(const char& JOBVS, const char& SORT, int (*ptr2func)(std::complex<double>*), const int& n, std::complex<double>* A, const int& lda, int* sdim, std::complex<double>* W, std::complex<double>* VS, const int& ldvs, std::complex<double>* WORK, const int& lwork, double* RWORK, int* BWORK, int* info) const
1672 {
1673 ZGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(SORT), ptr2func, &n, A, &lda, sdim, W, VS, &ldvs, WORK, &lwork, RWORK, BWORK, info);
1674 }
1675
1676
1677 void LAPACK<int, std::complex<double> >::GEES(const char& JOBVS, const int& n, std::complex<double>* A, const int& lda, int* sdim, double* WR, double* WI, std::complex<double>* VS, const int& ldvs, std::complex<double>* WORK, const int& lwork, double* RWORK, int* BWORK, int* info) const
1678 {
1679 int (*nullfptr)(std::complex<double>*) = NULL;
1680 std::vector< std::complex<double> > W(n);
1681 const char sort = 'N';
1682 ZGEES_F77(CHAR_MACRO(JOBVS), CHAR_MACRO(sort), nullfptr, &n, A, &lda, sdim, &W[0], VS, &ldvs, WORK, &lwork, RWORK, BWORK, info);
1683 for (int i=0; i<n; i++) {
1684 WR[i] = W[i].real();
1685 WI[i] = W[i].imag();
1686 }
1687 }
1688
1689
1690 void LAPACK<int, std::complex<double> >::GEEV(const char& JOBVL, const char& JOBVR, const int& n, std::complex<double>* A, const int& lda, std::complex<double>* W, std::complex<double>* VL, const int& ldvl, std::complex<double>* VR, const int& ldvr, std::complex<double>* WORK, const int& lwork, double* RWORK, int* info) const
1691 {
1692 ZGEEV_F77(CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), &n, A, &lda, W, VL, &ldvl, VR, &ldvr, WORK, &lwork, RWORK, info);
1693 }
1694
1695
1696 void LAPACK<int, std::complex<double> >::GEEV(const char& JOBVL, const char& JOBVR, const int& n, std::complex<double>* A, const int& lda, double* WR, double* WI, std::complex<double>* VL, const int& ldvl, std::complex<double>* VR, const int& ldvr, std::complex<double>* WORK, const int& lwork, double* RWORK, int* info) const
1697 {
1698 std::vector<std::complex<double> > w (n);
1699 std::complex<double>* w_rawPtr = (n == 0) ? NULL : &w[0];
1700 GEEV (JOBVL, JOBVR, n, A, lda, w_rawPtr, VL, ldvl, VR, ldvr, WORK, lwork, RWORK, info);
1701 if (*info == 0) {
1702 // The eigenvalues are only valid on output if INFO is zero.
1703 // Otherwise, we shouldn't even write to WR or WI.
1704 for (int k = 0; k < n; ++k) {
1705 WR[k] = w[k].real ();
1706 WI[k] = w[k].imag ();
1707 }
1708 }
1709 }
1710
1711
1712 void LAPACK<int, std::complex<double> >::GESVD(const char& JOBU, const char& JOBVT, const int& m, const int& n, std::complex<double> * A, const int& lda, double* S, std::complex<double> * U, const int& ldu, std::complex<double> * V, const int& ldv, std::complex<double> * WORK, const int& lwork, double* RWORK, int* info) const {
1713 ZGESVD_F77(CHAR_MACRO(JOBU), CHAR_MACRO(JOBVT), &m, &n, A, &lda, S, U, &ldu, V, &ldv, WORK, &lwork, RWORK, info);
1714 }
1715
1716 void LAPACK<int, std::complex<double> >::GEEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, std::complex<double>* A, const int& lda, std::complex<double>* W, std::complex<double>* VL, const int& ldvl, std::complex<double>* VR, const int& ldvr, int* ilo, int* ihi, double* SCALE, double* abnrm, double* RCONDE, double* RCONDV, std::complex<double>* WORK, const int& lwork, double* RWORK, int* info) const
1717 {
1718 ZGEEVX_F77(CHAR_MACRO(BALANC), CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), CHAR_MACRO(SENSE), &n, A, &lda, W, VL, &ldvl, VR, &ldvr, ilo, ihi, SCALE, abnrm, RCONDE, RCONDV, WORK, &lwork, RWORK, info);
1719 }
1720
1721 void LAPACK<int, std::complex<double> >::GGEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, std::complex<double>* A, const int& lda, std::complex<double>* B, const int& ldb, std::complex<double>* ALPHA, std::complex<double>* BETA, std::complex<double>* VL, const int& ldvl, std::complex<double>* VR, const int& ldvr, int* ilo, int* ihi, double* lscale, double* rscale, double* abnrm, double* bbnrm, double* RCONDE, double* RCONDV, std::complex<double>* WORK, const int& lwork, double* RWORK, int* IWORK, int* BWORK, int* info) const
1722 {
1723 ZGGEVX_F77(CHAR_MACRO(BALANC), CHAR_MACRO(JOBVL), CHAR_MACRO(JOBVR), CHAR_MACRO(SENSE), &n, A, &lda, B, &ldb, ALPHA, BETA, VL, &ldvl, VR, &ldvr, ilo, ihi, lscale, rscale, abnrm, bbnrm, RCONDE, RCONDV, WORK, &lwork, RWORK, IWORK, BWORK, info);
1724 }
1725
1726 void LAPACK<int, std::complex<double> >::GGEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, std::complex<double>* A, const int& lda, std::complex<double>* B, const int& ldb, double* ALPHAR, double* ALPHAI, std::complex<double>* BETA, std::complex<double>* VL, const int& ldvl, std::complex<double>* VR, const int& ldvr, int* ilo, int* ihi, double* lscale, double* rscale, double* abnrm, double* bbnrm, double* RCONDE, double* RCONDV, std::complex<double>* WORK, const int& lwork, double* RWORK, int* IWORK, int* BWORK, int* info) const
1727 {
1728 std::vector<std::complex<double> > w (n);
1729 std::complex<double>* w_rawPtr = (n == 0) ? NULL : &w[0];
1730 GGEVX(BALANC, JOBVL, JOBVR, SENSE, n, A, lda, B, ldb, w_rawPtr, BETA, VL, ldvl, VR, ldvr, ilo, ihi, lscale, rscale, abnrm, bbnrm, RCONDE, RCONDV, WORK, lwork, RWORK, IWORK, BWORK, info);
1731 if (*info == 0) {
1732 // The eigenvalues are only valid on output if INFO is zero.
1733 // Otherwise, we shouldn't even write to WR or WI.
1734 for (int k = 0; k < n; ++k) {
1735 ALPHAR[k] = w[k].real ();
1736 ALPHAI[k] = w[k].imag ();
1737 }
1738 }
1739 }
1740
1741 void LAPACK<int, std::complex<double> >::TREVC(const char& SIDE, const char& HOWMNY, int* select, const int& n, const std::complex<double>* T, const int& ldt, std::complex<double>* VL, const int& ldvl, std::complex<double>* VR, const int& ldvr, const int& mm, int* m, std::complex<double>* WORK, double* RWORK, int* info) const
1742 {
1743 ZTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(HOWMNY), select, &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, RWORK, info);
1744 }
1745
1746
1747 void LAPACK<int, std::complex<double> >::TREVC(const char& SIDE, const int& n, const std::complex<double>* T, const int& ldt, std::complex<double>* VL, const int& ldvl, std::complex<double>* VR, const int& ldvr, const int& mm, int* m, std::complex<double>* WORK, double* RWORK, int* info) const
1748 {
1749 std::vector<int> select(1);
1750 const char& whch = 'A';
1751 ZTREVC_F77(CHAR_MACRO(SIDE), CHAR_MACRO(whch), &select[0], &n, T, &ldt, VL, &ldvl, VR, &ldvr, &mm, m, WORK, RWORK, info);
1752 }
1753
1754 void LAPACK<int, std::complex<double> >::TREXC(const char& COMPQ, const int& n, std::complex<double>* T, const int& ldt, std::complex<double>* Q, const int& ldq, int* ifst, int* ilst, std::complex<double>* WORK, int* info) const
1755 {
1756 ZTREXC_F77(CHAR_MACRO(COMPQ), &n, T, &ldt, Q, &ldq, ifst, ilst, info);
1757 }
1758
1759 void LAPACK<int, std::complex<double> >::LARTG( const std::complex<double> f, const std::complex<double> g, double* c, std::complex<double>* s, std::complex<double>* r ) const
1760 {
1761 ZLARTG_F77(&f, &g, c, s, r);
1762 }
1763
1764
1765 void LAPACK<int, std::complex<double> >::LARFG( const int& n, std::complex<double>* alpha, std::complex<double>* x, const int& incx, std::complex<double>* tau ) const
1766 {
1767 ZLARFG_F77(&n, alpha, x, &incx, tau);
1768 }
1769
1770 void LAPACK<int, std::complex<double> >::GEBAL(const char& JOBZ, const int& n, std::complex<double>* A, const int& lda, int* ilo, int* ihi, double* scale, int* info) const
1771 {
1772 ZGEBAL_F77(CHAR_MACRO(JOBZ),&n, A, &lda, ilo, ihi, scale, info);
1773 }
1774
1775
1776 void LAPACK<int, std::complex<double> >::GEBAK(const char& JOBZ, const char& SIDE, const int& n, const int& ilo, const int& ihi, const double* scale, const int& m, std::complex<double>* V, const int& ldv, int* info) const
1777 {
1778 ZGEBAK_F77(CHAR_MACRO(JOBZ), CHAR_MACRO(SIDE), &n, &ilo, &ihi, scale, &m, V, &ldv, info);
1779 }
1780
1781
1782#ifdef HAVE_TEUCHOS_LAPACKLARND
1783 std::complex<double> LAPACK<int, std::complex<double> >::LARND( const int& idist, int* seed ) const
1784 {
1785 return(ZLARND_F77(&idist, seed));
1786 }
1787#endif
1788
1789 void LAPACK<int, std::complex<double> >::LARNV( const int& idist, int* seed, const int& n, std::complex<double>* v ) const
1790 {
1791 ZLARNV_F77(&idist, seed, &n, v);
1792 }
1793
1794
1795 int LAPACK<int, std::complex<double> >::ILAENV( const int& ispec, const std::string& NAME, const std::string& OPTS, const int& N1, const int& N2, const int& N3, const int& N4 ) const
1796 {
1797 unsigned int opts_length = OPTS.length();
1798 std::string temp_NAME = "z" + NAME;
1799 unsigned int name_length = temp_NAME.length();
1800 return ilaenv_wrapper(&ispec, &temp_NAME[0], name_length, &OPTS[0], opts_length, &N1, &N2, &N3, &N4);
1801 }
1802
1803 // END INT, COMPLEX<DOUBLE> SPECIALIZATION IMPLEMENTATION //
1804
1805#endif // HAVE_TEUCHOS_COMPLEX
1806
1807
1808#ifdef HAVE_TEUCHOSCORE_QUADMATH
1809
1810 // BEGIN int, __float128 SPECIALIZATION IMPLEMENTATION //
1811
1813 GEQRF(const int& m, const int& n, __float128* A, const int& lda, __float128* TAU, __float128* WORK, const int& lwork, int* info) const
1814 {
1815 Teuchos::Details::Lapack128 lapack;
1816 lapack.GEQRF (m, n, A, lda, TAU, WORK, lwork, info);
1817 }
1818
1820 GEQR2(const int& m, const int& n, __float128 A[], const int& lda, __float128 TAU[], __float128 WORK[], int* const info) const
1821 {
1822 Teuchos::Details::Lapack128 lapack;
1823 lapack.GEQR2 (m, n, A, lda, TAU, WORK, info);
1824 }
1825
1827 GETRF(const int& m, const int& n, __float128* A, const int& lda, int* IPIV, int* info) const
1828 {
1829 Teuchos::Details::Lapack128 lapack;
1830 lapack.GETRF (m, n, A, lda, IPIV, info);
1831 }
1832
1834 GETRS(const char& TRANS, const int& n, const int& nrhs, const __float128* A, const int& lda, const int* IPIV, __float128* B, const int& ldb, int* info) const
1835 {
1836 Teuchos::Details::Lapack128 lapack;
1837 lapack.GETRS (TRANS, n, nrhs, A, lda, IPIV, B, ldb, info);
1838 }
1839
1841 GETRI (const int& n, __float128* A, const int& lda, const int* IPIV, __float128* WORK, const int& lwork, int* info) const
1842 {
1843 Teuchos::Details::Lapack128 lapack;
1844 lapack.GETRI (n, A, lda, const_cast<int*> (IPIV), WORK, lwork, info);
1845 }
1846
1848 LASWP (const int& N, __float128 A[], const int& LDA, const int& K1, const int& K2, const int IPIV[], const int& INCX) const
1849 {
1850 Teuchos::Details::Lapack128 lapack;
1851 lapack.LASWP (N, A, LDA, K1, K2, IPIV, INCX);
1852 }
1853
1855 ORM2R(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& k, const __float128 A[], const int& lda, const __float128 TAU[], __float128 C[], const int& ldc, __float128 WORK[], int* const info) const
1856 {
1857 Teuchos::Details::Lapack128 lapack;
1858 lapack.ORM2R (SIDE, TRANS, m, n, k, A, lda, TAU, C, ldc, WORK, info);
1859 }
1860
1862 ORGQR(const int& m, const int& n, const int& k, __float128* A, const int& lda, const __float128* TAU, __float128* WORK, const int& lwork, int* info) const
1863 {
1864 Teuchos::Details::Lapack128 lapack;
1865 lapack.ORGQR (m, n, k, A, lda, TAU, WORK, lwork, info);
1866 }
1867
1869 UNGQR(const int& m, const int& n, const int& k, __float128* A, const int& lda, const __float128* TAU, __float128* WORK, const int& lwork, int* info) const
1870 {
1871 Teuchos::Details::Lapack128 lapack;
1872 lapack.UNGQR (m, n, k, A, lda, TAU, WORK, lwork, info);
1873 }
1874
1876 LARFG( const int& n, __float128* alpha, __float128* x, const int& incx, __float128* tau ) const
1877 {
1878 Teuchos::Details::Lapack128 lapack;
1879 lapack.LARFG (n, alpha, x, incx, tau);
1880 }
1881
1882 __float128 LAPACK<int, __float128>::
1883 LAPY2 (const __float128 x, const __float128 y) const
1884 {
1885 Teuchos::Details::Lapack128 lapack;
1886 return lapack.LAPY2 (x, y);
1887 }
1888
1890 GBTRF (const int& m, const int& n, const int& kl, const int& ku,
1891 __float128* A, const int& lda, int* IPIV, int* info) const
1892 {
1893 Teuchos::Details::Lapack128 lapack;
1894 return lapack.GBTRF (m, n, kl, ku, A, lda, IPIV, info);
1895 }
1896
1898 GBTRS (const char& TRANS, const int& n, const int& kl, const int& ku,
1899 const int& nrhs, const __float128* A, const int& lda, const int* IPIV,
1900 __float128* B, const int& ldb, int* info) const
1901 {
1902 Teuchos::Details::Lapack128 lapack;
1903 return lapack.GBTRS (TRANS, n, kl, ku, nrhs, A, lda, IPIV, B, ldb, info);
1904 }
1905
1907 LASCL (const char& TYPE, const int& kl, const int& ku, const __float128 cfrom,
1908 const __float128 cto, const int& m, const int& n, __float128* A,
1909 const int& lda, int* info) const
1910 {
1911 Teuchos::Details::Lapack128 lapack;
1912 return lapack.LASCL (TYPE, kl, ku, cfrom, cto, m, n, A, lda, info);
1913 }
1914
1915 // END int, __float128 SPECIALIZATION IMPLEMENTATION //
1916
1917#endif // HAVE_TEUCHOSCORE_QUADMATH
1918
1919#ifdef HAVE_TEUCHOS_LONG_DOUBLE
1920
1921 // BEGIN int, long double SPECIALIZATION IMPLEMENTATION //
1922
1924 GESV(const int& n, const int& nrhs, long double* A, const int& lda, int* IPIV, long double* B, const int& ldb, int* info) const
1925 {
1926 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: GESV not implemented for long double scalar type!");
1927 }
1929 GTTRS(const char& TRANS, const int& n, const int& nrhs, const long double* dl, const long double* d, const long double* du, const long double* du2, const int* IPIV, long double* B, const int& ldb, int* info) const
1930 {
1931 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: GTTRS not implemented for long double scalar type!");
1932 }
1934 GTTRF(const int& n, long double* dl, long double* d, long double* du, long double* du2, int* IPIV, int* info) const
1935 {
1936 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: GTTRF not implemented for long double scalar type!");
1937 }
1939 SYGV(const int& itype, const char& JOBZ, const char& UPLO, const int& n, long double* A, const int& lda, long double* B, const int& ldb, long double* W, long double* WORK, const int& lwork, int* info) const
1940 {
1941 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: GESV not implemented for long double scalar type!");
1942 }
1944 GEEV(const char& JOBVL, const char& JOBVR, const int& n, long double* A, const int& lda, long double* WR, long double* WI, long double* VL, const int& ldvl, long double* VR, const int& ldvr, long double* WORK, const int& lwork, int* info) const
1945 {
1946 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: GEEV not implemented for long double scalar type!");
1947 }
1949 GEEV(const char& JOBVL, const char& JOBVR, const int& n, long double* A, const int& lda, long double* WR, long double* WI, long double* VL, const int& ldvl, long double* VR, const int& ldvr, long double* WORK, const int& lwork, long double* /* RWORK */, int* info) const
1950 {
1951 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: GEEV not implemented for long double scalar type!");
1952 }
1954 GGEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, long double* A, const int& lda, long double* B, const int& ldb, long double* ALPHAR, long double* ALPHAI, long double* BETA, long double* VL, const int& ldvl, long double* VR, const int& ldvr, int* ilo, int* ihi, long double* lscale, long double* rscale, long double* abnrm, long double* bbnrm, long double* RCONDE, long double* RCONDV, long double* WORK, const int& lwork, int* IWORK, int* BWORK, int* info) const
1955 {
1956 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: GGEVX not implemented for long double scalar type!");
1957 }
1959 GGEVX(const char& BALANC, const char& JOBVL, const char& JOBVR, const char& SENSE, const int& n, long double* A, const int& lda, long double* B, const int& ldb, long double* ALPHAR, long double* ALPHAI, long double* BETA, long double* VL, const int& ldvl, long double* VR, const int& ldvr, int* ilo, int* ihi, long double* lscale, long double* rscale, long double* abnrm, long double* bbnrm, long double* RCONDE, long double* RCONDV, long double* WORK, const int& lwork, long double* /* RWORK */, int* IWORK, int* BWORK, int* info) const
1960 {
1961 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: GGEVX not implemented for long double scalar type!");
1962 }
1964 PORFS(const char& UPLO, const int& n, const int& nrhs, long double* A, const int& lda, const long double* AF, const int& ldaf, const long double* B, const int& ldb, long double* X, const int& ldx, long double* FERR, long double* BERR, long double* WORK, int* IWORK, int* info) const
1965 {
1966 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: PORFS not implemented for long double scalar type!");
1967 }
1969 PTEQR(const char& COMPZ, const int& n, long double* D, long double* E, long double* Z, const int& ldz, long double* WORK, int* info) const
1970 {
1971 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: PTEQR not implemented for long double scalar type!");
1972 }
1974 POTRF(const char& UPLO, const int& n, long double* A, const int& lda, int* info) const
1975 {
1976 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: POTRF not implemented for long double scalar type!");
1977 }
1979 POTRS(const char& UPLO, const int& n, const int& nrhs, const long double* A, const int& lda, long double* B, const int& ldb, int* info) const
1980 {
1981 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: POTRS not implemented for long double scalar type!");
1982 }
1984 POEQU(const int& n, const long double* A, const int& lda, long double* S, long double* scond, long double* amax, int* info) const
1985 {
1986 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ERROR in Teuchos::LAPACK: POEQU not implemented for long double scalar type!");
1987 }
1989 GEQRF(const int& m, const int& n, long double* A, const int& lda, long double* TAU, long double* WORK, const int& lwork, int* info) const
1990 {
1991 Teuchos::Details::LapackLongDouble lapack;
1992 lapack.GEQRF (m, n, A, lda, TAU, WORK, lwork, info);
1993 }
1994
1996 GEQR2(const int& m, const int& n, long double A[], const int& lda, long double TAU[], long double WORK[], int* const info) const
1997 {
1998 Teuchos::Details::LapackLongDouble lapack;
1999 lapack.GEQR2 (m, n, A, lda, TAU, WORK, info);
2000 }
2001
2003 GETRF(const int& m, const int& n, long double* A, const int& lda, int* IPIV, int* info) const
2004 {
2005 Teuchos::Details::LapackLongDouble lapack;
2006 lapack.GETRF (m, n, A, lda, IPIV, info);
2007 }
2008
2010 GETRS(const char& TRANS, const int& n, const int& nrhs, const long double* A, const int& lda, const int* IPIV, long double* B, const int& ldb, int* info) const
2011 {
2012 Teuchos::Details::LapackLongDouble lapack;
2013 lapack.GETRS (TRANS, n, nrhs, A, lda, IPIV, B, ldb, info);
2014 }
2015
2017 GETRI (const int& n, long double* A, const int& lda, const int* IPIV, long double* WORK, const int& lwork, int* info) const
2018 {
2019 Teuchos::Details::LapackLongDouble lapack;
2020 lapack.GETRI (n, A, lda, const_cast<int*> (IPIV), WORK, lwork, info);
2021 }
2022
2024 LASWP (const int& N, long double A[], const int& LDA, const int& K1, const int& K2, const int IPIV[], const int& INCX) const
2025 {
2026 Teuchos::Details::LapackLongDouble lapack;
2027 lapack.LASWP (N, A, LDA, K1, K2, IPIV, INCX);
2028 }
2029
2031 ORM2R(const char& SIDE, const char& TRANS, const int& m, const int& n, const int& k, const long double A[], const int& lda, const long double TAU[], long double C[], const int& ldc, long double WORK[], int* const info) const
2032 {
2033 Teuchos::Details::LapackLongDouble lapack;
2034 lapack.ORM2R (SIDE, TRANS, m, n, k, A, lda, TAU, C, ldc, WORK, info);
2035 }
2036
2038 ORGQR(const int& m, const int& n, const int& k, long double* A, const int& lda, const long double* TAU, long double* WORK, const int& lwork, int* info) const
2039 {
2040 Teuchos::Details::LapackLongDouble lapack;
2041 lapack.ORGQR (m, n, k, A, lda, TAU, WORK, lwork, info);
2042 }
2043
2045 UNGQR(const int& m, const int& n, const int& k, long double* A, const int& lda, const long double* TAU, long double* WORK, const int& lwork, int* info) const
2046 {
2047 Teuchos::Details::LapackLongDouble lapack;
2048 lapack.UNGQR (m, n, k, A, lda, TAU, WORK, lwork, info);
2049 }
2050
2052 LARFG( const int& n, long double* alpha, long double* x, const int& incx, long double* tau ) const
2053 {
2054 Teuchos::Details::LapackLongDouble lapack;
2055 lapack.LARFG (n, alpha, x, incx, tau);
2056 }
2057
2058 long double LAPACK<int, long double>::
2059 LAPY2 (const long double x, const long double y) const
2060 {
2061 Teuchos::Details::LapackLongDouble lapack;
2062 return lapack.LAPY2 (x, y);
2063 }
2064
2066 GBTRF (const int& m, const int& n, const int& kl, const int& ku,
2067 long double* A, const int& lda, int* IPIV, int* info) const
2068 {
2069 Teuchos::Details::LapackLongDouble lapack;
2070 return lapack.GBTRF (m, n, kl, ku, A, lda, IPIV, info);
2071 }
2072
2074 GBTRS (const char& TRANS, const int& n, const int& kl, const int& ku,
2075 const int& nrhs, const long double* A, const int& lda, const int* IPIV,
2076 long double* B, const int& ldb, int* info) const
2077 {
2078 Teuchos::Details::LapackLongDouble lapack;
2079 return lapack.GBTRS (TRANS, n, kl, ku, nrhs, A, lda, IPIV, B, ldb, info);
2080 }
2081
2083 LASCL (const char& TYPE, const int& kl, const int& ku, const long double cfrom,
2084 const long double cto, const int& m, const int& n, long double* A,
2085 const int& lda, int* info) const
2086 {
2087 Teuchos::Details::LapackLongDouble lapack;
2088 return lapack.LASCL (TYPE, kl, ku, cfrom, cto, m, n, A, lda, info);
2089 }
2090
2091 // END int, long double SPECIALIZATION IMPLEMENTATION //
2092
2093#endif // HAVE_TEUCHOS_LONG_DOUBLE
2094
2095
2096} // namespace Teuchos
Teuchos header file which uses auto-configuration information to include necessary C++ headers.
Declaration and definition of Teuchos::Details::Lapack128, a partial implementation of Teuchos::LAPAC...
Declaration and definition of Teuchos::Details::LapackLongDouble, a partial implementation of Teuchos...
Templated interface class to LAPACK routines.
The Templated LAPACK wrappers.
Defines basic traits for the scalar field type.
Standard test and throw macros.
void LARFG(const OrdinalType &n, ScalarType *alpha, ScalarType *x, const OrdinalType &incx, ScalarType *tau) const
Generates an elementary reflector of order n that zeros out the last n-1 components of the input vect...
void GEEV(const char &JOBVL, const char &JOBVR, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, MagnitudeType *WR, MagnitudeType *WI, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const
Computes for an n by n real nonsymmetric matrix A, the eigenvalues and, optionally,...
void TGSEN(const OrdinalType &ijob, const OrdinalType &wantq, const OrdinalType &wantz, const OrdinalType *SELECT, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, MagnitudeType *BETA, ScalarType *Q, const OrdinalType &ldq, ScalarType *Z, const OrdinalType &ldz, OrdinalType *M, MagnitudeType *PL, MagnitudeType *PR, MagnitudeType *DIF, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *IWORK, const OrdinalType &liwork, OrdinalType *info) const
void TREXC(const char &COMPQ, const OrdinalType &n, ScalarType *T, const OrdinalType &ldt, ScalarType *Q, const OrdinalType &ldq, OrdinalType *ifst, OrdinalType *ilst, ScalarType *WORK, OrdinalType *info) const
void SYGV(const OrdinalType &itype, const char &JOBZ, const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, ScalarType *W, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Computes all the eigenvalues and, optionally, eigenvectors of a symmetric n by n matrix pencil {A,...
void GEBAL(const char &JOBZ, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *ilo, OrdinalType *ihi, MagnitudeType *scale, OrdinalType *info) const
Balances a general matrix A, through similarity transformations to make the rows and columns as close...
void GBRFS(const char &TRANS, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const ScalarType *AF, const OrdinalType &ldaf, const OrdinalType *IPIV, const ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Improves the computed solution to a banded system of linear equations and provides error bounds and b...
ScalarTraits< ScalarType >::magnitudeType LANGB(const char &NORM, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const ScalarType *A, const OrdinalType &lda, MagnitudeType *WORK) const
Returns the value of the one norm, or the Frobenius norm, or the infinity norm, or the element of lar...
void TRTRS(const char &UPLO, const char &TRANS, const char &DIAG, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Solves a triangular linear system of the form A*X=B or A**T*X=B, where A is a triangular matrix.
void HEGV(const OrdinalType &itype, const char &JOBZ, const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, MagnitudeType *W, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const
Computes all the eigenvalues and, optionally, eigenvectors of a generalized Hermitian-definite n by n...
void STEQR(const char &COMPZ, const OrdinalType &n, MagnitudeType *D, MagnitudeType *E, ScalarType *Z, const OrdinalType &ldz, MagnitudeType *WORK, OrdinalType *info) const
Computes the eigenvalues and, optionally, eigenvectors of a symmetric tridiagonal n by n matrix A usi...
void GGEV(const char &JOBVL, const char &JOBVR, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, ScalarType *BETA, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
void ORMHR(const char &SIDE, const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &ilo, const OrdinalType &ihi, const ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *C, const OrdinalType &ldc, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Overwrites the general real m by n matrix C with the product of C and Q, which is a product of ihi-il...
void UNGQR(const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Compute explicit QR factor from QR factorization (GEQRF) (complex case).
void UNM2R(const char &SIDE, const char &TRANS, const OrdinalType &M, const OrdinalType &N, const OrdinalType &K, const ScalarType A[], const OrdinalType &LDA, const ScalarType TAU[], ScalarType C[], const OrdinalType &LDC, ScalarType WORK[], OrdinalType *const INFO) const
BLAS 2 version of UNMQR; known workspace size.
ScalarType LAMCH(const char &CMACH) const
Determines machine parameters for floating point characteristics.
void ORGHR(const OrdinalType &n, const OrdinalType &ilo, const OrdinalType &ihi, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Generates a real orthogonal matrix Q which is the product of ihi-ilo elementary reflectors of order n...
void PTTRS(const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *d, const ScalarType *e, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Solves a tridiagonal system A*X=B using the \L*D*L' factorization of A computed by PTTRF.
void GTTRF(const OrdinalType &n, ScalarType *dl, ScalarType *d, ScalarType *du, ScalarType *du2, OrdinalType *IPIV, OrdinalType *info) const
Computes an LU factorization of a n by n tridiagonal matrix A using partial pivoting with row interch...
void GESVX(const char &FACT, const char &TRANS, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *AF, const OrdinalType &ldaf, OrdinalType *IPIV, char *EQUED, ScalarType *R, ScalarType *C, ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *rcond, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Uses the LU factorization to compute the solution to a real system of linear equations A*X=B,...
void PTEQR(const char &COMPZ, const OrdinalType &n, MagnitudeType *D, MagnitudeType *E, ScalarType *Z, const OrdinalType &ldz, MagnitudeType *WORK, OrdinalType *info) const
Computes the eigenvalues and, optionally, eigenvectors of a symmetric positive-definite tridiagonal n...
void GETRI(const OrdinalType &n, ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Computes the inverse of a matrix A using the LU factorization computed by GETRF.
void GGLSE(const OrdinalType &m, const OrdinalType &n, const OrdinalType &p, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, ScalarType *C, ScalarType *D, ScalarType *X, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Solves the linear equality-constrained least squares (LSE) problem where A is an m by n matrix,...
void GEQR2(const OrdinalType &m, const OrdinalType &n, ScalarType A[], const OrdinalType &lda, ScalarType TAU[], ScalarType WORK[], OrdinalType *const info) const
BLAS 2 version of GEQRF, with known workspace size.
void GBTRF(const OrdinalType &m, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, OrdinalType *info) const
Computes an LU factorization of a general banded m by n matrix A using partial pivoting with row inte...
void PORFS(const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const ScalarType *AF, const OrdinalType &ldaf, const ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Improves the computed solution to a system of linear equations when the coefficient matrix is symmetr...
void TGEVC(const char &SIDE, const char &HOWMNY, const OrdinalType *SELECT, const OrdinalType &n, ScalarType *S, const OrdinalType &lds, ScalarType *P, const OrdinalType &ldp, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, const OrdinalType &mm, OrdinalType *M, ScalarType *WORK, OrdinalType *info) const
void GELSS(const OrdinalType &m, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, MagnitudeType *S, const MagnitudeType rcond, OrdinalType *rank, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const
Use the SVD to solve a possibly rank-deficient linear least-squares problem.
void GEHRD(const OrdinalType &n, const OrdinalType &ilo, const OrdinalType &ihi, ScalarType *A, const OrdinalType &lda, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Reduces a real general matrix A to upper Hessenberg form by orthogonal similarity transformations.
void POSVX(const char &FACT, const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *AF, const OrdinalType &ldaf, char *EQUED, ScalarType *S, ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *rcond, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Uses the Cholesky factorization to compute the solution to a real system of linear equations A*X=B,...
void ORGQR(const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Compute explicit Q factor from QR factorization (GEQRF) (real case).
void GEES(const char &JOBVS, const char &SORT, OrdinalType &(*ptr2func)(ScalarType *, ScalarType *), const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *sdim, ScalarType *WR, ScalarType *WI, ScalarType *VS, const OrdinalType &ldvs, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *BWORK, OrdinalType *info) const
void LATRS(const char &UPLO, const char &TRANS, const char &DIAG, const char &NORMIN, const OrdinalType &N, ScalarType *A, const OrdinalType &LDA, ScalarType *X, MagnitudeType *SCALE, MagnitudeType *CNORM, OrdinalType *INFO) const
Robustly solve a possibly singular triangular linear system.
void GBEQU(const OrdinalType &m, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const ScalarType *A, const OrdinalType &lda, MagnitudeType *R, MagnitudeType *C, MagnitudeType *rowcond, MagnitudeType *colcond, MagnitudeType *amax, OrdinalType *info) const
Computes row and column scalings intended to equilibrate an m by n banded matrix A and reduce its con...
void GESVD(const char &JOBU, const char &JOBVT, const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, MagnitudeType *S, ScalarType *U, const OrdinalType &ldu, ScalarType *V, const OrdinalType &ldv, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const
Computes the singular values (and optionally, vectors) of a real matrix A.
void POCON(const char &UPLO, const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, const ScalarType &anorm, ScalarType *rcond, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Estimates the reciprocal of the condition number (1-norm) of a real symmetric positive definite matri...
void GEEVX(const char &BALANC, const char &JOBVL, const char &JOBVR, const char &SENSE, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *WR, ScalarType *WI, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, OrdinalType *ilo, OrdinalType *ihi, MagnitudeType *SCALE, MagnitudeType *abnrm, MagnitudeType *RCONDE, MagnitudeType *RCONDV, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *IWORK, OrdinalType *info) const
void HSEQR(const char &JOB, const char &COMPZ, const OrdinalType &n, const OrdinalType &ilo, const OrdinalType &ihi, ScalarType *H, const OrdinalType &ldh, ScalarType *WR, ScalarType *WI, ScalarType *Z, const OrdinalType &ldz, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Computes the eigenvalues of a real upper Hessenberg matrix H and, optionally, the matrices T and Z fr...
void LARTG(const ScalarType &f, const ScalarType &g, MagnitudeType *c, ScalarType *s, ScalarType *r) const
Gnerates a plane rotation that zeros out the second component of the input vector.
void GELS(const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Solves an over/underdetermined real m by n linear system A using QR or LQ factorization of A.
void HEEV(const char &JOBZ, const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, MagnitudeType *W, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const
Computes all the eigenvalues and, optionally, eigenvectors of a Hermitian n by n matrix A.
void POTRS(const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Solves a system of linear equations A*X=B, where A is a symmetric positive definite matrix factored b...
void GGEVX(const char &BALANC, const char &JOBVL, const char &JOBVR, const char &SENSE, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, ScalarType *BETA, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, OrdinalType *ilo, OrdinalType *ihi, MagnitudeType *lscale, MagnitudeType *rscale, MagnitudeType *abnrm, MagnitudeType *bbnrm, MagnitudeType *RCONDE, MagnitudeType *RCONDV, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *IWORK, OrdinalType *BWORK, OrdinalType *info) const
ScalarType LARND(const OrdinalType &idist, OrdinalType *seed) const
Returns a random number from a uniform or normal distribution.
void UNMQR(const char &SIDE, const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *C, const OrdinalType &ldc, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Apply Householder reflectors (complex case).
ScalarType LAPY2(const ScalarType &x, const ScalarType &y) const
Computes x^2 + y^2 safely, to avoid overflow.
void ORM2R(const char &SIDE, const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const ScalarType A[], const OrdinalType &lda, const ScalarType TAU[], ScalarType C[], const OrdinalType &ldc, ScalarType WORK[], OrdinalType *const info) const
BLAS 2 version of ORMQR; known workspace size.
void POEQU(const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, MagnitudeType *S, MagnitudeType *scond, MagnitudeType *amax, OrdinalType *info) const
Computes row and column scalings intended to equilibrate a symmetric positive definite matrix A and r...
void GEEQU(const OrdinalType &m, const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, ScalarType *R, ScalarType *C, ScalarType *rowcond, ScalarType *colcond, ScalarType *amax, OrdinalType *info) const
Computes row and column scalings intended to equilibrate an m by n matrix A and reduce its condition ...
void TRTRI(const char &UPLO, const char &DIAG, const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, OrdinalType *info) const
Computes the inverse of an upper or lower triangular matrix A.
void POTRI(const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *info) const
Computes the inverse of a real symmetric positive definite matrix A using the Cholesky factorization ...
void TRSEN(const char &JOB, const char &COMPQ, const OrdinalType *SELECT, const OrdinalType &n, ScalarType *T, const OrdinalType &ldt, ScalarType *Q, const OrdinalType &ldq, MagnitudeType *WR, MagnitudeType *WI, OrdinalType *M, ScalarType *S, MagnitudeType *SEP, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *IWORK, const OrdinalType &liwork, OrdinalType *info) const
void PTTRF(const OrdinalType &n, ScalarType *d, ScalarType *e, OrdinalType *info) const
Computes the L*D*L' factorization of a Hermitian/symmetric positive definite tridiagonal matrix A.
void GESV(const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Computes the solution to a real system of linear equations A*X=B, where A is factored through GETRF a...
void GBTRS(const char &TRANS, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Solves a system of linear equations A*X=B or A'*X=B with a general banded n by n matrix A using the L...
void LARNV(const OrdinalType &idist, OrdinalType *seed, const OrdinalType &n, ScalarType *v) const
Returns a vector of random numbers from a chosen distribution.
void GGES(const char &JOBVL, const char &JOBVR, const char &SORT, OrdinalType &(*ptr2func)(ScalarType *, ScalarType *, ScalarType *), const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *sdim, MagnitudeType *ALPHAR, MagnitudeType *ALPHAI, MagnitudeType *BETA, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *BWORK, OrdinalType *info) const
void GBCON(const char &NORM, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, const ScalarType &anorm, ScalarType *rcond, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Estimates the reciprocal of the condition number of a general banded real matrix A,...
void POSV(const char &UPLO, const OrdinalType &n, const OrdinalType &nrhs, ScalarType *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Computes the solution to a real system of linear equations A*X=B, where A is a symmetric positive def...
OrdinalType ILAENV(const OrdinalType &ispec, const std::string &NAME, const std::string &OPTS, const OrdinalType &N1=-1, const OrdinalType &N2=-1, const OrdinalType &N3=-1, const OrdinalType &N4=-1) const
Chooses problem-dependent parameters for the local environment.
void GETRF(const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, OrdinalType *info) const
Computes an LU factorization of a general m by n matrix A using partial pivoting with row interchange...
void LASCL(const char &TYPE, const OrdinalType &kl, const OrdinalType &ku, const MagnitudeType cfrom, const MagnitudeType cto, const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *info) const
Multiplies the m by n matrix A by the real scalar cto/cfrom.
void GECON(const char &NORM, const OrdinalType &n, const ScalarType *A, const OrdinalType &lda, const ScalarType &anorm, ScalarType *rcond, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Estimates the reciprocal of the condition number of a general real matrix A, in either the 1-norm or ...
void GERFS(const char &TRANS, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const ScalarType *AF, const OrdinalType &ldaf, const OrdinalType *IPIV, const ScalarType *B, const OrdinalType &ldb, ScalarType *X, const OrdinalType &ldx, ScalarType *FERR, ScalarType *BERR, ScalarType *WORK, OrdinalType *IWORK, OrdinalType *info) const
Improves the computed solution to a system of linear equations and provides error bounds and backward...
void GTTRS(const char &TRANS, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *dl, const ScalarType *d, const ScalarType *du, const ScalarType *du2, const OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Solves a system of linear equations A*X=B or A'*X=B or A^H*X=B with a tridiagonal matrix A using the ...
void SYTRD(const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *D, ScalarType *E, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Reduces a real symmetric matrix A to tridiagonal form by orthogonal similarity transformations.
void GEBAK(const char &JOBZ, const char &SIDE, const OrdinalType &n, const OrdinalType &ilo, const OrdinalType &ihi, const MagnitudeType *scale, const OrdinalType &m, ScalarType *V, const OrdinalType &ldv, OrdinalType *info) const
Forms the left or right eigenvectors of a general matrix that has been balanced by GEBAL by backward ...
void GEQP3(const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *jpvt, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, MagnitudeType *RWORK, OrdinalType *info) const
Computes a QR factorization with column pivoting of a matrix A: A*P = Q*R using Level 3 BLAS.
void LASWP(const OrdinalType &N, ScalarType A[], const OrdinalType &LDA, const OrdinalType &K1, const OrdinalType &K2, const OrdinalType IPIV[], const OrdinalType &INCX) const
Apply a series of row interchanges to the matrix A.
void SYEV(const char &JOBZ, const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *W, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Computes all the eigenvalues and, optionally, eigenvectors of a symmetric n by n matrix A.
void GEQRF(const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Computes a QR factorization of a general m by n matrix A.
void POTRF(const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *info) const
Computes Cholesky factorization of a real symmetric positive definite matrix A.
void GETRS(const char &TRANS, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
Solves a system of linear equations A*X=B or A'*X=B with a general n by n matrix A using the LU facto...
void ORMQR(const char &SIDE, const char &TRANS, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *C, const OrdinalType &ldc, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
void TREVC(const char &SIDE, const char &HOWMNY, OrdinalType *select, const OrdinalType &n, const ScalarType *T, const OrdinalType &ldt, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, const OrdinalType &mm, OrdinalType *m, ScalarType *WORK, OrdinalType *info) const
void SPEV(const char &JOBZ, const char &UPLO, const OrdinalType &n, ScalarType *AP, ScalarType *W, ScalarType *Z, const OrdinalType &ldz, ScalarType *WORK, OrdinalType *info) const
Computes the eigenvalues and, optionally, eigenvectors of a symmetric n by n matrix A in packed stora...
Smart reference counting pointer class for automatic garbage collection.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Macro for throwing an exception with breakpointing to ease debugging.
The Teuchos namespace contains all of the classes, structs and enums used by Teuchos,...