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