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