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,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
11 *
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
16 *
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
21 *
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 *
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
34 *
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
37 */
38 #include "gmxpre.h"
39
40 #include "do_fit.h"
41
42 #include <cmath>
43 #include <cstdio>
44
45 #include "gromacs/linearalgebra/nrjac.h"
46 #include "gromacs/math/functions.h"
47 #include "gromacs/math/utilities.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/smalloc.h"
51
calc_similar_ind(gmx_bool bRho,int nind,const int * index,const real mass[],rvec x[],rvec xp[])52 real calc_similar_ind(gmx_bool bRho, int nind, const int* index, const real mass[], rvec x[], rvec xp[])
53 {
54 int i, j, d;
55 real m, tm, xs, xd, rs, rd;
56
57 tm = 0;
58 rs = 0;
59 rd = 0;
60 for (j = 0; j < nind; j++)
61 {
62 if (index)
63 {
64 i = index[j];
65 }
66 else
67 {
68 i = j;
69 }
70 m = mass[i];
71 tm += m;
72 for (d = 0; d < DIM; d++)
73 {
74 xd = x[i][d] - xp[i][d];
75 rd += m * gmx::square(xd);
76 if (bRho)
77 {
78 xs = x[i][d] + xp[i][d];
79 rs += m * gmx::square(xs);
80 }
81 }
82 }
83 if (bRho)
84 {
85 return 2 * std::sqrt(rd / rs);
86 }
87 else
88 {
89 return std::sqrt(rd / tm);
90 }
91 }
92
rmsdev_ind(int nind,int index[],real mass[],rvec x[],rvec xp[])93 real rmsdev_ind(int nind, int index[], real mass[], rvec x[], rvec xp[])
94 {
95 return calc_similar_ind(FALSE, nind, index, mass, x, xp);
96 }
97
rmsdev(int natoms,real mass[],rvec x[],rvec xp[])98 real rmsdev(int natoms, real mass[], rvec x[], rvec xp[])
99 {
100 return calc_similar_ind(FALSE, natoms, nullptr, mass, x, xp);
101 }
102
rhodev_ind(int nind,int index[],real mass[],rvec x[],rvec xp[])103 real rhodev_ind(int nind, int index[], real mass[], rvec x[], rvec xp[])
104 {
105 return calc_similar_ind(TRUE, nind, index, mass, x, xp);
106 }
107
rhodev(int natoms,real mass[],rvec x[],rvec xp[])108 real rhodev(int natoms, real mass[], rvec x[], rvec xp[])
109 {
110 return calc_similar_ind(TRUE, natoms, nullptr, mass, x, xp);
111 }
112
calc_fit_R(int ndim,int natoms,const real * w_rls,const rvec * xp,rvec * x,matrix R)113 void calc_fit_R(int ndim, int natoms, const real* w_rls, const rvec* xp, rvec* x, matrix R)
114 {
115 int c, r, n, j, i, irot, s;
116 double **omega, **om;
117 double d[2 * DIM], xnr, xpc;
118 matrix vh, vk, u;
119 real mn;
120 int index;
121 real max_d;
122
123 if (ndim != 3 && ndim != 2)
124 {
125 gmx_fatal(FARGS, "calc_fit_R called with ndim=%d instead of 3 or 2", ndim);
126 }
127
128 snew(omega, 2 * ndim);
129 snew(om, 2 * ndim);
130 for (i = 0; i < 2 * ndim; i++)
131 {
132 snew(omega[i], 2 * ndim);
133 snew(om[i], 2 * ndim);
134 }
135
136 for (i = 0; i < 2 * ndim; i++)
137 {
138 d[i] = 0;
139 for (j = 0; j < 2 * ndim; j++)
140 {
141 omega[i][j] = 0;
142 om[i][j] = 0;
143 }
144 }
145
146 /*calculate the matrix U*/
147 clear_mat(u);
148 for (n = 0; (n < natoms); n++)
149 {
150 if ((mn = w_rls[n]) != 0.0)
151 {
152 for (c = 0; (c < ndim); c++)
153 {
154 xpc = xp[n][c];
155 for (r = 0; (r < ndim); r++)
156 {
157 xnr = x[n][r];
158 u[c][r] += mn * xnr * xpc;
159 }
160 }
161 }
162 }
163
164 /*construct omega*/
165 /*omega is symmetric -> omega==omega' */
166 for (r = 0; r < 2 * ndim; r++)
167 {
168 for (c = 0; c <= r; c++)
169 {
170 if (r >= ndim && c < ndim)
171 {
172 omega[r][c] = u[r - ndim][c];
173 omega[c][r] = u[r - ndim][c];
174 }
175 else
176 {
177 omega[r][c] = 0;
178 omega[c][r] = 0;
179 }
180 }
181 }
182
183 /*determine h and k*/
184 jacobi(omega, 2 * ndim, d, om, &irot);
185 /*real **omega = input matrix a[0..n-1][0..n-1] must be symmetric
186 * int natoms = number of rows and columns
187 * real NULL = d[0]..d[n-1] are the eigenvalues of a[][]
188 * real **v = v[0..n-1][0..n-1] contains the vectors in columns
189 * int *irot = number of jacobi rotations
190 */
191
192 if (debug && irot == 0)
193 {
194 fprintf(debug, "IROT=0\n");
195 }
196
197 index = 0; /* For the compiler only */
198
199 /* Copy only the first ndim-1 eigenvectors */
200 for (j = 0; j < ndim - 1; j++)
201 {
202 max_d = -1000;
203 for (i = 0; i < 2 * ndim; i++)
204 {
205 if (d[i] > max_d)
206 {
207 max_d = d[i];
208 index = i;
209 }
210 }
211 d[index] = -10000;
212 for (i = 0; i < ndim; i++)
213 {
214 vh[j][i] = M_SQRT2 * om[i][index];
215 vk[j][i] = M_SQRT2 * om[i + ndim][index];
216 }
217 }
218 if (ndim == 3)
219 {
220 /* Calculate the last eigenvector as the outer-product of the first two.
221 * This insures that the conformation is not mirrored and
222 * prevents problems with completely flat reference structures.
223 */
224 cprod(vh[0], vh[1], vh[2]);
225 cprod(vk[0], vk[1], vk[2]);
226 }
227 else if (ndim == 2)
228 {
229 /* Calculate the last eigenvector from the first one */
230 vh[1][XX] = -vh[0][YY];
231 vh[1][YY] = vh[0][XX];
232 vk[1][XX] = -vk[0][YY];
233 vk[1][YY] = vk[0][XX];
234 }
235
236 /* determine R */
237 clear_mat(R);
238 for (r = 0; r < ndim; r++)
239 {
240 for (c = 0; c < ndim; c++)
241 {
242 for (s = 0; s < ndim; s++)
243 {
244 R[r][c] += vk[s][r] * vh[s][c];
245 }
246 }
247 }
248 for (r = ndim; r < DIM; r++)
249 {
250 R[r][r] = 1;
251 }
252
253 for (i = 0; i < 2 * ndim; i++)
254 {
255 sfree(omega[i]);
256 sfree(om[i]);
257 }
258 sfree(omega);
259 sfree(om);
260 }
261
do_fit_ndim(int ndim,int natoms,real * w_rls,const rvec * xp,rvec * x)262 void do_fit_ndim(int ndim, int natoms, real* w_rls, const rvec* xp, rvec* x)
263 {
264 int j, m, r, c;
265 matrix R;
266 rvec x_old;
267
268 /* Calculate the rotation matrix R */
269 calc_fit_R(ndim, natoms, w_rls, xp, x, R);
270
271 /*rotate X*/
272 for (j = 0; j < natoms; j++)
273 {
274 for (m = 0; m < DIM; m++)
275 {
276 x_old[m] = x[j][m];
277 }
278 for (r = 0; r < DIM; r++)
279 {
280 x[j][r] = 0;
281 for (c = 0; c < DIM; c++)
282 {
283 x[j][r] += R[r][c] * x_old[c];
284 }
285 }
286 }
287 }
288
do_fit(int natoms,real * w_rls,const rvec * xp,rvec * x)289 void do_fit(int natoms, real* w_rls, const rvec* xp, rvec* x)
290 {
291 do_fit_ndim(3, natoms, w_rls, xp, x);
292 }
293
reset_x_ndim(int ndim,int ncm,const int * ind_cm,int nreset,const int * ind_reset,rvec x[],const real mass[])294 void reset_x_ndim(int ndim, int ncm, const int* ind_cm, int nreset, const int* ind_reset, rvec x[], const real mass[])
295 {
296 int i, m, ai;
297 rvec xcm;
298 real tm, mm;
299
300 if (ndim > DIM)
301 {
302 gmx_incons("More than 3 dimensions not supported.");
303 }
304 tm = 0.0;
305 clear_rvec(xcm);
306 if (ind_cm != nullptr)
307 {
308 for (i = 0; i < ncm; i++)
309 {
310 ai = ind_cm[i];
311 mm = mass[ai];
312 for (m = 0; m < ndim; m++)
313 {
314 xcm[m] += mm * x[ai][m];
315 }
316 tm += mm;
317 }
318 }
319 else
320 {
321 for (i = 0; i < ncm; i++)
322 {
323 mm = mass[i];
324 for (m = 0; m < ndim; m++)
325 {
326 xcm[m] += mm * x[i][m];
327 }
328 tm += mm;
329 }
330 }
331 for (m = 0; m < ndim; m++)
332 {
333 xcm[m] /= tm;
334 }
335
336 if (ind_reset != nullptr)
337 {
338 for (i = 0; i < nreset; i++)
339 {
340 rvec_dec(x[ind_reset[i]], xcm);
341 }
342 }
343 else
344 {
345 for (i = 0; i < nreset; i++)
346 {
347 rvec_dec(x[i], xcm);
348 }
349 }
350 }
351
reset_x(int ncm,const int * ind_cm,int nreset,const int * ind_reset,rvec x[],const real mass[])352 void reset_x(int ncm, const int* ind_cm, int nreset, const int* ind_reset, rvec x[], const real mass[])
353 {
354 reset_x_ndim(3, ncm, ind_cm, nreset, ind_reset, x, mass);
355 }
356