Merge branch release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_nmens.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <math.h>
40 #include <string.h>
41
42 #include "statutil.h"
43 #include "sysstuff.h"
44 #include "typedefs.h"
45 #include "smalloc.h"
46 #include "macros.h"
47 #include "gmx_fatal.h"
48 #include "vec.h"
49 #include "futil.h"
50 #include "statutil.h"
51 #include "index.h"
52 #include "pdbio.h"
53 #include "tpxio.h"
54 #include "txtdump.h"
55 #include "physics.h"
56 #include "random.h"
57 #include "eigio.h"
58 #include "gmx_ana.h"
59
60
61 int gmx_nmens(int argc, char *argv[])
62 {
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."
71     };
72     static int  nstruct = 100, first = 7, last = -1, seed = -1;
73     static real temp    = 300.0;
74     t_pargs     pa[]    = {
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)" }
85     };
86 #define NPA asize(pa)
87
88     t_trxstatus        *out;
89     int                 status, trjout;
90     t_topology          top;
91     int                 ePBC;
92     t_atoms            *atoms;
93     rvec               *xtop, *xref, *xav, *xout1, *xout2;
94     gmx_bool            bDMR, bDMA, bFit;
95     int                 nvec, *eignr = NULL;
96     rvec              **eigvec = NULL;
97     matrix              box;
98     real               *eigval, totmass, *invsqrtm, t, disp;
99     int                 natoms, neigval;
100     char               *grpname, title[STRLEN];
101     const char         *indexfile;
102     int                 i, j, d, s, v;
103     int                 nout, *iout, noutvec, *outvec;
104     atom_id            *index;
105     real                rfac, invfr, rhalf, jr;
106     int          *      eigvalnr;
107     output_env_t        oenv;
108
109     unsigned long       jran;
110     const unsigned long im = 0xffff;
111     const unsigned long ia = 1093;
112     const unsigned long ic = 18257;
113
114
115     t_filenm fnm[] = {
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 }
121     };
122 #define NFILE asize(fnm)
123
124     parse_common_args(&argc, argv, PCA_BE_NICE,
125                       NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv);
126
127     indexfile = ftp2fn_null(efNDX, NFILE, fnm);
128
129     read_eigenvectors(opt2fn("-v", NFILE, fnm), &natoms, &bFit,
130                       &xref, &bDMR, &xav, &bDMA, &nvec, &eignr, &eigvec, &eigval);
131
132     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, bDMA);
133     atoms = &top.atoms;
134
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);
137     if (i != natoms)
138     {
139         gmx_fatal(FARGS, "you selected a group with %d elements instead of %d",
140                   i, natoms);
141     }
142     printf("\n");
143
144     snew(invsqrtm, natoms);
145     if (bDMA)
146     {
147         for (i = 0; (i < natoms); i++)
148         {
149             invsqrtm[i] = gmx_invsqrt(atoms->atom[index[i]].m);
150         }
151     }
152     else
153     {
154         for (i = 0; (i < natoms); i++)
155         {
156             invsqrtm[i] = 1.0;
157         }
158     }
159
160     if (last == -1)
161     {
162         last = natoms*DIM;
163     }
164     if (first > -1)
165     {
166         /* make an index from first to last */
167         nout = last-first+1;
168         snew(iout, nout);
169         for (i = 0; i < nout; i++)
170         {
171             iout[i] = first-1+i;
172         }
173     }
174     else
175     {
176         printf("Select eigenvectors for output, end your selection with 0\n");
177         nout = -1;
178         iout = NULL;
179         do
180         {
181             nout++;
182             srenew(iout, nout+1);
183             if (1 != scanf("%d", &iout[nout]))
184             {
185                 gmx_fatal(FARGS, "Error reading user input");
186             }
187             iout[nout]--;
188         }
189         while (iout[nout] >= 0);
190         printf("\n");
191     }
192
193     /* make an index of the eigenvectors which are present */
194     snew(outvec, nout);
195     noutvec = 0;
196     for (i = 0; i < nout; i++)
197     {
198         j = 0;
199         while ((j < nvec) && (eignr[j] != iout[i]))
200         {
201             j++;
202         }
203         if ((j < nvec) && (eignr[j] == iout[i]))
204         {
205             outvec[noutvec] = j;
206             iout[noutvec]   = iout[i];
207             noutvec++;
208         }
209     }
210
211     fprintf(stderr, "%d eigenvectors selected for output\n", noutvec);
212
213     if (seed == -1)
214     {
215         seed = make_seed();
216     }
217     fprintf(stderr, "Using seed %d and a temperature of %g K\n", seed, temp);
218
219     snew(xout1, natoms);
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++)
224     {
225         for (i = 0; i < natoms; i++)
226         {
227             copy_rvec(xav[i], xout1[i]);
228         }
229         for (j = 0; j < noutvec; j++)
230         {
231             v = outvec[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 */
234
235             rfac  = sqrt(3.0 * BOLTZ*temp/eigval[iout[j]]);
236             rhalf = 2.0*rfac;
237             rfac  = rfac/(real)im;
238
239             jran = (jran*ia+ic) & im;
240             jr   = (real)jran;
241             jran = (jran*ia+ic) & im;
242             jr  += (real)jran;
243             jran = (jran*ia+ic) & im;
244             jr  += (real)jran;
245             jran = (jran*ia+ic) & im;
246             jr  += (real)jran;
247             disp = rfac * jr - rhalf;
248
249             for (i = 0; i < natoms; i++)
250             {
251                 for (d = 0; d < DIM; d++)
252                 {
253                     xout1[i][d] += disp*eigvec[v][i][d]*invsqrtm[i];
254                 }
255             }
256         }
257         for (i = 0; i < natoms; i++)
258         {
259             copy_rvec(xout1[i], xout2[index[i]]);
260         }
261         t = s+1;
262         write_trx(out, natoms, index, atoms, 0, t, box, xout2, NULL, NULL);
263         fprintf(stderr, "\rGenerated %d structures", s+1);
264     }
265     fprintf(stderr, "\n");
266     close_trx(out);
267
268     return 0;
269 }