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