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