Update instructions in containers.rst
[gromacs.git] / src / gromacs / linearalgebra / gmx_arpack.h
blobe7b9511f0c59e2183fc367bcb4657b4faf052e37
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2004 David van der Spoel, Erik Lindahl, University of Groningen.
5 * Copyright (c) 2012,2013,2014,2018,2019, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
36 /*! \internal \file
37 * \brief
38 * Selected routines from ARPACK
40 * This file contains a subset of ARPACK functions to perform
41 * diagonalization and SVD for sparse matrices in Gromacs.
43 * Consult the main ARPACK site for detailed documentation:
44 * http://www.caam.rice.edu/software/ARPACK/
46 * Below, we just list the options and any specific differences
47 * from ARPACK. The code is essentially the same, but the routines
48 * have been made thread-safe by using extra workspace arrays.
50 #ifndef GMX_ARPACK_H
51 #define GMX_ARPACK_H
53 #include "config.h"
54 /*! \brief Implicitly Restarted Arnoldi Iteration, double precision.
56 * Reverse communication interface for the Implicitly Restarted Arnoldi
57 * Iteration. For symmetric problems this reduces to a variant of the
58 * Lanczos method. See the ARPACK site for details.
60 * \param ido Reverse communication flag. Set to 0 first time.
61 * Upon return with ido=-1 or ido=1 you should calculate
62 * Y=A*X and recall the routine. Return with ido=2 means
63 * Y=B*X should be calculated. ipntr[0] is the pointer in
64 * workd for X, ipntr[1] is the index for Y.
65 * Return with ido=99 means it finished.
66 * \param bmat 'I' for standard eigenproblem, 'G' for generalized.
67 * \param n Order of eigenproblem.
68 * \param which Which eigenvalues to calculate. 'LA' for largest
69 * algebraic, 'SA' for smallest algebraic, 'LM' for largest
70 * magnitude, 'SM' for smallest magnitude, and finally
71 * 'BE' (both ends) to calculate half from each end of
72 * the spectrum.
73 * \param nev Number of eigenvalues to calculate. 0<nev<n.
74 * \param tol Tolerance. Machine precision of it is 0.
75 * \param resid Optional starting residual vector at input if info=1,
76 * otherwise a random one is used. Final residual vector on
77 * return.
78 * \param ncv Number of columns in matrix v.
79 * \param v N*NCV matrix. V contain the Lanczos basis vectors.
80 * \param ldv Leading dimension of v.
81 * \param iparam Integer array, size 11. Same contents as arpack.
82 * \param ipntr Integer array, size 11. Points to starting locations
83 * in the workd/workl arrays. Same contents as arpack.
84 * \param workd Double precision work array, length 3*n+4.
85 * Provide the same array for all calls, and don't touch it.
86 * IMPORTANT: This is 4 units larger than standard ARPACK!
87 * \param iwork Integer work array, size 80.
88 * Provide the same array for all calls, and don't touch it.
89 * IMPORTANT: New argument compared to standard ARPACK!
90 * \param workl Double precision work array, length lwork.
91 * \param lworkl Length of the work array workl. Must be at least ncv*(ncv+8)
92 * \param info Set info to 0 to use random initial residual vector,
93 * or to 1 if you provide a one. On output, info=0 means
94 * normal exit, 1 that max number of iterations was reached,
95 * and 3 that no shifts could be applied. Negative numbers
96 * correspond to errors in the arguments provided.
98 void F77_FUNC(dsaupd, DSAUPD)(int* ido,
99 const char* bmat,
100 int* n,
101 const char* which,
102 int* nev,
103 double* tol,
104 double* resid,
105 int* ncv,
106 double* v,
107 int* ldv,
108 int* iparam,
109 int* ipntr,
110 double* workd,
111 int* iwork,
112 double* workl,
113 int* lworkl,
114 int* info);
117 /*! \brief Get eigenvalues/vectors after Arnoldi iteration, double prec.
119 * See the ARPACK site for details. You must have finished the interative
120 * part with dsaupd() before calling this function.
122 * \param rvec 1 if you want eigenvectors, 0 if not.
123 * \param howmny 'A' if you want all nvec vectors, 'S' if you
124 * provide a subset selection in select[].
125 * \param select Integer array, dimension nev. Indices of the
126 * eigenvectors to calculate. Fortran code means we
127 * start counting on 1. This array must be given even in
128 * howmny is 'A'. (Arpack documentation is wrong on this).
129 * \param d Double precision array, length nev. Eigenvalues.
130 * \param z Double precision array, n*nev. Eigenvectors.
131 * \param ldz Leading dimension of z. Normally n.
132 * \param sigma Shift if iparam[6] is 3,4, or 5. Ignored otherwise.
133 * \param bmat Provide the same argument as you did to dsaupd()
134 * \param n Provide the same argument as you did to dsaupd()
135 * \param which Provide the same argument as you did to dsaupd()
136 * \param nev Provide the same argument as you did to dsaupd()
137 * \param tol Provide the same argument as you did to dsaupd()
138 * \param resid Provide the same argument as you did to dsaupd()
139 * The array must not be touched between the two function calls!
140 * \param ncv Provide the same argument as you did to dsaupd()
141 * \param v Provide the same argument as you did to dsaupd()
142 * The array must not be touched between the two function calls!
143 * \param ldv Provide the same argument as you did to dsaupd()
144 * \param iparam Provide the same argument as you did to dsaupd()
145 * The array must not be touched between the two function calls!
146 * \param ipntr Provide the same argument as you did to dsaupd()
147 * The array must not be touched between the two function calls!
148 * \param workd Provide the same argument as you did to dsaupd()
149 * The array must not be touched between the two function calls!
150 * \param workl Double precision work array, length lwork.
151 * The array must not be touched between the two function calls!
152 * \param lworkl Provide the same argument as you did to dsaupd()
153 * \param info Provide the same argument as you did to dsaupd()
155 void F77_FUNC(dseupd, DSEUPD)(int* rvec,
156 const char* howmny,
157 int* select,
158 double* d,
159 double* z,
160 int* ldz,
161 double* sigma,
162 const char* bmat,
163 int* n,
164 const char* which,
165 int* nev,
166 double* tol,
167 double* resid,
168 int* ncv,
169 double* v,
170 int* ldv,
171 int* iparam,
172 int* ipntr,
173 double* workd,
174 double* workl,
175 int* lworkl,
176 int* info);
179 /*! \brief Implicitly Restarted Arnoldi Iteration, single precision.
181 * Reverse communication interface for the Implicitly Restarted Arnoldi
182 * Iteration. For symmetric problems this reduces to a variant of the
183 * Lanczos method. See the ARPACK site for details.
185 * \param ido Reverse communication flag. Set to 0 first time.
186 * Upon return with ido=-1 or ido=1 you should calculate
187 * Y=A*X and recall the routine. Return with ido=2 means
188 * Y=B*X should be calculated. ipntr[0] is the pointer in
189 * workd for X, ipntr[1] is the index for Y.
190 * Return with ido=99 means it finished.
191 * \param bmat 'I' for standard eigenproblem, 'G' for generalized.
192 * \param n Order of eigenproblem.
193 * \param which Which eigenvalues to calculate. 'LA' for largest
194 * algebraic, 'SA' for smallest algebraic, 'LM' for largest
195 * magnitude, 'SM' for smallest magnitude, and finally
196 * 'BE' (both ends) to calculate half from each end of
197 * the spectrum.
198 * \param nev Number of eigenvalues to calculate. 0<nev<n.
199 * \param tol Tolerance. Machine precision of it is 0.
200 * \param resid Optional starting residual vector at input if info=1,
201 * otherwise a random one is used. Final residual vector on
202 * return.
203 * \param ncv Number of columns in matrix v.
204 * \param v N*NCV matrix. V contain the Lanczos basis vectors.
205 * \param ldv Leading dimension of v.
206 * \param iparam Integer array, size 11. Same contents as arpack.
207 * \param ipntr Integer array, size 11. Points to starting locations
208 * in the workd/workl arrays. Same contents as arpack.
209 * \param workd Single precision work array, length 3*n+4.
210 * Provide the same array for all calls, and don't touch it.
211 * IMPORTANT: This is 4 units larger than standard ARPACK!
212 * \param iwork Integer work array, size 80.
213 * Provide the same array for all calls, and don't touch it.
214 * IMPORTANT: New argument compared to standard ARPACK!
215 * \param workl Single precision work array, length lwork.
216 * \param lworkl Length of the work array workl. Must be at least ncv*(ncv+8)
217 * \param info Set info to 0 to use random initial residual vector,
218 * or to 1 if you provide a one. On output, info=0 means
219 * normal exit, 1 that max number of iterations was reached,
220 * and 3 that no shifts could be applied. Negative numbers
221 * correspond to errors in the arguments provided.
223 void F77_FUNC(ssaupd, SSAUPD)(int* ido,
224 const char* bmat,
225 int* n,
226 const char* which,
227 int* nev,
228 float* tol,
229 float* resid,
230 int* ncv,
231 float* v,
232 int* ldv,
233 int* iparam,
234 int* ipntr,
235 float* workd,
236 int* iwork,
237 float* workl,
238 int* lworkl,
239 int* info);
242 /*! \brief Get eigenvalues/vectors after Arnoldi iteration, single prec.
244 * See the ARPACK site for details. You must have finished the interative
245 * part with ssaupd() before calling this function.
247 * \param rvec 1 if you want eigenvectors, 0 if not.
248 * \param howmny 'A' if you want all nvec vectors, 'S' if you
249 * provide a subset selection in select[].
250 * \param select Integer array, dimension nev. Indices of the
251 * eigenvectors to calculate. Fortran code means we
252 * start counting on 1. This array must be given even in
253 * howmny is 'A'. (Arpack documentation is wrong on this).
254 * \param d Single precision array, length nev. Eigenvalues.
255 * \param z Single precision array, n*nev. Eigenvectors.
256 * \param ldz Leading dimension of z. Normally n.
257 * \param sigma Shift if iparam[6] is 3,4, or 5. Ignored otherwise.
258 * \param bmat Provide the same argument as you did to ssaupd()
259 * \param n Provide the same argument as you did to ssaupd()
260 * \param which Provide the same argument as you did to ssaupd()
261 * \param nev Provide the same argument as you did to ssaupd()
262 * \param tol Provide the same argument as you did to ssaupd()
263 * \param resid Provide the same argument as you did to ssaupd()
264 * The array must not be touched between the two function calls!
265 * \param ncv Provide the same argument as you did to ssaupd()
266 * \param v Provide the same argument as you did to ssaupd()
267 * The array must not be touched between the two function calls!
268 * \param ldv Provide the same argument as you did to ssaupd()
269 * \param iparam Provide the same argument as you did to ssaupd()
270 * The array must not be touched between the two function calls!
271 * \param ipntr Provide the same argument as you did to ssaupd()
272 * The array must not be touched between the two function calls!
273 * \param workd Provide the same argument as you did to ssaupd()
274 * The array must not be touched between the two function calls!
275 * \param workl Single precision work array, length lwork.
276 * The array must not be touched between the two function calls!
277 * \param lworkl Provide the same argument as you did to ssaupd()
278 * \param info Provide the same argument as you did to ssaupd()
280 void F77_FUNC(sseupd, SSEUPD)(int* rvec,
281 const char* howmny,
282 int* select,
283 float* d,
284 float* z,
285 int* ldz,
286 float* sigma,
287 const char* bmat,
288 int* n,
289 const char* which,
290 int* nev,
291 float* tol,
292 float* resid,
293 int* ncv,
294 float* v,
295 int* ldv,
296 int* iparam,
297 int* ipntr,
298 float* workd,
299 float* workl,
300 int* lworkl,
301 int* info);
303 #endif