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, 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_lapack.h"
47 #include "gmx_arpack.h"
59 int il, iu, m, iw0, info;
76 /* Make jobz point to the character "V" if eigenvectors
77 * should be calculated, otherwise "N" (only eigenvalues).
79 jobz = (eigenvectors != NULL) ? "V" : "N";
81 /* allocate lapack stuff */
85 /* First time we ask the routine how much workspace it needs */
90 /* Convert indices to fortran standard */
94 /* Call LAPACK routine using fortran interface. Note that we use upper storage,
95 * but this corresponds to lower storage ("L") in Fortran.
98 F77_FUNC(dsyevr, DSYEVR) (jobz, "I", "L", &n, a, &n, &vl, &vu, &index_lower, &index_upper,
99 &abstol, &m, eigenvalues, eigenvectors, &n,
100 isuppz, &w0, &lwork, &iw0, &liwork, &info);
102 F77_FUNC(ssyevr, SSYEVR) (jobz, "I", "L", &n, a, &n, &vl, &vu, &index_lower, &index_upper,
103 &abstol, &m, eigenvalues, eigenvectors, &n,
104 isuppz, &w0, &lwork, &iw0, &liwork, &info);
110 gmx_fatal(FARGS, "Internal errror in LAPACK diagonalization.");
122 F77_FUNC(dsyevr, DSYEVR) (jobz, "I", "L", &n, a, &n, &vl, &vu, &index_lower, &index_upper,
123 &abstol, &m, eigenvalues, eigenvectors, &n,
124 isuppz, work, &lwork, iwork, &liwork, &info);
126 F77_FUNC(ssyevr, SSYEVR) (jobz, "I", "L", &n, a, &n, &vl, &vu, &index_lower, &index_upper,
127 &abstol, &m, eigenvalues, eigenvectors, &n,
128 isuppz, work, &lwork, iwork, &liwork, &info);
137 gmx_fatal(FARGS, "Internal errror in LAPACK diagonalization.");
145 sparse_parallel_eigensolver(gmx_sparsematrix_t * A,
159 int ido, info, lworkl, i, ncv, dovec;
165 MPI_Comm_size( MPI_COMM_WORLD, &nnodes );
166 MPI_Comm_rank( MPI_COMM_WORLD, &rank );
168 if (eigenvectors != NULL)
185 for (i = 0; i < 11; i++)
187 iparam[i] = ipntr[i] = 0;
190 iparam[0] = 1; /* Don't use explicit shifts */
191 iparam[2] = maxiter; /* Max number of iterations */
192 iparam[6] = 1; /* Standard symmetric eigenproblem */
194 lworkl = ncv*(8+ncv);
196 snew(workd, (3*n+4));
201 /* Use machine tolerance - roughly 1e-16 in double precision */
205 fprintf(stderr, "Calculation Ritz values and Lanczos vectors, max %d iterations...\n", maxiter);
211 F77_FUNC(pdsaupd, PDSAUPD) (&ido, "I", &n, "SA", &neig, &abstol,
212 resid, &ncv, v, &n, iparam, ipntr,
213 workd, iwork, workl, &lworkl, &info);
215 F77_FUNC(pssaupd, PSSAUPD) (&ido, "I", &n, "SA", &neig, &abstol,
216 resid, &ncv, v, &n, iparam, ipntr,
217 workd, iwork, workl, &lworkl, &info);
219 if (ido == -1 || ido == 1)
221 gmx_sparsematrix_vector_multiply(A, workd+ipntr[0]-1, workd+ipntr[1]-1);
224 fprintf(stderr, "\rIteration %4d: %3d out of %3d Ritz values converged.", iter++, iparam[4], neig);
226 while (info == 0 && (ido == -1 || ido == 1));
228 fprintf(stderr, "\n");
232 "Maximum number of iterations (%d) reached in Arnoldi\n"
233 "diagonalization, but only %d of %d eigenvectors converged.\n",
234 maxiter, iparam[4], neig);
238 gmx_fatal(FARGS, "Unspecified error from Arnoldi diagonalization:%d\n", info);
242 /* Extract eigenvalues and vectors from data */
243 fprintf(stderr, "Calculating eigenvalues and eigenvectors...\n");
246 F77_FUNC(pdseupd, PDSEUPD) (&dovec, "A", select, eigenvalues, eigenvectors,
247 &n, NULL, "I", &n, "SA", &neig, &abstol,
248 resid, &ncv, v, &n, iparam, ipntr,
249 workd, workl, &lworkl, &info);
251 F77_FUNC(psseupd, PSSEUPD) (&dovec, "A", select, eigenvalues, eigenvectors,
252 &n, NULL, "I", &n, "SA", &neig, &abstol,
253 resid, &ncv, v, &n, iparam, ipntr,
254 workd, workl, &lworkl, &info);
267 sparse_eigensolver(gmx_sparsematrix_t * A,
281 int ido, info, lworkl, i, ncv, dovec;
287 MPI_Comm_size( MPI_COMM_WORLD, &n );
290 sparse_parallel_eigensolver(A, neig, eigenvalues, eigenvectors, maxiter);
295 if (eigenvectors != NULL)
312 for (i = 0; i < 11; i++)
314 iparam[i] = ipntr[i] = 0;
317 iparam[0] = 1; /* Don't use explicit shifts */
318 iparam[2] = maxiter; /* Max number of iterations */
319 iparam[6] = 1; /* Standard symmetric eigenproblem */
321 lworkl = ncv*(8+ncv);
323 snew(workd, (3*n+4));
328 /* Use machine tolerance - roughly 1e-16 in double precision */
332 fprintf(stderr, "Calculation Ritz values and Lanczos vectors, max %d iterations...\n", maxiter);
338 F77_FUNC(dsaupd, DSAUPD) (&ido, "I", &n, "SA", &neig, &abstol,
339 resid, &ncv, v, &n, iparam, ipntr,
340 workd, iwork, workl, &lworkl, &info);
342 F77_FUNC(ssaupd, SSAUPD) (&ido, "I", &n, "SA", &neig, &abstol,
343 resid, &ncv, v, &n, iparam, ipntr,
344 workd, iwork, workl, &lworkl, &info);
346 if (ido == -1 || ido == 1)
348 gmx_sparsematrix_vector_multiply(A, workd+ipntr[0]-1, workd+ipntr[1]-1);
351 fprintf(stderr, "\rIteration %4d: %3d out of %3d Ritz values converged.", iter++, iparam[4], neig);
353 while (info == 0 && (ido == -1 || ido == 1));
355 fprintf(stderr, "\n");
359 "Maximum number of iterations (%d) reached in Arnoldi\n"
360 "diagonalization, but only %d of %d eigenvectors converged.\n",
361 maxiter, iparam[4], neig);
365 gmx_fatal(FARGS, "Unspecified error from Arnoldi diagonalization:%d\n", info);
369 /* Extract eigenvalues and vectors from data */
370 fprintf(stderr, "Calculating eigenvalues and eigenvectors...\n");
373 F77_FUNC(dseupd, DSEUPD) (&dovec, "A", select, eigenvalues, eigenvectors,
374 &n, NULL, "I", &n, "SA", &neig, &abstol,
375 resid, &ncv, v, &n, iparam, ipntr,
376 workd, workl, &lworkl, &info);
378 F77_FUNC(sseupd, SSEUPD) (&dovec, "A", select, eigenvalues, eigenvectors,
379 &n, NULL, "I", &n, "SA", &neig, &abstol,
380 resid, &ncv, v, &n, iparam, ipntr,
381 workd, workl, &lworkl, &info);