Merge branch release-2016
[alexxy/gromacs.git] / src / gromacs / mdlib / sim_util.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, 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 "sim_util.h"
40
41 #include "config.h"
42
43 #include <assert.h>
44 #include <math.h>
45 #include <stdio.h>
46 #include <string.h>
47
48 #include <cstdint>
49
50 #include <array>
51
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/essentialdynamics/edsam.h"
55 #include "gromacs/ewald/pme.h"
56 #include "gromacs/gmxlib/chargegroup.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/gmxlib/nrnb.h"
59 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
60 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
61 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
62 #include "gromacs/imd/imd.h"
63 #include "gromacs/listed-forces/bonded.h"
64 #include "gromacs/listed-forces/disre.h"
65 #include "gromacs/listed-forces/orires.h"
66 #include "gromacs/math/functions.h"
67 #include "gromacs/math/units.h"
68 #include "gromacs/math/vec.h"
69 #include "gromacs/math/vecdump.h"
70 #include "gromacs/mdlib/calcmu.h"
71 #include "gromacs/mdlib/constr.h"
72 #include "gromacs/mdlib/force.h"
73 #include "gromacs/mdlib/forcerec.h"
74 #include "gromacs/mdlib/genborn.h"
75 #include "gromacs/mdlib/gmx_omp_nthreads.h"
76 #include "gromacs/mdlib/mdrun.h"
77 #include "gromacs/mdlib/nb_verlet.h"
78 #include "gromacs/mdlib/nbnxn_atomdata.h"
79 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
80 #include "gromacs/mdlib/nbnxn_grid.h"
81 #include "gromacs/mdlib/nbnxn_search.h"
82 #include "gromacs/mdlib/qmmm.h"
83 #include "gromacs/mdlib/update.h"
84 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
85 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref.h"
86 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
87 #include "gromacs/mdlib/nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
88 #include "gromacs/mdtypes/commrec.h"
89 #include "gromacs/mdtypes/inputrec.h"
90 #include "gromacs/mdtypes/md_enums.h"
91 #include "gromacs/mdtypes/state.h"
92 #include "gromacs/pbcutil/ishift.h"
93 #include "gromacs/pbcutil/mshift.h"
94 #include "gromacs/pbcutil/pbc.h"
95 #include "gromacs/pulling/pull.h"
96 #include "gromacs/pulling/pull_rotation.h"
97 #include "gromacs/timing/cyclecounter.h"
98 #include "gromacs/timing/gpu_timing.h"
99 #include "gromacs/timing/wallcycle.h"
100 #include "gromacs/timing/wallcyclereporting.h"
101 #include "gromacs/timing/walltime_accounting.h"
102 #include "gromacs/utility/arrayref.h"
103 #include "gromacs/utility/basedefinitions.h"
104 #include "gromacs/utility/cstringutil.h"
105 #include "gromacs/utility/exceptions.h"
106 #include "gromacs/utility/fatalerror.h"
107 #include "gromacs/utility/gmxassert.h"
108 #include "gromacs/utility/gmxmpi.h"
109 #include "gromacs/utility/logger.h"
110 #include "gromacs/utility/pleasecite.h"
111 #include "gromacs/utility/smalloc.h"
112 #include "gromacs/utility/sysinfo.h"
113
114 #include "nbnxn_gpu.h"
115
116 void print_time(FILE                     *out,
117                 gmx_walltime_accounting_t walltime_accounting,
118                 gmx_int64_t               step,
119                 t_inputrec               *ir,
120                 t_commrec gmx_unused     *cr)
121 {
122     time_t finish;
123     char   timebuf[STRLEN];
124     double dt, elapsed_seconds, time_per_step;
125     char   buf[48];
126
127 #if !GMX_THREAD_MPI
128     if (!PAR(cr))
129 #endif
130     {
131         fprintf(out, "\r");
132     }
133     fprintf(out, "step %s", gmx_step_str(step, buf));
134     fflush(out);
135
136     if ((step >= ir->nstlist))
137     {
138         double seconds_since_epoch = gmx_gettime();
139         elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
140         time_per_step   = elapsed_seconds/(step - ir->init_step + 1);
141         dt              = (ir->nsteps + ir->init_step - step) * time_per_step;
142
143         if (ir->nsteps >= 0)
144         {
145             if (dt >= 300)
146             {
147                 finish = (time_t) (seconds_since_epoch + dt);
148                 gmx_ctime_r(&finish, timebuf, STRLEN);
149                 sprintf(buf, "%s", timebuf);
150                 buf[strlen(buf)-1] = '\0';
151                 fprintf(out, ", will finish %s", buf);
152             }
153             else
154             {
155                 fprintf(out, ", remaining wall clock time: %5d s          ", (int)dt);
156             }
157         }
158         else
159         {
160             fprintf(out, " performance: %.1f ns/day    ",
161                     ir->delta_t/1000*24*60*60/time_per_step);
162         }
163     }
164 #if !GMX_THREAD_MPI
165     if (PAR(cr))
166     {
167         fprintf(out, "\n");
168     }
169 #endif
170
171     fflush(out);
172 }
173
174 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
175                          double the_time)
176 {
177     char   time_string[STRLEN];
178
179     if (!fplog)
180     {
181         return;
182     }
183
184     {
185         int    i;
186         char   timebuf[STRLEN];
187         time_t temp_time = (time_t) the_time;
188
189         gmx_ctime_r(&temp_time, timebuf, STRLEN);
190         for (i = 0; timebuf[i] >= ' '; i++)
191         {
192             time_string[i] = timebuf[i];
193         }
194         time_string[i] = '\0';
195     }
196
197     fprintf(fplog, "%s on rank %d %s\n", title, nodeid, time_string);
198 }
199
200 void print_start(FILE *fplog, t_commrec *cr,
201                  gmx_walltime_accounting_t walltime_accounting,
202                  const char *name)
203 {
204     char buf[STRLEN];
205
206     sprintf(buf, "Started %s", name);
207     print_date_and_time(fplog, cr->nodeid, buf,
208                         walltime_accounting_get_start_time_stamp(walltime_accounting));
209 }
210
211 static void sum_forces(rvec f[], const PaddedRVecVector *forceToAdd)
212 {
213     /* TODO: remove this - 1 when padding is properly implemented */
214     int         end  = forceToAdd->size() - 1;
215     const rvec *fAdd = as_rvec_array(forceToAdd->data());
216
217     // cppcheck-suppress unreadVariable
218     int gmx_unused nt = gmx_omp_nthreads_get(emntDefault);
219 #pragma omp parallel for num_threads(nt) schedule(static)
220     for (int i = 0; i < end; i++)
221     {
222         rvec_inc(f[i], fAdd[i]);
223     }
224 }
225
226 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
227                         tensor vir_part, t_graph *graph, matrix box,
228                         t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
229 {
230     int    i;
231
232     /* The short-range virial from surrounding boxes */
233     clear_mat(vir_part);
234     calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
235     inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
236
237     /* Calculate partial virial, for local atoms only, based on short range.
238      * Total virial is computed in global_stat, called from do_md
239      */
240     f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
241     inc_nrnb(nrnb, eNR_VIRIAL, homenr);
242
243     /* Add position restraint contribution */
244     for (i = 0; i < DIM; i++)
245     {
246         vir_part[i][i] += fr->vir_diag_posres[i];
247     }
248
249     /* Add wall contribution */
250     for (i = 0; i < DIM; i++)
251     {
252         vir_part[i][ZZ] += fr->vir_wall_z[i];
253     }
254
255     if (debug)
256     {
257         pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
258     }
259 }
260
261 static void pull_potential_wrapper(t_commrec *cr,
262                                    t_inputrec *ir,
263                                    matrix box, rvec x[],
264                                    rvec f[],
265                                    tensor vir_force,
266                                    t_mdatoms *mdatoms,
267                                    gmx_enerdata_t *enerd,
268                                    real *lambda,
269                                    double t,
270                                    gmx_wallcycle_t wcycle)
271 {
272     t_pbc  pbc;
273     real   dvdl;
274
275     /* Calculate the center of mass forces, this requires communication,
276      * which is why pull_potential is called close to other communication.
277      * The virial contribution is calculated directly,
278      * which is why we call pull_potential after calc_virial.
279      */
280     wallcycle_start(wcycle, ewcPULLPOT);
281     set_pbc(&pbc, ir->ePBC, box);
282     dvdl                     = 0;
283     enerd->term[F_COM_PULL] +=
284         pull_potential(ir->pull_work, mdatoms, &pbc,
285                        cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
286     enerd->dvdl_lin[efptRESTRAINT] += dvdl;
287     wallcycle_stop(wcycle, ewcPULLPOT);
288 }
289
290 static void pme_receive_force_ener(t_commrec      *cr,
291                                    gmx_wallcycle_t wcycle,
292                                    gmx_enerdata_t *enerd,
293                                    t_forcerec     *fr)
294 {
295     real   e_q, e_lj, dvdl_q, dvdl_lj;
296     float  cycles_ppdpme, cycles_seppme;
297
298     cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
299     dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
300
301     /* In case of node-splitting, the PP nodes receive the long-range
302      * forces, virial and energy from the PME nodes here.
303      */
304     wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
305     dvdl_q  = 0;
306     dvdl_lj = 0;
307     gmx_pme_receive_f(cr, as_rvec_array(fr->f_novirsum->data()), fr->vir_el_recip, &e_q,
308                       fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
309                       &cycles_seppme);
310     enerd->term[F_COUL_RECIP] += e_q;
311     enerd->term[F_LJ_RECIP]   += e_lj;
312     enerd->dvdl_lin[efptCOUL] += dvdl_q;
313     enerd->dvdl_lin[efptVDW]  += dvdl_lj;
314
315     if (wcycle)
316     {
317         dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
318     }
319     wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
320 }
321
322 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
323                                gmx_int64_t step, real forceTolerance,
324                                const rvec *x, const rvec *f)
325 {
326     real           force2Tolerance = gmx::square(forceTolerance);
327     std::uintmax_t numNonFinite    = 0;
328     for (int i = 0; i < md->homenr; i++)
329     {
330         real force2    = norm2(f[i]);
331         bool nonFinite = !std::isfinite(force2);
332         if (force2 >= force2Tolerance || nonFinite)
333         {
334             fprintf(fp, "step %" GMX_PRId64 " atom %6d  x %8.3f %8.3f %8.3f  force %12.5e\n",
335                     step,
336                     ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], std::sqrt(force2));
337         }
338         if (nonFinite)
339         {
340             numNonFinite++;
341         }
342     }
343     if (numNonFinite > 0)
344     {
345         /* Note that with MPI this fatal call on one rank might interrupt
346          * the printing on other ranks. But we can only avoid that with
347          * an expensive MPI barrier that we would need at each step.
348          */
349         gmx_fatal(FARGS, "At step %" GMX_PRId64 " detected non-finite forces on %ju atoms", step, numNonFinite);
350     }
351 }
352
353 static void post_process_forces(t_commrec *cr,
354                                 gmx_int64_t step,
355                                 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
356                                 gmx_localtop_t *top,
357                                 matrix box, rvec x[],
358                                 rvec f[],
359                                 tensor vir_force,
360                                 t_mdatoms *mdatoms,
361                                 t_graph *graph,
362                                 t_forcerec *fr, gmx_vsite_t *vsite,
363                                 int flags)
364 {
365     if (fr->bF_NoVirSum)
366     {
367         if (vsite)
368         {
369             /* Spread the mesh force on virtual sites to the other particles...
370              * This is parallellized. MPI communication is performed
371              * if the constructing atoms aren't local.
372              */
373             wallcycle_start(wcycle, ewcVSITESPREAD);
374             spread_vsite_f(vsite, x, as_rvec_array(fr->f_novirsum->data()), nullptr,
375                            (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
376                            nrnb,
377                            &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
378             wallcycle_stop(wcycle, ewcVSITESPREAD);
379         }
380         if (flags & GMX_FORCE_VIRIAL)
381         {
382             /* Now add the forces, this is local */
383             sum_forces(f, fr->f_novirsum);
384
385             if (EEL_FULL(fr->eeltype))
386             {
387                 /* Add the mesh contribution to the virial */
388                 m_add(vir_force, fr->vir_el_recip, vir_force);
389             }
390             if (EVDW_PME(fr->vdwtype))
391             {
392                 /* Add the mesh contribution to the virial */
393                 m_add(vir_force, fr->vir_lj_recip, vir_force);
394             }
395             if (debug)
396             {
397                 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
398             }
399         }
400     }
401
402     if (fr->print_force >= 0)
403     {
404         print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
405     }
406 }
407
408 static void do_nb_verlet(t_forcerec *fr,
409                          interaction_const_t *ic,
410                          gmx_enerdata_t *enerd,
411                          int flags, int ilocality,
412                          int clearF,
413                          t_nrnb *nrnb,
414                          gmx_wallcycle_t wcycle)
415 {
416     int                        enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
417     nonbonded_verlet_group_t  *nbvg;
418     gmx_bool                   bUsingGpuKernels;
419
420     if (!(flags & GMX_FORCE_NONBONDED))
421     {
422         /* skip non-bonded calculation */
423         return;
424     }
425
426     nbvg = &fr->nbv->grp[ilocality];
427
428     /* GPU kernel launch overhead is already timed separately */
429     if (fr->cutoff_scheme != ecutsVERLET)
430     {
431         gmx_incons("Invalid cut-off scheme passed!");
432     }
433
434     bUsingGpuKernels = (nbvg->kernel_type == nbnxnk8x8x8_GPU);
435
436     if (!bUsingGpuKernels)
437     {
438         wallcycle_sub_start(wcycle, ewcsNONBONDED);
439     }
440     switch (nbvg->kernel_type)
441     {
442         case nbnxnk4x4_PlainC:
443             nbnxn_kernel_ref(&nbvg->nbl_lists,
444                              nbvg->nbat, ic,
445                              fr->shift_vec,
446                              flags,
447                              clearF,
448                              fr->fshift[0],
449                              enerd->grpp.ener[egCOULSR],
450                              fr->bBHAM ?
451                              enerd->grpp.ener[egBHAMSR] :
452                              enerd->grpp.ener[egLJSR]);
453             break;
454
455         case nbnxnk4xN_SIMD_4xN:
456             nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
457                                   nbvg->nbat, ic,
458                                   nbvg->ewald_excl,
459                                   fr->shift_vec,
460                                   flags,
461                                   clearF,
462                                   fr->fshift[0],
463                                   enerd->grpp.ener[egCOULSR],
464                                   fr->bBHAM ?
465                                   enerd->grpp.ener[egBHAMSR] :
466                                   enerd->grpp.ener[egLJSR]);
467             break;
468         case nbnxnk4xN_SIMD_2xNN:
469             nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
470                                    nbvg->nbat, ic,
471                                    nbvg->ewald_excl,
472                                    fr->shift_vec,
473                                    flags,
474                                    clearF,
475                                    fr->fshift[0],
476                                    enerd->grpp.ener[egCOULSR],
477                                    fr->bBHAM ?
478                                    enerd->grpp.ener[egBHAMSR] :
479                                    enerd->grpp.ener[egLJSR]);
480             break;
481
482         case nbnxnk8x8x8_GPU:
483             nbnxn_gpu_launch_kernel(fr->nbv->gpu_nbv, nbvg->nbat, flags, ilocality);
484             break;
485
486         case nbnxnk8x8x8_PlainC:
487             nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
488                                  nbvg->nbat, ic,
489                                  fr->shift_vec,
490                                  flags,
491                                  clearF,
492                                  nbvg->nbat->out[0].f,
493                                  fr->fshift[0],
494                                  enerd->grpp.ener[egCOULSR],
495                                  fr->bBHAM ?
496                                  enerd->grpp.ener[egBHAMSR] :
497                                  enerd->grpp.ener[egLJSR]);
498             break;
499
500         default:
501             gmx_incons("Invalid nonbonded kernel type passed!");
502
503     }
504     if (!bUsingGpuKernels)
505     {
506         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
507     }
508
509     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
510     {
511         enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
512     }
513     else if ((!bUsingGpuKernels && nbvg->ewald_excl == ewaldexclAnalytical) ||
514              (bUsingGpuKernels && nbnxn_gpu_is_kernel_ewald_analytical(fr->nbv->gpu_nbv)))
515     {
516         enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
517     }
518     else
519     {
520         enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
521     }
522     enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
523     if (flags & GMX_FORCE_ENERGY)
524     {
525         /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
526         enr_nbnxn_kernel_ljc += 1;
527         enr_nbnxn_kernel_lj  += 1;
528     }
529
530     inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
531              nbvg->nbl_lists.natpair_ljq);
532     inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
533              nbvg->nbl_lists.natpair_lj);
534     /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
535     inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
536              nbvg->nbl_lists.natpair_q);
537
538     if (ic->vdw_modifier == eintmodFORCESWITCH)
539     {
540         /* We add up the switch cost separately */
541         inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
542                  nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
543     }
544     if (ic->vdw_modifier == eintmodPOTSWITCH)
545     {
546         /* We add up the switch cost separately */
547         inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
548                  nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
549     }
550     if (ic->vdwtype == evdwPME)
551     {
552         /* We add up the LJ Ewald cost separately */
553         inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
554                  nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
555     }
556 }
557
558 static void do_nb_verlet_fep(nbnxn_pairlist_set_t *nbl_lists,
559                              t_forcerec           *fr,
560                              rvec                  x[],
561                              rvec                  f[],
562                              t_mdatoms            *mdatoms,
563                              t_lambda             *fepvals,
564                              real                 *lambda,
565                              gmx_enerdata_t       *enerd,
566                              int                   flags,
567                              t_nrnb               *nrnb,
568                              gmx_wallcycle_t       wcycle)
569 {
570     int              donb_flags;
571     nb_kernel_data_t kernel_data;
572     real             lam_i[efptNR];
573     real             dvdl_nb[efptNR];
574     int              th;
575     int              i, j;
576
577     donb_flags = 0;
578     /* Add short-range interactions */
579     donb_flags |= GMX_NONBONDED_DO_SR;
580
581     /* Currently all group scheme kernels always calculate (shift-)forces */
582     if (flags & GMX_FORCE_FORCES)
583     {
584         donb_flags |= GMX_NONBONDED_DO_FORCE;
585     }
586     if (flags & GMX_FORCE_VIRIAL)
587     {
588         donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
589     }
590     if (flags & GMX_FORCE_ENERGY)
591     {
592         donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
593     }
594
595     kernel_data.flags  = donb_flags;
596     kernel_data.lambda = lambda;
597     kernel_data.dvdl   = dvdl_nb;
598
599     kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
600     kernel_data.energygrp_vdw  = enerd->grpp.ener[egLJSR];
601
602     /* reset free energy components */
603     for (i = 0; i < efptNR; i++)
604     {
605         dvdl_nb[i]  = 0;
606     }
607
608     assert(gmx_omp_nthreads_get(emntNonbonded) == nbl_lists->nnbl);
609
610     wallcycle_sub_start(wcycle, ewcsNONBONDED);
611 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
612     for (th = 0; th < nbl_lists->nnbl; th++)
613     {
614         try
615         {
616             gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
617                                       x, f, fr, mdatoms, &kernel_data, nrnb);
618         }
619         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
620     }
621
622     if (fepvals->sc_alpha != 0)
623     {
624         enerd->dvdl_nonlin[efptVDW]  += dvdl_nb[efptVDW];
625         enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
626     }
627     else
628     {
629         enerd->dvdl_lin[efptVDW]  += dvdl_nb[efptVDW];
630         enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
631     }
632
633     /* If we do foreign lambda and we have soft-core interactions
634      * we have to recalculate the (non-linear) energies contributions.
635      */
636     if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
637     {
638         kernel_data.flags          = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
639         kernel_data.lambda         = lam_i;
640         kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
641         kernel_data.energygrp_vdw  = enerd->foreign_grpp.ener[egLJSR];
642         /* Note that we add to kernel_data.dvdl, but ignore the result */
643
644         for (i = 0; i < enerd->n_lambda; i++)
645         {
646             for (j = 0; j < efptNR; j++)
647             {
648                 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
649             }
650             reset_foreign_enerdata(enerd);
651 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
652             for (th = 0; th < nbl_lists->nnbl; th++)
653             {
654                 try
655                 {
656                     gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
657                                               x, f, fr, mdatoms, &kernel_data, nrnb);
658                 }
659                 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
660             }
661
662             sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
663             enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
664         }
665     }
666
667     wallcycle_sub_stop(wcycle, ewcsNONBONDED);
668 }
669
670 gmx_bool use_GPU(const nonbonded_verlet_t *nbv)
671 {
672     return nbv != nullptr && nbv->bUseGPU;
673 }
674
675 static gmx_inline void clear_rvecs_omp(int n, rvec v[])
676 {
677     int nth = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, n);
678
679     /* Note that we would like to avoid this conditional by putting it
680      * into the omp pragma instead, but then we still take the full
681      * omp parallel for overhead (at least with gcc5).
682      */
683     if (nth == 1)
684     {
685         for (int i = 0; i < n; i++)
686         {
687             clear_rvec(v[i]);
688         }
689     }
690     else
691     {
692 #pragma omp parallel for num_threads(nth) schedule(static)
693         for (int i = 0; i < n; i++)
694         {
695             clear_rvec(v[i]);
696         }
697     }
698 }
699
700 /*! \brief  This routine checks if the potential energy is finite.
701  *
702  * Note that passing this check does not guarantee finite forces,
703  * since those use slightly different arithmetics. But in most cases
704  * there is just a narrow coordinate range where forces are not finite
705  * and energies are finite.
706  *
707  * \param[in] enerd  The energy data; the non-bonded group energies need to be added in here before calling this routine
708  */
709 static void checkPotentialEnergyValidity(const gmx_enerdata_t *enerd)
710 {
711     if (!std::isfinite(enerd->term[F_EPOT]))
712     {
713         gmx_fatal(FARGS, "The total potential energy is %g, which is not finite. The LJ and electrostatic contributions to the energy are %g and %g, respectively. A non-finite potential energy can be caused by overlapping interactions in bonded interactions or very large or NaN coordinate values. Usually this is caused by a badly or non-equilibrated initial configuration or incorrect interactions or parameters in the topology.",
714                   enerd->term[F_EPOT],
715                   enerd->term[F_LJ],
716                   enerd->term[F_COUL_SR]);
717     }
718 }
719
720 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
721                          t_inputrec *inputrec,
722                          gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
723                          gmx_localtop_t *top,
724                          gmx_groups_t gmx_unused *groups,
725                          matrix box, rvec x[], history_t *hist,
726                          PaddedRVecVector *force,
727                          tensor vir_force,
728                          t_mdatoms *mdatoms,
729                          gmx_enerdata_t *enerd, t_fcdata *fcd,
730                          real *lambda, t_graph *graph,
731                          t_forcerec *fr, interaction_const_t *ic,
732                          gmx_vsite_t *vsite, rvec mu_tot,
733                          double t, gmx_edsam_t ed,
734                          gmx_bool bBornRadii,
735                          int flags)
736 {
737     int                 cg1, i, j;
738     double              mu[2*DIM];
739     gmx_bool            bStateChanged, bNS, bFillGrid, bCalcCGCM;
740     gmx_bool            bDoForces, bUseGPU, bUseOrEmulGPU;
741     gmx_bool            bDiffKernels = FALSE;
742     rvec                vzero, box_diag;
743     float               cycles_pme, cycles_force, cycles_wait_gpu;
744     /* TODO To avoid loss of precision, float can't be used for a
745      * cycle count. Build an object that can do this right and perhaps
746      * also be used by gmx_wallcycle_t */
747     gmx_cycles_t        cycleCountBeforeLocalWorkCompletes = 0;
748     nonbonded_verlet_t *nbv;
749
750     cycles_force    = 0;
751     cycles_wait_gpu = 0;
752     nbv             = fr->nbv;
753
754     const int start  = 0;
755     const int homenr = mdatoms->homenr;
756
757     clear_mat(vir_force);
758
759     if (DOMAINDECOMP(cr))
760     {
761         cg1 = cr->dd->ncg_tot;
762     }
763     else
764     {
765         cg1 = top->cgs.nr;
766     }
767     if (fr->n_tpi > 0)
768     {
769         cg1--;
770     }
771
772     bStateChanged = (flags & GMX_FORCE_STATECHANGED);
773     bNS           = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
774     bFillGrid     = (bNS && bStateChanged);
775     bCalcCGCM     = (bFillGrid && !DOMAINDECOMP(cr));
776     bDoForces     = (flags & GMX_FORCE_FORCES);
777     bUseGPU       = fr->nbv->bUseGPU;
778     bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
779
780     if (bStateChanged)
781     {
782         update_forcerec(fr, box);
783
784         if (inputrecNeedMutot(inputrec))
785         {
786             /* Calculate total (local) dipole moment in a temporary common array.
787              * This makes it possible to sum them over nodes faster.
788              */
789             calc_mu(start, homenr,
790                     x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
791                     mu, mu+DIM);
792         }
793     }
794
795     if (fr->ePBC != epbcNONE)
796     {
797         /* Compute shift vectors every step,
798          * because of pressure coupling or box deformation!
799          */
800         if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
801         {
802             calc_shifts(box, fr->shift_vec);
803         }
804
805         if (bCalcCGCM)
806         {
807             put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
808             inc_nrnb(nrnb, eNR_SHIFTX, homenr);
809         }
810         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
811         {
812             unshift_self(graph, box, x);
813         }
814     }
815
816     nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
817                                  fr->shift_vec, nbv->grp[0].nbat);
818
819 #if GMX_MPI
820     if (!(cr->duty & DUTY_PME))
821     {
822         gmx_bool bBS;
823         matrix   boxs;
824
825         /* Send particle coordinates to the pme nodes.
826          * Since this is only implemented for domain decomposition
827          * and domain decomposition does not use the graph,
828          * we do not need to worry about shifting.
829          */
830
831         wallcycle_start(wcycle, ewcPP_PMESENDX);
832
833         bBS = (inputrec->nwall == 2);
834         if (bBS)
835         {
836             copy_mat(box, boxs);
837             svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
838         }
839
840         gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
841                                  lambda[efptCOUL], lambda[efptVDW],
842                                  (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
843                                  step);
844
845         wallcycle_stop(wcycle, ewcPP_PMESENDX);
846     }
847 #endif /* GMX_MPI */
848
849     /* do gridding for pair search */
850     if (bNS)
851     {
852         if (graph && bStateChanged)
853         {
854             /* Calculate intramolecular shift vectors to make molecules whole */
855             mk_mshift(fplog, graph, fr->ePBC, box, x);
856         }
857
858         clear_rvec(vzero);
859         box_diag[XX] = box[XX][XX];
860         box_diag[YY] = box[YY][YY];
861         box_diag[ZZ] = box[ZZ][ZZ];
862
863         wallcycle_start(wcycle, ewcNS);
864         if (!fr->bDomDec)
865         {
866             wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
867             nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
868                               0, vzero, box_diag,
869                               0, mdatoms->homenr, -1, fr->cginfo, x,
870                               0, nullptr,
871                               nbv->grp[eintLocal].kernel_type,
872                               nbv->grp[eintLocal].nbat);
873             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
874         }
875         else
876         {
877             wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
878             nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
879                                        fr->cginfo, x,
880                                        nbv->grp[eintNonlocal].kernel_type,
881                                        nbv->grp[eintNonlocal].nbat);
882             wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
883         }
884
885         if (nbv->ngrp == 1 ||
886             nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
887         {
888             nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
889                                nbv->nbs, mdatoms, fr->cginfo);
890         }
891         else
892         {
893             nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
894                                nbv->nbs, mdatoms, fr->cginfo);
895             nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
896                                nbv->nbs, mdatoms, fr->cginfo);
897         }
898         wallcycle_stop(wcycle, ewcNS);
899     }
900
901     /* initialize the GPU atom data and copy shift vector */
902     if (bUseGPU)
903     {
904         if (bNS)
905         {
906             wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
907             nbnxn_gpu_init_atomdata(nbv->gpu_nbv, nbv->grp[eintLocal].nbat);
908             wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
909         }
910
911         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
912         nbnxn_gpu_upload_shiftvec(nbv->gpu_nbv, nbv->grp[eintLocal].nbat);
913         wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
914     }
915
916     /* do local pair search */
917     if (bNS)
918     {
919         wallcycle_start_nocount(wcycle, ewcNS);
920         wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
921         nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
922                             &top->excls,
923                             ic->rlist,
924                             nbv->min_ci_balanced,
925                             &nbv->grp[eintLocal].nbl_lists,
926                             eintLocal,
927                             nbv->grp[eintLocal].kernel_type,
928                             nrnb);
929         wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
930
931         if (bUseGPU)
932         {
933             /* initialize local pair-list on the GPU */
934             nbnxn_gpu_init_pairlist(nbv->gpu_nbv,
935                                     nbv->grp[eintLocal].nbl_lists.nbl[0],
936                                     eintLocal);
937         }
938         wallcycle_stop(wcycle, ewcNS);
939     }
940     else
941     {
942         wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
943         wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
944         nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
945                                         nbv->grp[eintLocal].nbat);
946         wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
947         wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
948     }
949
950     if (bUseGPU)
951     {
952         wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
953         /* launch local nonbonded F on GPU */
954         do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
955                      nrnb, wcycle);
956         wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
957     }
958
959     /* Communicate coordinates and sum dipole if necessary +
960        do non-local pair search */
961     if (DOMAINDECOMP(cr))
962     {
963         bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
964                         nbv->grp[eintLocal].kernel_type);
965
966         if (bDiffKernels)
967         {
968             /* With GPU+CPU non-bonded calculations we need to copy
969              * the local coordinates to the non-local nbat struct
970              * (in CPU format) as the non-local kernel call also
971              * calculates the local - non-local interactions.
972              */
973             wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
974             wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
975             nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
976                                             nbv->grp[eintNonlocal].nbat);
977             wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
978             wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
979         }
980
981         if (bNS)
982         {
983             wallcycle_start_nocount(wcycle, ewcNS);
984             wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
985
986             if (bDiffKernels)
987             {
988                 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
989             }
990
991             nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
992                                 &top->excls,
993                                 ic->rlist,
994                                 nbv->min_ci_balanced,
995                                 &nbv->grp[eintNonlocal].nbl_lists,
996                                 eintNonlocal,
997                                 nbv->grp[eintNonlocal].kernel_type,
998                                 nrnb);
999
1000             wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1001
1002             if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_GPU)
1003             {
1004                 /* initialize non-local pair-list on the GPU */
1005                 nbnxn_gpu_init_pairlist(nbv->gpu_nbv,
1006                                         nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1007                                         eintNonlocal);
1008             }
1009             wallcycle_stop(wcycle, ewcNS);
1010         }
1011         else
1012         {
1013             wallcycle_start(wcycle, ewcMOVEX);
1014             dd_move_x(cr->dd, box, x);
1015             wallcycle_stop(wcycle, ewcMOVEX);
1016
1017             wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1018             wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1019             nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1020                                             nbv->grp[eintNonlocal].nbat);
1021             wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1022             cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1023         }
1024
1025         if (bUseGPU && !bDiffKernels)
1026         {
1027             wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1028             /* launch non-local nonbonded F on GPU */
1029             do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1030                          nrnb, wcycle);
1031             cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1032         }
1033     }
1034
1035     if (bUseGPU)
1036     {
1037         /* launch D2H copy-back F */
1038         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1039         if (DOMAINDECOMP(cr) && !bDiffKernels)
1040         {
1041             nbnxn_gpu_launch_cpyback(nbv->gpu_nbv, nbv->grp[eintNonlocal].nbat,
1042                                      flags, eatNonlocal);
1043         }
1044         nbnxn_gpu_launch_cpyback(nbv->gpu_nbv, nbv->grp[eintLocal].nbat,
1045                                  flags, eatLocal);
1046         cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1047     }
1048
1049     if (bStateChanged && inputrecNeedMutot(inputrec))
1050     {
1051         if (PAR(cr))
1052         {
1053             gmx_sumd(2*DIM, mu, cr);
1054         }
1055
1056         for (i = 0; i < 2; i++)
1057         {
1058             for (j = 0; j < DIM; j++)
1059             {
1060                 fr->mu_tot[i][j] = mu[i*DIM + j];
1061             }
1062         }
1063     }
1064     if (fr->efep == efepNO)
1065     {
1066         copy_rvec(fr->mu_tot[0], mu_tot);
1067     }
1068     else
1069     {
1070         for (j = 0; j < DIM; j++)
1071         {
1072             mu_tot[j] =
1073                 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1074                 lambda[efptCOUL]*fr->mu_tot[1][j];
1075         }
1076     }
1077
1078     /* Reset energies */
1079     reset_enerdata(enerd);
1080     clear_rvecs(SHIFTS, fr->fshift);
1081
1082     if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1083     {
1084         wallcycle_start(wcycle, ewcPPDURINGPME);
1085         dd_force_flop_start(cr->dd, nrnb);
1086     }
1087
1088     if (inputrec->bRot)
1089     {
1090         /* Enforced rotation has its own cycle counter that starts after the collective
1091          * coordinates have been communicated. It is added to ddCyclF to allow
1092          * for proper load-balancing */
1093         wallcycle_start(wcycle, ewcROT);
1094         do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1095         wallcycle_stop(wcycle, ewcROT);
1096     }
1097
1098     /* Temporary solution until all routines take PaddedRVecVector */
1099     rvec *f = as_rvec_array(force->data());
1100
1101     /* Start the force cycle counter.
1102      * This counter is stopped after do_force_lowlevel.
1103      * No parallel communication should occur while this counter is running,
1104      * since that will interfere with the dynamic load balancing.
1105      */
1106     wallcycle_start(wcycle, ewcFORCE);
1107     if (bDoForces)
1108     {
1109         /* Reset forces for which the virial is calculated separately:
1110          * PME/Ewald forces if necessary */
1111         if (fr->bF_NoVirSum)
1112         {
1113             if (flags & GMX_FORCE_VIRIAL)
1114             {
1115                 fr->f_novirsum = fr->forceBufferNoVirialSummation;
1116             }
1117             else
1118             {
1119                 /* We are not calculating the pressure so we do not need
1120                  * a separate array for forces that do not contribute
1121                  * to the pressure.
1122                  */
1123                 fr->f_novirsum = force;
1124             }
1125         }
1126
1127         if (fr->bF_NoVirSum)
1128         {
1129             if (flags & GMX_FORCE_VIRIAL)
1130             {
1131                 /* TODO: remove this - 1 when padding is properly implemented */
1132                 clear_rvecs_omp(fr->f_novirsum->size() - 1,
1133                                 as_rvec_array(fr->f_novirsum->data()));
1134             }
1135         }
1136         /* Clear the short- and long-range forces */
1137         clear_rvecs_omp(fr->natoms_force_constr, f);
1138
1139         clear_rvec(fr->vir_diag_posres);
1140     }
1141
1142     if (inputrec->bPull && pull_have_constraint(inputrec->pull_work))
1143     {
1144         clear_pull_forces(inputrec->pull_work);
1145     }
1146
1147     /* We calculate the non-bonded forces, when done on the CPU, here.
1148      * We do this before calling do_force_lowlevel, because in that
1149      * function, the listed forces are calculated before PME, which
1150      * does communication.  With this order, non-bonded and listed
1151      * force calculation imbalance can be balanced out by the domain
1152      * decomposition load balancing.
1153      */
1154
1155     if (!bUseOrEmulGPU)
1156     {
1157         /* Maybe we should move this into do_force_lowlevel */
1158         do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1159                      nrnb, wcycle);
1160     }
1161
1162     if (fr->efep != efepNO)
1163     {
1164         /* Calculate the local and non-local free energy interactions here.
1165          * Happens here on the CPU both with and without GPU.
1166          */
1167         if (fr->nbv->grp[eintLocal].nbl_lists.nbl_fep[0]->nrj > 0)
1168         {
1169             do_nb_verlet_fep(&fr->nbv->grp[eintLocal].nbl_lists,
1170                              fr, x, f, mdatoms,
1171                              inputrec->fepvals, lambda,
1172                              enerd, flags, nrnb, wcycle);
1173         }
1174
1175         if (DOMAINDECOMP(cr) &&
1176             fr->nbv->grp[eintNonlocal].nbl_lists.nbl_fep[0]->nrj > 0)
1177         {
1178             do_nb_verlet_fep(&fr->nbv->grp[eintNonlocal].nbl_lists,
1179                              fr, x, f, mdatoms,
1180                              inputrec->fepvals, lambda,
1181                              enerd, flags, nrnb, wcycle);
1182         }
1183     }
1184
1185     if (!bUseOrEmulGPU || bDiffKernels)
1186     {
1187         int aloc;
1188
1189         if (DOMAINDECOMP(cr))
1190         {
1191             do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1192                          bDiffKernels ? enbvClearFYes : enbvClearFNo,
1193                          nrnb, wcycle);
1194         }
1195
1196         if (!bUseOrEmulGPU)
1197         {
1198             aloc = eintLocal;
1199         }
1200         else
1201         {
1202             aloc = eintNonlocal;
1203         }
1204
1205         /* Add all the non-bonded force to the normal force array.
1206          * This can be split into a local and a non-local part when overlapping
1207          * communication with calculation with domain decomposition.
1208          */
1209         cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1210         wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1211         wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1212         nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1213         wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1214         cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1215         wallcycle_start_nocount(wcycle, ewcFORCE);
1216
1217         /* if there are multiple fshift output buffers reduce them */
1218         if ((flags & GMX_FORCE_VIRIAL) &&
1219             nbv->grp[aloc].nbl_lists.nnbl > 1)
1220         {
1221             /* This is not in a subcounter because it takes a
1222                negligible and constant-sized amount of time */
1223             nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1224                                                      fr->fshift);
1225         }
1226     }
1227
1228     /* update QMMMrec, if necessary */
1229     if (fr->bQMMM)
1230     {
1231         update_QMMMrec(cr, fr, x, mdatoms, box, top);
1232     }
1233
1234     /* Compute the bonded and non-bonded energies and optionally forces */
1235     do_force_lowlevel(fr, inputrec, &(top->idef),
1236                       cr, nrnb, wcycle, mdatoms,
1237                       x, hist, f, enerd, fcd, top, fr->born,
1238                       bBornRadii, box,
1239                       inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1240                       flags, &cycles_pme);
1241
1242     cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1243
1244     if (ed)
1245     {
1246         do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1247     }
1248
1249     if (bUseOrEmulGPU && !bDiffKernels)
1250     {
1251         /* wait for non-local forces (or calculate in emulation mode) */
1252         if (DOMAINDECOMP(cr))
1253         {
1254             if (bUseGPU)
1255             {
1256                 float cycles_tmp;
1257
1258                 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1259                 nbnxn_gpu_wait_for_gpu(nbv->gpu_nbv,
1260                                        flags, eatNonlocal,
1261                                        enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1262                                        fr->fshift);
1263                 cycles_tmp       = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1264                 cycles_wait_gpu += cycles_tmp;
1265                 cycles_force    += cycles_tmp;
1266             }
1267             else
1268             {
1269                 wallcycle_start_nocount(wcycle, ewcFORCE);
1270                 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1271                              nrnb, wcycle);
1272                 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1273             }
1274             wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1275             wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1276             /* skip the reduction if there was no non-local work to do */
1277             if (nbv->grp[eintNonlocal].nbl_lists.nbl[0]->nsci > 0)
1278             {
1279                 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1280                                                nbv->grp[eintNonlocal].nbat, f);
1281             }
1282             wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1283             cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1284         }
1285     }
1286
1287     if (bDoForces && DOMAINDECOMP(cr))
1288     {
1289         if (bUseGPU)
1290         {
1291             /* We are done with the CPU compute, but the GPU local non-bonded
1292              * kernel can still be running while we communicate the forces.
1293              * We start a counter here, so we can, hopefully, time the rest
1294              * of the GPU kernel execution and data transfer.
1295              */
1296             cycleCountBeforeLocalWorkCompletes = gmx_cycles_read();
1297         }
1298
1299         /* Communicate the forces */
1300         wallcycle_start(wcycle, ewcMOVEF);
1301         dd_move_f(cr->dd, f, fr->fshift);
1302         wallcycle_stop(wcycle, ewcMOVEF);
1303     }
1304
1305     if (bUseOrEmulGPU)
1306     {
1307         /* wait for local forces (or calculate in emulation mode) */
1308         if (bUseGPU)
1309         {
1310             float       cycles_tmp, cycles_wait_est;
1311             /* Measured overhead on CUDA and OpenCL with(out) GPU sharing
1312              * is between 0.5 and 1.5 Mcycles. So 2 MCycles is an overestimate,
1313              * but even with a step of 0.1 ms the difference is less than 1%
1314              * of the step time.
1315              */
1316             const float gpuWaitApiOverheadMargin = 2e6f; /* cycles */
1317
1318             wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1319             nbnxn_gpu_wait_for_gpu(nbv->gpu_nbv,
1320                                    flags, eatLocal,
1321                                    enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1322                                    fr->fshift);
1323             cycles_tmp      = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1324
1325             if (bDoForces && DOMAINDECOMP(cr))
1326             {
1327                 cycles_wait_est = gmx_cycles_read() - cycleCountBeforeLocalWorkCompletes;
1328
1329                 if (cycles_tmp < gpuWaitApiOverheadMargin)
1330                 {
1331                     /* We measured few cycles, it could be that the kernel
1332                      * and transfer finished earlier and there was no actual
1333                      * wait time, only API call overhead.
1334                      * Then the actual time could be anywhere between 0 and
1335                      * cycles_wait_est. As a compromise, we use half the time.
1336                      */
1337                     cycles_wait_est *= 0.5f;
1338                 }
1339             }
1340             else
1341             {
1342                 /* No force communication so we actually timed the wait */
1343                 cycles_wait_est = cycles_tmp;
1344             }
1345             /* Even though this is after dd_move_f, the actual task we are
1346              * waiting for runs asynchronously with dd_move_f and we usually
1347              * have nothing to balance it with, so we can and should add
1348              * the time to the force time for load balancing.
1349              */
1350             cycles_force    += cycles_wait_est;
1351             cycles_wait_gpu += cycles_wait_est;
1352
1353             /* now clear the GPU outputs while we finish the step on the CPU */
1354             wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1355             nbnxn_gpu_clear_outputs(nbv->gpu_nbv, flags);
1356             wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1357         }
1358         else
1359         {
1360             wallcycle_start_nocount(wcycle, ewcFORCE);
1361             do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1362                          DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1363                          nrnb, wcycle);
1364             wallcycle_stop(wcycle, ewcFORCE);
1365         }
1366         wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1367         wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1368         nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1369                                        nbv->grp[eintLocal].nbat, f);
1370         wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1371         wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1372     }
1373
1374     if (DOMAINDECOMP(cr))
1375     {
1376         dd_force_flop_stop(cr->dd, nrnb);
1377         if (wcycle)
1378         {
1379             dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1380             if (bUseGPU)
1381             {
1382                 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1383             }
1384         }
1385     }
1386
1387     if (bDoForces)
1388     {
1389         /* Compute forces due to electric field */
1390         if (fr->efield != nullptr)
1391         {
1392             fr->efield->calculateForces(cr, mdatoms, fr->f_novirsum, t);
1393         }
1394
1395         /* If we have NoVirSum forces, but we do not calculate the virial,
1396          * we sum fr->f_novirsum=f later.
1397          */
1398         if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1399         {
1400             wallcycle_start(wcycle, ewcVSITESPREAD);
1401             spread_vsite_f(vsite, x, f, fr->fshift, FALSE, nullptr, nrnb,
1402                            &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1403             wallcycle_stop(wcycle, ewcVSITESPREAD);
1404         }
1405
1406         if (flags & GMX_FORCE_VIRIAL)
1407         {
1408             /* Calculation of the virial must be done after vsites! */
1409             calc_virial(0, mdatoms->homenr, x, f,
1410                         vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1411         }
1412     }
1413
1414     if (inputrec->bPull && pull_have_potential(inputrec->pull_work))
1415     {
1416         /* Since the COM pulling is always done mass-weighted, no forces are
1417          * applied to vsites and this call can be done after vsite spreading.
1418          */
1419         pull_potential_wrapper(cr, inputrec, box, x,
1420                                f, vir_force, mdatoms, enerd, lambda, t,
1421                                wcycle);
1422     }
1423
1424     /* Add the forces from enforced rotation potentials (if any) */
1425     if (inputrec->bRot)
1426     {
1427         wallcycle_start(wcycle, ewcROTadd);
1428         enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1429         wallcycle_stop(wcycle, ewcROTadd);
1430     }
1431
1432     /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1433     IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1434
1435     if (PAR(cr) && !(cr->duty & DUTY_PME))
1436     {
1437         /* In case of node-splitting, the PP nodes receive the long-range
1438          * forces, virial and energy from the PME nodes here.
1439          */
1440         pme_receive_force_ener(cr, wcycle, enerd, fr);
1441     }
1442
1443     if (bDoForces)
1444     {
1445         post_process_forces(cr, step, nrnb, wcycle,
1446                             top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1447                             flags);
1448     }
1449
1450     if (flags & GMX_FORCE_ENERGY)
1451     {
1452         /* Sum the potential energy terms from group contributions */
1453         sum_epot(&(enerd->grpp), enerd->term);
1454
1455         if (!EI_TPI(inputrec->eI))
1456         {
1457             checkPotentialEnergyValidity(enerd);
1458         }
1459     }
1460 }
1461
1462 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1463                         t_inputrec *inputrec,
1464                         gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1465                         gmx_localtop_t *top,
1466                         gmx_groups_t *groups,
1467                         matrix box, rvec x[], history_t *hist,
1468                         PaddedRVecVector *force,
1469                         tensor vir_force,
1470                         t_mdatoms *mdatoms,
1471                         gmx_enerdata_t *enerd, t_fcdata *fcd,
1472                         real *lambda, t_graph *graph,
1473                         t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1474                         double t, gmx_edsam_t ed,
1475                         gmx_bool bBornRadii,
1476                         int flags)
1477 {
1478     int        cg0, cg1, i, j;
1479     double     mu[2*DIM];
1480     gmx_bool   bStateChanged, bNS, bFillGrid, bCalcCGCM;
1481     gmx_bool   bDoForces;
1482     float      cycles_pme, cycles_force;
1483
1484     const int  start  = 0;
1485     const int  homenr = mdatoms->homenr;
1486
1487     clear_mat(vir_force);
1488
1489     cg0 = 0;
1490     if (DOMAINDECOMP(cr))
1491     {
1492         cg1 = cr->dd->ncg_tot;
1493     }
1494     else
1495     {
1496         cg1 = top->cgs.nr;
1497     }
1498     if (fr->n_tpi > 0)
1499     {
1500         cg1--;
1501     }
1502
1503     bStateChanged  = (flags & GMX_FORCE_STATECHANGED);
1504     bNS            = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1505     /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1506     bFillGrid      = (bNS && bStateChanged);
1507     bCalcCGCM      = (bFillGrid && !DOMAINDECOMP(cr));
1508     bDoForces      = (flags & GMX_FORCE_FORCES);
1509
1510     if (bStateChanged)
1511     {
1512         update_forcerec(fr, box);
1513
1514         if (inputrecNeedMutot(inputrec))
1515         {
1516             /* Calculate total (local) dipole moment in a temporary common array.
1517              * This makes it possible to sum them over nodes faster.
1518              */
1519             calc_mu(start, homenr,
1520                     x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1521                     mu, mu+DIM);
1522         }
1523     }
1524
1525     if (fr->ePBC != epbcNONE)
1526     {
1527         /* Compute shift vectors every step,
1528          * because of pressure coupling or box deformation!
1529          */
1530         if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1531         {
1532             calc_shifts(box, fr->shift_vec);
1533         }
1534
1535         if (bCalcCGCM)
1536         {
1537             put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1538                                      &(top->cgs), x, fr->cg_cm);
1539             inc_nrnb(nrnb, eNR_CGCM, homenr);
1540             inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1541         }
1542         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1543         {
1544             unshift_self(graph, box, x);
1545         }
1546     }
1547     else if (bCalcCGCM)
1548     {
1549         calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1550         inc_nrnb(nrnb, eNR_CGCM, homenr);
1551     }
1552
1553     if (bCalcCGCM && gmx_debug_at)
1554     {
1555         pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1556     }
1557
1558 #if GMX_MPI
1559     if (!(cr->duty & DUTY_PME))
1560     {
1561         gmx_bool bBS;
1562         matrix   boxs;
1563
1564         /* Send particle coordinates to the pme nodes.
1565          * Since this is only implemented for domain decomposition
1566          * and domain decomposition does not use the graph,
1567          * we do not need to worry about shifting.
1568          */
1569
1570         wallcycle_start(wcycle, ewcPP_PMESENDX);
1571
1572         bBS = (inputrec->nwall == 2);
1573         if (bBS)
1574         {
1575             copy_mat(box, boxs);
1576             svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1577         }
1578
1579         gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1580                                  lambda[efptCOUL], lambda[efptVDW],
1581                                  (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1582                                  step);
1583
1584         wallcycle_stop(wcycle, ewcPP_PMESENDX);
1585     }
1586 #endif /* GMX_MPI */
1587
1588     /* Communicate coordinates and sum dipole if necessary */
1589     if (DOMAINDECOMP(cr))
1590     {
1591         wallcycle_start(wcycle, ewcMOVEX);
1592         dd_move_x(cr->dd, box, x);
1593         wallcycle_stop(wcycle, ewcMOVEX);
1594     }
1595
1596     if (inputrecNeedMutot(inputrec))
1597     {
1598         if (bStateChanged)
1599         {
1600             if (PAR(cr))
1601             {
1602                 gmx_sumd(2*DIM, mu, cr);
1603             }
1604             for (i = 0; i < 2; i++)
1605             {
1606                 for (j = 0; j < DIM; j++)
1607                 {
1608                     fr->mu_tot[i][j] = mu[i*DIM + j];
1609                 }
1610             }
1611         }
1612         if (fr->efep == efepNO)
1613         {
1614             copy_rvec(fr->mu_tot[0], mu_tot);
1615         }
1616         else
1617         {
1618             for (j = 0; j < DIM; j++)
1619             {
1620                 mu_tot[j] =
1621                     (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1622             }
1623         }
1624     }
1625
1626     /* Reset energies */
1627     reset_enerdata(enerd);
1628     clear_rvecs(SHIFTS, fr->fshift);
1629
1630     if (bNS)
1631     {
1632         wallcycle_start(wcycle, ewcNS);
1633
1634         if (graph && bStateChanged)
1635         {
1636             /* Calculate intramolecular shift vectors to make molecules whole */
1637             mk_mshift(fplog, graph, fr->ePBC, box, x);
1638         }
1639
1640         /* Do the actual neighbour searching */
1641         ns(fplog, fr, box,
1642            groups, top, mdatoms,
1643            cr, nrnb, bFillGrid);
1644
1645         wallcycle_stop(wcycle, ewcNS);
1646     }
1647
1648     if (inputrec->implicit_solvent && bNS)
1649     {
1650         make_gb_nblist(cr, inputrec->gb_algorithm,
1651                        x, box, fr, &top->idef, graph, fr->born);
1652     }
1653
1654     if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1655     {
1656         wallcycle_start(wcycle, ewcPPDURINGPME);
1657         dd_force_flop_start(cr->dd, nrnb);
1658     }
1659
1660     if (inputrec->bRot)
1661     {
1662         /* Enforced rotation has its own cycle counter that starts after the collective
1663          * coordinates have been communicated. It is added to ddCyclF to allow
1664          * for proper load-balancing */
1665         wallcycle_start(wcycle, ewcROT);
1666         do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1667         wallcycle_stop(wcycle, ewcROT);
1668     }
1669
1670     /* Temporary solution until all routines take PaddedRVecVector */
1671     rvec *f = as_rvec_array(force->data());
1672
1673     /* Start the force cycle counter.
1674      * This counter is stopped after do_force_lowlevel.
1675      * No parallel communication should occur while this counter is running,
1676      * since that will interfere with the dynamic load balancing.
1677      */
1678     wallcycle_start(wcycle, ewcFORCE);
1679
1680     if (bDoForces)
1681     {
1682         /* Reset forces for which the virial is calculated separately:
1683          * PME/Ewald forces if necessary */
1684         if (fr->bF_NoVirSum)
1685         {
1686             if (flags & GMX_FORCE_VIRIAL)
1687             {
1688                 fr->f_novirsum = fr->forceBufferNoVirialSummation;
1689                 /* TODO: remove this - 1 when padding is properly implemented */
1690                 clear_rvecs(fr->f_novirsum->size() - 1,
1691                             as_rvec_array(fr->f_novirsum->data()));
1692             }
1693             else
1694             {
1695                 /* We are not calculating the pressure so we do not need
1696                  * a separate array for forces that do not contribute
1697                  * to the pressure.
1698                  */
1699                 fr->f_novirsum = force;
1700             }
1701         }
1702
1703         /* Clear the short- and long-range forces */
1704         clear_rvecs(fr->natoms_force_constr, f);
1705
1706         clear_rvec(fr->vir_diag_posres);
1707     }
1708     if (inputrec->bPull && pull_have_constraint(inputrec->pull_work))
1709     {
1710         clear_pull_forces(inputrec->pull_work);
1711     }
1712
1713     /* update QMMMrec, if necessary */
1714     if (fr->bQMMM)
1715     {
1716         update_QMMMrec(cr, fr, x, mdatoms, box, top);
1717     }
1718
1719     /* Compute the bonded and non-bonded energies and optionally forces */
1720     do_force_lowlevel(fr, inputrec, &(top->idef),
1721                       cr, nrnb, wcycle, mdatoms,
1722                       x, hist, f, enerd, fcd, top, fr->born,
1723                       bBornRadii, box,
1724                       inputrec->fepvals, lambda,
1725                       graph, &(top->excls), fr->mu_tot,
1726                       flags,
1727                       &cycles_pme);
1728
1729     cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1730
1731     if (ed)
1732     {
1733         do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1734     }
1735
1736     if (DOMAINDECOMP(cr))
1737     {
1738         dd_force_flop_stop(cr->dd, nrnb);
1739         if (wcycle)
1740         {
1741             dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1742         }
1743     }
1744
1745     if (bDoForces)
1746     {
1747         /* Compute forces due to electric field */
1748         if (fr->efield != nullptr)
1749         {
1750             fr->efield->calculateForces(cr, mdatoms, fr->f_novirsum, t);
1751         }
1752
1753         /* Communicate the forces */
1754         if (DOMAINDECOMP(cr))
1755         {
1756             wallcycle_start(wcycle, ewcMOVEF);
1757             dd_move_f(cr->dd, f, fr->fshift);
1758             /* Do we need to communicate the separate force array
1759              * for terms that do not contribute to the single sum virial?
1760              * Position restraints and electric fields do not introduce
1761              * inter-cg forces, only full electrostatics methods do.
1762              * When we do not calculate the virial, fr->f_novirsum = f,
1763              * so we have already communicated these forces.
1764              */
1765             if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1766                 (flags & GMX_FORCE_VIRIAL))
1767             {
1768                 dd_move_f(cr->dd, as_rvec_array(fr->f_novirsum->data()), nullptr);
1769             }
1770             wallcycle_stop(wcycle, ewcMOVEF);
1771         }
1772
1773         /* If we have NoVirSum forces, but we do not calculate the virial,
1774          * we sum fr->f_novirsum=f later.
1775          */
1776         if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1777         {
1778             wallcycle_start(wcycle, ewcVSITESPREAD);
1779             spread_vsite_f(vsite, x, f, fr->fshift, FALSE, nullptr, nrnb,
1780                            &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1781             wallcycle_stop(wcycle, ewcVSITESPREAD);
1782         }
1783
1784         if (flags & GMX_FORCE_VIRIAL)
1785         {
1786             /* Calculation of the virial must be done after vsites! */
1787             calc_virial(0, mdatoms->homenr, x, f,
1788                         vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1789         }
1790     }
1791
1792     if (inputrec->bPull && pull_have_potential(inputrec->pull_work))
1793     {
1794         pull_potential_wrapper(cr, inputrec, box, x,
1795                                f, vir_force, mdatoms, enerd, lambda, t,
1796                                wcycle);
1797     }
1798
1799     /* Add the forces from enforced rotation potentials (if any) */
1800     if (inputrec->bRot)
1801     {
1802         wallcycle_start(wcycle, ewcROTadd);
1803         enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1804         wallcycle_stop(wcycle, ewcROTadd);
1805     }
1806
1807     /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1808     IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1809
1810     if (PAR(cr) && !(cr->duty & DUTY_PME))
1811     {
1812         /* In case of node-splitting, the PP nodes receive the long-range
1813          * forces, virial and energy from the PME nodes here.
1814          */
1815         pme_receive_force_ener(cr, wcycle, enerd, fr);
1816     }
1817
1818     if (bDoForces)
1819     {
1820         post_process_forces(cr, step, nrnb, wcycle,
1821                             top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1822                             flags);
1823     }
1824
1825     if (flags & GMX_FORCE_ENERGY)
1826     {
1827         /* Sum the potential energy terms from group contributions */
1828         sum_epot(&(enerd->grpp), enerd->term);
1829
1830         if (!EI_TPI(inputrec->eI))
1831         {
1832             checkPotentialEnergyValidity(enerd);
1833         }
1834     }
1835
1836 }
1837
1838 void do_force(FILE *fplog, t_commrec *cr,
1839               t_inputrec *inputrec,
1840               gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1841               gmx_localtop_t *top,
1842               gmx_groups_t *groups,
1843               matrix box, PaddedRVecVector *coordinates, history_t *hist,
1844               PaddedRVecVector *force,
1845               tensor vir_force,
1846               t_mdatoms *mdatoms,
1847               gmx_enerdata_t *enerd, t_fcdata *fcd,
1848               gmx::ArrayRef<real> lambda, t_graph *graph,
1849               t_forcerec *fr,
1850               gmx_vsite_t *vsite, rvec mu_tot,
1851               double t, gmx_edsam_t ed,
1852               gmx_bool bBornRadii,
1853               int flags)
1854 {
1855     /* modify force flag if not doing nonbonded */
1856     if (!fr->bNonbonded)
1857     {
1858         flags &= ~GMX_FORCE_NONBONDED;
1859     }
1860
1861     GMX_ASSERT(coordinates->size() >= static_cast<unsigned int>(fr->natoms_force + 1), "We might need 1 element extra for SIMD");
1862     GMX_ASSERT(force->size() >= static_cast<unsigned int>(fr->natoms_force + 1), "We might need 1 element extra for SIMD");
1863
1864     rvec *x = as_rvec_array(coordinates->data());
1865
1866     switch (inputrec->cutoff_scheme)
1867     {
1868         case ecutsVERLET:
1869             do_force_cutsVERLET(fplog, cr, inputrec,
1870                                 step, nrnb, wcycle,
1871                                 top,
1872                                 groups,
1873                                 box, x, hist,
1874                                 force, vir_force,
1875                                 mdatoms,
1876                                 enerd, fcd,
1877                                 lambda.data(), graph,
1878                                 fr, fr->ic,
1879                                 vsite, mu_tot,
1880                                 t, ed,
1881                                 bBornRadii,
1882                                 flags);
1883             break;
1884         case ecutsGROUP:
1885             do_force_cutsGROUP(fplog, cr, inputrec,
1886                                step, nrnb, wcycle,
1887                                top,
1888                                groups,
1889                                box, x, hist,
1890                                force, vir_force,
1891                                mdatoms,
1892                                enerd, fcd,
1893                                lambda.data(), graph,
1894                                fr, vsite, mu_tot,
1895                                t, ed,
1896                                bBornRadii,
1897                                flags);
1898             break;
1899         default:
1900             gmx_incons("Invalid cut-off scheme passed!");
1901     }
1902 }
1903
1904
1905 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
1906                         t_inputrec *ir, t_mdatoms *md,
1907                         t_state *state, t_commrec *cr, t_nrnb *nrnb,
1908                         t_forcerec *fr, gmx_localtop_t *top)
1909 {
1910     int             i, m, start, end;
1911     gmx_int64_t     step;
1912     real            dt = ir->delta_t;
1913     real            dvdl_dum;
1914     rvec           *savex;
1915
1916     /* We need to allocate one element extra, since we might use
1917      * (unaligned) 4-wide SIMD loads to access rvec entries.
1918      */
1919     snew(savex, state->natoms + 1);
1920
1921     start = 0;
1922     end   = md->homenr;
1923
1924     if (debug)
1925     {
1926         fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
1927                 start, md->homenr, end);
1928     }
1929     /* Do a first constrain to reset particles... */
1930     step = ir->init_step;
1931     if (fplog)
1932     {
1933         char buf[STEPSTRSIZE];
1934         fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
1935                 gmx_step_str(step, buf));
1936     }
1937     dvdl_dum = 0;
1938
1939     /* constrain the current position */
1940     constrain(nullptr, TRUE, FALSE, constr, &(top->idef),
1941               ir, cr, step, 0, 1.0, md,
1942               as_rvec_array(state->x.data()), as_rvec_array(state->x.data()), nullptr,
1943               fr->bMolPBC, state->box,
1944               state->lambda[efptBONDED], &dvdl_dum,
1945               nullptr, nullptr, nrnb, econqCoord);
1946     if (EI_VV(ir->eI))
1947     {
1948         /* constrain the inital velocity, and save it */
1949         /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1950         constrain(nullptr, TRUE, FALSE, constr, &(top->idef),
1951                   ir, cr, step, 0, 1.0, md,
1952                   as_rvec_array(state->x.data()), as_rvec_array(state->v.data()), as_rvec_array(state->v.data()),
1953                   fr->bMolPBC, state->box,
1954                   state->lambda[efptBONDED], &dvdl_dum,
1955                   nullptr, nullptr, nrnb, econqVeloc);
1956     }
1957     /* constrain the inital velocities at t-dt/2 */
1958     if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
1959     {
1960         for (i = start; (i < end); i++)
1961         {
1962             for (m = 0; (m < DIM); m++)
1963             {
1964                 /* Reverse the velocity */
1965                 state->v[i][m] = -state->v[i][m];
1966                 /* Store the position at t-dt in buf */
1967                 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
1968             }
1969         }
1970         /* Shake the positions at t=-dt with the positions at t=0
1971          * as reference coordinates.
1972          */
1973         if (fplog)
1974         {
1975             char buf[STEPSTRSIZE];
1976             fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
1977                     gmx_step_str(step, buf));
1978         }
1979         dvdl_dum = 0;
1980         constrain(nullptr, TRUE, FALSE, constr, &(top->idef),
1981                   ir, cr, step, -1, 1.0, md,
1982                   as_rvec_array(state->x.data()), savex, nullptr,
1983                   fr->bMolPBC, state->box,
1984                   state->lambda[efptBONDED], &dvdl_dum,
1985                   as_rvec_array(state->v.data()), nullptr, nrnb, econqCoord);
1986
1987         for (i = start; i < end; i++)
1988         {
1989             for (m = 0; m < DIM; m++)
1990             {
1991                 /* Re-reverse the velocities */
1992                 state->v[i][m] = -state->v[i][m];
1993             }
1994         }
1995     }
1996     sfree(savex);
1997 }
1998
1999
2000 static void
2001 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2002                 double *enerout, double *virout)
2003 {
2004     double enersum, virsum;
2005     double invscale, invscale2, invscale3;
2006     double r, ea, eb, ec, pa, pb, pc, pd;
2007     double y0, f, g, h;
2008     int    ri, offset;
2009     double tabfactor;
2010
2011     invscale  = 1.0/scale;
2012     invscale2 = invscale*invscale;
2013     invscale3 = invscale*invscale2;
2014
2015     /* Following summation derived from cubic spline definition,
2016      * Numerical Recipies in C, second edition, p. 113-116.  Exact for
2017      * the cubic spline.  We first calculate the negative of the
2018      * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2019      * add the more standard, abrupt cutoff correction to that result,
2020      * yielding the long-range correction for a switched function.  We
2021      * perform both the pressure and energy loops at the same time for
2022      * simplicity, as the computational cost is low. */
2023
2024     if (offstart == 0)
2025     {
2026         /* Since the dispersion table has been scaled down a factor
2027          * 6.0 and the repulsion a factor 12.0 to compensate for the
2028          * c6/c12 parameters inside nbfp[] being scaled up (to save
2029          * flops in kernels), we need to correct for this.
2030          */
2031         tabfactor = 6.0;
2032     }
2033     else
2034     {
2035         tabfactor = 12.0;
2036     }
2037
2038     enersum = 0.0;
2039     virsum  = 0.0;
2040     for (ri = rstart; ri < rend; ++ri)
2041     {
2042         r  = ri*invscale;
2043         ea = invscale3;
2044         eb = 2.0*invscale2*r;
2045         ec = invscale*r*r;
2046
2047         pa = invscale3;
2048         pb = 3.0*invscale2*r;
2049         pc = 3.0*invscale*r*r;
2050         pd = r*r*r;
2051
2052         /* this "8" is from the packing in the vdwtab array - perhaps
2053            should be defined? */
2054
2055         offset = 8*ri + offstart;
2056         y0     = vdwtab[offset];
2057         f      = vdwtab[offset+1];
2058         g      = vdwtab[offset+2];
2059         h      = vdwtab[offset+3];
2060
2061         enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2) + g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
2062         virsum  +=  f*(pa/4 + pb/3 + pc/2 + pd) + 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
2063     }
2064     *enerout = 4.0*M_PI*enersum*tabfactor;
2065     *virout  = 4.0*M_PI*virsum*tabfactor;
2066 }
2067
2068 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2069 {
2070     double   eners[2], virs[2], enersum, virsum;
2071     double   r0, rc3, rc9;
2072     int      ri0, ri1, i;
2073     real     scale, *vdwtab;
2074
2075     fr->enershiftsix    = 0;
2076     fr->enershifttwelve = 0;
2077     fr->enerdiffsix     = 0;
2078     fr->enerdifftwelve  = 0;
2079     fr->virdiffsix      = 0;
2080     fr->virdifftwelve   = 0;
2081
2082     if (eDispCorr != edispcNO)
2083     {
2084         for (i = 0; i < 2; i++)
2085         {
2086             eners[i] = 0;
2087             virs[i]  = 0;
2088         }
2089         if ((fr->vdw_modifier == eintmodPOTSHIFT) ||
2090             (fr->vdw_modifier == eintmodPOTSWITCH) ||
2091             (fr->vdw_modifier == eintmodFORCESWITCH) ||
2092             (fr->vdwtype == evdwSHIFT) ||
2093             (fr->vdwtype == evdwSWITCH))
2094         {
2095             if (((fr->vdw_modifier == eintmodPOTSWITCH) ||
2096                  (fr->vdw_modifier == eintmodFORCESWITCH) ||
2097                  (fr->vdwtype == evdwSWITCH)) && fr->rvdw_switch == 0)
2098             {
2099                 gmx_fatal(FARGS,
2100                           "With dispersion correction rvdw-switch can not be zero "
2101                           "for vdw-type = %s", evdw_names[fr->vdwtype]);
2102             }
2103
2104             /* TODO This code depends on the logic in tables.c that
2105                constructs the table layout, which should be made
2106                explicit in future cleanup. */
2107             GMX_ASSERT(fr->dispersionCorrectionTable->interaction == GMX_TABLE_INTERACTION_VDWREP_VDWDISP,
2108                        "Dispersion-correction code needs a table with both repulsion and dispersion terms");
2109             scale  = fr->dispersionCorrectionTable->scale;
2110             vdwtab = fr->dispersionCorrectionTable->data;
2111
2112             /* Round the cut-offs to exact table values for precision */
2113             ri0  = static_cast<int>(floor(fr->rvdw_switch*scale));
2114             ri1  = static_cast<int>(ceil(fr->rvdw*scale));
2115
2116             /* The code below has some support for handling force-switching, i.e.
2117              * when the force (instead of potential) is switched over a limited
2118              * region. This leads to a constant shift in the potential inside the
2119              * switching region, which we can handle by adding a constant energy
2120              * term in the force-switch case just like when we do potential-shift.
2121              *
2122              * For now this is not enabled, but to keep the functionality in the
2123              * code we check separately for switch and shift. When we do force-switch
2124              * the shifting point is rvdw_switch, while it is the cutoff when we
2125              * have a classical potential-shift.
2126              *
2127              * For a pure potential-shift the potential has a constant shift
2128              * all the way out to the cutoff, and that is it. For other forms
2129              * we need to calculate the constant shift up to the point where we
2130              * start modifying the potential.
2131              */
2132             ri0  = (fr->vdw_modifier == eintmodPOTSHIFT) ? ri1 : ri0;
2133
2134             r0   = ri0/scale;
2135             rc3  = r0*r0*r0;
2136             rc9  = rc3*rc3*rc3;
2137
2138             if ((fr->vdw_modifier == eintmodFORCESWITCH) ||
2139                 (fr->vdwtype == evdwSHIFT))
2140             {
2141                 /* Determine the constant energy shift below rvdw_switch.
2142                  * Table has a scale factor since we have scaled it down to compensate
2143                  * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2144                  */
2145                 fr->enershiftsix    = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2146                 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2147             }
2148             else if (fr->vdw_modifier == eintmodPOTSHIFT)
2149             {
2150                 fr->enershiftsix    = (real)(-1.0/(rc3*rc3));
2151                 fr->enershifttwelve = (real)( 1.0/(rc9*rc3));
2152             }
2153
2154             /* Add the constant part from 0 to rvdw_switch.
2155              * This integration from 0 to rvdw_switch overcounts the number
2156              * of interactions by 1, as it also counts the self interaction.
2157              * We will correct for this later.
2158              */
2159             eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2160             eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2161
2162             /* Calculate the contribution in the range [r0,r1] where we
2163              * modify the potential. For a pure potential-shift modifier we will
2164              * have ri0==ri1, and there will not be any contribution here.
2165              */
2166             for (i = 0; i < 2; i++)
2167             {
2168                 enersum = 0;
2169                 virsum  = 0;
2170                 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2171                 eners[i] -= enersum;
2172                 virs[i]  -= virsum;
2173             }
2174
2175             /* Alright: Above we compensated by REMOVING the parts outside r0
2176              * corresponding to the ideal VdW 1/r6 and /r12 potentials.
2177              *
2178              * Regardless of whether r0 is the point where we start switching,
2179              * or the cutoff where we calculated the constant shift, we include
2180              * all the parts we are missing out to infinity from r0 by
2181              * calculating the analytical dispersion correction.
2182              */
2183             eners[0] += -4.0*M_PI/(3.0*rc3);
2184             eners[1] +=  4.0*M_PI/(9.0*rc9);
2185             virs[0]  +=  8.0*M_PI/rc3;
2186             virs[1]  += -16.0*M_PI/(3.0*rc9);
2187         }
2188         else if (fr->vdwtype == evdwCUT ||
2189                  EVDW_PME(fr->vdwtype) ||
2190                  fr->vdwtype == evdwUSER)
2191         {
2192             if (fr->vdwtype == evdwUSER && fplog)
2193             {
2194                 fprintf(fplog,
2195                         "WARNING: using dispersion correction with user tables\n");
2196             }
2197
2198             /* Note that with LJ-PME, the dispersion correction is multiplied
2199              * by the difference between the actual C6 and the value of C6
2200              * that would produce the combination rule.
2201              * This means the normal energy and virial difference formulas
2202              * can be used here.
2203              */
2204
2205             rc3  = fr->rvdw*fr->rvdw*fr->rvdw;
2206             rc9  = rc3*rc3*rc3;
2207             /* Contribution beyond the cut-off */
2208             eners[0] += -4.0*M_PI/(3.0*rc3);
2209             eners[1] +=  4.0*M_PI/(9.0*rc9);
2210             if (fr->vdw_modifier == eintmodPOTSHIFT)
2211             {
2212                 /* Contribution within the cut-off */
2213                 eners[0] += -4.0*M_PI/(3.0*rc3);
2214                 eners[1] +=  4.0*M_PI/(3.0*rc9);
2215             }
2216             /* Contribution beyond the cut-off */
2217             virs[0]  +=  8.0*M_PI/rc3;
2218             virs[1]  += -16.0*M_PI/(3.0*rc9);
2219         }
2220         else
2221         {
2222             gmx_fatal(FARGS,
2223                       "Dispersion correction is not implemented for vdw-type = %s",
2224                       evdw_names[fr->vdwtype]);
2225         }
2226
2227         /* When we deprecate the group kernels the code below can go too */
2228         if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2229         {
2230             /* Calculate self-interaction coefficient (assuming that
2231              * the reciprocal-space contribution is constant in the
2232              * region that contributes to the self-interaction).
2233              */
2234             fr->enershiftsix = gmx::power6(fr->ewaldcoeff_lj) / 6.0;
2235
2236             eners[0] += -gmx::power3(std::sqrt(M_PI)*fr->ewaldcoeff_lj)/3.0;
2237             virs[0]  +=  gmx::power3(std::sqrt(M_PI)*fr->ewaldcoeff_lj);
2238         }
2239
2240         fr->enerdiffsix    = eners[0];
2241         fr->enerdifftwelve = eners[1];
2242         /* The 0.5 is due to the Gromacs definition of the virial */
2243         fr->virdiffsix     = 0.5*virs[0];
2244         fr->virdifftwelve  = 0.5*virs[1];
2245     }
2246 }
2247
2248 void calc_dispcorr(t_inputrec *ir, t_forcerec *fr,
2249                    matrix box, real lambda, tensor pres, tensor virial,
2250                    real *prescorr, real *enercorr, real *dvdlcorr)
2251 {
2252     gmx_bool bCorrAll, bCorrPres;
2253     real     dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2254     int      m;
2255
2256     *prescorr = 0;
2257     *enercorr = 0;
2258     *dvdlcorr = 0;
2259
2260     clear_mat(virial);
2261     clear_mat(pres);
2262
2263     if (ir->eDispCorr != edispcNO)
2264     {
2265         bCorrAll  = (ir->eDispCorr == edispcAllEner ||
2266                      ir->eDispCorr == edispcAllEnerPres);
2267         bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2268                      ir->eDispCorr == edispcAllEnerPres);
2269
2270         invvol = 1/det(box);
2271         if (fr->n_tpi)
2272         {
2273             /* Only correct for the interactions with the inserted molecule */
2274             dens   = (fr->numAtomsForDispersionCorrection - fr->n_tpi)*invvol;
2275             ninter = fr->n_tpi;
2276         }
2277         else
2278         {
2279             dens   = fr->numAtomsForDispersionCorrection*invvol;
2280             ninter = 0.5*fr->numAtomsForDispersionCorrection;
2281         }
2282
2283         if (ir->efep == efepNO)
2284         {
2285             avcsix    = fr->avcsix[0];
2286             avctwelve = fr->avctwelve[0];
2287         }
2288         else
2289         {
2290             avcsix    = (1 - lambda)*fr->avcsix[0]    + lambda*fr->avcsix[1];
2291             avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2292         }
2293
2294         enerdiff   = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2295         *enercorr += avcsix*enerdiff;
2296         dvdlambda  = 0.0;
2297         if (ir->efep != efepNO)
2298         {
2299             dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2300         }
2301         if (bCorrAll)
2302         {
2303             enerdiff   = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2304             *enercorr += avctwelve*enerdiff;
2305             if (fr->efep != efepNO)
2306             {
2307                 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2308             }
2309         }
2310
2311         if (bCorrPres)
2312         {
2313             svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2314             if (ir->eDispCorr == edispcAllEnerPres)
2315             {
2316                 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2317             }
2318             /* The factor 2 is because of the Gromacs virial definition */
2319             spres = -2.0*invvol*svir*PRESFAC;
2320
2321             for (m = 0; m < DIM; m++)
2322             {
2323                 virial[m][m] += svir;
2324                 pres[m][m]   += spres;
2325             }
2326             *prescorr += spres;
2327         }
2328
2329         /* Can't currently control when it prints, for now, just print when degugging */
2330         if (debug)
2331         {
2332             if (bCorrAll)
2333             {
2334                 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2335                         avcsix, avctwelve);
2336             }
2337             if (bCorrPres)
2338             {
2339                 fprintf(debug,
2340                         "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2341                         *enercorr, spres, svir);
2342             }
2343             else
2344             {
2345                 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2346             }
2347         }
2348
2349         if (fr->efep != efepNO)
2350         {
2351             *dvdlcorr += dvdlambda;
2352         }
2353     }
2354 }
2355
2356 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2357                   t_graph *graph, rvec x[])
2358 {
2359     if (fplog)
2360     {
2361         fprintf(fplog, "Removing pbc first time\n");
2362     }
2363     calc_shifts(box, fr->shift_vec);
2364     if (graph)
2365     {
2366         mk_mshift(fplog, graph, fr->ePBC, box, x);
2367         if (gmx_debug_at)
2368         {
2369             p_graph(debug, "do_pbc_first 1", graph);
2370         }
2371         shift_self(graph, box, x);
2372         /* By doing an extra mk_mshift the molecules that are broken
2373          * because they were e.g. imported from another software
2374          * will be made whole again. Such are the healing powers
2375          * of GROMACS.
2376          */
2377         mk_mshift(fplog, graph, fr->ePBC, box, x);
2378         if (gmx_debug_at)
2379         {
2380             p_graph(debug, "do_pbc_first 2", graph);
2381         }
2382     }
2383     if (fplog)
2384     {
2385         fprintf(fplog, "Done rmpbc\n");
2386     }
2387 }
2388
2389 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2390                             const gmx_mtop_t *mtop, rvec x[],
2391                             gmx_bool bFirst)
2392 {
2393     t_graph        *graph;
2394     int             mb, as, mol;
2395     gmx_molblock_t *molb;
2396
2397     if (bFirst && fplog)
2398     {
2399         fprintf(fplog, "Removing pbc first time\n");
2400     }
2401
2402     snew(graph, 1);
2403     as = 0;
2404     for (mb = 0; mb < mtop->nmolblock; mb++)
2405     {
2406         molb = &mtop->molblock[mb];
2407         if (molb->natoms_mol == 1 ||
2408             (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2409         {
2410             /* Just one atom or charge group in the molecule, no PBC required */
2411             as += molb->nmol*molb->natoms_mol;
2412         }
2413         else
2414         {
2415             /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2416             mk_graph_ilist(nullptr, mtop->moltype[molb->type].ilist,
2417                            0, molb->natoms_mol, FALSE, FALSE, graph);
2418
2419             for (mol = 0; mol < molb->nmol; mol++)
2420             {
2421                 mk_mshift(fplog, graph, ePBC, box, x+as);
2422
2423                 shift_self(graph, box, x+as);
2424                 /* The molecule is whole now.
2425                  * We don't need the second mk_mshift call as in do_pbc_first,
2426                  * since we no longer need this graph.
2427                  */
2428
2429                 as += molb->natoms_mol;
2430             }
2431             done_graph(graph);
2432         }
2433     }
2434     sfree(graph);
2435 }
2436
2437 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2438                        const gmx_mtop_t *mtop, rvec x[])
2439 {
2440     low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2441 }
2442
2443 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2444                  gmx_mtop_t *mtop, rvec x[])
2445 {
2446     low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2447 }
2448
2449 void put_atoms_in_box_omp(int ePBC, const matrix box, int natoms, rvec x[])
2450 {
2451     int t, nth;
2452     nth = gmx_omp_nthreads_get(emntDefault);
2453
2454 #pragma omp parallel for num_threads(nth) schedule(static)
2455     for (t = 0; t < nth; t++)
2456     {
2457         try
2458         {
2459             int offset, len;
2460
2461             offset = (natoms*t    )/nth;
2462             len    = (natoms*(t + 1))/nth - offset;
2463             put_atoms_in_box(ePBC, box, len, x + offset);
2464         }
2465         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2466     }
2467 }
2468
2469 // TODO This can be cleaned up a lot, and move back to runner.cpp
2470 void finish_run(FILE *fplog, const gmx::MDLogger &mdlog, t_commrec *cr,
2471                 t_inputrec *inputrec,
2472                 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2473                 gmx_walltime_accounting_t walltime_accounting,
2474                 nonbonded_verlet_t *nbv,
2475                 gmx_bool bWriteStat)
2476 {
2477     t_nrnb *nrnb_tot = nullptr;
2478     double  delta_t  = 0;
2479     double  nbfs     = 0, mflop = 0;
2480     double  elapsed_time,
2481             elapsed_time_over_all_ranks,
2482             elapsed_time_over_all_threads,
2483             elapsed_time_over_all_threads_over_all_ranks;
2484
2485     if (!walltime_accounting_get_valid_finish(walltime_accounting))
2486     {
2487         GMX_LOG(mdlog.warning).asParagraph().appendText("Simulation ended prematurely, no performance report will be written.");
2488         return;
2489     }
2490
2491     if (cr->nnodes > 1)
2492     {
2493         snew(nrnb_tot, 1);
2494 #if GMX_MPI
2495         MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2496                       cr->mpi_comm_mysim);
2497 #endif
2498     }
2499     else
2500     {
2501         nrnb_tot = nrnb;
2502     }
2503
2504     elapsed_time                                 = walltime_accounting_get_elapsed_time(walltime_accounting);
2505     elapsed_time_over_all_ranks                  = elapsed_time;
2506     elapsed_time_over_all_threads                = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2507     elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2508 #if GMX_MPI
2509     if (cr->nnodes > 1)
2510     {
2511         /* reduce elapsed_time over all MPI ranks in the current simulation */
2512         MPI_Allreduce(&elapsed_time,
2513                       &elapsed_time_over_all_ranks,
2514                       1, MPI_DOUBLE, MPI_SUM,
2515                       cr->mpi_comm_mysim);
2516         elapsed_time_over_all_ranks /= cr->nnodes;
2517         /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2518          * current simulation. */
2519         MPI_Allreduce(&elapsed_time_over_all_threads,
2520                       &elapsed_time_over_all_threads_over_all_ranks,
2521                       1, MPI_DOUBLE, MPI_SUM,
2522                       cr->mpi_comm_mysim);
2523     }
2524 #endif
2525
2526     if (SIMMASTER(cr))
2527     {
2528         print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2529     }
2530     if (cr->nnodes > 1)
2531     {
2532         sfree(nrnb_tot);
2533     }
2534
2535     if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2536     {
2537         print_dd_statistics(cr, inputrec, fplog);
2538     }
2539
2540     /* TODO Move the responsibility for any scaling by thread counts
2541      * to the code that handled the thread region, so that there's a
2542      * mechanism to keep cycle counting working during the transition
2543      * to task parallelism. */
2544     int nthreads_pp  = gmx_omp_nthreads_get(emntNonbonded);
2545     int nthreads_pme = gmx_omp_nthreads_get(emntPME);
2546     wallcycle_scale_by_num_threads(wcycle, cr->duty == DUTY_PME, nthreads_pp, nthreads_pme);
2547     auto cycle_sum(wallcycle_sum(cr, wcycle));
2548
2549     if (SIMMASTER(cr))
2550     {
2551         struct gmx_wallclock_gpu_t* gputimes = use_GPU(nbv) ? nbnxn_gpu_get_timings(nbv->gpu_nbv) : nullptr;
2552
2553         wallcycle_print(fplog, mdlog, cr->nnodes, cr->npmenodes, nthreads_pp, nthreads_pme,
2554                         elapsed_time_over_all_ranks,
2555                         wcycle, cycle_sum, gputimes);
2556
2557         if (EI_DYNAMICS(inputrec->eI))
2558         {
2559             delta_t = inputrec->delta_t;
2560         }
2561
2562         if (fplog)
2563         {
2564             print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2565                        elapsed_time_over_all_ranks,
2566                        walltime_accounting_get_nsteps_done(walltime_accounting),
2567                        delta_t, nbfs, mflop);
2568         }
2569         if (bWriteStat)
2570         {
2571             print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2572                        elapsed_time_over_all_ranks,
2573                        walltime_accounting_get_nsteps_done(walltime_accounting),
2574                        delta_t, nbfs, mflop);
2575         }
2576     }
2577 }
2578
2579 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, gmx::ArrayRef<real> lambda, double *lam0)
2580 {
2581     /* this function works, but could probably use a logic rewrite to keep all the different
2582        types of efep straight. */
2583
2584     if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2585     {
2586         return;
2587     }
2588
2589     t_lambda *fep = ir->fepvals;
2590     *fep_state    = fep->init_fep_state; /* this might overwrite the checkpoint
2591                                             if checkpoint is set -- a kludge is in for now
2592                                             to prevent this.*/
2593
2594     for (int i = 0; i < efptNR; i++)
2595     {
2596         /* overwrite lambda state with init_lambda for now for backwards compatibility */
2597         if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2598         {
2599             lambda[i] = fep->init_lambda;
2600             if (lam0)
2601             {
2602                 lam0[i] = lambda[i];
2603             }
2604         }
2605         else
2606         {
2607             lambda[i] = fep->all_lambda[i][*fep_state];
2608             if (lam0)
2609             {
2610                 lam0[i] = lambda[i];
2611             }
2612         }
2613     }
2614     if (ir->bSimTemp)
2615     {
2616         /* need to rescale control temperatures to match current state */
2617         for (int i = 0; i < ir->opts.ngtc; i++)
2618         {
2619             if (ir->opts.ref_t[i] > 0)
2620             {
2621                 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2622             }
2623         }
2624     }
2625
2626     /* Send to the log the information on the current lambdas */
2627     if (fplog != nullptr)
2628     {
2629         fprintf(fplog, "Initial vector of lambda components:[ ");
2630         for (const auto &l : lambda)
2631         {
2632             fprintf(fplog, "%10.4f ", l);
2633         }
2634         fprintf(fplog, "]\n");
2635     }
2636     return;
2637 }
2638
2639
2640 void init_md(FILE *fplog,
2641              t_commrec *cr, t_inputrec *ir, const gmx_output_env_t *oenv,
2642              double *t, double *t0,
2643              gmx::ArrayRef<real> lambda, int *fep_state, double *lam0,
2644              t_nrnb *nrnb, gmx_mtop_t *mtop,
2645              gmx_update_t **upd,
2646              int nfile, const t_filenm fnm[],
2647              gmx_mdoutf_t *outf, t_mdebin **mdebin,
2648              tensor force_vir, tensor shake_vir, rvec mu_tot,
2649              gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags,
2650              gmx_wallcycle_t wcycle)
2651 {
2652     int  i;
2653
2654     /* Initial values */
2655     *t = *t0       = ir->init_t;
2656
2657     *bSimAnn = FALSE;
2658     for (i = 0; i < ir->opts.ngtc; i++)
2659     {
2660         /* set bSimAnn if any group is being annealed */
2661         if (ir->opts.annealing[i] != eannNO)
2662         {
2663             *bSimAnn = TRUE;
2664         }
2665     }
2666
2667     /* Initialize lambda variables */
2668     initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2669
2670     // TODO upd is never NULL in practice, but the analysers don't know that
2671     if (upd)
2672     {
2673         *upd = init_update(ir);
2674     }
2675     if (*bSimAnn)
2676     {
2677         update_annealing_target_temp(ir, ir->init_t, upd ? *upd : nullptr);
2678     }
2679
2680     if (vcm != nullptr)
2681     {
2682         *vcm = init_vcm(fplog, &mtop->groups, ir);
2683     }
2684
2685     if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2686     {
2687         if (ir->etc == etcBERENDSEN)
2688         {
2689             please_cite(fplog, "Berendsen84a");
2690         }
2691         if (ir->etc == etcVRESCALE)
2692         {
2693             please_cite(fplog, "Bussi2007a");
2694         }
2695         if (ir->eI == eiSD1)
2696         {
2697             please_cite(fplog, "Goga2012");
2698         }
2699     }
2700     init_nrnb(nrnb);
2701
2702     if (nfile != -1)
2703     {
2704         *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv, wcycle);
2705
2706         *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? nullptr : mdoutf_get_fp_ene(*outf),
2707                               mtop, ir, mdoutf_get_fp_dhdl(*outf));
2708     }
2709
2710     /* Initiate variables */
2711     clear_mat(force_vir);
2712     clear_mat(shake_vir);
2713     clear_rvec(mu_tot);
2714 }