ergo
template_lapack_gesv.h
Go to the documentation of this file.
1 /* Ergo, version 3.3, a program for linear scaling electronic structure
2  * calculations.
3  * Copyright (C) 2013 Elias Rudberg, Emanuel H. Rubensson, and Pawel Salek.
4  *
5  * This program is free software: you can redistribute it and/or modify
6  * it under the terms of the GNU General Public License as published by
7  * the Free Software Foundation, either version 3 of the License, or
8  * (at your option) any later version.
9  *
10  * This program is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  * GNU General Public License for more details.
14  *
15  * You should have received a copy of the GNU General Public License
16  * along with this program. If not, see <http://www.gnu.org/licenses/>.
17  *
18  * Primary academic reference:
19  * Kohn−Sham Density Functional Theory Electronic Structure Calculations
20  * with Linearly Scaling Computational Time and Memory Usage,
21  * Elias Rudberg, Emanuel H. Rubensson, and Pawel Salek,
22  * J. Chem. Theory Comput. 7, 340 (2011),
23  * <http://dx.doi.org/10.1021/ct100611z>
24  *
25  * For further information about Ergo, see <http://www.ergoscf.org>.
26  */
27 
28  /* This file belongs to the template_lapack part of the Ergo source
29  * code. The source files in the template_lapack directory are modified
30  * versions of files originally distributed as CLAPACK, see the
31  * Copyright/license notice in the file template_lapack/COPYING.
32  */
33 
34 
35 #ifndef TEMPLATE_LAPACK_GESV_HEADER
36 #define TEMPLATE_LAPACK_GESV_HEADER
37 
38 
39 template<class Treal>
40 int template_lapack_gesv(const integer *n, const integer *nrhs, Treal *a, const integer
41  *lda, integer *ipiv, Treal *b, const integer *ldb, integer *info)
42 {
43 /* -- LAPACK driver routine (version 3.0) --
44  Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
45  Courant Institute, Argonne National Lab, and Rice University
46  March 31, 1993
47 
48 
49  Purpose
50  =======
51 
52  DGESV computes the solution to a real system of linear equations
53  A * X = B,
54  where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
55 
56  The LU decomposition with partial pivoting and row interchanges is
57  used to factor A as
58  A = P * L * U,
59  where P is a permutation matrix, L is unit lower triangular, and U is
60  upper triangular. The factored form of A is then used to solve the
61  system of equations A * X = B.
62 
63  Arguments
64  =========
65 
66  N (input) INTEGER
67  The number of linear equations, i.e., the order of the
68  matrix A. N >= 0.
69 
70  NRHS (input) INTEGER
71  The number of right hand sides, i.e., the number of columns
72  of the matrix B. NRHS >= 0.
73 
74  A (input/output) DOUBLE PRECISION array, dimension (LDA,N)
75  On entry, the N-by-N coefficient matrix A.
76  On exit, the factors L and U from the factorization
77  A = P*L*U; the unit diagonal elements of L are not stored.
78 
79  LDA (input) INTEGER
80  The leading dimension of the array A. LDA >= max(1,N).
81 
82  IPIV (output) INTEGER array, dimension (N)
83  The pivot indices that define the permutation matrix P;
84  row i of the matrix was interchanged with row IPIV(i).
85 
86  B (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
87  On entry, the N-by-NRHS matrix of right hand side matrix B.
88  On exit, if INFO = 0, the N-by-NRHS solution matrix X.
89 
90  LDB (input) INTEGER
91  The leading dimension of the array B. LDB >= max(1,N).
92 
93  INFO (output) INTEGER
94  = 0: successful exit
95  < 0: if INFO = -i, the i-th argument had an illegal value
96  > 0: if INFO = i, U(i,i) is exactly zero. The factorization
97  has been completed, but the factor U is exactly
98  singular, so the solution could not be computed.
99 
100  =====================================================================
101 
102 
103  Test the input parameters.
104 
105  Parameter adjustments */
106  /* System generated locals */
107  integer a_dim1, a_offset, b_dim1, b_offset, i__1;
108  /* Local variables */
109 
110  a_dim1 = *lda;
111  a_offset = 1 + a_dim1 * 1;
112  a -= a_offset;
113  --ipiv;
114  b_dim1 = *ldb;
115  b_offset = 1 + b_dim1 * 1;
116  b -= b_offset;
117 
118  /* Function Body */
119  *info = 0;
120  if (*n < 0) {
121  *info = -1;
122  } else if (*nrhs < 0) {
123  *info = -2;
124  } else if (*lda < maxMACRO(1,*n)) {
125  *info = -4;
126  } else if (*ldb < maxMACRO(1,*n)) {
127  *info = -7;
128  }
129  if (*info != 0) {
130  i__1 = -(*info);
131  template_blas_erbla("GESV ", &i__1);
132  return 0;
133  }
134 
135 /* Compute the LU factorization of A. */
136 
137  template_lapack_getrf(n, n, &a[a_offset], lda, &ipiv[1], info);
138  if (*info == 0) {
139 
140 /* Solve the system A*X = B, overwriting B with X. */
141 
142  template_lapack_getrs("No transpose", n, nrhs, &a[a_offset], lda, &ipiv[1], &b[
143  b_offset], ldb, info);
144  }
145  return 0;
146 
147 /* End of DGESV */
148 
149 } /* dgesv_ */
150 
151 #endif
int integer
Definition: template_blas_common.h:38
#define maxMACRO(a, b)
Definition: template_blas_common.h:43
int template_blas_erbla(const char *srname, integer *info)
Definition: template_blas_common.cc:144
int template_lapack_getrs(const char *trans, const integer *n, const integer *nrhs, const Treal *a, const integer *lda, const integer *ipiv, Treal *b, const integer *ldb, integer *info)
Definition: template_lapack_getrs.h:40
int template_lapack_getrf(const integer *m, const integer *n, Treal *a, const integer *lda, integer *ipiv, integer *info)
Definition: template_lapack_getrf.h:40
int template_lapack_gesv(const integer *n, const integer *nrhs, Treal *a, const integer *lda, integer *ipiv, Treal *b, const integer *ldb, integer *info)
Definition: template_lapack_gesv.h:40