Merge 'origin/release-4-6' into master
[alexxy/gromacs.git] / src / gromacs / mdlib / force.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  * GROwing Monsters And Cloning Shrimps
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <math.h>
41 #include <string.h>
42 #include "sysstuff.h"
43 #include "typedefs.h"
44 #include "macros.h"
45 #include "smalloc.h"
46 #include "macros.h"
47 #include "physics.h"
48 #include "force.h"
49 #include "nonbonded.h"
50 #include "names.h"
51 #include "network.h"
52 #include "pbc.h"
53 #include "ns.h"
54 #include "nrnb.h"
55 #include "bondf.h"
56 #include "mshift.h"
57 #include "txtdump.h"
58 #include "coulomb.h"
59 #include "pme.h"
60 #include "mdrun.h"
61 #include "domdec.h"
62 #include "partdec.h"
63 #include "qmmm.h"
64
65
66 void ns(FILE *fp,
67         t_forcerec *fr,
68         rvec       x[],
69         matrix     box,
70         gmx_groups_t *groups,
71         t_grpopts  *opts,
72         gmx_localtop_t *top,
73         t_mdatoms  *md,
74         t_commrec  *cr,
75         t_nrnb     *nrnb,
76         real       lambda,
77         real       *dvdlambda,
78         gmx_grppairener_t *grppener,
79         gmx_bool       bFillGrid,
80         gmx_bool       bDoLongRange,
81         gmx_bool       bDoForces,
82         rvec       *f)
83 {
84   char   *ptr;
85   int    nsearch;
86
87
88   if (!fr->ns.nblist_initialized)
89   {
90       init_neighbor_list(fp, fr, md->homenr);
91   }
92     
93   if (fr->bTwinRange) 
94     fr->nlr=0;
95
96     nsearch = search_neighbours(fp,fr,x,box,top,groups,cr,nrnb,md,
97                                 lambda,dvdlambda,grppener,
98                                 bFillGrid,bDoLongRange,
99                                 bDoForces,f);
100   if (debug)
101     fprintf(debug,"nsearch = %d\n",nsearch);
102     
103   /* Check whether we have to do dynamic load balancing */
104   /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
105     count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
106     &(top->idef),opts->ngener);
107   */
108   if (fr->ns.dump_nl > 0)
109     dump_nblist(fp,cr,fr,fr->ns.dump_nl);
110 }
111
112 void do_force_lowlevel(FILE       *fplog,   gmx_large_int_t step,
113                        t_forcerec *fr,      t_inputrec *ir,
114                        t_idef     *idef,    t_commrec  *cr,
115                        t_nrnb     *nrnb,    gmx_wallcycle_t wcycle,
116                        t_mdatoms  *md,
117                        t_grpopts  *opts,
118                        rvec       x[],      history_t  *hist,
119                        rvec       f[],
120                        gmx_enerdata_t *enerd,
121                        t_fcdata   *fcd,
122                        gmx_mtop_t     *mtop,
123                        gmx_localtop_t *top,
124                        gmx_genborn_t *born,
125                        t_atomtypes *atype,
126                        gmx_bool       bBornRadii,
127                        matrix     box,
128                        real       lambda,  
129                        t_graph    *graph,
130                        t_blocka   *excl,    
131                        rvec       mu_tot[],
132                        int        flags,
133                        float      *cycles_pme)
134 {
135     int     i,status;
136     int     donb_flags;
137     gmx_bool    bDoEpot,bSepDVDL,bSB;
138     int     pme_flags;
139     matrix  boxs;
140     rvec    box_size;
141     real    dvdlambda,Vsr,Vlr,Vcorr=0,vdip,vcharge;
142     t_pbc   pbc;
143     real    dvdgb;
144     char    buf[22];
145     gmx_enerdata_t ed_lam;
146     double  lam_i;
147     real    dvdl_dum;
148
149 #ifdef GMX_MPI
150     double  t0=0.0,t1,t2,t3; /* time measurement for coarse load balancing */
151 #endif
152     
153 #define PRINT_SEPDVDL(s,v,dvdl) if (bSepDVDL) fprintf(fplog,sepdvdlformat,s,v,dvdl);
154     
155
156     set_pbc(&pbc,fr->ePBC,box);
157     
158     /* Reset box */
159     for(i=0; (i<DIM); i++)
160     {
161         box_size[i]=box[i][i];
162     }
163     
164     bSepDVDL=(fr->bSepDVDL && do_per_step(step,ir->nstlog));
165     debug_gmx();
166     
167     /* do QMMM first if requested */
168     if(fr->bQMMM)
169     {
170         enerd->term[F_EQM] = calculate_QMMM(cr,x,f,fr,md);
171     }
172     
173     if (bSepDVDL)
174     {
175         fprintf(fplog,"Step %s: non-bonded V and dVdl for node %d:\n",
176                 gmx_step_str(step,buf),cr->nodeid);
177     }
178     
179     /* Call the short range functions all in one go. */
180     
181     dvdlambda = 0;
182     
183 #ifdef GMX_MPI
184     /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
185 #define TAKETIME FALSE
186     if (TAKETIME)
187     {
188         MPI_Barrier(cr->mpi_comm_mygroup);
189         t0=MPI_Wtime();
190     }
191 #endif
192     
193     if (ir->nwall)
194     {
195         dvdlambda = do_walls(ir,fr,box,md,x,f,lambda,
196                              enerd->grpp.ener[egLJSR],nrnb);
197         PRINT_SEPDVDL("Walls",0.0,dvdlambda);
198         enerd->dvdl_lin += dvdlambda;
199     }
200                 
201         /* If doing GB, reset dvda and calculate the Born radii */
202         if (ir->implicit_solvent)
203         {
204                 /* wallcycle_start(wcycle,ewcGB); */
205                 
206                 for(i=0;i<born->nr;i++)
207                 {
208                         fr->dvda[i]=0;
209                 }
210                 
211                 if(bBornRadii)
212                 {
213                         calc_gb_rad(cr,fr,ir,top,atype,x,&(fr->gblist),born,md,nrnb);
214                 }
215                 
216                 /* wallcycle_stop(wcycle, ewcGB); */
217         }
218         
219     where();
220     donb_flags = 0;
221     if (flags & GMX_FORCE_FORCES)
222     {
223         donb_flags |= GMX_DONB_FORCES;
224     }
225     do_nonbonded(cr,fr,x,f,md,excl,
226                  fr->bBHAM ?
227                  enerd->grpp.ener[egBHAMSR] :
228                  enerd->grpp.ener[egLJSR],
229                  enerd->grpp.ener[egCOULSR],
230                                  enerd->grpp.ener[egGB],box_size,nrnb,
231                  lambda,&dvdlambda,-1,-1,donb_flags);
232     /* If we do foreign lambda and we have soft-core interactions
233      * we have to recalculate the (non-linear) energies contributions.
234      */
235     if (ir->n_flambda > 0 && (flags & GMX_FORCE_DHDL) && ir->sc_alpha != 0)
236     {
237         init_enerdata(mtop->groups.grps[egcENER].nr,ir->n_flambda,&ed_lam);
238         
239         for(i=0; i<enerd->n_lambda; i++)
240         {
241             lam_i = (i==0 ? lambda : ir->flambda[i-1]);
242             dvdl_dum = 0;
243             reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
244             do_nonbonded(cr,fr,x,f,md,excl,
245                          fr->bBHAM ?
246                          ed_lam.grpp.ener[egBHAMSR] :
247                          ed_lam.grpp.ener[egLJSR],
248                          ed_lam.grpp.ener[egCOULSR],
249                          enerd->grpp.ener[egGB], box_size,nrnb,
250                          lam_i,&dvdl_dum,-1,-1,
251                          GMX_DONB_FOREIGNLAMBDA);
252             sum_epot(&ir->opts,&ed_lam);
253             enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
254         }
255         destroy_enerdata(&ed_lam);
256     }
257     where();
258         
259         /* If we are doing GB, calculate bonded forces and apply corrections 
260          * to the solvation forces */
261         if (ir->implicit_solvent)  {
262                 calc_gb_forces(cr,md,born,top,atype,x,f,fr,idef,
263                        ir->gb_algorithm,ir->sa_algorithm,nrnb,bBornRadii,&pbc,graph,enerd);
264     }
265
266 #ifdef GMX_MPI
267     if (TAKETIME)
268     {
269         t1=MPI_Wtime();
270         fr->t_fnbf += t1-t0;
271     }
272 #endif
273     
274     if (ir->sc_alpha != 0)
275     {
276         enerd->dvdl_nonlin += dvdlambda;
277     }
278     else
279     {
280         enerd->dvdl_lin    += dvdlambda;
281     }
282     Vsr = 0;
283     if (bSepDVDL)
284     {
285         for(i=0; i<enerd->grpp.nener; i++)
286         {
287             Vsr +=
288                 (fr->bBHAM ?
289                  enerd->grpp.ener[egBHAMSR][i] :
290                  enerd->grpp.ener[egLJSR][i])
291                 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
292         }
293     }
294     PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",Vsr,dvdlambda);
295     debug_gmx();
296     
297     if (debug)
298     {
299         pr_rvecs(debug,0,"fshift after SR",fr->fshift,SHIFTS);
300     }
301     
302     /* Shift the coordinates. Must be done before bonded forces and PPPM, 
303      * but is also necessary for SHAKE and update, therefore it can NOT 
304      * go when no bonded forces have to be evaluated.
305      */
306     
307     /* Here sometimes we would not need to shift with NBFonly,
308      * but we do so anyhow for consistency of the returned coordinates.
309      */
310     if (graph)
311     {
312         shift_self(graph,box,x);
313         if (TRICLINIC(box))
314         {
315             inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
316         }
317         else
318         {
319             inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
320         }
321     }
322     /* Check whether we need to do bondeds or correct for exclusions */
323     if (fr->bMolPBC &&
324         ((flags & GMX_FORCE_BONDED)
325          || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
326     {
327         /* Since all atoms are in the rectangular or triclinic unit-cell,
328          * only single box vector shifts (2 in x) are required.
329          */
330         set_pbc_dd(&pbc,fr->ePBC,cr->dd,TRUE,box);
331     }
332     debug_gmx();
333     
334     if (flags & GMX_FORCE_BONDED)
335     {
336         calc_bonds(fplog,cr->ms,
337                    idef,x,hist,f,fr,&pbc,graph,enerd,nrnb,lambda,md,fcd,
338                    DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
339                    fr->bSepDVDL && do_per_step(step,ir->nstlog),step);
340         
341         /* Check if we have to determine energy differences
342          * at foreign lambda's.
343          */
344         if (ir->n_flambda > 0 && (flags & GMX_FORCE_DHDL) &&
345             idef->ilsort != ilsortNO_FE)
346         {
347             if (idef->ilsort != ilsortFE_SORTED)
348             {
349                 gmx_incons("The bonded interactions are not sorted for free energy");
350             }
351             init_enerdata(mtop->groups.grps[egcENER].nr,ir->n_flambda,&ed_lam);
352             
353             for(i=0; i<enerd->n_lambda; i++)
354             {
355                 lam_i = (i==0 ? lambda : ir->flambda[i-1]);
356                 dvdl_dum = 0;
357                 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
358                 calc_bonds_lambda(fplog,
359                                   idef,x,fr,&pbc,graph,&ed_lam,nrnb,lam_i,md,
360                                   fcd,
361                                   DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
362                 sum_epot(&ir->opts,&ed_lam);
363                 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
364             }
365             destroy_enerdata(&ed_lam);
366         }
367         debug_gmx();
368     }
369
370     where();
371
372     *cycles_pme = 0;
373     if (EEL_FULL(fr->eeltype))
374     {
375         bSB = (ir->nwall == 2);
376         if (bSB)
377         {
378             copy_mat(box,boxs);
379             svmul(ir->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
380             box_size[ZZ] *= ir->wall_ewald_zfac;
381         }
382         
383         clear_mat(fr->vir_el_recip);    
384         
385         if (fr->bEwald)
386         {
387             if (fr->n_tpi == 0)
388             {
389                 dvdlambda = 0;
390                 Vcorr = ewald_LRcorrection(fplog,md->start,md->start+md->homenr,
391                                            cr,fr,
392                                            md->chargeA,
393                                            md->nChargePerturbed ? md->chargeB : NULL,
394                                            excl,x,bSB ? boxs : box,mu_tot,
395                                            ir->ewald_geometry,
396                                            ir->epsilon_surface,
397                                            lambda,&dvdlambda,&vdip,&vcharge);
398                 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.",Vcorr,dvdlambda);
399                 enerd->dvdl_lin += dvdlambda;
400             }
401             else
402             {
403                 if (ir->ewald_geometry != eewg3D || ir->epsilon_surface != 0)
404                 {
405                     gmx_fatal(FARGS,"TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
406                 }
407                 /* The TPI molecule does not have exclusions with the rest
408                  * of the system and no intra-molecular PME grid contributions
409                  * will be calculated in gmx_pme_calc_energy.
410                  */
411                 Vcorr = 0;
412             }
413         }
414         
415         dvdlambda = 0;
416         status = 0;
417         switch (fr->eeltype)
418         {
419         case eelPME:
420         case eelPMESWITCH:
421         case eelPMEUSER:
422         case eelPMEUSERSWITCH:
423         case eelP3M_AD:
424             if (cr->duty & DUTY_PME)
425             {
426                 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
427                 {
428                     pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
429                     if (flags & GMX_FORCE_FORCES)
430                     {
431                         pme_flags |= GMX_PME_CALC_F;
432                     }
433                     if (flags & GMX_FORCE_VIRIAL)
434                     {
435                         pme_flags |= GMX_PME_CALC_ENER_VIR;
436                     }
437                     if (fr->n_tpi > 0)
438                     {
439                         /* We don't calculate f, but we do want the potential */
440                         pme_flags |= GMX_PME_CALC_POT;
441                     }
442                     wallcycle_start(wcycle,ewcPMEMESH);
443                     status = gmx_pme_do(fr->pmedata,
444                                         md->start,md->homenr - fr->n_tpi,
445                                         x,fr->f_novirsum,
446                                         md->chargeA,md->chargeB,
447                                         bSB ? boxs : box,cr,
448                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
449                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
450                                         nrnb,wcycle,
451                                         fr->vir_el_recip,fr->ewaldcoeff,
452                                         &Vlr,lambda,&dvdlambda,
453                                         pme_flags);
454                     *cycles_pme = wallcycle_stop(wcycle,ewcPMEMESH);
455
456                     /* We should try to do as little computation after
457                      * this as possible, because parallel PME synchronizes
458                      * the nodes, so we want all load imbalance of the rest
459                      * of the force calculation to be before the PME call.
460                      * DD load balancing is done on the whole time of
461                      * the force call (without PME).
462                      */
463                 }
464                 if (fr->n_tpi > 0)
465                 {
466                     /* Determine the PME grid energy of the test molecule
467                      * with the PME grid potential of the other charges.
468                      */
469                     gmx_pme_calc_energy(fr->pmedata,fr->n_tpi,
470                                         x + md->homenr - fr->n_tpi,
471                                         md->chargeA + md->homenr - fr->n_tpi,
472                                         &Vlr);
473                 }
474                 PRINT_SEPDVDL("PME mesh",Vlr,dvdlambda);
475             } 
476             else
477             {
478                 /* Energies and virial are obtained later from the PME nodes */
479                 /* but values have to be zeroed out here */
480                 Vlr=0.0;
481             }
482             break;
483         case eelEWALD:
484             Vlr = do_ewald(fplog,FALSE,ir,x,fr->f_novirsum,
485                            md->chargeA,md->chargeB,
486                            box_size,cr,md->homenr,
487                            fr->vir_el_recip,fr->ewaldcoeff,
488                            lambda,&dvdlambda,fr->ewald_table);
489             PRINT_SEPDVDL("Ewald long-range",Vlr,dvdlambda);
490             break;
491         default:
492             Vlr = 0;
493             gmx_fatal(FARGS,"No such electrostatics method implemented %s",
494                       eel_names[fr->eeltype]);
495         }
496         if (status != 0)
497         {
498             gmx_fatal(FARGS,"Error %d in long range electrostatics routine %s",
499                       status,EELTYPE(fr->eeltype));
500                 }
501         enerd->dvdl_lin += dvdlambda;
502         enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
503         if (debug)
504         {
505             fprintf(debug,"Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
506                     Vlr,Vcorr,enerd->term[F_COUL_RECIP]);
507             pr_rvecs(debug,0,"vir_el_recip after corr",fr->vir_el_recip,DIM);
508             pr_rvecs(debug,0,"fshift after LR Corrections",fr->fshift,SHIFTS);
509         }
510     }
511     else
512     {
513         if (EEL_RF(fr->eeltype))
514         {
515             dvdlambda = 0;
516             
517             if (fr->eeltype != eelRF_NEC)
518             {
519                 enerd->term[F_RF_EXCL] =
520                     RF_excl_correction(fplog,fr,graph,md,excl,x,f,
521                                        fr->fshift,&pbc,lambda,&dvdlambda);
522             }
523             
524             enerd->dvdl_lin += dvdlambda;
525             PRINT_SEPDVDL("RF exclusion correction",
526                           enerd->term[F_RF_EXCL],dvdlambda);
527         }
528     }
529     where();
530     debug_gmx();
531         
532     if (debug)
533     {
534         print_nrnb(debug,nrnb); 
535     }
536     debug_gmx();
537     
538 #ifdef GMX_MPI
539     if (TAKETIME)
540     {
541         t2=MPI_Wtime();
542         MPI_Barrier(cr->mpi_comm_mygroup);
543         t3=MPI_Wtime();
544         fr->t_wait += t3-t2;
545         if (fr->timesteps == 11)
546         {
547             fprintf(stderr,"* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n", 
548                     cr->nodeid, gmx_step_str(fr->timesteps,buf), 
549                     100*fr->t_wait/(fr->t_wait+fr->t_fnbf), 
550                     (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
551         }         
552         fr->timesteps++;
553     }
554 #endif
555     
556     if (debug)
557     {
558         pr_rvecs(debug,0,"fshift after bondeds",fr->fshift,SHIFTS);
559     }
560     
561 }
562
563 void init_enerdata(int ngener,int n_flambda,gmx_enerdata_t *enerd)
564 {
565     int i,n2;
566     
567     for(i=0; i<F_NRE; i++)
568     {
569         enerd->term[i] = 0;
570     }
571     
572     n2=ngener*ngener;
573     if (debug)
574     {
575         fprintf(debug,"Creating %d sized group matrix for energies\n",n2);
576     }
577     enerd->grpp.nener = n2;
578     for(i=0; (i<egNR); i++)
579     {
580         snew(enerd->grpp.ener[i],n2);
581     }
582
583     if (n_flambda)
584     {
585         enerd->n_lambda = 1 + n_flambda;
586         snew(enerd->enerpart_lambda,enerd->n_lambda);
587     }
588     else
589     {
590         enerd->n_lambda = 0;
591     }
592 }
593
594 void destroy_enerdata(gmx_enerdata_t *enerd)
595 {
596     int i;
597
598     for(i=0; (i<egNR); i++)
599     {
600         sfree(enerd->grpp.ener[i]);
601     }
602
603     if (enerd->n_lambda)
604     {
605         sfree(enerd->enerpart_lambda);
606     }
607 }
608
609 static real sum_v(int n,real v[])
610 {
611   real t;
612   int  i;
613   
614   t = 0.0;
615   for(i=0; (i<n); i++)
616     t = t + v[i];
617     
618   return t;
619 }
620
621 void sum_epot(t_grpopts *opts,gmx_enerdata_t *enerd)
622 {
623   gmx_grppairener_t *grpp;
624   real *epot;
625   int i;
626   
627   grpp = &enerd->grpp;
628   epot = enerd->term;
629
630   /* Accumulate energies */
631   epot[F_COUL_SR]  = sum_v(grpp->nener,grpp->ener[egCOULSR]);
632   epot[F_LJ]       = sum_v(grpp->nener,grpp->ener[egLJSR]);
633   epot[F_LJ14]     = sum_v(grpp->nener,grpp->ener[egLJ14]);
634   epot[F_COUL14]   = sum_v(grpp->nener,grpp->ener[egCOUL14]);
635   epot[F_COUL_LR]  = sum_v(grpp->nener,grpp->ener[egCOULLR]);
636   epot[F_LJ_LR]    = sum_v(grpp->nener,grpp->ener[egLJLR]);
637   /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
638   epot[F_GBPOL]   += sum_v(grpp->nener,grpp->ener[egGB]);
639     
640 /* lattice part of LR doesnt belong to any group
641  * and has been added earlier
642  */
643   epot[F_BHAM]     = sum_v(grpp->nener,grpp->ener[egBHAMSR]);
644   epot[F_BHAM_LR]  = sum_v(grpp->nener,grpp->ener[egBHAMLR]);
645
646   epot[F_EPOT] = 0;
647   for(i=0; (i<F_EPOT); i++)
648     if (i != F_DISRESVIOL && i != F_ORIRESDEV && i != F_DIHRESVIOL)
649       epot[F_EPOT] += epot[i];
650 }
651
652 void sum_dhdl(gmx_enerdata_t *enerd,double lambda,t_inputrec *ir)
653 {
654     int i;
655     double dlam,dhdl_lin;
656
657     enerd->term[F_DVDL] = enerd->dvdl_lin + enerd->dvdl_nonlin;
658
659     if (debug)
660     {
661         fprintf(debug,"dvdl: %f, non-linear %f + linear %f\n",
662                 enerd->term[F_DVDL],enerd->dvdl_nonlin,enerd->dvdl_lin);
663     }
664
665     /* Notes on the foreign lambda free energy difference evaluation:
666      * Adding the potential and ekin terms that depend linearly on lambda
667      * as delta lam * dvdl to the energy differences is exact.
668      * For the constraint dvdl this is not exact, but we have no other option.
669      * For the non-bonded LR term we assume that the soft-core (if present)
670      * no longer affects the energy beyond the short-range cut-off,
671      * which is a very good approximation (except for exotic settings).
672      */
673     for(i=1; i<enerd->n_lambda; i++)
674     {
675         dlam = (ir->flambda[i-1] - lambda);
676         dhdl_lin =
677             enerd->dvdl_lin + enerd->term[F_DKDL] + enerd->term[F_DHDL_CON];
678         if (debug)
679         {
680             fprintf(debug,"enerdiff lam %g: non-linear %f linear %f*%f\n",
681                     ir->flambda[i-1],
682                     enerd->enerpart_lambda[i] - enerd->enerpart_lambda[0],
683                     dlam,dhdl_lin);
684         }
685         enerd->enerpart_lambda[i] += dlam*dhdl_lin;
686
687     }
688 }
689
690 void reset_enerdata(t_grpopts *opts,
691                     t_forcerec *fr,gmx_bool bNS,
692                     gmx_enerdata_t *enerd,
693                     gmx_bool bMaster)
694 {
695   gmx_bool bKeepLR;
696   int  i,j;
697   
698   /* First reset all energy components, except for the long range terms
699    * on the master at non neighbor search steps, since the long range
700    * terms have already been summed at the last neighbor search step.
701    */
702   bKeepLR = (fr->bTwinRange && !bNS);
703   for(i=0; (i<egNR); i++) {
704     if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR))) {
705       for(j=0; (j<enerd->grpp.nener); j++)
706         enerd->grpp.ener[i][j] = 0.0;
707     }
708   }
709   enerd->dvdl_lin    = 0.0;
710   enerd->dvdl_nonlin = 0.0;
711
712   /* Normal potential energy components */
713   for(i=0; (i<=F_EPOT); i++) {
714     enerd->term[i] = 0.0;
715   }
716   /* Initialize the dVdlambda term with the long range contribution */
717   enerd->term[F_DVDL]     = 0.0;
718   enerd->term[F_DKDL]     = 0.0;
719   enerd->term[F_DHDL_CON] = 0.0;
720   if (enerd->n_lambda > 0)
721   {
722       for(i=0; i<enerd->n_lambda; i++)
723       {
724           enerd->enerpart_lambda[i] = 0.0;
725       }
726   }
727 }