Move vec.h to math/
[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 "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 "qmmm.h"
64 #include "gmx_omp_nthreads.h"
65
66 #include "gromacs/math/vec.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     if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
452     {
453         real Vlr             = 0, Vcorr = 0;
454         real dvdl_long_range = 0;
455         int  status          = 0;
456
457         bSB = (ir->nwall == 2);
458         if (bSB)
459         {
460             copy_mat(box, boxs);
461             svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
462             box_size[ZZ] *= ir->wall_ewald_zfac;
463         }
464     }
465
466     /* Do long-range electrostatics and/or LJ-PME, including related short-range
467      * corrections.
468      */
469
470     clear_mat(fr->vir_el_recip);
471     clear_mat(fr->vir_lj_recip);
472
473     if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
474     {
475         real Vlr_q             = 0, Vlr_lj = 0, Vcorr_q = 0, Vcorr_lj = 0;
476         real dvdl_long_range_q = 0, dvdl_long_range_lj = 0;
477         int  status            = 0;
478
479         if (EEL_PME_EWALD(fr->eeltype) || EVDW_PME(fr->vdwtype))
480         {
481             real dvdl_long_range_correction_q   = 0;
482             real dvdl_long_range_correction_lj  = 0;
483             /* With the Verlet scheme exclusion forces are calculated
484              * in the non-bonded kernel.
485              */
486             /* The TPI molecule does not have exclusions with the rest
487              * of the system and no intra-molecular PME grid
488              * contributions will be calculated in
489              * gmx_pme_calc_energy.
490              */
491             if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
492                 ir->ewald_geometry != eewg3D ||
493                 ir->epsilon_surface != 0)
494             {
495                 int nthreads, t;
496
497                 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
498
499                 if (fr->n_tpi > 0)
500                 {
501                     gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
502                 }
503
504                 nthreads = gmx_omp_nthreads_get(emntBonded);
505 #pragma omp parallel for num_threads(nthreads) schedule(static)
506                 for (t = 0; t < nthreads; t++)
507                 {
508                     int     s, e, i;
509                     rvec   *fnv;
510                     tensor *vir_q, *vir_lj;
511                     real   *Vcorrt_q, *Vcorrt_lj, *dvdlt_q, *dvdlt_lj;
512                     if (t == 0)
513                     {
514                         fnv       = fr->f_novirsum;
515                         vir_q     = &fr->vir_el_recip;
516                         vir_lj    = &fr->vir_lj_recip;
517                         Vcorrt_q  = &Vcorr_q;
518                         Vcorrt_lj = &Vcorr_lj;
519                         dvdlt_q   = &dvdl_long_range_correction_q;
520                         dvdlt_lj  = &dvdl_long_range_correction_lj;
521                     }
522                     else
523                     {
524                         fnv       = fr->f_t[t].f;
525                         vir_q     = &fr->f_t[t].vir_q;
526                         vir_lj    = &fr->f_t[t].vir_lj;
527                         Vcorrt_q  = &fr->f_t[t].Vcorr_q;
528                         Vcorrt_lj = &fr->f_t[t].Vcorr_lj;
529                         dvdlt_q   = &fr->f_t[t].dvdl[efptCOUL];
530                         dvdlt_lj  = &fr->f_t[t].dvdl[efptVDW];
531                         for (i = 0; i < fr->natoms_force; i++)
532                         {
533                             clear_rvec(fnv[i]);
534                         }
535                         clear_mat(*vir_q);
536                         clear_mat(*vir_lj);
537                     }
538                     *dvdlt_q  = 0;
539                     *dvdlt_lj = 0;
540
541                     ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
542                                        cr, t, fr,
543                                        md->chargeA,
544                                        md->nChargePerturbed ? md->chargeB : NULL,
545                                        md->sqrt_c6A,
546                                        md->nTypePerturbed ? md->sqrt_c6B : NULL,
547                                        md->sigmaA,
548                                        md->nTypePerturbed ? md->sigmaB : NULL,
549                                        md->sigma3A,
550                                        md->nTypePerturbed ? md->sigma3B : NULL,
551                                        ir->cutoff_scheme != ecutsVERLET,
552                                        excl, x, bSB ? boxs : box, mu_tot,
553                                        ir->ewald_geometry,
554                                        ir->epsilon_surface,
555                                        fnv, *vir_q, *vir_lj,
556                                        Vcorrt_q, Vcorrt_lj,
557                                        lambda[efptCOUL], lambda[efptVDW],
558                                        dvdlt_q, dvdlt_lj);
559                 }
560                 if (nthreads > 1)
561                 {
562                     reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
563                                          fr->vir_el_recip, fr->vir_lj_recip,
564                                          &Vcorr_q, &Vcorr_lj,
565                                          &dvdl_long_range_correction_q,
566                                          &dvdl_long_range_correction_lj,
567                                          nthreads, fr->f_t);
568                 }
569                 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
570             }
571
572             if (EEL_PME_EWALD(fr->eeltype) && fr->n_tpi == 0)
573             {
574                 Vcorr_q += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
575                                                    &dvdl_long_range_correction_q,
576                                                    fr->vir_el_recip);
577             }
578
579             PRINT_SEPDVDL("Ewald excl./charge/dip. corr.", Vcorr_q, dvdl_long_range_correction_q);
580             PRINT_SEPDVDL("Ewald excl. corr. LJ", Vcorr_lj, dvdl_long_range_correction_lj);
581             enerd->dvdl_lin[efptCOUL] += dvdl_long_range_correction_q;
582             enerd->dvdl_lin[efptVDW]  += dvdl_long_range_correction_lj;
583         }
584
585         if ((EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype)))
586         {
587             if (cr->duty & DUTY_PME)
588             {
589                 /* Do reciprocal PME for Coulomb and/or LJ. */
590                 assert(fr->n_tpi >= 0);
591                 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
592                 {
593                     pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
594                     if (EEL_PME(fr->eeltype))
595                     {
596                         pme_flags     |= GMX_PME_DO_COULOMB;
597                     }
598                     if (EVDW_PME(fr->vdwtype))
599                     {
600                         pme_flags |= GMX_PME_DO_LJ;
601                     }
602                     if (flags & GMX_FORCE_FORCES)
603                     {
604                         pme_flags |= GMX_PME_CALC_F;
605                     }
606                     if (flags & GMX_FORCE_VIRIAL)
607                     {
608                         pme_flags |= GMX_PME_CALC_ENER_VIR;
609                     }
610                     if (fr->n_tpi > 0)
611                     {
612                         /* We don't calculate f, but we do want the potential */
613                         pme_flags |= GMX_PME_CALC_POT;
614                     }
615                     wallcycle_start(wcycle, ewcPMEMESH);
616                     status = gmx_pme_do(fr->pmedata,
617                                         0, md->homenr - fr->n_tpi,
618                                         x, fr->f_novirsum,
619                                         md->chargeA, md->chargeB,
620                                         md->sqrt_c6A, md->sqrt_c6B,
621                                         md->sigmaA, md->sigmaB,
622                                         bSB ? boxs : box, cr,
623                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
624                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
625                                         nrnb, wcycle,
626                                         fr->vir_el_recip, fr->ewaldcoeff_q,
627                                         fr->vir_lj_recip, fr->ewaldcoeff_lj,
628                                         &Vlr_q, &Vlr_lj,
629                                         lambda[efptCOUL], lambda[efptVDW],
630                                         &dvdl_long_range_q, &dvdl_long_range_lj, pme_flags);
631                     *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
632                     if (status != 0)
633                     {
634                         gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
635                     }
636                     /* We should try to do as little computation after
637                      * this as possible, because parallel PME synchronizes
638                      * the nodes, so we want all load imbalance of the
639                      * rest of the force calculation to be before the PME
640                      * call.  DD load balancing is done on the whole time
641                      * of the force call (without PME).
642                      */
643                 }
644                 if (fr->n_tpi > 0)
645                 {
646                     if (EVDW_PME(ir->vdwtype))
647                     {
648
649                         gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
650                     }
651                     /* Determine the PME grid energy of the test molecule
652                      * with the PME grid potential of the other charges.
653                      */
654                     gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
655                                         x + md->homenr - fr->n_tpi,
656                                         md->chargeA + md->homenr - fr->n_tpi,
657                                         &Vlr_q);
658                 }
659                 PRINT_SEPDVDL("PME mesh", Vlr_q + Vlr_lj, dvdl_long_range_q+dvdl_long_range_lj);
660             }
661         }
662
663         if (!EEL_PME(fr->eeltype) && EEL_PME_EWALD(fr->eeltype))
664         {
665             Vlr_q = do_ewald(ir, x, fr->f_novirsum,
666                              md->chargeA, md->chargeB,
667                              box_size, cr, md->homenr,
668                              fr->vir_el_recip, fr->ewaldcoeff_q,
669                              lambda[efptCOUL], &dvdl_long_range_q, fr->ewald_table);
670             PRINT_SEPDVDL("Ewald long-range", Vlr_q, dvdl_long_range_q);
671         }
672
673         /* Note that with separate PME nodes we get the real energies later */
674         enerd->dvdl_lin[efptCOUL] += dvdl_long_range_q;
675         enerd->dvdl_lin[efptVDW]  += dvdl_long_range_lj;
676         enerd->term[F_COUL_RECIP]  = Vlr_q + Vcorr_q;
677         enerd->term[F_LJ_RECIP]    = Vlr_lj + Vcorr_lj;
678         if (debug)
679         {
680             fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
681                     Vlr_q, Vcorr_q, enerd->term[F_COUL_RECIP]);
682             pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
683             pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
684             fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
685                     Vlr_lj, Vcorr_lj, enerd->term[F_LJ_RECIP]);
686             pr_rvecs(debug, 0, "vir_lj_recip after corr", fr->vir_lj_recip, DIM);
687         }
688     }
689     else
690     {
691         /* Is there a reaction-field exclusion correction needed? */
692         if (EEL_RF(fr->eeltype) && eelRF_NEC != fr->eeltype)
693         {
694             /* With the Verlet scheme, exclusion forces are calculated
695              * in the non-bonded kernel.
696              */
697             if (ir->cutoff_scheme != ecutsVERLET)
698             {
699                 real dvdl_rf_excl      = 0;
700                 enerd->term[F_RF_EXCL] =
701                     RF_excl_correction(fr, graph, md, excl, x, f,
702                                        fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
703
704                 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
705                 PRINT_SEPDVDL("RF exclusion correction",
706                               enerd->term[F_RF_EXCL], dvdl_rf_excl);
707             }
708         }
709     }
710     where();
711     debug_gmx();
712
713     if (debug)
714     {
715         print_nrnb(debug, nrnb);
716     }
717     debug_gmx();
718
719 #ifdef GMX_MPI
720     if (TAKETIME)
721     {
722         t2 = MPI_Wtime();
723         MPI_Barrier(cr->mpi_comm_mygroup);
724         t3          = MPI_Wtime();
725         fr->t_wait += t3-t2;
726         if (fr->timesteps == 11)
727         {
728             fprintf(stderr, "* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
729                     cr->nodeid, gmx_step_str(fr->timesteps, buf),
730                     100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
731                     (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
732         }
733         fr->timesteps++;
734     }
735 #endif
736
737     if (debug)
738     {
739         pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
740     }
741
742 }
743
744 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
745 {
746     int i, n2;
747
748     for (i = 0; i < F_NRE; i++)
749     {
750         enerd->term[i]         = 0;
751         enerd->foreign_term[i] = 0;
752     }
753
754
755     for (i = 0; i < efptNR; i++)
756     {
757         enerd->dvdl_lin[i]     = 0;
758         enerd->dvdl_nonlin[i]  = 0;
759     }
760
761     n2 = ngener*ngener;
762     if (debug)
763     {
764         fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
765     }
766     enerd->grpp.nener         = n2;
767     enerd->foreign_grpp.nener = n2;
768     for (i = 0; (i < egNR); i++)
769     {
770         snew(enerd->grpp.ener[i], n2);
771         snew(enerd->foreign_grpp.ener[i], n2);
772     }
773
774     if (n_lambda)
775     {
776         enerd->n_lambda = 1 + n_lambda;
777         snew(enerd->enerpart_lambda, enerd->n_lambda);
778     }
779     else
780     {
781         enerd->n_lambda = 0;
782     }
783 }
784
785 void destroy_enerdata(gmx_enerdata_t *enerd)
786 {
787     int i;
788
789     for (i = 0; (i < egNR); i++)
790     {
791         sfree(enerd->grpp.ener[i]);
792     }
793
794     for (i = 0; (i < egNR); i++)
795     {
796         sfree(enerd->foreign_grpp.ener[i]);
797     }
798
799     if (enerd->n_lambda)
800     {
801         sfree(enerd->enerpart_lambda);
802     }
803 }
804
805 static real sum_v(int n, real v[])
806 {
807     real t;
808     int  i;
809
810     t = 0.0;
811     for (i = 0; (i < n); i++)
812     {
813         t = t + v[i];
814     }
815
816     return t;
817 }
818
819 void sum_epot(gmx_grppairener_t *grpp, real *epot)
820 {
821     int i;
822
823     /* Accumulate energies */
824     epot[F_COUL_SR]  = sum_v(grpp->nener, grpp->ener[egCOULSR]);
825     epot[F_LJ]       = sum_v(grpp->nener, grpp->ener[egLJSR]);
826     epot[F_LJ14]     = sum_v(grpp->nener, grpp->ener[egLJ14]);
827     epot[F_COUL14]   = sum_v(grpp->nener, grpp->ener[egCOUL14]);
828     epot[F_COUL_LR]  = sum_v(grpp->nener, grpp->ener[egCOULLR]);
829     epot[F_LJ_LR]    = sum_v(grpp->nener, grpp->ener[egLJLR]);
830     /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
831     epot[F_GBPOL]   += sum_v(grpp->nener, grpp->ener[egGB]);
832
833 /* lattice part of LR doesnt belong to any group
834  * and has been added earlier
835  */
836     epot[F_BHAM]     = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
837     epot[F_BHAM_LR]  = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
838
839     epot[F_EPOT] = 0;
840     for (i = 0; (i < F_EPOT); i++)
841     {
842         if (i != F_DISRESVIOL && i != F_ORIRESDEV)
843         {
844             epot[F_EPOT] += epot[i];
845         }
846     }
847 }
848
849 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
850 {
851     int    i, j, index;
852     double dlam;
853
854     enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW];  /* include dispersion correction */
855     enerd->term[F_DVDL]       = 0.0;
856     for (i = 0; i < efptNR; i++)
857     {
858         if (fepvals->separate_dvdl[i])
859         {
860             /* could this be done more readably/compactly? */
861             switch (i)
862             {
863                 case (efptMASS):
864                     index = F_DKDL;
865                     break;
866                 case (efptCOUL):
867                     index = F_DVDL_COUL;
868                     break;
869                 case (efptVDW):
870                     index = F_DVDL_VDW;
871                     break;
872                 case (efptBONDED):
873                     index = F_DVDL_BONDED;
874                     break;
875                 case (efptRESTRAINT):
876                     index = F_DVDL_RESTRAINT;
877                     break;
878                 default:
879                     index = F_DVDL;
880                     break;
881             }
882             enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
883             if (debug)
884             {
885                 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
886                         efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
887             }
888         }
889         else
890         {
891             enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
892             if (debug)
893             {
894                 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
895                         efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
896             }
897         }
898     }
899
900     /* Notes on the foreign lambda free energy difference evaluation:
901      * Adding the potential and ekin terms that depend linearly on lambda
902      * as delta lam * dvdl to the energy differences is exact.
903      * For the constraints this is not exact, but we have no other option
904      * without literally changing the lengths and reevaluating the energies at each step.
905      * (try to remedy this post 4.6 - MRS)
906      * For the non-bonded LR term we assume that the soft-core (if present)
907      * no longer affects the energy beyond the short-range cut-off,
908      * which is a very good approximation (except for exotic settings).
909      * (investigate how to overcome this post 4.6 - MRS)
910      */
911     if (fepvals->separate_dvdl[efptBONDED])
912     {
913         enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
914     }
915     else
916     {
917         enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
918     }
919     enerd->term[F_DVDL_CONSTR] = 0;
920
921     for (i = 0; i < fepvals->n_lambda; i++)
922     {
923         /* note we are iterating over fepvals here!
924            For the current lam, dlam = 0 automatically,
925            so we don't need to add anything to the
926            enerd->enerpart_lambda[0] */
927
928         /* we don't need to worry about dvdl_lin contributions to dE at
929            current lambda, because the contributions to the current
930            lambda are automatically zeroed */
931
932         for (j = 0; j < efptNR; j++)
933         {
934             /* Note that this loop is over all dhdl components, not just the separated ones */
935             dlam = (fepvals->all_lambda[j][i]-lambda[j]);
936             enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
937             if (debug)
938             {
939                 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
940                         fepvals->all_lambda[j][i], efpt_names[j],
941                         (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
942                         dlam, enerd->dvdl_lin[j]);
943             }
944         }
945     }
946 }
947
948
949 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
950 {
951     int  i, j;
952
953     /* First reset all foreign energy components.  Foreign energies always called on
954        neighbor search steps */
955     for (i = 0; (i < egNR); i++)
956     {
957         for (j = 0; (j < enerd->grpp.nener); j++)
958         {
959             enerd->foreign_grpp.ener[i][j] = 0.0;
960         }
961     }
962
963     /* potential energy components */
964     for (i = 0; (i <= F_EPOT); i++)
965     {
966         enerd->foreign_term[i] = 0.0;
967     }
968 }
969
970 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
971                     gmx_enerdata_t *enerd,
972                     gmx_bool bMaster)
973 {
974     gmx_bool bKeepLR;
975     int      i, j;
976
977     /* First reset all energy components, except for the long range terms
978      * on the master at non neighbor search steps, since the long range
979      * terms have already been summed at the last neighbor search step.
980      */
981     bKeepLR = (fr->bTwinRange && !bNS);
982     for (i = 0; (i < egNR); i++)
983     {
984         if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
985         {
986             for (j = 0; (j < enerd->grpp.nener); j++)
987             {
988                 enerd->grpp.ener[i][j] = 0.0;
989             }
990         }
991     }
992     for (i = 0; i < efptNR; i++)
993     {
994         enerd->dvdl_lin[i]    = 0.0;
995         enerd->dvdl_nonlin[i] = 0.0;
996     }
997
998     /* Normal potential energy components */
999     for (i = 0; (i <= F_EPOT); i++)
1000     {
1001         enerd->term[i] = 0.0;
1002     }
1003     /* Initialize the dVdlambda term with the long range contribution */
1004     /* Initialize the dvdl term with the long range contribution */
1005     enerd->term[F_DVDL]            = 0.0;
1006     enerd->term[F_DVDL_COUL]       = 0.0;
1007     enerd->term[F_DVDL_VDW]        = 0.0;
1008     enerd->term[F_DVDL_BONDED]     = 0.0;
1009     enerd->term[F_DVDL_RESTRAINT]  = 0.0;
1010     enerd->term[F_DKDL]            = 0.0;
1011     if (enerd->n_lambda > 0)
1012     {
1013         for (i = 0; i < enerd->n_lambda; i++)
1014         {
1015             enerd->enerpart_lambda[i] = 0.0;
1016         }
1017     }
1018     /* reset foreign energy data - separate function since we also call it elsewhere */
1019     reset_foreign_enerdata(enerd);
1020 }