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