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