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