Remove continuation from convert_tpr
[gromacs.git] / src / gromacs / linearalgebra / nrjac.cpp
blob37de48e9510aefc62ac6885aafdee6645d53463b
1 /*
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,2014,2015, 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.
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
40 #include "nrjac.h"
42 #include <cmath>
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/smalloc.h"
47 static gmx_inline
48 void do_rotate(double **a, int i, int j, int k, int l, double tau, double s)
50 double g, h;
51 g = a[i][j];
52 h = a[k][l];
53 a[i][j] = g - s * (h + g * tau);
54 a[k][l] = h + s * (g - h * tau);
57 void jacobi(double **a, int n, double d[], double **v, int *nrot)
59 int j, i;
60 int iq, ip;
61 double tresh, theta, tau, t, sm, s, h, g, c, *b, *z;
63 snew(b, n);
64 snew(z, n);
65 for (ip = 0; ip < n; ip++)
67 for (iq = 0; iq < n; iq++)
69 v[ip][iq] = 0.0;
71 v[ip][ip] = 1.0;
73 for (ip = 0; ip < n; ip++)
75 b[ip] = d[ip] = a[ip][ip];
76 z[ip] = 0.0;
78 *nrot = 0;
79 for (i = 1; i <= 50; i++)
81 sm = 0.0;
82 for (ip = 0; ip < n-1; ip++)
84 for (iq = ip+1; iq < n; iq++)
86 sm += std::abs(a[ip][iq]);
89 if (sm == 0.0)
91 sfree(z);
92 sfree(b);
93 return;
95 if (i < 4)
97 tresh = 0.2*sm/(n*n);
99 else
101 tresh = 0.0;
103 for (ip = 0; ip < n-1; ip++)
105 for (iq = ip+1; iq < n; iq++)
107 g = 100.0*std::abs(a[ip][iq]);
108 if (i > 4 && std::abs(d[ip])+g == std::abs(d[ip])
109 && std::abs(d[iq])+g == std::abs(d[iq]))
111 a[ip][iq] = 0.0;
113 else if (std::abs(a[ip][iq]) > tresh)
115 h = d[iq]-d[ip];
116 if (std::abs(h)+g == std::abs(h))
118 t = (a[ip][iq])/h;
120 else
122 theta = 0.5*h/(a[ip][iq]);
123 t = 1.0/(std::abs(theta)+std::sqrt(1.0+theta*theta));
124 if (theta < 0.0)
126 t = -t;
129 c = 1.0/std::sqrt(1+t*t);
130 s = t*c;
131 tau = s/(1.0+c);
132 h = t*a[ip][iq];
133 z[ip] -= h;
134 z[iq] += h;
135 d[ip] -= h;
136 d[iq] += h;
137 a[ip][iq] = 0.0;
138 for (j = 0; j < ip; j++)
140 do_rotate(a, j, ip, j, iq, tau, s);
142 for (j = ip+1; j < iq; j++)
144 do_rotate(a, ip, j, j, iq, tau, s);
146 for (j = iq+1; j < n; j++)
148 do_rotate(a, ip, j, iq, j, tau, s);
150 for (j = 0; j < n; j++)
152 do_rotate(v, j, ip, j, iq, tau, s);
154 ++(*nrot);
158 for (ip = 0; ip < n; ip++)
160 b[ip] += z[ip];
161 d[ip] = b[ip];
162 z[ip] = 0.0;
165 gmx_fatal(FARGS, "Error: Too many iterations in routine JACOBI\n");
168 int m_inv_gen(real **m, int n, real **minv)
170 double **md, **v, *eig, tol, s;
171 int nzero, i, j, k, nrot;
173 snew(md, n);
174 for (i = 0; i < n; i++)
176 snew(md[i], n);
178 snew(v, n);
179 for (i = 0; i < n; i++)
181 snew(v[i], n);
183 snew(eig, n);
184 for (i = 0; i < n; i++)
186 for (j = 0; j < n; j++)
188 md[i][j] = m[i][j];
192 tol = 0;
193 for (i = 0; i < n; i++)
195 tol += std::abs(md[i][i]);
197 tol = 1e-6*tol/n;
199 jacobi(md, n, eig, v, &nrot);
201 nzero = 0;
202 for (i = 0; i < n; i++)
204 if (std::abs(eig[i]) < tol)
206 eig[i] = 0;
207 nzero++;
209 else
211 eig[i] = 1.0/eig[i];
215 for (i = 0; i < n; i++)
217 for (j = 0; j < n; j++)
219 s = 0;
220 for (k = 0; k < n; k++)
222 s += eig[k]*v[i][k]*v[j][k];
224 minv[i][j] = s;
228 sfree(eig);
229 for (i = 0; i < n; i++)
231 sfree(v[i]);
233 sfree(v);
234 for (i = 0; i < n; i++)
236 sfree(md[i]);
238 sfree(md);
240 return nzero;