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