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.
43 #ifdef HAVE_SYS_TIME_H
48 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "chargegroup.h"
53 #include "gromacs/math/vec.h"
58 #include "gromacs/math/units.h"
71 #include "nbnxn_atomdata.h"
72 #include "nbnxn_search.h"
73 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
74 #include "nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
75 #include "nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
76 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
77 #include "nonbonded.h"
78 #include "../gmxlib/nonbonded/nb_kernel.h"
79 #include "../gmxlib/nonbonded/nb_free_energy.h"
81 #include "gromacs/legacyheaders/types/commrec.h"
82 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_data_mgmt.h"
83 #include "gromacs/pbcutil/ishift.h"
84 #include "gromacs/pbcutil/mshift.h"
85 #include "gromacs/timing/wallcycle.h"
86 #include "gromacs/timing/walltime_accounting.h"
87 #include "gromacs/utility/gmxmpi.h"
88 #include "gromacs/utility/smalloc.h"
89 #include "gromacs/essentialdynamics/edsam.h"
90 #include "gromacs/pulling/pull.h"
91 #include "gromacs/pulling/pull_rotation.h"
92 #include "gromacs/imd/imd.h"
96 #include "gmx_omp_nthreads.h"
98 #include "nbnxn_cuda/nbnxn_cuda.h"
100 #include "nb_verlet.h"
102 void print_time(FILE *out,
103 gmx_walltime_accounting_t walltime_accounting,
106 t_commrec gmx_unused *cr)
109 char timebuf[STRLEN];
110 double dt, elapsed_seconds, time_per_step;
113 #ifndef GMX_THREAD_MPI
119 fprintf(out, "step %s", gmx_step_str(step, buf));
120 if ((step >= ir->nstlist))
122 double seconds_since_epoch = gmx_gettime();
123 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
124 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
125 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
131 finish = (time_t) (seconds_since_epoch + dt);
132 gmx_ctime_r(&finish, timebuf, STRLEN);
133 sprintf(buf, "%s", timebuf);
134 buf[strlen(buf)-1] = '\0';
135 fprintf(out, ", will finish %s", buf);
139 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
144 fprintf(out, " performance: %.1f ns/day ",
145 ir->delta_t/1000*24*60*60/time_per_step);
148 #ifndef GMX_THREAD_MPI
158 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
161 char time_string[STRLEN];
170 char timebuf[STRLEN];
171 time_t temp_time = (time_t) the_time;
173 gmx_ctime_r(&temp_time, timebuf, STRLEN);
174 for (i = 0; timebuf[i] >= ' '; i++)
176 time_string[i] = timebuf[i];
178 time_string[i] = '\0';
181 fprintf(fplog, "%s on rank %d %s\n", title, nodeid, time_string);
184 void print_start(FILE *fplog, t_commrec *cr,
185 gmx_walltime_accounting_t walltime_accounting,
190 sprintf(buf, "Started %s", name);
191 print_date_and_time(fplog, cr->nodeid, buf,
192 walltime_accounting_get_start_time_stamp(walltime_accounting));
195 static void sum_forces(int start, int end, rvec f[], rvec flr[])
201 pr_rvecs(debug, 0, "fsr", f+start, end-start);
202 pr_rvecs(debug, 0, "flr", flr+start, end-start);
204 for (i = start; (i < end); i++)
206 rvec_inc(f[i], flr[i]);
211 * calc_f_el calculates forces due to an electric field.
213 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
215 * Et[] contains the parameters for the time dependent
216 * part of the field (not yet used).
217 * Ex[] contains the parameters for
218 * the spatial dependent part of the field. You can have cool periodic
219 * fields in principle, but only a constant field is supported
221 * The function should return the energy due to the electric field
222 * (if any) but for now returns 0.
225 * There can be problems with the virial.
226 * Since the field is not self-consistent this is unavoidable.
227 * For neutral molecules the virial is correct within this approximation.
228 * For neutral systems with many charged molecules the error is small.
229 * But for systems with a net charge or a few charged molecules
230 * the error can be significant when the field is high.
231 * Solution: implement a self-consitent electric field into PME.
233 static void calc_f_el(FILE *fp, int start, int homenr,
234 real charge[], rvec f[],
235 t_cosines Ex[], t_cosines Et[], double t)
241 for (m = 0; (m < DIM); m++)
248 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
252 Ext[m] = cos(Et[m].a[0]*t);
261 /* Convert the field strength from V/nm to MD-units */
262 Ext[m] *= Ex[m].a[0]*FIELDFAC;
263 for (i = start; (i < start+homenr); i++)
265 f[i][m] += charge[i]*Ext[m];
275 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
276 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
280 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
281 tensor vir_part, t_graph *graph, matrix box,
282 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
287 /* The short-range virial from surrounding boxes */
289 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
290 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
292 /* Calculate partial virial, for local atoms only, based on short range.
293 * Total virial is computed in global_stat, called from do_md
295 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
296 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
298 /* Add position restraint contribution */
299 for (i = 0; i < DIM; i++)
301 vir_part[i][i] += fr->vir_diag_posres[i];
304 /* Add wall contribution */
305 for (i = 0; i < DIM; i++)
307 vir_part[i][ZZ] += fr->vir_wall_z[i];
312 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
316 static void posres_wrapper(FILE *fplog,
322 matrix box, rvec x[],
323 gmx_enerdata_t *enerd,
331 /* Position restraints always require full pbc */
332 set_pbc(&pbc, ir->ePBC, box);
334 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
335 top->idef.iparams_posres,
336 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
337 ir->ePBC == epbcNONE ? NULL : &pbc,
338 lambda[efptRESTRAINT], &dvdl,
339 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
342 gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
344 enerd->term[F_POSRES] += v;
345 /* If just the force constant changes, the FEP term is linear,
346 * but if k changes, it is not.
348 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
349 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
351 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
353 for (i = 0; i < enerd->n_lambda; i++)
355 real dvdl_dum, lambda_dum;
357 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
358 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
359 top->idef.iparams_posres,
360 (const rvec*)x, NULL, NULL,
361 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
362 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
363 enerd->enerpart_lambda[i] += v;
368 static void fbposres_wrapper(t_inputrec *ir,
371 matrix box, rvec x[],
372 gmx_enerdata_t *enerd,
378 /* Flat-bottomed position restraints always require full pbc */
379 set_pbc(&pbc, ir->ePBC, box);
380 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
381 top->idef.iparams_fbposres,
382 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
383 ir->ePBC == epbcNONE ? NULL : &pbc,
384 fr->rc_scaling, fr->ePBC, fr->posres_com);
385 enerd->term[F_FBPOSRES] += v;
386 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
389 static void pull_potential_wrapper(FILE *fplog,
393 matrix box, rvec x[],
397 gmx_enerdata_t *enerd,
404 /* Calculate the center of mass forces, this requires communication,
405 * which is why pull_potential is called close to other communication.
406 * The virial contribution is calculated directly,
407 * which is why we call pull_potential after calc_virial.
409 set_pbc(&pbc, ir->ePBC, box);
411 enerd->term[F_COM_PULL] +=
412 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
413 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
416 gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
418 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
421 static void pme_receive_force_ener(FILE *fplog,
424 gmx_wallcycle_t wcycle,
425 gmx_enerdata_t *enerd,
428 real e_q, e_lj, v, dvdl_q, dvdl_lj;
429 float cycles_ppdpme, cycles_seppme;
431 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
432 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
434 /* In case of node-splitting, the PP nodes receive the long-range
435 * forces, virial and energy from the PME nodes here.
437 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
440 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
441 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
445 gmx_print_sepdvdl(fplog, "Electrostatic PME mesh", e_q, dvdl_q);
446 gmx_print_sepdvdl(fplog, "Lennard-Jones PME mesh", e_lj, dvdl_lj);
448 enerd->term[F_COUL_RECIP] += e_q;
449 enerd->term[F_LJ_RECIP] += e_lj;
450 enerd->dvdl_lin[efptCOUL] += dvdl_q;
451 enerd->dvdl_lin[efptVDW] += dvdl_lj;
455 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
457 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
460 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
461 gmx_int64_t step, real pforce, rvec *x, rvec *f)
465 char buf[STEPSTRSIZE];
468 for (i = 0; i < md->homenr; i++)
471 /* We also catch NAN, if the compiler does not optimize this away. */
472 if (fn2 >= pf2 || fn2 != fn2)
474 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
475 gmx_step_str(step, buf),
476 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
481 static void post_process_forces(t_commrec *cr,
483 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
485 matrix box, rvec x[],
490 t_forcerec *fr, gmx_vsite_t *vsite,
497 /* Spread the mesh force on virtual sites to the other particles...
498 * This is parallellized. MPI communication is performed
499 * if the constructing atoms aren't local.
501 wallcycle_start(wcycle, ewcVSITESPREAD);
502 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
503 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
505 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
506 wallcycle_stop(wcycle, ewcVSITESPREAD);
508 if (flags & GMX_FORCE_VIRIAL)
510 /* Now add the forces, this is local */
513 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
517 sum_forces(0, mdatoms->homenr,
520 if (EEL_FULL(fr->eeltype))
522 /* Add the mesh contribution to the virial */
523 m_add(vir_force, fr->vir_el_recip, vir_force);
525 if (EVDW_PME(fr->vdwtype))
527 /* Add the mesh contribution to the virial */
528 m_add(vir_force, fr->vir_lj_recip, vir_force);
532 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
537 if (fr->print_force >= 0)
539 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
543 static void do_nb_verlet(t_forcerec *fr,
544 interaction_const_t *ic,
545 gmx_enerdata_t *enerd,
546 int flags, int ilocality,
549 gmx_wallcycle_t wcycle)
551 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
553 nonbonded_verlet_group_t *nbvg;
556 if (!(flags & GMX_FORCE_NONBONDED))
558 /* skip non-bonded calculation */
562 nbvg = &fr->nbv->grp[ilocality];
564 /* CUDA kernel launch overhead is already timed separately */
565 if (fr->cutoff_scheme != ecutsVERLET)
567 gmx_incons("Invalid cut-off scheme passed!");
570 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
574 wallcycle_sub_start(wcycle, ewcsNONBONDED);
576 switch (nbvg->kernel_type)
578 case nbnxnk4x4_PlainC:
579 nbnxn_kernel_ref(&nbvg->nbl_lists,
585 enerd->grpp.ener[egCOULSR],
587 enerd->grpp.ener[egBHAMSR] :
588 enerd->grpp.ener[egLJSR]);
591 case nbnxnk4xN_SIMD_4xN:
592 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
599 enerd->grpp.ener[egCOULSR],
601 enerd->grpp.ener[egBHAMSR] :
602 enerd->grpp.ener[egLJSR]);
604 case nbnxnk4xN_SIMD_2xNN:
605 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
612 enerd->grpp.ener[egCOULSR],
614 enerd->grpp.ener[egBHAMSR] :
615 enerd->grpp.ener[egLJSR]);
618 case nbnxnk8x8x8_CUDA:
619 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
622 case nbnxnk8x8x8_PlainC:
623 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
628 nbvg->nbat->out[0].f,
630 enerd->grpp.ener[egCOULSR],
632 enerd->grpp.ener[egBHAMSR] :
633 enerd->grpp.ener[egLJSR]);
637 gmx_incons("Invalid nonbonded kernel type passed!");
642 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
645 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
647 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
649 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
650 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
652 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
656 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
658 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
659 if (flags & GMX_FORCE_ENERGY)
661 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
662 enr_nbnxn_kernel_ljc += 1;
663 enr_nbnxn_kernel_lj += 1;
666 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
667 nbvg->nbl_lists.natpair_ljq);
668 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
669 nbvg->nbl_lists.natpair_lj);
670 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
671 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
672 nbvg->nbl_lists.natpair_q);
674 if (ic->vdw_modifier == eintmodFORCESWITCH)
676 /* We add up the switch cost separately */
677 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
678 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
680 if (ic->vdw_modifier == eintmodPOTSWITCH)
682 /* We add up the switch cost separately */
683 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
684 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
686 if (ic->vdwtype == evdwPME)
688 /* We add up the LJ Ewald cost separately */
689 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
690 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
694 static void do_nb_verlet_fep(nbnxn_pairlist_set_t *nbl_lists,
701 gmx_enerdata_t *enerd,
704 gmx_wallcycle_t wcycle)
707 nb_kernel_data_t kernel_data;
709 real dvdl_nb[efptNR];
714 /* Add short-range interactions */
715 donb_flags |= GMX_NONBONDED_DO_SR;
717 /* Currently all group scheme kernels always calculate (shift-)forces */
718 if (flags & GMX_FORCE_FORCES)
720 donb_flags |= GMX_NONBONDED_DO_FORCE;
722 if (flags & GMX_FORCE_VIRIAL)
724 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
726 if (flags & GMX_FORCE_ENERGY)
728 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
730 if (flags & GMX_FORCE_DO_LR)
732 donb_flags |= GMX_NONBONDED_DO_LR;
735 kernel_data.flags = donb_flags;
736 kernel_data.lambda = lambda;
737 kernel_data.dvdl = dvdl_nb;
739 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
740 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
742 /* reset free energy components */
743 for (i = 0; i < efptNR; i++)
748 assert(gmx_omp_nthreads_get(emntNonbonded) == nbl_lists->nnbl);
750 wallcycle_sub_start(wcycle, ewcsNONBONDED);
751 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
752 for (th = 0; th < nbl_lists->nnbl; th++)
754 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
755 x, f, fr, mdatoms, &kernel_data, nrnb);
758 if (fepvals->sc_alpha != 0)
760 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
761 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
765 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
766 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
769 /* If we do foreign lambda and we have soft-core interactions
770 * we have to recalculate the (non-linear) energies contributions.
772 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
774 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
775 kernel_data.lambda = lam_i;
776 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
777 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
778 /* Note that we add to kernel_data.dvdl, but ignore the result */
780 for (i = 0; i < enerd->n_lambda; i++)
782 for (j = 0; j < efptNR; j++)
784 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
786 reset_foreign_enerdata(enerd);
787 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
788 for (th = 0; th < nbl_lists->nnbl; th++)
790 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
791 x, f, fr, mdatoms, &kernel_data, nrnb);
794 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
795 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
799 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
802 gmx_bool use_GPU(const nonbonded_verlet_t *nbv)
804 return nbv != NULL && nbv->bUseGPU;
807 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
808 t_inputrec *inputrec,
809 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
811 gmx_groups_t gmx_unused *groups,
812 matrix box, rvec x[], history_t *hist,
816 gmx_enerdata_t *enerd, t_fcdata *fcd,
817 real *lambda, t_graph *graph,
818 t_forcerec *fr, interaction_const_t *ic,
819 gmx_vsite_t *vsite, rvec mu_tot,
820 double t, FILE *field, gmx_edsam_t ed,
828 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
829 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
830 gmx_bool bDiffKernels = FALSE;
832 rvec vzero, box_diag;
834 float cycles_pme, cycles_force, cycles_wait_gpu;
835 nonbonded_verlet_t *nbv;
840 nb_kernel_type = fr->nbv->grp[0].kernel_type;
843 homenr = mdatoms->homenr;
845 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
847 clear_mat(vir_force);
850 if (DOMAINDECOMP(cr))
852 cg1 = cr->dd->ncg_tot;
863 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
864 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
865 bFillGrid = (bNS && bStateChanged);
866 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
867 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
868 bDoForces = (flags & GMX_FORCE_FORCES);
869 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
870 bUseGPU = fr->nbv->bUseGPU;
871 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
875 update_forcerec(fr, box);
877 if (NEED_MUTOT(*inputrec))
879 /* Calculate total (local) dipole moment in a temporary common array.
880 * This makes it possible to sum them over nodes faster.
882 calc_mu(start, homenr,
883 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
888 if (fr->ePBC != epbcNONE)
890 /* Compute shift vectors every step,
891 * because of pressure coupling or box deformation!
893 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
895 calc_shifts(box, fr->shift_vec);
900 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
901 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
903 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
905 unshift_self(graph, box, x);
909 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
910 fr->shift_vec, nbv->grp[0].nbat);
913 if (!(cr->duty & DUTY_PME))
915 /* Send particle coordinates to the pme nodes.
916 * Since this is only implemented for domain decomposition
917 * and domain decomposition does not use the graph,
918 * we do not need to worry about shifting.
923 wallcycle_start(wcycle, ewcPP_PMESENDX);
925 bBS = (inputrec->nwall == 2);
929 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
932 if (EEL_PME(fr->eeltype))
934 pme_flags |= GMX_PME_DO_COULOMB;
937 if (EVDW_PME(fr->vdwtype))
939 pme_flags |= GMX_PME_DO_LJ;
942 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
943 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
944 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
947 wallcycle_stop(wcycle, ewcPP_PMESENDX);
951 /* do gridding for pair search */
954 if (graph && bStateChanged)
956 /* Calculate intramolecular shift vectors to make molecules whole */
957 mk_mshift(fplog, graph, fr->ePBC, box, x);
961 box_diag[XX] = box[XX][XX];
962 box_diag[YY] = box[YY][YY];
963 box_diag[ZZ] = box[ZZ][ZZ];
965 wallcycle_start(wcycle, ewcNS);
968 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
969 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
971 0, mdatoms->homenr, -1, fr->cginfo, x,
973 nbv->grp[eintLocal].kernel_type,
974 nbv->grp[eintLocal].nbat);
975 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
979 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
980 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
982 nbv->grp[eintNonlocal].kernel_type,
983 nbv->grp[eintNonlocal].nbat);
984 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
987 if (nbv->ngrp == 1 ||
988 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
990 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
991 nbv->nbs, mdatoms, fr->cginfo);
995 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
996 nbv->nbs, mdatoms, fr->cginfo);
997 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
998 nbv->nbs, mdatoms, fr->cginfo);
1000 wallcycle_stop(wcycle, ewcNS);
1003 /* initialize the GPU atom data and copy shift vector */
1008 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1009 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
1010 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1013 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1014 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
1015 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1018 /* do local pair search */
1021 wallcycle_start_nocount(wcycle, ewcNS);
1022 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1023 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
1026 nbv->min_ci_balanced,
1027 &nbv->grp[eintLocal].nbl_lists,
1029 nbv->grp[eintLocal].kernel_type,
1031 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1035 /* initialize local pair-list on the GPU */
1036 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1037 nbv->grp[eintLocal].nbl_lists.nbl[0],
1040 wallcycle_stop(wcycle, ewcNS);
1044 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1045 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1046 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
1047 nbv->grp[eintLocal].nbat);
1048 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1049 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1054 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1055 /* launch local nonbonded F on GPU */
1056 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
1058 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1061 /* Communicate coordinates and sum dipole if necessary +
1062 do non-local pair search */
1063 if (DOMAINDECOMP(cr))
1065 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
1066 nbv->grp[eintLocal].kernel_type);
1070 /* With GPU+CPU non-bonded calculations we need to copy
1071 * the local coordinates to the non-local nbat struct
1072 * (in CPU format) as the non-local kernel call also
1073 * calculates the local - non-local interactions.
1075 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1076 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1077 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
1078 nbv->grp[eintNonlocal].nbat);
1079 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1080 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1085 wallcycle_start_nocount(wcycle, ewcNS);
1086 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1090 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
1093 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
1096 nbv->min_ci_balanced,
1097 &nbv->grp[eintNonlocal].nbl_lists,
1099 nbv->grp[eintNonlocal].kernel_type,
1102 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1104 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1106 /* initialize non-local pair-list on the GPU */
1107 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1108 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1111 wallcycle_stop(wcycle, ewcNS);
1115 wallcycle_start(wcycle, ewcMOVEX);
1116 dd_move_x(cr->dd, box, x);
1118 /* When we don't need the total dipole we sum it in global_stat */
1119 if (bStateChanged && NEED_MUTOT(*inputrec))
1121 gmx_sumd(2*DIM, mu, cr);
1123 wallcycle_stop(wcycle, ewcMOVEX);
1125 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1126 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1127 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1128 nbv->grp[eintNonlocal].nbat);
1129 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1130 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1133 if (bUseGPU && !bDiffKernels)
1135 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1136 /* launch non-local nonbonded F on GPU */
1137 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1139 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1145 /* launch D2H copy-back F */
1146 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1147 if (DOMAINDECOMP(cr) && !bDiffKernels)
1149 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1150 flags, eatNonlocal);
1152 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1154 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1157 if (bStateChanged && NEED_MUTOT(*inputrec))
1161 gmx_sumd(2*DIM, mu, cr);
1164 for (i = 0; i < 2; i++)
1166 for (j = 0; j < DIM; j++)
1168 fr->mu_tot[i][j] = mu[i*DIM + j];
1172 if (fr->efep == efepNO)
1174 copy_rvec(fr->mu_tot[0], mu_tot);
1178 for (j = 0; j < DIM; j++)
1181 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1182 lambda[efptCOUL]*fr->mu_tot[1][j];
1186 /* Reset energies */
1187 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1188 clear_rvecs(SHIFTS, fr->fshift);
1190 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1192 wallcycle_start(wcycle, ewcPPDURINGPME);
1193 dd_force_flop_start(cr->dd, nrnb);
1198 /* Enforced rotation has its own cycle counter that starts after the collective
1199 * coordinates have been communicated. It is added to ddCyclF to allow
1200 * for proper load-balancing */
1201 wallcycle_start(wcycle, ewcROT);
1202 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1203 wallcycle_stop(wcycle, ewcROT);
1206 /* Start the force cycle counter.
1207 * This counter is stopped in do_forcelow_level.
1208 * No parallel communication should occur while this counter is running,
1209 * since that will interfere with the dynamic load balancing.
1211 wallcycle_start(wcycle, ewcFORCE);
1214 /* Reset forces for which the virial is calculated separately:
1215 * PME/Ewald forces if necessary */
1216 if (fr->bF_NoVirSum)
1218 if (flags & GMX_FORCE_VIRIAL)
1220 fr->f_novirsum = fr->f_novirsum_alloc;
1223 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1227 clear_rvecs(homenr, fr->f_novirsum+start);
1232 /* We are not calculating the pressure so we do not need
1233 * a separate array for forces that do not contribute
1240 /* Clear the short- and long-range forces */
1241 clear_rvecs(fr->natoms_force_constr, f);
1242 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1244 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1247 clear_rvec(fr->vir_diag_posres);
1250 if (inputrec->ePull == epullCONSTRAINT)
1252 clear_pull_forces(inputrec->pull);
1255 /* We calculate the non-bonded forces, when done on the CPU, here.
1256 * We do this before calling do_force_lowlevel, as in there bondeds
1257 * forces are calculated before PME, which does communication.
1258 * With this order, non-bonded and bonded force calculation imbalance
1259 * can be balanced out by the domain decomposition load balancing.
1264 /* Maybe we should move this into do_force_lowlevel */
1265 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1269 if (fr->efep != efepNO)
1271 /* Calculate the local and non-local free energy interactions here.
1272 * Happens here on the CPU both with and without GPU.
1274 if (fr->nbv->grp[eintLocal].nbl_lists.nbl_fep[0]->nrj > 0)
1276 do_nb_verlet_fep(&fr->nbv->grp[eintLocal].nbl_lists,
1278 inputrec->fepvals, lambda,
1279 enerd, flags, nrnb, wcycle);
1282 if (DOMAINDECOMP(cr) &&
1283 fr->nbv->grp[eintNonlocal].nbl_lists.nbl_fep[0]->nrj > 0)
1285 do_nb_verlet_fep(&fr->nbv->grp[eintNonlocal].nbl_lists,
1287 inputrec->fepvals, lambda,
1288 enerd, flags, nrnb, wcycle);
1292 if (!bUseOrEmulGPU || bDiffKernels)
1296 if (DOMAINDECOMP(cr))
1298 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1299 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1309 aloc = eintNonlocal;
1312 /* Add all the non-bonded force to the normal force array.
1313 * This can be split into a local a non-local part when overlapping
1314 * communication with calculation with domain decomposition.
1316 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1317 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1318 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1319 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1320 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1321 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1322 wallcycle_start_nocount(wcycle, ewcFORCE);
1324 /* if there are multiple fshift output buffers reduce them */
1325 if ((flags & GMX_FORCE_VIRIAL) &&
1326 nbv->grp[aloc].nbl_lists.nnbl > 1)
1328 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1333 /* update QMMMrec, if necessary */
1336 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1339 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1341 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1345 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1347 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1350 /* Compute the bonded and non-bonded energies and optionally forces */
1351 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1352 cr, nrnb, wcycle, mdatoms,
1353 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1354 &(top->atomtypes), bBornRadii, box,
1355 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1356 flags, &cycles_pme);
1360 if (do_per_step(step, inputrec->nstcalclr))
1362 /* Add the long range forces to the short range forces */
1363 for (i = 0; i < fr->natoms_force_constr; i++)
1365 rvec_add(fr->f_twin[i], f[i], f[i]);
1370 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1374 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1377 if (bUseOrEmulGPU && !bDiffKernels)
1379 /* wait for non-local forces (or calculate in emulation mode) */
1380 if (DOMAINDECOMP(cr))
1386 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1387 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1388 nbv->grp[eintNonlocal].nbat,
1390 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1392 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1393 cycles_wait_gpu += cycles_tmp;
1394 cycles_force += cycles_tmp;
1398 wallcycle_start_nocount(wcycle, ewcFORCE);
1399 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1401 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1403 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1404 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1405 /* skip the reduction if there was no non-local work to do */
1406 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1408 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1409 nbv->grp[eintNonlocal].nbat, f);
1411 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1412 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1416 if (bDoForces && DOMAINDECOMP(cr))
1418 /* Communicate the forces */
1419 wallcycle_start(wcycle, ewcMOVEF);
1420 dd_move_f(cr->dd, f, fr->fshift);
1421 /* Do we need to communicate the separate force array
1422 * for terms that do not contribute to the single sum virial?
1423 * Position restraints and electric fields do not introduce
1424 * inter-cg forces, only full electrostatics methods do.
1425 * When we do not calculate the virial, fr->f_novirsum = f,
1426 * so we have already communicated these forces.
1428 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1429 (flags & GMX_FORCE_VIRIAL))
1431 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1435 /* We should not update the shift forces here,
1436 * since f_twin is already included in f.
1438 dd_move_f(cr->dd, fr->f_twin, NULL);
1440 wallcycle_stop(wcycle, ewcMOVEF);
1445 /* wait for local forces (or calculate in emulation mode) */
1448 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1449 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1450 nbv->grp[eintLocal].nbat,
1452 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1454 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1456 /* now clear the GPU outputs while we finish the step on the CPU */
1458 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1459 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1460 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1464 wallcycle_start_nocount(wcycle, ewcFORCE);
1465 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1466 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1468 wallcycle_stop(wcycle, ewcFORCE);
1470 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1471 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1472 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1474 /* skip the reduction if there was no non-local work to do */
1475 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1476 nbv->grp[eintLocal].nbat, f);
1478 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1479 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1482 if (DOMAINDECOMP(cr))
1484 dd_force_flop_stop(cr->dd, nrnb);
1487 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1490 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1497 if (IR_ELEC_FIELD(*inputrec))
1499 /* Compute forces due to electric field */
1500 calc_f_el(MASTER(cr) ? field : NULL,
1501 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1502 inputrec->ex, inputrec->et, t);
1505 /* If we have NoVirSum forces, but we do not calculate the virial,
1506 * we sum fr->f_novirum=f later.
1508 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1510 wallcycle_start(wcycle, ewcVSITESPREAD);
1511 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1512 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1513 wallcycle_stop(wcycle, ewcVSITESPREAD);
1517 wallcycle_start(wcycle, ewcVSITESPREAD);
1518 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1520 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1521 wallcycle_stop(wcycle, ewcVSITESPREAD);
1525 if (flags & GMX_FORCE_VIRIAL)
1527 /* Calculation of the virial must be done after vsites! */
1528 calc_virial(0, mdatoms->homenr, x, f,
1529 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1533 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1535 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1536 f, vir_force, mdatoms, enerd, lambda, t);
1539 /* Add the forces from enforced rotation potentials (if any) */
1542 wallcycle_start(wcycle, ewcROTadd);
1543 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1544 wallcycle_stop(wcycle, ewcROTadd);
1547 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1548 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1550 if (PAR(cr) && !(cr->duty & DUTY_PME))
1552 /* In case of node-splitting, the PP nodes receive the long-range
1553 * forces, virial and energy from the PME nodes here.
1555 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1560 post_process_forces(cr, step, nrnb, wcycle,
1561 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1565 /* Sum the potential energy terms from group contributions */
1566 sum_epot(&(enerd->grpp), enerd->term);
1569 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1570 t_inputrec *inputrec,
1571 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1572 gmx_localtop_t *top,
1573 gmx_groups_t *groups,
1574 matrix box, rvec x[], history_t *hist,
1578 gmx_enerdata_t *enerd, t_fcdata *fcd,
1579 real *lambda, t_graph *graph,
1580 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1581 double t, FILE *field, gmx_edsam_t ed,
1582 gmx_bool bBornRadii,
1588 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1589 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1590 gmx_bool bDoAdressWF;
1592 rvec vzero, box_diag;
1593 real e, v, dvdlambda[efptNR];
1595 float cycles_pme, cycles_force;
1598 homenr = mdatoms->homenr;
1600 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1602 clear_mat(vir_force);
1605 if (DOMAINDECOMP(cr))
1607 cg1 = cr->dd->ncg_tot;
1618 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1619 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1620 /* Should we update the long-range neighborlists at this step? */
1621 bDoLongRangeNS = fr->bTwinRange && bNS;
1622 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1623 bFillGrid = (bNS && bStateChanged);
1624 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1625 bDoForces = (flags & GMX_FORCE_FORCES);
1626 bDoPotential = (flags & GMX_FORCE_ENERGY);
1627 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1628 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1630 /* should probably move this to the forcerec since it doesn't change */
1631 bDoAdressWF = ((fr->adress_type != eAdressOff));
1635 update_forcerec(fr, box);
1637 if (NEED_MUTOT(*inputrec))
1639 /* Calculate total (local) dipole moment in a temporary common array.
1640 * This makes it possible to sum them over nodes faster.
1642 calc_mu(start, homenr,
1643 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1648 if (fr->ePBC != epbcNONE)
1650 /* Compute shift vectors every step,
1651 * because of pressure coupling or box deformation!
1653 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1655 calc_shifts(box, fr->shift_vec);
1660 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1661 &(top->cgs), x, fr->cg_cm);
1662 inc_nrnb(nrnb, eNR_CGCM, homenr);
1663 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1665 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1667 unshift_self(graph, box, x);
1672 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1673 inc_nrnb(nrnb, eNR_CGCM, homenr);
1676 if (bCalcCGCM && gmx_debug_at)
1678 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1682 if (!(cr->duty & DUTY_PME))
1684 /* Send particle coordinates to the pme nodes.
1685 * Since this is only implemented for domain decomposition
1686 * and domain decomposition does not use the graph,
1687 * we do not need to worry about shifting.
1692 wallcycle_start(wcycle, ewcPP_PMESENDX);
1694 bBS = (inputrec->nwall == 2);
1697 copy_mat(box, boxs);
1698 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1701 if (EEL_PME(fr->eeltype))
1703 pme_flags |= GMX_PME_DO_COULOMB;
1706 if (EVDW_PME(fr->vdwtype))
1708 pme_flags |= GMX_PME_DO_LJ;
1711 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1712 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1713 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1716 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1718 #endif /* GMX_MPI */
1720 /* Communicate coordinates and sum dipole if necessary */
1721 if (DOMAINDECOMP(cr))
1723 wallcycle_start(wcycle, ewcMOVEX);
1724 dd_move_x(cr->dd, box, x);
1725 wallcycle_stop(wcycle, ewcMOVEX);
1728 /* update adress weight beforehand */
1729 if (bStateChanged && bDoAdressWF)
1731 /* need pbc for adress weight calculation with pbc_dx */
1732 set_pbc(&pbc, inputrec->ePBC, box);
1733 if (fr->adress_site == eAdressSITEcog)
1735 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1736 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1738 else if (fr->adress_site == eAdressSITEcom)
1740 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1741 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1743 else if (fr->adress_site == eAdressSITEatomatom)
1745 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1746 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1750 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1751 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1755 if (NEED_MUTOT(*inputrec))
1762 gmx_sumd(2*DIM, mu, cr);
1764 for (i = 0; i < 2; i++)
1766 for (j = 0; j < DIM; j++)
1768 fr->mu_tot[i][j] = mu[i*DIM + j];
1772 if (fr->efep == efepNO)
1774 copy_rvec(fr->mu_tot[0], mu_tot);
1778 for (j = 0; j < DIM; j++)
1781 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1786 /* Reset energies */
1787 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1788 clear_rvecs(SHIFTS, fr->fshift);
1792 wallcycle_start(wcycle, ewcNS);
1794 if (graph && bStateChanged)
1796 /* Calculate intramolecular shift vectors to make molecules whole */
1797 mk_mshift(fplog, graph, fr->ePBC, box, x);
1800 /* Do the actual neighbour searching */
1802 groups, top, mdatoms,
1803 cr, nrnb, bFillGrid,
1806 wallcycle_stop(wcycle, ewcNS);
1809 if (inputrec->implicit_solvent && bNS)
1811 make_gb_nblist(cr, inputrec->gb_algorithm,
1812 x, box, fr, &top->idef, graph, fr->born);
1815 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1817 wallcycle_start(wcycle, ewcPPDURINGPME);
1818 dd_force_flop_start(cr->dd, nrnb);
1823 /* Enforced rotation has its own cycle counter that starts after the collective
1824 * coordinates have been communicated. It is added to ddCyclF to allow
1825 * for proper load-balancing */
1826 wallcycle_start(wcycle, ewcROT);
1827 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1828 wallcycle_stop(wcycle, ewcROT);
1831 /* Start the force cycle counter.
1832 * This counter is stopped in do_forcelow_level.
1833 * No parallel communication should occur while this counter is running,
1834 * since that will interfere with the dynamic load balancing.
1836 wallcycle_start(wcycle, ewcFORCE);
1840 /* Reset forces for which the virial is calculated separately:
1841 * PME/Ewald forces if necessary */
1842 if (fr->bF_NoVirSum)
1844 if (flags & GMX_FORCE_VIRIAL)
1846 fr->f_novirsum = fr->f_novirsum_alloc;
1849 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1853 clear_rvecs(homenr, fr->f_novirsum+start);
1858 /* We are not calculating the pressure so we do not need
1859 * a separate array for forces that do not contribute
1866 /* Clear the short- and long-range forces */
1867 clear_rvecs(fr->natoms_force_constr, f);
1868 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1870 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1873 clear_rvec(fr->vir_diag_posres);
1875 if (inputrec->ePull == epullCONSTRAINT)
1877 clear_pull_forces(inputrec->pull);
1880 /* update QMMMrec, if necessary */
1883 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1886 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1888 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1892 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1894 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1897 /* Compute the bonded and non-bonded energies and optionally forces */
1898 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1899 cr, nrnb, wcycle, mdatoms,
1900 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1901 &(top->atomtypes), bBornRadii, box,
1902 inputrec->fepvals, lambda,
1903 graph, &(top->excls), fr->mu_tot,
1909 if (do_per_step(step, inputrec->nstcalclr))
1911 /* Add the long range forces to the short range forces */
1912 for (i = 0; i < fr->natoms_force_constr; i++)
1914 rvec_add(fr->f_twin[i], f[i], f[i]);
1919 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1923 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1926 if (DOMAINDECOMP(cr))
1928 dd_force_flop_stop(cr->dd, nrnb);
1931 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1937 if (IR_ELEC_FIELD(*inputrec))
1939 /* Compute forces due to electric field */
1940 calc_f_el(MASTER(cr) ? field : NULL,
1941 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1942 inputrec->ex, inputrec->et, t);
1945 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1947 /* Compute thermodynamic force in hybrid AdResS region */
1948 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1949 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1952 /* Communicate the forces */
1953 if (DOMAINDECOMP(cr))
1955 wallcycle_start(wcycle, ewcMOVEF);
1956 dd_move_f(cr->dd, f, fr->fshift);
1957 /* Do we need to communicate the separate force array
1958 * for terms that do not contribute to the single sum virial?
1959 * Position restraints and electric fields do not introduce
1960 * inter-cg forces, only full electrostatics methods do.
1961 * When we do not calculate the virial, fr->f_novirsum = f,
1962 * so we have already communicated these forces.
1964 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1965 (flags & GMX_FORCE_VIRIAL))
1967 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1971 /* We should not update the shift forces here,
1972 * since f_twin is already included in f.
1974 dd_move_f(cr->dd, fr->f_twin, NULL);
1976 wallcycle_stop(wcycle, ewcMOVEF);
1979 /* If we have NoVirSum forces, but we do not calculate the virial,
1980 * we sum fr->f_novirum=f later.
1982 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1984 wallcycle_start(wcycle, ewcVSITESPREAD);
1985 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1986 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1987 wallcycle_stop(wcycle, ewcVSITESPREAD);
1991 wallcycle_start(wcycle, ewcVSITESPREAD);
1992 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1994 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1995 wallcycle_stop(wcycle, ewcVSITESPREAD);
1999 if (flags & GMX_FORCE_VIRIAL)
2001 /* Calculation of the virial must be done after vsites! */
2002 calc_virial(0, mdatoms->homenr, x, f,
2003 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
2007 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
2009 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
2010 f, vir_force, mdatoms, enerd, lambda, t);
2013 /* Add the forces from enforced rotation potentials (if any) */
2016 wallcycle_start(wcycle, ewcROTadd);
2017 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
2018 wallcycle_stop(wcycle, ewcROTadd);
2021 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
2022 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
2024 if (PAR(cr) && !(cr->duty & DUTY_PME))
2026 /* In case of node-splitting, the PP nodes receive the long-range
2027 * forces, virial and energy from the PME nodes here.
2029 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
2034 post_process_forces(cr, step, nrnb, wcycle,
2035 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
2039 /* Sum the potential energy terms from group contributions */
2040 sum_epot(&(enerd->grpp), enerd->term);
2043 void do_force(FILE *fplog, t_commrec *cr,
2044 t_inputrec *inputrec,
2045 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2046 gmx_localtop_t *top,
2047 gmx_groups_t *groups,
2048 matrix box, rvec x[], history_t *hist,
2052 gmx_enerdata_t *enerd, t_fcdata *fcd,
2053 real *lambda, t_graph *graph,
2055 gmx_vsite_t *vsite, rvec mu_tot,
2056 double t, FILE *field, gmx_edsam_t ed,
2057 gmx_bool bBornRadii,
2060 /* modify force flag if not doing nonbonded */
2061 if (!fr->bNonbonded)
2063 flags &= ~GMX_FORCE_NONBONDED;
2066 switch (inputrec->cutoff_scheme)
2069 do_force_cutsVERLET(fplog, cr, inputrec,
2085 do_force_cutsGROUP(fplog, cr, inputrec,
2100 gmx_incons("Invalid cut-off scheme passed!");
2105 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2106 t_inputrec *ir, t_mdatoms *md,
2107 t_state *state, t_commrec *cr, t_nrnb *nrnb,
2108 t_forcerec *fr, gmx_localtop_t *top)
2110 int i, m, start, end;
2112 real dt = ir->delta_t;
2116 snew(savex, state->natoms);
2123 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2124 start, md->homenr, end);
2126 /* Do a first constrain to reset particles... */
2127 step = ir->init_step;
2130 char buf[STEPSTRSIZE];
2131 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2132 gmx_step_str(step, buf));
2136 /* constrain the current position */
2137 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2138 ir, NULL, cr, step, 0, 1.0, md,
2139 state->x, state->x, NULL,
2140 fr->bMolPBC, state->box,
2141 state->lambda[efptBONDED], &dvdl_dum,
2142 NULL, NULL, nrnb, econqCoord,
2143 ir->epc == epcMTTK, state->veta, state->veta);
2146 /* constrain the inital velocity, and save it */
2147 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2148 /* might not yet treat veta correctly */
2149 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2150 ir, NULL, cr, step, 0, 1.0, md,
2151 state->x, state->v, state->v,
2152 fr->bMolPBC, state->box,
2153 state->lambda[efptBONDED], &dvdl_dum,
2154 NULL, NULL, nrnb, econqVeloc,
2155 ir->epc == epcMTTK, state->veta, state->veta);
2157 /* constrain the inital velocities at t-dt/2 */
2158 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2160 for (i = start; (i < end); i++)
2162 for (m = 0; (m < DIM); m++)
2164 /* Reverse the velocity */
2165 state->v[i][m] = -state->v[i][m];
2166 /* Store the position at t-dt in buf */
2167 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2170 /* Shake the positions at t=-dt with the positions at t=0
2171 * as reference coordinates.
2175 char buf[STEPSTRSIZE];
2176 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2177 gmx_step_str(step, buf));
2180 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2181 ir, NULL, cr, step, -1, 1.0, md,
2182 state->x, savex, NULL,
2183 fr->bMolPBC, state->box,
2184 state->lambda[efptBONDED], &dvdl_dum,
2185 state->v, NULL, nrnb, econqCoord,
2186 ir->epc == epcMTTK, state->veta, state->veta);
2188 for (i = start; i < end; i++)
2190 for (m = 0; m < DIM; m++)
2192 /* Re-reverse the velocities */
2193 state->v[i][m] = -state->v[i][m];
2202 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2203 double *enerout, double *virout)
2205 double enersum, virsum;
2206 double invscale, invscale2, invscale3;
2207 double r, ea, eb, ec, pa, pb, pc, pd;
2209 int ri, offset, tabfactor;
2211 invscale = 1.0/scale;
2212 invscale2 = invscale*invscale;
2213 invscale3 = invscale*invscale2;
2215 /* Following summation derived from cubic spline definition,
2216 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2217 * the cubic spline. We first calculate the negative of the
2218 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2219 * add the more standard, abrupt cutoff correction to that result,
2220 * yielding the long-range correction for a switched function. We
2221 * perform both the pressure and energy loops at the same time for
2222 * simplicity, as the computational cost is low. */
2226 /* Since the dispersion table has been scaled down a factor
2227 * 6.0 and the repulsion a factor 12.0 to compensate for the
2228 * c6/c12 parameters inside nbfp[] being scaled up (to save
2229 * flops in kernels), we need to correct for this.
2240 for (ri = rstart; ri < rend; ++ri)
2244 eb = 2.0*invscale2*r;
2248 pb = 3.0*invscale2*r;
2249 pc = 3.0*invscale*r*r;
2252 /* this "8" is from the packing in the vdwtab array - perhaps
2253 should be defined? */
2255 offset = 8*ri + offstart;
2256 y0 = vdwtab[offset];
2257 f = vdwtab[offset+1];
2258 g = vdwtab[offset+2];
2259 h = vdwtab[offset+3];
2261 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);
2262 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);
2264 *enerout = 4.0*M_PI*enersum*tabfactor;
2265 *virout = 4.0*M_PI*virsum*tabfactor;
2268 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2270 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2271 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2272 double invscale, invscale2, invscale3;
2273 int ri0, ri1, ri, i, offstart, offset;
2274 real scale, *vdwtab, tabfactor, tmp;
2276 fr->enershiftsix = 0;
2277 fr->enershifttwelve = 0;
2278 fr->enerdiffsix = 0;
2279 fr->enerdifftwelve = 0;
2281 fr->virdifftwelve = 0;
2283 if (eDispCorr != edispcNO)
2285 for (i = 0; i < 2; i++)
2290 if ((fr->vdw_modifier == eintmodPOTSHIFT) ||
2291 (fr->vdw_modifier == eintmodPOTSWITCH) ||
2292 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2293 (fr->vdwtype == evdwSHIFT) ||
2294 (fr->vdwtype == evdwSWITCH))
2296 if (((fr->vdw_modifier == eintmodPOTSWITCH) ||
2297 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2298 (fr->vdwtype == evdwSWITCH)) && fr->rvdw_switch == 0)
2301 "With dispersion correction rvdw-switch can not be zero "
2302 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2305 scale = fr->nblists[0].table_vdw.scale;
2306 vdwtab = fr->nblists[0].table_vdw.data;
2308 /* Round the cut-offs to exact table values for precision */
2309 ri0 = floor(fr->rvdw_switch*scale);
2310 ri1 = ceil(fr->rvdw*scale);
2312 /* The code below has some support for handling force-switching, i.e.
2313 * when the force (instead of potential) is switched over a limited
2314 * region. This leads to a constant shift in the potential inside the
2315 * switching region, which we can handle by adding a constant energy
2316 * term in the force-switch case just like when we do potential-shift.
2318 * For now this is not enabled, but to keep the functionality in the
2319 * code we check separately for switch and shift. When we do force-switch
2320 * the shifting point is rvdw_switch, while it is the cutoff when we
2321 * have a classical potential-shift.
2323 * For a pure potential-shift the potential has a constant shift
2324 * all the way out to the cutoff, and that is it. For other forms
2325 * we need to calculate the constant shift up to the point where we
2326 * start modifying the potential.
2328 ri0 = (fr->vdw_modifier == eintmodPOTSHIFT) ? ri1 : ri0;
2335 if ((fr->vdw_modifier == eintmodFORCESWITCH) ||
2336 (fr->vdwtype == evdwSHIFT))
2338 /* Determine the constant energy shift below rvdw_switch.
2339 * Table has a scale factor since we have scaled it down to compensate
2340 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2342 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2343 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2345 else if (fr->vdw_modifier == eintmodPOTSHIFT)
2347 fr->enershiftsix = (real)(-1.0/(rc3*rc3));
2348 fr->enershifttwelve = (real)( 1.0/(rc9*rc3));
2351 /* Add the constant part from 0 to rvdw_switch.
2352 * This integration from 0 to rvdw_switch overcounts the number
2353 * of interactions by 1, as it also counts the self interaction.
2354 * We will correct for this later.
2356 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2357 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2359 /* Calculate the contribution in the range [r0,r1] where we
2360 * modify the potential. For a pure potential-shift modifier we will
2361 * have ri0==ri1, and there will not be any contribution here.
2363 for (i = 0; i < 2; i++)
2367 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2368 eners[i] -= enersum;
2372 /* Alright: Above we compensated by REMOVING the parts outside r0
2373 * corresponding to the ideal VdW 1/r6 and /r12 potentials.
2375 * Regardless of whether r0 is the point where we start switching,
2376 * or the cutoff where we calculated the constant shift, we include
2377 * all the parts we are missing out to infinity from r0 by
2378 * calculating the analytical dispersion correction.
2380 eners[0] += -4.0*M_PI/(3.0*rc3);
2381 eners[1] += 4.0*M_PI/(9.0*rc9);
2382 virs[0] += 8.0*M_PI/rc3;
2383 virs[1] += -16.0*M_PI/(3.0*rc9);
2385 else if (fr->vdwtype == evdwCUT ||
2386 EVDW_PME(fr->vdwtype) ||
2387 fr->vdwtype == evdwUSER)
2389 if (fr->vdwtype == evdwUSER && fplog)
2392 "WARNING: using dispersion correction with user tables\n");
2395 /* Note that with LJ-PME, the dispersion correction is multiplied
2396 * by the difference between the actual C6 and the value of C6
2397 * that would produce the combination rule.
2398 * This means the normal energy and virial difference formulas
2402 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2404 /* Contribution beyond the cut-off */
2405 eners[0] += -4.0*M_PI/(3.0*rc3);
2406 eners[1] += 4.0*M_PI/(9.0*rc9);
2407 if (fr->vdw_modifier == eintmodPOTSHIFT)
2409 /* Contribution within the cut-off */
2410 eners[0] += -4.0*M_PI/(3.0*rc3);
2411 eners[1] += 4.0*M_PI/(3.0*rc9);
2413 /* Contribution beyond the cut-off */
2414 virs[0] += 8.0*M_PI/rc3;
2415 virs[1] += -16.0*M_PI/(3.0*rc9);
2420 "Dispersion correction is not implemented for vdw-type = %s",
2421 evdw_names[fr->vdwtype]);
2424 /* When we deprecate the group kernels the code below can go too */
2425 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2427 /* Calculate self-interaction coefficient (assuming that
2428 * the reciprocal-space contribution is constant in the
2429 * region that contributes to the self-interaction).
2431 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2433 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2434 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2437 fr->enerdiffsix = eners[0];
2438 fr->enerdifftwelve = eners[1];
2439 /* The 0.5 is due to the Gromacs definition of the virial */
2440 fr->virdiffsix = 0.5*virs[0];
2441 fr->virdifftwelve = 0.5*virs[1];
2445 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2446 gmx_int64_t step, int natoms,
2447 matrix box, real lambda, tensor pres, tensor virial,
2448 real *prescorr, real *enercorr, real *dvdlcorr)
2450 gmx_bool bCorrAll, bCorrPres;
2451 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2461 if (ir->eDispCorr != edispcNO)
2463 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2464 ir->eDispCorr == edispcAllEnerPres);
2465 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2466 ir->eDispCorr == edispcAllEnerPres);
2468 invvol = 1/det(box);
2471 /* Only correct for the interactions with the inserted molecule */
2472 dens = (natoms - fr->n_tpi)*invvol;
2477 dens = natoms*invvol;
2478 ninter = 0.5*natoms;
2481 if (ir->efep == efepNO)
2483 avcsix = fr->avcsix[0];
2484 avctwelve = fr->avctwelve[0];
2488 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2489 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2492 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2493 *enercorr += avcsix*enerdiff;
2495 if (ir->efep != efepNO)
2497 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2501 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2502 *enercorr += avctwelve*enerdiff;
2503 if (fr->efep != efepNO)
2505 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2511 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2512 if (ir->eDispCorr == edispcAllEnerPres)
2514 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2516 /* The factor 2 is because of the Gromacs virial definition */
2517 spres = -2.0*invvol*svir*PRESFAC;
2519 for (m = 0; m < DIM; m++)
2521 virial[m][m] += svir;
2522 pres[m][m] += spres;
2527 /* Can't currently control when it prints, for now, just print when degugging */
2532 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2538 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2539 *enercorr, spres, svir);
2543 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2547 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2549 gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
2551 if (fr->efep != efepNO)
2553 *dvdlcorr += dvdlambda;
2558 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2559 t_graph *graph, rvec x[])
2563 fprintf(fplog, "Removing pbc first time\n");
2565 calc_shifts(box, fr->shift_vec);
2568 mk_mshift(fplog, graph, fr->ePBC, box, x);
2571 p_graph(debug, "do_pbc_first 1", graph);
2573 shift_self(graph, box, x);
2574 /* By doing an extra mk_mshift the molecules that are broken
2575 * because they were e.g. imported from another software
2576 * will be made whole again. Such are the healing powers
2579 mk_mshift(fplog, graph, fr->ePBC, box, x);
2582 p_graph(debug, "do_pbc_first 2", graph);
2587 fprintf(fplog, "Done rmpbc\n");
2591 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2592 gmx_mtop_t *mtop, rvec x[],
2597 gmx_molblock_t *molb;
2599 if (bFirst && fplog)
2601 fprintf(fplog, "Removing pbc first time\n");
2606 for (mb = 0; mb < mtop->nmolblock; mb++)
2608 molb = &mtop->molblock[mb];
2609 if (molb->natoms_mol == 1 ||
2610 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2612 /* Just one atom or charge group in the molecule, no PBC required */
2613 as += molb->nmol*molb->natoms_mol;
2617 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2618 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2619 0, molb->natoms_mol, FALSE, FALSE, graph);
2621 for (mol = 0; mol < molb->nmol; mol++)
2623 mk_mshift(fplog, graph, ePBC, box, x+as);
2625 shift_self(graph, box, x+as);
2626 /* The molecule is whole now.
2627 * We don't need the second mk_mshift call as in do_pbc_first,
2628 * since we no longer need this graph.
2631 as += molb->natoms_mol;
2639 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2640 gmx_mtop_t *mtop, rvec x[])
2642 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2645 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2646 gmx_mtop_t *mtop, rvec x[])
2648 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2651 void finish_run(FILE *fplog, t_commrec *cr,
2652 t_inputrec *inputrec,
2653 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2654 gmx_walltime_accounting_t walltime_accounting,
2655 nonbonded_verlet_t *nbv,
2656 gmx_bool bWriteStat)
2659 t_nrnb *nrnb_tot = NULL;
2662 double elapsed_time,
2663 elapsed_time_over_all_ranks,
2664 elapsed_time_over_all_threads,
2665 elapsed_time_over_all_threads_over_all_ranks;
2666 wallcycle_sum(cr, wcycle);
2672 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2673 cr->mpi_comm_mysim);
2681 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2682 elapsed_time_over_all_ranks = elapsed_time;
2683 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2684 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2688 /* reduce elapsed_time over all MPI ranks in the current simulation */
2689 MPI_Allreduce(&elapsed_time,
2690 &elapsed_time_over_all_ranks,
2691 1, MPI_DOUBLE, MPI_SUM,
2692 cr->mpi_comm_mysim);
2693 elapsed_time_over_all_ranks /= cr->nnodes;
2694 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2695 * current simulation. */
2696 MPI_Allreduce(&elapsed_time_over_all_threads,
2697 &elapsed_time_over_all_threads_over_all_ranks,
2698 1, MPI_DOUBLE, MPI_SUM,
2699 cr->mpi_comm_mysim);
2705 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2712 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2714 print_dd_statistics(cr, inputrec, fplog);
2719 wallclock_gpu_t* gputimes = use_GPU(nbv) ?
2720 nbnxn_cuda_get_timings(nbv->cu_nbv) : NULL;
2721 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2722 elapsed_time_over_all_ranks,
2725 if (EI_DYNAMICS(inputrec->eI))
2727 delta_t = inputrec->delta_t;
2736 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2737 elapsed_time_over_all_ranks,
2738 walltime_accounting_get_nsteps_done(walltime_accounting),
2739 delta_t, nbfs, mflop);
2743 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2744 elapsed_time_over_all_ranks,
2745 walltime_accounting_get_nsteps_done(walltime_accounting),
2746 delta_t, nbfs, mflop);
2751 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2753 /* this function works, but could probably use a logic rewrite to keep all the different
2754 types of efep straight. */
2757 t_lambda *fep = ir->fepvals;
2759 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2761 for (i = 0; i < efptNR; i++)
2773 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2774 if checkpoint is set -- a kludge is in for now
2776 for (i = 0; i < efptNR; i++)
2778 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2779 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2781 lambda[i] = fep->init_lambda;
2784 lam0[i] = lambda[i];
2789 lambda[i] = fep->all_lambda[i][*fep_state];
2792 lam0[i] = lambda[i];
2798 /* need to rescale control temperatures to match current state */
2799 for (i = 0; i < ir->opts.ngtc; i++)
2801 if (ir->opts.ref_t[i] > 0)
2803 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2809 /* Send to the log the information on the current lambdas */
2812 fprintf(fplog, "Initial vector of lambda components:[ ");
2813 for (i = 0; i < efptNR; i++)
2815 fprintf(fplog, "%10.4f ", lambda[i]);
2817 fprintf(fplog, "]\n");
2823 void init_md(FILE *fplog,
2824 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2825 double *t, double *t0,
2826 real *lambda, int *fep_state, double *lam0,
2827 t_nrnb *nrnb, gmx_mtop_t *mtop,
2829 int nfile, const t_filenm fnm[],
2830 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2831 tensor force_vir, tensor shake_vir, rvec mu_tot,
2832 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags)
2837 /* Initial values */
2838 *t = *t0 = ir->init_t;
2841 for (i = 0; i < ir->opts.ngtc; i++)
2843 /* set bSimAnn if any group is being annealed */
2844 if (ir->opts.annealing[i] != eannNO)
2851 update_annealing_target_temp(&(ir->opts), ir->init_t);
2854 /* Initialize lambda variables */
2855 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2859 *upd = init_update(ir);
2865 *vcm = init_vcm(fplog, &mtop->groups, ir);
2868 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2870 if (ir->etc == etcBERENDSEN)
2872 please_cite(fplog, "Berendsen84a");
2874 if (ir->etc == etcVRESCALE)
2876 please_cite(fplog, "Bussi2007a");
2878 if (ir->eI == eiSD1)
2880 please_cite(fplog, "Goga2012");
2888 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv);
2890 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2891 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2896 please_cite(fplog, "Fritsch12");
2897 please_cite(fplog, "Junghans10");
2899 /* Initiate variables */
2900 clear_mat(force_vir);
2901 clear_mat(shake_vir);