2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016, 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.
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.
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.
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.
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.
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.
39 #include "gromacs/legacyheaders/sim_util.h"
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/essentialdynamics/edsam.h"
50 #include "gromacs/ewald/pme.h"
51 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
52 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
53 #include "gromacs/imd/imd.h"
54 #include "gromacs/legacyheaders/calcmu.h"
55 #include "gromacs/legacyheaders/chargegroup.h"
56 #include "gromacs/legacyheaders/constr.h"
57 #include "gromacs/legacyheaders/copyrite.h"
58 #include "gromacs/legacyheaders/disre.h"
59 #include "gromacs/legacyheaders/force.h"
60 #include "gromacs/legacyheaders/genborn.h"
61 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
62 #include "gromacs/legacyheaders/mdatoms.h"
63 #include "gromacs/legacyheaders/mdrun.h"
64 #include "gromacs/legacyheaders/names.h"
65 #include "gromacs/legacyheaders/network.h"
66 #include "gromacs/legacyheaders/nonbonded.h"
67 #include "gromacs/legacyheaders/nrnb.h"
68 #include "gromacs/legacyheaders/orires.h"
69 #include "gromacs/legacyheaders/qmmm.h"
70 #include "gromacs/legacyheaders/txtdump.h"
71 #include "gromacs/legacyheaders/typedefs.h"
72 #include "gromacs/legacyheaders/update.h"
73 #include "gromacs/legacyheaders/types/commrec.h"
74 #include "gromacs/listed-forces/bonded.h"
75 #include "gromacs/math/units.h"
76 #include "gromacs/math/vec.h"
77 #include "gromacs/mdlib/nb_verlet.h"
78 #include "gromacs/mdlib/nbnxn_atomdata.h"
79 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
80 #include "gromacs/mdlib/nbnxn_search.h"
81 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
82 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref.h"
83 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
84 #include "gromacs/mdlib/nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
85 #include "gromacs/pbcutil/ishift.h"
86 #include "gromacs/pbcutil/mshift.h"
87 #include "gromacs/pbcutil/pbc.h"
88 #include "gromacs/pulling/pull.h"
89 #include "gromacs/pulling/pull_rotation.h"
90 #include "gromacs/timing/gpu_timing.h"
91 #include "gromacs/timing/wallcycle.h"
92 #include "gromacs/timing/walltime_accounting.h"
93 #include "gromacs/utility/cstringutil.h"
94 #include "gromacs/utility/gmxmpi.h"
95 #include "gromacs/utility/smalloc.h"
96 #include "gromacs/utility/sysinfo.h"
99 #include "nbnxn_gpu.h"
101 void print_time(FILE *out,
102 gmx_walltime_accounting_t walltime_accounting,
105 t_commrec gmx_unused *cr)
108 char timebuf[STRLEN];
109 double dt, elapsed_seconds, time_per_step;
112 #ifndef GMX_THREAD_MPI
118 fprintf(out, "step %s", gmx_step_str(step, buf));
119 if ((step >= ir->nstlist))
121 double seconds_since_epoch = gmx_gettime();
122 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
123 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
124 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
130 finish = (time_t) (seconds_since_epoch + dt);
131 gmx_ctime_r(&finish, timebuf, STRLEN);
132 sprintf(buf, "%s", timebuf);
133 buf[strlen(buf)-1] = '\0';
134 fprintf(out, ", will finish %s", buf);
138 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
143 fprintf(out, " performance: %.1f ns/day ",
144 ir->delta_t/1000*24*60*60/time_per_step);
147 #ifndef GMX_THREAD_MPI
157 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
160 char time_string[STRLEN];
169 char timebuf[STRLEN];
170 time_t temp_time = (time_t) the_time;
172 gmx_ctime_r(&temp_time, timebuf, STRLEN);
173 for (i = 0; timebuf[i] >= ' '; i++)
175 time_string[i] = timebuf[i];
177 time_string[i] = '\0';
180 fprintf(fplog, "%s on rank %d %s\n", title, nodeid, time_string);
183 void print_start(FILE *fplog, t_commrec *cr,
184 gmx_walltime_accounting_t walltime_accounting,
189 sprintf(buf, "Started %s", name);
190 print_date_and_time(fplog, cr->nodeid, buf,
191 walltime_accounting_get_start_time_stamp(walltime_accounting));
194 static void sum_forces(int start, int end, rvec f[], rvec flr[])
200 pr_rvecs(debug, 0, "fsr", f+start, end-start);
201 pr_rvecs(debug, 0, "flr", flr+start, end-start);
203 for (i = start; (i < end); i++)
205 rvec_inc(f[i], flr[i]);
210 * calc_f_el calculates forces due to an electric field.
212 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
214 * Et[] contains the parameters for the time dependent
216 * Ex[] contains the parameters for
217 * the spatial dependent part of the field.
218 * The function should return the energy due to the electric field
219 * (if any) but for now returns 0.
222 * There can be problems with the virial.
223 * Since the field is not self-consistent this is unavoidable.
224 * For neutral molecules the virial is correct within this approximation.
225 * For neutral systems with many charged molecules the error is small.
226 * But for systems with a net charge or a few charged molecules
227 * the error can be significant when the field is high.
228 * Solution: implement a self-consistent electric field into PME.
230 static void calc_f_el(FILE *fp, int start, int homenr,
231 real charge[], rvec f[],
232 t_cosines Ex[], t_cosines Et[], double t)
238 for (m = 0; (m < DIM); m++)
245 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
249 Ext[m] = cos(Et[m].a[0]*t);
258 /* Convert the field strength from V/nm to MD-units */
259 Ext[m] *= Ex[m].a[0]*FIELDFAC;
260 for (i = start; (i < start+homenr); i++)
262 f[i][m] += charge[i]*Ext[m];
272 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
273 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
277 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
278 tensor vir_part, t_graph *graph, matrix box,
279 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
283 /* The short-range virial from surrounding boxes */
285 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
286 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
288 /* Calculate partial virial, for local atoms only, based on short range.
289 * Total virial is computed in global_stat, called from do_md
291 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
292 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
294 /* Add position restraint contribution */
295 for (i = 0; i < DIM; i++)
297 vir_part[i][i] += fr->vir_diag_posres[i];
300 /* Add wall contribution */
301 for (i = 0; i < DIM; i++)
303 vir_part[i][ZZ] += fr->vir_wall_z[i];
308 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
312 static void pull_potential_wrapper(t_commrec *cr,
314 matrix box, rvec x[],
318 gmx_enerdata_t *enerd,
321 gmx_wallcycle_t wcycle)
326 /* Calculate the center of mass forces, this requires communication,
327 * which is why pull_potential is called close to other communication.
328 * The virial contribution is calculated directly,
329 * which is why we call pull_potential after calc_virial.
331 wallcycle_start(wcycle, ewcPULLPOT);
332 set_pbc(&pbc, ir->ePBC, box);
334 enerd->term[F_COM_PULL] +=
335 pull_potential(ir->pull_work, mdatoms, &pbc,
336 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
337 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
338 wallcycle_stop(wcycle, ewcPULLPOT);
341 static void pme_receive_force_ener(t_commrec *cr,
342 gmx_wallcycle_t wcycle,
343 gmx_enerdata_t *enerd,
346 real e_q, e_lj, dvdl_q, dvdl_lj;
347 float cycles_ppdpme, cycles_seppme;
349 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
350 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
352 /* In case of node-splitting, the PP nodes receive the long-range
353 * forces, virial and energy from the PME nodes here.
355 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
358 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
359 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
361 enerd->term[F_COUL_RECIP] += e_q;
362 enerd->term[F_LJ_RECIP] += e_lj;
363 enerd->dvdl_lin[efptCOUL] += dvdl_q;
364 enerd->dvdl_lin[efptVDW] += dvdl_lj;
368 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
370 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
373 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
374 gmx_int64_t step, real pforce, rvec *x, rvec *f)
378 char buf[STEPSTRSIZE];
381 for (i = 0; i < md->homenr; i++)
384 /* We also catch NAN, if the compiler does not optimize this away. */
385 if (fn2 >= pf2 || fn2 != fn2)
387 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
388 gmx_step_str(step, buf),
389 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
394 static void post_process_forces(t_commrec *cr,
396 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
398 matrix box, rvec x[],
403 t_forcerec *fr, gmx_vsite_t *vsite,
410 /* Spread the mesh force on virtual sites to the other particles...
411 * This is parallellized. MPI communication is performed
412 * if the constructing atoms aren't local.
414 wallcycle_start(wcycle, ewcVSITESPREAD);
415 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
416 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
418 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
419 wallcycle_stop(wcycle, ewcVSITESPREAD);
421 if (flags & GMX_FORCE_VIRIAL)
423 /* Now add the forces, this is local */
426 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
430 sum_forces(0, mdatoms->homenr,
433 if (EEL_FULL(fr->eeltype))
435 /* Add the mesh contribution to the virial */
436 m_add(vir_force, fr->vir_el_recip, vir_force);
438 if (EVDW_PME(fr->vdwtype))
440 /* Add the mesh contribution to the virial */
441 m_add(vir_force, fr->vir_lj_recip, vir_force);
445 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
450 if (fr->print_force >= 0)
452 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
456 static void do_nb_verlet(t_forcerec *fr,
457 interaction_const_t *ic,
458 gmx_enerdata_t *enerd,
459 int flags, int ilocality,
462 gmx_wallcycle_t wcycle)
464 int enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
465 nonbonded_verlet_group_t *nbvg;
466 gmx_bool bUsingGpuKernels;
468 if (!(flags & GMX_FORCE_NONBONDED))
470 /* skip non-bonded calculation */
474 nbvg = &fr->nbv->grp[ilocality];
476 /* GPU kernel launch overhead is already timed separately */
477 if (fr->cutoff_scheme != ecutsVERLET)
479 gmx_incons("Invalid cut-off scheme passed!");
482 bUsingGpuKernels = (nbvg->kernel_type == nbnxnk8x8x8_GPU);
484 if (!bUsingGpuKernels)
486 wallcycle_sub_start(wcycle, ewcsNONBONDED);
488 switch (nbvg->kernel_type)
490 case nbnxnk4x4_PlainC:
491 nbnxn_kernel_ref(&nbvg->nbl_lists,
497 enerd->grpp.ener[egCOULSR],
499 enerd->grpp.ener[egBHAMSR] :
500 enerd->grpp.ener[egLJSR]);
503 case nbnxnk4xN_SIMD_4xN:
504 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
511 enerd->grpp.ener[egCOULSR],
513 enerd->grpp.ener[egBHAMSR] :
514 enerd->grpp.ener[egLJSR]);
516 case nbnxnk4xN_SIMD_2xNN:
517 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
524 enerd->grpp.ener[egCOULSR],
526 enerd->grpp.ener[egBHAMSR] :
527 enerd->grpp.ener[egLJSR]);
530 case nbnxnk8x8x8_GPU:
531 nbnxn_gpu_launch_kernel(fr->nbv->gpu_nbv, nbvg->nbat, flags, ilocality);
534 case nbnxnk8x8x8_PlainC:
535 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
540 nbvg->nbat->out[0].f,
542 enerd->grpp.ener[egCOULSR],
544 enerd->grpp.ener[egBHAMSR] :
545 enerd->grpp.ener[egLJSR]);
549 gmx_incons("Invalid nonbonded kernel type passed!");
552 if (!bUsingGpuKernels)
554 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
557 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
559 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
561 else if ((!bUsingGpuKernels && nbvg->ewald_excl == ewaldexclAnalytical) ||
562 (bUsingGpuKernels && nbnxn_gpu_is_kernel_ewald_analytical(fr->nbv->gpu_nbv)))
564 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
568 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
570 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
571 if (flags & GMX_FORCE_ENERGY)
573 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
574 enr_nbnxn_kernel_ljc += 1;
575 enr_nbnxn_kernel_lj += 1;
578 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
579 nbvg->nbl_lists.natpair_ljq);
580 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
581 nbvg->nbl_lists.natpair_lj);
582 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
583 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
584 nbvg->nbl_lists.natpair_q);
586 if (ic->vdw_modifier == eintmodFORCESWITCH)
588 /* We add up the switch cost separately */
589 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
590 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
592 if (ic->vdw_modifier == eintmodPOTSWITCH)
594 /* We add up the switch cost separately */
595 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
596 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
598 if (ic->vdwtype == evdwPME)
600 /* We add up the LJ Ewald cost separately */
601 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
602 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
606 static void do_nb_verlet_fep(nbnxn_pairlist_set_t *nbl_lists,
613 gmx_enerdata_t *enerd,
616 gmx_wallcycle_t wcycle)
619 nb_kernel_data_t kernel_data;
621 real dvdl_nb[efptNR];
626 /* Add short-range interactions */
627 donb_flags |= GMX_NONBONDED_DO_SR;
629 /* Currently all group scheme kernels always calculate (shift-)forces */
630 if (flags & GMX_FORCE_FORCES)
632 donb_flags |= GMX_NONBONDED_DO_FORCE;
634 if (flags & GMX_FORCE_VIRIAL)
636 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
638 if (flags & GMX_FORCE_ENERGY)
640 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
642 if (flags & GMX_FORCE_DO_LR)
644 donb_flags |= GMX_NONBONDED_DO_LR;
647 kernel_data.flags = donb_flags;
648 kernel_data.lambda = lambda;
649 kernel_data.dvdl = dvdl_nb;
651 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
652 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
654 /* reset free energy components */
655 for (i = 0; i < efptNR; i++)
660 assert(gmx_omp_nthreads_get(emntNonbonded) == nbl_lists->nnbl);
662 wallcycle_sub_start(wcycle, ewcsNONBONDED);
663 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
664 for (th = 0; th < nbl_lists->nnbl; th++)
666 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
667 x, f, fr, mdatoms, &kernel_data, nrnb);
670 if (fepvals->sc_alpha != 0)
672 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
673 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
677 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
678 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
681 /* If we do foreign lambda and we have soft-core interactions
682 * we have to recalculate the (non-linear) energies contributions.
684 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
686 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
687 kernel_data.lambda = lam_i;
688 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
689 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
690 /* Note that we add to kernel_data.dvdl, but ignore the result */
692 for (i = 0; i < enerd->n_lambda; i++)
694 for (j = 0; j < efptNR; j++)
696 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
698 reset_foreign_enerdata(enerd);
699 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
700 for (th = 0; th < nbl_lists->nnbl; th++)
702 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
703 x, f, fr, mdatoms, &kernel_data, nrnb);
706 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
707 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
711 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
714 gmx_bool use_GPU(const nonbonded_verlet_t *nbv)
716 return nbv != NULL && nbv->bUseGPU;
719 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
720 t_inputrec *inputrec,
721 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
723 gmx_groups_t gmx_unused *groups,
724 matrix box, rvec x[], history_t *hist,
728 gmx_enerdata_t *enerd, t_fcdata *fcd,
729 real *lambda, t_graph *graph,
730 t_forcerec *fr, interaction_const_t *ic,
731 gmx_vsite_t *vsite, rvec mu_tot,
732 double t, FILE *field, gmx_edsam_t ed,
739 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
740 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
741 gmx_bool bDiffKernels = FALSE;
742 rvec vzero, box_diag;
743 float cycles_pme, cycles_force, cycles_wait_gpu;
744 nonbonded_verlet_t *nbv;
751 homenr = mdatoms->homenr;
753 clear_mat(vir_force);
755 if (DOMAINDECOMP(cr))
757 cg1 = cr->dd->ncg_tot;
768 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
769 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
770 bFillGrid = (bNS && bStateChanged);
771 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
772 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
773 bDoForces = (flags & GMX_FORCE_FORCES);
774 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
775 bUseGPU = fr->nbv->bUseGPU;
776 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
780 update_forcerec(fr, box);
782 if (NEED_MUTOT(*inputrec))
784 /* Calculate total (local) dipole moment in a temporary common array.
785 * This makes it possible to sum them over nodes faster.
787 calc_mu(start, homenr,
788 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
793 if (fr->ePBC != epbcNONE)
795 /* Compute shift vectors every step,
796 * because of pressure coupling or box deformation!
798 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
800 calc_shifts(box, fr->shift_vec);
805 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
806 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
808 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
810 unshift_self(graph, box, x);
814 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
815 fr->shift_vec, nbv->grp[0].nbat);
818 if (!(cr->duty & DUTY_PME))
823 /* Send particle coordinates to the pme nodes.
824 * Since this is only implemented for domain decomposition
825 * and domain decomposition does not use the graph,
826 * we do not need to worry about shifting.
831 wallcycle_start(wcycle, ewcPP_PMESENDX);
833 bBS = (inputrec->nwall == 2);
837 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
840 if (EEL_PME(fr->eeltype))
842 pme_flags |= GMX_PME_DO_COULOMB;
845 if (EVDW_PME(fr->vdwtype))
847 pme_flags |= GMX_PME_DO_LJ;
850 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
851 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
852 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
855 wallcycle_stop(wcycle, ewcPP_PMESENDX);
859 /* do gridding for pair search */
862 if (graph && bStateChanged)
864 /* Calculate intramolecular shift vectors to make molecules whole */
865 mk_mshift(fplog, graph, fr->ePBC, box, x);
869 box_diag[XX] = box[XX][XX];
870 box_diag[YY] = box[YY][YY];
871 box_diag[ZZ] = box[ZZ][ZZ];
873 wallcycle_start(wcycle, ewcNS);
876 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
877 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
879 0, mdatoms->homenr, -1, fr->cginfo, x,
881 nbv->grp[eintLocal].kernel_type,
882 nbv->grp[eintLocal].nbat);
883 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
887 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
888 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
890 nbv->grp[eintNonlocal].kernel_type,
891 nbv->grp[eintNonlocal].nbat);
892 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
895 if (nbv->ngrp == 1 ||
896 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
898 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
899 nbv->nbs, mdatoms, fr->cginfo);
903 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
904 nbv->nbs, mdatoms, fr->cginfo);
905 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
906 nbv->nbs, mdatoms, fr->cginfo);
908 wallcycle_stop(wcycle, ewcNS);
911 /* initialize the GPU atom data and copy shift vector */
916 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
917 nbnxn_gpu_init_atomdata(nbv->gpu_nbv, nbv->grp[eintLocal].nbat);
918 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
921 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
922 nbnxn_gpu_upload_shiftvec(nbv->gpu_nbv, nbv->grp[eintLocal].nbat);
923 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
926 /* do local pair search */
929 wallcycle_start_nocount(wcycle, ewcNS);
930 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
931 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
934 nbv->min_ci_balanced,
935 &nbv->grp[eintLocal].nbl_lists,
937 nbv->grp[eintLocal].kernel_type,
939 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
943 /* initialize local pair-list on the GPU */
944 nbnxn_gpu_init_pairlist(nbv->gpu_nbv,
945 nbv->grp[eintLocal].nbl_lists.nbl[0],
948 wallcycle_stop(wcycle, ewcNS);
952 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
953 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
954 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
955 nbv->grp[eintLocal].nbat);
956 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
957 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
962 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
963 /* launch local nonbonded F on GPU */
964 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
966 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
969 /* Communicate coordinates and sum dipole if necessary +
970 do non-local pair search */
971 if (DOMAINDECOMP(cr))
973 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
974 nbv->grp[eintLocal].kernel_type);
978 /* With GPU+CPU non-bonded calculations we need to copy
979 * the local coordinates to the non-local nbat struct
980 * (in CPU format) as the non-local kernel call also
981 * calculates the local - non-local interactions.
983 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
984 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
985 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
986 nbv->grp[eintNonlocal].nbat);
987 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
988 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
993 wallcycle_start_nocount(wcycle, ewcNS);
994 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
998 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
1001 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
1004 nbv->min_ci_balanced,
1005 &nbv->grp[eintNonlocal].nbl_lists,
1007 nbv->grp[eintNonlocal].kernel_type,
1010 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1012 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_GPU)
1014 /* initialize non-local pair-list on the GPU */
1015 nbnxn_gpu_init_pairlist(nbv->gpu_nbv,
1016 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1019 wallcycle_stop(wcycle, ewcNS);
1023 wallcycle_start(wcycle, ewcMOVEX);
1024 dd_move_x(cr->dd, box, x);
1025 wallcycle_stop(wcycle, ewcMOVEX);
1027 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1028 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1029 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1030 nbv->grp[eintNonlocal].nbat);
1031 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1032 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1035 if (bUseGPU && !bDiffKernels)
1037 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1038 /* launch non-local nonbonded F on GPU */
1039 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1041 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1047 /* launch D2H copy-back F */
1048 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1049 if (DOMAINDECOMP(cr) && !bDiffKernels)
1051 nbnxn_gpu_launch_cpyback(nbv->gpu_nbv, nbv->grp[eintNonlocal].nbat,
1052 flags, eatNonlocal);
1054 nbnxn_gpu_launch_cpyback(nbv->gpu_nbv, nbv->grp[eintLocal].nbat,
1056 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1059 if (bStateChanged && NEED_MUTOT(*inputrec))
1063 gmx_sumd(2*DIM, mu, cr);
1066 for (i = 0; i < 2; i++)
1068 for (j = 0; j < DIM; j++)
1070 fr->mu_tot[i][j] = mu[i*DIM + j];
1074 if (fr->efep == efepNO)
1076 copy_rvec(fr->mu_tot[0], mu_tot);
1080 for (j = 0; j < DIM; j++)
1083 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1084 lambda[efptCOUL]*fr->mu_tot[1][j];
1088 /* Reset energies */
1089 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1090 clear_rvecs(SHIFTS, fr->fshift);
1092 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1094 wallcycle_start(wcycle, ewcPPDURINGPME);
1095 dd_force_flop_start(cr->dd, nrnb);
1100 /* Enforced rotation has its own cycle counter that starts after the collective
1101 * coordinates have been communicated. It is added to ddCyclF to allow
1102 * for proper load-balancing */
1103 wallcycle_start(wcycle, ewcROT);
1104 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1105 wallcycle_stop(wcycle, ewcROT);
1108 /* Start the force cycle counter.
1109 * This counter is stopped after do_force_lowlevel.
1110 * No parallel communication should occur while this counter is running,
1111 * since that will interfere with the dynamic load balancing.
1113 wallcycle_start(wcycle, ewcFORCE);
1116 /* Reset forces for which the virial is calculated separately:
1117 * PME/Ewald forces if necessary */
1118 if (fr->bF_NoVirSum)
1120 if (flags & GMX_FORCE_VIRIAL)
1122 fr->f_novirsum = fr->f_novirsum_alloc;
1125 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1129 clear_rvecs(homenr, fr->f_novirsum+start);
1134 /* We are not calculating the pressure so we do not need
1135 * a separate array for forces that do not contribute
1142 /* Clear the short- and long-range forces */
1143 clear_rvecs(fr->natoms_force_constr, f);
1144 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1146 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1149 clear_rvec(fr->vir_diag_posres);
1152 if (inputrec->bPull && pull_have_constraint(inputrec->pull_work))
1154 clear_pull_forces(inputrec->pull_work);
1157 /* We calculate the non-bonded forces, when done on the CPU, here.
1158 * We do this before calling do_force_lowlevel, because in that
1159 * function, the listed forces are calculated before PME, which
1160 * does communication. With this order, non-bonded and listed
1161 * force calculation imbalance can be balanced out by the domain
1162 * decomposition load balancing.
1167 /* Maybe we should move this into do_force_lowlevel */
1168 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1172 if (fr->efep != efepNO)
1174 /* Calculate the local and non-local free energy interactions here.
1175 * Happens here on the CPU both with and without GPU.
1177 if (fr->nbv->grp[eintLocal].nbl_lists.nbl_fep[0]->nrj > 0)
1179 do_nb_verlet_fep(&fr->nbv->grp[eintLocal].nbl_lists,
1181 inputrec->fepvals, lambda,
1182 enerd, flags, nrnb, wcycle);
1185 if (DOMAINDECOMP(cr) &&
1186 fr->nbv->grp[eintNonlocal].nbl_lists.nbl_fep[0]->nrj > 0)
1188 do_nb_verlet_fep(&fr->nbv->grp[eintNonlocal].nbl_lists,
1190 inputrec->fepvals, lambda,
1191 enerd, flags, nrnb, wcycle);
1195 if (!bUseOrEmulGPU || bDiffKernels)
1199 if (DOMAINDECOMP(cr))
1201 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1202 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1212 aloc = eintNonlocal;
1215 /* Add all the non-bonded force to the normal force array.
1216 * This can be split into a local and a non-local part when overlapping
1217 * communication with calculation with domain decomposition.
1219 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1220 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1221 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1222 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1223 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1224 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1225 wallcycle_start_nocount(wcycle, ewcFORCE);
1227 /* if there are multiple fshift output buffers reduce them */
1228 if ((flags & GMX_FORCE_VIRIAL) &&
1229 nbv->grp[aloc].nbl_lists.nnbl > 1)
1231 /* This is not in a subcounter because it takes a
1232 negligible and constant-sized amount of time */
1233 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1238 /* update QMMMrec, if necessary */
1241 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1244 /* Compute the bonded and non-bonded energies and optionally forces */
1245 do_force_lowlevel(fr, inputrec, &(top->idef),
1246 cr, nrnb, wcycle, mdatoms,
1247 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1249 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1250 flags, &cycles_pme);
1254 if (do_per_step(step, inputrec->nstcalclr))
1256 /* Add the long range forces to the short range forces */
1257 for (i = 0; i < fr->natoms_force_constr; i++)
1259 rvec_add(fr->f_twin[i], f[i], f[i]);
1264 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1268 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1271 if (bUseOrEmulGPU && !bDiffKernels)
1273 /* wait for non-local forces (or calculate in emulation mode) */
1274 if (DOMAINDECOMP(cr))
1280 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1281 nbnxn_gpu_wait_for_gpu(nbv->gpu_nbv,
1282 nbv->grp[eintNonlocal].nbat,
1284 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1286 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1287 cycles_wait_gpu += cycles_tmp;
1288 cycles_force += cycles_tmp;
1292 wallcycle_start_nocount(wcycle, ewcFORCE);
1293 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1295 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1297 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1298 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1299 /* skip the reduction if there was no non-local work to do */
1300 if (nbv->grp[eintNonlocal].nbl_lists.nbl[0]->nsci > 0)
1302 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1303 nbv->grp[eintNonlocal].nbat, f);
1305 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1306 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1310 if (bDoForces && DOMAINDECOMP(cr))
1314 /* We are done with the CPU compute, but the GPU local non-bonded
1315 * kernel can still be running while we communicate the forces.
1316 * We start a counter here, so we can, hopefully, time the rest
1317 * of the GPU kernel execution and data transfer.
1319 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L_EST);
1322 /* Communicate the forces */
1323 wallcycle_start(wcycle, ewcMOVEF);
1324 dd_move_f(cr->dd, f, fr->fshift);
1327 /* We should not update the shift forces here,
1328 * since f_twin is already included in f.
1330 dd_move_f(cr->dd, fr->f_twin, NULL);
1332 wallcycle_stop(wcycle, ewcMOVEF);
1337 /* wait for local forces (or calculate in emulation mode) */
1340 #if defined(GMX_GPU) && !defined(GMX_USE_OPENCL)
1341 float cycles_tmp, cycles_wait_est;
1342 const float cuda_api_overhead_margin = 50000.0f; /* cycles */
1344 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1345 nbnxn_gpu_wait_for_gpu(nbv->gpu_nbv,
1346 nbv->grp[eintLocal].nbat,
1348 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1350 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1352 if (bDoForces && DOMAINDECOMP(cr))
1354 cycles_wait_est = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L_EST);
1356 if (cycles_tmp < cuda_api_overhead_margin)
1358 /* We measured few cycles, it could be that the kernel
1359 * and transfer finished earlier and there was no actual
1360 * wait time, only API call overhead.
1361 * Then the actual time could be anywhere between 0 and
1362 * cycles_wait_est. As a compromise, we use half the time.
1364 cycles_wait_est *= 0.5f;
1369 /* No force communication so we actually timed the wait */
1370 cycles_wait_est = cycles_tmp;
1372 /* Even though this is after dd_move_f, the actual task we are
1373 * waiting for runs asynchronously with dd_move_f and we usually
1374 * have nothing to balance it with, so we can and should add
1375 * the time to the force time for load balancing.
1377 cycles_force += cycles_wait_est;
1378 cycles_wait_gpu += cycles_wait_est;
1380 #elif defined(GMX_GPU) && defined(GMX_USE_OPENCL)
1382 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1383 nbnxn_gpu_wait_for_gpu(nbv->gpu_nbv,
1384 nbv->grp[eintLocal].nbat,
1386 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1388 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1391 /* now clear the GPU outputs while we finish the step on the CPU */
1392 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1393 nbnxn_gpu_clear_outputs(nbv->gpu_nbv, flags);
1394 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1398 wallcycle_start_nocount(wcycle, ewcFORCE);
1399 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1400 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1402 wallcycle_stop(wcycle, ewcFORCE);
1404 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1405 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1406 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1407 nbv->grp[eintLocal].nbat, f);
1408 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1409 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1412 if (DOMAINDECOMP(cr))
1414 dd_force_flop_stop(cr->dd, nrnb);
1417 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1420 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1427 if (IR_ELEC_FIELD(*inputrec))
1429 /* Compute forces due to electric field */
1430 calc_f_el(MASTER(cr) ? field : NULL,
1431 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1432 inputrec->ex, inputrec->et, t);
1435 /* If we have NoVirSum forces, but we do not calculate the virial,
1436 * we sum fr->f_novirsum=f later.
1438 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1440 wallcycle_start(wcycle, ewcVSITESPREAD);
1441 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1442 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1443 wallcycle_stop(wcycle, ewcVSITESPREAD);
1447 wallcycle_start(wcycle, ewcVSITESPREAD);
1448 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1450 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1451 wallcycle_stop(wcycle, ewcVSITESPREAD);
1455 if (flags & GMX_FORCE_VIRIAL)
1457 /* Calculation of the virial must be done after vsites! */
1458 calc_virial(0, mdatoms->homenr, x, f,
1459 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1463 if (inputrec->bPull && pull_have_potential(inputrec->pull_work))
1465 /* Since the COM pulling is always done mass-weighted, no forces are
1466 * applied to vsites and this call can be done after vsite spreading.
1468 pull_potential_wrapper(cr, inputrec, box, x,
1469 f, vir_force, mdatoms, enerd, lambda, t,
1473 /* Add the forces from enforced rotation potentials (if any) */
1476 wallcycle_start(wcycle, ewcROTadd);
1477 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1478 wallcycle_stop(wcycle, ewcROTadd);
1481 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1482 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1484 if (PAR(cr) && !(cr->duty & DUTY_PME))
1486 /* In case of node-splitting, the PP nodes receive the long-range
1487 * forces, virial and energy from the PME nodes here.
1489 pme_receive_force_ener(cr, wcycle, enerd, fr);
1494 post_process_forces(cr, step, nrnb, wcycle,
1495 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1499 /* Sum the potential energy terms from group contributions */
1500 sum_epot(&(enerd->grpp), enerd->term);
1503 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1504 t_inputrec *inputrec,
1505 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1506 gmx_localtop_t *top,
1507 gmx_groups_t *groups,
1508 matrix box, rvec x[], history_t *hist,
1512 gmx_enerdata_t *enerd, t_fcdata *fcd,
1513 real *lambda, t_graph *graph,
1514 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1515 double t, FILE *field, gmx_edsam_t ed,
1516 gmx_bool bBornRadii,
1522 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
1523 gmx_bool bDoLongRangeNS, bDoForces, bSepLRF;
1524 gmx_bool bDoAdressWF;
1526 float cycles_pme, cycles_force;
1529 homenr = mdatoms->homenr;
1531 clear_mat(vir_force);
1534 if (DOMAINDECOMP(cr))
1536 cg1 = cr->dd->ncg_tot;
1547 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1548 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1549 /* Should we update the long-range neighborlists at this step? */
1550 bDoLongRangeNS = fr->bTwinRange && bNS;
1551 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1552 bFillGrid = (bNS && bStateChanged);
1553 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1554 bDoForces = (flags & GMX_FORCE_FORCES);
1555 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1556 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1558 /* should probably move this to the forcerec since it doesn't change */
1559 bDoAdressWF = ((fr->adress_type != eAdressOff));
1563 update_forcerec(fr, box);
1565 if (NEED_MUTOT(*inputrec))
1567 /* Calculate total (local) dipole moment in a temporary common array.
1568 * This makes it possible to sum them over nodes faster.
1570 calc_mu(start, homenr,
1571 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1576 if (fr->ePBC != epbcNONE)
1578 /* Compute shift vectors every step,
1579 * because of pressure coupling or box deformation!
1581 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1583 calc_shifts(box, fr->shift_vec);
1588 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1589 &(top->cgs), x, fr->cg_cm);
1590 inc_nrnb(nrnb, eNR_CGCM, homenr);
1591 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1593 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1595 unshift_self(graph, box, x);
1600 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1601 inc_nrnb(nrnb, eNR_CGCM, homenr);
1604 if (bCalcCGCM && gmx_debug_at)
1606 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1610 if (!(cr->duty & DUTY_PME))
1615 /* Send particle coordinates to the pme nodes.
1616 * Since this is only implemented for domain decomposition
1617 * and domain decomposition does not use the graph,
1618 * we do not need to worry about shifting.
1623 wallcycle_start(wcycle, ewcPP_PMESENDX);
1625 bBS = (inputrec->nwall == 2);
1628 copy_mat(box, boxs);
1629 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1632 if (EEL_PME(fr->eeltype))
1634 pme_flags |= GMX_PME_DO_COULOMB;
1637 if (EVDW_PME(fr->vdwtype))
1639 pme_flags |= GMX_PME_DO_LJ;
1642 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1643 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1644 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1647 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1649 #endif /* GMX_MPI */
1651 /* Communicate coordinates and sum dipole if necessary */
1652 if (DOMAINDECOMP(cr))
1654 wallcycle_start(wcycle, ewcMOVEX);
1655 dd_move_x(cr->dd, box, x);
1656 wallcycle_stop(wcycle, ewcMOVEX);
1659 /* update adress weight beforehand */
1660 if (bStateChanged && bDoAdressWF)
1662 /* need pbc for adress weight calculation with pbc_dx */
1663 set_pbc(&pbc, inputrec->ePBC, box);
1664 if (fr->adress_site == eAdressSITEcog)
1666 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1667 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1669 else if (fr->adress_site == eAdressSITEcom)
1671 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1672 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1674 else if (fr->adress_site == eAdressSITEatomatom)
1676 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1677 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1681 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1682 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1686 if (NEED_MUTOT(*inputrec))
1693 gmx_sumd(2*DIM, mu, cr);
1695 for (i = 0; i < 2; i++)
1697 for (j = 0; j < DIM; j++)
1699 fr->mu_tot[i][j] = mu[i*DIM + j];
1703 if (fr->efep == efepNO)
1705 copy_rvec(fr->mu_tot[0], mu_tot);
1709 for (j = 0; j < DIM; j++)
1712 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1717 /* Reset energies */
1718 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1719 clear_rvecs(SHIFTS, fr->fshift);
1723 wallcycle_start(wcycle, ewcNS);
1725 if (graph && bStateChanged)
1727 /* Calculate intramolecular shift vectors to make molecules whole */
1728 mk_mshift(fplog, graph, fr->ePBC, box, x);
1731 /* Do the actual neighbour searching */
1733 groups, top, mdatoms,
1734 cr, nrnb, bFillGrid,
1737 wallcycle_stop(wcycle, ewcNS);
1740 if (inputrec->implicit_solvent && bNS)
1742 make_gb_nblist(cr, inputrec->gb_algorithm,
1743 x, box, fr, &top->idef, graph, fr->born);
1746 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1748 wallcycle_start(wcycle, ewcPPDURINGPME);
1749 dd_force_flop_start(cr->dd, nrnb);
1754 /* Enforced rotation has its own cycle counter that starts after the collective
1755 * coordinates have been communicated. It is added to ddCyclF to allow
1756 * for proper load-balancing */
1757 wallcycle_start(wcycle, ewcROT);
1758 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1759 wallcycle_stop(wcycle, ewcROT);
1762 /* Start the force cycle counter.
1763 * This counter is stopped after do_force_lowlevel.
1764 * No parallel communication should occur while this counter is running,
1765 * since that will interfere with the dynamic load balancing.
1767 wallcycle_start(wcycle, ewcFORCE);
1771 /* Reset forces for which the virial is calculated separately:
1772 * PME/Ewald forces if necessary */
1773 if (fr->bF_NoVirSum)
1775 if (flags & GMX_FORCE_VIRIAL)
1777 fr->f_novirsum = fr->f_novirsum_alloc;
1780 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1784 clear_rvecs(homenr, fr->f_novirsum+start);
1789 /* We are not calculating the pressure so we do not need
1790 * a separate array for forces that do not contribute
1797 /* Clear the short- and long-range forces */
1798 clear_rvecs(fr->natoms_force_constr, f);
1799 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1801 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1804 clear_rvec(fr->vir_diag_posres);
1806 if (inputrec->bPull && pull_have_constraint(inputrec->pull_work))
1808 clear_pull_forces(inputrec->pull_work);
1811 /* update QMMMrec, if necessary */
1814 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1817 /* Compute the bonded and non-bonded energies and optionally forces */
1818 do_force_lowlevel(fr, inputrec, &(top->idef),
1819 cr, nrnb, wcycle, mdatoms,
1820 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1822 inputrec->fepvals, lambda,
1823 graph, &(top->excls), fr->mu_tot,
1829 if (do_per_step(step, inputrec->nstcalclr))
1831 /* Add the long range forces to the short range forces */
1832 for (i = 0; i < fr->natoms_force_constr; i++)
1834 rvec_add(fr->f_twin[i], f[i], f[i]);
1839 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1843 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1846 if (DOMAINDECOMP(cr))
1848 dd_force_flop_stop(cr->dd, nrnb);
1851 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1857 if (IR_ELEC_FIELD(*inputrec))
1859 /* Compute forces due to electric field */
1860 calc_f_el(MASTER(cr) ? field : NULL,
1861 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1862 inputrec->ex, inputrec->et, t);
1865 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1867 /* Compute thermodynamic force in hybrid AdResS region */
1868 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1869 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1872 /* Communicate the forces */
1873 if (DOMAINDECOMP(cr))
1875 wallcycle_start(wcycle, ewcMOVEF);
1876 dd_move_f(cr->dd, f, fr->fshift);
1877 /* Do we need to communicate the separate force array
1878 * for terms that do not contribute to the single sum virial?
1879 * Position restraints and electric fields do not introduce
1880 * inter-cg forces, only full electrostatics methods do.
1881 * When we do not calculate the virial, fr->f_novirsum = f,
1882 * so we have already communicated these forces.
1884 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1885 (flags & GMX_FORCE_VIRIAL))
1887 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1891 /* We should not update the shift forces here,
1892 * since f_twin is already included in f.
1894 dd_move_f(cr->dd, fr->f_twin, NULL);
1896 wallcycle_stop(wcycle, ewcMOVEF);
1899 /* If we have NoVirSum forces, but we do not calculate the virial,
1900 * we sum fr->f_novirsum=f later.
1902 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1904 wallcycle_start(wcycle, ewcVSITESPREAD);
1905 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1906 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1907 wallcycle_stop(wcycle, ewcVSITESPREAD);
1911 wallcycle_start(wcycle, ewcVSITESPREAD);
1912 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1914 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1915 wallcycle_stop(wcycle, ewcVSITESPREAD);
1919 if (flags & GMX_FORCE_VIRIAL)
1921 /* Calculation of the virial must be done after vsites! */
1922 calc_virial(0, mdatoms->homenr, x, f,
1923 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1927 if (inputrec->bPull && pull_have_potential(inputrec->pull_work))
1929 pull_potential_wrapper(cr, inputrec, box, x,
1930 f, vir_force, mdatoms, enerd, lambda, t,
1934 /* Add the forces from enforced rotation potentials (if any) */
1937 wallcycle_start(wcycle, ewcROTadd);
1938 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1939 wallcycle_stop(wcycle, ewcROTadd);
1942 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1943 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1945 if (PAR(cr) && !(cr->duty & DUTY_PME))
1947 /* In case of node-splitting, the PP nodes receive the long-range
1948 * forces, virial and energy from the PME nodes here.
1950 pme_receive_force_ener(cr, wcycle, enerd, fr);
1955 post_process_forces(cr, step, nrnb, wcycle,
1956 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1960 /* Sum the potential energy terms from group contributions */
1961 sum_epot(&(enerd->grpp), enerd->term);
1964 void do_force(FILE *fplog, t_commrec *cr,
1965 t_inputrec *inputrec,
1966 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1967 gmx_localtop_t *top,
1968 gmx_groups_t *groups,
1969 matrix box, rvec x[], history_t *hist,
1973 gmx_enerdata_t *enerd, t_fcdata *fcd,
1974 real *lambda, t_graph *graph,
1976 gmx_vsite_t *vsite, rvec mu_tot,
1977 double t, FILE *field, gmx_edsam_t ed,
1978 gmx_bool bBornRadii,
1981 /* modify force flag if not doing nonbonded */
1982 if (!fr->bNonbonded)
1984 flags &= ~GMX_FORCE_NONBONDED;
1987 switch (inputrec->cutoff_scheme)
1990 do_force_cutsVERLET(fplog, cr, inputrec,
2006 do_force_cutsGROUP(fplog, cr, inputrec,
2021 gmx_incons("Invalid cut-off scheme passed!");
2026 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2027 t_inputrec *ir, t_mdatoms *md,
2028 t_state *state, t_commrec *cr, t_nrnb *nrnb,
2029 t_forcerec *fr, gmx_localtop_t *top)
2031 int i, m, start, end;
2033 real dt = ir->delta_t;
2037 snew(savex, state->natoms);
2044 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2045 start, md->homenr, end);
2047 /* Do a first constrain to reset particles... */
2048 step = ir->init_step;
2051 char buf[STEPSTRSIZE];
2052 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2053 gmx_step_str(step, buf));
2057 /* constrain the current position */
2058 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2059 ir, cr, step, 0, 1.0, md,
2060 state->x, state->x, NULL,
2061 fr->bMolPBC, state->box,
2062 state->lambda[efptBONDED], &dvdl_dum,
2063 NULL, NULL, nrnb, econqCoord);
2066 /* constrain the inital velocity, and save it */
2067 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2068 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2069 ir, cr, step, 0, 1.0, md,
2070 state->x, state->v, state->v,
2071 fr->bMolPBC, state->box,
2072 state->lambda[efptBONDED], &dvdl_dum,
2073 NULL, NULL, nrnb, econqVeloc);
2075 /* constrain the inital velocities at t-dt/2 */
2076 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2078 for (i = start; (i < end); i++)
2080 for (m = 0; (m < DIM); m++)
2082 /* Reverse the velocity */
2083 state->v[i][m] = -state->v[i][m];
2084 /* Store the position at t-dt in buf */
2085 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2088 /* Shake the positions at t=-dt with the positions at t=0
2089 * as reference coordinates.
2093 char buf[STEPSTRSIZE];
2094 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2095 gmx_step_str(step, buf));
2098 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2099 ir, cr, step, -1, 1.0, md,
2100 state->x, savex, NULL,
2101 fr->bMolPBC, state->box,
2102 state->lambda[efptBONDED], &dvdl_dum,
2103 state->v, NULL, nrnb, econqCoord);
2105 for (i = start; i < end; i++)
2107 for (m = 0; m < DIM; m++)
2109 /* Re-reverse the velocities */
2110 state->v[i][m] = -state->v[i][m];
2119 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2120 double *enerout, double *virout)
2122 double enersum, virsum;
2123 double invscale, invscale2, invscale3;
2124 double r, ea, eb, ec, pa, pb, pc, pd;
2129 invscale = 1.0/scale;
2130 invscale2 = invscale*invscale;
2131 invscale3 = invscale*invscale2;
2133 /* Following summation derived from cubic spline definition,
2134 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2135 * the cubic spline. We first calculate the negative of the
2136 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2137 * add the more standard, abrupt cutoff correction to that result,
2138 * yielding the long-range correction for a switched function. We
2139 * perform both the pressure and energy loops at the same time for
2140 * simplicity, as the computational cost is low. */
2144 /* Since the dispersion table has been scaled down a factor
2145 * 6.0 and the repulsion a factor 12.0 to compensate for the
2146 * c6/c12 parameters inside nbfp[] being scaled up (to save
2147 * flops in kernels), we need to correct for this.
2158 for (ri = rstart; ri < rend; ++ri)
2162 eb = 2.0*invscale2*r;
2166 pb = 3.0*invscale2*r;
2167 pc = 3.0*invscale*r*r;
2170 /* this "8" is from the packing in the vdwtab array - perhaps
2171 should be defined? */
2173 offset = 8*ri + offstart;
2174 y0 = vdwtab[offset];
2175 f = vdwtab[offset+1];
2176 g = vdwtab[offset+2];
2177 h = vdwtab[offset+3];
2179 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);
2180 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);
2182 *enerout = 4.0*M_PI*enersum*tabfactor;
2183 *virout = 4.0*M_PI*virsum*tabfactor;
2186 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2188 double eners[2], virs[2], enersum, virsum;
2189 double r0, rc3, rc9;
2191 real scale, *vdwtab;
2193 fr->enershiftsix = 0;
2194 fr->enershifttwelve = 0;
2195 fr->enerdiffsix = 0;
2196 fr->enerdifftwelve = 0;
2198 fr->virdifftwelve = 0;
2200 if (eDispCorr != edispcNO)
2202 for (i = 0; i < 2; i++)
2207 if ((fr->vdw_modifier == eintmodPOTSHIFT) ||
2208 (fr->vdw_modifier == eintmodPOTSWITCH) ||
2209 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2210 (fr->vdwtype == evdwSHIFT) ||
2211 (fr->vdwtype == evdwSWITCH))
2213 if (((fr->vdw_modifier == eintmodPOTSWITCH) ||
2214 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2215 (fr->vdwtype == evdwSWITCH)) && fr->rvdw_switch == 0)
2218 "With dispersion correction rvdw-switch can not be zero "
2219 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2222 scale = fr->nblists[0].table_vdw.scale;
2223 vdwtab = fr->nblists[0].table_vdw.data;
2225 /* Round the cut-offs to exact table values for precision */
2226 ri0 = static_cast<int>(floor(fr->rvdw_switch*scale));
2227 ri1 = static_cast<int>(ceil(fr->rvdw*scale));
2229 /* The code below has some support for handling force-switching, i.e.
2230 * when the force (instead of potential) is switched over a limited
2231 * region. This leads to a constant shift in the potential inside the
2232 * switching region, which we can handle by adding a constant energy
2233 * term in the force-switch case just like when we do potential-shift.
2235 * For now this is not enabled, but to keep the functionality in the
2236 * code we check separately for switch and shift. When we do force-switch
2237 * the shifting point is rvdw_switch, while it is the cutoff when we
2238 * have a classical potential-shift.
2240 * For a pure potential-shift the potential has a constant shift
2241 * all the way out to the cutoff, and that is it. For other forms
2242 * we need to calculate the constant shift up to the point where we
2243 * start modifying the potential.
2245 ri0 = (fr->vdw_modifier == eintmodPOTSHIFT) ? ri1 : ri0;
2251 if ((fr->vdw_modifier == eintmodFORCESWITCH) ||
2252 (fr->vdwtype == evdwSHIFT))
2254 /* Determine the constant energy shift below rvdw_switch.
2255 * Table has a scale factor since we have scaled it down to compensate
2256 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2258 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2259 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2261 else if (fr->vdw_modifier == eintmodPOTSHIFT)
2263 fr->enershiftsix = (real)(-1.0/(rc3*rc3));
2264 fr->enershifttwelve = (real)( 1.0/(rc9*rc3));
2267 /* Add the constant part from 0 to rvdw_switch.
2268 * This integration from 0 to rvdw_switch overcounts the number
2269 * of interactions by 1, as it also counts the self interaction.
2270 * We will correct for this later.
2272 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2273 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2275 /* Calculate the contribution in the range [r0,r1] where we
2276 * modify the potential. For a pure potential-shift modifier we will
2277 * have ri0==ri1, and there will not be any contribution here.
2279 for (i = 0; i < 2; i++)
2283 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2284 eners[i] -= enersum;
2288 /* Alright: Above we compensated by REMOVING the parts outside r0
2289 * corresponding to the ideal VdW 1/r6 and /r12 potentials.
2291 * Regardless of whether r0 is the point where we start switching,
2292 * or the cutoff where we calculated the constant shift, we include
2293 * all the parts we are missing out to infinity from r0 by
2294 * calculating the analytical dispersion correction.
2296 eners[0] += -4.0*M_PI/(3.0*rc3);
2297 eners[1] += 4.0*M_PI/(9.0*rc9);
2298 virs[0] += 8.0*M_PI/rc3;
2299 virs[1] += -16.0*M_PI/(3.0*rc9);
2301 else if (fr->vdwtype == evdwCUT ||
2302 EVDW_PME(fr->vdwtype) ||
2303 fr->vdwtype == evdwUSER)
2305 if (fr->vdwtype == evdwUSER && fplog)
2308 "WARNING: using dispersion correction with user tables\n");
2311 /* Note that with LJ-PME, the dispersion correction is multiplied
2312 * by the difference between the actual C6 and the value of C6
2313 * that would produce the combination rule.
2314 * This means the normal energy and virial difference formulas
2318 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2320 /* Contribution beyond the cut-off */
2321 eners[0] += -4.0*M_PI/(3.0*rc3);
2322 eners[1] += 4.0*M_PI/(9.0*rc9);
2323 if (fr->vdw_modifier == eintmodPOTSHIFT)
2325 /* Contribution within the cut-off */
2326 eners[0] += -4.0*M_PI/(3.0*rc3);
2327 eners[1] += 4.0*M_PI/(3.0*rc9);
2329 /* Contribution beyond the cut-off */
2330 virs[0] += 8.0*M_PI/rc3;
2331 virs[1] += -16.0*M_PI/(3.0*rc9);
2336 "Dispersion correction is not implemented for vdw-type = %s",
2337 evdw_names[fr->vdwtype]);
2340 /* When we deprecate the group kernels the code below can go too */
2341 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2343 /* Calculate self-interaction coefficient (assuming that
2344 * the reciprocal-space contribution is constant in the
2345 * region that contributes to the self-interaction).
2347 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2349 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2350 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2353 fr->enerdiffsix = eners[0];
2354 fr->enerdifftwelve = eners[1];
2355 /* The 0.5 is due to the Gromacs definition of the virial */
2356 fr->virdiffsix = 0.5*virs[0];
2357 fr->virdifftwelve = 0.5*virs[1];
2361 void calc_dispcorr(t_inputrec *ir, t_forcerec *fr,
2363 matrix box, real lambda, tensor pres, tensor virial,
2364 real *prescorr, real *enercorr, real *dvdlcorr)
2366 gmx_bool bCorrAll, bCorrPres;
2367 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2377 if (ir->eDispCorr != edispcNO)
2379 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2380 ir->eDispCorr == edispcAllEnerPres);
2381 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2382 ir->eDispCorr == edispcAllEnerPres);
2384 invvol = 1/det(box);
2387 /* Only correct for the interactions with the inserted molecule */
2388 dens = (natoms - fr->n_tpi)*invvol;
2393 dens = natoms*invvol;
2394 ninter = 0.5*natoms;
2397 if (ir->efep == efepNO)
2399 avcsix = fr->avcsix[0];
2400 avctwelve = fr->avctwelve[0];
2404 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2405 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2408 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2409 *enercorr += avcsix*enerdiff;
2411 if (ir->efep != efepNO)
2413 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2417 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2418 *enercorr += avctwelve*enerdiff;
2419 if (fr->efep != efepNO)
2421 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2427 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2428 if (ir->eDispCorr == edispcAllEnerPres)
2430 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2432 /* The factor 2 is because of the Gromacs virial definition */
2433 spres = -2.0*invvol*svir*PRESFAC;
2435 for (m = 0; m < DIM; m++)
2437 virial[m][m] += svir;
2438 pres[m][m] += spres;
2443 /* Can't currently control when it prints, for now, just print when degugging */
2448 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2454 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2455 *enercorr, spres, svir);
2459 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2463 if (fr->efep != efepNO)
2465 *dvdlcorr += dvdlambda;
2470 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2471 t_graph *graph, rvec x[])
2475 fprintf(fplog, "Removing pbc first time\n");
2477 calc_shifts(box, fr->shift_vec);
2480 mk_mshift(fplog, graph, fr->ePBC, box, x);
2483 p_graph(debug, "do_pbc_first 1", graph);
2485 shift_self(graph, box, x);
2486 /* By doing an extra mk_mshift the molecules that are broken
2487 * because they were e.g. imported from another software
2488 * will be made whole again. Such are the healing powers
2491 mk_mshift(fplog, graph, fr->ePBC, box, x);
2494 p_graph(debug, "do_pbc_first 2", graph);
2499 fprintf(fplog, "Done rmpbc\n");
2503 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2504 gmx_mtop_t *mtop, rvec x[],
2509 gmx_molblock_t *molb;
2511 if (bFirst && fplog)
2513 fprintf(fplog, "Removing pbc first time\n");
2518 for (mb = 0; mb < mtop->nmolblock; mb++)
2520 molb = &mtop->molblock[mb];
2521 if (molb->natoms_mol == 1 ||
2522 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2524 /* Just one atom or charge group in the molecule, no PBC required */
2525 as += molb->nmol*molb->natoms_mol;
2529 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2530 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2531 0, molb->natoms_mol, FALSE, FALSE, graph);
2533 for (mol = 0; mol < molb->nmol; mol++)
2535 mk_mshift(fplog, graph, ePBC, box, x+as);
2537 shift_self(graph, box, x+as);
2538 /* The molecule is whole now.
2539 * We don't need the second mk_mshift call as in do_pbc_first,
2540 * since we no longer need this graph.
2543 as += molb->natoms_mol;
2551 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2552 gmx_mtop_t *mtop, rvec x[])
2554 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2557 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2558 gmx_mtop_t *mtop, rvec x[])
2560 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2563 void finish_run(FILE *fplog, t_commrec *cr,
2564 t_inputrec *inputrec,
2565 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2566 gmx_walltime_accounting_t walltime_accounting,
2567 nonbonded_verlet_t *nbv,
2568 gmx_bool bWriteStat)
2570 t_nrnb *nrnb_tot = NULL;
2572 double nbfs = 0, mflop = 0;
2573 double elapsed_time,
2574 elapsed_time_over_all_ranks,
2575 elapsed_time_over_all_threads,
2576 elapsed_time_over_all_threads_over_all_ranks;
2577 wallcycle_sum(cr, wcycle);
2583 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2584 cr->mpi_comm_mysim);
2592 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2593 elapsed_time_over_all_ranks = elapsed_time;
2594 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2595 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2599 /* reduce elapsed_time over all MPI ranks in the current simulation */
2600 MPI_Allreduce(&elapsed_time,
2601 &elapsed_time_over_all_ranks,
2602 1, MPI_DOUBLE, MPI_SUM,
2603 cr->mpi_comm_mysim);
2604 elapsed_time_over_all_ranks /= cr->nnodes;
2605 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2606 * current simulation. */
2607 MPI_Allreduce(&elapsed_time_over_all_threads,
2608 &elapsed_time_over_all_threads_over_all_ranks,
2609 1, MPI_DOUBLE, MPI_SUM,
2610 cr->mpi_comm_mysim);
2616 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2623 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2625 print_dd_statistics(cr, inputrec, fplog);
2630 struct gmx_wallclock_gpu_t* gputimes = use_GPU(nbv) ? nbnxn_gpu_get_timings(nbv->gpu_nbv) : NULL;
2632 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2633 elapsed_time_over_all_ranks,
2636 if (EI_DYNAMICS(inputrec->eI))
2638 delta_t = inputrec->delta_t;
2643 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2644 elapsed_time_over_all_ranks,
2645 walltime_accounting_get_nsteps_done(walltime_accounting),
2646 delta_t, nbfs, mflop);
2650 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2651 elapsed_time_over_all_ranks,
2652 walltime_accounting_get_nsteps_done(walltime_accounting),
2653 delta_t, nbfs, mflop);
2658 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2660 /* this function works, but could probably use a logic rewrite to keep all the different
2661 types of efep straight. */
2664 t_lambda *fep = ir->fepvals;
2666 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2668 for (i = 0; i < efptNR; i++)
2680 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2681 if checkpoint is set -- a kludge is in for now
2683 for (i = 0; i < efptNR; i++)
2685 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2686 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2688 lambda[i] = fep->init_lambda;
2691 lam0[i] = lambda[i];
2696 lambda[i] = fep->all_lambda[i][*fep_state];
2699 lam0[i] = lambda[i];
2705 /* need to rescale control temperatures to match current state */
2706 for (i = 0; i < ir->opts.ngtc; i++)
2708 if (ir->opts.ref_t[i] > 0)
2710 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2716 /* Send to the log the information on the current lambdas */
2719 fprintf(fplog, "Initial vector of lambda components:[ ");
2720 for (i = 0; i < efptNR; i++)
2722 fprintf(fplog, "%10.4f ", lambda[i]);
2724 fprintf(fplog, "]\n");
2730 void init_md(FILE *fplog,
2731 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2732 double *t, double *t0,
2733 real *lambda, int *fep_state, double *lam0,
2734 t_nrnb *nrnb, gmx_mtop_t *mtop,
2736 int nfile, const t_filenm fnm[],
2737 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2738 tensor force_vir, tensor shake_vir, rvec mu_tot,
2739 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags,
2740 gmx_wallcycle_t wcycle)
2744 /* Initial values */
2745 *t = *t0 = ir->init_t;
2748 for (i = 0; i < ir->opts.ngtc; i++)
2750 /* set bSimAnn if any group is being annealed */
2751 if (ir->opts.annealing[i] != eannNO)
2758 update_annealing_target_temp(&(ir->opts), ir->init_t);
2761 /* Initialize lambda variables */
2762 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2766 *upd = init_update(ir);
2772 *vcm = init_vcm(fplog, &mtop->groups, ir);
2775 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2777 if (ir->etc == etcBERENDSEN)
2779 please_cite(fplog, "Berendsen84a");
2781 if (ir->etc == etcVRESCALE)
2783 please_cite(fplog, "Bussi2007a");
2785 if (ir->eI == eiSD1)
2787 please_cite(fplog, "Goga2012");
2790 if ((ir->et[XX].n > 0) || (ir->et[YY].n > 0) || (ir->et[ZZ].n > 0))
2792 please_cite(fplog, "Caleman2008a");
2798 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv, wcycle);
2800 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2801 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2806 please_cite(fplog, "Fritsch12");
2807 please_cite(fplog, "Junghans10");
2809 /* Initiate variables */
2810 clear_mat(force_vir);
2811 clear_mat(shake_vir);