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