2432b763c4368c31191e67b418d4b712582de5ed
[alexxy/gromacs.git] / src / 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 <assert.h>
43 #include "sysstuff.h"
44 #include "typedefs.h"
45 #include "macros.h"
46 #include "smalloc.h"
47 #include "macros.h"
48 #include "physics.h"
49 #include "force.h"
50 #include "nonbonded.h"
51 #include "names.h"
52 #include "network.h"
53 #include "pbc.h"
54 #include "ns.h"
55 #include "nrnb.h"
56 #include "bondf.h"
57 #include "mshift.h"
58 #include "txtdump.h"
59 #include "coulomb.h"
60 #include "pme.h"
61 #include "mdrun.h"
62 #include "domdec.h"
63 #include "partdec.h"
64 #include "qmmm.h"
65 #include "mpelogging.h"
66 #include "gmx_omp_nthreads.h"
67
68 void ns(FILE *fp,
69         t_forcerec *fr,
70         rvec       x[],
71         matrix     box,
72         gmx_groups_t *groups,
73         t_grpopts  *opts,
74         gmx_localtop_t *top,
75         t_mdatoms  *md,
76         t_commrec  *cr,
77         t_nrnb     *nrnb,
78         real       *lambda,
79         real       *dvdlambda,
80         gmx_grppairener_t *grppener,
81         gmx_bool       bFillGrid,
82         gmx_bool       bDoLongRangeNS)
83 {
84   char   *ptr;
85   int    nsearch;
86
87   GMX_MPE_LOG(ev_ns_start);
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,bDoLongRangeNS,TRUE);
99   if (debug)
100     fprintf(debug,"nsearch = %d\n",nsearch);
101
102   /* Check whether we have to do dynamic load balancing */
103   /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
104     count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
105     &(top->idef),opts->ngener);
106   */
107   if (fr->ns.dump_nl > 0)
108     dump_nblist(fp,cr,fr,fr->ns.dump_nl);
109
110   GMX_MPE_LOG(ev_ns_finish);
111 }
112
113 static void reduce_thread_forces(int n,rvec *f,
114                                  tensor vir,
115                                  real *Vcorr,
116                                  int efpt_ind,real *dvdl,
117                                  int nthreads,f_thread_t *f_t)
118 {
119     int t,i;
120
121     /* This reduction can run over any number of threads */
122 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntBonded)) private(t) schedule(static)
123     for(i=0; i<n; i++)
124     {
125         for(t=1; t<nthreads; t++)
126         {
127             rvec_inc(f[i],f_t[t].f[i]);
128         }
129     }
130     for(t=1; t<nthreads; t++)
131     {
132         *Vcorr += f_t[t].Vcorr;
133         *dvdl  += f_t[t].dvdl[efpt_ind];
134         m_add(vir,f_t[t].vir,vir);
135     }
136 }
137
138 void do_force_lowlevel(FILE       *fplog,   gmx_large_int_t step,
139                        t_forcerec *fr,      t_inputrec *ir,
140                        t_idef     *idef,    t_commrec  *cr,
141                        t_nrnb     *nrnb,    gmx_wallcycle_t wcycle,
142                        t_mdatoms  *md,
143                        t_grpopts  *opts,
144                        rvec       x[],      history_t  *hist,
145                        rvec       f[],
146                        rvec       f_longrange[],
147                        gmx_enerdata_t *enerd,
148                        t_fcdata   *fcd,
149                        gmx_mtop_t     *mtop,
150                        gmx_localtop_t *top,
151                        gmx_genborn_t *born,
152                        t_atomtypes *atype,
153                        gmx_bool       bBornRadii,
154                        matrix     box,
155                        t_lambda   *fepvals,
156                        real       *lambda,
157                        t_graph    *graph,
158                        t_blocka   *excl,
159                        rvec       mu_tot[],
160                        int        flags,
161                        float      *cycles_pme)
162 {
163     int     i,j,status;
164     int     donb_flags;
165     gmx_bool    bDoEpot,bSepDVDL,bSB;
166     int     pme_flags;
167     matrix  boxs;
168     rvec    box_size;
169     real    Vsr,Vlr,Vcorr=0;
170     t_pbc   pbc;
171     real    dvdgb;
172     char    buf[22];
173     double  clam_i,vlam_i;
174     real    dvdl_dum[efptNR], dvdl, dvdl_nb[efptNR], lam_i[efptNR];
175     real    dvdlsum;
176
177 #ifdef GMX_MPI
178     double  t0=0.0,t1,t2,t3; /* time measurement for coarse load balancing */
179 #endif
180
181 #define PRINT_SEPDVDL(s,v,dvdlambda) if (bSepDVDL) fprintf(fplog,sepdvdlformat,s,v,dvdlambda);
182
183     GMX_MPE_LOG(ev_force_start);
184     set_pbc(&pbc,fr->ePBC,box);
185
186     /* reset free energy components */
187     for (i=0;i<efptNR;i++)
188     {
189         dvdl_nb[i]  = 0;
190         dvdl_dum[i] = 0;
191     }
192
193     /* Reset box */
194     for(i=0; (i<DIM); i++)
195     {
196         box_size[i]=box[i][i];
197     }
198
199     bSepDVDL=(fr->bSepDVDL && do_per_step(step,ir->nstlog));
200     debug_gmx();
201
202     /* do QMMM first if requested */
203     if(fr->bQMMM)
204     {
205         enerd->term[F_EQM] = calculate_QMMM(cr,x,f,fr,md);
206     }
207
208     if (bSepDVDL)
209     {
210         fprintf(fplog,"Step %s: non-bonded V and dVdl for node %d:\n",
211                 gmx_step_str(step,buf),cr->nodeid);
212     }
213
214     /* Call the short range functions all in one go. */
215     GMX_MPE_LOG(ev_do_fnbf_start);
216
217 #ifdef GMX_MPI
218     /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
219 #define TAKETIME FALSE
220     if (TAKETIME)
221     {
222         MPI_Barrier(cr->mpi_comm_mygroup);
223         t0=MPI_Wtime();
224     }
225 #endif
226
227     if (ir->nwall)
228     {
229         /* foreign lambda component for walls */
230         dvdl = do_walls(ir,fr,box,md,x,f,lambda[efptVDW],
231                         enerd->grpp.ener[egLJSR],nrnb);
232         PRINT_SEPDVDL("Walls",0.0,dvdl);
233         enerd->dvdl_lin[efptVDW] += dvdl;
234     }
235
236         /* If doing GB, reset dvda and calculate the Born radii */
237         if (ir->implicit_solvent)
238         {
239         wallcycle_sub_start(wcycle, ewcsNONBONDED);
240
241                 for(i=0;i<born->nr;i++)
242                 {
243                         fr->dvda[i]=0;
244                 }
245
246                 if(bBornRadii)
247                 {
248                         calc_gb_rad(cr,fr,ir,top,atype,x,&(fr->gblist),born,md,nrnb);
249                 }
250
251         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
252         }
253
254     where();
255     /* We only do non-bonded calculation with group scheme here, the verlet
256      * calls are done from do_force_cutsVERLET(). */
257     if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
258     {
259         donb_flags = 0;
260         /* Add short-range interactions */
261         donb_flags |= GMX_NONBONDED_DO_SR;
262
263         if (flags & GMX_FORCE_FORCES)
264         {
265             donb_flags |= GMX_NONBONDED_DO_FORCE;
266         }
267         if (flags & GMX_FORCE_ENERGY)
268         {
269             donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
270         }
271         if (flags & GMX_FORCE_DO_LR)
272         {
273             donb_flags |= GMX_NONBONDED_DO_LR;
274         }
275
276         wallcycle_sub_start(wcycle, ewcsNONBONDED);
277         do_nonbonded(cr,fr,x,f,f_longrange,md,excl,
278                     &enerd->grpp,box_size,nrnb,
279                     lambda,dvdl_nb,-1,-1,donb_flags);
280
281         /* If we do foreign lambda and we have soft-core interactions
282          * we have to recalculate the (non-linear) energies contributions.
283          */
284         if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
285         {
286             for(i=0; i<enerd->n_lambda; i++)
287             {
288                 for (j=0;j<efptNR;j++)
289                 {
290                     lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
291                 }
292                 reset_foreign_enerdata(enerd);
293                 do_nonbonded(cr,fr,x,f,f_longrange,md,excl,
294                              &(enerd->foreign_grpp),box_size,nrnb,
295                              lam_i,dvdl_dum,-1,-1,
296                              (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
297                 sum_epot(&ir->opts,&(enerd->foreign_grpp),enerd->foreign_term);
298                 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
299             }
300         }
301         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
302         where();
303     }
304
305         /* If we are doing GB, calculate bonded forces and apply corrections
306          * to the solvation forces */
307     /* MRS: Eventually, many need to include free energy contribution here! */
308         if (ir->implicit_solvent)
309     {
310         wallcycle_sub_start(wcycle, ewcsBONDED);
311                 calc_gb_forces(cr,md,born,top,atype,x,f,fr,idef,
312                        ir->gb_algorithm,ir->sa_algorithm,nrnb,bBornRadii,&pbc,graph,enerd);
313         wallcycle_sub_stop(wcycle, ewcsBONDED);
314     }
315
316 #ifdef GMX_MPI
317     if (TAKETIME)
318     {
319         t1=MPI_Wtime();
320         fr->t_fnbf += t1-t0;
321     }
322 #endif
323
324     if (fepvals->sc_alpha!=0)
325     {
326         enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
327     }
328     else
329     {
330         enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
331     }
332
333     if (fepvals->sc_alpha!=0)
334
335         /* even though coulomb part is linear, we already added it, beacuse we
336            need to go through the vdw calculation anyway */
337     {
338         enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
339     }
340     else
341     {
342         enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
343     }
344
345     Vsr = 0;
346     if (bSepDVDL)
347     {
348         for(i=0; i<enerd->grpp.nener; i++)
349         {
350             Vsr +=
351                 (fr->bBHAM ?
352                  enerd->grpp.ener[egBHAMSR][i] :
353                  enerd->grpp.ener[egLJSR][i])
354                 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
355         }
356         dvdlsum = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
357         PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",Vsr,dvdlsum);
358     }
359     debug_gmx();
360
361     GMX_MPE_LOG(ev_do_fnbf_finish);
362
363     if (debug)
364     {
365         pr_rvecs(debug,0,"fshift after SR",fr->fshift,SHIFTS);
366     }
367
368     /* Shift the coordinates. Must be done before bonded forces and PPPM,
369      * but is also necessary for SHAKE and update, therefore it can NOT
370      * go when no bonded forces have to be evaluated.
371      */
372
373     /* Here sometimes we would not need to shift with NBFonly,
374      * but we do so anyhow for consistency of the returned coordinates.
375      */
376     if (graph)
377     {
378         shift_self(graph,box,x);
379         if (TRICLINIC(box))
380         {
381             inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
382         }
383         else
384         {
385             inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
386         }
387     }
388     /* Check whether we need to do bondeds or correct for exclusions */
389     if (fr->bMolPBC &&
390         ((flags & GMX_FORCE_BONDED)
391          || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
392     {
393         /* Since all atoms are in the rectangular or triclinic unit-cell,
394          * only single box vector shifts (2 in x) are required.
395          */
396         set_pbc_dd(&pbc,fr->ePBC,cr->dd,TRUE,box);
397     }
398     debug_gmx();
399
400     if (flags & GMX_FORCE_BONDED)
401     {
402         GMX_MPE_LOG(ev_calc_bonds_start);
403
404         wallcycle_sub_start(wcycle, ewcsBONDED);
405         calc_bonds(fplog,cr->ms,
406                    idef,x,hist,f,fr,&pbc,graph,enerd,nrnb,lambda,md,fcd,
407                    DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
408                    flags,
409                    fr->bSepDVDL && do_per_step(step,ir->nstlog),step);
410
411         /* Check if we have to determine energy differences
412          * at foreign lambda's.
413          */
414         if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
415             idef->ilsort != ilsortNO_FE)
416         {
417             if (idef->ilsort != ilsortFE_SORTED)
418             {
419                 gmx_incons("The bonded interactions are not sorted for free energy");
420             }
421             for(i=0; i<enerd->n_lambda; i++)
422             {
423                 reset_foreign_enerdata(enerd);
424                 for (j=0;j<efptNR;j++)
425                 {
426                     lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
427                 }
428                 calc_bonds_lambda(fplog,idef,x,fr,&pbc,graph,&(enerd->foreign_grpp),enerd->foreign_term,nrnb,lam_i,md,
429                                   fcd,DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
430                 sum_epot(&ir->opts,&(enerd->foreign_grpp),enerd->foreign_term);
431                 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
432             }
433         }
434         debug_gmx();
435         GMX_MPE_LOG(ev_calc_bonds_finish);
436         wallcycle_sub_stop(wcycle, ewcsBONDED);
437     }
438
439     where();
440
441     *cycles_pme = 0;
442     if (EEL_FULL(fr->eeltype))
443     {
444         bSB = (ir->nwall == 2);
445         if (bSB)
446         {
447             copy_mat(box,boxs);
448             svmul(ir->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
449             box_size[ZZ] *= ir->wall_ewald_zfac;
450         }
451
452         clear_mat(fr->vir_el_recip);
453
454         if (fr->bEwald)
455         {
456             Vcorr = 0;
457             dvdl  = 0;
458
459             /* With the Verlet scheme exclusion forces are calculated
460              * in the non-bonded kernel.
461              */
462             /* The TPI molecule does not have exclusions with the rest
463              * of the system and no intra-molecular PME grid contributions
464              * will be calculated in gmx_pme_calc_energy.
465              */
466             if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
467                 ir->ewald_geometry != eewg3D ||
468                 ir->epsilon_surface != 0)
469             {
470                 int nthreads,t;
471
472                 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
473
474                 if (fr->n_tpi > 0)
475                 {
476                     gmx_fatal(FARGS,"TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
477                 }
478
479                 nthreads = gmx_omp_nthreads_get(emntBonded);
480 #pragma omp parallel for num_threads(nthreads) schedule(static)
481                 for(t=0; t<nthreads; t++)
482                 {
483                     int s,e,i;
484                     rvec *fnv;
485                     tensor *vir;
486                     real *Vcorrt,*dvdlt;
487                     if (t == 0)
488                     {
489                         fnv    = fr->f_novirsum;
490                         vir    = &fr->vir_el_recip;
491                         Vcorrt = &Vcorr;
492                         dvdlt  = &dvdl;
493                     }
494                     else
495                     {
496                         fnv    = fr->f_t[t].f;
497                         vir    = &fr->f_t[t].vir;
498                         Vcorrt = &fr->f_t[t].Vcorr;
499                         dvdlt  = &fr->f_t[t].dvdl[efptCOUL];
500                         for(i=0; i<fr->natoms_force; i++)
501                         {
502                             clear_rvec(fnv[i]);
503                         }
504                         clear_mat(*vir);
505                     }
506                     *dvdlt = 0;
507                     *Vcorrt =
508                         ewald_LRcorrection(fplog,
509                                            fr->excl_load[t],fr->excl_load[t+1],
510                                            cr,t,fr,
511                                            md->chargeA,
512                                            md->nChargePerturbed ? md->chargeB : NULL,
513                                            ir->cutoff_scheme != ecutsVERLET,
514                                            excl,x,bSB ? boxs : box,mu_tot,
515                                            ir->ewald_geometry,
516                                            ir->epsilon_surface,
517                                            fnv,*vir,
518                                            lambda[efptCOUL],dvdlt);
519                 }
520                 if (nthreads > 1)
521                 {
522                     reduce_thread_forces(fr->natoms_force,fr->f_novirsum,
523                                          fr->vir_el_recip,
524                                          &Vcorr,efptCOUL,&dvdl,
525                                          nthreads,fr->f_t);
526                 }
527
528                 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
529             }
530
531             if (fr->n_tpi == 0)
532             {
533                 Vcorr += ewald_charge_correction(cr,fr,lambda[efptCOUL],box,
534                                                  &dvdl,fr->vir_el_recip);
535             }
536
537             PRINT_SEPDVDL("Ewald excl./charge/dip. corr.",Vcorr,dvdl);
538             enerd->dvdl_lin[efptCOUL] += dvdl;
539         }
540
541         status = 0;
542         Vlr  = 0;
543         dvdl = 0;
544         switch (fr->eeltype)
545         {
546         case eelPME:
547         case eelPMESWITCH:
548         case eelPMEUSER:
549         case eelPMEUSERSWITCH:
550         case eelP3M_AD:
551             if (cr->duty & DUTY_PME)
552             {
553                 assert(fr->n_tpi >= 0);
554                 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
555                 {
556                     pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
557                     if (flags & GMX_FORCE_FORCES)
558                     {
559                         pme_flags |= GMX_PME_CALC_F;
560                     }
561                     if (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY))
562                     {
563                         pme_flags |= GMX_PME_CALC_ENER_VIR;
564                     }
565                     if (fr->n_tpi > 0)
566                     {
567                         /* We don't calculate f, but we do want the potential */
568                         pme_flags |= GMX_PME_CALC_POT;
569                     }
570                     wallcycle_start(wcycle,ewcPMEMESH);
571                     status = gmx_pme_do(fr->pmedata,
572                                         md->start,md->homenr - fr->n_tpi,
573                                         x,fr->f_novirsum,
574                                         md->chargeA,md->chargeB,
575                                         bSB ? boxs : box,cr,
576                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
577                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
578                                         nrnb,wcycle,
579                                         fr->vir_el_recip,fr->ewaldcoeff,
580                                         &Vlr,lambda[efptCOUL],&dvdl,
581                                         pme_flags);
582                     *cycles_pme = wallcycle_stop(wcycle,ewcPMEMESH);
583
584                     /* We should try to do as little computation after
585                      * this as possible, because parallel PME synchronizes
586                      * the nodes, so we want all load imbalance of the rest
587                      * of the force calculation to be before the PME call.
588                      * DD load balancing is done on the whole time of
589                      * the force call (without PME).
590                      */
591                 }
592                 if (fr->n_tpi > 0)
593                 {
594                     /* Determine the PME grid energy of the test molecule
595                      * with the PME grid potential of the other charges.
596                      */
597                     gmx_pme_calc_energy(fr->pmedata,fr->n_tpi,
598                                         x + md->homenr - fr->n_tpi,
599                                         md->chargeA + md->homenr - fr->n_tpi,
600                                         &Vlr);
601                 }
602                 PRINT_SEPDVDL("PME mesh",Vlr,dvdl);
603             }
604             break;
605         case eelEWALD:
606             Vlr = do_ewald(fplog,FALSE,ir,x,fr->f_novirsum,
607                            md->chargeA,md->chargeB,
608                            box_size,cr,md->homenr,
609                            fr->vir_el_recip,fr->ewaldcoeff,
610                            lambda[efptCOUL],&dvdl,fr->ewald_table);
611             PRINT_SEPDVDL("Ewald long-range",Vlr,dvdl);
612             break;
613         default:
614             gmx_fatal(FARGS,"No such electrostatics method implemented %s",
615                       eel_names[fr->eeltype]);
616         }
617         if (status != 0)
618         {
619             gmx_fatal(FARGS,"Error %d in long range electrostatics routine %s",
620                       status,EELTYPE(fr->eeltype));
621                 }
622         /* Note that with separate PME nodes we get the real energies later */
623         enerd->dvdl_lin[efptCOUL] += dvdl;
624         enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
625         if (debug)
626         {
627             fprintf(debug,"Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
628                     Vlr,Vcorr,enerd->term[F_COUL_RECIP]);
629             pr_rvecs(debug,0,"vir_el_recip after corr",fr->vir_el_recip,DIM);
630             pr_rvecs(debug,0,"fshift after LR Corrections",fr->fshift,SHIFTS);
631         }
632     }
633     else
634     {
635         if (EEL_RF(fr->eeltype))
636         {
637             /* With the Verlet scheme exclusion forces are calculated
638              * in the non-bonded kernel.
639              */
640             if (ir->cutoff_scheme != ecutsVERLET && fr->eeltype != eelRF_NEC)
641             {
642                 dvdl = 0;
643                 enerd->term[F_RF_EXCL] =
644                     RF_excl_correction(fplog,fr,graph,md,excl,x,f,
645                                        fr->fshift,&pbc,lambda[efptCOUL],&dvdl);
646             }
647
648             enerd->dvdl_lin[efptCOUL] += dvdl;
649             PRINT_SEPDVDL("RF exclusion correction",
650                           enerd->term[F_RF_EXCL],dvdl);
651         }
652     }
653     where();
654     debug_gmx();
655
656     if (debug)
657     {
658         print_nrnb(debug,nrnb);
659     }
660     debug_gmx();
661
662 #ifdef GMX_MPI
663     if (TAKETIME)
664     {
665         t2=MPI_Wtime();
666         MPI_Barrier(cr->mpi_comm_mygroup);
667         t3=MPI_Wtime();
668         fr->t_wait += t3-t2;
669         if (fr->timesteps == 11)
670         {
671             fprintf(stderr,"* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
672                     cr->nodeid, gmx_step_str(fr->timesteps,buf),
673                     100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
674                     (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
675         }
676         fr->timesteps++;
677     }
678 #endif
679
680     if (debug)
681     {
682         pr_rvecs(debug,0,"fshift after bondeds",fr->fshift,SHIFTS);
683     }
684
685     GMX_MPE_LOG(ev_force_finish);
686
687 }
688
689 void init_enerdata(int ngener,int n_lambda,gmx_enerdata_t *enerd)
690 {
691     int i,n2;
692
693     for(i=0; i<F_NRE; i++)
694     {
695         enerd->term[i] = 0;
696         enerd->foreign_term[i] = 0;
697     }
698
699
700     for(i=0; i<efptNR; i++) {
701         enerd->dvdl_lin[i]  = 0;
702         enerd->dvdl_nonlin[i]  = 0;
703     }
704
705     n2=ngener*ngener;
706     if (debug)
707     {
708         fprintf(debug,"Creating %d sized group matrix for energies\n",n2);
709     }
710     enerd->grpp.nener = n2;
711     enerd->foreign_grpp.nener = n2;
712     for(i=0; (i<egNR); i++)
713     {
714         snew(enerd->grpp.ener[i],n2);
715         snew(enerd->foreign_grpp.ener[i],n2);
716     }
717
718     if (n_lambda)
719     {
720         enerd->n_lambda = 1 + n_lambda;
721         snew(enerd->enerpart_lambda,enerd->n_lambda);
722     }
723     else
724     {
725         enerd->n_lambda = 0;
726     }
727 }
728
729 void destroy_enerdata(gmx_enerdata_t *enerd)
730 {
731     int i;
732
733     for(i=0; (i<egNR); i++)
734     {
735         sfree(enerd->grpp.ener[i]);
736     }
737
738     for(i=0; (i<egNR); i++)
739     {
740         sfree(enerd->foreign_grpp.ener[i]);
741     }
742
743     if (enerd->n_lambda)
744     {
745         sfree(enerd->enerpart_lambda);
746     }
747 }
748
749 static real sum_v(int n,real v[])
750 {
751   real t;
752   int  i;
753
754   t = 0.0;
755   for(i=0; (i<n); i++)
756     t = t + v[i];
757
758   return t;
759 }
760
761 void sum_epot(t_grpopts *opts, gmx_grppairener_t *grpp, real *epot)
762 {
763   int i;
764
765   /* Accumulate energies */
766   epot[F_COUL_SR]  = sum_v(grpp->nener,grpp->ener[egCOULSR]);
767   epot[F_LJ]       = sum_v(grpp->nener,grpp->ener[egLJSR]);
768   epot[F_LJ14]     = sum_v(grpp->nener,grpp->ener[egLJ14]);
769   epot[F_COUL14]   = sum_v(grpp->nener,grpp->ener[egCOUL14]);
770   epot[F_COUL_LR]  = sum_v(grpp->nener,grpp->ener[egCOULLR]);
771   epot[F_LJ_LR]    = sum_v(grpp->nener,grpp->ener[egLJLR]);
772   /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
773   epot[F_GBPOL]   += sum_v(grpp->nener,grpp->ener[egGB]);
774
775 /* lattice part of LR doesnt belong to any group
776  * and has been added earlier
777  */
778   epot[F_BHAM]     = sum_v(grpp->nener,grpp->ener[egBHAMSR]);
779   epot[F_BHAM_LR]  = sum_v(grpp->nener,grpp->ener[egBHAMLR]);
780
781   epot[F_EPOT] = 0;
782   for(i=0; (i<F_EPOT); i++)
783   {
784       if (i != F_DISRESVIOL && i != F_ORIRESDEV)
785       {
786           epot[F_EPOT] += epot[i];
787       }
788   }
789 }
790
791 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
792 {
793     int i,j,index;
794     double dlam;
795
796     enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW];  /* include dispersion correction */
797     enerd->term[F_DVDL] = 0.0;
798     for (i=0;i<efptNR;i++)
799     {
800         if (fepvals->separate_dvdl[i])
801         {
802             /* could this be done more readably/compactly? */
803             switch (i) {
804             case (efptCOUL):
805                 index = F_DVDL_COUL;
806                 break;
807             case (efptVDW):
808                 index = F_DVDL_VDW;
809                 break;
810             case (efptBONDED):
811                 index = F_DVDL_BONDED;
812                 break;
813             case (efptRESTRAINT):
814                 index = F_DVDL_RESTRAINT;
815                 break;
816             case (efptMASS):
817                 index = F_DKDL;
818                 break;
819             default:
820                 index = F_DVDL;
821                 break;
822             }
823             enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
824             if (debug)
825             {
826                 fprintf(debug,"dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
827                         efpt_names[i],i,enerd->term[index],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
828             }
829         }
830         else
831         {
832             enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
833             if (debug)
834             {
835                 fprintf(debug,"dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
836                         efpt_names[0],i,enerd->term[F_DVDL],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
837             }
838         }
839     }
840
841     /* Notes on the foreign lambda free energy difference evaluation:
842      * Adding the potential and ekin terms that depend linearly on lambda
843      * as delta lam * dvdl to the energy differences is exact.
844      * For the constraints this is not exact, but we have no other option
845      * without literally changing the lengths and reevaluating the energies at each step.
846      * (try to remedy this post 4.6 - MRS)
847      * For the non-bonded LR term we assume that the soft-core (if present)
848      * no longer affects the energy beyond the short-range cut-off,
849      * which is a very good approximation (except for exotic settings).
850      * (investigate how to overcome this post 4.6 - MRS)
851      */
852
853     for(i=0; i<fepvals->n_lambda; i++)
854     {                                         /* note we are iterating over fepvals here!
855                                                  For the current lam, dlam = 0 automatically,
856                                                  so we don't need to add anything to the
857                                                  enerd->enerpart_lambda[0] */
858
859         /* we don't need to worry about dvdl contributions to the current lambda, because
860            it's automatically zero */
861
862         /* first kinetic energy term */
863         dlam = (fepvals->all_lambda[efptMASS][i] - lambda[efptMASS]);
864
865         enerd->enerpart_lambda[i+1] += enerd->term[F_DKDL]*dlam;
866
867         for (j=0;j<efptNR;j++)
868         {
869             if (j==efptMASS) {continue;} /* no other mass term to worry about */
870
871             dlam = (fepvals->all_lambda[j][i]-lambda[j]);
872             enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
873             if (debug)
874             {
875                 fprintf(debug,"enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
876                         fepvals->all_lambda[j][i],efpt_names[j],
877                         (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
878                         dlam,enerd->dvdl_lin[j]);
879             }
880         }
881     }
882 }
883
884
885 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
886 {
887     int  i,j;
888
889     /* First reset all foreign energy components.  Foreign energies always called on
890        neighbor search steps */
891     for(i=0; (i<egNR); i++)
892     {
893         for(j=0; (j<enerd->grpp.nener); j++)
894         {
895             enerd->foreign_grpp.ener[i][j] = 0.0;
896         }
897     }
898
899     /* potential energy components */
900     for(i=0; (i<=F_EPOT); i++)
901     {
902         enerd->foreign_term[i] = 0.0;
903     }
904 }
905
906 void reset_enerdata(t_grpopts *opts,
907                     t_forcerec *fr,gmx_bool bNS,
908                     gmx_enerdata_t *enerd,
909                     gmx_bool bMaster)
910 {
911     gmx_bool bKeepLR;
912     int  i,j;
913
914     /* First reset all energy components, except for the long range terms
915      * on the master at non neighbor search steps, since the long range
916      * terms have already been summed at the last neighbor search step.
917      */
918     bKeepLR = (fr->bTwinRange && !bNS);
919     for(i=0; (i<egNR); i++) {
920         if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR))) {
921             for(j=0; (j<enerd->grpp.nener); j++)
922                 enerd->grpp.ener[i][j] = 0.0;
923         }
924     }
925     for (i=0;i<efptNR;i++)
926     {
927         enerd->dvdl_lin[i]    = 0.0;
928         enerd->dvdl_nonlin[i] = 0.0;
929     }
930
931     /* Normal potential energy components */
932     for(i=0; (i<=F_EPOT); i++) {
933         enerd->term[i] = 0.0;
934     }
935     /* Initialize the dVdlambda term with the long range contribution */
936     /* Initialize the dvdl term with the long range contribution */
937     enerd->term[F_DVDL]            = 0.0;
938     enerd->term[F_DVDL_COUL]       = 0.0;
939     enerd->term[F_DVDL_VDW]        = 0.0;
940     enerd->term[F_DVDL_BONDED]     = 0.0;
941     enerd->term[F_DVDL_RESTRAINT]  = 0.0;
942     enerd->term[F_DKDL]            = 0.0;
943     if (enerd->n_lambda > 0)
944     {
945         for(i=0; i<enerd->n_lambda; i++)
946         {
947             enerd->enerpart_lambda[i] = 0.0;
948         }
949     }
950     /* reset foreign energy data - separate function since we also call it elsewhere */
951     reset_foreign_enerdata(enerd);
952 }