Merge branch release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / linearalgebra / nrjac.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2012,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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "nrjac.h"
39
40 #include <math.h>
41
42 #include "gmx_fatal.h"
43 #include "gromacs/utility/smalloc.h"
44
45 static gmx_inline
46 void do_rotate(double **a, int i, int j, int k, int l, double tau, double s)
47 {
48     double g, h;
49     g       = a[i][j];
50     h       = a[k][l];
51     a[i][j] = g - s * (h + g * tau);
52     a[k][l] = h + s * (g - h * tau);
53 }
54
55 void jacobi(double **a, int n, double d[], double **v, int *nrot)
56 {
57     int    j, i;
58     int    iq, ip;
59     double tresh, theta, tau, t, sm, s, h, g, c, *b, *z;
60
61     snew(b, n);
62     snew(z, n);
63     for (ip = 0; ip < n; ip++)
64     {
65         for (iq = 0; iq < n; iq++)
66         {
67             v[ip][iq] = 0.0;
68         }
69         v[ip][ip] = 1.0;
70     }
71     for (ip = 0; ip < n; ip++)
72     {
73         b[ip] = d[ip] = a[ip][ip];
74         z[ip] = 0.0;
75     }
76     *nrot = 0;
77     for (i = 1; i <= 50; i++)
78     {
79         sm = 0.0;
80         for (ip = 0; ip < n-1; ip++)
81         {
82             for (iq = ip+1; iq < n; iq++)
83             {
84                 sm += fabs(a[ip][iq]);
85             }
86         }
87         if (sm == 0.0)
88         {
89             sfree(z);
90             sfree(b);
91             return;
92         }
93         if (i < 4)
94         {
95             tresh = 0.2*sm/(n*n);
96         }
97         else
98         {
99             tresh = 0.0;
100         }
101         for (ip = 0; ip < n-1; ip++)
102         {
103             for (iq = ip+1; iq < n; iq++)
104             {
105                 g = 100.0*fabs(a[ip][iq]);
106                 if (i > 4 && fabs(d[ip])+g == fabs(d[ip])
107                     && fabs(d[iq])+g == fabs(d[iq]))
108                 {
109                     a[ip][iq] = 0.0;
110                 }
111                 else if (fabs(a[ip][iq]) > tresh)
112                 {
113                     h = d[iq]-d[ip];
114                     if (fabs(h)+g == fabs(h))
115                     {
116                         t = (a[ip][iq])/h;
117                     }
118                     else
119                     {
120                         theta = 0.5*h/(a[ip][iq]);
121                         t     = 1.0/(fabs(theta)+sqrt(1.0+theta*theta));
122                         if (theta < 0.0)
123                         {
124                             t = -t;
125                         }
126                     }
127                     c         = 1.0/sqrt(1+t*t);
128                     s         = t*c;
129                     tau       = s/(1.0+c);
130                     h         = t*a[ip][iq];
131                     z[ip]    -= h;
132                     z[iq]    += h;
133                     d[ip]    -= h;
134                     d[iq]    += h;
135                     a[ip][iq] = 0.0;
136                     for (j = 0; j < ip; j++)
137                     {
138                         do_rotate(a, j, ip, j, iq, tau, s);
139                     }
140                     for (j = ip+1; j < iq; j++)
141                     {
142                         do_rotate(a, ip, j, j, iq, tau, s);
143                     }
144                     for (j = iq+1; j < n; j++)
145                     {
146                         do_rotate(a, ip, j, iq, j, tau, s);
147                     }
148                     for (j = 0; j < n; j++)
149                     {
150                         do_rotate(v, j, ip, j, iq, tau, s);
151                     }
152                     ++(*nrot);
153                 }
154             }
155         }
156         for (ip = 0; ip < n; ip++)
157         {
158             b[ip] +=  z[ip];
159             d[ip]  =  b[ip];
160             z[ip]  =  0.0;
161         }
162     }
163     gmx_fatal(FARGS, "Error: Too many iterations in routine JACOBI\n");
164 }
165
166 int m_inv_gen(real **m, int n, real **minv)
167 {
168     double **md, **v, *eig, tol, s;
169     int      nzero, i, j, k, nrot;
170
171     snew(md, n);
172     for (i = 0; i < n; i++)
173     {
174         snew(md[i], n);
175     }
176     snew(v, n);
177     for (i = 0; i < n; i++)
178     {
179         snew(v[i], n);
180     }
181     snew(eig, n);
182     for (i = 0; i < n; i++)
183     {
184         for (j = 0; j < n; j++)
185         {
186             md[i][j] = m[i][j];
187         }
188     }
189
190     tol = 0;
191     for (i = 0; i < n; i++)
192     {
193         tol += fabs(md[i][i]);
194     }
195     tol = 1e-6*tol/n;
196
197     jacobi(md, n, eig, v, &nrot);
198
199     nzero = 0;
200     for (i = 0; i < n; i++)
201     {
202         if (fabs(eig[i]) < tol)
203         {
204             eig[i] = 0;
205             nzero++;
206         }
207         else
208         {
209             eig[i] = 1.0/eig[i];
210         }
211     }
212
213     for (i = 0; i < n; i++)
214     {
215         for (j = 0; j < n; j++)
216         {
217             s = 0;
218             for (k = 0; k < n; k++)
219             {
220                 s += eig[k]*v[i][k]*v[j][k];
221             }
222             minv[i][j] = s;
223         }
224     }
225
226     sfree(eig);
227     for (i = 0; i < n; i++)
228     {
229         sfree(v[i]);
230     }
231     sfree(v);
232     for (i = 0; i < n; i++)
233     {
234         sfree(md[i]);
235     }
236     sfree(md);
237
238     return nzero;
239 }