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