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