Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / tools / 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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43 #include <string.h>
44
45 #include "statutil.h"
46 #include "sysstuff.h"
47 #include "typedefs.h"
48 #include "smalloc.h"
49 #include "macros.h"
50 #include "gmx_fatal.h"
51 #include "vec.h"
52 #include "copyrite.h"
53 #include "futil.h"
54 #include "statutil.h"
55 #include "index.h"
56 #include "pdbio.h"
57 #include "tpxio.h"
58 #include "txtdump.h"
59 #include "physics.h"
60 #include "random.h"
61 #include "eigio.h"
62 #include "gmx_ana.h"
63
64
65 int gmx_nmens(int argc,char *argv[])
66 {
67   const char *desc[] = {
68     "[TT]g_nmens[tt] generates an ensemble around an average structure",
69     "in a subspace that is defined by a set of normal modes (eigenvectors).",
70     "The eigenvectors are assumed to be mass-weighted.",
71     "The position along each eigenvector is randomly taken from a Gaussian",
72     "distribution with variance kT/eigenvalue.[PAR]",
73     "By default the starting eigenvector is set to 7, since the first six",
74     "normal modes are the translational and rotational degrees of freedom." 
75   };
76   static int  nstruct=100,first=7,last=-1,seed=-1;
77   static real temp=300.0;
78   t_pargs pa[] = {
79     { "-temp",  FALSE, etREAL, {&temp}, 
80       "Temperature in Kelvin" },
81     { "-seed", FALSE, etINT, {&seed},     
82       "Random seed, -1 generates a seed from time and pid" },
83     { "-num", FALSE, etINT, {&nstruct},     
84       "Number of structures to generate" },
85     { "-first", FALSE, etINT, {&first},     
86       "First eigenvector to use (-1 is select)" },
87     { "-last",  FALSE, etINT, {&last}, 
88       "Last eigenvector to use (-1 is till the last)" }
89   };
90 #define NPA asize(pa)
91   
92   t_trxstatus *out;
93   int        status,trjout;
94   t_topology top;
95   int        ePBC;
96   t_atoms    *atoms;
97   rvec       *xtop,*xref,*xav,*xout1,*xout2;
98   gmx_bool       bDMR,bDMA,bFit;
99   int        nvec,*eignr=NULL;
100   rvec       **eigvec=NULL;
101   matrix     box;
102   real       *eigval,totmass,*invsqrtm,t,disp;
103   int        natoms,neigval;
104   char       *grpname,title[STRLEN];
105   const char *indexfile;
106   int        i,j,d,s,v;
107   int        nout,*iout,noutvec,*outvec;
108   atom_id    *index;
109   real       rfac,invfr,rhalf,jr;
110   int *      eigvalnr;
111   output_env_t oenv;
112   
113   unsigned long      jran;
114   const unsigned long im = 0xffff;
115   const unsigned long ia = 1093;
116   const unsigned long ic = 18257;
117
118
119   t_filenm fnm[] = { 
120     { efTRN, "-v",    "eigenvec",    ffREAD  },
121     { efXVG, "-e",    "eigenval",    ffREAD  },
122     { efTPS, NULL,    NULL,          ffREAD },
123     { efNDX, NULL,    NULL,          ffOPTRD },
124     { efTRO, "-o",    "ensemble",    ffWRITE }
125   }; 
126 #define NFILE asize(fnm) 
127
128   CopyRight(stderr,argv[0]); 
129   parse_common_args(&argc,argv,PCA_BE_NICE,
130                     NFILE,fnm,NPA,pa,asize(desc),desc,0,NULL,&oenv); 
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     gmx_fatal(FARGS,"you selected a group with %d elements instead of %d",
144                 i,natoms);
145   printf("\n");
146   
147   snew(invsqrtm,natoms);
148   if (bDMA) {
149     for(i=0; (i<natoms); i++)
150       invsqrtm[i] = gmx_invsqrt(atoms->atom[index[i]].m);
151   } else {
152     for(i=0; (i<natoms); i++)
153       invsqrtm[i]=1.0;
154   }
155   
156   if (last==-1)
157       last=natoms*DIM;
158   if (first>-1) 
159   {
160       /* make an index from first to last */
161       nout=last-first+1;
162       snew(iout,nout);
163       for(i=0; i<nout; i++)
164           iout[i]=first-1+i;
165   }
166   else 
167   {
168       printf("Select eigenvectors for output, end your selection with 0\n");
169       nout=-1;
170       iout=NULL;
171       do {
172           nout++;
173           srenew(iout,nout+1);
174           if(1 != scanf("%d",&iout[nout]))
175           {
176               gmx_fatal(FARGS,"Error reading user input");
177           }
178           iout[nout]--;
179       } while (iout[nout]>=0);
180       printf("\n");
181   }
182   
183   /* make an index of the eigenvectors which are present */
184   snew(outvec,nout);
185   noutvec=0;
186   for(i=0; i<nout; i++)
187   {
188       j=0;
189       while ((j<nvec) && (eignr[j]!=iout[i]))
190           j++;
191       if ((j<nvec) && (eignr[j]==iout[i]))
192       {
193           outvec[noutvec] = j;
194           iout[noutvec] = iout[i];
195           noutvec++;
196       }
197   }
198   
199   fprintf(stderr,"%d eigenvectors selected for output\n",noutvec);
200
201   if (seed == -1)
202     seed = make_seed();
203   fprintf(stderr,"Using seed %d and a temperature of %g K\n",seed,temp);
204
205   snew(xout1,natoms);
206   snew(xout2,atoms->nr);
207   out=open_trx(ftp2fn(efTRO,NFILE,fnm),"w");
208   jran = (unsigned long)((real)im*rando(&seed));
209   for(s=0; s<nstruct; s++) {
210     for(i=0; i<natoms; i++)
211       copy_rvec(xav[i],xout1[i]);
212     for(j=0; j<noutvec; j++) {
213       v = outvec[j];
214       /* (r-0.5) n times:  var_n = n * var_1 = n/12
215          n=4:  var_n = 1/3, so multiply with 3 */
216       
217       rfac  = sqrt(3.0 * BOLTZ*temp/eigval[iout[j]]);
218       rhalf = 2.0*rfac; 
219       rfac  = rfac/(real)im;
220
221       jran = (jran*ia+ic) & im;
222       jr = (real)jran;
223       jran = (jran*ia+ic) & im;
224       jr += (real)jran;
225       jran = (jran*ia+ic) & im;
226       jr += (real)jran;
227       jran = (jran*ia+ic) & im;
228       jr += (real)jran;
229       disp = rfac * jr - rhalf;
230       
231       for(i=0; i<natoms; i++)
232           for(d=0; d<DIM; d++)
233               xout1[i][d] += disp*eigvec[v][i][d]*invsqrtm[i];
234     }
235     for(i=0; i<natoms; i++)
236         copy_rvec(xout1[i],xout2[index[i]]);
237     t = s+1;
238     write_trx(out,natoms,index,atoms,0,t,box,xout2,NULL,NULL);
239     fprintf(stderr,"\rGenerated %d structures",s+1);
240   }
241   fprintf(stderr,"\n");
242   close_trx(out);
243   
244   return 0;
245 }
246