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