Redefine the default boolean type to gmx_bool.
[alexxy/gromacs.git] / src / mdlib / stat.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  *                This source code is part of
5  * 
6  *                 G   R   O   M   A   C   S
7  * 
8  *          GROningen MAchine for Chemical Simulations
9  * 
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  * 
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  * 
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * GROningen Mixture of Alchemy and Childrens' Stories
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <string.h>
41 #include <stdio.h>
42 #include "typedefs.h"
43 #include "sysstuff.h"
44 #include "gmx_fatal.h"
45 #include "network.h"
46 #include "txtdump.h"
47 #include "names.h"
48 #include "physics.h"
49 #include "vec.h"
50 #include "maths.h"
51 #include "mvdata.h"
52 #include "main.h"
53 #include "force.h"
54 #include "vcm.h"
55 #include "smalloc.h"
56 #include "futil.h"
57 #include "network.h"
58 #include "rbin.h"
59 #include "tgroup.h"
60 #include "xtcio.h"
61 #include "gmxfio.h"
62 #include "trnio.h"
63 #include "statutil.h"
64 #include "domdec.h"
65 #include "partdec.h"
66 #include "constr.h"
67 #include "checkpoint.h"
68 #include "mdrun.h"
69 #include "xvgr.h"
70
71 typedef struct gmx_global_stat
72 {
73     t_bin *rb;
74     int   *itc0;
75     int   *itc1;
76 } t_gmx_global_stat;
77
78 gmx_global_stat_t global_stat_init(t_inputrec *ir)
79 {
80     gmx_global_stat_t gs;
81
82     snew(gs,1);
83     
84     gs->rb = mk_bin();
85     snew(gs->itc0,ir->opts.ngtc);
86     snew(gs->itc1,ir->opts.ngtc);
87
88     return gs;
89 }
90
91 void global_stat_destroy(gmx_global_stat_t gs)
92 {
93     destroy_bin(gs->rb);
94     sfree(gs->itc0);
95     sfree(gs->itc1);
96     sfree(gs);
97 }
98
99 static int filter_enerdterm(real *afrom, gmx_bool bToBuffer, real *ato,
100                             gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner) {
101     int i,to,from;
102
103     from = 0;
104     to   = 0;
105     for (i=0;i<F_NRE;i++)
106     {
107         if (bToBuffer)
108         {
109             from = i;
110         }
111         else
112         {
113             to = i;
114         }
115         switch (i) {
116         case F_EKIN:
117         case F_TEMP:
118         case F_DKDL:
119             if (bTemp)
120             {
121                 ato[to++] = afrom[from++];
122             }
123             break;
124         case F_PRES:    
125         case F_PDISPCORR:
126         case F_VTEMP:
127             if (bPres)
128             {
129                 ato[to++] = afrom[from++];
130             }
131             break;
132         default:
133             if (bEner)
134             {
135                 ato[to++] = afrom[from++];
136             }
137             break;
138         }
139     }
140
141     return to;
142 }
143
144 void global_stat(FILE *fplog,gmx_global_stat_t gs,
145                  t_commrec *cr,gmx_enerdata_t *enerd,
146                  tensor fvir,tensor svir,rvec mu_tot,
147                  t_inputrec *inputrec,
148                  gmx_ekindata_t *ekind,gmx_constr_t constr,
149                  t_vcm *vcm,
150                  int nsig,real *sig,
151                  gmx_mtop_t *top_global, t_state *state_local, 
152                  gmx_bool bSumEkinhOld, int flags)
153 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
154 {
155   t_bin  *rb;
156   int    *itc0,*itc1;
157   int    ie=0,ifv=0,isv=0,irmsd=0,imu=0;
158   int    idedl=0,idvdll=0,idvdlnl=0,iepl=0,icm=0,imass=0,ica=0,inb=0;
159   int    isig=-1;
160   int    icj=-1,ici=-1,icx=-1;
161   int    inn[egNR];
162   real   copyenerd[F_NRE];
163   int    nener,j;
164   real   *rmsd_data=NULL;
165   double nb;
166   gmx_bool   bVV,bTemp,bEner,bPres,bConstrVir,bEkinAveVel,bFirstIterate,bReadEkin;
167
168   bVV           = EI_VV(inputrec->eI);
169   bTemp         = flags & CGLO_TEMPERATURE;
170   bEner         = flags & CGLO_ENERGY;
171 /* FIX ME after 4.5 */
172 /* temporary hack because we are using gmx_bool (unsigned char) */
173   bPres         = (flags & CGLO_PRESSURE) != 0; 
174   bConstrVir    = (flags & CGLO_CONSTRAINT) != 0;
175   bFirstIterate = (flags & CGLO_FIRSTITERATE) != 0;
176   bEkinAveVel   = (inputrec->eI==eiVV || (inputrec->eI==eiVVAK && bPres));
177   bReadEkin     = (flags & CGLO_READEKIN) != 0;
178
179   rb   = gs->rb;
180   itc0 = gs->itc0;
181   itc1 = gs->itc1;
182   
183
184   reset_bin(rb);
185   /* This routine copies all the data to be summed to one big buffer
186    * using the t_bin struct. 
187    */
188
189   /* First, we neeed to identify which enerd->term should be
190      communicated.  Temperature and pressure terms should only be
191      communicated and summed when they need to be, to avoid repeating
192      the sums and overcounting. */
193
194   nener = filter_enerdterm(enerd->term,TRUE,copyenerd,bTemp,bPres,bEner);
195   
196   /* First, the data that needs to be communicated with velocity verlet every time
197      This is just the constraint virial.*/
198   if (bConstrVir) {
199       isv = add_binr(rb,DIM*DIM,svir[0]);
200       where();
201   }
202   
203 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
204   if (bTemp || !bVV)
205   {
206       if (ekind) 
207       {
208           for(j=0; (j<inputrec->opts.ngtc); j++) 
209           {
210               if (bSumEkinhOld) 
211               {
212                   itc0[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinh_old[0]);
213               }
214               if (bEkinAveVel && !bReadEkin) 
215               {
216                   itc1[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinf[0]);
217               } 
218               else if (!bReadEkin)
219               {
220                   itc1[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinh[0]);
221               }
222           }
223           /* these probably need to be put into one of these categories */
224           where();
225           idedl = add_binr(rb,1,&(ekind->dekindl));
226           where();
227           ica   = add_binr(rb,1,&(ekind->cosacc.mvcos));
228           where();
229       }  
230   }      
231   where();
232   
233   if ((bPres || !bVV) && bFirstIterate)
234   {
235       ifv = add_binr(rb,DIM*DIM,fvir[0]);
236   }
237
238
239   if (bEner) 
240   { 
241       where();
242       if (bFirstIterate) 
243       {
244           ie  = add_binr(rb,nener,copyenerd);
245       }
246       where();
247       if (constr) 
248       {
249           rmsd_data = constr_rmsd_data(constr);
250           if (rmsd_data) 
251           {
252               irmsd = add_binr(rb,inputrec->eI==eiSD2 ? 3 : 2,rmsd_data);
253           }
254       } 
255       if (!NEED_MUTOT(*inputrec)) 
256       {
257           imu = add_binr(rb,DIM,mu_tot);
258           where();
259       }
260       
261       if (bFirstIterate) 
262       {
263           for(j=0; (j<egNR); j++)
264           {
265               inn[j]=add_binr(rb,enerd->grpp.nener,enerd->grpp.ener[j]);
266           }
267           where();
268           if (inputrec->efep != efepNO) 
269           {
270               idvdll  = add_bind(rb,1,&enerd->dvdl_lin);
271               idvdlnl = add_bind(rb,1,&enerd->dvdl_nonlin);
272               if (enerd->n_lambda > 0) 
273               {
274                   iepl = add_bind(rb,enerd->n_lambda,enerd->enerpart_lambda);
275               }
276           }
277       }
278
279       if (vcm) 
280       {
281           icm   = add_binr(rb,DIM*vcm->nr,vcm->group_p[0]);
282           where();
283           imass = add_binr(rb,vcm->nr,vcm->group_mass);
284           where();
285           if (vcm->mode == ecmANGULAR) 
286           {
287               icj   = add_binr(rb,DIM*vcm->nr,vcm->group_j[0]);
288               where();
289               icx   = add_binr(rb,DIM*vcm->nr,vcm->group_x[0]);
290               where();
291               ici   = add_binr(rb,DIM*DIM*vcm->nr,vcm->group_i[0][0]);
292               where();
293           }
294       }
295   }
296   if (DOMAINDECOMP(cr)) 
297   {
298       nb = cr->dd->nbonded_local;
299       inb = add_bind(rb,1,&nb);
300       }
301   where();
302   if (nsig > 0) 
303   {
304       isig = add_binr(rb,nsig,sig);
305   }
306
307   /* Global sum it all */
308   if (debug)
309   {
310       fprintf(debug,"Summing %d energies\n",rb->maxreal);
311   }
312   sum_bin(rb,cr);
313   where();
314
315   /* Extract all the data locally */
316
317   if (bConstrVir) 
318   {
319       extract_binr(rb,isv ,DIM*DIM,svir[0]);
320   }
321
322   /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
323   if (bTemp || !bVV)
324   {
325       if (ekind) 
326       {
327           for(j=0; (j<inputrec->opts.ngtc); j++) 
328           {
329               if (bSumEkinhOld)
330               {
331                   extract_binr(rb,itc0[j],DIM*DIM,ekind->tcstat[j].ekinh_old[0]);
332               }
333               if (bEkinAveVel && !bReadEkin) {
334                   extract_binr(rb,itc1[j],DIM*DIM,ekind->tcstat[j].ekinf[0]);
335               }
336               else if (!bReadEkin)
337               {
338                   extract_binr(rb,itc1[j],DIM*DIM,ekind->tcstat[j].ekinh[0]);              
339               }
340           }
341           extract_binr(rb,idedl,1,&(ekind->dekindl));
342           extract_binr(rb,ica,1,&(ekind->cosacc.mvcos));
343           where();
344       }
345   }
346   if ((bPres || !bVV) && bFirstIterate)
347   {
348       extract_binr(rb,ifv ,DIM*DIM,fvir[0]);
349   }
350
351   if (bEner) 
352   {
353       if (bFirstIterate) 
354       {
355           extract_binr(rb,ie,nener,copyenerd);
356           if (rmsd_data) 
357           {
358               extract_binr(rb,irmsd,inputrec->eI==eiSD2 ? 3 : 2,rmsd_data);
359           }
360           if (!NEED_MUTOT(*inputrec))
361           {
362               extract_binr(rb,imu,DIM,mu_tot);
363           }
364
365           for(j=0; (j<egNR); j++)
366           {
367               extract_binr(rb,inn[j],enerd->grpp.nener,enerd->grpp.ener[j]);
368           }
369           if (inputrec->efep != efepNO) 
370           {
371               extract_bind(rb,idvdll ,1,&enerd->dvdl_lin);
372               extract_bind(rb,idvdlnl,1,&enerd->dvdl_nonlin);
373               if (enerd->n_lambda > 0) 
374               {
375                   extract_bind(rb,iepl,enerd->n_lambda,enerd->enerpart_lambda);
376               }
377           }
378           /* should this be here, or with ekin?*/
379           if (vcm) 
380           {
381               extract_binr(rb,icm,DIM*vcm->nr,vcm->group_p[0]);
382               where();
383               extract_binr(rb,imass,vcm->nr,vcm->group_mass);
384               where();
385               if (vcm->mode == ecmANGULAR) 
386               {
387                   extract_binr(rb,icj,DIM*vcm->nr,vcm->group_j[0]);
388                   where();
389                   extract_binr(rb,icx,DIM*vcm->nr,vcm->group_x[0]);
390                   where();
391                   extract_binr(rb,ici,DIM*DIM*vcm->nr,vcm->group_i[0][0]);
392                   where();
393               }
394           }
395           if (DOMAINDECOMP(cr)) 
396           {
397               extract_bind(rb,inb,1,&nb);
398               if ((int)(nb + 0.5) != cr->dd->nbonded_global) 
399               {
400                   dd_print_missing_interactions(fplog,cr,(int)(nb + 0.5),top_global,state_local);
401               }
402           }
403           where();
404
405           filter_enerdterm(copyenerd,FALSE,enerd->term,bTemp,bPres,bEner);    
406 /* Small hack for temp only - not entirely clear if still needed?*/
407           /* enerd->term[F_TEMP] /= (cr->nnodes - cr->npmenodes); */
408       }
409   }
410
411   if (nsig > 0) 
412   {
413       extract_binr(rb,isig,nsig,sig);
414   }
415   where();
416 }
417
418 int do_per_step(gmx_large_int_t step,gmx_large_int_t nstep)
419 {
420   if (nstep != 0) 
421     return ((step % nstep)==0); 
422   else 
423     return 0;
424 }
425
426 static void moveit(t_commrec *cr,
427                    int left,int right,const char *s,rvec xx[])
428 {
429   if (!xx) 
430     return;
431
432   move_rvecs(cr,FALSE,FALSE,left,right,
433              xx,NULL,(cr->nnodes-cr->npmenodes)-1,NULL);
434 }
435
436 gmx_mdoutf_t *init_mdoutf(int nfile,const t_filenm fnm[],int mdrun_flags,
437                           const t_commrec *cr,const t_inputrec *ir,
438                           const output_env_t oenv)
439 {
440     gmx_mdoutf_t *of;
441     char filemode[3];
442     gmx_bool bAppendFiles;
443
444     snew(of,1);
445
446     of->fp_trn   = NULL;
447     of->fp_ene   = NULL;
448     of->fp_xtc   = NULL;
449     of->fp_dhdl  = NULL;
450     of->fp_field = NULL;
451     
452     of->eIntegrator     = ir->eI;
453     of->simulation_part = ir->simulation_part;
454
455     if (MASTER(cr))
456     {
457         bAppendFiles = (mdrun_flags & MD_APPENDFILES) != 0;
458
459         of->bKeepAndNumCPT = (mdrun_flags & MD_KEEPANDNUMCPT) != 0;
460
461         sprintf(filemode, bAppendFiles ? "a+" : "w+");  
462         
463         if (ir->eI != eiNM 
464 #ifndef GMX_FAHCORE
465             &&
466             !(EI_DYNAMICS(ir->eI) &&
467               ir->nstxout == 0 &&
468               ir->nstvout == 0 &&
469               ir->nstfout == 0)
470 #endif
471             )
472         {
473             of->fp_trn = open_trn(ftp2fn(efTRN,nfile,fnm), filemode);
474         }
475         if (!EI_ENERGY_MINIMIZATION(ir->eI) &&
476             ir->nstxtcout > 0)
477         {
478             of->fp_xtc = open_xtc(ftp2fn(efXTC,nfile,fnm), filemode);
479             of->xtc_prec = ir->xtcprec;
480         }
481         of->fp_ene = open_enx(ftp2fn(efEDR,nfile,fnm), filemode);
482         of->fn_cpt = opt2fn("-cpo",nfile,fnm);
483         
484         if (ir->efep != efepNO && ir->nstdhdl > 0 &&
485             (ir->separate_dhdl_file == sepdhdlfileYES ) && 
486             !EI_ENERGY_MINIMIZATION(ir->eI))
487         {
488             if (bAppendFiles)
489             {
490                 of->fp_dhdl = gmx_fio_fopen(opt2fn("-dhdl",nfile,fnm),filemode);
491             }
492             else
493             {
494                 of->fp_dhdl = open_dhdl(opt2fn("-dhdl",nfile,fnm),ir,oenv);
495             }
496         }
497         
498         if (opt2bSet("-field",nfile,fnm) &&
499             (ir->ex[XX].n || ir->ex[YY].n || ir->ex[ZZ].n))
500         {
501             if (bAppendFiles)
502             {
503                 of->fp_dhdl = gmx_fio_fopen(opt2fn("-field",nfile,fnm),
504                                             filemode);
505             }
506             else
507             {                             
508                 of->fp_field = xvgropen(opt2fn("-field",nfile,fnm),
509                                         "Applied electric field","Time (ps)",
510                                         "E (V/nm)",oenv);
511             }
512         }
513     }
514
515     return of;
516 }
517
518 void done_mdoutf(gmx_mdoutf_t *of)
519 {
520     if (of->fp_ene != NULL)
521     {
522         close_enx(of->fp_ene);
523     }
524     if (of->fp_xtc)
525     {
526         close_xtc(of->fp_xtc);
527     }
528     if (of->fp_trn)
529     {
530         close_trn(of->fp_trn);
531     }
532     if (of->fp_dhdl != NULL)
533     {
534         gmx_fio_fclose(of->fp_dhdl);
535     }
536     if (of->fp_field != NULL)
537     {
538         gmx_fio_fclose(of->fp_field);
539     }
540
541     sfree(of);
542 }
543
544 void write_traj(FILE *fplog,t_commrec *cr,
545                 gmx_mdoutf_t *of,
546                 int mdof_flags,
547                 gmx_mtop_t *top_global,
548                 gmx_large_int_t step,double t,
549                 t_state *state_local,t_state *state_global,
550                 rvec *f_local,rvec *f_global,
551                 int *n_xtc,rvec **x_xtc)
552 {
553     int     i,j;
554     gmx_groups_t *groups;
555     rvec    *xxtc;
556     rvec *local_v;
557     rvec *global_v;
558     
559 #define MX(xvf) moveit(cr,GMX_LEFT,GMX_RIGHT,#xvf,xvf)
560
561     /* MRS -- defining these variables is to manage the difference
562      * between half step and full step velocities, but there must be a better way . . . */
563
564     local_v  = state_local->v;
565     global_v = state_global->v;
566     
567     if (DOMAINDECOMP(cr))
568     {
569         if (mdof_flags & MDOF_CPT)
570         {
571             dd_collect_state(cr->dd,state_local,state_global);
572         }
573         else
574         {
575             if (mdof_flags & (MDOF_X | MDOF_XTC))
576             {
577                 dd_collect_vec(cr->dd,state_local,state_local->x,
578                                state_global->x);
579             }
580             if (mdof_flags & MDOF_V)
581             {
582                 dd_collect_vec(cr->dd,state_local,local_v,
583                                global_v);
584             }
585         }
586         if (mdof_flags & MDOF_F)
587         {
588             dd_collect_vec(cr->dd,state_local,f_local,f_global);
589         }
590     }
591     else
592     {
593         if (mdof_flags & MDOF_CPT)
594         {
595             /* All pointers in state_local are equal to state_global,
596              * but we need to copy the non-pointer entries.
597              */
598             state_global->lambda = state_local->lambda;
599             state_global->veta = state_local->veta;
600             state_global->vol0 = state_local->vol0;
601             copy_mat(state_local->box,state_global->box);
602             copy_mat(state_local->boxv,state_global->boxv);
603             copy_mat(state_local->svir_prev,state_global->svir_prev);
604             copy_mat(state_local->fvir_prev,state_global->fvir_prev);
605             copy_mat(state_local->pres_prev,state_global->pres_prev);
606         }
607         if (cr->nnodes > 1)
608         {
609             /* Particle decomposition, collect the data on the master node */
610             if (mdof_flags & MDOF_CPT)
611             {
612                 if (state_local->flags & (1<<estX))   MX(state_global->x);
613                 if (state_local->flags & (1<<estV))   MX(state_global->v);
614                 if (state_local->flags & (1<<estSDX)) MX(state_global->sd_X);
615                 if (state_global->nrngi > 1) {
616                     if (state_local->flags & (1<<estLD_RNG)) {
617 #ifdef GMX_MPI
618                         MPI_Gather(state_local->ld_rng ,
619                                    state_local->nrng*sizeof(state_local->ld_rng[0]),MPI_BYTE,
620                                    state_global->ld_rng,
621                                    state_local->nrng*sizeof(state_local->ld_rng[0]),MPI_BYTE,
622                                    MASTERRANK(cr),cr->mpi_comm_mygroup);
623 #endif
624                     }
625                     if (state_local->flags & (1<<estLD_RNGI))
626                     {
627 #ifdef GMX_MPI
628                         MPI_Gather(state_local->ld_rngi,
629                                    sizeof(state_local->ld_rngi[0]),MPI_BYTE,
630                                    state_global->ld_rngi,
631                                    sizeof(state_local->ld_rngi[0]),MPI_BYTE,
632                                    MASTERRANK(cr),cr->mpi_comm_mygroup);
633 #endif
634                     }
635                 }
636             }
637             else
638             {
639                 if (mdof_flags & (MDOF_X | MDOF_XTC)) MX(state_global->x);
640                 if (mdof_flags & MDOF_V)              MX(global_v);
641             }
642             if (mdof_flags & MDOF_F) MX(f_global);
643          }
644      }
645
646      if (MASTER(cr))
647      {
648          if (mdof_flags & MDOF_CPT)
649          {
650              write_checkpoint(of->fn_cpt,of->bKeepAndNumCPT,
651                               fplog,cr,of->eIntegrator,
652                               of->simulation_part,step,t,state_global);
653          }
654
655          if (mdof_flags & (MDOF_X | MDOF_V | MDOF_F))
656          {
657             fwrite_trn(of->fp_trn,step,t,state_local->lambda,
658                        state_local->box,top_global->natoms,
659                        (mdof_flags & MDOF_X) ? state_global->x : NULL,
660                        (mdof_flags & MDOF_V) ? global_v : NULL,
661                        (mdof_flags & MDOF_F) ? f_global : NULL);
662             if (gmx_fio_flush(of->fp_trn) != 0)
663             {
664                 gmx_file("Cannot write trajectory; maybe you are out of quota?");
665             }
666             gmx_fio_check_file_position(of->fp_trn);
667         }      
668         if (mdof_flags & MDOF_XTC) {
669             groups = &top_global->groups;
670             if (*n_xtc == -1)
671             {
672                 *n_xtc = 0;
673                 for(i=0; (i<top_global->natoms); i++)
674                 {
675                     if (ggrpnr(groups,egcXTC,i) == 0)
676                     {
677                         (*n_xtc)++;
678                     }
679                 }
680                 if (*n_xtc != top_global->natoms)
681                 {
682                     snew(*x_xtc,*n_xtc);
683                 }
684             }
685             if (*n_xtc == top_global->natoms)
686             {
687                 xxtc = state_global->x;
688             }
689             else
690             {
691                 xxtc = *x_xtc;
692                 j = 0;
693                 for(i=0; (i<top_global->natoms); i++)
694                 {
695                     if (ggrpnr(groups,egcXTC,i) == 0)
696                     {
697                         copy_rvec(state_global->x[i],xxtc[j++]);
698                     }
699                 }
700             }
701             if (write_xtc(of->fp_xtc,*n_xtc,step,t,
702                           state_local->box,xxtc,of->xtc_prec) == 0)
703             {
704                 gmx_fatal(FARGS,"XTC error - maybe you are out of quota?");
705             }
706             gmx_fio_check_file_position(of->fp_xtc);
707         }
708     }
709 }
710