Add conserved quantity for Berendsen P-couple
[gromacs.git] / src / gromacs / linearalgebra / gmx_arpack.h
blobdcf556fe0c06c3c85ccdd336c0447a17bc629228
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, 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 #ifdef __cplusplus
55 extern "C" {
56 #endif
58 /*! \brief Implicitly Restarted Arnoldi Iteration, double precision.
60 * Reverse communication interface for the Implicitly Restarted Arnoldi
61 * Iteration. For symmetric problems this reduces to a variant of the
62 * Lanczos method. See the ARPACK site for details.
64 * \param ido Reverse communication flag. Set to 0 first time.
65 * Upon return with ido=-1 or ido=1 you should calculate
66 * Y=A*X and recall the routine. Return with ido=2 means
67 * Y=B*X should be calculated. ipntr[0] is the pointer in
68 * workd for X, ipntr[1] is the index for Y.
69 * Return with ido=99 means it finished.
70 * \param bmat 'I' for standard eigenproblem, 'G' for generalized.
71 * \param n Order of eigenproblem.
72 * \param which Which eigenvalues to calculate. 'LA' for largest
73 * algebraic, 'SA' for smallest algebraic, 'LM' for largest
74 * magnitude, 'SM' for smallest magnitude, and finally
75 * 'BE' (both ends) to calculate half from each end of
76 * the spectrum.
77 * \param nev Number of eigenvalues to calculate. 0<nev<n.
78 * \param tol Tolerance. Machine precision of it is 0.
79 * \param resid Optional starting residual vector at input if info=1,
80 * otherwise a random one is used. Final residual vector on
81 * return.
82 * \param ncv Number of columns in matrix v.
83 * \param v N*NCV matrix. V contain the Lanczos basis vectors.
84 * \param ldv Leading dimension of v.
85 * \param iparam Integer array, size 11. Same contents as arpack.
86 * \param ipntr Integer array, size 11. Points to starting locations
87 * in the workd/workl arrays. Same contents as arpack.
88 * \param workd Double precision work array, length 3*n+4.
89 * Provide the same array for all calls, and don't touch it.
90 * IMPORTANT: This is 4 units larger than standard ARPACK!
91 * \param iwork Integer work array, size 80.
92 * Provide the same array for all calls, and don't touch it.
93 * IMPORTANT: New argument compared to standard ARPACK!
94 * \param workl Double precision work array, length lwork.
95 * \param lworkl Length of the work array workl. Must be at least ncv*(ncv+8)
96 * \param info Set info to 0 to use random initial residual vector,
97 * or to 1 if you provide a one. On output, info=0 means
98 * normal exit, 1 that max number of iterations was reached,
99 * and 3 that no shifts could be applied. Negative numbers
100 * correspond to errors in the arguments provided.
102 void
103 F77_FUNC(dsaupd, DSAUPD) (int * ido,
104 const char * bmat,
105 int * n,
106 const char * which,
107 int * nev,
108 double * tol,
109 double * resid,
110 int * ncv,
111 double * v,
112 int * ldv,
113 int * iparam,
114 int * ipntr,
115 double * workd,
116 int * iwork,
117 double * workl,
118 int * lworkl,
119 int * info);
123 /*! \brief Get eigenvalues/vectors after Arnoldi iteration, double prec.
125 * See the ARPACK site for details. You must have finished the interative
126 * part with dsaupd() before calling this function.
128 * \param rvec 1 if you want eigenvectors, 0 if not.
129 * \param howmny 'A' if you want all nvec vectors, 'S' if you
130 * provide a subset selection in select[].
131 * \param select Integer array, dimension nev. Indices of the
132 * eigenvectors to calculate. Fortran code means we
133 * start counting on 1. This array must be given even in
134 * howmny is 'A'. (Arpack documentation is wrong on this).
135 * \param d Double precision array, length nev. Eigenvalues.
136 * \param z Double precision array, n*nev. Eigenvectors.
137 * \param ldz Leading dimension of z. Normally n.
138 * \param sigma Shift if iparam[6] is 3,4, or 5. Ignored otherwise.
139 * \param bmat Provide the same argument as you did to dsaupd()
140 * \param n Provide the same argument as you did to dsaupd()
141 * \param which Provide the same argument as you did to dsaupd()
142 * \param nev Provide the same argument as you did to dsaupd()
143 * \param tol Provide the same argument as you did to dsaupd()
144 * \param resid Provide the same argument as you did to dsaupd()
145 * The array must not be touched between the two function calls!
146 * \param ncv Provide the same argument as you did to dsaupd()
147 * \param v Provide the same argument as you did to dsaupd()
148 * The array must not be touched between the two function calls!
149 * \param ldv Provide the same argument as you did to dsaupd()
150 * \param iparam Provide the same argument as you did to dsaupd()
151 * The array must not be touched between the two function calls!
152 * \param ipntr Provide the same argument as you did to dsaupd()
153 * The array must not be touched between the two function calls!
154 * \param workd Provide the same argument as you did to dsaupd()
155 * The array must not be touched between the two function calls!
156 * \param workl Double precision work array, length lwork.
157 * The array must not be touched between the two function calls!
158 * \param lworkl Provide the same argument as you did to dsaupd()
159 * \param info Provide the same argument as you did to dsaupd()
161 void
162 F77_FUNC(dseupd, DSEUPD) (int * rvec,
163 const char * howmny,
164 int * select,
165 double * d,
166 double * z,
167 int * ldz,
168 double * sigma,
169 const char * bmat,
170 int * n,
171 const char * which,
172 int * nev,
173 double * tol,
174 double * resid,
175 int * ncv,
176 double * v,
177 int * ldv,
178 int * iparam,
179 int * ipntr,
180 double * workd,
181 double * workl,
182 int * lworkl,
183 int * info);
189 /*! \brief Implicitly Restarted Arnoldi Iteration, single precision.
191 * Reverse communication interface for the Implicitly Restarted Arnoldi
192 * Iteration. For symmetric problems this reduces to a variant of the
193 * Lanczos method. See the ARPACK site for details.
195 * \param ido Reverse communication flag. Set to 0 first time.
196 * Upon return with ido=-1 or ido=1 you should calculate
197 * Y=A*X and recall the routine. Return with ido=2 means
198 * Y=B*X should be calculated. ipntr[0] is the pointer in
199 * workd for X, ipntr[1] is the index for Y.
200 * Return with ido=99 means it finished.
201 * \param bmat 'I' for standard eigenproblem, 'G' for generalized.
202 * \param n Order of eigenproblem.
203 * \param which Which eigenvalues to calculate. 'LA' for largest
204 * algebraic, 'SA' for smallest algebraic, 'LM' for largest
205 * magnitude, 'SM' for smallest magnitude, and finally
206 * 'BE' (both ends) to calculate half from each end of
207 * the spectrum.
208 * \param nev Number of eigenvalues to calculate. 0<nev<n.
209 * \param tol Tolerance. Machine precision of it is 0.
210 * \param resid Optional starting residual vector at input if info=1,
211 * otherwise a random one is used. Final residual vector on
212 * return.
213 * \param ncv Number of columns in matrix v.
214 * \param v N*NCV matrix. V contain the Lanczos basis vectors.
215 * \param ldv Leading dimension of v.
216 * \param iparam Integer array, size 11. Same contents as arpack.
217 * \param ipntr Integer array, size 11. Points to starting locations
218 * in the workd/workl arrays. Same contents as arpack.
219 * \param workd Single precision work array, length 3*n+4.
220 * Provide the same array for all calls, and don't touch it.
221 * IMPORTANT: This is 4 units larger than standard ARPACK!
222 * \param iwork Integer work array, size 80.
223 * Provide the same array for all calls, and don't touch it.
224 * IMPORTANT: New argument compared to standard ARPACK!
225 * \param workl Single precision work array, length lwork.
226 * \param lworkl Length of the work array workl. Must be at least ncv*(ncv+8)
227 * \param info Set info to 0 to use random initial residual vector,
228 * or to 1 if you provide a one. On output, info=0 means
229 * normal exit, 1 that max number of iterations was reached,
230 * and 3 that no shifts could be applied. Negative numbers
231 * correspond to errors in the arguments provided.
233 void
234 F77_FUNC(ssaupd, SSAUPD) (int * ido,
235 const char * bmat,
236 int * n,
237 const char * which,
238 int * nev,
239 float * tol,
240 float * resid,
241 int * ncv,
242 float * v,
243 int * ldv,
244 int * iparam,
245 int * ipntr,
246 float * workd,
247 int * iwork,
248 float * workl,
249 int * lworkl,
250 int * info);
256 /*! \brief Get eigenvalues/vectors after Arnoldi iteration, single prec.
258 * See the ARPACK site for details. You must have finished the interative
259 * part with ssaupd() before calling this function.
261 * \param rvec 1 if you want eigenvectors, 0 if not.
262 * \param howmny 'A' if you want all nvec vectors, 'S' if you
263 * provide a subset selection in select[].
264 * \param select Integer array, dimension nev. Indices of the
265 * eigenvectors to calculate. Fortran code means we
266 * start counting on 1. This array must be given even in
267 * howmny is 'A'. (Arpack documentation is wrong on this).
268 * \param d Single precision array, length nev. Eigenvalues.
269 * \param z Single precision array, n*nev. Eigenvectors.
270 * \param ldz Leading dimension of z. Normally n.
271 * \param sigma Shift if iparam[6] is 3,4, or 5. Ignored otherwise.
272 * \param bmat Provide the same argument as you did to ssaupd()
273 * \param n Provide the same argument as you did to ssaupd()
274 * \param which Provide the same argument as you did to ssaupd()
275 * \param nev Provide the same argument as you did to ssaupd()
276 * \param tol Provide the same argument as you did to ssaupd()
277 * \param resid Provide the same argument as you did to ssaupd()
278 * The array must not be touched between the two function calls!
279 * \param ncv Provide the same argument as you did to ssaupd()
280 * \param v Provide the same argument as you did to ssaupd()
281 * The array must not be touched between the two function calls!
282 * \param ldv Provide the same argument as you did to ssaupd()
283 * \param iparam Provide the same argument as you did to ssaupd()
284 * The array must not be touched between the two function calls!
285 * \param ipntr Provide the same argument as you did to ssaupd()
286 * The array must not be touched between the two function calls!
287 * \param workd Provide the same argument as you did to ssaupd()
288 * The array must not be touched between the two function calls!
289 * \param workl Single precision work array, length lwork.
290 * The array must not be touched between the two function calls!
291 * \param lworkl Provide the same argument as you did to ssaupd()
292 * \param info Provide the same argument as you did to ssaupd()
294 void
295 F77_FUNC(sseupd, SSEUPD) (int * rvec,
296 const char * howmny,
297 int * select,
298 float * d,
299 float * z,
300 int * ldz,
301 float * sigma,
302 const char * bmat,
303 int * n,
304 const char * which,
305 int * nev,
306 float * tol,
307 float * resid,
308 int * ncv,
309 float * v,
310 int * ldv,
311 int * iparam,
312 int * ipntr,
313 float * workd,
314 float * workl,
315 int * lworkl,
316 int * info);
318 #ifdef __cplusplus
320 #endif
322 #endif