Create fileio module
[alexxy/gromacs.git] / src / contrib / ehanal.c
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.3.3
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-2008, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
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.
19  * 
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.
26  * 
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.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Groningen Machine for Chemical Simulation
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <stdlib.h>
40 #include <math.h>
41 #include <string.h>
42 #include "typedefs.h"
43 #include "smalloc.h"
44 #include "macros.h"
45 #include "statutil.h"
46 #include "gmx_fatal.h"
47 #include "random.h"
48 #include "gromacs/fileio/pdbio.h"
49 #include "gromacs/fileio/futil.h"
50 #include "physics.h"
51 #include "xvgr.h"
52 #include "vec.h"
53 #include "names.h"
54 #include "ehdata.h"
55 #include "gromacs/fileio/pdbio.h"
56
57 t_histo *init_histo(int np,real minx,real maxx)
58 {
59   t_histo *h;
60   
61   snew(h,1);
62   snew(h->y,np+1);
63   snew(h->nh,np+1);
64   h->np   = np;
65   if (maxx <= minx)
66     gmx_fatal(FARGS,"minx (%f) should be less than maxx (%f) in init_histo",minx,maxx);
67   h->minx = minx;
68   h->maxx = maxx;
69   h->dx_1 = np/(maxx-minx);
70   h->dx   = (maxx-minx)/np;
71   
72   return h;
73 }
74
75 void done_histo(t_histo *h)
76 {
77   sfree(h->y);
78   sfree(h->nh);
79   h->np = 0;
80 }
81
82 void add_histo(t_histo *h,real x,real y)
83 {
84   int n;
85   
86   n = (x-h->minx)*h->dx_1;
87   if ((n < 0) || (n > h->np)) 
88     gmx_fatal(FARGS,"Invalid x (%f) in add_histo. Should be in %f - %f",x,h->minx,h->maxx);
89   h->y[n] += y;
90   h->nh[n]++;
91 }
92
93 void dump_histo(t_histo *h,char *fn,char *title,char *xaxis,char *yaxis,
94                 int enorm,real norm_fac)
95 {
96   FILE *fp;
97   int  i,nn;
98   
99   for(nn=h->np; (nn > 0); nn--)
100     if (h->nh[nn] != 0) 
101       break;
102   for(i=0; (i<nn); i++)
103     if (h->nh[i] > 0)
104       break;
105   fp = xvgropen(fn,title,xaxis,yaxis);
106   for(  ; (i<nn); i++) {
107     switch (enorm) {
108     case enormNO:
109       fprintf(fp,"%12f  %12f  %d\n",h->minx+h->dx*i,h->y[i],h->nh[i]);
110       break;
111     case enormFAC:
112       fprintf(fp,"%12f  %12f  %d\n",h->minx+h->dx*i,h->y[i]*norm_fac,h->nh[i]);
113       break;
114     case enormNP:
115       if (h->nh[i] > 0)
116         fprintf(fp,"%12f  %12f  %d\n",
117                 h->minx+h->dx*i,h->y[i]*norm_fac/h->nh[i],h->nh[i]);
118       break;
119     default:
120       gmx_fatal(FARGS,"Wrong value for enorm (%d)",enorm);
121     }
122   }
123   ffclose(fp);
124 }
125
126 /*******************************************************************
127  *
128  * Functions to analyse and monitor scattering
129  *
130  *******************************************************************/   
131
132 void add_scatter_event(t_ana_scat *scatter,rvec pos,gmx_bool bInel,
133                        real t,real ekin)
134 {
135   int np = scatter->np;
136   
137   if (np == scatter->maxp) {
138     scatter->maxp += 32;
139     srenew(scatter->time,scatter->maxp);
140     srenew(scatter->ekin,scatter->maxp);
141     srenew(scatter->bInel,scatter->maxp);
142     srenew(scatter->pos,scatter->maxp);
143   }
144   scatter->time[np]  = t;
145   scatter->bInel[np] = np;
146   scatter->ekin[np]  = ekin;
147   copy_rvec(pos,scatter->pos[np]);
148   scatter->np++;
149 }
150
151 void reset_ana_scat(t_ana_scat *scatter)
152 {
153   scatter->np = 0;
154 }
155
156 void done_scatter(t_ana_scat *scatter)
157 {
158   scatter->np = 0;
159   sfree(scatter->time);
160   sfree(scatter->ekin);
161   sfree(scatter->bInel);
162   sfree(scatter->pos);
163 }
164
165 void analyse_scatter(t_ana_scat *scatter,t_histo *hmfp)
166 {
167   int   i,n;
168   rvec  dx;
169   
170   if (scatter->np > 1) {
171     for(i=1; (i<scatter->np); i++) {
172       rvec_sub(scatter->pos[i],scatter->pos[i-1],dx);
173       add_histo(hmfp,scatter->ekin[i],norm(dx));
174     }
175   }
176 }
177
178 /*******************************************************************
179  *
180  * Functions to analyse structural changes
181  *
182  *******************************************************************/   
183
184 t_ana_struct *init_ana_struct(int nstep,int nsave,real timestep,
185                               int maxparticle)
186 {
187   t_ana_struct *anal;
188   
189   snew(anal,1);
190   anal->nanal = 1.2*((nstep / nsave)+1);
191   anal->index = 0;
192   anal->dt    = nsave*timestep;
193   snew(anal->t,anal->nanal);
194   snew(anal->maxdist,anal->nanal);
195   snew(anal->d2_com,anal->nanal);
196   snew(anal->d2_origin,anal->nanal);
197   snew(anal->nion,anal->nanal);
198   anal->nstruct   = 1;
199   anal->nparticle = 1;
200   anal->maxparticle = maxparticle;
201   snew(anal->q,1);
202   snew(anal->x,1);
203   snew(anal->x[0],maxparticle);
204   
205   return anal;
206 }
207
208 void done_ana_struct(t_ana_struct *anal)
209 {
210   int i;
211   
212   sfree(anal->t);
213   sfree(anal->maxdist);
214   sfree(anal->d2_com);
215   sfree(anal->d2_origin);
216   sfree(anal->nion);
217   sfree(anal->q);
218   for(i=0; (i<anal->nstruct); i++)
219     sfree(anal->x[i]);
220   sfree(anal->x);
221 }
222
223 void reset_ana_struct(t_ana_struct *anal)
224 {
225   int i;
226   
227   for(i=0; (i<anal->nanal); i++) {
228     anal->t[i] = 0;
229     anal->maxdist[i] = 0;
230     clear_rvec(anal->d2_com[i]);
231     clear_rvec(anal->d2_origin[i]);
232     anal->nion[i] = 0;
233   }
234   anal->index = 0;
235 }
236
237 void add_ana_struct(t_ana_struct *total,t_ana_struct *add)
238 {
239   int i,m;
240   
241   if (total->index == 0)
242     total->index = add->index;
243   else if (total->index != add->index)
244     gmx_fatal(FARGS,"Analysis incompatible (total: %d, add: %d) %s, %d",
245                 total->index,add->index,__FILE__,__LINE__);
246   for(i=0; (i<total->index); i++) {
247     if (total->t[i] == 0)
248       total->t[i] = add->t[i];
249     else if (total->t[i] != add->t[i])
250       gmx_fatal(FARGS,"Inconsistent times in analysis (%f-%f) %s, %d",
251                   total->t[i],add->t[i],__FILE__,__LINE__);
252     if (add->maxdist[i] > total->maxdist[i])
253       total->maxdist[i]  = add->maxdist[i];
254     for(m=0; (m<DIM); m++) {
255       total->d2_com[i][m]    += add->d2_com[i][m]/add->nion[i];
256       total->d2_origin[i][m] += add->d2_origin[i][m]/add->nion[i];
257     }
258     total->nion[i]     += add->nion[i];
259   }
260 }
261
262 static void do_add_struct(t_ana_struct *anal,int nparticle,rvec x[])
263 {
264   int i,j;
265   
266   if (nparticle > anal->nparticle) {
267     for(i=0; (i<anal->nstruct); i++) {
268       for(j=anal->nparticle; (j<nparticle); j++)
269         copy_rvec(x[j],anal->x[i][j]);
270     }
271   }
272   anal->nparticle=nparticle;
273   srenew(anal->x,anal->nstruct+1);
274   snew(anal->x[anal->nstruct],anal->maxparticle);
275   for(j=0; (j<nparticle); j++)
276     copy_rvec(x[j],anal->x[anal->nstruct][j]);
277   anal->nstruct++;
278 }
279
280 void analyse_structure(t_ana_struct *anal,real t,rvec center,
281                        rvec x[],int nparticle,real charge[])
282 {
283   int  i,j,m,nel,n=0;
284   rvec dx,com;
285   real dx2,dx1;
286   
287   j = anal->index;
288   if (j >= anal->nanal)
289     gmx_fatal(FARGS,"Too many points in analyse_structure");
290   anal->t[j]       = t;
291   anal->maxdist[j] = 0;
292   clear_rvec(com);
293   nel = 0;
294   for(i=0; (i<nparticle); i++) {
295     if (charge[i] < 0) {
296       rvec_inc(com,x[i]);
297       nel++;
298     }
299   }
300   if (nel > 0)
301     for(m=0; (m<3); m++)
302       com[m] /= nel;
303   for(i=0; (i<nparticle); i++) {
304     if (charge[i] < 0) {
305       rvec_sub(x[i],com,dx);
306       for(m=0; (m<DIM); m++) {
307         anal->d2_com[j][m]    += sqr(dx[m]);
308         anal->d2_origin[j][m] += sqr(x[i][m]);
309       }
310       dx2 = iprod(x[i],x[i]);
311       dx1 = sqrt(dx2);
312       if (dx1 > anal->maxdist[j])
313         anal->maxdist[j] = dx1;
314       n++;
315     }
316   }
317   do_add_struct(anal,nparticle,x);
318   anal->nion[j] = n;
319   anal->index++;
320 }
321
322 void dump_ana_struct(char *rmax,char *nion,char *gyr_com,char *gyr_origin,
323                      t_ana_struct *anal,int nsim)
324 {
325   FILE *fp,*gp,*hp,*kp;
326   int  i,j;
327   real t,d2;
328   char *legend[] = { "Rg", "RgX", "RgY", "RgZ" };
329   
330   fp = xvgropen(rmax,"rmax","Time (fs)","r (nm)");
331   gp = xvgropen(nion,"N ion","Time (fs)","N ions");
332   hp = xvgropen(gyr_com,"Radius of gyration wrt. C.O.M.",
333                 "Time (fs)","Rg (nm)");
334   xvgr_legend(hp,asize(legend),legend);
335   kp = xvgropen(gyr_origin,"Radius of gyration wrt. Origin",
336                 "Time (fs)","Rg (nm)");
337   xvgr_legend(kp,asize(legend),legend);
338   for(i=0; (i<anal->index); i++) {
339     t = 1000*anal->t[i];
340     fprintf(fp,"%12g  %10.3f\n",t,anal->maxdist[i]);
341     fprintf(gp,"%12g  %10.3f\n",t,(1.0*anal->nion[i])/nsim-1);
342     d2 = anal->d2_com[i][XX] + anal->d2_com[i][YY] + anal->d2_com[i][ZZ];
343     fprintf(hp,"%12g  %10.3f  %10.3f  %10.3f  %10.3f\n",
344             t,sqrt(d2/nsim),
345             sqrt(anal->d2_com[i][XX]/nsim),
346             sqrt(anal->d2_com[i][YY]/nsim),
347             sqrt(anal->d2_com[i][ZZ]/nsim));
348     d2 = anal->d2_origin[i][XX] + anal->d2_origin[i][YY] + anal->d2_origin[i][ZZ];
349     fprintf(kp,"%12g  %10.3f  %10.3f  %10.3f  %10.3f\n",
350             t,sqrt(d2/nsim),
351             sqrt(anal->d2_origin[i][XX]/nsim),
352             sqrt(anal->d2_origin[i][YY]/nsim),
353             sqrt(anal->d2_origin[i][ZZ]/nsim));
354   }
355   ffclose(hp);
356   ffclose(gp);
357   ffclose(fp);
358   ffclose(kp);
359 }
360
361 void dump_as_pdb(char *pdb,t_ana_struct *anal)
362 {
363   FILE *kp;
364   int  i,j;
365   real t;
366   
367   kp = ffopen(pdb,"w");
368   for(i=0; (i<anal->nstruct); i++) {
369     t = 1000*anal->t[i];
370     fprintf(kp,"MODEL  %d  time %g fs\n",i+1,t);
371     for(j=0; (j<anal->nparticle); j++) {
372       fprintf(kp,get_pdbformat(),"ATOM",i+1,(j < anal->nion[i]) ? "O" : "N",
373               "PLS",' ',1,
374               anal->x[i][j][XX]/100,
375               anal->x[i][j][YY]/100,
376               anal->x[i][j][ZZ]/100);
377       fprintf(kp,"\n");
378     }
379     fprintf(kp,"ENDMDL\n");
380   }
381   ffclose(kp);
382 }
383
384 char *enms[eNR] = {
385   "Coulomb", "Repulsion", "Potential",
386   "EkHole",  "EkElectron", "EkLattice", "Kinetic",
387   "Total"
388 };
389
390 void add_ana_ener(t_ana_ener *ae,int nn,real e[])
391 {
392   int i;
393  
394   /* First time around we are constantly increasing the array size */ 
395   if (nn >= ae->nx) {
396     if (ae->nx == ae->maxx) {
397       ae->maxx += 1024;
398       srenew(ae->e,ae->maxx);
399     }
400     for(i=0; (i<eNR); i++)
401       ae->e[ae->nx][i] = e[i];
402     ae->nx++;
403   }
404   else {
405     for(i=0; (i<eNR); i++)
406       ae->e[nn][i] += e[i];
407   }
408 }
409
410 void dump_ana_ener(t_ana_ener *ae,int nsim,real dt,char *edump,
411                    t_ana_struct *total)
412 {
413   FILE *fp;
414   int  i,j;
415   real fac;
416   
417   fac = 1.0/(nsim*ELECTRONVOLT);
418   fp=xvgropen(edump,"Energies","Time (fs)","E (eV)");
419   xvgr_legend(fp,eNR,enms);
420   fprintf(fp,"@ s%d legend \"Ek/Nelec\"\n",eNR);
421   fprintf(fp,"@ type nxy\n");
422   for(i=0; (i<ae->nx); i++) {
423     fprintf(fp,"%10f",1000.0*dt*i);
424     for(j=0; (j<eNR); j++)
425       fprintf(fp,"  %8.3f",ae->e[i][j]*fac);
426     fprintf(fp,"  %8.3f\n",ae->e[i][eELECTRON]/(ELECTRONVOLT*total->nion[i]));
427   }    
428   fprintf(fp,"&\n");
429   ffclose(fp);
430 }
431