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, 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_search.h"
80 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda.h"
81 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_data_mgmt.h"
82 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
83 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref.h"
84 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
85 #include "gromacs/mdlib/nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
86 #include "gromacs/pbcutil/ishift.h"
87 #include "gromacs/pbcutil/mshift.h"
88 #include "gromacs/pbcutil/pbc.h"
89 #include "gromacs/pulling/pull.h"
90 #include "gromacs/pulling/pull_rotation.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"
100 void print_time(FILE *out,
101 gmx_walltime_accounting_t walltime_accounting,
104 t_commrec gmx_unused *cr)
107 char timebuf[STRLEN];
108 double dt, elapsed_seconds, time_per_step;
111 #ifndef GMX_THREAD_MPI
117 fprintf(out, "step %s", gmx_step_str(step, buf));
118 if ((step >= ir->nstlist))
120 double seconds_since_epoch = gmx_gettime();
121 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
122 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
123 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
129 finish = (time_t) (seconds_since_epoch + dt);
130 gmx_ctime_r(&finish, timebuf, STRLEN);
131 sprintf(buf, "%s", timebuf);
132 buf[strlen(buf)-1] = '\0';
133 fprintf(out, ", will finish %s", buf);
137 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
142 fprintf(out, " performance: %.1f ns/day ",
143 ir->delta_t/1000*24*60*60/time_per_step);
146 #ifndef GMX_THREAD_MPI
156 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
159 char time_string[STRLEN];
168 char timebuf[STRLEN];
169 time_t temp_time = (time_t) the_time;
171 gmx_ctime_r(&temp_time, timebuf, STRLEN);
172 for (i = 0; timebuf[i] >= ' '; i++)
174 time_string[i] = timebuf[i];
176 time_string[i] = '\0';
179 fprintf(fplog, "%s on rank %d %s\n", title, nodeid, time_string);
182 void print_start(FILE *fplog, t_commrec *cr,
183 gmx_walltime_accounting_t walltime_accounting,
188 sprintf(buf, "Started %s", name);
189 print_date_and_time(fplog, cr->nodeid, buf,
190 walltime_accounting_get_start_time_stamp(walltime_accounting));
193 static void sum_forces(int start, int end, rvec f[], rvec flr[])
199 pr_rvecs(debug, 0, "fsr", f+start, end-start);
200 pr_rvecs(debug, 0, "flr", flr+start, end-start);
202 for (i = start; (i < end); i++)
204 rvec_inc(f[i], flr[i]);
209 * calc_f_el calculates forces due to an electric field.
211 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
213 * Et[] contains the parameters for the time dependent
214 * part of the field (not yet used).
215 * Ex[] contains the parameters for
216 * the spatial dependent part of the field. You can have cool periodic
217 * fields in principle, but only a constant field is supported
219 * The function should return the energy due to the electric field
220 * (if any) but for now returns 0.
223 * There can be problems with the virial.
224 * Since the field is not self-consistent this is unavoidable.
225 * For neutral molecules the virial is correct within this approximation.
226 * For neutral systems with many charged molecules the error is small.
227 * But for systems with a net charge or a few charged molecules
228 * the error can be significant when the field is high.
229 * Solution: implement a self-consitent electric field into PME.
231 static void calc_f_el(FILE *fp, int start, int homenr,
232 real charge[], rvec f[],
233 t_cosines Ex[], t_cosines Et[], double t)
239 for (m = 0; (m < DIM); m++)
246 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
250 Ext[m] = cos(Et[m].a[0]*t);
259 /* Convert the field strength from V/nm to MD-units */
260 Ext[m] *= Ex[m].a[0]*FIELDFAC;
261 for (i = start; (i < start+homenr); i++)
263 f[i][m] += charge[i]*Ext[m];
273 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
274 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
278 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
279 tensor vir_part, t_graph *graph, matrix box,
280 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
284 /* The short-range virial from surrounding boxes */
286 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
287 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
289 /* Calculate partial virial, for local atoms only, based on short range.
290 * Total virial is computed in global_stat, called from do_md
292 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
293 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
295 /* Add position restraint contribution */
296 for (i = 0; i < DIM; i++)
298 vir_part[i][i] += fr->vir_diag_posres[i];
301 /* Add wall contribution */
302 for (i = 0; i < DIM; i++)
304 vir_part[i][ZZ] += fr->vir_wall_z[i];
309 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
313 static void pull_potential_wrapper(t_commrec *cr,
315 matrix box, rvec x[],
319 gmx_enerdata_t *enerd,
322 gmx_wallcycle_t wcycle)
327 /* Calculate the center of mass forces, this requires communication,
328 * which is why pull_potential is called close to other communication.
329 * The virial contribution is calculated directly,
330 * which is why we call pull_potential after calc_virial.
332 wallcycle_start(wcycle, ewcPULLPOT);
333 set_pbc(&pbc, ir->ePBC, box);
335 enerd->term[F_COM_PULL] +=
336 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
337 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
338 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
339 wallcycle_stop(wcycle, ewcPULLPOT);
342 static void pme_receive_force_ener(t_commrec *cr,
343 gmx_wallcycle_t wcycle,
344 gmx_enerdata_t *enerd,
347 real e_q, e_lj, dvdl_q, dvdl_lj;
348 float cycles_ppdpme, cycles_seppme;
350 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
351 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
353 /* In case of node-splitting, the PP nodes receive the long-range
354 * forces, virial and energy from the PME nodes here.
356 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
359 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
360 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
362 enerd->term[F_COUL_RECIP] += e_q;
363 enerd->term[F_LJ_RECIP] += e_lj;
364 enerd->dvdl_lin[efptCOUL] += dvdl_q;
365 enerd->dvdl_lin[efptVDW] += dvdl_lj;
369 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
371 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
374 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
375 gmx_int64_t step, real pforce, rvec *x, rvec *f)
379 char buf[STEPSTRSIZE];
382 for (i = 0; i < md->homenr; i++)
385 /* We also catch NAN, if the compiler does not optimize this away. */
386 if (fn2 >= pf2 || fn2 != fn2)
388 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
389 gmx_step_str(step, buf),
390 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
395 static void post_process_forces(t_commrec *cr,
397 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
399 matrix box, rvec x[],
404 t_forcerec *fr, gmx_vsite_t *vsite,
411 /* Spread the mesh force on virtual sites to the other particles...
412 * This is parallellized. MPI communication is performed
413 * if the constructing atoms aren't local.
415 wallcycle_start(wcycle, ewcVSITESPREAD);
416 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
417 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
419 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
420 wallcycle_stop(wcycle, ewcVSITESPREAD);
422 if (flags & GMX_FORCE_VIRIAL)
424 /* Now add the forces, this is local */
427 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
431 sum_forces(0, mdatoms->homenr,
434 if (EEL_FULL(fr->eeltype))
436 /* Add the mesh contribution to the virial */
437 m_add(vir_force, fr->vir_el_recip, vir_force);
439 if (EVDW_PME(fr->vdwtype))
441 /* Add the mesh contribution to the virial */
442 m_add(vir_force, fr->vir_lj_recip, vir_force);
446 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
451 if (fr->print_force >= 0)
453 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
457 static void do_nb_verlet(t_forcerec *fr,
458 interaction_const_t *ic,
459 gmx_enerdata_t *enerd,
460 int flags, int ilocality,
463 gmx_wallcycle_t wcycle)
465 int enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
466 nonbonded_verlet_group_t *nbvg;
469 if (!(flags & GMX_FORCE_NONBONDED))
471 /* skip non-bonded calculation */
475 nbvg = &fr->nbv->grp[ilocality];
477 /* CUDA kernel launch overhead is already timed separately */
478 if (fr->cutoff_scheme != ecutsVERLET)
480 gmx_incons("Invalid cut-off scheme passed!");
483 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
487 wallcycle_sub_start(wcycle, ewcsNONBONDED);
489 switch (nbvg->kernel_type)
491 case nbnxnk4x4_PlainC:
492 nbnxn_kernel_ref(&nbvg->nbl_lists,
498 enerd->grpp.ener[egCOULSR],
500 enerd->grpp.ener[egBHAMSR] :
501 enerd->grpp.ener[egLJSR]);
504 case nbnxnk4xN_SIMD_4xN:
505 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
512 enerd->grpp.ener[egCOULSR],
514 enerd->grpp.ener[egBHAMSR] :
515 enerd->grpp.ener[egLJSR]);
517 case nbnxnk4xN_SIMD_2xNN:
518 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
525 enerd->grpp.ener[egCOULSR],
527 enerd->grpp.ener[egBHAMSR] :
528 enerd->grpp.ener[egLJSR]);
531 case nbnxnk8x8x8_CUDA:
532 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
535 case nbnxnk8x8x8_PlainC:
536 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
541 nbvg->nbat->out[0].f,
543 enerd->grpp.ener[egCOULSR],
545 enerd->grpp.ener[egBHAMSR] :
546 enerd->grpp.ener[egLJSR]);
550 gmx_incons("Invalid nonbonded kernel type passed!");
555 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
558 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
560 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
562 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
563 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
565 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
569 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
571 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
572 if (flags & GMX_FORCE_ENERGY)
574 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
575 enr_nbnxn_kernel_ljc += 1;
576 enr_nbnxn_kernel_lj += 1;
579 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
580 nbvg->nbl_lists.natpair_ljq);
581 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
582 nbvg->nbl_lists.natpair_lj);
583 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
584 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
585 nbvg->nbl_lists.natpair_q);
587 if (ic->vdw_modifier == eintmodFORCESWITCH)
589 /* We add up the switch cost separately */
590 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
591 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
593 if (ic->vdw_modifier == eintmodPOTSWITCH)
595 /* We add up the switch cost separately */
596 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
597 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
599 if (ic->vdwtype == evdwPME)
601 /* We add up the LJ Ewald cost separately */
602 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
603 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
607 static void do_nb_verlet_fep(nbnxn_pairlist_set_t *nbl_lists,
614 gmx_enerdata_t *enerd,
617 gmx_wallcycle_t wcycle)
620 nb_kernel_data_t kernel_data;
622 real dvdl_nb[efptNR];
627 /* Add short-range interactions */
628 donb_flags |= GMX_NONBONDED_DO_SR;
630 /* Currently all group scheme kernels always calculate (shift-)forces */
631 if (flags & GMX_FORCE_FORCES)
633 donb_flags |= GMX_NONBONDED_DO_FORCE;
635 if (flags & GMX_FORCE_VIRIAL)
637 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
639 if (flags & GMX_FORCE_ENERGY)
641 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
643 if (flags & GMX_FORCE_DO_LR)
645 donb_flags |= GMX_NONBONDED_DO_LR;
648 kernel_data.flags = donb_flags;
649 kernel_data.lambda = lambda;
650 kernel_data.dvdl = dvdl_nb;
652 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
653 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
655 /* reset free energy components */
656 for (i = 0; i < efptNR; i++)
661 assert(gmx_omp_nthreads_get(emntNonbonded) == nbl_lists->nnbl);
663 wallcycle_sub_start(wcycle, ewcsNONBONDED);
664 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
665 for (th = 0; th < nbl_lists->nnbl; th++)
667 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
668 x, f, fr, mdatoms, &kernel_data, nrnb);
671 if (fepvals->sc_alpha != 0)
673 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
674 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
678 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
679 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
682 /* If we do foreign lambda and we have soft-core interactions
683 * we have to recalculate the (non-linear) energies contributions.
685 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
687 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
688 kernel_data.lambda = lam_i;
689 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
690 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
691 /* Note that we add to kernel_data.dvdl, but ignore the result */
693 for (i = 0; i < enerd->n_lambda; i++)
695 for (j = 0; j < efptNR; j++)
697 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
699 reset_foreign_enerdata(enerd);
700 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
701 for (th = 0; th < nbl_lists->nnbl; th++)
703 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
704 x, f, fr, mdatoms, &kernel_data, nrnb);
707 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
708 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
712 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
715 gmx_bool use_GPU(const nonbonded_verlet_t *nbv)
717 return nbv != NULL && nbv->bUseGPU;
720 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
721 t_inputrec *inputrec,
722 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
724 gmx_groups_t gmx_unused *groups,
725 matrix box, rvec x[], history_t *hist,
729 gmx_enerdata_t *enerd, t_fcdata *fcd,
730 real *lambda, t_graph *graph,
731 t_forcerec *fr, interaction_const_t *ic,
732 gmx_vsite_t *vsite, rvec mu_tot,
733 double t, FILE *field, gmx_edsam_t ed,
740 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
741 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
742 gmx_bool bDiffKernels = FALSE;
743 rvec vzero, box_diag;
744 float cycles_pme, cycles_force, cycles_wait_gpu;
745 nonbonded_verlet_t *nbv;
752 homenr = mdatoms->homenr;
754 clear_mat(vir_force);
756 if (DOMAINDECOMP(cr))
758 cg1 = cr->dd->ncg_tot;
769 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
770 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
771 bFillGrid = (bNS && bStateChanged);
772 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
773 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
774 bDoForces = (flags & GMX_FORCE_FORCES);
775 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
776 bUseGPU = fr->nbv->bUseGPU;
777 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
781 update_forcerec(fr, box);
783 if (NEED_MUTOT(*inputrec))
785 /* Calculate total (local) dipole moment in a temporary common array.
786 * This makes it possible to sum them over nodes faster.
788 calc_mu(start, homenr,
789 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
794 if (fr->ePBC != epbcNONE)
796 /* Compute shift vectors every step,
797 * because of pressure coupling or box deformation!
799 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
801 calc_shifts(box, fr->shift_vec);
806 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
807 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
809 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
811 unshift_self(graph, box, x);
815 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
816 fr->shift_vec, nbv->grp[0].nbat);
819 if (!(cr->duty & DUTY_PME))
824 /* Send particle coordinates to the pme nodes.
825 * Since this is only implemented for domain decomposition
826 * and domain decomposition does not use the graph,
827 * we do not need to worry about shifting.
832 wallcycle_start(wcycle, ewcPP_PMESENDX);
834 bBS = (inputrec->nwall == 2);
838 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
841 if (EEL_PME(fr->eeltype))
843 pme_flags |= GMX_PME_DO_COULOMB;
846 if (EVDW_PME(fr->vdwtype))
848 pme_flags |= GMX_PME_DO_LJ;
851 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
852 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
853 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
856 wallcycle_stop(wcycle, ewcPP_PMESENDX);
860 /* do gridding for pair search */
863 if (graph && bStateChanged)
865 /* Calculate intramolecular shift vectors to make molecules whole */
866 mk_mshift(fplog, graph, fr->ePBC, box, x);
870 box_diag[XX] = box[XX][XX];
871 box_diag[YY] = box[YY][YY];
872 box_diag[ZZ] = box[ZZ][ZZ];
874 wallcycle_start(wcycle, ewcNS);
877 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
878 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
880 0, mdatoms->homenr, -1, fr->cginfo, x,
882 nbv->grp[eintLocal].kernel_type,
883 nbv->grp[eintLocal].nbat);
884 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
888 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
889 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
891 nbv->grp[eintNonlocal].kernel_type,
892 nbv->grp[eintNonlocal].nbat);
893 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
896 if (nbv->ngrp == 1 ||
897 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
899 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
900 nbv->nbs, mdatoms, fr->cginfo);
904 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
905 nbv->nbs, mdatoms, fr->cginfo);
906 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
907 nbv->nbs, mdatoms, fr->cginfo);
909 wallcycle_stop(wcycle, ewcNS);
912 /* initialize the GPU atom data and copy shift vector */
917 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
918 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
919 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
922 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
923 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
924 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
927 /* do local pair search */
930 wallcycle_start_nocount(wcycle, ewcNS);
931 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
932 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
935 nbv->min_ci_balanced,
936 &nbv->grp[eintLocal].nbl_lists,
938 nbv->grp[eintLocal].kernel_type,
940 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
944 /* initialize local pair-list on the GPU */
945 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
946 nbv->grp[eintLocal].nbl_lists.nbl[0],
949 wallcycle_stop(wcycle, ewcNS);
953 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
954 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
955 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
956 nbv->grp[eintLocal].nbat);
957 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
958 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
963 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
964 /* launch local nonbonded F on GPU */
965 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
967 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
970 /* Communicate coordinates and sum dipole if necessary +
971 do non-local pair search */
972 if (DOMAINDECOMP(cr))
974 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
975 nbv->grp[eintLocal].kernel_type);
979 /* With GPU+CPU non-bonded calculations we need to copy
980 * the local coordinates to the non-local nbat struct
981 * (in CPU format) as the non-local kernel call also
982 * calculates the local - non-local interactions.
984 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
985 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
986 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
987 nbv->grp[eintNonlocal].nbat);
988 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
989 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
994 wallcycle_start_nocount(wcycle, ewcNS);
995 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
999 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
1002 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
1005 nbv->min_ci_balanced,
1006 &nbv->grp[eintNonlocal].nbl_lists,
1008 nbv->grp[eintNonlocal].kernel_type,
1011 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1013 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1015 /* initialize non-local pair-list on the GPU */
1016 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1017 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1020 wallcycle_stop(wcycle, ewcNS);
1024 wallcycle_start(wcycle, ewcMOVEX);
1025 dd_move_x(cr->dd, box, x);
1027 /* When we don't need the total dipole we sum it in global_stat */
1028 if (bStateChanged && NEED_MUTOT(*inputrec))
1030 gmx_sumd(2*DIM, mu, cr);
1032 wallcycle_stop(wcycle, ewcMOVEX);
1034 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1035 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1036 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1037 nbv->grp[eintNonlocal].nbat);
1038 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1039 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1042 if (bUseGPU && !bDiffKernels)
1044 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1045 /* launch non-local nonbonded F on GPU */
1046 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1048 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1054 /* launch D2H copy-back F */
1055 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1056 if (DOMAINDECOMP(cr) && !bDiffKernels)
1058 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1059 flags, eatNonlocal);
1061 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1063 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1066 if (bStateChanged && NEED_MUTOT(*inputrec))
1070 gmx_sumd(2*DIM, mu, cr);
1073 for (i = 0; i < 2; i++)
1075 for (j = 0; j < DIM; j++)
1077 fr->mu_tot[i][j] = mu[i*DIM + j];
1081 if (fr->efep == efepNO)
1083 copy_rvec(fr->mu_tot[0], mu_tot);
1087 for (j = 0; j < DIM; j++)
1090 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1091 lambda[efptCOUL]*fr->mu_tot[1][j];
1095 /* Reset energies */
1096 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1097 clear_rvecs(SHIFTS, fr->fshift);
1099 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1101 wallcycle_start(wcycle, ewcPPDURINGPME);
1102 dd_force_flop_start(cr->dd, nrnb);
1107 /* Enforced rotation has its own cycle counter that starts after the collective
1108 * coordinates have been communicated. It is added to ddCyclF to allow
1109 * for proper load-balancing */
1110 wallcycle_start(wcycle, ewcROT);
1111 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1112 wallcycle_stop(wcycle, ewcROT);
1115 /* Start the force cycle counter.
1116 * This counter is stopped in do_forcelow_level.
1117 * No parallel communication should occur while this counter is running,
1118 * since that will interfere with the dynamic load balancing.
1120 wallcycle_start(wcycle, ewcFORCE);
1123 /* Reset forces for which the virial is calculated separately:
1124 * PME/Ewald forces if necessary */
1125 if (fr->bF_NoVirSum)
1127 if (flags & GMX_FORCE_VIRIAL)
1129 fr->f_novirsum = fr->f_novirsum_alloc;
1132 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1136 clear_rvecs(homenr, fr->f_novirsum+start);
1141 /* We are not calculating the pressure so we do not need
1142 * a separate array for forces that do not contribute
1149 /* Clear the short- and long-range forces */
1150 clear_rvecs(fr->natoms_force_constr, f);
1151 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1153 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1156 clear_rvec(fr->vir_diag_posres);
1159 if (inputrec->ePull == epullCONSTRAINT)
1161 clear_pull_forces(inputrec->pull);
1164 /* We calculate the non-bonded forces, when done on the CPU, here.
1165 * We do this before calling do_force_lowlevel, because in that
1166 * function, the listed forces are calculated before PME, which
1167 * does communication. With this order, non-bonded and listed
1168 * force calculation imbalance can be balanced out by the domain
1169 * decomposition load balancing.
1174 /* Maybe we should move this into do_force_lowlevel */
1175 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1179 if (fr->efep != efepNO)
1181 /* Calculate the local and non-local free energy interactions here.
1182 * Happens here on the CPU both with and without GPU.
1184 if (fr->nbv->grp[eintLocal].nbl_lists.nbl_fep[0]->nrj > 0)
1186 do_nb_verlet_fep(&fr->nbv->grp[eintLocal].nbl_lists,
1188 inputrec->fepvals, lambda,
1189 enerd, flags, nrnb, wcycle);
1192 if (DOMAINDECOMP(cr) &&
1193 fr->nbv->grp[eintNonlocal].nbl_lists.nbl_fep[0]->nrj > 0)
1195 do_nb_verlet_fep(&fr->nbv->grp[eintNonlocal].nbl_lists,
1197 inputrec->fepvals, lambda,
1198 enerd, flags, nrnb, wcycle);
1202 if (!bUseOrEmulGPU || bDiffKernels)
1206 if (DOMAINDECOMP(cr))
1208 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1209 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1219 aloc = eintNonlocal;
1222 /* Add all the non-bonded force to the normal force array.
1223 * This can be split into a local a non-local part when overlapping
1224 * communication with calculation with domain decomposition.
1226 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1227 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1228 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1229 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1230 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1231 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1232 wallcycle_start_nocount(wcycle, ewcFORCE);
1234 /* if there are multiple fshift output buffers reduce them */
1235 if ((flags & GMX_FORCE_VIRIAL) &&
1236 nbv->grp[aloc].nbl_lists.nnbl > 1)
1238 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1243 /* update QMMMrec, if necessary */
1246 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1249 /* Compute the bonded and non-bonded energies and optionally forces */
1250 do_force_lowlevel(fr, inputrec, &(top->idef),
1251 cr, nrnb, wcycle, mdatoms,
1252 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1254 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1255 flags, &cycles_pme);
1259 if (do_per_step(step, inputrec->nstcalclr))
1261 /* Add the long range forces to the short range forces */
1262 for (i = 0; i < fr->natoms_force_constr; i++)
1264 rvec_add(fr->f_twin[i], f[i], f[i]);
1269 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1273 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1276 if (bUseOrEmulGPU && !bDiffKernels)
1278 /* wait for non-local forces (or calculate in emulation mode) */
1279 if (DOMAINDECOMP(cr))
1285 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1286 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1287 nbv->grp[eintNonlocal].nbat,
1289 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1291 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1292 cycles_wait_gpu += cycles_tmp;
1293 cycles_force += cycles_tmp;
1297 wallcycle_start_nocount(wcycle, ewcFORCE);
1298 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1300 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1302 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1303 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1304 /* skip the reduction if there was no non-local work to do */
1305 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1307 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1308 nbv->grp[eintNonlocal].nbat, f);
1310 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1311 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1315 if (bDoForces && DOMAINDECOMP(cr))
1319 /* We are done with the CPU compute, but the GPU local non-bonded
1320 * kernel can still be running while we communicate the forces.
1321 * We start a counter here, so we can, hopefully, time the rest
1322 * of the GPU kernel execution and data transfer.
1324 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L_EST);
1327 /* Communicate the forces */
1328 wallcycle_start(wcycle, ewcMOVEF);
1329 dd_move_f(cr->dd, f, fr->fshift);
1332 /* We should not update the shift forces here,
1333 * since f_twin is already included in f.
1335 dd_move_f(cr->dd, fr->f_twin, NULL);
1337 wallcycle_stop(wcycle, ewcMOVEF);
1342 /* wait for local forces (or calculate in emulation mode) */
1345 float cycles_tmp, cycles_wait_est;
1346 const float cuda_api_overhead_margin = 50000.0f; /* cycles */
1348 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1349 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1350 nbv->grp[eintLocal].nbat,
1352 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1354 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1356 if (bDoForces && DOMAINDECOMP(cr))
1358 cycles_wait_est = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L_EST);
1360 if (cycles_tmp < cuda_api_overhead_margin)
1362 /* We measured few cycles, it could be that the kernel
1363 * and transfer finished earlier and there was no actual
1364 * wait time, only API call overhead.
1365 * Then the actual time could be anywhere between 0 and
1366 * cycles_wait_est. As a compromise, we use half the time.
1368 cycles_wait_est *= 0.5f;
1373 /* No force communication so we actually timed the wait */
1374 cycles_wait_est = cycles_tmp;
1376 /* Even though this is after dd_move_f, the actual task we are
1377 * waiting for runs asynchronously with dd_move_f and we usually
1378 * have nothing to balance it with, so we can and should add
1379 * the time to the force time for load balancing.
1381 cycles_force += cycles_wait_est;
1382 cycles_wait_gpu += cycles_wait_est;
1384 /* now clear the GPU outputs while we finish the step on the CPU */
1386 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1387 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1388 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1392 wallcycle_start_nocount(wcycle, ewcFORCE);
1393 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1394 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1396 wallcycle_stop(wcycle, ewcFORCE);
1398 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1399 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1400 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1402 /* skip the reduction if there was no non-local work to do */
1403 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1404 nbv->grp[eintLocal].nbat, f);
1406 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1407 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1410 if (DOMAINDECOMP(cr))
1412 dd_force_flop_stop(cr->dd, nrnb);
1415 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1418 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1425 if (IR_ELEC_FIELD(*inputrec))
1427 /* Compute forces due to electric field */
1428 calc_f_el(MASTER(cr) ? field : NULL,
1429 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1430 inputrec->ex, inputrec->et, t);
1433 /* If we have NoVirSum forces, but we do not calculate the virial,
1434 * we sum fr->f_novirsum=f later.
1436 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1438 wallcycle_start(wcycle, ewcVSITESPREAD);
1439 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1440 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1441 wallcycle_stop(wcycle, ewcVSITESPREAD);
1445 wallcycle_start(wcycle, ewcVSITESPREAD);
1446 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1448 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1449 wallcycle_stop(wcycle, ewcVSITESPREAD);
1453 if (flags & GMX_FORCE_VIRIAL)
1455 /* Calculation of the virial must be done after vsites! */
1456 calc_virial(0, mdatoms->homenr, x, f,
1457 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1461 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1463 /* Since the COM pulling is always done mass-weighted, no forces are
1464 * applied to vsites and this call can be done after vsite spreading.
1466 pull_potential_wrapper(cr, inputrec, box, x,
1467 f, vir_force, mdatoms, enerd, lambda, t,
1471 /* Add the forces from enforced rotation potentials (if any) */
1474 wallcycle_start(wcycle, ewcROTadd);
1475 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1476 wallcycle_stop(wcycle, ewcROTadd);
1479 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1480 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1482 if (PAR(cr) && !(cr->duty & DUTY_PME))
1484 /* In case of node-splitting, the PP nodes receive the long-range
1485 * forces, virial and energy from the PME nodes here.
1487 pme_receive_force_ener(cr, wcycle, enerd, fr);
1492 post_process_forces(cr, step, nrnb, wcycle,
1493 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1497 /* Sum the potential energy terms from group contributions */
1498 sum_epot(&(enerd->grpp), enerd->term);
1501 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1502 t_inputrec *inputrec,
1503 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1504 gmx_localtop_t *top,
1505 gmx_groups_t *groups,
1506 matrix box, rvec x[], history_t *hist,
1510 gmx_enerdata_t *enerd, t_fcdata *fcd,
1511 real *lambda, t_graph *graph,
1512 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1513 double t, FILE *field, gmx_edsam_t ed,
1514 gmx_bool bBornRadii,
1520 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
1521 gmx_bool bDoLongRangeNS, bDoForces, bSepLRF;
1522 gmx_bool bDoAdressWF;
1524 float cycles_pme, cycles_force;
1527 homenr = mdatoms->homenr;
1529 clear_mat(vir_force);
1532 if (DOMAINDECOMP(cr))
1534 cg1 = cr->dd->ncg_tot;
1545 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1546 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1547 /* Should we update the long-range neighborlists at this step? */
1548 bDoLongRangeNS = fr->bTwinRange && bNS;
1549 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1550 bFillGrid = (bNS && bStateChanged);
1551 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1552 bDoForces = (flags & GMX_FORCE_FORCES);
1553 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1554 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1556 /* should probably move this to the forcerec since it doesn't change */
1557 bDoAdressWF = ((fr->adress_type != eAdressOff));
1561 update_forcerec(fr, box);
1563 if (NEED_MUTOT(*inputrec))
1565 /* Calculate total (local) dipole moment in a temporary common array.
1566 * This makes it possible to sum them over nodes faster.
1568 calc_mu(start, homenr,
1569 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1574 if (fr->ePBC != epbcNONE)
1576 /* Compute shift vectors every step,
1577 * because of pressure coupling or box deformation!
1579 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1581 calc_shifts(box, fr->shift_vec);
1586 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1587 &(top->cgs), x, fr->cg_cm);
1588 inc_nrnb(nrnb, eNR_CGCM, homenr);
1589 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1591 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1593 unshift_self(graph, box, x);
1598 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1599 inc_nrnb(nrnb, eNR_CGCM, homenr);
1602 if (bCalcCGCM && gmx_debug_at)
1604 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1608 if (!(cr->duty & DUTY_PME))
1613 /* Send particle coordinates to the pme nodes.
1614 * Since this is only implemented for domain decomposition
1615 * and domain decomposition does not use the graph,
1616 * we do not need to worry about shifting.
1621 wallcycle_start(wcycle, ewcPP_PMESENDX);
1623 bBS = (inputrec->nwall == 2);
1626 copy_mat(box, boxs);
1627 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1630 if (EEL_PME(fr->eeltype))
1632 pme_flags |= GMX_PME_DO_COULOMB;
1635 if (EVDW_PME(fr->vdwtype))
1637 pme_flags |= GMX_PME_DO_LJ;
1640 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1641 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1642 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1645 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1647 #endif /* GMX_MPI */
1649 /* Communicate coordinates and sum dipole if necessary */
1650 if (DOMAINDECOMP(cr))
1652 wallcycle_start(wcycle, ewcMOVEX);
1653 dd_move_x(cr->dd, box, x);
1654 wallcycle_stop(wcycle, ewcMOVEX);
1657 /* update adress weight beforehand */
1658 if (bStateChanged && bDoAdressWF)
1660 /* need pbc for adress weight calculation with pbc_dx */
1661 set_pbc(&pbc, inputrec->ePBC, box);
1662 if (fr->adress_site == eAdressSITEcog)
1664 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1665 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1667 else if (fr->adress_site == eAdressSITEcom)
1669 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1670 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1672 else if (fr->adress_site == eAdressSITEatomatom)
1674 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1675 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1679 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1680 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1684 if (NEED_MUTOT(*inputrec))
1691 gmx_sumd(2*DIM, mu, cr);
1693 for (i = 0; i < 2; i++)
1695 for (j = 0; j < DIM; j++)
1697 fr->mu_tot[i][j] = mu[i*DIM + j];
1701 if (fr->efep == efepNO)
1703 copy_rvec(fr->mu_tot[0], mu_tot);
1707 for (j = 0; j < DIM; j++)
1710 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1715 /* Reset energies */
1716 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1717 clear_rvecs(SHIFTS, fr->fshift);
1721 wallcycle_start(wcycle, ewcNS);
1723 if (graph && bStateChanged)
1725 /* Calculate intramolecular shift vectors to make molecules whole */
1726 mk_mshift(fplog, graph, fr->ePBC, box, x);
1729 /* Do the actual neighbour searching */
1731 groups, top, mdatoms,
1732 cr, nrnb, bFillGrid,
1735 wallcycle_stop(wcycle, ewcNS);
1738 if (inputrec->implicit_solvent && bNS)
1740 make_gb_nblist(cr, inputrec->gb_algorithm,
1741 x, box, fr, &top->idef, graph, fr->born);
1744 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1746 wallcycle_start(wcycle, ewcPPDURINGPME);
1747 dd_force_flop_start(cr->dd, nrnb);
1752 /* Enforced rotation has its own cycle counter that starts after the collective
1753 * coordinates have been communicated. It is added to ddCyclF to allow
1754 * for proper load-balancing */
1755 wallcycle_start(wcycle, ewcROT);
1756 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1757 wallcycle_stop(wcycle, ewcROT);
1760 /* Start the force cycle counter.
1761 * This counter is stopped in do_forcelow_level.
1762 * No parallel communication should occur while this counter is running,
1763 * since that will interfere with the dynamic load balancing.
1765 wallcycle_start(wcycle, ewcFORCE);
1769 /* Reset forces for which the virial is calculated separately:
1770 * PME/Ewald forces if necessary */
1771 if (fr->bF_NoVirSum)
1773 if (flags & GMX_FORCE_VIRIAL)
1775 fr->f_novirsum = fr->f_novirsum_alloc;
1778 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1782 clear_rvecs(homenr, fr->f_novirsum+start);
1787 /* We are not calculating the pressure so we do not need
1788 * a separate array for forces that do not contribute
1795 /* Clear the short- and long-range forces */
1796 clear_rvecs(fr->natoms_force_constr, f);
1797 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1799 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1802 clear_rvec(fr->vir_diag_posres);
1804 if (inputrec->ePull == epullCONSTRAINT)
1806 clear_pull_forces(inputrec->pull);
1809 /* update QMMMrec, if necessary */
1812 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1815 /* Compute the bonded and non-bonded energies and optionally forces */
1816 do_force_lowlevel(fr, inputrec, &(top->idef),
1817 cr, nrnb, wcycle, mdatoms,
1818 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1820 inputrec->fepvals, lambda,
1821 graph, &(top->excls), fr->mu_tot,
1827 if (do_per_step(step, inputrec->nstcalclr))
1829 /* Add the long range forces to the short range forces */
1830 for (i = 0; i < fr->natoms_force_constr; i++)
1832 rvec_add(fr->f_twin[i], f[i], f[i]);
1837 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1841 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1844 if (DOMAINDECOMP(cr))
1846 dd_force_flop_stop(cr->dd, nrnb);
1849 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1855 if (IR_ELEC_FIELD(*inputrec))
1857 /* Compute forces due to electric field */
1858 calc_f_el(MASTER(cr) ? field : NULL,
1859 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1860 inputrec->ex, inputrec->et, t);
1863 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1865 /* Compute thermodynamic force in hybrid AdResS region */
1866 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1867 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1870 /* Communicate the forces */
1871 if (DOMAINDECOMP(cr))
1873 wallcycle_start(wcycle, ewcMOVEF);
1874 dd_move_f(cr->dd, f, fr->fshift);
1875 /* Do we need to communicate the separate force array
1876 * for terms that do not contribute to the single sum virial?
1877 * Position restraints and electric fields do not introduce
1878 * inter-cg forces, only full electrostatics methods do.
1879 * When we do not calculate the virial, fr->f_novirsum = f,
1880 * so we have already communicated these forces.
1882 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1883 (flags & GMX_FORCE_VIRIAL))
1885 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1889 /* We should not update the shift forces here,
1890 * since f_twin is already included in f.
1892 dd_move_f(cr->dd, fr->f_twin, NULL);
1894 wallcycle_stop(wcycle, ewcMOVEF);
1897 /* If we have NoVirSum forces, but we do not calculate the virial,
1898 * we sum fr->f_novirsum=f later.
1900 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1902 wallcycle_start(wcycle, ewcVSITESPREAD);
1903 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1904 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1905 wallcycle_stop(wcycle, ewcVSITESPREAD);
1909 wallcycle_start(wcycle, ewcVSITESPREAD);
1910 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1912 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1913 wallcycle_stop(wcycle, ewcVSITESPREAD);
1917 if (flags & GMX_FORCE_VIRIAL)
1919 /* Calculation of the virial must be done after vsites! */
1920 calc_virial(0, mdatoms->homenr, x, f,
1921 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1925 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1927 pull_potential_wrapper(cr, inputrec, box, x,
1928 f, vir_force, mdatoms, enerd, lambda, t,
1932 /* Add the forces from enforced rotation potentials (if any) */
1935 wallcycle_start(wcycle, ewcROTadd);
1936 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1937 wallcycle_stop(wcycle, ewcROTadd);
1940 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1941 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1943 if (PAR(cr) && !(cr->duty & DUTY_PME))
1945 /* In case of node-splitting, the PP nodes receive the long-range
1946 * forces, virial and energy from the PME nodes here.
1948 pme_receive_force_ener(cr, wcycle, enerd, fr);
1953 post_process_forces(cr, step, nrnb, wcycle,
1954 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1958 /* Sum the potential energy terms from group contributions */
1959 sum_epot(&(enerd->grpp), enerd->term);
1962 void do_force(FILE *fplog, t_commrec *cr,
1963 t_inputrec *inputrec,
1964 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1965 gmx_localtop_t *top,
1966 gmx_groups_t *groups,
1967 matrix box, rvec x[], history_t *hist,
1971 gmx_enerdata_t *enerd, t_fcdata *fcd,
1972 real *lambda, t_graph *graph,
1974 gmx_vsite_t *vsite, rvec mu_tot,
1975 double t, FILE *field, gmx_edsam_t ed,
1976 gmx_bool bBornRadii,
1979 /* modify force flag if not doing nonbonded */
1980 if (!fr->bNonbonded)
1982 flags &= ~GMX_FORCE_NONBONDED;
1985 switch (inputrec->cutoff_scheme)
1988 do_force_cutsVERLET(fplog, cr, inputrec,
2004 do_force_cutsGROUP(fplog, cr, inputrec,
2019 gmx_incons("Invalid cut-off scheme passed!");
2024 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2025 t_inputrec *ir, t_mdatoms *md,
2026 t_state *state, t_commrec *cr, t_nrnb *nrnb,
2027 t_forcerec *fr, gmx_localtop_t *top)
2029 int i, m, start, end;
2031 real dt = ir->delta_t;
2035 snew(savex, state->natoms);
2042 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2043 start, md->homenr, end);
2045 /* Do a first constrain to reset particles... */
2046 step = ir->init_step;
2049 char buf[STEPSTRSIZE];
2050 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2051 gmx_step_str(step, buf));
2055 /* constrain the current position */
2056 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2057 ir, cr, step, 0, 1.0, md,
2058 state->x, state->x, NULL,
2059 fr->bMolPBC, state->box,
2060 state->lambda[efptBONDED], &dvdl_dum,
2061 NULL, NULL, nrnb, econqCoord);
2064 /* constrain the inital velocity, and save it */
2065 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2066 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2067 ir, cr, step, 0, 1.0, md,
2068 state->x, state->v, state->v,
2069 fr->bMolPBC, state->box,
2070 state->lambda[efptBONDED], &dvdl_dum,
2071 NULL, NULL, nrnb, econqVeloc);
2073 /* constrain the inital velocities at t-dt/2 */
2074 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2076 for (i = start; (i < end); i++)
2078 for (m = 0; (m < DIM); m++)
2080 /* Reverse the velocity */
2081 state->v[i][m] = -state->v[i][m];
2082 /* Store the position at t-dt in buf */
2083 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2086 /* Shake the positions at t=-dt with the positions at t=0
2087 * as reference coordinates.
2091 char buf[STEPSTRSIZE];
2092 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2093 gmx_step_str(step, buf));
2096 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2097 ir, cr, step, -1, 1.0, md,
2098 state->x, savex, NULL,
2099 fr->bMolPBC, state->box,
2100 state->lambda[efptBONDED], &dvdl_dum,
2101 state->v, NULL, nrnb, econqCoord);
2103 for (i = start; i < end; i++)
2105 for (m = 0; m < DIM; m++)
2107 /* Re-reverse the velocities */
2108 state->v[i][m] = -state->v[i][m];
2117 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2118 double *enerout, double *virout)
2120 double enersum, virsum;
2121 double invscale, invscale2, invscale3;
2122 double r, ea, eb, ec, pa, pb, pc, pd;
2127 invscale = 1.0/scale;
2128 invscale2 = invscale*invscale;
2129 invscale3 = invscale*invscale2;
2131 /* Following summation derived from cubic spline definition,
2132 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2133 * the cubic spline. We first calculate the negative of the
2134 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2135 * add the more standard, abrupt cutoff correction to that result,
2136 * yielding the long-range correction for a switched function. We
2137 * perform both the pressure and energy loops at the same time for
2138 * simplicity, as the computational cost is low. */
2142 /* Since the dispersion table has been scaled down a factor
2143 * 6.0 and the repulsion a factor 12.0 to compensate for the
2144 * c6/c12 parameters inside nbfp[] being scaled up (to save
2145 * flops in kernels), we need to correct for this.
2156 for (ri = rstart; ri < rend; ++ri)
2160 eb = 2.0*invscale2*r;
2164 pb = 3.0*invscale2*r;
2165 pc = 3.0*invscale*r*r;
2168 /* this "8" is from the packing in the vdwtab array - perhaps
2169 should be defined? */
2171 offset = 8*ri + offstart;
2172 y0 = vdwtab[offset];
2173 f = vdwtab[offset+1];
2174 g = vdwtab[offset+2];
2175 h = vdwtab[offset+3];
2177 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);
2178 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);
2180 *enerout = 4.0*M_PI*enersum*tabfactor;
2181 *virout = 4.0*M_PI*virsum*tabfactor;
2184 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2186 double eners[2], virs[2], enersum, virsum;
2187 double r0, rc3, rc9;
2189 real scale, *vdwtab;
2191 fr->enershiftsix = 0;
2192 fr->enershifttwelve = 0;
2193 fr->enerdiffsix = 0;
2194 fr->enerdifftwelve = 0;
2196 fr->virdifftwelve = 0;
2198 if (eDispCorr != edispcNO)
2200 for (i = 0; i < 2; i++)
2205 if ((fr->vdw_modifier == eintmodPOTSHIFT) ||
2206 (fr->vdw_modifier == eintmodPOTSWITCH) ||
2207 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2208 (fr->vdwtype == evdwSHIFT) ||
2209 (fr->vdwtype == evdwSWITCH))
2211 if (((fr->vdw_modifier == eintmodPOTSWITCH) ||
2212 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2213 (fr->vdwtype == evdwSWITCH)) && fr->rvdw_switch == 0)
2216 "With dispersion correction rvdw-switch can not be zero "
2217 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2220 scale = fr->nblists[0].table_vdw.scale;
2221 vdwtab = fr->nblists[0].table_vdw.data;
2223 /* Round the cut-offs to exact table values for precision */
2224 ri0 = static_cast<int>(floor(fr->rvdw_switch*scale));
2225 ri1 = static_cast<int>(ceil(fr->rvdw*scale));
2227 /* The code below has some support for handling force-switching, i.e.
2228 * when the force (instead of potential) is switched over a limited
2229 * region. This leads to a constant shift in the potential inside the
2230 * switching region, which we can handle by adding a constant energy
2231 * term in the force-switch case just like when we do potential-shift.
2233 * For now this is not enabled, but to keep the functionality in the
2234 * code we check separately for switch and shift. When we do force-switch
2235 * the shifting point is rvdw_switch, while it is the cutoff when we
2236 * have a classical potential-shift.
2238 * For a pure potential-shift the potential has a constant shift
2239 * all the way out to the cutoff, and that is it. For other forms
2240 * we need to calculate the constant shift up to the point where we
2241 * start modifying the potential.
2243 ri0 = (fr->vdw_modifier == eintmodPOTSHIFT) ? ri1 : ri0;
2249 if ((fr->vdw_modifier == eintmodFORCESWITCH) ||
2250 (fr->vdwtype == evdwSHIFT))
2252 /* Determine the constant energy shift below rvdw_switch.
2253 * Table has a scale factor since we have scaled it down to compensate
2254 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2256 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2257 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2259 else if (fr->vdw_modifier == eintmodPOTSHIFT)
2261 fr->enershiftsix = (real)(-1.0/(rc3*rc3));
2262 fr->enershifttwelve = (real)( 1.0/(rc9*rc3));
2265 /* Add the constant part from 0 to rvdw_switch.
2266 * This integration from 0 to rvdw_switch overcounts the number
2267 * of interactions by 1, as it also counts the self interaction.
2268 * We will correct for this later.
2270 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2271 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2273 /* Calculate the contribution in the range [r0,r1] where we
2274 * modify the potential. For a pure potential-shift modifier we will
2275 * have ri0==ri1, and there will not be any contribution here.
2277 for (i = 0; i < 2; i++)
2281 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2282 eners[i] -= enersum;
2286 /* Alright: Above we compensated by REMOVING the parts outside r0
2287 * corresponding to the ideal VdW 1/r6 and /r12 potentials.
2289 * Regardless of whether r0 is the point where we start switching,
2290 * or the cutoff where we calculated the constant shift, we include
2291 * all the parts we are missing out to infinity from r0 by
2292 * calculating the analytical dispersion correction.
2294 eners[0] += -4.0*M_PI/(3.0*rc3);
2295 eners[1] += 4.0*M_PI/(9.0*rc9);
2296 virs[0] += 8.0*M_PI/rc3;
2297 virs[1] += -16.0*M_PI/(3.0*rc9);
2299 else if (fr->vdwtype == evdwCUT ||
2300 EVDW_PME(fr->vdwtype) ||
2301 fr->vdwtype == evdwUSER)
2303 if (fr->vdwtype == evdwUSER && fplog)
2306 "WARNING: using dispersion correction with user tables\n");
2309 /* Note that with LJ-PME, the dispersion correction is multiplied
2310 * by the difference between the actual C6 and the value of C6
2311 * that would produce the combination rule.
2312 * This means the normal energy and virial difference formulas
2316 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2318 /* Contribution beyond the cut-off */
2319 eners[0] += -4.0*M_PI/(3.0*rc3);
2320 eners[1] += 4.0*M_PI/(9.0*rc9);
2321 if (fr->vdw_modifier == eintmodPOTSHIFT)
2323 /* Contribution within the cut-off */
2324 eners[0] += -4.0*M_PI/(3.0*rc3);
2325 eners[1] += 4.0*M_PI/(3.0*rc9);
2327 /* Contribution beyond the cut-off */
2328 virs[0] += 8.0*M_PI/rc3;
2329 virs[1] += -16.0*M_PI/(3.0*rc9);
2334 "Dispersion correction is not implemented for vdw-type = %s",
2335 evdw_names[fr->vdwtype]);
2338 /* When we deprecate the group kernels the code below can go too */
2339 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2341 /* Calculate self-interaction coefficient (assuming that
2342 * the reciprocal-space contribution is constant in the
2343 * region that contributes to the self-interaction).
2345 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2347 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2348 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2351 fr->enerdiffsix = eners[0];
2352 fr->enerdifftwelve = eners[1];
2353 /* The 0.5 is due to the Gromacs definition of the virial */
2354 fr->virdiffsix = 0.5*virs[0];
2355 fr->virdifftwelve = 0.5*virs[1];
2359 void calc_dispcorr(t_inputrec *ir, t_forcerec *fr,
2361 matrix box, real lambda, tensor pres, tensor virial,
2362 real *prescorr, real *enercorr, real *dvdlcorr)
2364 gmx_bool bCorrAll, bCorrPres;
2365 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2375 if (ir->eDispCorr != edispcNO)
2377 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2378 ir->eDispCorr == edispcAllEnerPres);
2379 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2380 ir->eDispCorr == edispcAllEnerPres);
2382 invvol = 1/det(box);
2385 /* Only correct for the interactions with the inserted molecule */
2386 dens = (natoms - fr->n_tpi)*invvol;
2391 dens = natoms*invvol;
2392 ninter = 0.5*natoms;
2395 if (ir->efep == efepNO)
2397 avcsix = fr->avcsix[0];
2398 avctwelve = fr->avctwelve[0];
2402 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2403 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2406 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2407 *enercorr += avcsix*enerdiff;
2409 if (ir->efep != efepNO)
2411 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2415 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2416 *enercorr += avctwelve*enerdiff;
2417 if (fr->efep != efepNO)
2419 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2425 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2426 if (ir->eDispCorr == edispcAllEnerPres)
2428 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2430 /* The factor 2 is because of the Gromacs virial definition */
2431 spres = -2.0*invvol*svir*PRESFAC;
2433 for (m = 0; m < DIM; m++)
2435 virial[m][m] += svir;
2436 pres[m][m] += spres;
2441 /* Can't currently control when it prints, for now, just print when degugging */
2446 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2452 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2453 *enercorr, spres, svir);
2457 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2461 if (fr->efep != efepNO)
2463 *dvdlcorr += dvdlambda;
2468 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2469 t_graph *graph, rvec x[])
2473 fprintf(fplog, "Removing pbc first time\n");
2475 calc_shifts(box, fr->shift_vec);
2478 mk_mshift(fplog, graph, fr->ePBC, box, x);
2481 p_graph(debug, "do_pbc_first 1", graph);
2483 shift_self(graph, box, x);
2484 /* By doing an extra mk_mshift the molecules that are broken
2485 * because they were e.g. imported from another software
2486 * will be made whole again. Such are the healing powers
2489 mk_mshift(fplog, graph, fr->ePBC, box, x);
2492 p_graph(debug, "do_pbc_first 2", graph);
2497 fprintf(fplog, "Done rmpbc\n");
2501 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2502 gmx_mtop_t *mtop, rvec x[],
2507 gmx_molblock_t *molb;
2509 if (bFirst && fplog)
2511 fprintf(fplog, "Removing pbc first time\n");
2516 for (mb = 0; mb < mtop->nmolblock; mb++)
2518 molb = &mtop->molblock[mb];
2519 if (molb->natoms_mol == 1 ||
2520 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2522 /* Just one atom or charge group in the molecule, no PBC required */
2523 as += molb->nmol*molb->natoms_mol;
2527 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2528 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2529 0, molb->natoms_mol, FALSE, FALSE, graph);
2531 for (mol = 0; mol < molb->nmol; mol++)
2533 mk_mshift(fplog, graph, ePBC, box, x+as);
2535 shift_self(graph, box, x+as);
2536 /* The molecule is whole now.
2537 * We don't need the second mk_mshift call as in do_pbc_first,
2538 * since we no longer need this graph.
2541 as += molb->natoms_mol;
2549 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2550 gmx_mtop_t *mtop, rvec x[])
2552 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2555 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2556 gmx_mtop_t *mtop, rvec x[])
2558 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2561 void finish_run(FILE *fplog, t_commrec *cr,
2562 t_inputrec *inputrec,
2563 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2564 gmx_walltime_accounting_t walltime_accounting,
2565 nonbonded_verlet_t *nbv,
2566 gmx_bool bWriteStat)
2568 t_nrnb *nrnb_tot = NULL;
2570 double nbfs = 0, mflop = 0;
2571 double elapsed_time,
2572 elapsed_time_over_all_ranks,
2573 elapsed_time_over_all_threads,
2574 elapsed_time_over_all_threads_over_all_ranks;
2575 wallcycle_sum(cr, wcycle);
2581 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2582 cr->mpi_comm_mysim);
2590 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2591 elapsed_time_over_all_ranks = elapsed_time;
2592 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2593 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2597 /* reduce elapsed_time over all MPI ranks in the current simulation */
2598 MPI_Allreduce(&elapsed_time,
2599 &elapsed_time_over_all_ranks,
2600 1, MPI_DOUBLE, MPI_SUM,
2601 cr->mpi_comm_mysim);
2602 elapsed_time_over_all_ranks /= cr->nnodes;
2603 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2604 * current simulation. */
2605 MPI_Allreduce(&elapsed_time_over_all_threads,
2606 &elapsed_time_over_all_threads_over_all_ranks,
2607 1, MPI_DOUBLE, MPI_SUM,
2608 cr->mpi_comm_mysim);
2614 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2621 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2623 print_dd_statistics(cr, inputrec, fplog);
2628 wallclock_gpu_t* gputimes = use_GPU(nbv) ?
2629 nbnxn_cuda_get_timings(nbv->cu_nbv) : NULL;
2630 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2631 elapsed_time_over_all_ranks,
2634 if (EI_DYNAMICS(inputrec->eI))
2636 delta_t = inputrec->delta_t;
2641 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2642 elapsed_time_over_all_ranks,
2643 walltime_accounting_get_nsteps_done(walltime_accounting),
2644 delta_t, nbfs, mflop);
2648 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2649 elapsed_time_over_all_ranks,
2650 walltime_accounting_get_nsteps_done(walltime_accounting),
2651 delta_t, nbfs, mflop);
2656 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2658 /* this function works, but could probably use a logic rewrite to keep all the different
2659 types of efep straight. */
2662 t_lambda *fep = ir->fepvals;
2664 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2666 for (i = 0; i < efptNR; i++)
2678 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2679 if checkpoint is set -- a kludge is in for now
2681 for (i = 0; i < efptNR; i++)
2683 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2684 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2686 lambda[i] = fep->init_lambda;
2689 lam0[i] = lambda[i];
2694 lambda[i] = fep->all_lambda[i][*fep_state];
2697 lam0[i] = lambda[i];
2703 /* need to rescale control temperatures to match current state */
2704 for (i = 0; i < ir->opts.ngtc; i++)
2706 if (ir->opts.ref_t[i] > 0)
2708 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2714 /* Send to the log the information on the current lambdas */
2717 fprintf(fplog, "Initial vector of lambda components:[ ");
2718 for (i = 0; i < efptNR; i++)
2720 fprintf(fplog, "%10.4f ", lambda[i]);
2722 fprintf(fplog, "]\n");
2728 void init_md(FILE *fplog,
2729 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2730 double *t, double *t0,
2731 real *lambda, int *fep_state, double *lam0,
2732 t_nrnb *nrnb, gmx_mtop_t *mtop,
2734 int nfile, const t_filenm fnm[],
2735 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2736 tensor force_vir, tensor shake_vir, rvec mu_tot,
2737 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags,
2738 gmx_wallcycle_t wcycle)
2742 /* Initial values */
2743 *t = *t0 = ir->init_t;
2746 for (i = 0; i < ir->opts.ngtc; i++)
2748 /* set bSimAnn if any group is being annealed */
2749 if (ir->opts.annealing[i] != eannNO)
2756 update_annealing_target_temp(&(ir->opts), ir->init_t);
2759 /* Initialize lambda variables */
2760 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2764 *upd = init_update(ir);
2770 *vcm = init_vcm(fplog, &mtop->groups, ir);
2773 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2775 if (ir->etc == etcBERENDSEN)
2777 please_cite(fplog, "Berendsen84a");
2779 if (ir->etc == etcVRESCALE)
2781 please_cite(fplog, "Bussi2007a");
2783 if (ir->eI == eiSD1)
2785 please_cite(fplog, "Goga2012");
2793 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv, wcycle);
2795 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2796 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2801 please_cite(fplog, "Fritsch12");
2802 please_cite(fplog, "Junghans10");
2804 /* Initiate variables */
2805 clear_mat(force_vir);
2806 clear_mat(shake_vir);