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 "copyrite.h"
50 #include "futil.h"
51 #include "statutil.h"
52 #include "index.h"
53 #include "pdbio.h"
54 #include "tpxio.h"
55 #include "txtdump.h"
56 #include "physics.h"
57 #include "random.h"
58 #include "eigio.h"
59 #include "gmx_ana.h"
60
61
62 int gmx_nmens(int argc, char *argv[])
63 {
64     const char *desc[] = {
65         "[TT]g_nmens[tt] 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."
72     };
73     static int  nstruct = 100, first = 7, last = -1, seed = -1;
74     static real temp    = 300.0;
75     t_pargs     pa[]    = {
76         { "-temp",  FALSE, etREAL, {&temp},
77           "Temperature in Kelvin" },
78         { "-seed", FALSE, etINT, {&seed},
79           "Random seed, -1 generates a seed from time and pid" },
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)" }
86     };
87 #define NPA asize(pa)
88
89     t_trxstatus        *out;
90     int                 status, trjout;
91     t_topology          top;
92     int                 ePBC;
93     t_atoms            *atoms;
94     rvec               *xtop, *xref, *xav, *xout1, *xout2;
95     gmx_bool            bDMR, bDMA, bFit;
96     int                 nvec, *eignr = NULL;
97     rvec              **eigvec = NULL;
98     matrix              box;
99     real               *eigval, totmass, *invsqrtm, t, disp;
100     int                 natoms, neigval;
101     char               *grpname, title[STRLEN];
102     const char         *indexfile;
103     int                 i, j, d, s, v;
104     int                 nout, *iout, noutvec, *outvec;
105     atom_id            *index;
106     real                rfac, invfr, rhalf, jr;
107     int          *      eigvalnr;
108     output_env_t        oenv;
109
110     unsigned long       jran;
111     const unsigned long im = 0xffff;
112     const unsigned long ia = 1093;
113     const unsigned long ic = 18257;
114
115
116     t_filenm fnm[] = {
117         { efTRN, "-v",    "eigenvec",    ffREAD  },
118         { efXVG, "-e",    "eigenval",    ffREAD  },
119         { efTPS, NULL,    NULL,          ffREAD },
120         { efNDX, NULL,    NULL,          ffOPTRD },
121         { efTRO, "-o",    "ensemble",    ffWRITE }
122     };
123 #define NFILE asize(fnm)
124
125     parse_common_args(&argc, argv, PCA_BE_NICE,
126                       NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv);
127
128     indexfile = ftp2fn_null(efNDX, NFILE, fnm);
129
130     read_eigenvectors(opt2fn("-v", NFILE, fnm), &natoms, &bFit,
131                       &xref, &bDMR, &xav, &bDMA, &nvec, &eignr, &eigvec, &eigval);
132
133     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, bDMA);
134     atoms = &top.atoms;
135
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);
138     if (i != natoms)
139     {
140         gmx_fatal(FARGS, "you selected a group with %d elements instead of %d",
141                   i, natoms);
142     }
143     printf("\n");
144
145     snew(invsqrtm, natoms);
146     if (bDMA)
147     {
148         for (i = 0; (i < natoms); i++)
149         {
150             invsqrtm[i] = gmx_invsqrt(atoms->atom[index[i]].m);
151         }
152     }
153     else
154     {
155         for (i = 0; (i < natoms); i++)
156         {
157             invsqrtm[i] = 1.0;
158         }
159     }
160
161     if (last == -1)
162     {
163         last = natoms*DIM;
164     }
165     if (first > -1)
166     {
167         /* make an index from first to last */
168         nout = last-first+1;
169         snew(iout, nout);
170         for (i = 0; i < nout; i++)
171         {
172             iout[i] = first-1+i;
173         }
174     }
175     else
176     {
177         printf("Select eigenvectors for output, end your selection with 0\n");
178         nout = -1;
179         iout = NULL;
180         do
181         {
182             nout++;
183             srenew(iout, nout+1);
184             if (1 != scanf("%d", &iout[nout]))
185             {
186                 gmx_fatal(FARGS, "Error reading user input");
187             }
188             iout[nout]--;
189         }
190         while (iout[nout] >= 0);
191         printf("\n");
192     }
193
194     /* make an index of the eigenvectors which are present */
195     snew(outvec, nout);
196     noutvec = 0;
197     for (i = 0; i < nout; i++)
198     {
199         j = 0;
200         while ((j < nvec) && (eignr[j] != iout[i]))
201         {
202             j++;
203         }
204         if ((j < nvec) && (eignr[j] == iout[i]))
205         {
206             outvec[noutvec] = j;
207             iout[noutvec]   = iout[i];
208             noutvec++;
209         }
210     }
211
212     fprintf(stderr, "%d eigenvectors selected for output\n", noutvec);
213
214     if (seed == -1)
215     {
216         seed = make_seed();
217     }
218     fprintf(stderr, "Using seed %d and a temperature of %g K\n", seed, temp);
219
220     snew(xout1, natoms);
221     snew(xout2, atoms->nr);
222     out  = open_trx(ftp2fn(efTRO, NFILE, fnm), "w");
223     jran = (unsigned long)((real)im*rando(&seed));
224     for (s = 0; s < nstruct; s++)
225     {
226         for (i = 0; i < natoms; i++)
227         {
228             copy_rvec(xav[i], xout1[i]);
229         }
230         for (j = 0; j < noutvec; j++)
231         {
232             v = outvec[j];
233             /* (r-0.5) n times:  var_n = n * var_1 = n/12
234                n=4:  var_n = 1/3, so multiply with 3 */
235
236             rfac  = sqrt(3.0 * BOLTZ*temp/eigval[iout[j]]);
237             rhalf = 2.0*rfac;
238             rfac  = rfac/(real)im;
239
240             jran = (jran*ia+ic) & im;
241             jr   = (real)jran;
242             jran = (jran*ia+ic) & im;
243             jr  += (real)jran;
244             jran = (jran*ia+ic) & im;
245             jr  += (real)jran;
246             jran = (jran*ia+ic) & im;
247             jr  += (real)jran;
248             disp = rfac * jr - rhalf;
249
250             for (i = 0; i < natoms; i++)
251             {
252                 for (d = 0; d < DIM; d++)
253                 {
254                     xout1[i][d] += disp*eigvec[v][i][d]*invsqrtm[i];
255                 }
256             }
257         }
258         for (i = 0; i < natoms; i++)
259         {
260             copy_rvec(xout1[i], xout2[index[i]]);
261         }
262         t = s+1;
263         write_trx(out, natoms, index, atoms, 0, t, box, xout2, NULL, NULL);
264         fprintf(stderr, "\rGenerated %d structures", s+1);
265     }
266     fprintf(stderr, "\n");
267     close_trx(out);
268
269     return 0;
270 }