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