Refactor gmx_group_t to SimulationAtomGroups
[alexxy/gromacs.git] / src / gromacs / mdlib / force.cpp
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,2015,2016,2017,2018,2019, 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 #include "gmxpre.h"
38
39 #include "force.h"
40
41 #include "config.h"
42
43 #include <cassert>
44 #include <cmath>
45 #include <cstring>
46
47 #include "gromacs/domdec/dlbtiming.h"
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/domdec/domdec_struct.h"
50 #include "gromacs/ewald/ewald.h"
51 #include "gromacs/ewald/long_range_correction.h"
52 #include "gromacs/ewald/pme.h"
53 #include "gromacs/gmxlib/network.h"
54 #include "gromacs/gmxlib/nrnb.h"
55 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
56 #include "gromacs/listed_forces/listed_forces.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/math/vecdump.h"
59 #include "gromacs/mdlib/force_flags.h"
60 #include "gromacs/mdlib/forcerec_threading.h"
61 #include "gromacs/mdlib/ns.h"
62 #include "gromacs/mdlib/qmmm.h"
63 #include "gromacs/mdlib/rf_util.h"
64 #include "gromacs/mdlib/wall.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/enerdata.h"
67 #include "gromacs/mdtypes/forceoutput.h"
68 #include "gromacs/mdtypes/forcerec.h"
69 #include "gromacs/mdtypes/inputrec.h"
70 #include "gromacs/mdtypes/md_enums.h"
71 #include "gromacs/pbcutil/ishift.h"
72 #include "gromacs/pbcutil/mshift.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/timing/wallcycle.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/exceptions.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/smalloc.h"
79
80 void ns(FILE                      *fp,
81         t_forcerec                *fr,
82         matrix                     box,
83         const SimulationGroups    *groups,
84         gmx_localtop_t            *top,
85         const t_mdatoms           *md,
86         const t_commrec           *cr,
87         t_nrnb                    *nrnb,
88         gmx_bool                   bFillGrid)
89 {
90     int     nsearch;
91
92
93     if (!fr->ns->nblist_initialized)
94     {
95         init_neighbor_list(fp, fr, md->homenr);
96     }
97
98     nsearch = search_neighbours(fp, fr, box, top, groups, cr, nrnb, md,
99                                 bFillGrid);
100     if (debug)
101     {
102         fprintf(debug, "nsearch = %d\n", nsearch);
103     }
104
105     /* Check whether we have to do dynamic load balancing */
106     /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
107        count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
108        &(top->idef),opts->ngener);
109      */
110     if (fr->ns->dump_nl > 0)
111     {
112         dump_nblist(fp, cr, fr, fr->ns->dump_nl);
113     }
114 }
115
116 static void clearEwaldThreadOutput(ewald_corr_thread_t *ewc_t)
117 {
118     ewc_t->Vcorr_q        = 0;
119     ewc_t->Vcorr_lj       = 0;
120     ewc_t->dvdl[efptCOUL] = 0;
121     ewc_t->dvdl[efptVDW]  = 0;
122     clear_mat(ewc_t->vir_q);
123     clear_mat(ewc_t->vir_lj);
124 }
125
126 static void reduceEwaldThreadOuput(int nthreads, ewald_corr_thread_t *ewc_t)
127 {
128     ewald_corr_thread_t &dest = ewc_t[0];
129
130     for (int t = 1; t < nthreads; t++)
131     {
132         dest.Vcorr_q        += ewc_t[t].Vcorr_q;
133         dest.Vcorr_lj       += ewc_t[t].Vcorr_lj;
134         dest.dvdl[efptCOUL] += ewc_t[t].dvdl[efptCOUL];
135         dest.dvdl[efptVDW]  += ewc_t[t].dvdl[efptVDW];
136         m_add(dest.vir_q,  ewc_t[t].vir_q,  dest.vir_q);
137         m_add(dest.vir_lj, ewc_t[t].vir_lj, dest.vir_lj);
138     }
139 }
140
141 void do_force_lowlevel(t_forcerec                   *fr,
142                        const t_inputrec             *ir,
143                        const t_idef                 *idef,
144                        const t_commrec              *cr,
145                        const gmx_multisim_t         *ms,
146                        t_nrnb                       *nrnb,
147                        gmx_wallcycle_t               wcycle,
148                        const t_mdatoms              *md,
149                        rvec                          x[],
150                        history_t                    *hist,
151                        rvec                         *forceForUseWithShiftForces,
152                        gmx::ForceWithVirial         *forceWithVirial,
153                        gmx_enerdata_t               *enerd,
154                        t_fcdata                     *fcd,
155                        matrix                        box,
156                        t_lambda                     *fepvals,
157                        real                         *lambda,
158                        const t_graph                *graph,
159                        const t_blocka               *excl,
160                        rvec                          mu_tot[],
161                        int                           flags,
162                        const DDBalanceRegionHandler &ddBalanceRegionHandler)
163 {
164     int         i, j;
165     int         donb_flags;
166     int         pme_flags;
167     real        dvdl_dum[efptNR], dvdl_nb[efptNR];
168
169 #if GMX_MPI
170     double  t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
171 #endif
172
173     /* reset free energy components */
174     for (i = 0; i < efptNR; i++)
175     {
176         dvdl_nb[i]  = 0;
177         dvdl_dum[i] = 0;
178     }
179
180     /* do QMMM first if requested */
181     if (fr->bQMMM)
182     {
183         enerd->term[F_EQM] = calculate_QMMM(cr, forceForUseWithShiftForces, fr);
184     }
185
186     /* Call the short range functions all in one go. */
187
188 #if GMX_MPI
189     /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
190 #define TAKETIME FALSE
191     if (TAKETIME)
192     {
193         MPI_Barrier(cr->mpi_comm_mygroup);
194         t0 = MPI_Wtime();
195     }
196 #endif
197
198     if (ir->nwall)
199     {
200         /* foreign lambda component for walls */
201         real dvdl_walls = do_walls(*ir, *fr, box, *md, x,
202                                    forceWithVirial, lambda[efptVDW],
203                                    enerd->grpp.ener[egLJSR], nrnb);
204         enerd->dvdl_lin[efptVDW] += dvdl_walls;
205     }
206
207     /* We only do non-bonded calculation with group scheme here, the verlet
208      * calls are done from do_force_cutsVERLET(). */
209     if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
210     {
211         donb_flags = 0;
212         /* Add short-range interactions */
213         donb_flags |= GMX_NONBONDED_DO_SR;
214
215         /* Currently all group scheme kernels always calculate (shift-)forces */
216         if (flags & GMX_FORCE_FORCES)
217         {
218             donb_flags |= GMX_NONBONDED_DO_FORCE;
219         }
220         if (flags & GMX_FORCE_VIRIAL)
221         {
222             donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
223         }
224         if (flags & GMX_FORCE_ENERGY)
225         {
226             donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
227         }
228
229         wallcycle_sub_start(wcycle, ewcsNONBONDED);
230         do_nonbonded(fr, x, forceForUseWithShiftForces, md, excl,
231                      &enerd->grpp, nrnb,
232                      lambda, dvdl_nb, -1, -1, donb_flags);
233
234         /* If we do foreign lambda and we have soft-core interactions
235          * we have to recalculate the (non-linear) energies contributions.
236          */
237         if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
238         {
239             for (i = 0; i < enerd->n_lambda; i++)
240             {
241                 real lam_i[efptNR];
242
243                 for (j = 0; j < efptNR; j++)
244                 {
245                     lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
246                 }
247                 reset_foreign_enerdata(enerd);
248                 do_nonbonded(fr, x, forceForUseWithShiftForces, md, excl,
249                              &(enerd->foreign_grpp), nrnb,
250                              lam_i, dvdl_dum, -1, -1,
251                              (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
252                 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
253                 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
254             }
255         }
256         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
257     }
258
259 #if GMX_MPI
260     if (TAKETIME)
261     {
262         t1          = MPI_Wtime();
263         fr->t_fnbf += t1-t0;
264     }
265 #endif
266
267     if (fepvals->sc_alpha != 0)
268     {
269         enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
270     }
271     else
272     {
273         enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
274     }
275
276     if (fepvals->sc_alpha != 0)
277
278     /* even though coulomb part is linear, we already added it, beacuse we
279        need to go through the vdw calculation anyway */
280     {
281         enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
282     }
283     else
284     {
285         enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
286     }
287
288     if (debug)
289     {
290         pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
291     }
292
293     /* Shift the coordinates. Must be done before listed forces and PPPM,
294      * but is also necessary for SHAKE and update, therefore it can NOT
295      * go when no listed forces have to be evaluated.
296      *
297      * The shifting and PBC code is deliberately not timed, since with
298      * the Verlet scheme it only takes non-zero time with triclinic
299      * boxes, and even then the time is around a factor of 100 less
300      * than the next smallest counter.
301      */
302
303
304     /* Here sometimes we would not need to shift with NBFonly,
305      * but we do so anyhow for consistency of the returned coordinates.
306      */
307     if (graph)
308     {
309         shift_self(graph, box, x);
310         if (TRICLINIC(box))
311         {
312             inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
313         }
314         else
315         {
316             inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
317         }
318     }
319
320     t_pbc      pbc;
321     const bool useRfWithGroupScheme = (fr->cutoff_scheme == ecutsGROUP) && EEL_RF(fr->ic->eeltype);
322     /* Check whether we need to take into account PBC in the following force tasks:
323      * listed interactions or when correcting for exclusions (for the Group scheme with RF) */
324     if (fr->bMolPBC)
325     {
326         const auto needPbcForListedForces = bool(flags & GMX_FORCE_LISTED) && haveCpuListedForces(*fr, *idef, *fcd);
327         if (needPbcForListedForces || useRfWithGroupScheme)
328         {
329             /* Since all atoms are in the rectangular or triclinic unit-cell,
330              * only single box vector shifts (2 in x) are required.
331              */
332             set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd->nc : nullptr,
333                        TRUE, box);
334         }
335     }
336
337     do_force_listed(wcycle, box, ir->fepvals, cr, ms,
338                     idef, x, hist,
339                     forceForUseWithShiftForces, forceWithVirial,
340                     fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
341                     DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr,
342                     flags);
343
344
345     /* Do long-range electrostatics and/or LJ-PME, including related short-range
346      * corrections.
347      */
348     if (EEL_FULL(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
349     {
350         int  status            = 0;
351         real Vlr_q             = 0, Vlr_lj = 0;
352
353         /* We reduce all virial, dV/dlambda and energy contributions, except
354          * for the reciprocal energies (Vlr_q, Vlr_lj) into the same struct.
355          */
356         ewald_corr_thread_t &ewaldOutput = fr->ewc_t[0];
357         clearEwaldThreadOutput(&ewaldOutput);
358
359         if (EEL_PME_EWALD(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
360         {
361             /* With the Verlet scheme exclusion forces are calculated
362              * in the non-bonded kernel.
363              */
364             /* The TPI molecule does not have exclusions with the rest
365              * of the system and no intra-molecular PME grid
366              * contributions will be calculated in
367              * gmx_pme_calc_energy.
368              */
369             if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
370                 ir->ewald_geometry != eewg3D ||
371                 ir->epsilon_surface != 0)
372             {
373                 int nthreads, t;
374
375                 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
376
377                 if (fr->n_tpi > 0)
378                 {
379                     gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
380                 }
381
382                 nthreads = fr->nthread_ewc;
383 #pragma omp parallel for num_threads(nthreads) schedule(static)
384                 for (t = 0; t < nthreads; t++)
385                 {
386                     try
387                     {
388                         ewald_corr_thread_t &ewc_t = fr->ewc_t[t];
389                         if (t > 0)
390                         {
391                             clearEwaldThreadOutput(&ewc_t);
392                         }
393
394                         /* Threading is only supported with the Verlet cut-off
395                          * scheme and then only single particle forces (no
396                          * exclusion forces) are calculated, so we can store
397                          * the forces in the normal, single forceWithVirial->force_ array.
398                          */
399                         ewald_LRcorrection(md->homenr, cr, nthreads, t, fr, ir,
400                                            md->chargeA, md->chargeB,
401                                            md->sqrt_c6A, md->sqrt_c6B,
402                                            md->sigmaA, md->sigmaB,
403                                            md->sigma3A, md->sigma3B,
404                                            (md->nChargePerturbed != 0) || (md->nTypePerturbed != 0),
405                                            ir->cutoff_scheme != ecutsVERLET,
406                                            excl, x, box, mu_tot,
407                                            ir->ewald_geometry,
408                                            ir->epsilon_surface,
409                                            as_rvec_array(forceWithVirial->force_.data()),
410                                            ewc_t.vir_q, ewc_t.vir_lj,
411                                            &ewc_t.Vcorr_q, &ewc_t.Vcorr_lj,
412                                            lambda[efptCOUL], lambda[efptVDW],
413                                            &ewc_t.dvdl[efptCOUL], &ewc_t.dvdl[efptVDW]);
414                     }
415                     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
416                 }
417                 if (nthreads > 1)
418                 {
419                     reduceEwaldThreadOuput(nthreads, fr->ewc_t);
420                 }
421                 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
422             }
423
424             if (EEL_PME_EWALD(fr->ic->eeltype) && fr->n_tpi == 0)
425             {
426                 /* This is not in a subcounter because it takes a
427                    negligible and constant-sized amount of time */
428                 ewaldOutput.Vcorr_q +=
429                     ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
430                                             &ewaldOutput.dvdl[efptCOUL],
431                                             ewaldOutput.vir_q);
432             }
433
434             if ((EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype)) &&
435                 thisRankHasDuty(cr, DUTY_PME) && (pme_run_mode(fr->pmedata) == PmeRunMode::CPU))
436             {
437                 /* Do reciprocal PME for Coulomb and/or LJ. */
438                 assert(fr->n_tpi >= 0);
439                 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
440                 {
441                     pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
442
443                     if (flags & GMX_FORCE_FORCES)
444                     {
445                         pme_flags |= GMX_PME_CALC_F;
446                     }
447                     if (flags & GMX_FORCE_VIRIAL)
448                     {
449                         pme_flags |= GMX_PME_CALC_ENER_VIR;
450                     }
451                     if (fr->n_tpi > 0)
452                     {
453                         /* We don't calculate f, but we do want the potential */
454                         pme_flags |= GMX_PME_CALC_POT;
455                     }
456
457                     /* With domain decomposition we close the CPU side load
458                      * balancing region here, because PME does global
459                      * communication that acts as a global barrier.
460                      */
461                     ddBalanceRegionHandler.closeAfterForceComputationCpu();
462
463                     wallcycle_start(wcycle, ewcPMEMESH);
464                     status = gmx_pme_do(fr->pmedata,
465                                         0, md->homenr - fr->n_tpi,
466                                         x,
467                                         as_rvec_array(forceWithVirial->force_.data()),
468                                         md->chargeA, md->chargeB,
469                                         md->sqrt_c6A, md->sqrt_c6B,
470                                         md->sigmaA, md->sigmaB,
471                                         box, cr,
472                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
473                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
474                                         nrnb, wcycle,
475                                         ewaldOutput.vir_q, ewaldOutput.vir_lj,
476                                         &Vlr_q, &Vlr_lj,
477                                         lambda[efptCOUL], lambda[efptVDW],
478                                         &ewaldOutput.dvdl[efptCOUL],
479                                         &ewaldOutput.dvdl[efptVDW],
480                                         pme_flags);
481                     wallcycle_stop(wcycle, ewcPMEMESH);
482                     if (status != 0)
483                     {
484                         gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
485                     }
486
487                     /* We should try to do as little computation after
488                      * this as possible, because parallel PME synchronizes
489                      * the nodes, so we want all load imbalance of the
490                      * rest of the force calculation to be before the PME
491                      * call.  DD load balancing is done on the whole time
492                      * of the force call (without PME).
493                      */
494                 }
495                 if (fr->n_tpi > 0)
496                 {
497                     if (EVDW_PME(ir->vdwtype))
498                     {
499
500                         gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
501                     }
502                     /* Determine the PME grid energy of the test molecule
503                      * with the PME grid potential of the other charges.
504                      */
505                     gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
506                                         x + md->homenr - fr->n_tpi,
507                                         md->chargeA + md->homenr - fr->n_tpi,
508                                         &Vlr_q);
509                 }
510             }
511         }
512
513         if (!EEL_PME(fr->ic->eeltype) && EEL_PME_EWALD(fr->ic->eeltype))
514         {
515             Vlr_q = do_ewald(ir, x, as_rvec_array(forceWithVirial->force_.data()),
516                              md->chargeA, md->chargeB,
517                              box, cr, md->homenr,
518                              ewaldOutput.vir_q, fr->ic->ewaldcoeff_q,
519                              lambda[efptCOUL], &ewaldOutput.dvdl[efptCOUL],
520                              fr->ewald_table);
521         }
522
523         /* Note that with separate PME nodes we get the real energies later */
524         // TODO it would be simpler if we just accumulated a single
525         // long-range virial contribution.
526         forceWithVirial->addVirialContribution(ewaldOutput.vir_q);
527         forceWithVirial->addVirialContribution(ewaldOutput.vir_lj);
528         enerd->dvdl_lin[efptCOUL] += ewaldOutput.dvdl[efptCOUL];
529         enerd->dvdl_lin[efptVDW]  += ewaldOutput.dvdl[efptVDW];
530         enerd->term[F_COUL_RECIP]  = Vlr_q + ewaldOutput.Vcorr_q;
531         enerd->term[F_LJ_RECIP]    = Vlr_lj + ewaldOutput.Vcorr_lj;
532
533         if (debug)
534         {
535             fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
536                     Vlr_q, ewaldOutput.Vcorr_q, enerd->term[F_COUL_RECIP]);
537             pr_rvecs(debug, 0, "vir_el_recip after corr", ewaldOutput.vir_q, DIM);
538             pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
539             fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
540                     Vlr_lj, ewaldOutput.Vcorr_lj, enerd->term[F_LJ_RECIP]);
541             pr_rvecs(debug, 0, "vir_lj_recip after corr", ewaldOutput.vir_lj, DIM);
542         }
543     }
544     else
545     {
546         /* Is there a reaction-field exclusion correction needed?
547          * With the Verlet scheme, exclusion forces are calculated
548          * in the non-bonded kernel.
549          */
550         if (useRfWithGroupScheme)
551         {
552             real dvdl_rf_excl      = 0;
553             enerd->term[F_RF_EXCL] =
554                 RF_excl_correction(fr, graph, md, excl, DOMAINDECOMP(cr),
555                                    x, forceForUseWithShiftForces,
556                                    fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
557
558             enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
559         }
560     }
561
562     if (debug)
563     {
564         print_nrnb(debug, nrnb);
565     }
566
567 #if GMX_MPI
568     if (TAKETIME)
569     {
570         t2 = MPI_Wtime();
571         MPI_Barrier(cr->mpi_comm_mygroup);
572         t3          = MPI_Wtime();
573         fr->t_wait += t3-t2;
574         if (fr->timesteps == 11)
575         {
576             char buf[22];
577             fprintf(stderr, "* PP load balancing info: rank %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
578                     cr->nodeid, gmx_step_str(fr->timesteps, buf),
579                     100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
580                     (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
581         }
582         fr->timesteps++;
583     }
584 #endif
585
586     if (debug)
587     {
588         pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
589     }
590
591 }
592
593 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
594 {
595     int i, n2;
596
597     for (i = 0; i < F_NRE; i++)
598     {
599         enerd->term[i]         = 0;
600         enerd->foreign_term[i] = 0;
601     }
602
603
604     for (i = 0; i < efptNR; i++)
605     {
606         enerd->dvdl_lin[i]     = 0;
607         enerd->dvdl_nonlin[i]  = 0;
608     }
609
610     n2 = ngener*ngener;
611     if (debug)
612     {
613         fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
614     }
615     enerd->grpp.nener         = n2;
616     enerd->foreign_grpp.nener = n2;
617     for (i = 0; (i < egNR); i++)
618     {
619         snew(enerd->grpp.ener[i], n2);
620         snew(enerd->foreign_grpp.ener[i], n2);
621     }
622
623     if (n_lambda)
624     {
625         enerd->n_lambda = 1 + n_lambda;
626         snew(enerd->enerpart_lambda, enerd->n_lambda);
627     }
628     else
629     {
630         enerd->n_lambda = 0;
631     }
632 }
633
634 void destroy_enerdata(gmx_enerdata_t *enerd)
635 {
636     int i;
637
638     for (i = 0; (i < egNR); i++)
639     {
640         sfree(enerd->grpp.ener[i]);
641     }
642
643     for (i = 0; (i < egNR); i++)
644     {
645         sfree(enerd->foreign_grpp.ener[i]);
646     }
647
648     if (enerd->n_lambda)
649     {
650         sfree(enerd->enerpart_lambda);
651     }
652 }
653
654 static real sum_v(int n, const real v[])
655 {
656     real t;
657     int  i;
658
659     t = 0.0;
660     for (i = 0; (i < n); i++)
661     {
662         t = t + v[i];
663     }
664
665     return t;
666 }
667
668 void sum_epot(gmx_grppairener_t *grpp, real *epot)
669 {
670     int i;
671
672     /* Accumulate energies */
673     epot[F_COUL_SR]  = sum_v(grpp->nener, grpp->ener[egCOULSR]);
674     epot[F_LJ]       = sum_v(grpp->nener, grpp->ener[egLJSR]);
675     epot[F_LJ14]     = sum_v(grpp->nener, grpp->ener[egLJ14]);
676     epot[F_COUL14]   = sum_v(grpp->nener, grpp->ener[egCOUL14]);
677
678 /* lattice part of LR doesnt belong to any group
679  * and has been added earlier
680  */
681     epot[F_BHAM]     = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
682
683     epot[F_EPOT] = 0;
684     for (i = 0; (i < F_EPOT); i++)
685     {
686         if (i != F_DISRESVIOL && i != F_ORIRESDEV)
687         {
688             epot[F_EPOT] += epot[i];
689         }
690     }
691 }
692
693 void sum_dhdl(gmx_enerdata_t *enerd, gmx::ArrayRef<const real> lambda, t_lambda *fepvals)
694 {
695     int    index;
696
697     enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW];  /* include dispersion correction */
698     enerd->term[F_DVDL]       = 0.0;
699     for (int i = 0; i < efptNR; i++)
700     {
701         if (fepvals->separate_dvdl[i])
702         {
703             /* could this be done more readably/compactly? */
704             switch (i)
705             {
706                 case (efptMASS):
707                     index = F_DKDL;
708                     break;
709                 case (efptCOUL):
710                     index = F_DVDL_COUL;
711                     break;
712                 case (efptVDW):
713                     index = F_DVDL_VDW;
714                     break;
715                 case (efptBONDED):
716                     index = F_DVDL_BONDED;
717                     break;
718                 case (efptRESTRAINT):
719                     index = F_DVDL_RESTRAINT;
720                     break;
721                 default:
722                     index = F_DVDL;
723                     break;
724             }
725             enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
726             if (debug)
727             {
728                 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
729                         efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
730             }
731         }
732         else
733         {
734             enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
735             if (debug)
736             {
737                 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
738                         efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
739             }
740         }
741     }
742
743     if (fepvals->separate_dvdl[efptBONDED])
744     {
745         enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
746     }
747     else
748     {
749         enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
750     }
751
752     for (int i = 0; i < fepvals->n_lambda; i++)
753     {
754         /* note we are iterating over fepvals here!
755            For the current lam, dlam = 0 automatically,
756            so we don't need to add anything to the
757            enerd->enerpart_lambda[0] */
758
759         /* we don't need to worry about dvdl_lin contributions to dE at
760            current lambda, because the contributions to the current
761            lambda are automatically zeroed */
762
763         double &enerpart_lambda = enerd->enerpart_lambda[i + 1];
764
765         for (gmx::index j = 0; j < lambda.ssize(); j++)
766         {
767             /* Note that this loop is over all dhdl components, not just the separated ones */
768             const double dlam  = fepvals->all_lambda[j][i] - lambda[j];
769
770             enerpart_lambda   += dlam*enerd->dvdl_lin[j];
771
772             /* Constraints can not be evaluated at foreign lambdas, so we add
773              * a linear extrapolation. This is an approximation, but usually
774              * quite accurate since constraints change little between lambdas.
775              */
776             if ((j == efptBONDED && fepvals->separate_dvdl[efptBONDED]) ||
777                 (j == efptFEP && !fepvals->separate_dvdl[efptBONDED]))
778             {
779                 enerpart_lambda += dlam*enerd->term[F_DVDL_CONSTR];
780             }
781
782             if (j == efptMASS && !fepvals->separate_dvdl[j])
783             {
784                 enerpart_lambda += dlam*enerd->term[F_DKDL];
785             }
786
787             if (debug)
788             {
789                 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
790                         fepvals->all_lambda[j][i], efpt_names[j],
791                         enerpart_lambda - enerd->enerpart_lambda[0],
792                         dlam, enerd->dvdl_lin[j]);
793             }
794         }
795     }
796
797     /* The constrain contribution is now included in other terms, so clear it */
798     enerd->term[F_DVDL_CONSTR] = 0;
799 }
800
801
802 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
803 {
804     int  i, j;
805
806     /* First reset all foreign energy components.  Foreign energies always called on
807        neighbor search steps */
808     for (i = 0; (i < egNR); i++)
809     {
810         for (j = 0; (j < enerd->grpp.nener); j++)
811         {
812             enerd->foreign_grpp.ener[i][j] = 0.0;
813         }
814     }
815
816     /* potential energy components */
817     for (i = 0; (i <= F_EPOT); i++)
818     {
819         enerd->foreign_term[i] = 0.0;
820     }
821 }
822
823 void reset_enerdata(gmx_enerdata_t *enerd)
824 {
825     int      i, j;
826
827     /* First reset all energy components. */
828     for (i = 0; (i < egNR); i++)
829     {
830         for (j = 0; (j < enerd->grpp.nener); j++)
831         {
832             enerd->grpp.ener[i][j] = 0.0;
833         }
834     }
835     for (i = 0; i < efptNR; i++)
836     {
837         enerd->dvdl_lin[i]    = 0.0;
838         enerd->dvdl_nonlin[i] = 0.0;
839     }
840
841     /* Normal potential energy components */
842     for (i = 0; (i <= F_EPOT); i++)
843     {
844         enerd->term[i] = 0.0;
845     }
846     enerd->term[F_DVDL]            = 0.0;
847     enerd->term[F_DVDL_COUL]       = 0.0;
848     enerd->term[F_DVDL_VDW]        = 0.0;
849     enerd->term[F_DVDL_BONDED]     = 0.0;
850     enerd->term[F_DVDL_RESTRAINT]  = 0.0;
851     enerd->term[F_DKDL]            = 0.0;
852     if (enerd->n_lambda > 0)
853     {
854         for (i = 0; i < enerd->n_lambda; i++)
855         {
856             enerd->enerpart_lambda[i] = 0.0;
857         }
858     }
859     /* reset foreign energy data - separate function since we also call it elsewhere */
860     reset_foreign_enerdata(enerd);
861 }