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
47 #include "gromacs/legacyheaders/typedefs.h"
48 #include "gromacs/utility/cstringutil.h"
49 #include "gromacs/legacyheaders/names.h"
50 #include "gromacs/legacyheaders/txtdump.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "gromacs/legacyheaders/chargegroup.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/legacyheaders/nrnb.h"
55 #include "gromacs/legacyheaders/mdrun.h"
56 #include "gromacs/legacyheaders/sim_util.h"
57 #include "gromacs/legacyheaders/update.h"
58 #include "gromacs/math/units.h"
59 #include "gromacs/legacyheaders/mdatoms.h"
60 #include "gromacs/legacyheaders/force.h"
61 #include "gromacs/legacyheaders/bondf.h"
62 #include "gromacs/legacyheaders/pme.h"
63 #include "gromacs/legacyheaders/disre.h"
64 #include "gromacs/legacyheaders/orires.h"
65 #include "gromacs/legacyheaders/network.h"
66 #include "gromacs/legacyheaders/calcmu.h"
67 #include "gromacs/legacyheaders/constr.h"
68 #include "gromacs/legacyheaders/copyrite.h"
69 #include "gromacs/legacyheaders/domdec.h"
70 #include "gromacs/legacyheaders/genborn.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 "gromacs/legacyheaders/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"
94 #include "gromacs/legacyheaders/qmmm.h"
96 #include "gromacs/legacyheaders/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(int flags,
320 matrix box, rvec x[],
321 gmx_enerdata_t *enerd,
329 /* Position restraints always require full pbc */
330 set_pbc(&pbc, ir->ePBC, box);
332 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
333 top->idef.iparams_posres,
334 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
335 ir->ePBC == epbcNONE ? NULL : &pbc,
336 lambda[efptRESTRAINT], &dvdl,
337 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
338 enerd->term[F_POSRES] += v;
339 /* If just the force constant changes, the FEP term is linear,
340 * but if k changes, it is not.
342 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
343 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
345 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
347 for (i = 0; i < enerd->n_lambda; i++)
349 real dvdl_dum, lambda_dum;
351 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
352 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
353 top->idef.iparams_posres,
354 (const rvec*)x, NULL, NULL,
355 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
356 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
357 enerd->enerpart_lambda[i] += v;
362 static void fbposres_wrapper(t_inputrec *ir,
365 matrix box, rvec x[],
366 gmx_enerdata_t *enerd,
372 /* Flat-bottomed position restraints always require full pbc */
373 set_pbc(&pbc, ir->ePBC, box);
374 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
375 top->idef.iparams_fbposres,
376 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
377 ir->ePBC == epbcNONE ? NULL : &pbc,
378 fr->rc_scaling, fr->ePBC, fr->posres_com);
379 enerd->term[F_FBPOSRES] += v;
380 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
383 static void pull_potential_wrapper(t_commrec *cr,
385 matrix box, rvec x[],
389 gmx_enerdata_t *enerd,
396 /* Calculate the center of mass forces, this requires communication,
397 * which is why pull_potential is called close to other communication.
398 * The virial contribution is calculated directly,
399 * which is why we call pull_potential after calc_virial.
401 set_pbc(&pbc, ir->ePBC, box);
403 enerd->term[F_COM_PULL] +=
404 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
405 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
406 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
409 static void pme_receive_force_ener(t_commrec *cr,
410 gmx_wallcycle_t wcycle,
411 gmx_enerdata_t *enerd,
414 real e_q, e_lj, v, dvdl_q, dvdl_lj;
415 float cycles_ppdpme, cycles_seppme;
417 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
418 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
420 /* In case of node-splitting, the PP nodes receive the long-range
421 * forces, virial and energy from the PME nodes here.
423 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
426 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
427 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
429 enerd->term[F_COUL_RECIP] += e_q;
430 enerd->term[F_LJ_RECIP] += e_lj;
431 enerd->dvdl_lin[efptCOUL] += dvdl_q;
432 enerd->dvdl_lin[efptVDW] += dvdl_lj;
436 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
438 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
441 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
442 gmx_int64_t step, real pforce, rvec *x, rvec *f)
446 char buf[STEPSTRSIZE];
449 for (i = 0; i < md->homenr; i++)
452 /* We also catch NAN, if the compiler does not optimize this away. */
453 if (fn2 >= pf2 || fn2 != fn2)
455 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
456 gmx_step_str(step, buf),
457 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
462 static void post_process_forces(t_commrec *cr,
464 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
466 matrix box, rvec x[],
471 t_forcerec *fr, gmx_vsite_t *vsite,
478 /* Spread the mesh force on virtual sites to the other particles...
479 * This is parallellized. MPI communication is performed
480 * if the constructing atoms aren't local.
482 wallcycle_start(wcycle, ewcVSITESPREAD);
483 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
484 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
486 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
487 wallcycle_stop(wcycle, ewcVSITESPREAD);
489 if (flags & GMX_FORCE_VIRIAL)
491 /* Now add the forces, this is local */
494 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
498 sum_forces(0, mdatoms->homenr,
501 if (EEL_FULL(fr->eeltype))
503 /* Add the mesh contribution to the virial */
504 m_add(vir_force, fr->vir_el_recip, vir_force);
506 if (EVDW_PME(fr->vdwtype))
508 /* Add the mesh contribution to the virial */
509 m_add(vir_force, fr->vir_lj_recip, vir_force);
513 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
518 if (fr->print_force >= 0)
520 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
524 static void do_nb_verlet(t_forcerec *fr,
525 interaction_const_t *ic,
526 gmx_enerdata_t *enerd,
527 int flags, int ilocality,
530 gmx_wallcycle_t wcycle)
532 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
534 nonbonded_verlet_group_t *nbvg;
537 if (!(flags & GMX_FORCE_NONBONDED))
539 /* skip non-bonded calculation */
543 nbvg = &fr->nbv->grp[ilocality];
545 /* CUDA kernel launch overhead is already timed separately */
546 if (fr->cutoff_scheme != ecutsVERLET)
548 gmx_incons("Invalid cut-off scheme passed!");
551 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
555 wallcycle_sub_start(wcycle, ewcsNONBONDED);
557 switch (nbvg->kernel_type)
559 case nbnxnk4x4_PlainC:
560 nbnxn_kernel_ref(&nbvg->nbl_lists,
566 enerd->grpp.ener[egCOULSR],
568 enerd->grpp.ener[egBHAMSR] :
569 enerd->grpp.ener[egLJSR]);
572 case nbnxnk4xN_SIMD_4xN:
573 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
580 enerd->grpp.ener[egCOULSR],
582 enerd->grpp.ener[egBHAMSR] :
583 enerd->grpp.ener[egLJSR]);
585 case nbnxnk4xN_SIMD_2xNN:
586 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
593 enerd->grpp.ener[egCOULSR],
595 enerd->grpp.ener[egBHAMSR] :
596 enerd->grpp.ener[egLJSR]);
599 case nbnxnk8x8x8_CUDA:
600 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
603 case nbnxnk8x8x8_PlainC:
604 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
609 nbvg->nbat->out[0].f,
611 enerd->grpp.ener[egCOULSR],
613 enerd->grpp.ener[egBHAMSR] :
614 enerd->grpp.ener[egLJSR]);
618 gmx_incons("Invalid nonbonded kernel type passed!");
623 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
626 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
628 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
630 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
631 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
633 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
637 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
639 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
640 if (flags & GMX_FORCE_ENERGY)
642 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
643 enr_nbnxn_kernel_ljc += 1;
644 enr_nbnxn_kernel_lj += 1;
647 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
648 nbvg->nbl_lists.natpair_ljq);
649 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
650 nbvg->nbl_lists.natpair_lj);
651 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
652 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
653 nbvg->nbl_lists.natpair_q);
655 if (ic->vdw_modifier == eintmodFORCESWITCH)
657 /* We add up the switch cost separately */
658 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
659 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
661 if (ic->vdw_modifier == eintmodPOTSWITCH)
663 /* We add up the switch cost separately */
664 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
665 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
667 if (ic->vdwtype == evdwPME)
669 /* We add up the LJ Ewald cost separately */
670 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
671 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
675 static void do_nb_verlet_fep(nbnxn_pairlist_set_t *nbl_lists,
682 gmx_enerdata_t *enerd,
685 gmx_wallcycle_t wcycle)
688 nb_kernel_data_t kernel_data;
690 real dvdl_nb[efptNR];
695 /* Add short-range interactions */
696 donb_flags |= GMX_NONBONDED_DO_SR;
698 /* Currently all group scheme kernels always calculate (shift-)forces */
699 if (flags & GMX_FORCE_FORCES)
701 donb_flags |= GMX_NONBONDED_DO_FORCE;
703 if (flags & GMX_FORCE_VIRIAL)
705 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
707 if (flags & GMX_FORCE_ENERGY)
709 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
711 if (flags & GMX_FORCE_DO_LR)
713 donb_flags |= GMX_NONBONDED_DO_LR;
716 kernel_data.flags = donb_flags;
717 kernel_data.lambda = lambda;
718 kernel_data.dvdl = dvdl_nb;
720 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
721 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
723 /* reset free energy components */
724 for (i = 0; i < efptNR; i++)
729 assert(gmx_omp_nthreads_get(emntNonbonded) == nbl_lists->nnbl);
731 wallcycle_sub_start(wcycle, ewcsNONBONDED);
732 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
733 for (th = 0; th < nbl_lists->nnbl; th++)
735 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
736 x, f, fr, mdatoms, &kernel_data, nrnb);
739 if (fepvals->sc_alpha != 0)
741 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
742 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
746 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
747 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
750 /* If we do foreign lambda and we have soft-core interactions
751 * we have to recalculate the (non-linear) energies contributions.
753 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
755 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
756 kernel_data.lambda = lam_i;
757 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
758 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
759 /* Note that we add to kernel_data.dvdl, but ignore the result */
761 for (i = 0; i < enerd->n_lambda; i++)
763 for (j = 0; j < efptNR; j++)
765 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
767 reset_foreign_enerdata(enerd);
768 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
769 for (th = 0; th < nbl_lists->nnbl; th++)
771 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
772 x, f, fr, mdatoms, &kernel_data, nrnb);
775 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
776 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
780 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
783 gmx_bool use_GPU(const nonbonded_verlet_t *nbv)
785 return nbv != NULL && nbv->bUseGPU;
788 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
789 t_inputrec *inputrec,
790 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
792 gmx_groups_t gmx_unused *groups,
793 matrix box, rvec x[], history_t *hist,
797 gmx_enerdata_t *enerd, t_fcdata *fcd,
798 real *lambda, t_graph *graph,
799 t_forcerec *fr, interaction_const_t *ic,
800 gmx_vsite_t *vsite, rvec mu_tot,
801 double t, FILE *field, gmx_edsam_t ed,
809 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
810 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
811 gmx_bool bDiffKernels = FALSE;
813 rvec vzero, box_diag;
815 float cycles_pme, cycles_force, cycles_wait_gpu;
816 nonbonded_verlet_t *nbv;
821 nb_kernel_type = fr->nbv->grp[0].kernel_type;
824 homenr = mdatoms->homenr;
826 clear_mat(vir_force);
829 if (DOMAINDECOMP(cr))
831 cg1 = cr->dd->ncg_tot;
842 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
843 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
844 bFillGrid = (bNS && bStateChanged);
845 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
846 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
847 bDoForces = (flags & GMX_FORCE_FORCES);
848 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
849 bUseGPU = fr->nbv->bUseGPU;
850 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
854 update_forcerec(fr, box);
856 if (NEED_MUTOT(*inputrec))
858 /* Calculate total (local) dipole moment in a temporary common array.
859 * This makes it possible to sum them over nodes faster.
861 calc_mu(start, homenr,
862 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
867 if (fr->ePBC != epbcNONE)
869 /* Compute shift vectors every step,
870 * because of pressure coupling or box deformation!
872 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
874 calc_shifts(box, fr->shift_vec);
879 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
880 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
882 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
884 unshift_self(graph, box, x);
888 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
889 fr->shift_vec, nbv->grp[0].nbat);
892 if (!(cr->duty & DUTY_PME))
894 /* Send particle coordinates to the pme nodes.
895 * Since this is only implemented for domain decomposition
896 * and domain decomposition does not use the graph,
897 * we do not need to worry about shifting.
902 wallcycle_start(wcycle, ewcPP_PMESENDX);
904 bBS = (inputrec->nwall == 2);
908 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
911 if (EEL_PME(fr->eeltype))
913 pme_flags |= GMX_PME_DO_COULOMB;
916 if (EVDW_PME(fr->vdwtype))
918 pme_flags |= GMX_PME_DO_LJ;
921 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
922 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
923 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
926 wallcycle_stop(wcycle, ewcPP_PMESENDX);
930 /* do gridding for pair search */
933 if (graph && bStateChanged)
935 /* Calculate intramolecular shift vectors to make molecules whole */
936 mk_mshift(fplog, graph, fr->ePBC, box, x);
940 box_diag[XX] = box[XX][XX];
941 box_diag[YY] = box[YY][YY];
942 box_diag[ZZ] = box[ZZ][ZZ];
944 wallcycle_start(wcycle, ewcNS);
947 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
948 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
950 0, mdatoms->homenr, -1, fr->cginfo, x,
952 nbv->grp[eintLocal].kernel_type,
953 nbv->grp[eintLocal].nbat);
954 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
958 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
959 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
961 nbv->grp[eintNonlocal].kernel_type,
962 nbv->grp[eintNonlocal].nbat);
963 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
966 if (nbv->ngrp == 1 ||
967 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
969 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
970 nbv->nbs, mdatoms, fr->cginfo);
974 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
975 nbv->nbs, mdatoms, fr->cginfo);
976 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
977 nbv->nbs, mdatoms, fr->cginfo);
979 wallcycle_stop(wcycle, ewcNS);
982 /* initialize the GPU atom data and copy shift vector */
987 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
988 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
989 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
992 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
993 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
994 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
997 /* do local pair search */
1000 wallcycle_start_nocount(wcycle, ewcNS);
1001 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1002 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
1005 nbv->min_ci_balanced,
1006 &nbv->grp[eintLocal].nbl_lists,
1008 nbv->grp[eintLocal].kernel_type,
1010 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1014 /* initialize local pair-list on the GPU */
1015 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1016 nbv->grp[eintLocal].nbl_lists.nbl[0],
1019 wallcycle_stop(wcycle, ewcNS);
1023 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1024 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1025 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
1026 nbv->grp[eintLocal].nbat);
1027 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1028 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1033 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1034 /* launch local nonbonded F on GPU */
1035 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
1037 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1040 /* Communicate coordinates and sum dipole if necessary +
1041 do non-local pair search */
1042 if (DOMAINDECOMP(cr))
1044 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
1045 nbv->grp[eintLocal].kernel_type);
1049 /* With GPU+CPU non-bonded calculations we need to copy
1050 * the local coordinates to the non-local nbat struct
1051 * (in CPU format) as the non-local kernel call also
1052 * calculates the local - non-local interactions.
1054 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1055 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1056 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
1057 nbv->grp[eintNonlocal].nbat);
1058 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1059 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1064 wallcycle_start_nocount(wcycle, ewcNS);
1065 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1069 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
1072 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
1075 nbv->min_ci_balanced,
1076 &nbv->grp[eintNonlocal].nbl_lists,
1078 nbv->grp[eintNonlocal].kernel_type,
1081 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1083 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1085 /* initialize non-local pair-list on the GPU */
1086 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1087 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1090 wallcycle_stop(wcycle, ewcNS);
1094 wallcycle_start(wcycle, ewcMOVEX);
1095 dd_move_x(cr->dd, box, x);
1097 /* When we don't need the total dipole we sum it in global_stat */
1098 if (bStateChanged && NEED_MUTOT(*inputrec))
1100 gmx_sumd(2*DIM, mu, cr);
1102 wallcycle_stop(wcycle, ewcMOVEX);
1104 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1105 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1106 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1107 nbv->grp[eintNonlocal].nbat);
1108 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1109 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1112 if (bUseGPU && !bDiffKernels)
1114 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1115 /* launch non-local nonbonded F on GPU */
1116 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1118 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1124 /* launch D2H copy-back F */
1125 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1126 if (DOMAINDECOMP(cr) && !bDiffKernels)
1128 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1129 flags, eatNonlocal);
1131 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1133 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1136 if (bStateChanged && NEED_MUTOT(*inputrec))
1140 gmx_sumd(2*DIM, mu, cr);
1143 for (i = 0; i < 2; i++)
1145 for (j = 0; j < DIM; j++)
1147 fr->mu_tot[i][j] = mu[i*DIM + j];
1151 if (fr->efep == efepNO)
1153 copy_rvec(fr->mu_tot[0], mu_tot);
1157 for (j = 0; j < DIM; j++)
1160 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1161 lambda[efptCOUL]*fr->mu_tot[1][j];
1165 /* Reset energies */
1166 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1167 clear_rvecs(SHIFTS, fr->fshift);
1169 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1171 wallcycle_start(wcycle, ewcPPDURINGPME);
1172 dd_force_flop_start(cr->dd, nrnb);
1177 /* Enforced rotation has its own cycle counter that starts after the collective
1178 * coordinates have been communicated. It is added to ddCyclF to allow
1179 * for proper load-balancing */
1180 wallcycle_start(wcycle, ewcROT);
1181 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1182 wallcycle_stop(wcycle, ewcROT);
1185 /* Start the force cycle counter.
1186 * This counter is stopped in do_forcelow_level.
1187 * No parallel communication should occur while this counter is running,
1188 * since that will interfere with the dynamic load balancing.
1190 wallcycle_start(wcycle, ewcFORCE);
1193 /* Reset forces for which the virial is calculated separately:
1194 * PME/Ewald forces if necessary */
1195 if (fr->bF_NoVirSum)
1197 if (flags & GMX_FORCE_VIRIAL)
1199 fr->f_novirsum = fr->f_novirsum_alloc;
1202 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1206 clear_rvecs(homenr, fr->f_novirsum+start);
1211 /* We are not calculating the pressure so we do not need
1212 * a separate array for forces that do not contribute
1219 /* Clear the short- and long-range forces */
1220 clear_rvecs(fr->natoms_force_constr, f);
1221 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1223 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1226 clear_rvec(fr->vir_diag_posres);
1229 if (inputrec->ePull == epullCONSTRAINT)
1231 clear_pull_forces(inputrec->pull);
1234 /* We calculate the non-bonded forces, when done on the CPU, here.
1235 * We do this before calling do_force_lowlevel, as in there bondeds
1236 * forces are calculated before PME, which does communication.
1237 * With this order, non-bonded and bonded force calculation imbalance
1238 * can be balanced out by the domain decomposition load balancing.
1243 /* Maybe we should move this into do_force_lowlevel */
1244 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1248 if (fr->efep != efepNO)
1250 /* Calculate the local and non-local free energy interactions here.
1251 * Happens here on the CPU both with and without GPU.
1253 if (fr->nbv->grp[eintLocal].nbl_lists.nbl_fep[0]->nrj > 0)
1255 do_nb_verlet_fep(&fr->nbv->grp[eintLocal].nbl_lists,
1257 inputrec->fepvals, lambda,
1258 enerd, flags, nrnb, wcycle);
1261 if (DOMAINDECOMP(cr) &&
1262 fr->nbv->grp[eintNonlocal].nbl_lists.nbl_fep[0]->nrj > 0)
1264 do_nb_verlet_fep(&fr->nbv->grp[eintNonlocal].nbl_lists,
1266 inputrec->fepvals, lambda,
1267 enerd, flags, nrnb, wcycle);
1271 if (!bUseOrEmulGPU || bDiffKernels)
1275 if (DOMAINDECOMP(cr))
1277 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1278 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1288 aloc = eintNonlocal;
1291 /* Add all the non-bonded force to the normal force array.
1292 * This can be split into a local a non-local part when overlapping
1293 * communication with calculation with domain decomposition.
1295 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1296 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1297 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1298 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1299 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1300 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1301 wallcycle_start_nocount(wcycle, ewcFORCE);
1303 /* if there are multiple fshift output buffers reduce them */
1304 if ((flags & GMX_FORCE_VIRIAL) &&
1305 nbv->grp[aloc].nbl_lists.nnbl > 1)
1307 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1312 /* update QMMMrec, if necessary */
1315 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1318 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1320 posres_wrapper(flags, inputrec, nrnb, top, box, x,
1324 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1326 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1329 /* Compute the bonded and non-bonded energies and optionally forces */
1330 do_force_lowlevel(fr, inputrec, &(top->idef),
1331 cr, nrnb, wcycle, mdatoms,
1332 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1333 &(top->atomtypes), bBornRadii, box,
1334 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1335 flags, &cycles_pme);
1339 if (do_per_step(step, inputrec->nstcalclr))
1341 /* Add the long range forces to the short range forces */
1342 for (i = 0; i < fr->natoms_force_constr; i++)
1344 rvec_add(fr->f_twin[i], f[i], f[i]);
1349 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1353 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1356 if (bUseOrEmulGPU && !bDiffKernels)
1358 /* wait for non-local forces (or calculate in emulation mode) */
1359 if (DOMAINDECOMP(cr))
1365 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1366 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1367 nbv->grp[eintNonlocal].nbat,
1369 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1371 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1372 cycles_wait_gpu += cycles_tmp;
1373 cycles_force += cycles_tmp;
1377 wallcycle_start_nocount(wcycle, ewcFORCE);
1378 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1380 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1382 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1383 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1384 /* skip the reduction if there was no non-local work to do */
1385 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1387 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1388 nbv->grp[eintNonlocal].nbat, f);
1390 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1391 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1395 if (bDoForces && DOMAINDECOMP(cr))
1397 /* Communicate the forces */
1398 wallcycle_start(wcycle, ewcMOVEF);
1399 dd_move_f(cr->dd, f, fr->fshift);
1400 /* Do we need to communicate the separate force array
1401 * for terms that do not contribute to the single sum virial?
1402 * Position restraints and electric fields do not introduce
1403 * inter-cg forces, only full electrostatics methods do.
1404 * When we do not calculate the virial, fr->f_novirsum = f,
1405 * so we have already communicated these forces.
1407 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1408 (flags & GMX_FORCE_VIRIAL))
1410 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1414 /* We should not update the shift forces here,
1415 * since f_twin is already included in f.
1417 dd_move_f(cr->dd, fr->f_twin, NULL);
1419 wallcycle_stop(wcycle, ewcMOVEF);
1424 /* wait for local forces (or calculate in emulation mode) */
1427 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1428 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1429 nbv->grp[eintLocal].nbat,
1431 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1433 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1435 /* now clear the GPU outputs while we finish the step on the CPU */
1437 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1438 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1439 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1443 wallcycle_start_nocount(wcycle, ewcFORCE);
1444 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1445 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1447 wallcycle_stop(wcycle, ewcFORCE);
1449 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1450 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1451 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1453 /* skip the reduction if there was no non-local work to do */
1454 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1455 nbv->grp[eintLocal].nbat, f);
1457 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1458 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1461 if (DOMAINDECOMP(cr))
1463 dd_force_flop_stop(cr->dd, nrnb);
1466 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1469 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1476 if (IR_ELEC_FIELD(*inputrec))
1478 /* Compute forces due to electric field */
1479 calc_f_el(MASTER(cr) ? field : NULL,
1480 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1481 inputrec->ex, inputrec->et, t);
1484 /* If we have NoVirSum forces, but we do not calculate the virial,
1485 * we sum fr->f_novirum=f later.
1487 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1489 wallcycle_start(wcycle, ewcVSITESPREAD);
1490 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1491 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1492 wallcycle_stop(wcycle, ewcVSITESPREAD);
1496 wallcycle_start(wcycle, ewcVSITESPREAD);
1497 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1499 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1500 wallcycle_stop(wcycle, ewcVSITESPREAD);
1504 if (flags & GMX_FORCE_VIRIAL)
1506 /* Calculation of the virial must be done after vsites! */
1507 calc_virial(0, mdatoms->homenr, x, f,
1508 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1512 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1514 pull_potential_wrapper(cr, inputrec, box, x,
1515 f, vir_force, mdatoms, enerd, lambda, t);
1518 /* Add the forces from enforced rotation potentials (if any) */
1521 wallcycle_start(wcycle, ewcROTadd);
1522 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1523 wallcycle_stop(wcycle, ewcROTadd);
1526 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1527 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1529 if (PAR(cr) && !(cr->duty & DUTY_PME))
1531 /* In case of node-splitting, the PP nodes receive the long-range
1532 * forces, virial and energy from the PME nodes here.
1534 pme_receive_force_ener(cr, wcycle, enerd, fr);
1539 post_process_forces(cr, step, nrnb, wcycle,
1540 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1544 /* Sum the potential energy terms from group contributions */
1545 sum_epot(&(enerd->grpp), enerd->term);
1548 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1549 t_inputrec *inputrec,
1550 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1551 gmx_localtop_t *top,
1552 gmx_groups_t *groups,
1553 matrix box, rvec x[], history_t *hist,
1557 gmx_enerdata_t *enerd, t_fcdata *fcd,
1558 real *lambda, t_graph *graph,
1559 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1560 double t, FILE *field, gmx_edsam_t ed,
1561 gmx_bool bBornRadii,
1567 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1568 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1569 gmx_bool bDoAdressWF;
1571 rvec vzero, box_diag;
1572 real e, v, dvdlambda[efptNR];
1574 float cycles_pme, cycles_force;
1577 homenr = mdatoms->homenr;
1579 clear_mat(vir_force);
1582 if (DOMAINDECOMP(cr))
1584 cg1 = cr->dd->ncg_tot;
1595 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1596 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1597 /* Should we update the long-range neighborlists at this step? */
1598 bDoLongRangeNS = fr->bTwinRange && bNS;
1599 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1600 bFillGrid = (bNS && bStateChanged);
1601 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1602 bDoForces = (flags & GMX_FORCE_FORCES);
1603 bDoPotential = (flags & GMX_FORCE_ENERGY);
1604 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1605 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1607 /* should probably move this to the forcerec since it doesn't change */
1608 bDoAdressWF = ((fr->adress_type != eAdressOff));
1612 update_forcerec(fr, box);
1614 if (NEED_MUTOT(*inputrec))
1616 /* Calculate total (local) dipole moment in a temporary common array.
1617 * This makes it possible to sum them over nodes faster.
1619 calc_mu(start, homenr,
1620 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1625 if (fr->ePBC != epbcNONE)
1627 /* Compute shift vectors every step,
1628 * because of pressure coupling or box deformation!
1630 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1632 calc_shifts(box, fr->shift_vec);
1637 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1638 &(top->cgs), x, fr->cg_cm);
1639 inc_nrnb(nrnb, eNR_CGCM, homenr);
1640 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1642 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1644 unshift_self(graph, box, x);
1649 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1650 inc_nrnb(nrnb, eNR_CGCM, homenr);
1653 if (bCalcCGCM && gmx_debug_at)
1655 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1659 if (!(cr->duty & DUTY_PME))
1661 /* Send particle coordinates to the pme nodes.
1662 * Since this is only implemented for domain decomposition
1663 * and domain decomposition does not use the graph,
1664 * we do not need to worry about shifting.
1669 wallcycle_start(wcycle, ewcPP_PMESENDX);
1671 bBS = (inputrec->nwall == 2);
1674 copy_mat(box, boxs);
1675 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1678 if (EEL_PME(fr->eeltype))
1680 pme_flags |= GMX_PME_DO_COULOMB;
1683 if (EVDW_PME(fr->vdwtype))
1685 pme_flags |= GMX_PME_DO_LJ;
1688 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1689 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1690 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1693 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1695 #endif /* GMX_MPI */
1697 /* Communicate coordinates and sum dipole if necessary */
1698 if (DOMAINDECOMP(cr))
1700 wallcycle_start(wcycle, ewcMOVEX);
1701 dd_move_x(cr->dd, box, x);
1702 wallcycle_stop(wcycle, ewcMOVEX);
1705 /* update adress weight beforehand */
1706 if (bStateChanged && bDoAdressWF)
1708 /* need pbc for adress weight calculation with pbc_dx */
1709 set_pbc(&pbc, inputrec->ePBC, box);
1710 if (fr->adress_site == eAdressSITEcog)
1712 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1713 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1715 else if (fr->adress_site == eAdressSITEcom)
1717 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1718 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1720 else if (fr->adress_site == eAdressSITEatomatom)
1722 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1723 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1727 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1728 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1732 if (NEED_MUTOT(*inputrec))
1739 gmx_sumd(2*DIM, mu, cr);
1741 for (i = 0; i < 2; i++)
1743 for (j = 0; j < DIM; j++)
1745 fr->mu_tot[i][j] = mu[i*DIM + j];
1749 if (fr->efep == efepNO)
1751 copy_rvec(fr->mu_tot[0], mu_tot);
1755 for (j = 0; j < DIM; j++)
1758 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1763 /* Reset energies */
1764 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1765 clear_rvecs(SHIFTS, fr->fshift);
1769 wallcycle_start(wcycle, ewcNS);
1771 if (graph && bStateChanged)
1773 /* Calculate intramolecular shift vectors to make molecules whole */
1774 mk_mshift(fplog, graph, fr->ePBC, box, x);
1777 /* Do the actual neighbour searching */
1779 groups, top, mdatoms,
1780 cr, nrnb, bFillGrid,
1783 wallcycle_stop(wcycle, ewcNS);
1786 if (inputrec->implicit_solvent && bNS)
1788 make_gb_nblist(cr, inputrec->gb_algorithm,
1789 x, box, fr, &top->idef, graph, fr->born);
1792 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1794 wallcycle_start(wcycle, ewcPPDURINGPME);
1795 dd_force_flop_start(cr->dd, nrnb);
1800 /* Enforced rotation has its own cycle counter that starts after the collective
1801 * coordinates have been communicated. It is added to ddCyclF to allow
1802 * for proper load-balancing */
1803 wallcycle_start(wcycle, ewcROT);
1804 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1805 wallcycle_stop(wcycle, ewcROT);
1808 /* Start the force cycle counter.
1809 * This counter is stopped in do_forcelow_level.
1810 * No parallel communication should occur while this counter is running,
1811 * since that will interfere with the dynamic load balancing.
1813 wallcycle_start(wcycle, ewcFORCE);
1817 /* Reset forces for which the virial is calculated separately:
1818 * PME/Ewald forces if necessary */
1819 if (fr->bF_NoVirSum)
1821 if (flags & GMX_FORCE_VIRIAL)
1823 fr->f_novirsum = fr->f_novirsum_alloc;
1826 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1830 clear_rvecs(homenr, fr->f_novirsum+start);
1835 /* We are not calculating the pressure so we do not need
1836 * a separate array for forces that do not contribute
1843 /* Clear the short- and long-range forces */
1844 clear_rvecs(fr->natoms_force_constr, f);
1845 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1847 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1850 clear_rvec(fr->vir_diag_posres);
1852 if (inputrec->ePull == epullCONSTRAINT)
1854 clear_pull_forces(inputrec->pull);
1857 /* update QMMMrec, if necessary */
1860 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1863 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1865 posres_wrapper(flags, inputrec, nrnb, top, box, x,
1869 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1871 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1874 /* Compute the bonded and non-bonded energies and optionally forces */
1875 do_force_lowlevel(fr, inputrec, &(top->idef),
1876 cr, nrnb, wcycle, mdatoms,
1877 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1878 &(top->atomtypes), bBornRadii, box,
1879 inputrec->fepvals, lambda,
1880 graph, &(top->excls), fr->mu_tot,
1886 if (do_per_step(step, inputrec->nstcalclr))
1888 /* Add the long range forces to the short range forces */
1889 for (i = 0; i < fr->natoms_force_constr; i++)
1891 rvec_add(fr->f_twin[i], f[i], f[i]);
1896 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1900 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1903 if (DOMAINDECOMP(cr))
1905 dd_force_flop_stop(cr->dd, nrnb);
1908 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1914 if (IR_ELEC_FIELD(*inputrec))
1916 /* Compute forces due to electric field */
1917 calc_f_el(MASTER(cr) ? field : NULL,
1918 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1919 inputrec->ex, inputrec->et, t);
1922 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1924 /* Compute thermodynamic force in hybrid AdResS region */
1925 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1926 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1929 /* Communicate the forces */
1930 if (DOMAINDECOMP(cr))
1932 wallcycle_start(wcycle, ewcMOVEF);
1933 dd_move_f(cr->dd, f, fr->fshift);
1934 /* Do we need to communicate the separate force array
1935 * for terms that do not contribute to the single sum virial?
1936 * Position restraints and electric fields do not introduce
1937 * inter-cg forces, only full electrostatics methods do.
1938 * When we do not calculate the virial, fr->f_novirsum = f,
1939 * so we have already communicated these forces.
1941 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1942 (flags & GMX_FORCE_VIRIAL))
1944 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1948 /* We should not update the shift forces here,
1949 * since f_twin is already included in f.
1951 dd_move_f(cr->dd, fr->f_twin, NULL);
1953 wallcycle_stop(wcycle, ewcMOVEF);
1956 /* If we have NoVirSum forces, but we do not calculate the virial,
1957 * we sum fr->f_novirum=f later.
1959 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1961 wallcycle_start(wcycle, ewcVSITESPREAD);
1962 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1963 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1964 wallcycle_stop(wcycle, ewcVSITESPREAD);
1968 wallcycle_start(wcycle, ewcVSITESPREAD);
1969 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1971 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1972 wallcycle_stop(wcycle, ewcVSITESPREAD);
1976 if (flags & GMX_FORCE_VIRIAL)
1978 /* Calculation of the virial must be done after vsites! */
1979 calc_virial(0, mdatoms->homenr, x, f,
1980 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1984 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1986 pull_potential_wrapper(cr, inputrec, box, x,
1987 f, vir_force, mdatoms, enerd, lambda, t);
1990 /* Add the forces from enforced rotation potentials (if any) */
1993 wallcycle_start(wcycle, ewcROTadd);
1994 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1995 wallcycle_stop(wcycle, ewcROTadd);
1998 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1999 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
2001 if (PAR(cr) && !(cr->duty & DUTY_PME))
2003 /* In case of node-splitting, the PP nodes receive the long-range
2004 * forces, virial and energy from the PME nodes here.
2006 pme_receive_force_ener(cr, wcycle, enerd, fr);
2011 post_process_forces(cr, step, nrnb, wcycle,
2012 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
2016 /* Sum the potential energy terms from group contributions */
2017 sum_epot(&(enerd->grpp), enerd->term);
2020 void do_force(FILE *fplog, t_commrec *cr,
2021 t_inputrec *inputrec,
2022 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2023 gmx_localtop_t *top,
2024 gmx_groups_t *groups,
2025 matrix box, rvec x[], history_t *hist,
2029 gmx_enerdata_t *enerd, t_fcdata *fcd,
2030 real *lambda, t_graph *graph,
2032 gmx_vsite_t *vsite, rvec mu_tot,
2033 double t, FILE *field, gmx_edsam_t ed,
2034 gmx_bool bBornRadii,
2037 /* modify force flag if not doing nonbonded */
2038 if (!fr->bNonbonded)
2040 flags &= ~GMX_FORCE_NONBONDED;
2043 switch (inputrec->cutoff_scheme)
2046 do_force_cutsVERLET(fplog, cr, inputrec,
2062 do_force_cutsGROUP(fplog, cr, inputrec,
2077 gmx_incons("Invalid cut-off scheme passed!");
2082 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2083 t_inputrec *ir, t_mdatoms *md,
2084 t_state *state, t_commrec *cr, t_nrnb *nrnb,
2085 t_forcerec *fr, gmx_localtop_t *top)
2087 int i, m, start, end;
2089 real dt = ir->delta_t;
2093 snew(savex, state->natoms);
2100 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2101 start, md->homenr, end);
2103 /* Do a first constrain to reset particles... */
2104 step = ir->init_step;
2107 char buf[STEPSTRSIZE];
2108 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2109 gmx_step_str(step, buf));
2113 /* constrain the current position */
2114 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2115 ir, NULL, cr, step, 0, 1.0, md,
2116 state->x, state->x, NULL,
2117 fr->bMolPBC, state->box,
2118 state->lambda[efptBONDED], &dvdl_dum,
2119 NULL, NULL, nrnb, econqCoord,
2120 ir->epc == epcMTTK, state->veta, state->veta);
2123 /* constrain the inital velocity, and save it */
2124 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2125 /* might not yet treat veta correctly */
2126 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2127 ir, NULL, cr, step, 0, 1.0, md,
2128 state->x, state->v, state->v,
2129 fr->bMolPBC, state->box,
2130 state->lambda[efptBONDED], &dvdl_dum,
2131 NULL, NULL, nrnb, econqVeloc,
2132 ir->epc == epcMTTK, state->veta, state->veta);
2134 /* constrain the inital velocities at t-dt/2 */
2135 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2137 for (i = start; (i < end); i++)
2139 for (m = 0; (m < DIM); m++)
2141 /* Reverse the velocity */
2142 state->v[i][m] = -state->v[i][m];
2143 /* Store the position at t-dt in buf */
2144 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2147 /* Shake the positions at t=-dt with the positions at t=0
2148 * as reference coordinates.
2152 char buf[STEPSTRSIZE];
2153 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2154 gmx_step_str(step, buf));
2157 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2158 ir, NULL, cr, step, -1, 1.0, md,
2159 state->x, savex, NULL,
2160 fr->bMolPBC, state->box,
2161 state->lambda[efptBONDED], &dvdl_dum,
2162 state->v, NULL, nrnb, econqCoord,
2163 ir->epc == epcMTTK, state->veta, state->veta);
2165 for (i = start; i < end; i++)
2167 for (m = 0; m < DIM; m++)
2169 /* Re-reverse the velocities */
2170 state->v[i][m] = -state->v[i][m];
2179 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2180 double *enerout, double *virout)
2182 double enersum, virsum;
2183 double invscale, invscale2, invscale3;
2184 double r, ea, eb, ec, pa, pb, pc, pd;
2186 int ri, offset, tabfactor;
2188 invscale = 1.0/scale;
2189 invscale2 = invscale*invscale;
2190 invscale3 = invscale*invscale2;
2192 /* Following summation derived from cubic spline definition,
2193 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2194 * the cubic spline. We first calculate the negative of the
2195 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2196 * add the more standard, abrupt cutoff correction to that result,
2197 * yielding the long-range correction for a switched function. We
2198 * perform both the pressure and energy loops at the same time for
2199 * simplicity, as the computational cost is low. */
2203 /* Since the dispersion table has been scaled down a factor
2204 * 6.0 and the repulsion a factor 12.0 to compensate for the
2205 * c6/c12 parameters inside nbfp[] being scaled up (to save
2206 * flops in kernels), we need to correct for this.
2217 for (ri = rstart; ri < rend; ++ri)
2221 eb = 2.0*invscale2*r;
2225 pb = 3.0*invscale2*r;
2226 pc = 3.0*invscale*r*r;
2229 /* this "8" is from the packing in the vdwtab array - perhaps
2230 should be defined? */
2232 offset = 8*ri + offstart;
2233 y0 = vdwtab[offset];
2234 f = vdwtab[offset+1];
2235 g = vdwtab[offset+2];
2236 h = vdwtab[offset+3];
2238 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);
2239 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);
2241 *enerout = 4.0*M_PI*enersum*tabfactor;
2242 *virout = 4.0*M_PI*virsum*tabfactor;
2245 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2247 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2248 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2249 double invscale, invscale2, invscale3;
2250 int ri0, ri1, ri, i, offstart, offset;
2251 real scale, *vdwtab, tabfactor, tmp;
2253 fr->enershiftsix = 0;
2254 fr->enershifttwelve = 0;
2255 fr->enerdiffsix = 0;
2256 fr->enerdifftwelve = 0;
2258 fr->virdifftwelve = 0;
2260 if (eDispCorr != edispcNO)
2262 for (i = 0; i < 2; i++)
2267 if ((fr->vdw_modifier == eintmodPOTSHIFT) ||
2268 (fr->vdw_modifier == eintmodPOTSWITCH) ||
2269 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2270 (fr->vdwtype == evdwSHIFT) ||
2271 (fr->vdwtype == evdwSWITCH))
2273 if (((fr->vdw_modifier == eintmodPOTSWITCH) ||
2274 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2275 (fr->vdwtype == evdwSWITCH)) && fr->rvdw_switch == 0)
2278 "With dispersion correction rvdw-switch can not be zero "
2279 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2282 scale = fr->nblists[0].table_vdw.scale;
2283 vdwtab = fr->nblists[0].table_vdw.data;
2285 /* Round the cut-offs to exact table values for precision */
2286 ri0 = floor(fr->rvdw_switch*scale);
2287 ri1 = ceil(fr->rvdw*scale);
2289 /* The code below has some support for handling force-switching, i.e.
2290 * when the force (instead of potential) is switched over a limited
2291 * region. This leads to a constant shift in the potential inside the
2292 * switching region, which we can handle by adding a constant energy
2293 * term in the force-switch case just like when we do potential-shift.
2295 * For now this is not enabled, but to keep the functionality in the
2296 * code we check separately for switch and shift. When we do force-switch
2297 * the shifting point is rvdw_switch, while it is the cutoff when we
2298 * have a classical potential-shift.
2300 * For a pure potential-shift the potential has a constant shift
2301 * all the way out to the cutoff, and that is it. For other forms
2302 * we need to calculate the constant shift up to the point where we
2303 * start modifying the potential.
2305 ri0 = (fr->vdw_modifier == eintmodPOTSHIFT) ? ri1 : ri0;
2312 if ((fr->vdw_modifier == eintmodFORCESWITCH) ||
2313 (fr->vdwtype == evdwSHIFT))
2315 /* Determine the constant energy shift below rvdw_switch.
2316 * Table has a scale factor since we have scaled it down to compensate
2317 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2319 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2320 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2322 else if (fr->vdw_modifier == eintmodPOTSHIFT)
2324 fr->enershiftsix = (real)(-1.0/(rc3*rc3));
2325 fr->enershifttwelve = (real)( 1.0/(rc9*rc3));
2328 /* Add the constant part from 0 to rvdw_switch.
2329 * This integration from 0 to rvdw_switch overcounts the number
2330 * of interactions by 1, as it also counts the self interaction.
2331 * We will correct for this later.
2333 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2334 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2336 /* Calculate the contribution in the range [r0,r1] where we
2337 * modify the potential. For a pure potential-shift modifier we will
2338 * have ri0==ri1, and there will not be any contribution here.
2340 for (i = 0; i < 2; i++)
2344 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2345 eners[i] -= enersum;
2349 /* Alright: Above we compensated by REMOVING the parts outside r0
2350 * corresponding to the ideal VdW 1/r6 and /r12 potentials.
2352 * Regardless of whether r0 is the point where we start switching,
2353 * or the cutoff where we calculated the constant shift, we include
2354 * all the parts we are missing out to infinity from r0 by
2355 * calculating the analytical dispersion correction.
2357 eners[0] += -4.0*M_PI/(3.0*rc3);
2358 eners[1] += 4.0*M_PI/(9.0*rc9);
2359 virs[0] += 8.0*M_PI/rc3;
2360 virs[1] += -16.0*M_PI/(3.0*rc9);
2362 else if (fr->vdwtype == evdwCUT ||
2363 EVDW_PME(fr->vdwtype) ||
2364 fr->vdwtype == evdwUSER)
2366 if (fr->vdwtype == evdwUSER && fplog)
2369 "WARNING: using dispersion correction with user tables\n");
2372 /* Note that with LJ-PME, the dispersion correction is multiplied
2373 * by the difference between the actual C6 and the value of C6
2374 * that would produce the combination rule.
2375 * This means the normal energy and virial difference formulas
2379 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2381 /* Contribution beyond the cut-off */
2382 eners[0] += -4.0*M_PI/(3.0*rc3);
2383 eners[1] += 4.0*M_PI/(9.0*rc9);
2384 if (fr->vdw_modifier == eintmodPOTSHIFT)
2386 /* Contribution within the cut-off */
2387 eners[0] += -4.0*M_PI/(3.0*rc3);
2388 eners[1] += 4.0*M_PI/(3.0*rc9);
2390 /* Contribution beyond the cut-off */
2391 virs[0] += 8.0*M_PI/rc3;
2392 virs[1] += -16.0*M_PI/(3.0*rc9);
2397 "Dispersion correction is not implemented for vdw-type = %s",
2398 evdw_names[fr->vdwtype]);
2401 /* When we deprecate the group kernels the code below can go too */
2402 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2404 /* Calculate self-interaction coefficient (assuming that
2405 * the reciprocal-space contribution is constant in the
2406 * region that contributes to the self-interaction).
2408 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2410 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2411 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2414 fr->enerdiffsix = eners[0];
2415 fr->enerdifftwelve = eners[1];
2416 /* The 0.5 is due to the Gromacs definition of the virial */
2417 fr->virdiffsix = 0.5*virs[0];
2418 fr->virdifftwelve = 0.5*virs[1];
2422 void calc_dispcorr(t_inputrec *ir, t_forcerec *fr,
2424 matrix box, real lambda, tensor pres, tensor virial,
2425 real *prescorr, real *enercorr, real *dvdlcorr)
2427 gmx_bool bCorrAll, bCorrPres;
2428 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2438 if (ir->eDispCorr != edispcNO)
2440 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2441 ir->eDispCorr == edispcAllEnerPres);
2442 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2443 ir->eDispCorr == edispcAllEnerPres);
2445 invvol = 1/det(box);
2448 /* Only correct for the interactions with the inserted molecule */
2449 dens = (natoms - fr->n_tpi)*invvol;
2454 dens = natoms*invvol;
2455 ninter = 0.5*natoms;
2458 if (ir->efep == efepNO)
2460 avcsix = fr->avcsix[0];
2461 avctwelve = fr->avctwelve[0];
2465 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2466 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2469 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2470 *enercorr += avcsix*enerdiff;
2472 if (ir->efep != efepNO)
2474 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2478 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2479 *enercorr += avctwelve*enerdiff;
2480 if (fr->efep != efepNO)
2482 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2488 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2489 if (ir->eDispCorr == edispcAllEnerPres)
2491 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2493 /* The factor 2 is because of the Gromacs virial definition */
2494 spres = -2.0*invvol*svir*PRESFAC;
2496 for (m = 0; m < DIM; m++)
2498 virial[m][m] += svir;
2499 pres[m][m] += spres;
2504 /* Can't currently control when it prints, for now, just print when degugging */
2509 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2515 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2516 *enercorr, spres, svir);
2520 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2524 if (fr->efep != efepNO)
2526 *dvdlcorr += dvdlambda;
2531 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2532 t_graph *graph, rvec x[])
2536 fprintf(fplog, "Removing pbc first time\n");
2538 calc_shifts(box, fr->shift_vec);
2541 mk_mshift(fplog, graph, fr->ePBC, box, x);
2544 p_graph(debug, "do_pbc_first 1", graph);
2546 shift_self(graph, box, x);
2547 /* By doing an extra mk_mshift the molecules that are broken
2548 * because they were e.g. imported from another software
2549 * will be made whole again. Such are the healing powers
2552 mk_mshift(fplog, graph, fr->ePBC, box, x);
2555 p_graph(debug, "do_pbc_first 2", graph);
2560 fprintf(fplog, "Done rmpbc\n");
2564 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2565 gmx_mtop_t *mtop, rvec x[],
2570 gmx_molblock_t *molb;
2572 if (bFirst && fplog)
2574 fprintf(fplog, "Removing pbc first time\n");
2579 for (mb = 0; mb < mtop->nmolblock; mb++)
2581 molb = &mtop->molblock[mb];
2582 if (molb->natoms_mol == 1 ||
2583 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2585 /* Just one atom or charge group in the molecule, no PBC required */
2586 as += molb->nmol*molb->natoms_mol;
2590 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2591 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2592 0, molb->natoms_mol, FALSE, FALSE, graph);
2594 for (mol = 0; mol < molb->nmol; mol++)
2596 mk_mshift(fplog, graph, ePBC, box, x+as);
2598 shift_self(graph, box, x+as);
2599 /* The molecule is whole now.
2600 * We don't need the second mk_mshift call as in do_pbc_first,
2601 * since we no longer need this graph.
2604 as += molb->natoms_mol;
2612 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2613 gmx_mtop_t *mtop, rvec x[])
2615 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2618 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2619 gmx_mtop_t *mtop, rvec x[])
2621 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2624 void finish_run(FILE *fplog, t_commrec *cr,
2625 t_inputrec *inputrec,
2626 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2627 gmx_walltime_accounting_t walltime_accounting,
2628 nonbonded_verlet_t *nbv,
2629 gmx_bool bWriteStat)
2632 t_nrnb *nrnb_tot = NULL;
2635 double elapsed_time,
2636 elapsed_time_over_all_ranks,
2637 elapsed_time_over_all_threads,
2638 elapsed_time_over_all_threads_over_all_ranks;
2639 wallcycle_sum(cr, wcycle);
2645 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2646 cr->mpi_comm_mysim);
2654 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2655 elapsed_time_over_all_ranks = elapsed_time;
2656 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2657 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2661 /* reduce elapsed_time over all MPI ranks in the current simulation */
2662 MPI_Allreduce(&elapsed_time,
2663 &elapsed_time_over_all_ranks,
2664 1, MPI_DOUBLE, MPI_SUM,
2665 cr->mpi_comm_mysim);
2666 elapsed_time_over_all_ranks /= cr->nnodes;
2667 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2668 * current simulation. */
2669 MPI_Allreduce(&elapsed_time_over_all_threads,
2670 &elapsed_time_over_all_threads_over_all_ranks,
2671 1, MPI_DOUBLE, MPI_SUM,
2672 cr->mpi_comm_mysim);
2678 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2685 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2687 print_dd_statistics(cr, inputrec, fplog);
2692 wallclock_gpu_t* gputimes = use_GPU(nbv) ?
2693 nbnxn_cuda_get_timings(nbv->cu_nbv) : NULL;
2694 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2695 elapsed_time_over_all_ranks,
2698 if (EI_DYNAMICS(inputrec->eI))
2700 delta_t = inputrec->delta_t;
2709 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2710 elapsed_time_over_all_ranks,
2711 walltime_accounting_get_nsteps_done(walltime_accounting),
2712 delta_t, nbfs, mflop);
2716 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2717 elapsed_time_over_all_ranks,
2718 walltime_accounting_get_nsteps_done(walltime_accounting),
2719 delta_t, nbfs, mflop);
2724 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2726 /* this function works, but could probably use a logic rewrite to keep all the different
2727 types of efep straight. */
2730 t_lambda *fep = ir->fepvals;
2732 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2734 for (i = 0; i < efptNR; i++)
2746 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2747 if checkpoint is set -- a kludge is in for now
2749 for (i = 0; i < efptNR; i++)
2751 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2752 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2754 lambda[i] = fep->init_lambda;
2757 lam0[i] = lambda[i];
2762 lambda[i] = fep->all_lambda[i][*fep_state];
2765 lam0[i] = lambda[i];
2771 /* need to rescale control temperatures to match current state */
2772 for (i = 0; i < ir->opts.ngtc; i++)
2774 if (ir->opts.ref_t[i] > 0)
2776 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2782 /* Send to the log the information on the current lambdas */
2785 fprintf(fplog, "Initial vector of lambda components:[ ");
2786 for (i = 0; i < efptNR; i++)
2788 fprintf(fplog, "%10.4f ", lambda[i]);
2790 fprintf(fplog, "]\n");
2796 void init_md(FILE *fplog,
2797 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2798 double *t, double *t0,
2799 real *lambda, int *fep_state, double *lam0,
2800 t_nrnb *nrnb, gmx_mtop_t *mtop,
2802 int nfile, const t_filenm fnm[],
2803 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2804 tensor force_vir, tensor shake_vir, rvec mu_tot,
2805 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags)
2810 /* Initial values */
2811 *t = *t0 = ir->init_t;
2814 for (i = 0; i < ir->opts.ngtc; i++)
2816 /* set bSimAnn if any group is being annealed */
2817 if (ir->opts.annealing[i] != eannNO)
2824 update_annealing_target_temp(&(ir->opts), ir->init_t);
2827 /* Initialize lambda variables */
2828 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2832 *upd = init_update(ir);
2838 *vcm = init_vcm(fplog, &mtop->groups, ir);
2841 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2843 if (ir->etc == etcBERENDSEN)
2845 please_cite(fplog, "Berendsen84a");
2847 if (ir->etc == etcVRESCALE)
2849 please_cite(fplog, "Bussi2007a");
2851 if (ir->eI == eiSD1)
2853 please_cite(fplog, "Goga2012");
2861 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv);
2863 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2864 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2869 please_cite(fplog, "Fritsch12");
2870 please_cite(fplog, "Junghans10");
2872 /* Initiate variables */
2873 clear_mat(force_vir);
2874 clear_mat(shake_vir);