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