Used IWYU to partially clean up some includes
[alexxy/gromacs.git] / src / gromacs / math / do_fit.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 #include "gmxpre.h"
38
39 #include "do_fit.h"
40
41 #include <math.h>
42 #include <stdio.h>
43
44 #include "gromacs/linearalgebra/nrjac.h"
45 #include "gromacs/math/utilities.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/smalloc.h"
49
50 real calc_similar_ind(gmx_bool bRho, int nind, atom_id *index, real mass[],
51                       rvec x[], rvec xp[])
52 {
53     int  i, j, d;
54     real m, tm, xs, xd, rs, rd;
55
56     tm = 0;
57     rs = 0;
58     rd = 0;
59     for (j = 0; j < nind; j++)
60     {
61         if (index)
62         {
63             i = index[j];
64         }
65         else
66         {
67             i = j;
68         }
69         m   = mass[i];
70         tm += m;
71         for (d = 0; d < DIM; d++)
72         {
73             xd  = x[i][d] - xp[i][d];
74             rd += m * sqr(xd);
75             if (bRho)
76             {
77                 xs  = x[i][d] + xp[i][d];
78                 rs += m * sqr(xs);
79             }
80         }
81     }
82     if (bRho)
83     {
84         return 2*sqrt(rd/rs);
85     }
86     else
87     {
88         return sqrt(rd/tm);
89     }
90 }
91
92 real rmsdev_ind(int nind, atom_id index[], real mass[], rvec x[], rvec xp[])
93 {
94     return calc_similar_ind(FALSE, nind, index, mass, x, xp);
95 }
96
97 real rmsdev(int natoms, real mass[], rvec x[], rvec xp[])
98 {
99     return calc_similar_ind(FALSE, natoms, NULL, mass, x, xp);
100 }
101
102 real rhodev_ind(int nind, atom_id index[], real mass[], rvec x[], rvec xp[])
103 {
104     return calc_similar_ind(TRUE, nind, index, mass, x, xp);
105 }
106
107 real rhodev(int natoms, real mass[], rvec x[], rvec xp[])
108 {
109     return calc_similar_ind(TRUE, natoms, NULL, mass, x, xp);
110 }
111
112 void calc_fit_R(int ndim, int natoms, real *w_rls, rvec *xp, rvec *x, matrix R)
113 {
114     int      c, r, n, j, m, i, irot, s;
115     double **omega, **om;
116     double   d[2*DIM], xnr, xpc;
117     matrix   vh, vk, u;
118     real     mn;
119     int      index;
120     real     max_d;
121
122     if (ndim != 3 && ndim != 2)
123     {
124         gmx_fatal(FARGS, "calc_fit_R called with ndim=%d instead of 3 or 2", ndim);
125     }
126
127     snew(omega, 2*ndim);
128     snew(om, 2*ndim);
129     for (i = 0; i < 2*ndim; i++)
130     {
131         snew(omega[i], 2*ndim);
132         snew(om[i], 2*ndim);
133     }
134
135     for (i = 0; i < 2*ndim; i++)
136     {
137         d[i] = 0;
138         for (j = 0; j < 2*ndim; j++)
139         {
140             omega[i][j] = 0;
141             om[i][j]    = 0;
142         }
143     }
144
145     /*calculate the matrix U*/
146     clear_mat(u);
147     for (n = 0; (n < natoms); n++)
148     {
149         if ((mn = w_rls[n]) != 0.0)
150         {
151             for (c = 0; (c < ndim); c++)
152             {
153                 xpc = xp[n][c];
154                 for (r = 0; (r < ndim); r++)
155                 {
156                     xnr      = x[n][r];
157                     u[c][r] += mn*xnr*xpc;
158                 }
159             }
160         }
161     }
162
163     /*construct omega*/
164     /*omega is symmetric -> omega==omega' */
165     for (r = 0; r < 2*ndim; r++)
166     {
167         for (c = 0; c <= r; c++)
168         {
169             if (r >= ndim && c < ndim)
170             {
171                 omega[r][c] = u[r-ndim][c];
172                 omega[c][r] = u[r-ndim][c];
173             }
174             else
175             {
176                 omega[r][c] = 0;
177                 omega[c][r] = 0;
178             }
179         }
180     }
181
182     /*determine h and k*/
183     jacobi(omega, 2*ndim, d, om, &irot);
184     /*real   **omega = input matrix a[0..n-1][0..n-1] must be symmetric
185      * int     natoms = number of rows and columns
186      * real      NULL = d[0]..d[n-1] are the eigenvalues of a[][]
187      * real       **v = v[0..n-1][0..n-1] contains the vectors in columns
188      * int      *irot = number of jacobi rotations
189      */
190
191     if (debug && irot == 0)
192     {
193         fprintf(debug, "IROT=0\n");
194     }
195
196     index = 0; /* For the compiler only */
197
198     /* Copy only the first ndim-1 eigenvectors */
199     for (j = 0; j < ndim-1; j++)
200     {
201         max_d = -1000;
202         for (i = 0; i < 2*ndim; i++)
203         {
204             if (d[i] > max_d)
205             {
206                 max_d = d[i];
207                 index = i;
208             }
209         }
210         d[index] = -10000;
211         for (i = 0; i < ndim; i++)
212         {
213             vh[j][i] = M_SQRT2*om[i][index];
214             vk[j][i] = M_SQRT2*om[i+ndim][index];
215         }
216     }
217     if (ndim == 3)
218     {
219         /* Calculate the last eigenvector as the outer-product of the first two.
220          * This insures that the conformation is not mirrored and
221          * prevents problems with completely flat reference structures.
222          */
223         cprod(vh[0], vh[1], vh[2]);
224         cprod(vk[0], vk[1], vk[2]);
225     }
226     else if (ndim == 2)
227     {
228         /* Calculate the last eigenvector from the first one */
229         vh[1][XX] = -vh[0][YY];
230         vh[1][YY] =  vh[0][XX];
231         vk[1][XX] = -vk[0][YY];
232         vk[1][YY] =  vk[0][XX];
233     }
234
235     /* determine R */
236     clear_mat(R);
237     for (r = 0; r < ndim; r++)
238     {
239         for (c = 0; c < ndim; c++)
240         {
241             for (s = 0; s < ndim; s++)
242             {
243                 R[r][c] += vk[s][r]*vh[s][c];
244             }
245         }
246     }
247     for (r = ndim; r < DIM; r++)
248     {
249         R[r][r] = 1;
250     }
251
252     for (i = 0; i < 2*ndim; i++)
253     {
254         sfree(omega[i]);
255         sfree(om[i]);
256     }
257     sfree(omega);
258     sfree(om);
259 }
260
261 void do_fit_ndim(int ndim, int natoms, real *w_rls, rvec *xp, rvec *x)
262 {
263     int    i, j, m, r, c;
264     matrix R;
265     rvec   x_old;
266
267     /* Calculate the rotation matrix R */
268     calc_fit_R(ndim, natoms, w_rls, xp, x, R);
269
270     /*rotate X*/
271     for (j = 0; j < natoms; j++)
272     {
273         for (m = 0; m < DIM; m++)
274         {
275             x_old[m] = x[j][m];
276         }
277         for (r = 0; r < DIM; r++)
278         {
279             x[j][r] = 0;
280             for (c = 0; c < DIM; c++)
281             {
282                 x[j][r] += R[r][c]*x_old[c];
283             }
284         }
285     }
286 }
287
288 void do_fit(int natoms, real *w_rls, rvec *xp, rvec *x)
289 {
290     do_fit_ndim(3, natoms, w_rls, xp, x);
291 }
292
293 void reset_x_ndim(int ndim, int ncm, const atom_id *ind_cm,
294                   int nreset, const atom_id *ind_reset,
295                   rvec x[], const real mass[])
296 {
297     int  i, m, ai;
298     rvec xcm;
299     real tm, mm;
300
301     if (ndim > DIM)
302     {
303         gmx_incons("More than 3 dimensions not supported.");
304     }
305     tm = 0.0;
306     clear_rvec(xcm);
307     if (ind_cm != NULL)
308     {
309         for (i = 0; i < ncm; i++)
310         {
311             ai = ind_cm[i];
312             mm = mass[ai];
313             for (m = 0; m < ndim; m++)
314             {
315                 xcm[m] += mm*x[ai][m];
316             }
317             tm += mm;
318         }
319     }
320     else
321     {
322         for (i = 0; i < ncm; i++)
323         {
324             mm = mass[i];
325             for (m = 0; m < ndim; m++)
326             {
327                 xcm[m] += mm*x[i][m];
328             }
329             tm += mm;
330         }
331     }
332     for (m = 0; m < ndim; m++)
333     {
334         xcm[m] /= tm;
335     }
336
337     if (ind_reset != NULL)
338     {
339         for (i = 0; i < nreset; i++)
340         {
341             rvec_dec(x[ind_reset[i]], xcm);
342         }
343     }
344     else
345     {
346         for (i = 0; i < nreset; i++)
347         {
348             rvec_dec(x[i], xcm);
349         }
350     }
351 }
352
353 void reset_x(int ncm, const atom_id *ind_cm,
354              int nreset, const atom_id *ind_reset,
355              rvec x[], const real mass[])
356 {
357     reset_x_ndim(3, ncm, ind_cm, nreset, ind_reset, x, mass);
358 }