Code beautification with uncrustify
[alexxy/gromacs.git] / src / gromacs / gmxlib / nrjac.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <stdio.h>
41 #include <stdlib.h>
42 #include <math.h>
43 #include "gstat.h"
44 #include "smalloc.h"
45 #include "gmx_fatal.h"
46 #include "nrjac.h"
47
48
49 static inline
50 void do_rotate(double **a, int i, int j, int k, int l, double tau, double s)
51 {
52     double g, h;
53     g       = a[i][j];
54     h       = a[k][l];
55     a[i][j] = g - s * (h + g * tau);
56     a[k][l] = h + s * (g - h * tau);
57 }
58
59 void jacobi(double **a, int n, double d[], double **v, int *nrot)
60 {
61     int    j, i;
62     int    iq, ip;
63     double tresh, theta, tau, t, sm, s, h, g, c, *b, *z;
64
65     snew(b, n);
66     snew(z, n);
67     for (ip = 0; ip < n; ip++)
68     {
69         for (iq = 0; iq < n; iq++)
70         {
71             v[ip][iq] = 0.0;
72         }
73         v[ip][ip] = 1.0;
74     }
75     for (ip = 0; ip < n; ip++)
76     {
77         b[ip] = d[ip] = a[ip][ip];
78         z[ip] = 0.0;
79     }
80     *nrot = 0;
81     for (i = 1; i <= 50; i++)
82     {
83         sm = 0.0;
84         for (ip = 0; ip < n-1; ip++)
85         {
86             for (iq = ip+1; iq < n; iq++)
87             {
88                 sm += fabs(a[ip][iq]);
89             }
90         }
91         if (sm == 0.0)
92         {
93             sfree(z);
94             sfree(b);
95             return;
96         }
97         if (i < 4)
98         {
99             tresh = 0.2*sm/(n*n);
100         }
101         else
102         {
103             tresh = 0.0;
104         }
105         for (ip = 0; ip < n-1; ip++)
106         {
107             for (iq = ip+1; iq < n; iq++)
108             {
109                 g = 100.0*fabs(a[ip][iq]);
110                 if (i > 4 && fabs(d[ip])+g == fabs(d[ip])
111                     && fabs(d[iq])+g == fabs(d[iq]))
112                 {
113                     a[ip][iq] = 0.0;
114                 }
115                 else if (fabs(a[ip][iq]) > tresh)
116                 {
117                     h = d[iq]-d[ip];
118                     if (fabs(h)+g == fabs(h))
119                     {
120                         t = (a[ip][iq])/h;
121                     }
122                     else
123                     {
124                         theta = 0.5*h/(a[ip][iq]);
125                         t     = 1.0/(fabs(theta)+sqrt(1.0+theta*theta));
126                         if (theta < 0.0)
127                         {
128                             t = -t;
129                         }
130                     }
131                     c         = 1.0/sqrt(1+t*t);
132                     s         = t*c;
133                     tau       = s/(1.0+c);
134                     h         = t*a[ip][iq];
135                     z[ip]    -= h;
136                     z[iq]    += h;
137                     d[ip]    -= h;
138                     d[iq]    += h;
139                     a[ip][iq] = 0.0;
140                     for (j = 0; j < ip; j++)
141                     {
142                         do_rotate(a, j, ip, j, iq, tau, s);
143                     }
144                     for (j = ip+1; j < iq; j++)
145                     {
146                         do_rotate(a, ip, j, j, iq, tau, s);
147                     }
148                     for (j = iq+1; j < n; j++)
149                     {
150                         do_rotate(a, ip, j, iq, j, tau, s);
151                     }
152                     for (j = 0; j < n; j++)
153                     {
154                         do_rotate(v, j, ip, j, iq, tau, s);
155                     }
156                     ++(*nrot);
157                 }
158             }
159         }
160         for (ip = 0; ip < n; ip++)
161         {
162             b[ip] +=  z[ip];
163             d[ip]  =  b[ip];
164             z[ip]  =  0.0;
165         }
166     }
167     gmx_fatal(FARGS, "Error: Too many iterations in routine JACOBI\n");
168 }
169
170 int m_inv_gen(real **m, int n, real **minv)
171 {
172     double **md, **v, *eig, tol, s;
173     int      nzero, i, j, k, nrot;
174
175     snew(md, n);
176     for (i = 0; i < n; i++)
177     {
178         snew(md[i], n);
179     }
180     snew(v, n);
181     for (i = 0; i < n; i++)
182     {
183         snew(v[i], n);
184     }
185     snew(eig, n);
186     for (i = 0; i < n; i++)
187     {
188         for (j = 0; j < n; j++)
189         {
190             md[i][j] = m[i][j];
191         }
192     }
193
194     tol = 0;
195     for (i = 0; i < n; i++)
196     {
197         tol += fabs(md[i][i]);
198     }
199     tol = 1e-6*tol/n;
200
201     jacobi(md, n, eig, v, &nrot);
202
203     nzero = 0;
204     for (i = 0; i < n; i++)
205     {
206         if (fabs(eig[i]) < tol)
207         {
208             eig[i] = 0;
209             nzero++;
210         }
211         else
212         {
213             eig[i] = 1.0/eig[i];
214         }
215     }
216
217     for (i = 0; i < n; i++)
218     {
219         for (j = 0; j < n; j++)
220         {
221             s = 0;
222             for (k = 0; k < n; k++)
223             {
224                 s += eig[k]*v[i][k]*v[j][k];
225             }
226             minv[i][j] = s;
227         }
228     }
229
230     sfree(eig);
231     for (i = 0; i < n; i++)
232     {
233         sfree(v[i]);
234     }
235     sfree(v);
236     for (i = 0; i < n; i++)
237     {
238         sfree(md[i]);
239     }
240     sfree(md);
241
242     return nzero;
243 }