3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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
33 * Green Red Orange Magenta Azure Cyan Skyblue
43 #include "gmx_fatal.h"
59 static int search_str2(int nstr
,char **str
,char *key
)
62 int keylen
= strlen(key
);
65 while( (n
<keylen
) && ((key
[n
]<'0') || (key
[n
]>'9')) )
67 for(i
=0; (i
<nstr
); i
++)
68 if (gmx_strncasecmp(str
[i
],key
,n
)==0)
74 int gmx_enemat(int argc
,char *argv
[])
76 const char *desc
[] = {
77 "g_enemat extracts an energy matrix from the energy file ([TT]-f[tt]).",
78 "With [TT]-groups[tt] a file must be supplied with on each",
79 "line a group of atoms to be used. For these groups matrix of",
80 "interaction energies will be extracted from the energy file",
81 "by looking for energy groups with names corresponding to pairs",
82 "of groups of atoms. E.g. if your [TT]-groups[tt] file contains:[BR]",
84 "[TT]Protein[tt][BR]",
86 "then energy groups with names like 'Coul-SR:Protein-SOL' and ",
87 "'LJ:Protein-SOL' are expected in the energy file (although",
88 "[TT]g_enemat[tt] is most useful if many groups are analyzed",
89 "simultaneously). Matrices for different energy types are written",
90 "out separately, as controlled by the",
91 "[TT]-[no]coul[tt], [TT]-[no]coulr[tt], [TT]-[no]coul14[tt], ",
92 "[TT]-[no]lj[tt], [TT]-[no]lj14[tt], ",
93 "[TT]-[no]bham[tt] and [TT]-[no]free[tt] options.",
94 "Finally, the total interaction energy energy per group can be ",
95 "calculated ([TT]-etot[tt]).[PAR]",
97 "An approximation of the free energy can be calculated using:",
98 "E(free) = E0 + kT log( <exp((E-E0)/kT)> ), where '<>'",
99 "stands for time-average. A file with reference free energies",
100 "can be supplied to calculate the free energy difference",
101 "with some reference state. Group names (e.g. residue names)",
102 "in the reference file should correspond to the group names",
103 "as used in the [TT]-groups[tt] file, but a appended number",
104 "(e.g. residue number) in the [TT]-groups[tt] will be ignored",
107 static gmx_bool bSum
=FALSE
;
108 static gmx_bool bMeanEmtx
=TRUE
;
109 static int skip
=0,nlevels
=20;
110 static real cutmax
=1e20
,cutmin
=-1e20
,reftemp
=300.0;
111 static gmx_bool bCoulSR
=TRUE
,bCoulLR
=FALSE
,bCoul14
=FALSE
;
112 static gmx_bool bLJSR
=TRUE
,bLJLR
=FALSE
,bLJ14
=FALSE
,bBhamSR
=FALSE
,bBhamLR
=FALSE
,
115 { "-sum", FALSE
, etBOOL
, {&bSum
},
116 "Sum the energy terms selected rather than display them all" },
117 { "-skip", FALSE
, etINT
, {&skip
},
118 "Skip number of frames between data points" },
119 { "-mean", FALSE
, etBOOL
, {&bMeanEmtx
},
120 "with -groups extracts matrix of mean energies instead of "
121 "matrix for each timestep" },
122 { "-nlevels", FALSE
, etINT
, {&nlevels
},"number of levels for matrix colors"},
123 { "-max",FALSE
, etREAL
, {&cutmax
},"max value for energies"},
124 { "-min",FALSE
, etREAL
, {&cutmin
},"min value for energies"},
125 { "-coul", FALSE
, etBOOL
, {&bCoulSR
},"extract Coulomb SR energies"},
126 { "-coulr", FALSE
, etBOOL
, {&bCoulLR
},"extract Coulomb LR energies"},
127 { "-coul14",FALSE
, etBOOL
, {&bCoul14
},"extract Coulomb 1-4 energies"},
128 { "-lj", FALSE
, etBOOL
, {&bLJSR
},"extract Lennard-Jones SR energies"},
129 { "-lj", FALSE
, etBOOL
, {&bLJLR
},"extract Lennard-Jones LR energies"},
130 { "-lj14",FALSE
, etBOOL
, {&bLJ14
},"extract Lennard-Jones 1-4 energies"},
131 { "-bhamsr",FALSE
, etBOOL
, {&bBhamSR
},"extract Buckingham SR energies"},
132 { "-bhamlr",FALSE
, etBOOL
, {&bBhamLR
},"extract Buckingham LR energies"},
133 { "-free",FALSE
, etBOOL
, {&bFree
},"calculate free energy"},
134 { "-temp",FALSE
, etREAL
, {&reftemp
},
135 "reference temperature for free energy calculation"}
137 /* We will define egSP more energy-groups:
138 egTotal (total energy) */
141 gmx_bool egrp_use
[egNR
+egSP
];
145 gmx_enxnm_t
*enm
=NULL
;
150 gmx_bool bCutmax
,bCutmin
;
151 real
**eneset
,*time
=NULL
;
152 int *set
,i
,j
,k
,prevk
,m
=0,n
,nre
,nset
,nenergy
;
153 char **groups
= NULL
;
154 char groupname
[255],fn
[255];
158 real
***emat
,**etot
,*groupnr
;
159 double beta
,expE
,**e
,*eaver
,*efree
=NULL
,edum
;
161 char **ereflines
,**erefres
=NULL
;
162 real
*eref
=NULL
,*edif
=NULL
;
167 { efEDR
, "-f", NULL
, ffOPTRD
},
168 { efDAT
, "-groups", "groups.dat", ffREAD
},
169 { efDAT
, "-eref", "eref.dat", ffOPTRD
},
170 { efXPM
, "-emat", "emat", ffWRITE
},
171 { efXVG
, "-etot", "energy", ffWRITE
}
173 #define NFILE asize(fnm)
175 CopyRight(stderr
,argv
[0]);
176 parse_common_args(&argc
,argv
,PCA_CAN_VIEW
| PCA_CAN_TIME
| PCA_BE_NICE
,
177 NFILE
,fnm
,asize(pa
),pa
,asize(desc
),desc
,0,NULL
,&oenv
);
179 egrp_use
[egCOULSR
]=bCoulSR
;
180 egrp_use
[egLJSR
]=bLJSR
;
181 egrp_use
[egBHAMSR
]=bBhamSR
;
182 egrp_use
[egCOULLR
]=bCoulLR
;
183 egrp_use
[egLJLR
]=bLJLR
;
184 egrp_use
[egBHAMLR
]=bBhamLR
;
185 egrp_use
[egCOUL14
]=bCoul14
;
186 egrp_use
[egLJ14
]=bLJ14
;
187 egrp_use
[egTotal
]=TRUE
;
189 bRef
=opt2bSet("-eref",NFILE
,fnm
);
190 in
=open_enx(ftp2fn(efEDR
,NFILE
,fnm
),"r");
191 do_enxnms(in
,&nre
,&enm
);
194 gmx_fatal(FARGS
,"No energies!\n");
196 bCutmax
=opt2parg_bSet("-max",asize(pa
),pa
);
197 bCutmin
=opt2parg_bSet("-min",asize(pa
),pa
);
201 /* Read groupnames from input file and construct selection of
202 energy groups from it*/
204 fprintf(stderr
,"Will read groupnames from inputfile\n");
205 ngroups
= get_lines(opt2fn("-groups",NFILE
,fnm
), &groups
);
206 fprintf(stderr
,"Read %d groups\n",ngroups
);
207 snew(set
,sqr(ngroups
)*egNR
/2);
210 for (i
=0; (i
<ngroups
); i
++) {
211 fprintf(stderr
,"\rgroup %d",i
);
212 for (j
=i
; (j
<ngroups
); j
++)
213 for(m
=0; (m
<egNR
); m
++)
215 sprintf(groupname
,"%s:%s-%s",egrp_nm
[m
],groups
[i
],groups
[j
]);
217 fprintf(stderr
,"\r%-15s %5d",groupname
,n
);
219 for(k
=prevk
; (k
<prevk
+nre
); k
++)
220 if (strcmp(enm
[k
%nre
].name
,groupname
) == 0) {
225 fprintf(stderr
,"WARNING! could not find group %s (%d,%d)"
226 "in energy file\n",groupname
,i
,j
);
231 fprintf(stderr
,"\n");
234 fprintf(stderr
,"Will select half-matrix of energies with %d elements\n",n
);
236 /* Start reading energy frames */
242 timecheck
=check_times(fr
->t
);
243 } while (bCont
&& (timecheck
< 0));
245 if (timecheck
== 0) {
246 #define DONTSKIP(cnt) (skip) ? ((cnt % skip) == 0) : TRUE
249 fprintf(stderr
,"\rRead frame: %d, Time: %.3f",teller
,fr
->t
);
251 if ((nenergy
% 1000) == 0) {
252 srenew(time
,nenergy
+1000);
253 for(i
=0; (i
<=nset
); i
++)
254 srenew(eneset
[i
],nenergy
+1000);
256 time
[nenergy
] = fr
->t
;
258 for(i
=0; (i
<nset
); i
++) {
259 eneset
[i
][nenergy
] = fr
->ener
[set
[i
]].e
;
260 sum
+= fr
->ener
[set
[i
]].e
;
263 eneset
[nset
][nenergy
] = sum
;
268 } while (bCont
&& (timecheck
== 0));
270 fprintf(stderr
,"\n");
272 fprintf(stderr
,"Will build energy half-matrix of %d groups, %d elements, "
273 "over %d frames\n",ngroups
,nset
,nenergy
);
275 snew(emat
,egNR
+egSP
);
276 for(j
=0; (j
<egNR
+egSP
); j
++)
278 snew(emat
[j
],ngroups
);
279 for (i
=0; (i
<ngroups
); i
++)
280 snew(emat
[j
][i
],ngroups
);
282 snew(groupnr
,ngroups
);
283 for (i
=0; (i
<ngroups
); i
++)
285 rlo
.r
= 1.0, rlo
.g
= 0.0, rlo
.b
= 0.0;
286 rmid
.r
= 1.0, rmid
.g
= 1.0, rmid
.b
= 1.0;
287 rhi
.r
= 0.0, rhi
.g
= 0.0, rhi
.b
= 1.0;
290 for (i
=0; (i
<ngroups
); i
++)
293 for (i
=0; (i
<ngroups
); i
++) {
294 for (j
=i
; (j
<ngroups
); j
++) {
295 for (m
=0; (m
<egNR
); m
++) {
297 for (k
=0; (k
<nenergy
); k
++) {
298 emat
[m
][i
][j
] += eneset
[n
][k
];
299 e
[i
][k
] += eneset
[n
][k
];/* *0.5; */
300 e
[j
][k
] += eneset
[n
][k
];/* *0.5; */
303 emat
[egTotal
][i
][j
]+=emat
[m
][i
][j
];
304 emat
[m
][i
][j
]/=nenergy
;
305 emat
[m
][j
][i
]=emat
[m
][i
][j
];
308 emat
[egTotal
][i
][j
]/=nenergy
;
309 emat
[egTotal
][j
][i
]=emat
[egTotal
][i
][j
];
314 fprintf(stderr
,"Will read reference energies from inputfile\n");
315 neref
= get_lines(opt2fn("-eref",NFILE
,fnm
), &ereflines
);
316 fprintf(stderr
,"Read %d reference energies\n",neref
);
318 snew(erefres
, neref
);
319 for(i
=0; (i
<neref
); i
++) {
321 sscanf(ereflines
[i
],"%s %lf",erefres
[i
],&edum
);
326 for (i
=0; (i
<ngroups
); i
++) {
327 for (k
=0; (k
<nenergy
); k
++)
331 beta
= 1.0/(BOLTZ
*reftemp
);
334 for (i
=0; (i
<ngroups
); i
++) {
336 for (k
=0; (k
<nenergy
); k
++) {
337 expE
+= exp(beta
*(e
[i
][k
]-eaver
[i
]));
339 efree
[i
] = log(expE
/nenergy
)/beta
+ eaver
[i
];
341 n
= search_str2(neref
,erefres
,groups
[i
]);
343 edif
[i
] = efree
[i
]-eref
[n
];
346 fprintf(stderr
,"WARNING: group %s not found "
347 "in reference energies.\n",groups
[i
]);
354 emid
= 0.0;/*(emin+emax)*0.5;*/
355 for(m
=0; (m
<egNR
); m
++)
356 egrp_nm
[m
]=egrp_nm
[m
];
357 egrp_nm
[egTotal
]="total";
358 for (m
=0; (m
<egNR
+egSP
); m
++)
362 for (i
=0; (i
<ngroups
); i
++) {
363 for (j
=i
; (j
<ngroups
); j
++) {
364 if (emat
[m
][i
][j
] > emax
)
366 else if (emat
[m
][i
][j
] < emin
)
371 fprintf(stderr
,"Matrix of %s energy is uniform at %f "
372 "(will not produce output).\n",egrp_nm
[m
],emax
);
374 fprintf(stderr
,"Matrix of %s energy ranges from %f to %f\n",
375 egrp_nm
[m
],emin
,emax
);
376 if ((bCutmax
) || (emax
>cutmax
)) emax
=cutmax
;
377 if ((bCutmin
) || (emin
<cutmin
)) emin
=cutmin
;
378 if ((emax
==cutmax
) || (emin
==cutmin
))
379 fprintf(stderr
,"Energy range adjusted: %f to %f\n",emin
,emax
);
381 sprintf(fn
,"%s%s",egrp_nm
[m
],ftp2fn(efXPM
,NFILE
,fnm
));
382 sprintf(label
,"%s Interaction Energies",egrp_nm
[m
]);
385 write_xpm(out
,0,label
,"Energy (kJ/mol)",
386 "Residue Index","Residue Index",
387 ngroups
,ngroups
,groupnr
,groupnr
,emat
[m
],
388 emid
,emax
,rmid
,rhi
,&nlevels
);
390 write_xpm(out
,0,label
,"Energy (kJ/mol)",
391 "Residue Index","Residue Index",
392 ngroups
,ngroups
,groupnr
,groupnr
,emat
[m
],
393 emin
,emid
,rlo
,rmid
,&nlevels
);
395 write_xpm3(out
,0,label
,"Energy (kJ/mol)",
396 "Residue Index","Residue Index",
397 ngroups
,ngroups
,groupnr
,groupnr
,emat
[m
],
398 emin
,emid
,emax
,rlo
,rmid
,rhi
,&nlevels
);
402 snew(etot
,egNR
+egSP
);
403 for (m
=0; (m
<egNR
+egSP
); m
++) {
404 snew(etot
[m
],ngroups
);
405 for (i
=0; (i
<ngroups
); i
++) {
406 for (j
=0; (j
<ngroups
); j
++)
407 etot
[m
][i
]+=emat
[m
][i
][j
];
411 out
=xvgropen(ftp2fn(efXVG
,NFILE
,fnm
),"Mean Energy","Residue","kJ/mol",
413 xvgr_legend(out
,0,NULL
, oenv
);
415 for (m
=0; (m
<egNR
+egSP
); m
++)
417 fprintf(out
,"@ legend string %d \"%s\"\n",j
++,egrp_nm
[m
]);
419 fprintf(out
,"@ legend string %d \"%s\"\n",j
++,"Free");
421 fprintf(out
,"@ legend string %d \"%s\"\n",j
++,"Diff");
422 fprintf(out
,"@TYPE xy\n");
423 fprintf(out
,"#%3s","grp");
424 for (m
=0; (m
<egNR
+egSP
); m
++)
426 fprintf(out
," %9s",egrp_nm
[m
]);
428 fprintf(out
," %9s","Free");
430 fprintf(out
," %9s","Diff");
432 for (i
=0; (i
<ngroups
); i
++) {
433 fprintf(out
,"%3.0f",groupnr
[i
]);
434 for (m
=0; (m
<egNR
+egSP
); m
++)
436 fprintf(out
," %9.5g",etot
[m
][i
]);
438 fprintf(out
," %9.5g",efree
[i
]);
440 fprintf(out
," %9.5g",edif
[i
]);
445 fprintf(stderr
,"While typing at your keyboard, suddenly...\n"
446 "...nothing happens.\nWARNING: Not Implemented Yet\n");
448 out=ftp2FILE(efMAT,NFILE,fnm,"w");
451 for (k=0; (k<nenergy); k++) {
452 for (i=0; (i<ngroups); i++)
453 for (j=i+1; (j<ngroups); j++)
454 emat[i][j]=eneset[n][k];
455 sprintf(label,"t=%.0f ps",time[k]);
456 write_matrix(out,ngroups,1,ngroups,groupnr,emat,label,emin,emax,nlevels);