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