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,2015,2016,2017, 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.
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/fileio/pdbio.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/gmxana/eigio.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/math/functions.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/random/threefry.h"
52 #include "gromacs/random/uniformintdistribution.h"
53 #include "gromacs/topology/index.h"
54 #include "gromacs/topology/topology.h"
55 #include "gromacs/utility/arraysize.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/smalloc.h"
62 int gmx_nmens(int argc, char *argv[])
64 const char *desc[] = {
65 "[THISMODULE] generates an ensemble around an average structure",
66 "in a subspace that is defined by a set of normal modes (eigenvectors).",
67 "The eigenvectors are assumed to be mass-weighted.",
68 "The position along each eigenvector is randomly taken from a Gaussian",
69 "distribution with variance kT/eigenvalue.[PAR]",
70 "By default the starting eigenvector is set to 7, since the first six",
71 "normal modes are the translational and rotational degrees of freedom."
73 static int nstruct = 100, first = 7, last = -1, seed = 0;
74 static real temp = 300.0;
76 { "-temp", FALSE, etREAL, {&temp},
77 "Temperature in Kelvin" },
78 { "-seed", FALSE, etINT, {&seed},
79 "Random seed (0 means generate)" },
80 { "-num", FALSE, etINT, {&nstruct},
81 "Number of structures to generate" },
82 { "-first", FALSE, etINT, {&first},
83 "First eigenvector to use (-1 is select)" },
84 { "-last", FALSE, etINT, {&last},
85 "Last eigenvector to use (-1 is till the last)" }
93 rvec *xtop, *xref, *xav, *xout1, *xout2;
94 gmx_bool bDMR, bDMA, bFit;
95 int nvec, *eignr = nullptr;
96 rvec **eigvec = nullptr;
98 real *eigval, *invsqrtm, t, disp;
101 const char *indexfile;
103 int nout, *iout, noutvec, *outvec;
105 real rfac, rhalf, jr;
106 gmx_output_env_t *oenv;
108 const unsigned long im = 0xffff;
109 const unsigned long ia = 1093;
110 const unsigned long ic = 18257;
114 { efTRN, "-v", "eigenvec", ffREAD },
115 { efXVG, "-e", "eigenval", ffREAD },
116 { efTPS, nullptr, nullptr, ffREAD },
117 { efNDX, nullptr, nullptr, ffOPTRD },
118 { efTRO, "-o", "ensemble", ffWRITE }
120 #define NFILE asize(fnm)
122 if (!parse_common_args(&argc, argv, 0,
123 NFILE, fnm, NPA, pa, asize(desc), desc, 0, nullptr, &oenv))
128 indexfile = ftp2fn_null(efNDX, NFILE, fnm);
130 read_eigenvectors(opt2fn("-v", NFILE, fnm), &natoms, &bFit,
131 &xref, &bDMR, &xav, &bDMA, &nvec, &eignr, &eigvec, &eigval);
133 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &xtop, nullptr, box, bDMA);
136 printf("\nSelect an index group of %d elements that corresponds to the eigenvectors\n", natoms);
137 get_index(atoms, indexfile, 1, &i, &index, &grpname);
140 gmx_fatal(FARGS, "you selected a group with %d elements instead of %d",
145 snew(invsqrtm, natoms);
148 for (i = 0; (i < natoms); i++)
150 invsqrtm[i] = gmx::invsqrt(atoms->atom[index[i]].m);
155 for (i = 0; (i < natoms); i++)
167 /* make an index from first to last */
170 for (i = 0; i < nout; i++)
177 printf("Select eigenvectors for output, end your selection with 0\n");
183 srenew(iout, nout+1);
184 if (1 != scanf("%d", &iout[nout]))
186 gmx_fatal(FARGS, "Error reading user input");
190 while (iout[nout] >= 0);
194 /* make an index of the eigenvectors which are present */
197 for (i = 0; i < nout; i++)
200 while ((j < nvec) && (eignr[j] != iout[i]))
204 if ((j < nvec) && (eignr[j] == iout[i]))
207 iout[noutvec] = iout[i];
212 fprintf(stderr, "%d eigenvectors selected for output\n", noutvec);
217 // Make do with 32 bits for now to avoid changing user input to hex
218 seed = static_cast<int>(gmx::makeRandomSeed());
221 gmx::DefaultRandomEngine rng(seed);
223 fprintf(stderr, "Using random seed %d and a temperature of %g K.\n", seed, temp);
225 gmx::UniformIntDistribution<int> dist(0, im-1);
229 snew(xout2, atoms->nr);
230 out = open_trx(ftp2fn(efTRO, NFILE, fnm), "w");
232 for (s = 0; s < nstruct; s++)
234 for (i = 0; i < natoms; i++)
236 copy_rvec(xav[i], xout1[i]);
238 for (j = 0; j < noutvec; j++)
241 /* (r-0.5) n times: var_n = n * var_1 = n/12
242 n=4: var_n = 1/3, so multiply with 3 */
244 rfac = std::sqrt(3.0 * BOLTZ*temp/eigval[iout[j]]);
248 jran = (jran*ia+ic) & im;
250 jran = (jran*ia+ic) & im;
252 jran = (jran*ia+ic) & im;
254 jran = (jran*ia+ic) & im;
256 disp = rfac * jr - rhalf;
258 for (i = 0; i < natoms; i++)
260 for (d = 0; d < DIM; d++)
262 xout1[i][d] += disp*eigvec[v][i][d]*invsqrtm[i];
266 for (i = 0; i < natoms; i++)
268 copy_rvec(xout1[i], xout2[index[i]]);
271 write_trx(out, natoms, index, atoms, 0, t, box, xout2, nullptr, nullptr);
272 fprintf(stderr, "\rGenerated %d structures", s+1);
275 fprintf(stderr, "\n");