Fixed make_edi.c
[gromacs/rigid-bodies.git] / src / tools / gmx_nmens.c
blob160847727c98af651861c8bc8313bc1747e088c5
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Green Red Orange Magenta Azure Cyan Skyblue
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include <math.h>
40 #include <string.h>
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"
62 int gmx_nmens(int argc,char *argv[])
64 const char *desc[] = {
65 "[TT]g_nmens[tt] generates an ensemble around an average structure",
66 "in a subspace which 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."
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)" }
87 #define NPA asize(pa)
89 int out;
90 int status,trjout;
91 t_topology top;
92 int ePBC;
93 t_atoms *atoms;
94 rvec *xtop,*xref,*xav,*xout1,*xout2;
95 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;
110 unsigned long jran;
111 const unsigned long im = 0xffff;
112 const unsigned long ia = 1093;
113 const unsigned long ic = 18257;
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 }
123 #define NFILE asize(fnm)
125 CopyRight(stderr,argv[0]);
126 parse_common_args(&argc,argv,PCA_BE_NICE,
127 NFILE,fnm,NPA,pa,asize(desc),desc,0,NULL,&oenv);
129 indexfile=ftp2fn_null(efNDX,NFILE,fnm);
131 read_eigenvectors(opt2fn("-v",NFILE,fnm),&natoms,&bFit,
132 &xref,&bDMR,&xav,&bDMA,&nvec,&eignr,&eigvec,&eigval);
134 read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,&xtop,NULL,box,bDMA);
135 atoms=&top.atoms;
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 gmx_fatal(FARGS,"you selected a group with %d elements instead of %d",
141 i,natoms);
142 printf("\n");
144 snew(invsqrtm,natoms);
145 if (bDMA) {
146 for(i=0; (i<natoms); i++)
147 invsqrtm[i] = gmx_invsqrt(atoms->atom[index[i]].m);
148 } else {
149 for(i=0; (i<natoms); i++)
150 invsqrtm[i]=1.0;
153 if (last==-1)
154 last=natoms*DIM;
155 if (first>-1)
157 /* make an index from first to last */
158 nout=last-first+1;
159 snew(iout,nout);
160 for(i=0; i<nout; i++)
161 iout[i]=first-1+i;
163 else
165 printf("Select eigenvectors for output, end your selection with 0\n");
166 nout=-1;
167 iout=NULL;
168 do {
169 nout++;
170 srenew(iout,nout+1);
171 if(1 != scanf("%d",&iout[nout]))
173 gmx_fatal(FARGS,"Error reading user input");
175 iout[nout]--;
176 } while (iout[nout]>=0);
177 printf("\n");
180 /* make an index of the eigenvectors which are present */
181 snew(outvec,nout);
182 noutvec=0;
183 for(i=0; i<nout; i++)
185 j=0;
186 while ((j<nvec) && (eignr[j]!=iout[i]))
187 j++;
188 if ((j<nvec) && (eignr[j]==iout[i]))
190 outvec[noutvec] = j;
191 iout[noutvec] = iout[i];
192 noutvec++;
196 fprintf(stderr,"%d eigenvectors selected for output\n",noutvec);
198 if (seed == -1)
199 seed = make_seed();
200 fprintf(stderr,"Using seed %d and a temperature of %g K\n",seed,temp);
202 snew(xout1,natoms);
203 snew(xout2,atoms->nr);
204 out=open_trx(ftp2fn(efTRO,NFILE,fnm),"w");
205 jran = (unsigned long)((real)im*rando(&seed));
206 for(s=0; s<nstruct; s++) {
207 for(i=0; i<natoms; i++)
208 copy_rvec(xav[i],xout1[i]);
209 for(j=0; j<noutvec; j++) {
210 v = outvec[j];
211 /* (r-0.5) n times: var_n = n * var_1 = n/12
212 n=4: var_n = 1/3, so multiply with 3 */
214 rfac = sqrt(3.0 * BOLTZ*temp/eigval[iout[j]]);
215 rhalf = 2.0*rfac;
216 rfac = rfac/(real)im;
218 jran = (jran*ia+ic) & im;
219 jr = (real)jran;
220 jran = (jran*ia+ic) & im;
221 jr += (real)jran;
222 jran = (jran*ia+ic) & im;
223 jr += (real)jran;
224 jran = (jran*ia+ic) & im;
225 jr += (real)jran;
226 disp = rfac * jr - rhalf;
228 for(i=0; i<natoms; i++)
229 for(d=0; d<DIM; d++)
230 xout1[i][d] += disp*eigvec[v][i][d]*invsqrtm[i];
232 for(i=0; i<natoms; i++)
233 copy_rvec(xout1[i],xout2[index[i]]);
234 t = s+1;
235 write_trx(out,natoms,index,atoms,0,t,box,xout2,NULL,NULL);
236 fprintf(stderr,"\rGenerated %d structures",s+1);
238 fprintf(stderr,"\n");
239 close_trx(out);
241 return 0;