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