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, 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);
1026 /* When we don't need the total dipole we sum it in global_stat */
1027 if (bStateChanged && NEED_MUTOT(*inputrec))
1029 gmx_sumd(2*DIM, mu, cr);
1031 wallcycle_stop(wcycle, ewcMOVEX);
1033 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1034 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1035 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1036 nbv->grp[eintNonlocal].nbat);
1037 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1038 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1041 if (bUseGPU && !bDiffKernels)
1043 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1044 /* launch non-local nonbonded F on GPU */
1045 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1047 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1053 /* launch D2H copy-back F */
1054 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1055 if (DOMAINDECOMP(cr) && !bDiffKernels)
1057 nbnxn_gpu_launch_cpyback(nbv->gpu_nbv, nbv->grp[eintNonlocal].nbat,
1058 flags, eatNonlocal);
1060 nbnxn_gpu_launch_cpyback(nbv->gpu_nbv, nbv->grp[eintLocal].nbat,
1062 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1065 if (bStateChanged && NEED_MUTOT(*inputrec))
1069 gmx_sumd(2*DIM, mu, cr);
1072 for (i = 0; i < 2; i++)
1074 for (j = 0; j < DIM; j++)
1076 fr->mu_tot[i][j] = mu[i*DIM + j];
1080 if (fr->efep == efepNO)
1082 copy_rvec(fr->mu_tot[0], mu_tot);
1086 for (j = 0; j < DIM; j++)
1089 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1090 lambda[efptCOUL]*fr->mu_tot[1][j];
1094 /* Reset energies */
1095 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1096 clear_rvecs(SHIFTS, fr->fshift);
1098 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1100 wallcycle_start(wcycle, ewcPPDURINGPME);
1101 dd_force_flop_start(cr->dd, nrnb);
1106 /* Enforced rotation has its own cycle counter that starts after the collective
1107 * coordinates have been communicated. It is added to ddCyclF to allow
1108 * for proper load-balancing */
1109 wallcycle_start(wcycle, ewcROT);
1110 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1111 wallcycle_stop(wcycle, ewcROT);
1114 /* Start the force cycle counter.
1115 * This counter is stopped after do_force_lowlevel.
1116 * No parallel communication should occur while this counter is running,
1117 * since that will interfere with the dynamic load balancing.
1119 wallcycle_start(wcycle, ewcFORCE);
1122 /* Reset forces for which the virial is calculated separately:
1123 * PME/Ewald forces if necessary */
1124 if (fr->bF_NoVirSum)
1126 if (flags & GMX_FORCE_VIRIAL)
1128 fr->f_novirsum = fr->f_novirsum_alloc;
1131 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1135 clear_rvecs(homenr, fr->f_novirsum+start);
1140 /* We are not calculating the pressure so we do not need
1141 * a separate array for forces that do not contribute
1148 /* Clear the short- and long-range forces */
1149 clear_rvecs(fr->natoms_force_constr, f);
1150 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1152 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1155 clear_rvec(fr->vir_diag_posres);
1158 if (inputrec->bPull && pull_have_constraint(inputrec->pull_work))
1160 clear_pull_forces(inputrec->pull_work);
1163 /* We calculate the non-bonded forces, when done on the CPU, here.
1164 * We do this before calling do_force_lowlevel, because in that
1165 * function, the listed forces are calculated before PME, which
1166 * does communication. With this order, non-bonded and listed
1167 * force calculation imbalance can be balanced out by the domain
1168 * decomposition load balancing.
1173 /* Maybe we should move this into do_force_lowlevel */
1174 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1178 if (fr->efep != efepNO)
1180 /* Calculate the local and non-local free energy interactions here.
1181 * Happens here on the CPU both with and without GPU.
1183 if (fr->nbv->grp[eintLocal].nbl_lists.nbl_fep[0]->nrj > 0)
1185 do_nb_verlet_fep(&fr->nbv->grp[eintLocal].nbl_lists,
1187 inputrec->fepvals, lambda,
1188 enerd, flags, nrnb, wcycle);
1191 if (DOMAINDECOMP(cr) &&
1192 fr->nbv->grp[eintNonlocal].nbl_lists.nbl_fep[0]->nrj > 0)
1194 do_nb_verlet_fep(&fr->nbv->grp[eintNonlocal].nbl_lists,
1196 inputrec->fepvals, lambda,
1197 enerd, flags, nrnb, wcycle);
1201 if (!bUseOrEmulGPU || bDiffKernels)
1205 if (DOMAINDECOMP(cr))
1207 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1208 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1218 aloc = eintNonlocal;
1221 /* Add all the non-bonded force to the normal force array.
1222 * This can be split into a local and a non-local part when overlapping
1223 * communication with calculation with domain decomposition.
1225 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1226 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1227 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1228 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1229 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1230 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1231 wallcycle_start_nocount(wcycle, ewcFORCE);
1233 /* if there are multiple fshift output buffers reduce them */
1234 if ((flags & GMX_FORCE_VIRIAL) &&
1235 nbv->grp[aloc].nbl_lists.nnbl > 1)
1237 /* This is not in a subcounter because it takes a
1238 negligible and constant-sized amount of time */
1239 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1244 /* update QMMMrec, if necessary */
1247 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1250 /* Compute the bonded and non-bonded energies and optionally forces */
1251 do_force_lowlevel(fr, inputrec, &(top->idef),
1252 cr, nrnb, wcycle, mdatoms,
1253 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1255 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1256 flags, &cycles_pme);
1260 if (do_per_step(step, inputrec->nstcalclr))
1262 /* Add the long range forces to the short range forces */
1263 for (i = 0; i < fr->natoms_force_constr; i++)
1265 rvec_add(fr->f_twin[i], f[i], f[i]);
1270 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1274 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1277 if (bUseOrEmulGPU && !bDiffKernels)
1279 /* wait for non-local forces (or calculate in emulation mode) */
1280 if (DOMAINDECOMP(cr))
1286 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1287 nbnxn_gpu_wait_for_gpu(nbv->gpu_nbv,
1288 nbv->grp[eintNonlocal].nbat,
1290 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1292 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1293 cycles_wait_gpu += cycles_tmp;
1294 cycles_force += cycles_tmp;
1298 wallcycle_start_nocount(wcycle, ewcFORCE);
1299 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1301 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1303 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1304 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1305 /* skip the reduction if there was no non-local work to do */
1306 if (nbv->grp[eintNonlocal].nbl_lists.nbl[0]->nsci > 0)
1308 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1309 nbv->grp[eintNonlocal].nbat, f);
1311 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1312 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1316 if (bDoForces && DOMAINDECOMP(cr))
1320 /* We are done with the CPU compute, but the GPU local non-bonded
1321 * kernel can still be running while we communicate the forces.
1322 * We start a counter here, so we can, hopefully, time the rest
1323 * of the GPU kernel execution and data transfer.
1325 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L_EST);
1328 /* Communicate the forces */
1329 wallcycle_start(wcycle, ewcMOVEF);
1330 dd_move_f(cr->dd, f, fr->fshift);
1333 /* We should not update the shift forces here,
1334 * since f_twin is already included in f.
1336 dd_move_f(cr->dd, fr->f_twin, NULL);
1338 wallcycle_stop(wcycle, ewcMOVEF);
1343 /* wait for local forces (or calculate in emulation mode) */
1346 #if defined(GMX_GPU) && !defined(GMX_USE_OPENCL)
1347 float cycles_tmp, cycles_wait_est;
1348 const float cuda_api_overhead_margin = 50000.0f; /* cycles */
1350 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1351 nbnxn_gpu_wait_for_gpu(nbv->gpu_nbv,
1352 nbv->grp[eintLocal].nbat,
1354 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1356 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1358 if (bDoForces && DOMAINDECOMP(cr))
1360 cycles_wait_est = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L_EST);
1362 if (cycles_tmp < cuda_api_overhead_margin)
1364 /* We measured few cycles, it could be that the kernel
1365 * and transfer finished earlier and there was no actual
1366 * wait time, only API call overhead.
1367 * Then the actual time could be anywhere between 0 and
1368 * cycles_wait_est. As a compromise, we use half the time.
1370 cycles_wait_est *= 0.5f;
1375 /* No force communication so we actually timed the wait */
1376 cycles_wait_est = cycles_tmp;
1378 /* Even though this is after dd_move_f, the actual task we are
1379 * waiting for runs asynchronously with dd_move_f and we usually
1380 * have nothing to balance it with, so we can and should add
1381 * the time to the force time for load balancing.
1383 cycles_force += cycles_wait_est;
1384 cycles_wait_gpu += cycles_wait_est;
1386 #elif defined(GMX_GPU) && defined(GMX_USE_OPENCL)
1388 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1389 nbnxn_gpu_wait_for_gpu(nbv->gpu_nbv,
1390 nbv->grp[eintLocal].nbat,
1392 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1394 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1397 /* now clear the GPU outputs while we finish the step on the CPU */
1398 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1399 nbnxn_gpu_clear_outputs(nbv->gpu_nbv, flags);
1400 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1404 wallcycle_start_nocount(wcycle, ewcFORCE);
1405 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1406 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1408 wallcycle_stop(wcycle, ewcFORCE);
1410 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1411 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1412 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1413 nbv->grp[eintLocal].nbat, f);
1414 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1415 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1418 if (DOMAINDECOMP(cr))
1420 dd_force_flop_stop(cr->dd, nrnb);
1423 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1426 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1433 if (IR_ELEC_FIELD(*inputrec))
1435 /* Compute forces due to electric field */
1436 calc_f_el(MASTER(cr) ? field : NULL,
1437 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1438 inputrec->ex, inputrec->et, t);
1441 /* If we have NoVirSum forces, but we do not calculate the virial,
1442 * we sum fr->f_novirsum=f later.
1444 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1446 wallcycle_start(wcycle, ewcVSITESPREAD);
1447 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1448 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1449 wallcycle_stop(wcycle, ewcVSITESPREAD);
1453 wallcycle_start(wcycle, ewcVSITESPREAD);
1454 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1456 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1457 wallcycle_stop(wcycle, ewcVSITESPREAD);
1461 if (flags & GMX_FORCE_VIRIAL)
1463 /* Calculation of the virial must be done after vsites! */
1464 calc_virial(0, mdatoms->homenr, x, f,
1465 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1469 if (inputrec->bPull && pull_have_potential(inputrec->pull_work))
1471 /* Since the COM pulling is always done mass-weighted, no forces are
1472 * applied to vsites and this call can be done after vsite spreading.
1474 pull_potential_wrapper(cr, inputrec, box, x,
1475 f, vir_force, mdatoms, enerd, lambda, t,
1479 /* Add the forces from enforced rotation potentials (if any) */
1482 wallcycle_start(wcycle, ewcROTadd);
1483 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1484 wallcycle_stop(wcycle, ewcROTadd);
1487 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1488 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1490 if (PAR(cr) && !(cr->duty & DUTY_PME))
1492 /* In case of node-splitting, the PP nodes receive the long-range
1493 * forces, virial and energy from the PME nodes here.
1495 pme_receive_force_ener(cr, wcycle, enerd, fr);
1500 post_process_forces(cr, step, nrnb, wcycle,
1501 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1505 /* Sum the potential energy terms from group contributions */
1506 sum_epot(&(enerd->grpp), enerd->term);
1509 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1510 t_inputrec *inputrec,
1511 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1512 gmx_localtop_t *top,
1513 gmx_groups_t *groups,
1514 matrix box, rvec x[], history_t *hist,
1518 gmx_enerdata_t *enerd, t_fcdata *fcd,
1519 real *lambda, t_graph *graph,
1520 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1521 double t, FILE *field, gmx_edsam_t ed,
1522 gmx_bool bBornRadii,
1528 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
1529 gmx_bool bDoLongRangeNS, bDoForces, bSepLRF;
1530 gmx_bool bDoAdressWF;
1532 float cycles_pme, cycles_force;
1535 homenr = mdatoms->homenr;
1537 clear_mat(vir_force);
1540 if (DOMAINDECOMP(cr))
1542 cg1 = cr->dd->ncg_tot;
1553 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1554 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1555 /* Should we update the long-range neighborlists at this step? */
1556 bDoLongRangeNS = fr->bTwinRange && bNS;
1557 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1558 bFillGrid = (bNS && bStateChanged);
1559 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1560 bDoForces = (flags & GMX_FORCE_FORCES);
1561 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1562 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1564 /* should probably move this to the forcerec since it doesn't change */
1565 bDoAdressWF = ((fr->adress_type != eAdressOff));
1569 update_forcerec(fr, box);
1571 if (NEED_MUTOT(*inputrec))
1573 /* Calculate total (local) dipole moment in a temporary common array.
1574 * This makes it possible to sum them over nodes faster.
1576 calc_mu(start, homenr,
1577 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1582 if (fr->ePBC != epbcNONE)
1584 /* Compute shift vectors every step,
1585 * because of pressure coupling or box deformation!
1587 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1589 calc_shifts(box, fr->shift_vec);
1594 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1595 &(top->cgs), x, fr->cg_cm);
1596 inc_nrnb(nrnb, eNR_CGCM, homenr);
1597 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1599 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1601 unshift_self(graph, box, x);
1606 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1607 inc_nrnb(nrnb, eNR_CGCM, homenr);
1610 if (bCalcCGCM && gmx_debug_at)
1612 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1616 if (!(cr->duty & DUTY_PME))
1621 /* Send particle coordinates to the pme nodes.
1622 * Since this is only implemented for domain decomposition
1623 * and domain decomposition does not use the graph,
1624 * we do not need to worry about shifting.
1629 wallcycle_start(wcycle, ewcPP_PMESENDX);
1631 bBS = (inputrec->nwall == 2);
1634 copy_mat(box, boxs);
1635 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1638 if (EEL_PME(fr->eeltype))
1640 pme_flags |= GMX_PME_DO_COULOMB;
1643 if (EVDW_PME(fr->vdwtype))
1645 pme_flags |= GMX_PME_DO_LJ;
1648 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1649 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1650 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1653 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1655 #endif /* GMX_MPI */
1657 /* Communicate coordinates and sum dipole if necessary */
1658 if (DOMAINDECOMP(cr))
1660 wallcycle_start(wcycle, ewcMOVEX);
1661 dd_move_x(cr->dd, box, x);
1662 wallcycle_stop(wcycle, ewcMOVEX);
1665 /* update adress weight beforehand */
1666 if (bStateChanged && bDoAdressWF)
1668 /* need pbc for adress weight calculation with pbc_dx */
1669 set_pbc(&pbc, inputrec->ePBC, box);
1670 if (fr->adress_site == eAdressSITEcog)
1672 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1673 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1675 else if (fr->adress_site == eAdressSITEcom)
1677 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1678 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1680 else if (fr->adress_site == eAdressSITEatomatom)
1682 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1683 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1687 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1688 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1692 if (NEED_MUTOT(*inputrec))
1699 gmx_sumd(2*DIM, mu, cr);
1701 for (i = 0; i < 2; i++)
1703 for (j = 0; j < DIM; j++)
1705 fr->mu_tot[i][j] = mu[i*DIM + j];
1709 if (fr->efep == efepNO)
1711 copy_rvec(fr->mu_tot[0], mu_tot);
1715 for (j = 0; j < DIM; j++)
1718 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1723 /* Reset energies */
1724 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1725 clear_rvecs(SHIFTS, fr->fshift);
1729 wallcycle_start(wcycle, ewcNS);
1731 if (graph && bStateChanged)
1733 /* Calculate intramolecular shift vectors to make molecules whole */
1734 mk_mshift(fplog, graph, fr->ePBC, box, x);
1737 /* Do the actual neighbour searching */
1739 groups, top, mdatoms,
1740 cr, nrnb, bFillGrid,
1743 wallcycle_stop(wcycle, ewcNS);
1746 if (inputrec->implicit_solvent && bNS)
1748 make_gb_nblist(cr, inputrec->gb_algorithm,
1749 x, box, fr, &top->idef, graph, fr->born);
1752 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1754 wallcycle_start(wcycle, ewcPPDURINGPME);
1755 dd_force_flop_start(cr->dd, nrnb);
1760 /* Enforced rotation has its own cycle counter that starts after the collective
1761 * coordinates have been communicated. It is added to ddCyclF to allow
1762 * for proper load-balancing */
1763 wallcycle_start(wcycle, ewcROT);
1764 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1765 wallcycle_stop(wcycle, ewcROT);
1768 /* Start the force cycle counter.
1769 * This counter is stopped after do_force_lowlevel.
1770 * No parallel communication should occur while this counter is running,
1771 * since that will interfere with the dynamic load balancing.
1773 wallcycle_start(wcycle, ewcFORCE);
1777 /* Reset forces for which the virial is calculated separately:
1778 * PME/Ewald forces if necessary */
1779 if (fr->bF_NoVirSum)
1781 if (flags & GMX_FORCE_VIRIAL)
1783 fr->f_novirsum = fr->f_novirsum_alloc;
1786 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1790 clear_rvecs(homenr, fr->f_novirsum+start);
1795 /* We are not calculating the pressure so we do not need
1796 * a separate array for forces that do not contribute
1803 /* Clear the short- and long-range forces */
1804 clear_rvecs(fr->natoms_force_constr, f);
1805 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1807 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1810 clear_rvec(fr->vir_diag_posres);
1812 if (inputrec->bPull && pull_have_constraint(inputrec->pull_work))
1814 clear_pull_forces(inputrec->pull_work);
1817 /* update QMMMrec, if necessary */
1820 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1823 /* Compute the bonded and non-bonded energies and optionally forces */
1824 do_force_lowlevel(fr, inputrec, &(top->idef),
1825 cr, nrnb, wcycle, mdatoms,
1826 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1828 inputrec->fepvals, lambda,
1829 graph, &(top->excls), fr->mu_tot,
1835 if (do_per_step(step, inputrec->nstcalclr))
1837 /* Add the long range forces to the short range forces */
1838 for (i = 0; i < fr->natoms_force_constr; i++)
1840 rvec_add(fr->f_twin[i], f[i], f[i]);
1845 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1849 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1852 if (DOMAINDECOMP(cr))
1854 dd_force_flop_stop(cr->dd, nrnb);
1857 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1863 if (IR_ELEC_FIELD(*inputrec))
1865 /* Compute forces due to electric field */
1866 calc_f_el(MASTER(cr) ? field : NULL,
1867 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1868 inputrec->ex, inputrec->et, t);
1871 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1873 /* Compute thermodynamic force in hybrid AdResS region */
1874 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1875 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1878 /* Communicate the forces */
1879 if (DOMAINDECOMP(cr))
1881 wallcycle_start(wcycle, ewcMOVEF);
1882 dd_move_f(cr->dd, f, fr->fshift);
1883 /* Do we need to communicate the separate force array
1884 * for terms that do not contribute to the single sum virial?
1885 * Position restraints and electric fields do not introduce
1886 * inter-cg forces, only full electrostatics methods do.
1887 * When we do not calculate the virial, fr->f_novirsum = f,
1888 * so we have already communicated these forces.
1890 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1891 (flags & GMX_FORCE_VIRIAL))
1893 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1897 /* We should not update the shift forces here,
1898 * since f_twin is already included in f.
1900 dd_move_f(cr->dd, fr->f_twin, NULL);
1902 wallcycle_stop(wcycle, ewcMOVEF);
1905 /* If we have NoVirSum forces, but we do not calculate the virial,
1906 * we sum fr->f_novirsum=f later.
1908 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1910 wallcycle_start(wcycle, ewcVSITESPREAD);
1911 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1912 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1913 wallcycle_stop(wcycle, ewcVSITESPREAD);
1917 wallcycle_start(wcycle, ewcVSITESPREAD);
1918 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1920 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1921 wallcycle_stop(wcycle, ewcVSITESPREAD);
1925 if (flags & GMX_FORCE_VIRIAL)
1927 /* Calculation of the virial must be done after vsites! */
1928 calc_virial(0, mdatoms->homenr, x, f,
1929 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1933 if (inputrec->bPull && pull_have_potential(inputrec->pull_work))
1935 pull_potential_wrapper(cr, inputrec, box, x,
1936 f, vir_force, mdatoms, enerd, lambda, t,
1940 /* Add the forces from enforced rotation potentials (if any) */
1943 wallcycle_start(wcycle, ewcROTadd);
1944 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1945 wallcycle_stop(wcycle, ewcROTadd);
1948 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1949 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1951 if (PAR(cr) && !(cr->duty & DUTY_PME))
1953 /* In case of node-splitting, the PP nodes receive the long-range
1954 * forces, virial and energy from the PME nodes here.
1956 pme_receive_force_ener(cr, wcycle, enerd, fr);
1961 post_process_forces(cr, step, nrnb, wcycle,
1962 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1966 /* Sum the potential energy terms from group contributions */
1967 sum_epot(&(enerd->grpp), enerd->term);
1970 void do_force(FILE *fplog, t_commrec *cr,
1971 t_inputrec *inputrec,
1972 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1973 gmx_localtop_t *top,
1974 gmx_groups_t *groups,
1975 matrix box, rvec x[], history_t *hist,
1979 gmx_enerdata_t *enerd, t_fcdata *fcd,
1980 real *lambda, t_graph *graph,
1982 gmx_vsite_t *vsite, rvec mu_tot,
1983 double t, FILE *field, gmx_edsam_t ed,
1984 gmx_bool bBornRadii,
1987 /* modify force flag if not doing nonbonded */
1988 if (!fr->bNonbonded)
1990 flags &= ~GMX_FORCE_NONBONDED;
1993 switch (inputrec->cutoff_scheme)
1996 do_force_cutsVERLET(fplog, cr, inputrec,
2012 do_force_cutsGROUP(fplog, cr, inputrec,
2027 gmx_incons("Invalid cut-off scheme passed!");
2032 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2033 t_inputrec *ir, t_mdatoms *md,
2034 t_state *state, t_commrec *cr, t_nrnb *nrnb,
2035 t_forcerec *fr, gmx_localtop_t *top)
2037 int i, m, start, end;
2039 real dt = ir->delta_t;
2043 snew(savex, state->natoms);
2050 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2051 start, md->homenr, end);
2053 /* Do a first constrain to reset particles... */
2054 step = ir->init_step;
2057 char buf[STEPSTRSIZE];
2058 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2059 gmx_step_str(step, buf));
2063 /* constrain the current position */
2064 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2065 ir, cr, step, 0, 1.0, md,
2066 state->x, state->x, NULL,
2067 fr->bMolPBC, state->box,
2068 state->lambda[efptBONDED], &dvdl_dum,
2069 NULL, NULL, nrnb, econqCoord);
2072 /* constrain the inital velocity, and save it */
2073 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2074 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2075 ir, cr, step, 0, 1.0, md,
2076 state->x, state->v, state->v,
2077 fr->bMolPBC, state->box,
2078 state->lambda[efptBONDED], &dvdl_dum,
2079 NULL, NULL, nrnb, econqVeloc);
2081 /* constrain the inital velocities at t-dt/2 */
2082 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2084 for (i = start; (i < end); i++)
2086 for (m = 0; (m < DIM); m++)
2088 /* Reverse the velocity */
2089 state->v[i][m] = -state->v[i][m];
2090 /* Store the position at t-dt in buf */
2091 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2094 /* Shake the positions at t=-dt with the positions at t=0
2095 * as reference coordinates.
2099 char buf[STEPSTRSIZE];
2100 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2101 gmx_step_str(step, buf));
2104 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2105 ir, cr, step, -1, 1.0, md,
2106 state->x, savex, NULL,
2107 fr->bMolPBC, state->box,
2108 state->lambda[efptBONDED], &dvdl_dum,
2109 state->v, NULL, nrnb, econqCoord);
2111 for (i = start; i < end; i++)
2113 for (m = 0; m < DIM; m++)
2115 /* Re-reverse the velocities */
2116 state->v[i][m] = -state->v[i][m];
2125 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2126 double *enerout, double *virout)
2128 double enersum, virsum;
2129 double invscale, invscale2, invscale3;
2130 double r, ea, eb, ec, pa, pb, pc, pd;
2135 invscale = 1.0/scale;
2136 invscale2 = invscale*invscale;
2137 invscale3 = invscale*invscale2;
2139 /* Following summation derived from cubic spline definition,
2140 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2141 * the cubic spline. We first calculate the negative of the
2142 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2143 * add the more standard, abrupt cutoff correction to that result,
2144 * yielding the long-range correction for a switched function. We
2145 * perform both the pressure and energy loops at the same time for
2146 * simplicity, as the computational cost is low. */
2150 /* Since the dispersion table has been scaled down a factor
2151 * 6.0 and the repulsion a factor 12.0 to compensate for the
2152 * c6/c12 parameters inside nbfp[] being scaled up (to save
2153 * flops in kernels), we need to correct for this.
2164 for (ri = rstart; ri < rend; ++ri)
2168 eb = 2.0*invscale2*r;
2172 pb = 3.0*invscale2*r;
2173 pc = 3.0*invscale*r*r;
2176 /* this "8" is from the packing in the vdwtab array - perhaps
2177 should be defined? */
2179 offset = 8*ri + offstart;
2180 y0 = vdwtab[offset];
2181 f = vdwtab[offset+1];
2182 g = vdwtab[offset+2];
2183 h = vdwtab[offset+3];
2185 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);
2186 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);
2188 *enerout = 4.0*M_PI*enersum*tabfactor;
2189 *virout = 4.0*M_PI*virsum*tabfactor;
2192 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2194 double eners[2], virs[2], enersum, virsum;
2195 double r0, rc3, rc9;
2197 real scale, *vdwtab;
2199 fr->enershiftsix = 0;
2200 fr->enershifttwelve = 0;
2201 fr->enerdiffsix = 0;
2202 fr->enerdifftwelve = 0;
2204 fr->virdifftwelve = 0;
2206 if (eDispCorr != edispcNO)
2208 for (i = 0; i < 2; i++)
2213 if ((fr->vdw_modifier == eintmodPOTSHIFT) ||
2214 (fr->vdw_modifier == eintmodPOTSWITCH) ||
2215 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2216 (fr->vdwtype == evdwSHIFT) ||
2217 (fr->vdwtype == evdwSWITCH))
2219 if (((fr->vdw_modifier == eintmodPOTSWITCH) ||
2220 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2221 (fr->vdwtype == evdwSWITCH)) && fr->rvdw_switch == 0)
2224 "With dispersion correction rvdw-switch can not be zero "
2225 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2228 scale = fr->nblists[0].table_vdw.scale;
2229 vdwtab = fr->nblists[0].table_vdw.data;
2231 /* Round the cut-offs to exact table values for precision */
2232 ri0 = static_cast<int>(floor(fr->rvdw_switch*scale));
2233 ri1 = static_cast<int>(ceil(fr->rvdw*scale));
2235 /* The code below has some support for handling force-switching, i.e.
2236 * when the force (instead of potential) is switched over a limited
2237 * region. This leads to a constant shift in the potential inside the
2238 * switching region, which we can handle by adding a constant energy
2239 * term in the force-switch case just like when we do potential-shift.
2241 * For now this is not enabled, but to keep the functionality in the
2242 * code we check separately for switch and shift. When we do force-switch
2243 * the shifting point is rvdw_switch, while it is the cutoff when we
2244 * have a classical potential-shift.
2246 * For a pure potential-shift the potential has a constant shift
2247 * all the way out to the cutoff, and that is it. For other forms
2248 * we need to calculate the constant shift up to the point where we
2249 * start modifying the potential.
2251 ri0 = (fr->vdw_modifier == eintmodPOTSHIFT) ? ri1 : ri0;
2257 if ((fr->vdw_modifier == eintmodFORCESWITCH) ||
2258 (fr->vdwtype == evdwSHIFT))
2260 /* Determine the constant energy shift below rvdw_switch.
2261 * Table has a scale factor since we have scaled it down to compensate
2262 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2264 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2265 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2267 else if (fr->vdw_modifier == eintmodPOTSHIFT)
2269 fr->enershiftsix = (real)(-1.0/(rc3*rc3));
2270 fr->enershifttwelve = (real)( 1.0/(rc9*rc3));
2273 /* Add the constant part from 0 to rvdw_switch.
2274 * This integration from 0 to rvdw_switch overcounts the number
2275 * of interactions by 1, as it also counts the self interaction.
2276 * We will correct for this later.
2278 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2279 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2281 /* Calculate the contribution in the range [r0,r1] where we
2282 * modify the potential. For a pure potential-shift modifier we will
2283 * have ri0==ri1, and there will not be any contribution here.
2285 for (i = 0; i < 2; i++)
2289 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2290 eners[i] -= enersum;
2294 /* Alright: Above we compensated by REMOVING the parts outside r0
2295 * corresponding to the ideal VdW 1/r6 and /r12 potentials.
2297 * Regardless of whether r0 is the point where we start switching,
2298 * or the cutoff where we calculated the constant shift, we include
2299 * all the parts we are missing out to infinity from r0 by
2300 * calculating the analytical dispersion correction.
2302 eners[0] += -4.0*M_PI/(3.0*rc3);
2303 eners[1] += 4.0*M_PI/(9.0*rc9);
2304 virs[0] += 8.0*M_PI/rc3;
2305 virs[1] += -16.0*M_PI/(3.0*rc9);
2307 else if (fr->vdwtype == evdwCUT ||
2308 EVDW_PME(fr->vdwtype) ||
2309 fr->vdwtype == evdwUSER)
2311 if (fr->vdwtype == evdwUSER && fplog)
2314 "WARNING: using dispersion correction with user tables\n");
2317 /* Note that with LJ-PME, the dispersion correction is multiplied
2318 * by the difference between the actual C6 and the value of C6
2319 * that would produce the combination rule.
2320 * This means the normal energy and virial difference formulas
2324 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2326 /* Contribution beyond the cut-off */
2327 eners[0] += -4.0*M_PI/(3.0*rc3);
2328 eners[1] += 4.0*M_PI/(9.0*rc9);
2329 if (fr->vdw_modifier == eintmodPOTSHIFT)
2331 /* Contribution within the cut-off */
2332 eners[0] += -4.0*M_PI/(3.0*rc3);
2333 eners[1] += 4.0*M_PI/(3.0*rc9);
2335 /* Contribution beyond the cut-off */
2336 virs[0] += 8.0*M_PI/rc3;
2337 virs[1] += -16.0*M_PI/(3.0*rc9);
2342 "Dispersion correction is not implemented for vdw-type = %s",
2343 evdw_names[fr->vdwtype]);
2346 /* When we deprecate the group kernels the code below can go too */
2347 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2349 /* Calculate self-interaction coefficient (assuming that
2350 * the reciprocal-space contribution is constant in the
2351 * region that contributes to the self-interaction).
2353 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2355 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2356 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2359 fr->enerdiffsix = eners[0];
2360 fr->enerdifftwelve = eners[1];
2361 /* The 0.5 is due to the Gromacs definition of the virial */
2362 fr->virdiffsix = 0.5*virs[0];
2363 fr->virdifftwelve = 0.5*virs[1];
2367 void calc_dispcorr(t_inputrec *ir, t_forcerec *fr,
2369 matrix box, real lambda, tensor pres, tensor virial,
2370 real *prescorr, real *enercorr, real *dvdlcorr)
2372 gmx_bool bCorrAll, bCorrPres;
2373 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2383 if (ir->eDispCorr != edispcNO)
2385 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2386 ir->eDispCorr == edispcAllEnerPres);
2387 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2388 ir->eDispCorr == edispcAllEnerPres);
2390 invvol = 1/det(box);
2393 /* Only correct for the interactions with the inserted molecule */
2394 dens = (natoms - fr->n_tpi)*invvol;
2399 dens = natoms*invvol;
2400 ninter = 0.5*natoms;
2403 if (ir->efep == efepNO)
2405 avcsix = fr->avcsix[0];
2406 avctwelve = fr->avctwelve[0];
2410 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2411 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2414 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2415 *enercorr += avcsix*enerdiff;
2417 if (ir->efep != efepNO)
2419 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2423 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2424 *enercorr += avctwelve*enerdiff;
2425 if (fr->efep != efepNO)
2427 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2433 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2434 if (ir->eDispCorr == edispcAllEnerPres)
2436 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2438 /* The factor 2 is because of the Gromacs virial definition */
2439 spres = -2.0*invvol*svir*PRESFAC;
2441 for (m = 0; m < DIM; m++)
2443 virial[m][m] += svir;
2444 pres[m][m] += spres;
2449 /* Can't currently control when it prints, for now, just print when degugging */
2454 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2460 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2461 *enercorr, spres, svir);
2465 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2469 if (fr->efep != efepNO)
2471 *dvdlcorr += dvdlambda;
2476 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2477 t_graph *graph, rvec x[])
2481 fprintf(fplog, "Removing pbc first time\n");
2483 calc_shifts(box, fr->shift_vec);
2486 mk_mshift(fplog, graph, fr->ePBC, box, x);
2489 p_graph(debug, "do_pbc_first 1", graph);
2491 shift_self(graph, box, x);
2492 /* By doing an extra mk_mshift the molecules that are broken
2493 * because they were e.g. imported from another software
2494 * will be made whole again. Such are the healing powers
2497 mk_mshift(fplog, graph, fr->ePBC, box, x);
2500 p_graph(debug, "do_pbc_first 2", graph);
2505 fprintf(fplog, "Done rmpbc\n");
2509 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2510 gmx_mtop_t *mtop, rvec x[],
2515 gmx_molblock_t *molb;
2517 if (bFirst && fplog)
2519 fprintf(fplog, "Removing pbc first time\n");
2524 for (mb = 0; mb < mtop->nmolblock; mb++)
2526 molb = &mtop->molblock[mb];
2527 if (molb->natoms_mol == 1 ||
2528 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2530 /* Just one atom or charge group in the molecule, no PBC required */
2531 as += molb->nmol*molb->natoms_mol;
2535 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2536 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2537 0, molb->natoms_mol, FALSE, FALSE, graph);
2539 for (mol = 0; mol < molb->nmol; mol++)
2541 mk_mshift(fplog, graph, ePBC, box, x+as);
2543 shift_self(graph, box, x+as);
2544 /* The molecule is whole now.
2545 * We don't need the second mk_mshift call as in do_pbc_first,
2546 * since we no longer need this graph.
2549 as += molb->natoms_mol;
2557 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2558 gmx_mtop_t *mtop, rvec x[])
2560 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2563 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2564 gmx_mtop_t *mtop, rvec x[])
2566 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2569 void finish_run(FILE *fplog, t_commrec *cr,
2570 t_inputrec *inputrec,
2571 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2572 gmx_walltime_accounting_t walltime_accounting,
2573 nonbonded_verlet_t *nbv,
2574 gmx_bool bWriteStat)
2576 t_nrnb *nrnb_tot = NULL;
2578 double nbfs = 0, mflop = 0;
2579 double elapsed_time,
2580 elapsed_time_over_all_ranks,
2581 elapsed_time_over_all_threads,
2582 elapsed_time_over_all_threads_over_all_ranks;
2583 wallcycle_sum(cr, wcycle);
2589 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2590 cr->mpi_comm_mysim);
2598 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2599 elapsed_time_over_all_ranks = elapsed_time;
2600 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2601 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2605 /* reduce elapsed_time over all MPI ranks in the current simulation */
2606 MPI_Allreduce(&elapsed_time,
2607 &elapsed_time_over_all_ranks,
2608 1, MPI_DOUBLE, MPI_SUM,
2609 cr->mpi_comm_mysim);
2610 elapsed_time_over_all_ranks /= cr->nnodes;
2611 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2612 * current simulation. */
2613 MPI_Allreduce(&elapsed_time_over_all_threads,
2614 &elapsed_time_over_all_threads_over_all_ranks,
2615 1, MPI_DOUBLE, MPI_SUM,
2616 cr->mpi_comm_mysim);
2622 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2629 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2631 print_dd_statistics(cr, inputrec, fplog);
2636 struct gmx_wallclock_gpu_t* gputimes = use_GPU(nbv) ? nbnxn_gpu_get_timings(nbv->gpu_nbv) : NULL;
2638 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2639 elapsed_time_over_all_ranks,
2642 if (EI_DYNAMICS(inputrec->eI))
2644 delta_t = inputrec->delta_t;
2649 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2650 elapsed_time_over_all_ranks,
2651 walltime_accounting_get_nsteps_done(walltime_accounting),
2652 delta_t, nbfs, mflop);
2656 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2657 elapsed_time_over_all_ranks,
2658 walltime_accounting_get_nsteps_done(walltime_accounting),
2659 delta_t, nbfs, mflop);
2664 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2666 /* this function works, but could probably use a logic rewrite to keep all the different
2667 types of efep straight. */
2670 t_lambda *fep = ir->fepvals;
2672 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2674 for (i = 0; i < efptNR; i++)
2686 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2687 if checkpoint is set -- a kludge is in for now
2689 for (i = 0; i < efptNR; i++)
2691 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2692 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2694 lambda[i] = fep->init_lambda;
2697 lam0[i] = lambda[i];
2702 lambda[i] = fep->all_lambda[i][*fep_state];
2705 lam0[i] = lambda[i];
2711 /* need to rescale control temperatures to match current state */
2712 for (i = 0; i < ir->opts.ngtc; i++)
2714 if (ir->opts.ref_t[i] > 0)
2716 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2722 /* Send to the log the information on the current lambdas */
2725 fprintf(fplog, "Initial vector of lambda components:[ ");
2726 for (i = 0; i < efptNR; i++)
2728 fprintf(fplog, "%10.4f ", lambda[i]);
2730 fprintf(fplog, "]\n");
2736 void init_md(FILE *fplog,
2737 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2738 double *t, double *t0,
2739 real *lambda, int *fep_state, double *lam0,
2740 t_nrnb *nrnb, gmx_mtop_t *mtop,
2742 int nfile, const t_filenm fnm[],
2743 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2744 tensor force_vir, tensor shake_vir, rvec mu_tot,
2745 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags,
2746 gmx_wallcycle_t wcycle)
2750 /* Initial values */
2751 *t = *t0 = ir->init_t;
2754 for (i = 0; i < ir->opts.ngtc; i++)
2756 /* set bSimAnn if any group is being annealed */
2757 if (ir->opts.annealing[i] != eannNO)
2764 update_annealing_target_temp(&(ir->opts), ir->init_t);
2767 /* Initialize lambda variables */
2768 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2772 *upd = init_update(ir);
2778 *vcm = init_vcm(fplog, &mtop->groups, ir);
2781 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2783 if (ir->etc == etcBERENDSEN)
2785 please_cite(fplog, "Berendsen84a");
2787 if (ir->etc == etcVRESCALE)
2789 please_cite(fplog, "Bussi2007a");
2791 if (ir->eI == eiSD1)
2793 please_cite(fplog, "Goga2012");
2796 if ((ir->et[XX].n > 0) || (ir->et[YY].n > 0) || (ir->et[ZZ].n > 0))
2798 please_cite(fplog, "Caleman2008a");
2804 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv, wcycle);
2806 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2807 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2812 please_cite(fplog, "Fritsch12");
2813 please_cite(fplog, "Junghans10");
2815 /* Initiate variables */
2816 clear_mat(force_vir);
2817 clear_mat(shake_vir);