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