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