2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2012,2013,2014,2015,2016,2017,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 #include "eigensolver.h"
41 #include "gromacs/linearalgebra/sparsematrix.h"
42 #include "gromacs/utility/fatalerror.h"
43 #include "gromacs/utility/real.h"
44 #include "gromacs/utility/smalloc.h"
46 #include "gmx_arpack.h"
47 #include "gmx_lapack.h"
49 void eigensolver(real* a, int n, int index_lower, int index_upper, real* eigenvalues, real* eigenvectors)
70 /* Make jobz point to the character "V" if eigenvectors
71 * should be calculated, otherwise "N" (only eigenvalues).
73 jobz = (eigenvectors != nullptr) ? "V" : "N";
75 /* allocate lapack stuff */
79 /* First time we ask the routine how much workspace it needs */
84 /* Convert indices to fortran standard */
88 /* Call LAPACK routine using fortran interface. Note that we use upper storage,
89 * but this corresponds to lower storage ("L") in Fortran.
92 F77_FUNC(dsyevr, DSYEVR)
93 (jobz, "I", "L", &n, a, &n, &vl, &vu, &index_lower, &index_upper, &abstol, &m, eigenvalues,
94 eigenvectors, &n, isuppz, &w0, &lwork, &iw0, &liwork, &info);
96 F77_FUNC(ssyevr, SSYEVR)
97 (jobz, "I", "L", &n, a, &n, &vl, &vu, &index_lower, &index_upper, &abstol, &m, eigenvalues,
98 eigenvectors, &n, isuppz, &w0, &lwork, &iw0, &liwork, &info);
104 gmx_fatal(FARGS, "Internal errror in LAPACK diagonalization.");
107 lwork = static_cast<int>(w0);
116 F77_FUNC(dsyevr, DSYEVR)
117 (jobz, "I", "L", &n, a, &n, &vl, &vu, &index_lower, &index_upper, &abstol, &m, eigenvalues,
118 eigenvectors, &n, isuppz, work, &lwork, iwork, &liwork, &info);
120 F77_FUNC(ssyevr, SSYEVR)
121 (jobz, "I", "L", &n, a, &n, &vl, &vu, &index_lower, &index_upper, &abstol, &m, eigenvalues,
122 eigenvectors, &n, isuppz, work, &lwork, iwork, &liwork, &info);
131 gmx_fatal(FARGS, "Internal errror in LAPACK diagonalization.");
137 void sparse_parallel_eigensolver(gmx_sparsematrix_t* A, int neig, real* eigenvalues, real* eigenvectors, int maxiter)
147 int ido, info, lworkl, i, ncv, dovec;
153 MPI_Comm_size(MPI_COMM_WORLD, &nnodes);
154 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
156 if (eigenvectors != NULL)
173 for (i = 0; i < 11; i++)
175 iparam[i] = ipntr[i] = 0;
178 iparam[0] = 1; /* Don't use explicit shifts */
179 iparam[2] = maxiter; /* Max number of iterations */
180 iparam[6] = 1; /* Standard symmetric eigenproblem */
182 lworkl = ncv * (8 + ncv);
184 snew(workd, (3 * n + 4));
189 /* Use machine tolerance - roughly 1e-16 in double precision */
193 fprintf(stderr, "Calculation Ritz values and Lanczos vectors, max %d iterations...\n", maxiter);
199 F77_FUNC(pdsaupd, PDSAUPD)
200 (&ido, "I", &n, "SA", &neig, &abstol, resid, &ncv, v, &n, iparam, ipntr, workd, iwork,
201 workl, &lworkl, &info);
203 F77_FUNC(pssaupd, PSSAUPD)
204 (&ido, "I", &n, "SA", &neig, &abstol, resid, &ncv, v, &n, iparam, ipntr, workd, iwork,
205 workl, &lworkl, &info);
207 if (ido == -1 || ido == 1)
209 gmx_sparsematrix_vector_multiply(A, workd + ipntr[0] - 1, workd + ipntr[1] - 1);
212 fprintf(stderr, "\rIteration %4d: %3d out of %3d Ritz values converged.", iter++, iparam[4], neig);
214 } while (info == 0 && (ido == -1 || ido == 1));
216 fprintf(stderr, "\n");
220 "Maximum number of iterations (%d) reached in Arnoldi\n"
221 "diagonalization, but only %d of %d eigenvectors converged.\n",
222 maxiter, iparam[4], neig);
226 gmx_fatal(FARGS, "Unspecified error from Arnoldi diagonalization:%d\n", info);
230 /* Extract eigenvalues and vectors from data */
231 fprintf(stderr, "Calculating eigenvalues and eigenvectors...\n");
234 F77_FUNC(pdseupd, PDSEUPD)
235 (&dovec, "A", select, eigenvalues, eigenvectors, &n, NULL, "I", &n, "SA", &neig, &abstol, resid,
236 &ncv, v, &n, iparam, ipntr, workd, workl, &lworkl, &info);
238 F77_FUNC(psseupd, PSSEUPD)
239 (&dovec, "A", select, eigenvalues, eigenvectors, &n, NULL, "I", &n, "SA", &neig, &abstol, resid,
240 &ncv, v, &n, iparam, ipntr, workd, workl, &lworkl, &info);
252 void sparse_eigensolver(gmx_sparsematrix_t* A, int neig, real* eigenvalues, real* eigenvectors, int maxiter)
262 int ido, info, lworkl, i, ncv, dovec;
268 MPI_Comm_size(MPI_COMM_WORLD, &n);
271 sparse_parallel_eigensolver(A, neig, eigenvalues, eigenvectors, maxiter);
276 if (eigenvectors != nullptr)
293 for (i = 0; i < 11; i++)
295 iparam[i] = ipntr[i] = 0;
298 iparam[0] = 1; /* Don't use explicit shifts */
299 iparam[2] = maxiter; /* Max number of iterations */
300 iparam[6] = 1; /* Standard symmetric eigenproblem */
302 lworkl = ncv * (8 + ncv);
304 snew(workd, (3 * n + 4));
309 /* Use machine tolerance - roughly 1e-16 in double precision */
313 fprintf(stderr, "Calculation Ritz values and Lanczos vectors, max %d iterations...\n", maxiter);
319 F77_FUNC(dsaupd, DSAUPD)
320 (&ido, "I", &n, "SA", &neig, &abstol, resid, &ncv, v, &n, iparam, ipntr, workd, iwork,
321 workl, &lworkl, &info);
323 F77_FUNC(ssaupd, SSAUPD)
324 (&ido, "I", &n, "SA", &neig, &abstol, resid, &ncv, v, &n, iparam, ipntr, workd, iwork,
325 workl, &lworkl, &info);
327 if (ido == -1 || ido == 1)
329 gmx_sparsematrix_vector_multiply(A, workd + ipntr[0] - 1, workd + ipntr[1] - 1);
332 fprintf(stderr, "\rIteration %4d: %3d out of %3d Ritz values converged.", iter++, iparam[4], neig);
334 } while (info == 0 && (ido == -1 || ido == 1));
336 fprintf(stderr, "\n");
340 "Maximum number of iterations (%d) reached in Arnoldi\n"
341 "diagonalization, but only %d of %d eigenvectors converged.\n",
342 maxiter, iparam[4], neig);
346 gmx_fatal(FARGS, "Unspecified error from Arnoldi diagonalization:%d\n", info);
350 /* Extract eigenvalues and vectors from data */
351 fprintf(stderr, "Calculating eigenvalues and eigenvectors...\n");
354 F77_FUNC(dseupd, DSEUPD)
355 (&dovec, "A", select, eigenvalues, eigenvectors, &n, nullptr, "I", &n, "SA", &neig, &abstol,
356 resid, &ncv, v, &n, iparam, ipntr, workd, workl, &lworkl, &info);
358 F77_FUNC(sseupd, SSEUPD)
359 (&dovec, "A", select, eigenvalues, eigenvectors, &n, nullptr, "I", &n, "SA", &neig, &abstol,
360 resid, &ncv, v, &n, iparam, ipntr, workd, workl, &lworkl, &info);