Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_nmens.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include <string.h>
43
44 #include "gromacs/commandline/pargs.h"
45 #include "sysstuff.h"
46 #include "typedefs.h"
47 #include "smalloc.h"
48 #include "macros.h"
49 #include "gmx_fatal.h"
50 #include "vec.h"
51 #include "gromacs/fileio/futil.h"
52 #include "index.h"
53 #include "gromacs/fileio/pdbio.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/fileio/trxio.h"
56 #include "txtdump.h"
57 #include "physics.h"
58 #include "gromacs/random/random.h"
59 #include "eigio.h"
60 #include "gmx_ana.h"
61
62
63 int gmx_nmens(int argc, char *argv[])
64 {
65     const char *desc[] = {
66         "[THISMODULE] generates an ensemble around an average structure",
67         "in a subspace that is defined by a set of normal modes (eigenvectors).",
68         "The eigenvectors are assumed to be mass-weighted.",
69         "The position along each eigenvector is randomly taken from a Gaussian",
70         "distribution with variance kT/eigenvalue.[PAR]",
71         "By default the starting eigenvector is set to 7, since the first six",
72         "normal modes are the translational and rotational degrees of freedom."
73     };
74     static int  nstruct = 100, first = 7, last = -1, seed = -1;
75     static real temp    = 300.0;
76     t_pargs     pa[]    = {
77         { "-temp",  FALSE, etREAL, {&temp},
78           "Temperature in Kelvin" },
79         { "-seed", FALSE, etINT, {&seed},
80           "Random seed, -1 generates a seed from time and pid" },
81         { "-num", FALSE, etINT, {&nstruct},
82           "Number of structures to generate" },
83         { "-first", FALSE, etINT, {&first},
84           "First eigenvector to use (-1 is select)" },
85         { "-last",  FALSE, etINT, {&last},
86           "Last eigenvector to use (-1 is till the last)" }
87     };
88 #define NPA asize(pa)
89
90     t_trxstatus        *out;
91     int                 status, trjout;
92     t_topology          top;
93     int                 ePBC;
94     t_atoms            *atoms;
95     rvec               *xtop, *xref, *xav, *xout1, *xout2;
96     gmx_bool            bDMR, bDMA, bFit;
97     int                 nvec, *eignr = NULL;
98     rvec              **eigvec = NULL;
99     matrix              box;
100     real               *eigval, totmass, *invsqrtm, t, disp;
101     int                 natoms, neigval;
102     char               *grpname, title[STRLEN];
103     const char         *indexfile;
104     int                 i, j, d, s, v;
105     int                 nout, *iout, noutvec, *outvec;
106     atom_id            *index;
107     real                rfac, invfr, rhalf, jr;
108     int          *      eigvalnr;
109     output_env_t        oenv;
110     gmx_rng_t           rng;
111     unsigned long       jran;
112     const unsigned long im = 0xffff;
113     const unsigned long ia = 1093;
114     const unsigned long ic = 18257;
115
116
117     t_filenm fnm[] = {
118         { efTRN, "-v",    "eigenvec",    ffREAD  },
119         { efXVG, "-e",    "eigenval",    ffREAD  },
120         { efTPS, NULL,    NULL,          ffREAD },
121         { efNDX, NULL,    NULL,          ffOPTRD },
122         { efTRO, "-o",    "ensemble",    ffWRITE }
123     };
124 #define NFILE asize(fnm)
125
126     if (!parse_common_args(&argc, argv, PCA_BE_NICE,
127                            NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
128     {
129         return 0;
130     }
131
132     indexfile = ftp2fn_null(efNDX, NFILE, fnm);
133
134     read_eigenvectors(opt2fn("-v", NFILE, fnm), &natoms, &bFit,
135                       &xref, &bDMR, &xav, &bDMA, &nvec, &eignr, &eigvec, &eigval);
136
137     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, bDMA);
138     atoms = &top.atoms;
139
140     printf("\nSelect an index group of %d elements that corresponds to the eigenvectors\n", natoms);
141     get_index(atoms, indexfile, 1, &i, &index, &grpname);
142     if (i != natoms)
143     {
144         gmx_fatal(FARGS, "you selected a group with %d elements instead of %d",
145                   i, natoms);
146     }
147     printf("\n");
148
149     snew(invsqrtm, natoms);
150     if (bDMA)
151     {
152         for (i = 0; (i < natoms); i++)
153         {
154             invsqrtm[i] = gmx_invsqrt(atoms->atom[index[i]].m);
155         }
156     }
157     else
158     {
159         for (i = 0; (i < natoms); i++)
160         {
161             invsqrtm[i] = 1.0;
162         }
163     }
164
165     if (last == -1)
166     {
167         last = natoms*DIM;
168     }
169     if (first > -1)
170     {
171         /* make an index from first to last */
172         nout = last-first+1;
173         snew(iout, nout);
174         for (i = 0; i < nout; i++)
175         {
176             iout[i] = first-1+i;
177         }
178     }
179     else
180     {
181         printf("Select eigenvectors for output, end your selection with 0\n");
182         nout = -1;
183         iout = NULL;
184         do
185         {
186             nout++;
187             srenew(iout, nout+1);
188             if (1 != scanf("%d", &iout[nout]))
189             {
190                 gmx_fatal(FARGS, "Error reading user input");
191             }
192             iout[nout]--;
193         }
194         while (iout[nout] >= 0);
195         printf("\n");
196     }
197
198     /* make an index of the eigenvectors which are present */
199     snew(outvec, nout);
200     noutvec = 0;
201     for (i = 0; i < nout; i++)
202     {
203         j = 0;
204         while ((j < nvec) && (eignr[j] != iout[i]))
205         {
206             j++;
207         }
208         if ((j < nvec) && (eignr[j] == iout[i]))
209         {
210             outvec[noutvec] = j;
211             iout[noutvec]   = iout[i];
212             noutvec++;
213         }
214     }
215
216     fprintf(stderr, "%d eigenvectors selected for output\n", noutvec);
217
218     if (seed == -1)
219     {
220         seed = (int)gmx_rng_make_seed();
221         rng  = gmx_rng_init(seed);
222     }
223     else
224     {
225         rng = gmx_rng_init(seed);
226     }
227     fprintf(stderr, "Using seed %d and a temperature of %g K\n", seed, temp);
228
229     snew(xout1, natoms);
230     snew(xout2, atoms->nr);
231     out  = open_trx(ftp2fn(efTRO, NFILE, fnm), "w");
232     jran = (unsigned long)((real)im*gmx_rng_uniform_real(rng));
233     gmx_rng_destroy(rng);
234     for (s = 0; s < nstruct; s++)
235     {
236         for (i = 0; i < natoms; i++)
237         {
238             copy_rvec(xav[i], xout1[i]);
239         }
240         for (j = 0; j < noutvec; j++)
241         {
242             v = outvec[j];
243             /* (r-0.5) n times:  var_n = n * var_1 = n/12
244                n=4:  var_n = 1/3, so multiply with 3 */
245
246             rfac  = sqrt(3.0 * BOLTZ*temp/eigval[iout[j]]);
247             rhalf = 2.0*rfac;
248             rfac  = rfac/(real)im;
249
250             jran = (jran*ia+ic) & im;
251             jr   = (real)jran;
252             jran = (jran*ia+ic) & im;
253             jr  += (real)jran;
254             jran = (jran*ia+ic) & im;
255             jr  += (real)jran;
256             jran = (jran*ia+ic) & im;
257             jr  += (real)jran;
258             disp = rfac * jr - rhalf;
259
260             for (i = 0; i < natoms; i++)
261             {
262                 for (d = 0; d < DIM; d++)
263                 {
264                     xout1[i][d] += disp*eigvec[v][i][d]*invsqrtm[i];
265                 }
266             }
267         }
268         for (i = 0; i < natoms; i++)
269         {
270             copy_rvec(xout1[i], xout2[index[i]]);
271         }
272         t = s+1;
273         write_trx(out, natoms, index, atoms, 0, t, box, xout2, NULL, NULL);
274         fprintf(stderr, "\rGenerated %d structures", s+1);
275     }
276     fprintf(stderr, "\n");
277     close_trx(out);
278
279     return 0;
280 }