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) 2013,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.
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.
44 #include "gromacs/commandline/pargs.h"
49 #include "gmx_fatal.h"
51 #include "gromacs/fileio/futil.h"
53 #include "gromacs/fileio/pdbio.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/fileio/trxio.h"
62 int gmx_nmtraj(int argc, char *argv[])
66 "[THISMODULE] generates an virtual trajectory from an eigenvector, ",
67 "corresponding to a harmonic Cartesian oscillation around the average ",
68 "structure. The eigenvectors should normally be mass-weighted, but you can ",
69 "use non-weighted eigenvectors to generate orthogonal motions. ",
70 "The output frames are written as a trajectory file covering an entire period, and ",
71 "the first frame is the average structure. If you write the trajectory in (or convert to) ",
72 "PDB format you can view it directly in PyMol and also render a photorealistic movie. ",
73 "Motion amplitudes are calculated from the eigenvalues and a preset temperature, ",
74 "assuming equipartition of the energy over all modes. To make the motion clearly visible ",
75 "in PyMol you might want to amplify it by setting an unrealistically high temperature. ",
76 "However, be aware that both the linear Cartesian displacements and mass weighting will ",
77 "lead to serious structure deformation for high amplitudes - this is is simply a limitation ",
78 "of the Cartesian normal mode model. By default the selected eigenvector is set to 7, since ",
79 " the first six normal modes are the translational and rotational degrees of freedom."
82 static real refamplitude = 0.25;
83 static int nframes = 30;
84 static real temp = 300.0;
85 static const char *eignrvec = "7";
86 static const char *phasevec = "0.0";
90 { "-eignr", FALSE, etSTR, {&eignrvec}, "String of eigenvectors to use (first is 1)" },
91 { "-phases", FALSE, etSTR, {&phasevec}, "String of phases (default is 0.0)" },
92 { "-temp", FALSE, etREAL, {&temp}, "Temperature (K)" },
93 { "-amplitude", FALSE, etREAL, {&refamplitude}, "Amplitude for modes with eigenvalue<=0" },
94 { "-nframes", FALSE, etINT, {&nframes}, "Number of frames to generate" }
103 rvec *xtop, *xref, *xav, *xout;
104 int nvec, *eignr = NULL;
106 rvec **eigvec = NULL;
109 int i, j, k, kmode, d, s, v;
110 gmx_bool bDMR, bDMA, bFit;
123 real omega, Ekin, sum, m, vel;
136 { efTPS, NULL, NULL, ffREAD },
137 { efTRN, "-v", "eigenvec", ffREAD },
138 { efTRO, "-o", "nmtraj", ffWRITE }
141 #define NFILE asize(fnm)
143 if (!parse_common_args(&argc, argv, PCA_BE_NICE,
144 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
149 read_eigenvectors(opt2fn("-v", NFILE, fnm), &natoms, &bFit,
150 &xref, &bDMR, &xav, &bDMA, &nvec, &eignr, &eigvec, &eigval);
152 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, bDMA);
154 /* Find vectors and phases */
156 /* first find number of args in string */
161 dum = strtod(p, &pe);
166 snew(imodes, nmodes);
168 for (i = 0; i < nmodes; i++)
170 /* C indices start on 0 */
171 imodes[i] = strtol(p, &pe, 10)-1;
175 /* Now read phases */
180 dum = strtod(p, &pe);
184 if (nphases > nmodes)
186 gmx_fatal(FARGS, "More phases than eigenvector indices specified.\n");
189 snew(phases, nmodes);
192 for (i = 0; i < nphases; i++)
194 phases[i] = strtod(p, &pe);
198 if (nmodes > nphases)
200 printf("Warning: Setting phase of last %d modes to zero...\n", nmodes-nphases);
203 for (i = nphases; i < nmodes; i++)
210 if (atoms->nr != natoms)
212 gmx_fatal(FARGS, "Different number of atoms in topology and eigenvectors.\n");
216 for (i = 0; i < natoms; i++)
221 /* Find the eigenvalue/vector to match our select one */
222 snew(out_eigidx, nmodes);
223 for (i = 0; i < nmodes; i++)
228 for (i = 0; i < nvec; i++)
230 for (j = 0; j < nmodes; j++)
232 if (imodes[j] == eignr[i])
238 for (i = 0; i < nmodes; i++)
240 if (out_eigidx[i] == -1)
242 gmx_fatal(FARGS, "Could not find mode %d in eigenvector file.\n", imodes[i]);
247 snew(invsqrtm, natoms);
251 for (i = 0; (i < natoms); i++)
253 invsqrtm[i] = gmx_invsqrt(atoms->atom[i].m);
258 for (i = 0; (i < natoms); i++)
265 snew(amplitude, nmodes);
267 printf("mode phases: %g %g\n", phases[0], phases[1]);
269 for (i = 0; i < nmodes; i++)
271 kmode = out_eigidx[i];
272 this_eigvec = eigvec[kmode];
274 if ( (kmode >= 6) && (eigval[kmode] > 0))
276 /* Derive amplitude from temperature and eigenvalue if we can */
278 /* Convert eigenvalue to angular frequency, in units s^(-1) */
279 omega = sqrt(eigval[kmode]*1.0E21/(AVOGADRO*AMU));
280 /* Harmonic motion will be x=x0 + A*sin(omega*t)*eigenvec.
281 * The velocity is thus:
283 * v = A*omega*cos(omega*t)*eigenvec.
285 * And the average kinetic energy the integral of mass*v*v/2 over a
288 * (1/4)*mass*A*omega*eigenvec
290 * For t =2*pi*n, all energy will be kinetic, and v=A*omega*eigenvec.
291 * The kinetic energy will be sum(0.5*mass*v*v) if we temporarily set A to 1,
292 * and the average over a period half of this.
296 for (k = 0; k < natoms; k++)
298 m = atoms->atom[k].m;
299 for (d = 0; d < DIM; d++)
301 vel = omega*this_eigvec[k][d];
302 Ekin += 0.5*0.5*m*vel*vel;
306 /* Convert Ekin from amu*(nm/s)^2 to J, i.e., kg*(m/s)^2
307 * This will also be proportional to A^2
311 /* Set the amplitude so the energy is kT/2 */
312 amplitude[i] = sqrt(0.5*BOLTZMANN*temp/Ekin);
316 amplitude[i] = refamplitude;
320 out = open_trx(ftp2fn(efTRO, NFILE, fnm), "w");
322 /* Write a sine oscillation around the average structure,
323 * modulated by the eigenvector with selected amplitude.
326 for (i = 0; i < nframes; i++)
328 fraction = (real)i/(real)nframes;
329 for (j = 0; j < natoms; j++)
331 copy_rvec(xav[j], xout[j]);
334 for (k = 0; k < nmodes; k++)
336 kmode = out_eigidx[k];
337 this_eigvec = eigvec[kmode];
339 for (j = 0; j < natoms; j++)
341 for (d = 0; d < DIM; d++)
343 xout[j][d] += amplitude[k]*sin(2*M_PI*(fraction+phases[k]/360.0))*this_eigvec[j][d];
347 write_trx(out, natoms, dummy, atoms, i, (real)i/(real)nframes, box, xout, NULL, NULL);
350 fprintf(stderr, "\n");