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