MADNESS  version 0.9
tensor_lapack.h
Go to the documentation of this file.
1 /*
2  This file is part of MADNESS.
3 
4  Copyright (C) 2007,2010 Oak Ridge National Laboratory
5 
6  This program is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2 of the License, or
9  (at your option) any later version.
10 
11  This program is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with this program; if not, write to the Free Software
18  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 
20  For more information please contact:
21 
22  Robert J. Harrison
23  Oak Ridge National Laboratory
24  One Bethel Valley Road
25  P.O. Box 2008, MS-6367
26 
27  email: harrisonrj@ornl.gov
28  tel: 865-241-3937
29  fax: 865-572-0680
30 
31 
32  $Id$
33 */
34 
35 
36 #ifndef MADNESS_LINALG_TENSOR_LAPACK_H__INCLUDED
37 #define MADNESS_LINALG_TENSOR_LAPACK_H__INCLUDED
38 
39 #include <madness/tensor/tensor.h>
40 #include <madness/fortran_ctypes.h>
41 
49 namespace madness {
50 
52 
54  template <typename T>
55  void svd(const Tensor<T>& a, Tensor<T>& U,
56  Tensor< typename Tensor<T>::scalar_type >& s, Tensor<T>& VT);
57 
59  template <typename T>
60  void svd_result(Tensor<T>& a, Tensor<T>& U,
61  Tensor< typename Tensor<T>::scalar_type >& s, Tensor<T>& VT, Tensor<T>& work);
62 
64 
66  template <typename T>
67  void gesv(const Tensor<T>& a, const Tensor<T>& b, Tensor<T>& x);
68 
70 
72  template <typename T>
73  void gelss(const Tensor<T>& a, const Tensor<T>& b, double rcond,
74  Tensor<T>& x, Tensor< typename Tensor<T>::scalar_type >& s,
75  long &rank, Tensor<typename Tensor<T>::scalar_type>& sumsq);
76 
78 
80  template <typename T>
81  void syev(const Tensor<T>& A,
82  Tensor<T>& V, Tensor< typename Tensor<T>::scalar_type >& e);
83 
85 
87  template <typename T>
88  void gesv(const Tensor<T>& a, const Tensor<T>& b, Tensor<T>& x);
89 
90 
92 
94  template <typename T>
95  void sygv(const Tensor<T>& A, const Tensor<T>& B, int itype,
96  Tensor<T>& V, Tensor< typename Tensor<T>::scalar_type >& e);
97 
98 #ifdef MADNESS_HAS_ELEMENTAL
99  class World; // UGH!
101 
103  template <typename T>
104  void sygvp(World& world, const Tensor<T>& A, const Tensor<T>& B, int itype,
105  Tensor<T>& V, Tensor< typename Tensor<T>::scalar_type >& e);
106 
108 
110  template <typename T>
111  void gesvp(World& world, const Tensor<T>& a, const Tensor<T>& b, Tensor<T>& x);
112 
113 #endif
114 
116 
118  template <typename T>
119  void cholesky(Tensor<T>& A);
120 
122  template<typename T>
123  void qr(Tensor<T>& A, Tensor<T>& R);
124 
126  template<typename T>
127  void lq(Tensor<T>& A, Tensor<T>& L);
128 
129  template <typename T>
130  void geqp3(Tensor<T>& A, Tensor<T>& tau, Tensor<integer>& jpvt);
131 
133 
138  template <typename T>
139  void orgqr(Tensor<T>& A, const Tensor<T>& tau);
140 
141 
143 
144 // /// \ingroup linalg
145 // template <typename T>
146 // void triangular_solve(const Tensor<T>& L, Tensor<T>& B,
147 // const char* side, const char* transa);
148 
150 
152  bool test_tensor_lapack();
153 
155 
157  void init_tensor_lapack();
158 }
159 #endif // MADNESS_LINALG_TENSOR_LAPACK_H__INCLUDED
Corresponding C and Fortran types.
Tensor< double > B
Definition: tdse1d.cc:167
void gesv(const Tensor< T > &a, const Tensor< T > &b, Tensor< T > &x)
Solve Ax = b for general A using the LAPACK *gesv routines.
Definition: lapack.cc:339
void qr(Tensor< T > &A, Tensor< T > &R)
compute the QR decomposition of the matrix A
Definition: lapack.cc:608
void gelss(const Tensor< T > &a, const Tensor< T > &b, double rcond, Tensor< T > &x, Tensor< typename Tensor< T >::scalar_type > &s, long &rank, Tensor< typename Tensor< T >::scalar_type > &sumsq)
Solve Ax = b for general A using the LAPACK *gelss routines.
Definition: lapack.cc:393
void lq(Tensor< T > &A, Tensor< T > &R)
compute the LQ decomposition of the matrix A = L Q
Definition: lapack.cc:630
double V(const Vector< double, 3 > &r)
Definition: apps/ii/hatom_energy.cc:46
void syev(const Tensor< T > &A, Tensor< T > &V, Tensor< typename Tensor< T >::scalar_type > &e)
Real-symmetric or complex-Hermitian eigenproblem.
Definition: lapack.cc:478
Defines and implements most of Tensor.
void sygvp(World &world, const Tensor< T > &a, const Tensor< T > &B, int itype, Tensor< T > &V, Tensor< typename Tensor< T >::scalar_type > &e)
Definition: elem.h:508
FLOAT a(int j, FLOAT z)
Definition: y1.cc:86
bool test_tensor_lapack()
Test the Tensor-LAPACK interface ... currently always returns true!
Definition: lapack.cc:908
void svd(const Tensor< T > &a, Tensor< T > &U, Tensor< typename Tensor< T >::scalar_type > &s, Tensor< T > &VT)
Compute the singluar value decomposition of an n-by-m matrix using *gesvd.
Definition: lapack.cc:273
void gesvp(World &world, const Tensor< T > &a, const Tensor< T > &b, Tensor< T > &x)
Definition: elem.h:516
void init_tensor_lapack()
World/MRA initialization calls this before going multithreaded due to static data in dlamch...
Definition: lapack.cc:893
void orgqr(Tensor< T > &A, const Tensor< T > &tau)
reconstruct the orthogonal matrix Q (e.g. from QR factorization)
Definition: lapack.cc:706
void cholesky(Tensor< T > &A)
Compute the Cholesky factorization.
Definition: lapack.cc:551
void svd_result(Tensor< T > &a, Tensor< T > &U, Tensor< typename Tensor< T >::scalar_type > &s, Tensor< T > &VT, Tensor< T > &work)
same as svd, but it optimizes away the tensor construction: a = U * diag(s) * VT
Definition: lapack.cc:308
void sygv(const Tensor< T > &A, const Tensor< T > &B, int itype, Tensor< T > &V, Tensor< typename Tensor< T >::scalar_type > &e)
Generalized real-symmetric or complex-Hermitian eigenproblem.
Definition: lapack.cc:519
void geqp3(Tensor< T > &A, Tensor< T > &tau, Tensor< integer > &jpvt)
Compute the QR factorization.
Definition: lapack.cc:570
Holds machinery to set up Functions/FuncImpls using various Factories and Interfaces.
Definition: chem/atomutil.cc:45
FLOAT b(int j, FLOAT z)
Definition: y1.cc:79