3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
47 #include "gmx_fatal.h"
61 int gmx_nmens(int argc, char *argv[])
63 const char *desc[] = {
64 "[TT]g_nmens[tt] 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 parse_common_args(&argc, argv, PCA_BE_NICE,
125 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv);
127 indexfile = ftp2fn_null(efNDX, NFILE, fnm);
129 read_eigenvectors(opt2fn("-v", NFILE, fnm), &natoms, &bFit,
130 &xref, &bDMR, &xav, &bDMA, &nvec, &eignr, &eigvec, &eigval);
132 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, bDMA);
135 printf("\nSelect an index group of %d elements that corresponds to the eigenvectors\n", natoms);
136 get_index(atoms, indexfile, 1, &i, &index, &grpname);
139 gmx_fatal(FARGS, "you selected a group with %d elements instead of %d",
144 snew(invsqrtm, natoms);
147 for (i = 0; (i < natoms); i++)
149 invsqrtm[i] = gmx_invsqrt(atoms->atom[index[i]].m);
154 for (i = 0; (i < natoms); i++)
166 /* make an index from first to last */
169 for (i = 0; i < nout; i++)
176 printf("Select eigenvectors for output, end your selection with 0\n");
182 srenew(iout, nout+1);
183 if (1 != scanf("%d", &iout[nout]))
185 gmx_fatal(FARGS, "Error reading user input");
189 while (iout[nout] >= 0);
193 /* make an index of the eigenvectors which are present */
196 for (i = 0; i < nout; i++)
199 while ((j < nvec) && (eignr[j] != iout[i]))
203 if ((j < nvec) && (eignr[j] == iout[i]))
206 iout[noutvec] = iout[i];
211 fprintf(stderr, "%d eigenvectors selected for output\n", noutvec);
217 fprintf(stderr, "Using seed %d and a temperature of %g K\n", seed, temp);
220 snew(xout2, atoms->nr);
221 out = open_trx(ftp2fn(efTRO, NFILE, fnm), "w");
222 jran = (unsigned long)((real)im*rando(&seed));
223 for (s = 0; s < nstruct; s++)
225 for (i = 0; i < natoms; i++)
227 copy_rvec(xav[i], xout1[i]);
229 for (j = 0; j < noutvec; j++)
232 /* (r-0.5) n times: var_n = n * var_1 = n/12
233 n=4: var_n = 1/3, so multiply with 3 */
235 rfac = sqrt(3.0 * BOLTZ*temp/eigval[iout[j]]);
237 rfac = rfac/(real)im;
239 jran = (jran*ia+ic) & im;
241 jran = (jran*ia+ic) & im;
243 jran = (jran*ia+ic) & im;
245 jran = (jran*ia+ic) & im;
247 disp = rfac * jr - rhalf;
249 for (i = 0; i < natoms; i++)
251 for (d = 0; d < DIM; d++)
253 xout1[i][d] += disp*eigvec[v][i][d]*invsqrtm[i];
257 for (i = 0; i < natoms; i++)
259 copy_rvec(xout1[i], xout2[index[i]]);
262 write_trx(out, natoms, index, atoms, 0, t, box, xout2, NULL, NULL);
263 fprintf(stderr, "\rGenerated %d structures", s+1);
265 fprintf(stderr, "\n");