Removed the unnecessary hacks from some boolean statements.
[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   bPres         = (flags & CGLO_PRESSURE); 
172   bConstrVir    = (flags & CGLO_CONSTRAINT);
173   bFirstIterate = (flags & CGLO_FIRSTITERATE);
174   bEkinAveVel   = (inputrec->eI==eiVV || (inputrec->eI==eiVVAK && bPres));
175   bReadEkin     = (flags & CGLO_READEKIN);
176
177   rb   = gs->rb;
178   itc0 = gs->itc0;
179   itc1 = gs->itc1;
180   
181
182   reset_bin(rb);
183   /* This routine copies all the data to be summed to one big buffer
184    * using the t_bin struct. 
185    */
186
187   /* First, we neeed to identify which enerd->term should be
188      communicated.  Temperature and pressure terms should only be
189      communicated and summed when they need to be, to avoid repeating
190      the sums and overcounting. */
191
192   nener = filter_enerdterm(enerd->term,TRUE,copyenerd,bTemp,bPres,bEner);
193   
194   /* First, the data that needs to be communicated with velocity verlet every time
195      This is just the constraint virial.*/
196   if (bConstrVir) {
197       isv = add_binr(rb,DIM*DIM,svir[0]);
198       where();
199   }
200   
201 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
202   if (bTemp || !bVV)
203   {
204       if (ekind) 
205       {
206           for(j=0; (j<inputrec->opts.ngtc); j++) 
207           {
208               if (bSumEkinhOld) 
209               {
210                   itc0[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinh_old[0]);
211               }
212               if (bEkinAveVel && !bReadEkin) 
213               {
214                   itc1[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinf[0]);
215               } 
216               else if (!bReadEkin)
217               {
218                   itc1[j]=add_binr(rb,DIM*DIM,ekind->tcstat[j].ekinh[0]);
219               }
220           }
221           /* these probably need to be put into one of these categories */
222           where();
223           idedl = add_binr(rb,1,&(ekind->dekindl));
224           where();
225           ica   = add_binr(rb,1,&(ekind->cosacc.mvcos));
226           where();
227       }  
228   }      
229   where();
230   
231   if ((bPres || !bVV) && bFirstIterate)
232   {
233       ifv = add_binr(rb,DIM*DIM,fvir[0]);
234   }
235
236
237   if (bEner) 
238   { 
239       where();
240       if (bFirstIterate) 
241       {
242           ie  = add_binr(rb,nener,copyenerd);
243       }
244       where();
245       if (constr) 
246       {
247           rmsd_data = constr_rmsd_data(constr);
248           if (rmsd_data) 
249           {
250               irmsd = add_binr(rb,inputrec->eI==eiSD2 ? 3 : 2,rmsd_data);
251           }
252       } 
253       if (!NEED_MUTOT(*inputrec)) 
254       {
255           imu = add_binr(rb,DIM,mu_tot);
256           where();
257       }
258       
259       if (bFirstIterate) 
260       {
261           for(j=0; (j<egNR); j++)
262           {
263               inn[j]=add_binr(rb,enerd->grpp.nener,enerd->grpp.ener[j]);
264           }
265           where();
266           if (inputrec->efep != efepNO) 
267           {
268               idvdll  = add_bind(rb,1,&enerd->dvdl_lin);
269               idvdlnl = add_bind(rb,1,&enerd->dvdl_nonlin);
270               if (enerd->n_lambda > 0) 
271               {
272                   iepl = add_bind(rb,enerd->n_lambda,enerd->enerpart_lambda);
273               }
274           }
275       }
276
277       if (vcm) 
278       {
279           icm   = add_binr(rb,DIM*vcm->nr,vcm->group_p[0]);
280           where();
281           imass = add_binr(rb,vcm->nr,vcm->group_mass);
282           where();
283           if (vcm->mode == ecmANGULAR) 
284           {
285               icj   = add_binr(rb,DIM*vcm->nr,vcm->group_j[0]);
286               where();
287               icx   = add_binr(rb,DIM*vcm->nr,vcm->group_x[0]);
288               where();
289               ici   = add_binr(rb,DIM*DIM*vcm->nr,vcm->group_i[0][0]);
290               where();
291           }
292       }
293   }
294   if (DOMAINDECOMP(cr)) 
295   {
296       nb = cr->dd->nbonded_local;
297       inb = add_bind(rb,1,&nb);
298       }
299   where();
300   if (nsig > 0) 
301   {
302       isig = add_binr(rb,nsig,sig);
303   }
304
305   /* Global sum it all */
306   if (debug)
307   {
308       fprintf(debug,"Summing %d energies\n",rb->maxreal);
309   }
310   sum_bin(rb,cr);
311   where();
312
313   /* Extract all the data locally */
314
315   if (bConstrVir) 
316   {
317       extract_binr(rb,isv ,DIM*DIM,svir[0]);
318   }
319
320   /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
321   if (bTemp || !bVV)
322   {
323       if (ekind) 
324       {
325           for(j=0; (j<inputrec->opts.ngtc); j++) 
326           {
327               if (bSumEkinhOld)
328               {
329                   extract_binr(rb,itc0[j],DIM*DIM,ekind->tcstat[j].ekinh_old[0]);
330               }
331               if (bEkinAveVel && !bReadEkin) {
332                   extract_binr(rb,itc1[j],DIM*DIM,ekind->tcstat[j].ekinf[0]);
333               }
334               else if (!bReadEkin)
335               {
336                   extract_binr(rb,itc1[j],DIM*DIM,ekind->tcstat[j].ekinh[0]);              
337               }
338           }
339           extract_binr(rb,idedl,1,&(ekind->dekindl));
340           extract_binr(rb,ica,1,&(ekind->cosacc.mvcos));
341           where();
342       }
343   }
344   if ((bPres || !bVV) && bFirstIterate)
345   {
346       extract_binr(rb,ifv ,DIM*DIM,fvir[0]);
347   }
348
349   if (bEner) 
350   {
351       if (bFirstIterate) 
352       {
353           extract_binr(rb,ie,nener,copyenerd);
354           if (rmsd_data) 
355           {
356               extract_binr(rb,irmsd,inputrec->eI==eiSD2 ? 3 : 2,rmsd_data);
357           }
358           if (!NEED_MUTOT(*inputrec))
359           {
360               extract_binr(rb,imu,DIM,mu_tot);
361           }
362
363           for(j=0; (j<egNR); j++)
364           {
365               extract_binr(rb,inn[j],enerd->grpp.nener,enerd->grpp.ener[j]);
366           }
367           if (inputrec->efep != efepNO) 
368           {
369               extract_bind(rb,idvdll ,1,&enerd->dvdl_lin);
370               extract_bind(rb,idvdlnl,1,&enerd->dvdl_nonlin);
371               if (enerd->n_lambda > 0) 
372               {
373                   extract_bind(rb,iepl,enerd->n_lambda,enerd->enerpart_lambda);
374               }
375           }
376           /* should this be here, or with ekin?*/
377           if (vcm) 
378           {
379               extract_binr(rb,icm,DIM*vcm->nr,vcm->group_p[0]);
380               where();
381               extract_binr(rb,imass,vcm->nr,vcm->group_mass);
382               where();
383               if (vcm->mode == ecmANGULAR) 
384               {
385                   extract_binr(rb,icj,DIM*vcm->nr,vcm->group_j[0]);
386                   where();
387                   extract_binr(rb,icx,DIM*vcm->nr,vcm->group_x[0]);
388                   where();
389                   extract_binr(rb,ici,DIM*DIM*vcm->nr,vcm->group_i[0][0]);
390                   where();
391               }
392           }
393           if (DOMAINDECOMP(cr)) 
394           {
395               extract_bind(rb,inb,1,&nb);
396               if ((int)(nb + 0.5) != cr->dd->nbonded_global) 
397               {
398                   dd_print_missing_interactions(fplog,cr,(int)(nb + 0.5),top_global,state_local);
399               }
400           }
401           where();
402
403           filter_enerdterm(copyenerd,FALSE,enerd->term,bTemp,bPres,bEner);    
404 /* Small hack for temp only - not entirely clear if still needed?*/
405           /* enerd->term[F_TEMP] /= (cr->nnodes - cr->npmenodes); */
406       }
407   }
408
409   if (nsig > 0) 
410   {
411       extract_binr(rb,isig,nsig,sig);
412   }
413   where();
414 }
415
416 int do_per_step(gmx_large_int_t step,gmx_large_int_t nstep)
417 {
418   if (nstep != 0) 
419     return ((step % nstep)==0); 
420   else 
421     return 0;
422 }
423
424 static void moveit(t_commrec *cr,
425                    int left,int right,const char *s,rvec xx[])
426 {
427   if (!xx) 
428     return;
429
430   move_rvecs(cr,FALSE,FALSE,left,right,
431              xx,NULL,(cr->nnodes-cr->npmenodes)-1,NULL);
432 }
433
434 gmx_mdoutf_t *init_mdoutf(int nfile,const t_filenm fnm[],int mdrun_flags,
435                           const t_commrec *cr,const t_inputrec *ir,
436                           const output_env_t oenv)
437 {
438     gmx_mdoutf_t *of;
439     char filemode[3];
440     gmx_bool bAppendFiles;
441
442     snew(of,1);
443
444     of->fp_trn   = NULL;
445     of->fp_ene   = NULL;
446     of->fp_xtc   = NULL;
447     of->fp_dhdl  = NULL;
448     of->fp_field = NULL;
449     
450     of->eIntegrator     = ir->eI;
451     of->simulation_part = ir->simulation_part;
452
453     if (MASTER(cr))
454     {
455         bAppendFiles = (mdrun_flags & MD_APPENDFILES);
456
457         of->bKeepAndNumCPT = (mdrun_flags & MD_KEEPANDNUMCPT);
458
459         sprintf(filemode, bAppendFiles ? "a+" : "w+");  
460         
461         if (ir->eI != eiNM 
462 #ifndef GMX_FAHCORE
463             &&
464             !(EI_DYNAMICS(ir->eI) &&
465               ir->nstxout == 0 &&
466               ir->nstvout == 0 &&
467               ir->nstfout == 0)
468 #endif
469             )
470         {
471             of->fp_trn = open_trn(ftp2fn(efTRN,nfile,fnm), filemode);
472         }
473         if (!EI_ENERGY_MINIMIZATION(ir->eI) &&
474             ir->nstxtcout > 0)
475         {
476             of->fp_xtc = open_xtc(ftp2fn(efXTC,nfile,fnm), filemode);
477             of->xtc_prec = ir->xtcprec;
478         }
479         of->fp_ene = open_enx(ftp2fn(efEDR,nfile,fnm), filemode);
480         of->fn_cpt = opt2fn("-cpo",nfile,fnm);
481         
482         if (ir->efep != efepNO && ir->nstdhdl > 0 &&
483             (ir->separate_dhdl_file == sepdhdlfileYES ) && 
484             !EI_ENERGY_MINIMIZATION(ir->eI))
485         {
486             if (bAppendFiles)
487             {
488                 of->fp_dhdl = gmx_fio_fopen(opt2fn("-dhdl",nfile,fnm),filemode);
489             }
490             else
491             {
492                 of->fp_dhdl = open_dhdl(opt2fn("-dhdl",nfile,fnm),ir,oenv);
493             }
494         }
495         
496         if (opt2bSet("-field",nfile,fnm) &&
497             (ir->ex[XX].n || ir->ex[YY].n || ir->ex[ZZ].n))
498         {
499             if (bAppendFiles)
500             {
501                 of->fp_dhdl = gmx_fio_fopen(opt2fn("-field",nfile,fnm),
502                                             filemode);
503             }
504             else
505             {                             
506                 of->fp_field = xvgropen(opt2fn("-field",nfile,fnm),
507                                         "Applied electric field","Time (ps)",
508                                         "E (V/nm)",oenv);
509             }
510         }
511     }
512
513     return of;
514 }
515
516 void done_mdoutf(gmx_mdoutf_t *of)
517 {
518     if (of->fp_ene != NULL)
519     {
520         close_enx(of->fp_ene);
521     }
522     if (of->fp_xtc)
523     {
524         close_xtc(of->fp_xtc);
525     }
526     if (of->fp_trn)
527     {
528         close_trn(of->fp_trn);
529     }
530     if (of->fp_dhdl != NULL)
531     {
532         gmx_fio_fclose(of->fp_dhdl);
533     }
534     if (of->fp_field != NULL)
535     {
536         gmx_fio_fclose(of->fp_field);
537     }
538
539     sfree(of);
540 }
541
542 void write_traj(FILE *fplog,t_commrec *cr,
543                 gmx_mdoutf_t *of,
544                 int mdof_flags,
545                 gmx_mtop_t *top_global,
546                 gmx_large_int_t step,double t,
547                 t_state *state_local,t_state *state_global,
548                 rvec *f_local,rvec *f_global,
549                 int *n_xtc,rvec **x_xtc)
550 {
551     int     i,j;
552     gmx_groups_t *groups;
553     rvec    *xxtc;
554     rvec *local_v;
555     rvec *global_v;
556     
557 #define MX(xvf) moveit(cr,GMX_LEFT,GMX_RIGHT,#xvf,xvf)
558
559     /* MRS -- defining these variables is to manage the difference
560      * between half step and full step velocities, but there must be a better way . . . */
561
562     local_v  = state_local->v;
563     global_v = state_global->v;
564     
565     if (DOMAINDECOMP(cr))
566     {
567         if (mdof_flags & MDOF_CPT)
568         {
569             dd_collect_state(cr->dd,state_local,state_global);
570         }
571         else
572         {
573             if (mdof_flags & (MDOF_X | MDOF_XTC))
574             {
575                 dd_collect_vec(cr->dd,state_local,state_local->x,
576                                state_global->x);
577             }
578             if (mdof_flags & MDOF_V)
579             {
580                 dd_collect_vec(cr->dd,state_local,local_v,
581                                global_v);
582             }
583         }
584         if (mdof_flags & MDOF_F)
585         {
586             dd_collect_vec(cr->dd,state_local,f_local,f_global);
587         }
588     }
589     else
590     {
591         if (mdof_flags & MDOF_CPT)
592         {
593             /* All pointers in state_local are equal to state_global,
594              * but we need to copy the non-pointer entries.
595              */
596             state_global->lambda = state_local->lambda;
597             state_global->veta = state_local->veta;
598             state_global->vol0 = state_local->vol0;
599             copy_mat(state_local->box,state_global->box);
600             copy_mat(state_local->boxv,state_global->boxv);
601             copy_mat(state_local->svir_prev,state_global->svir_prev);
602             copy_mat(state_local->fvir_prev,state_global->fvir_prev);
603             copy_mat(state_local->pres_prev,state_global->pres_prev);
604         }
605         if (cr->nnodes > 1)
606         {
607             /* Particle decomposition, collect the data on the master node */
608             if (mdof_flags & MDOF_CPT)
609             {
610                 if (state_local->flags & (1<<estX))   MX(state_global->x);
611                 if (state_local->flags & (1<<estV))   MX(state_global->v);
612                 if (state_local->flags & (1<<estSDX)) MX(state_global->sd_X);
613                 if (state_global->nrngi > 1) {
614                     if (state_local->flags & (1<<estLD_RNG)) {
615 #ifdef GMX_MPI
616                         MPI_Gather(state_local->ld_rng ,
617                                    state_local->nrng*sizeof(state_local->ld_rng[0]),MPI_BYTE,
618                                    state_global->ld_rng,
619                                    state_local->nrng*sizeof(state_local->ld_rng[0]),MPI_BYTE,
620                                    MASTERRANK(cr),cr->mpi_comm_mygroup);
621 #endif
622                     }
623                     if (state_local->flags & (1<<estLD_RNGI))
624                     {
625 #ifdef GMX_MPI
626                         MPI_Gather(state_local->ld_rngi,
627                                    sizeof(state_local->ld_rngi[0]),MPI_BYTE,
628                                    state_global->ld_rngi,
629                                    sizeof(state_local->ld_rngi[0]),MPI_BYTE,
630                                    MASTERRANK(cr),cr->mpi_comm_mygroup);
631 #endif
632                     }
633                 }
634             }
635             else
636             {
637                 if (mdof_flags & (MDOF_X | MDOF_XTC)) MX(state_global->x);
638                 if (mdof_flags & MDOF_V)              MX(global_v);
639             }
640             if (mdof_flags & MDOF_F) MX(f_global);
641          }
642      }
643
644      if (MASTER(cr))
645      {
646          if (mdof_flags & MDOF_CPT)
647          {
648              write_checkpoint(of->fn_cpt,of->bKeepAndNumCPT,
649                               fplog,cr,of->eIntegrator,
650                               of->simulation_part,step,t,state_global);
651          }
652
653          if (mdof_flags & (MDOF_X | MDOF_V | MDOF_F))
654          {
655             fwrite_trn(of->fp_trn,step,t,state_local->lambda,
656                        state_local->box,top_global->natoms,
657                        (mdof_flags & MDOF_X) ? state_global->x : NULL,
658                        (mdof_flags & MDOF_V) ? global_v : NULL,
659                        (mdof_flags & MDOF_F) ? f_global : NULL);
660             if (gmx_fio_flush(of->fp_trn) != 0)
661             {
662                 gmx_file("Cannot write trajectory; maybe you are out of quota?");
663             }
664             gmx_fio_check_file_position(of->fp_trn);
665         }      
666         if (mdof_flags & MDOF_XTC) {
667             groups = &top_global->groups;
668             if (*n_xtc == -1)
669             {
670                 *n_xtc = 0;
671                 for(i=0; (i<top_global->natoms); i++)
672                 {
673                     if (ggrpnr(groups,egcXTC,i) == 0)
674                     {
675                         (*n_xtc)++;
676                     }
677                 }
678                 if (*n_xtc != top_global->natoms)
679                 {
680                     snew(*x_xtc,*n_xtc);
681                 }
682             }
683             if (*n_xtc == top_global->natoms)
684             {
685                 xxtc = state_global->x;
686             }
687             else
688             {
689                 xxtc = *x_xtc;
690                 j = 0;
691                 for(i=0; (i<top_global->natoms); i++)
692                 {
693                     if (ggrpnr(groups,egcXTC,i) == 0)
694                     {
695                         copy_rvec(state_global->x[i],xxtc[j++]);
696                     }
697                 }
698             }
699             if (write_xtc(of->fp_xtc,*n_xtc,step,t,
700                           state_local->box,xxtc,of->xtc_prec) == 0)
701             {
702                 gmx_fatal(FARGS,"XTC error - maybe you are out of quota?");
703             }
704             gmx_fio_check_file_position(of->fp_xtc);
705         }
706     }
707 }
708