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, 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/legacyheaders/macros.h"
49 #include "gromacs/legacyheaders/txtdump.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/random/random.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
61 int gmx_nmens(int argc, char *argv[])
63 const char *desc[] = {
64 "[THISMODULE] generates an ensemble around an average structure",
65 "in a subspace that is defined by a set of normal modes (eigenvectors).",
66 "The eigenvectors are assumed to be mass-weighted.",
67 "The position along each eigenvector is randomly taken from a Gaussian",
68 "distribution with variance kT/eigenvalue.[PAR]",
69 "By default the starting eigenvector is set to 7, since the first six",
70 "normal modes are the translational and rotational degrees of freedom."
72 static int nstruct = 100, first = 7, last = -1, seed = -1;
73 static real temp = 300.0;
75 { "-temp", FALSE, etREAL, {&temp},
76 "Temperature in Kelvin" },
77 { "-seed", FALSE, etINT, {&seed},
78 "Random seed, -1 generates a seed from time and pid" },
79 { "-num", FALSE, etINT, {&nstruct},
80 "Number of structures to generate" },
81 { "-first", FALSE, etINT, {&first},
82 "First eigenvector to use (-1 is select)" },
83 { "-last", FALSE, etINT, {&last},
84 "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 = NULL;
98 real *eigval, totmass, *invsqrtm, t, disp;
100 char *grpname, title[STRLEN];
101 const char *indexfile;
103 int nout, *iout, noutvec, *outvec;
105 real rfac, invfr, rhalf, jr;
110 const unsigned long im = 0xffff;
111 const unsigned long ia = 1093;
112 const unsigned long ic = 18257;
116 { efTRN, "-v", "eigenvec", ffREAD },
117 { efXVG, "-e", "eigenval", ffREAD },
118 { efTPS, NULL, NULL, ffREAD },
119 { efNDX, NULL, NULL, ffOPTRD },
120 { efTRO, "-o", "ensemble", ffWRITE }
122 #define NFILE asize(fnm)
124 if (!parse_common_args(&argc, argv, 0,
125 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
130 indexfile = ftp2fn_null(efNDX, NFILE, fnm);
132 read_eigenvectors(opt2fn("-v", NFILE, fnm), &natoms, &bFit,
133 &xref, &bDMR, &xav, &bDMA, &nvec, &eignr, &eigvec, &eigval);
135 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, bDMA);
138 printf("\nSelect an index group of %d elements that corresponds to the eigenvectors\n", natoms);
139 get_index(atoms, indexfile, 1, &i, &index, &grpname);
142 gmx_fatal(FARGS, "you selected a group with %d elements instead of %d",
147 snew(invsqrtm, natoms);
150 for (i = 0; (i < natoms); i++)
152 invsqrtm[i] = gmx_invsqrt(atoms->atom[index[i]].m);
157 for (i = 0; (i < natoms); i++)
169 /* make an index from first to last */
172 for (i = 0; i < nout; i++)
179 printf("Select eigenvectors for output, end your selection with 0\n");
185 srenew(iout, nout+1);
186 if (1 != scanf("%d", &iout[nout]))
188 gmx_fatal(FARGS, "Error reading user input");
192 while (iout[nout] >= 0);
196 /* make an index of the eigenvectors which are present */
199 for (i = 0; i < nout; i++)
202 while ((j < nvec) && (eignr[j] != iout[i]))
206 if ((j < nvec) && (eignr[j] == iout[i]))
209 iout[noutvec] = iout[i];
214 fprintf(stderr, "%d eigenvectors selected for output\n", noutvec);
218 seed = (int)gmx_rng_make_seed();
219 rng = gmx_rng_init(seed);
223 rng = gmx_rng_init(seed);
225 fprintf(stderr, "Using seed %d and a temperature of %g K\n", seed, temp);
228 snew(xout2, atoms->nr);
229 out = open_trx(ftp2fn(efTRO, NFILE, fnm), "w");
230 jran = (unsigned long)((real)im*gmx_rng_uniform_real(rng));
231 gmx_rng_destroy(rng);
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 = sqrt(3.0 * BOLTZ*temp/eigval[iout[j]]);
246 rfac = rfac/(real)im;
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, NULL, NULL);
272 fprintf(stderr, "\rGenerated %d structures", s+1);
274 fprintf(stderr, "\n");