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