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.
45 #ifdef HAVE_SYS_TIME_H
49 #include "gromacs/legacyheaders/typedefs.h"
50 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/legacyheaders/txtdump.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/legacyheaders/chargegroup.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/legacyheaders/nrnb.h"
57 #include "gromacs/legacyheaders/mdrun.h"
58 #include "gromacs/legacyheaders/sim_util.h"
59 #include "gromacs/legacyheaders/update.h"
60 #include "gromacs/math/units.h"
61 #include "gromacs/legacyheaders/mdatoms.h"
62 #include "gromacs/legacyheaders/force.h"
63 #include "gromacs/legacyheaders/bondf.h"
64 #include "gromacs/legacyheaders/pme.h"
65 #include "gromacs/legacyheaders/disre.h"
66 #include "gromacs/legacyheaders/orires.h"
67 #include "gromacs/legacyheaders/network.h"
68 #include "gromacs/legacyheaders/calcmu.h"
69 #include "gromacs/legacyheaders/constr.h"
70 #include "gromacs/legacyheaders/copyrite.h"
71 #include "gromacs/legacyheaders/domdec.h"
72 #include "gromacs/legacyheaders/genborn.h"
73 #include "nbnxn_atomdata.h"
74 #include "nbnxn_search.h"
75 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
76 #include "nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
77 #include "nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
78 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
79 #include "gromacs/legacyheaders/nonbonded.h"
80 #include "../gmxlib/nonbonded/nb_kernel.h"
81 #include "../gmxlib/nonbonded/nb_free_energy.h"
83 #include "gromacs/legacyheaders/types/commrec.h"
84 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_data_mgmt.h"
85 #include "gromacs/pbcutil/ishift.h"
86 #include "gromacs/pbcutil/mshift.h"
87 #include "gromacs/timing/wallcycle.h"
88 #include "gromacs/timing/walltime_accounting.h"
89 #include "gromacs/utility/gmxmpi.h"
90 #include "gromacs/utility/smalloc.h"
91 #include "gromacs/essentialdynamics/edsam.h"
92 #include "gromacs/pulling/pull.h"
93 #include "gromacs/pulling/pull_rotation.h"
94 #include "gromacs/imd/imd.h"
96 #include "gromacs/legacyheaders/qmmm.h"
98 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
100 #include "nbnxn_cuda/nbnxn_cuda.h"
102 #include "nb_verlet.h"
104 void print_time(FILE *out,
105 gmx_walltime_accounting_t walltime_accounting,
108 t_commrec gmx_unused *cr)
111 char timebuf[STRLEN];
112 double dt, elapsed_seconds, time_per_step;
115 #ifndef GMX_THREAD_MPI
121 fprintf(out, "step %s", gmx_step_str(step, buf));
122 if ((step >= ir->nstlist))
124 double seconds_since_epoch = gmx_gettime();
125 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
126 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
127 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
133 finish = (time_t) (seconds_since_epoch + dt);
134 gmx_ctime_r(&finish, timebuf, STRLEN);
135 sprintf(buf, "%s", timebuf);
136 buf[strlen(buf)-1] = '\0';
137 fprintf(out, ", will finish %s", buf);
141 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
146 fprintf(out, " performance: %.1f ns/day ",
147 ir->delta_t/1000*24*60*60/time_per_step);
150 #ifndef GMX_THREAD_MPI
160 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
163 char time_string[STRLEN];
172 char timebuf[STRLEN];
173 time_t temp_time = (time_t) the_time;
175 gmx_ctime_r(&temp_time, timebuf, STRLEN);
176 for (i = 0; timebuf[i] >= ' '; i++)
178 time_string[i] = timebuf[i];
180 time_string[i] = '\0';
183 fprintf(fplog, "%s on rank %d %s\n", title, nodeid, time_string);
186 void print_start(FILE *fplog, t_commrec *cr,
187 gmx_walltime_accounting_t walltime_accounting,
192 sprintf(buf, "Started %s", name);
193 print_date_and_time(fplog, cr->nodeid, buf,
194 walltime_accounting_get_start_time_stamp(walltime_accounting));
197 static void sum_forces(int start, int end, rvec f[], rvec flr[])
203 pr_rvecs(debug, 0, "fsr", f+start, end-start);
204 pr_rvecs(debug, 0, "flr", flr+start, end-start);
206 for (i = start; (i < end); i++)
208 rvec_inc(f[i], flr[i]);
213 * calc_f_el calculates forces due to an electric field.
215 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
217 * Et[] contains the parameters for the time dependent
218 * part of the field (not yet used).
219 * Ex[] contains the parameters for
220 * the spatial dependent part of the field. You can have cool periodic
221 * fields in principle, but only a constant field is supported
223 * The function should return the energy due to the electric field
224 * (if any) but for now returns 0.
227 * There can be problems with the virial.
228 * Since the field is not self-consistent this is unavoidable.
229 * For neutral molecules the virial is correct within this approximation.
230 * For neutral systems with many charged molecules the error is small.
231 * But for systems with a net charge or a few charged molecules
232 * the error can be significant when the field is high.
233 * Solution: implement a self-consitent electric field into PME.
235 static void calc_f_el(FILE *fp, int start, int homenr,
236 real charge[], rvec f[],
237 t_cosines Ex[], t_cosines Et[], double t)
243 for (m = 0; (m < DIM); m++)
250 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
254 Ext[m] = cos(Et[m].a[0]*t);
263 /* Convert the field strength from V/nm to MD-units */
264 Ext[m] *= Ex[m].a[0]*FIELDFAC;
265 for (i = start; (i < start+homenr); i++)
267 f[i][m] += charge[i]*Ext[m];
277 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
278 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
282 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
283 tensor vir_part, t_graph *graph, matrix box,
284 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
289 /* The short-range virial from surrounding boxes */
291 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
292 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
294 /* Calculate partial virial, for local atoms only, based on short range.
295 * Total virial is computed in global_stat, called from do_md
297 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
298 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
300 /* Add position restraint contribution */
301 for (i = 0; i < DIM; i++)
303 vir_part[i][i] += fr->vir_diag_posres[i];
306 /* Add wall contribution */
307 for (i = 0; i < DIM; i++)
309 vir_part[i][ZZ] += fr->vir_wall_z[i];
314 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
318 static void posres_wrapper(int flags,
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);
340 enerd->term[F_POSRES] += v;
341 /* If just the force constant changes, the FEP term is linear,
342 * but if k changes, it is not.
344 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
345 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
347 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
349 for (i = 0; i < enerd->n_lambda; i++)
351 real dvdl_dum, lambda_dum;
353 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
354 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
355 top->idef.iparams_posres,
356 (const rvec*)x, NULL, NULL,
357 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
358 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
359 enerd->enerpart_lambda[i] += v;
364 static void fbposres_wrapper(t_inputrec *ir,
367 matrix box, rvec x[],
368 gmx_enerdata_t *enerd,
374 /* Flat-bottomed position restraints always require full pbc */
375 set_pbc(&pbc, ir->ePBC, box);
376 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
377 top->idef.iparams_fbposres,
378 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
379 ir->ePBC == epbcNONE ? NULL : &pbc,
380 fr->rc_scaling, fr->ePBC, fr->posres_com);
381 enerd->term[F_FBPOSRES] += v;
382 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
385 static void pull_potential_wrapper(t_commrec *cr,
387 matrix box, rvec x[],
391 gmx_enerdata_t *enerd,
398 /* Calculate the center of mass forces, this requires communication,
399 * which is why pull_potential is called close to other communication.
400 * The virial contribution is calculated directly,
401 * which is why we call pull_potential after calc_virial.
403 set_pbc(&pbc, ir->ePBC, box);
405 enerd->term[F_COM_PULL] +=
406 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
407 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
408 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
411 static void pme_receive_force_ener(t_commrec *cr,
412 gmx_wallcycle_t wcycle,
413 gmx_enerdata_t *enerd,
416 real e_q, e_lj, v, dvdl_q, dvdl_lj;
417 float cycles_ppdpme, cycles_seppme;
419 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
420 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
422 /* In case of node-splitting, the PP nodes receive the long-range
423 * forces, virial and energy from the PME nodes here.
425 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
428 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
429 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
431 enerd->term[F_COUL_RECIP] += e_q;
432 enerd->term[F_LJ_RECIP] += e_lj;
433 enerd->dvdl_lin[efptCOUL] += dvdl_q;
434 enerd->dvdl_lin[efptVDW] += dvdl_lj;
438 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
440 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
443 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
444 gmx_int64_t step, real pforce, rvec *x, rvec *f)
448 char buf[STEPSTRSIZE];
451 for (i = 0; i < md->homenr; i++)
454 /* We also catch NAN, if the compiler does not optimize this away. */
455 if (fn2 >= pf2 || fn2 != fn2)
457 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
458 gmx_step_str(step, buf),
459 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
464 static void post_process_forces(t_commrec *cr,
466 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
468 matrix box, rvec x[],
473 t_forcerec *fr, gmx_vsite_t *vsite,
480 /* Spread the mesh force on virtual sites to the other particles...
481 * This is parallellized. MPI communication is performed
482 * if the constructing atoms aren't local.
484 wallcycle_start(wcycle, ewcVSITESPREAD);
485 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
486 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
488 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
489 wallcycle_stop(wcycle, ewcVSITESPREAD);
491 if (flags & GMX_FORCE_VIRIAL)
493 /* Now add the forces, this is local */
496 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
500 sum_forces(0, mdatoms->homenr,
503 if (EEL_FULL(fr->eeltype))
505 /* Add the mesh contribution to the virial */
506 m_add(vir_force, fr->vir_el_recip, vir_force);
508 if (EVDW_PME(fr->vdwtype))
510 /* Add the mesh contribution to the virial */
511 m_add(vir_force, fr->vir_lj_recip, vir_force);
515 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
520 if (fr->print_force >= 0)
522 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
526 static void do_nb_verlet(t_forcerec *fr,
527 interaction_const_t *ic,
528 gmx_enerdata_t *enerd,
529 int flags, int ilocality,
532 gmx_wallcycle_t wcycle)
534 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
536 nonbonded_verlet_group_t *nbvg;
539 if (!(flags & GMX_FORCE_NONBONDED))
541 /* skip non-bonded calculation */
545 nbvg = &fr->nbv->grp[ilocality];
547 /* CUDA kernel launch overhead is already timed separately */
548 if (fr->cutoff_scheme != ecutsVERLET)
550 gmx_incons("Invalid cut-off scheme passed!");
553 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
557 wallcycle_sub_start(wcycle, ewcsNONBONDED);
559 switch (nbvg->kernel_type)
561 case nbnxnk4x4_PlainC:
562 nbnxn_kernel_ref(&nbvg->nbl_lists,
568 enerd->grpp.ener[egCOULSR],
570 enerd->grpp.ener[egBHAMSR] :
571 enerd->grpp.ener[egLJSR]);
574 case nbnxnk4xN_SIMD_4xN:
575 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
582 enerd->grpp.ener[egCOULSR],
584 enerd->grpp.ener[egBHAMSR] :
585 enerd->grpp.ener[egLJSR]);
587 case nbnxnk4xN_SIMD_2xNN:
588 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
595 enerd->grpp.ener[egCOULSR],
597 enerd->grpp.ener[egBHAMSR] :
598 enerd->grpp.ener[egLJSR]);
601 case nbnxnk8x8x8_CUDA:
602 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
605 case nbnxnk8x8x8_PlainC:
606 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
611 nbvg->nbat->out[0].f,
613 enerd->grpp.ener[egCOULSR],
615 enerd->grpp.ener[egBHAMSR] :
616 enerd->grpp.ener[egLJSR]);
620 gmx_incons("Invalid nonbonded kernel type passed!");
625 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
628 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
630 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
632 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
633 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
635 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
639 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
641 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
642 if (flags & GMX_FORCE_ENERGY)
644 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
645 enr_nbnxn_kernel_ljc += 1;
646 enr_nbnxn_kernel_lj += 1;
649 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
650 nbvg->nbl_lists.natpair_ljq);
651 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
652 nbvg->nbl_lists.natpair_lj);
653 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
654 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
655 nbvg->nbl_lists.natpair_q);
657 if (ic->vdw_modifier == eintmodFORCESWITCH)
659 /* We add up the switch cost separately */
660 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
661 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
663 if (ic->vdw_modifier == eintmodPOTSWITCH)
665 /* We add up the switch cost separately */
666 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
667 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
669 if (ic->vdwtype == evdwPME)
671 /* We add up the LJ Ewald cost separately */
672 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
673 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
677 static void do_nb_verlet_fep(nbnxn_pairlist_set_t *nbl_lists,
684 gmx_enerdata_t *enerd,
687 gmx_wallcycle_t wcycle)
690 nb_kernel_data_t kernel_data;
692 real dvdl_nb[efptNR];
697 /* Add short-range interactions */
698 donb_flags |= GMX_NONBONDED_DO_SR;
700 /* Currently all group scheme kernels always calculate (shift-)forces */
701 if (flags & GMX_FORCE_FORCES)
703 donb_flags |= GMX_NONBONDED_DO_FORCE;
705 if (flags & GMX_FORCE_VIRIAL)
707 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
709 if (flags & GMX_FORCE_ENERGY)
711 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
713 if (flags & GMX_FORCE_DO_LR)
715 donb_flags |= GMX_NONBONDED_DO_LR;
718 kernel_data.flags = donb_flags;
719 kernel_data.lambda = lambda;
720 kernel_data.dvdl = dvdl_nb;
722 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
723 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
725 /* reset free energy components */
726 for (i = 0; i < efptNR; i++)
731 assert(gmx_omp_nthreads_get(emntNonbonded) == nbl_lists->nnbl);
733 wallcycle_sub_start(wcycle, ewcsNONBONDED);
734 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
735 for (th = 0; th < nbl_lists->nnbl; th++)
737 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
738 x, f, fr, mdatoms, &kernel_data, nrnb);
741 if (fepvals->sc_alpha != 0)
743 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
744 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
748 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
749 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
752 /* If we do foreign lambda and we have soft-core interactions
753 * we have to recalculate the (non-linear) energies contributions.
755 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
757 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
758 kernel_data.lambda = lam_i;
759 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
760 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
761 /* Note that we add to kernel_data.dvdl, but ignore the result */
763 for (i = 0; i < enerd->n_lambda; i++)
765 for (j = 0; j < efptNR; j++)
767 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
769 reset_foreign_enerdata(enerd);
770 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
771 for (th = 0; th < nbl_lists->nnbl; th++)
773 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
774 x, f, fr, mdatoms, &kernel_data, nrnb);
777 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
778 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
782 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
785 gmx_bool use_GPU(const nonbonded_verlet_t *nbv)
787 return nbv != NULL && nbv->bUseGPU;
790 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
791 t_inputrec *inputrec,
792 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
794 gmx_groups_t gmx_unused *groups,
795 matrix box, rvec x[], history_t *hist,
799 gmx_enerdata_t *enerd, t_fcdata *fcd,
800 real *lambda, t_graph *graph,
801 t_forcerec *fr, interaction_const_t *ic,
802 gmx_vsite_t *vsite, rvec mu_tot,
803 double t, FILE *field, gmx_edsam_t ed,
811 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
812 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
813 gmx_bool bDiffKernels = FALSE;
815 rvec vzero, box_diag;
817 float cycles_pme, cycles_force, cycles_wait_gpu;
818 nonbonded_verlet_t *nbv;
823 nb_kernel_type = fr->nbv->grp[0].kernel_type;
826 homenr = mdatoms->homenr;
828 clear_mat(vir_force);
831 if (DOMAINDECOMP(cr))
833 cg1 = cr->dd->ncg_tot;
844 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
845 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
846 bFillGrid = (bNS && bStateChanged);
847 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
848 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
849 bDoForces = (flags & GMX_FORCE_FORCES);
850 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
851 bUseGPU = fr->nbv->bUseGPU;
852 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
856 update_forcerec(fr, box);
858 if (NEED_MUTOT(*inputrec))
860 /* Calculate total (local) dipole moment in a temporary common array.
861 * This makes it possible to sum them over nodes faster.
863 calc_mu(start, homenr,
864 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
869 if (fr->ePBC != epbcNONE)
871 /* Compute shift vectors every step,
872 * because of pressure coupling or box deformation!
874 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
876 calc_shifts(box, fr->shift_vec);
881 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
882 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
884 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
886 unshift_self(graph, box, x);
890 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
891 fr->shift_vec, nbv->grp[0].nbat);
894 if (!(cr->duty & DUTY_PME))
896 /* Send particle coordinates to the pme nodes.
897 * Since this is only implemented for domain decomposition
898 * and domain decomposition does not use the graph,
899 * we do not need to worry about shifting.
904 wallcycle_start(wcycle, ewcPP_PMESENDX);
906 bBS = (inputrec->nwall == 2);
910 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
913 if (EEL_PME(fr->eeltype))
915 pme_flags |= GMX_PME_DO_COULOMB;
918 if (EVDW_PME(fr->vdwtype))
920 pme_flags |= GMX_PME_DO_LJ;
923 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
924 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
925 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
928 wallcycle_stop(wcycle, ewcPP_PMESENDX);
932 /* do gridding for pair search */
935 if (graph && bStateChanged)
937 /* Calculate intramolecular shift vectors to make molecules whole */
938 mk_mshift(fplog, graph, fr->ePBC, box, x);
942 box_diag[XX] = box[XX][XX];
943 box_diag[YY] = box[YY][YY];
944 box_diag[ZZ] = box[ZZ][ZZ];
946 wallcycle_start(wcycle, ewcNS);
949 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
950 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
952 0, mdatoms->homenr, -1, fr->cginfo, x,
954 nbv->grp[eintLocal].kernel_type,
955 nbv->grp[eintLocal].nbat);
956 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
960 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
961 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
963 nbv->grp[eintNonlocal].kernel_type,
964 nbv->grp[eintNonlocal].nbat);
965 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
968 if (nbv->ngrp == 1 ||
969 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
971 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
972 nbv->nbs, mdatoms, fr->cginfo);
976 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
977 nbv->nbs, mdatoms, fr->cginfo);
978 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
979 nbv->nbs, mdatoms, fr->cginfo);
981 wallcycle_stop(wcycle, ewcNS);
984 /* initialize the GPU atom data and copy shift vector */
989 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
990 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
991 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
994 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
995 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
996 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
999 /* do local pair search */
1002 wallcycle_start_nocount(wcycle, ewcNS);
1003 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1004 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
1007 nbv->min_ci_balanced,
1008 &nbv->grp[eintLocal].nbl_lists,
1010 nbv->grp[eintLocal].kernel_type,
1012 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1016 /* initialize local pair-list on the GPU */
1017 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1018 nbv->grp[eintLocal].nbl_lists.nbl[0],
1021 wallcycle_stop(wcycle, ewcNS);
1025 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1026 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1027 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
1028 nbv->grp[eintLocal].nbat);
1029 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1030 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1035 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1036 /* launch local nonbonded F on GPU */
1037 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
1039 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1042 /* Communicate coordinates and sum dipole if necessary +
1043 do non-local pair search */
1044 if (DOMAINDECOMP(cr))
1046 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
1047 nbv->grp[eintLocal].kernel_type);
1051 /* With GPU+CPU non-bonded calculations we need to copy
1052 * the local coordinates to the non-local nbat struct
1053 * (in CPU format) as the non-local kernel call also
1054 * calculates the local - non-local interactions.
1056 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1057 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1058 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
1059 nbv->grp[eintNonlocal].nbat);
1060 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1061 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1066 wallcycle_start_nocount(wcycle, ewcNS);
1067 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1071 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
1074 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
1077 nbv->min_ci_balanced,
1078 &nbv->grp[eintNonlocal].nbl_lists,
1080 nbv->grp[eintNonlocal].kernel_type,
1083 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1085 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1087 /* initialize non-local pair-list on the GPU */
1088 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1089 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1092 wallcycle_stop(wcycle, ewcNS);
1096 wallcycle_start(wcycle, ewcMOVEX);
1097 dd_move_x(cr->dd, box, x);
1099 /* When we don't need the total dipole we sum it in global_stat */
1100 if (bStateChanged && NEED_MUTOT(*inputrec))
1102 gmx_sumd(2*DIM, mu, cr);
1104 wallcycle_stop(wcycle, ewcMOVEX);
1106 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1107 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1108 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1109 nbv->grp[eintNonlocal].nbat);
1110 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1111 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1114 if (bUseGPU && !bDiffKernels)
1116 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1117 /* launch non-local nonbonded F on GPU */
1118 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1120 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1126 /* launch D2H copy-back F */
1127 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1128 if (DOMAINDECOMP(cr) && !bDiffKernels)
1130 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1131 flags, eatNonlocal);
1133 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1135 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1138 if (bStateChanged && NEED_MUTOT(*inputrec))
1142 gmx_sumd(2*DIM, mu, cr);
1145 for (i = 0; i < 2; i++)
1147 for (j = 0; j < DIM; j++)
1149 fr->mu_tot[i][j] = mu[i*DIM + j];
1153 if (fr->efep == efepNO)
1155 copy_rvec(fr->mu_tot[0], mu_tot);
1159 for (j = 0; j < DIM; j++)
1162 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1163 lambda[efptCOUL]*fr->mu_tot[1][j];
1167 /* Reset energies */
1168 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1169 clear_rvecs(SHIFTS, fr->fshift);
1171 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1173 wallcycle_start(wcycle, ewcPPDURINGPME);
1174 dd_force_flop_start(cr->dd, nrnb);
1179 /* Enforced rotation has its own cycle counter that starts after the collective
1180 * coordinates have been communicated. It is added to ddCyclF to allow
1181 * for proper load-balancing */
1182 wallcycle_start(wcycle, ewcROT);
1183 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1184 wallcycle_stop(wcycle, ewcROT);
1187 /* Start the force cycle counter.
1188 * This counter is stopped in do_forcelow_level.
1189 * No parallel communication should occur while this counter is running,
1190 * since that will interfere with the dynamic load balancing.
1192 wallcycle_start(wcycle, ewcFORCE);
1195 /* Reset forces for which the virial is calculated separately:
1196 * PME/Ewald forces if necessary */
1197 if (fr->bF_NoVirSum)
1199 if (flags & GMX_FORCE_VIRIAL)
1201 fr->f_novirsum = fr->f_novirsum_alloc;
1204 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1208 clear_rvecs(homenr, fr->f_novirsum+start);
1213 /* We are not calculating the pressure so we do not need
1214 * a separate array for forces that do not contribute
1221 /* Clear the short- and long-range forces */
1222 clear_rvecs(fr->natoms_force_constr, f);
1223 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1225 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1228 clear_rvec(fr->vir_diag_posres);
1231 if (inputrec->ePull == epullCONSTRAINT)
1233 clear_pull_forces(inputrec->pull);
1236 /* We calculate the non-bonded forces, when done on the CPU, here.
1237 * We do this before calling do_force_lowlevel, as in there bondeds
1238 * forces are calculated before PME, which does communication.
1239 * With this order, non-bonded and bonded force calculation imbalance
1240 * can be balanced out by the domain decomposition load balancing.
1245 /* Maybe we should move this into do_force_lowlevel */
1246 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1250 if (fr->efep != efepNO)
1252 /* Calculate the local and non-local free energy interactions here.
1253 * Happens here on the CPU both with and without GPU.
1255 if (fr->nbv->grp[eintLocal].nbl_lists.nbl_fep[0]->nrj > 0)
1257 do_nb_verlet_fep(&fr->nbv->grp[eintLocal].nbl_lists,
1259 inputrec->fepvals, lambda,
1260 enerd, flags, nrnb, wcycle);
1263 if (DOMAINDECOMP(cr) &&
1264 fr->nbv->grp[eintNonlocal].nbl_lists.nbl_fep[0]->nrj > 0)
1266 do_nb_verlet_fep(&fr->nbv->grp[eintNonlocal].nbl_lists,
1268 inputrec->fepvals, lambda,
1269 enerd, flags, nrnb, wcycle);
1273 if (!bUseOrEmulGPU || bDiffKernels)
1277 if (DOMAINDECOMP(cr))
1279 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1280 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1290 aloc = eintNonlocal;
1293 /* Add all the non-bonded force to the normal force array.
1294 * This can be split into a local a non-local part when overlapping
1295 * communication with calculation with domain decomposition.
1297 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1298 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1299 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1300 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1301 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1302 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1303 wallcycle_start_nocount(wcycle, ewcFORCE);
1305 /* if there are multiple fshift output buffers reduce them */
1306 if ((flags & GMX_FORCE_VIRIAL) &&
1307 nbv->grp[aloc].nbl_lists.nnbl > 1)
1309 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1314 /* update QMMMrec, if necessary */
1317 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1320 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1322 posres_wrapper(flags, inputrec, nrnb, top, box, x,
1326 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1328 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1331 /* Compute the bonded and non-bonded energies and optionally forces */
1332 do_force_lowlevel(fr, inputrec, &(top->idef),
1333 cr, nrnb, wcycle, mdatoms,
1334 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1335 &(top->atomtypes), bBornRadii, box,
1336 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1337 flags, &cycles_pme);
1341 if (do_per_step(step, inputrec->nstcalclr))
1343 /* Add the long range forces to the short range forces */
1344 for (i = 0; i < fr->natoms_force_constr; i++)
1346 rvec_add(fr->f_twin[i], f[i], f[i]);
1351 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1355 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1358 if (bUseOrEmulGPU && !bDiffKernels)
1360 /* wait for non-local forces (or calculate in emulation mode) */
1361 if (DOMAINDECOMP(cr))
1367 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1368 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1369 nbv->grp[eintNonlocal].nbat,
1371 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1373 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1374 cycles_wait_gpu += cycles_tmp;
1375 cycles_force += cycles_tmp;
1379 wallcycle_start_nocount(wcycle, ewcFORCE);
1380 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1382 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1384 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1385 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1386 /* skip the reduction if there was no non-local work to do */
1387 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1389 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1390 nbv->grp[eintNonlocal].nbat, f);
1392 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1393 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1397 if (bDoForces && DOMAINDECOMP(cr))
1399 /* Communicate the forces */
1400 wallcycle_start(wcycle, ewcMOVEF);
1401 dd_move_f(cr->dd, f, fr->fshift);
1402 /* Do we need to communicate the separate force array
1403 * for terms that do not contribute to the single sum virial?
1404 * Position restraints and electric fields do not introduce
1405 * inter-cg forces, only full electrostatics methods do.
1406 * When we do not calculate the virial, fr->f_novirsum = f,
1407 * so we have already communicated these forces.
1409 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1410 (flags & GMX_FORCE_VIRIAL))
1412 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1416 /* We should not update the shift forces here,
1417 * since f_twin is already included in f.
1419 dd_move_f(cr->dd, fr->f_twin, NULL);
1421 wallcycle_stop(wcycle, ewcMOVEF);
1426 /* wait for local forces (or calculate in emulation mode) */
1429 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1430 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1431 nbv->grp[eintLocal].nbat,
1433 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1435 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1437 /* now clear the GPU outputs while we finish the step on the CPU */
1439 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1440 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1441 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1445 wallcycle_start_nocount(wcycle, ewcFORCE);
1446 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1447 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1449 wallcycle_stop(wcycle, ewcFORCE);
1451 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1452 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1453 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1455 /* skip the reduction if there was no non-local work to do */
1456 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1457 nbv->grp[eintLocal].nbat, f);
1459 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1460 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1463 if (DOMAINDECOMP(cr))
1465 dd_force_flop_stop(cr->dd, nrnb);
1468 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1471 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1478 if (IR_ELEC_FIELD(*inputrec))
1480 /* Compute forces due to electric field */
1481 calc_f_el(MASTER(cr) ? field : NULL,
1482 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1483 inputrec->ex, inputrec->et, t);
1486 /* If we have NoVirSum forces, but we do not calculate the virial,
1487 * we sum fr->f_novirum=f later.
1489 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1491 wallcycle_start(wcycle, ewcVSITESPREAD);
1492 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1493 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1494 wallcycle_stop(wcycle, ewcVSITESPREAD);
1498 wallcycle_start(wcycle, ewcVSITESPREAD);
1499 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1501 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1502 wallcycle_stop(wcycle, ewcVSITESPREAD);
1506 if (flags & GMX_FORCE_VIRIAL)
1508 /* Calculation of the virial must be done after vsites! */
1509 calc_virial(0, mdatoms->homenr, x, f,
1510 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1514 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1516 pull_potential_wrapper(cr, inputrec, box, x,
1517 f, vir_force, mdatoms, enerd, lambda, t);
1520 /* Add the forces from enforced rotation potentials (if any) */
1523 wallcycle_start(wcycle, ewcROTadd);
1524 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1525 wallcycle_stop(wcycle, ewcROTadd);
1528 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1529 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1531 if (PAR(cr) && !(cr->duty & DUTY_PME))
1533 /* In case of node-splitting, the PP nodes receive the long-range
1534 * forces, virial and energy from the PME nodes here.
1536 pme_receive_force_ener(cr, wcycle, enerd, fr);
1541 post_process_forces(cr, step, nrnb, wcycle,
1542 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1546 /* Sum the potential energy terms from group contributions */
1547 sum_epot(&(enerd->grpp), enerd->term);
1550 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1551 t_inputrec *inputrec,
1552 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1553 gmx_localtop_t *top,
1554 gmx_groups_t *groups,
1555 matrix box, rvec x[], history_t *hist,
1559 gmx_enerdata_t *enerd, t_fcdata *fcd,
1560 real *lambda, t_graph *graph,
1561 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1562 double t, FILE *field, gmx_edsam_t ed,
1563 gmx_bool bBornRadii,
1569 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1570 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1571 gmx_bool bDoAdressWF;
1573 rvec vzero, box_diag;
1574 real e, v, dvdlambda[efptNR];
1576 float cycles_pme, cycles_force;
1579 homenr = mdatoms->homenr;
1581 clear_mat(vir_force);
1584 if (DOMAINDECOMP(cr))
1586 cg1 = cr->dd->ncg_tot;
1597 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1598 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1599 /* Should we update the long-range neighborlists at this step? */
1600 bDoLongRangeNS = fr->bTwinRange && bNS;
1601 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1602 bFillGrid = (bNS && bStateChanged);
1603 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1604 bDoForces = (flags & GMX_FORCE_FORCES);
1605 bDoPotential = (flags & GMX_FORCE_ENERGY);
1606 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1607 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1609 /* should probably move this to the forcerec since it doesn't change */
1610 bDoAdressWF = ((fr->adress_type != eAdressOff));
1614 update_forcerec(fr, box);
1616 if (NEED_MUTOT(*inputrec))
1618 /* Calculate total (local) dipole moment in a temporary common array.
1619 * This makes it possible to sum them over nodes faster.
1621 calc_mu(start, homenr,
1622 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1627 if (fr->ePBC != epbcNONE)
1629 /* Compute shift vectors every step,
1630 * because of pressure coupling or box deformation!
1632 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1634 calc_shifts(box, fr->shift_vec);
1639 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1640 &(top->cgs), x, fr->cg_cm);
1641 inc_nrnb(nrnb, eNR_CGCM, homenr);
1642 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1644 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1646 unshift_self(graph, box, x);
1651 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1652 inc_nrnb(nrnb, eNR_CGCM, homenr);
1655 if (bCalcCGCM && gmx_debug_at)
1657 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1661 if (!(cr->duty & DUTY_PME))
1663 /* Send particle coordinates to the pme nodes.
1664 * Since this is only implemented for domain decomposition
1665 * and domain decomposition does not use the graph,
1666 * we do not need to worry about shifting.
1671 wallcycle_start(wcycle, ewcPP_PMESENDX);
1673 bBS = (inputrec->nwall == 2);
1676 copy_mat(box, boxs);
1677 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1680 if (EEL_PME(fr->eeltype))
1682 pme_flags |= GMX_PME_DO_COULOMB;
1685 if (EVDW_PME(fr->vdwtype))
1687 pme_flags |= GMX_PME_DO_LJ;
1690 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1691 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1692 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1695 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1697 #endif /* GMX_MPI */
1699 /* Communicate coordinates and sum dipole if necessary */
1700 if (DOMAINDECOMP(cr))
1702 wallcycle_start(wcycle, ewcMOVEX);
1703 dd_move_x(cr->dd, box, x);
1704 wallcycle_stop(wcycle, ewcMOVEX);
1707 /* update adress weight beforehand */
1708 if (bStateChanged && bDoAdressWF)
1710 /* need pbc for adress weight calculation with pbc_dx */
1711 set_pbc(&pbc, inputrec->ePBC, box);
1712 if (fr->adress_site == eAdressSITEcog)
1714 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1715 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1717 else if (fr->adress_site == eAdressSITEcom)
1719 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1720 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1722 else if (fr->adress_site == eAdressSITEatomatom)
1724 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1725 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1729 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1730 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1734 if (NEED_MUTOT(*inputrec))
1741 gmx_sumd(2*DIM, mu, cr);
1743 for (i = 0; i < 2; i++)
1745 for (j = 0; j < DIM; j++)
1747 fr->mu_tot[i][j] = mu[i*DIM + j];
1751 if (fr->efep == efepNO)
1753 copy_rvec(fr->mu_tot[0], mu_tot);
1757 for (j = 0; j < DIM; j++)
1760 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1765 /* Reset energies */
1766 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1767 clear_rvecs(SHIFTS, fr->fshift);
1771 wallcycle_start(wcycle, ewcNS);
1773 if (graph && bStateChanged)
1775 /* Calculate intramolecular shift vectors to make molecules whole */
1776 mk_mshift(fplog, graph, fr->ePBC, box, x);
1779 /* Do the actual neighbour searching */
1781 groups, top, mdatoms,
1782 cr, nrnb, bFillGrid,
1785 wallcycle_stop(wcycle, ewcNS);
1788 if (inputrec->implicit_solvent && bNS)
1790 make_gb_nblist(cr, inputrec->gb_algorithm,
1791 x, box, fr, &top->idef, graph, fr->born);
1794 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1796 wallcycle_start(wcycle, ewcPPDURINGPME);
1797 dd_force_flop_start(cr->dd, nrnb);
1802 /* Enforced rotation has its own cycle counter that starts after the collective
1803 * coordinates have been communicated. It is added to ddCyclF to allow
1804 * for proper load-balancing */
1805 wallcycle_start(wcycle, ewcROT);
1806 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1807 wallcycle_stop(wcycle, ewcROT);
1810 /* Start the force cycle counter.
1811 * This counter is stopped in do_forcelow_level.
1812 * No parallel communication should occur while this counter is running,
1813 * since that will interfere with the dynamic load balancing.
1815 wallcycle_start(wcycle, ewcFORCE);
1819 /* Reset forces for which the virial is calculated separately:
1820 * PME/Ewald forces if necessary */
1821 if (fr->bF_NoVirSum)
1823 if (flags & GMX_FORCE_VIRIAL)
1825 fr->f_novirsum = fr->f_novirsum_alloc;
1828 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1832 clear_rvecs(homenr, fr->f_novirsum+start);
1837 /* We are not calculating the pressure so we do not need
1838 * a separate array for forces that do not contribute
1845 /* Clear the short- and long-range forces */
1846 clear_rvecs(fr->natoms_force_constr, f);
1847 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1849 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1852 clear_rvec(fr->vir_diag_posres);
1854 if (inputrec->ePull == epullCONSTRAINT)
1856 clear_pull_forces(inputrec->pull);
1859 /* update QMMMrec, if necessary */
1862 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1865 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1867 posres_wrapper(flags, inputrec, nrnb, top, box, x,
1871 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1873 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1876 /* Compute the bonded and non-bonded energies and optionally forces */
1877 do_force_lowlevel(fr, inputrec, &(top->idef),
1878 cr, nrnb, wcycle, mdatoms,
1879 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1880 &(top->atomtypes), bBornRadii, box,
1881 inputrec->fepvals, lambda,
1882 graph, &(top->excls), fr->mu_tot,
1888 if (do_per_step(step, inputrec->nstcalclr))
1890 /* Add the long range forces to the short range forces */
1891 for (i = 0; i < fr->natoms_force_constr; i++)
1893 rvec_add(fr->f_twin[i], f[i], f[i]);
1898 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1902 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1905 if (DOMAINDECOMP(cr))
1907 dd_force_flop_stop(cr->dd, nrnb);
1910 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1916 if (IR_ELEC_FIELD(*inputrec))
1918 /* Compute forces due to electric field */
1919 calc_f_el(MASTER(cr) ? field : NULL,
1920 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1921 inputrec->ex, inputrec->et, t);
1924 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1926 /* Compute thermodynamic force in hybrid AdResS region */
1927 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1928 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1931 /* Communicate the forces */
1932 if (DOMAINDECOMP(cr))
1934 wallcycle_start(wcycle, ewcMOVEF);
1935 dd_move_f(cr->dd, f, fr->fshift);
1936 /* Do we need to communicate the separate force array
1937 * for terms that do not contribute to the single sum virial?
1938 * Position restraints and electric fields do not introduce
1939 * inter-cg forces, only full electrostatics methods do.
1940 * When we do not calculate the virial, fr->f_novirsum = f,
1941 * so we have already communicated these forces.
1943 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1944 (flags & GMX_FORCE_VIRIAL))
1946 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1950 /* We should not update the shift forces here,
1951 * since f_twin is already included in f.
1953 dd_move_f(cr->dd, fr->f_twin, NULL);
1955 wallcycle_stop(wcycle, ewcMOVEF);
1958 /* If we have NoVirSum forces, but we do not calculate the virial,
1959 * we sum fr->f_novirum=f later.
1961 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1963 wallcycle_start(wcycle, ewcVSITESPREAD);
1964 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1965 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1966 wallcycle_stop(wcycle, ewcVSITESPREAD);
1970 wallcycle_start(wcycle, ewcVSITESPREAD);
1971 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1973 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1974 wallcycle_stop(wcycle, ewcVSITESPREAD);
1978 if (flags & GMX_FORCE_VIRIAL)
1980 /* Calculation of the virial must be done after vsites! */
1981 calc_virial(0, mdatoms->homenr, x, f,
1982 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1986 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1988 pull_potential_wrapper(cr, inputrec, box, x,
1989 f, vir_force, mdatoms, enerd, lambda, t);
1992 /* Add the forces from enforced rotation potentials (if any) */
1995 wallcycle_start(wcycle, ewcROTadd);
1996 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1997 wallcycle_stop(wcycle, ewcROTadd);
2000 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
2001 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
2003 if (PAR(cr) && !(cr->duty & DUTY_PME))
2005 /* In case of node-splitting, the PP nodes receive the long-range
2006 * forces, virial and energy from the PME nodes here.
2008 pme_receive_force_ener(cr, wcycle, enerd, fr);
2013 post_process_forces(cr, step, nrnb, wcycle,
2014 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
2018 /* Sum the potential energy terms from group contributions */
2019 sum_epot(&(enerd->grpp), enerd->term);
2022 void do_force(FILE *fplog, t_commrec *cr,
2023 t_inputrec *inputrec,
2024 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2025 gmx_localtop_t *top,
2026 gmx_groups_t *groups,
2027 matrix box, rvec x[], history_t *hist,
2031 gmx_enerdata_t *enerd, t_fcdata *fcd,
2032 real *lambda, t_graph *graph,
2034 gmx_vsite_t *vsite, rvec mu_tot,
2035 double t, FILE *field, gmx_edsam_t ed,
2036 gmx_bool bBornRadii,
2039 /* modify force flag if not doing nonbonded */
2040 if (!fr->bNonbonded)
2042 flags &= ~GMX_FORCE_NONBONDED;
2045 switch (inputrec->cutoff_scheme)
2048 do_force_cutsVERLET(fplog, cr, inputrec,
2064 do_force_cutsGROUP(fplog, cr, inputrec,
2079 gmx_incons("Invalid cut-off scheme passed!");
2084 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2085 t_inputrec *ir, t_mdatoms *md,
2086 t_state *state, t_commrec *cr, t_nrnb *nrnb,
2087 t_forcerec *fr, gmx_localtop_t *top)
2089 int i, m, start, end;
2091 real dt = ir->delta_t;
2095 snew(savex, state->natoms);
2102 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2103 start, md->homenr, end);
2105 /* Do a first constrain to reset particles... */
2106 step = ir->init_step;
2109 char buf[STEPSTRSIZE];
2110 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2111 gmx_step_str(step, buf));
2115 /* constrain the current position */
2116 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2117 ir, NULL, cr, step, 0, 1.0, md,
2118 state->x, state->x, NULL,
2119 fr->bMolPBC, state->box,
2120 state->lambda[efptBONDED], &dvdl_dum,
2121 NULL, NULL, nrnb, econqCoord,
2122 ir->epc == epcMTTK, state->veta, state->veta);
2125 /* constrain the inital velocity, and save it */
2126 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2127 /* might not yet treat veta correctly */
2128 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2129 ir, NULL, cr, step, 0, 1.0, md,
2130 state->x, state->v, state->v,
2131 fr->bMolPBC, state->box,
2132 state->lambda[efptBONDED], &dvdl_dum,
2133 NULL, NULL, nrnb, econqVeloc,
2134 ir->epc == epcMTTK, state->veta, state->veta);
2136 /* constrain the inital velocities at t-dt/2 */
2137 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2139 for (i = start; (i < end); i++)
2141 for (m = 0; (m < DIM); m++)
2143 /* Reverse the velocity */
2144 state->v[i][m] = -state->v[i][m];
2145 /* Store the position at t-dt in buf */
2146 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2149 /* Shake the positions at t=-dt with the positions at t=0
2150 * as reference coordinates.
2154 char buf[STEPSTRSIZE];
2155 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2156 gmx_step_str(step, buf));
2159 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2160 ir, NULL, cr, step, -1, 1.0, md,
2161 state->x, savex, NULL,
2162 fr->bMolPBC, state->box,
2163 state->lambda[efptBONDED], &dvdl_dum,
2164 state->v, NULL, nrnb, econqCoord,
2165 ir->epc == epcMTTK, state->veta, state->veta);
2167 for (i = start; i < end; i++)
2169 for (m = 0; m < DIM; m++)
2171 /* Re-reverse the velocities */
2172 state->v[i][m] = -state->v[i][m];
2181 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2182 double *enerout, double *virout)
2184 double enersum, virsum;
2185 double invscale, invscale2, invscale3;
2186 double r, ea, eb, ec, pa, pb, pc, pd;
2188 int ri, offset, tabfactor;
2190 invscale = 1.0/scale;
2191 invscale2 = invscale*invscale;
2192 invscale3 = invscale*invscale2;
2194 /* Following summation derived from cubic spline definition,
2195 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2196 * the cubic spline. We first calculate the negative of the
2197 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2198 * add the more standard, abrupt cutoff correction to that result,
2199 * yielding the long-range correction for a switched function. We
2200 * perform both the pressure and energy loops at the same time for
2201 * simplicity, as the computational cost is low. */
2205 /* Since the dispersion table has been scaled down a factor
2206 * 6.0 and the repulsion a factor 12.0 to compensate for the
2207 * c6/c12 parameters inside nbfp[] being scaled up (to save
2208 * flops in kernels), we need to correct for this.
2219 for (ri = rstart; ri < rend; ++ri)
2223 eb = 2.0*invscale2*r;
2227 pb = 3.0*invscale2*r;
2228 pc = 3.0*invscale*r*r;
2231 /* this "8" is from the packing in the vdwtab array - perhaps
2232 should be defined? */
2234 offset = 8*ri + offstart;
2235 y0 = vdwtab[offset];
2236 f = vdwtab[offset+1];
2237 g = vdwtab[offset+2];
2238 h = vdwtab[offset+3];
2240 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);
2241 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);
2243 *enerout = 4.0*M_PI*enersum*tabfactor;
2244 *virout = 4.0*M_PI*virsum*tabfactor;
2247 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2249 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2250 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2251 double invscale, invscale2, invscale3;
2252 int ri0, ri1, ri, i, offstart, offset;
2253 real scale, *vdwtab, tabfactor, tmp;
2255 fr->enershiftsix = 0;
2256 fr->enershifttwelve = 0;
2257 fr->enerdiffsix = 0;
2258 fr->enerdifftwelve = 0;
2260 fr->virdifftwelve = 0;
2262 if (eDispCorr != edispcNO)
2264 for (i = 0; i < 2; i++)
2269 if ((fr->vdw_modifier == eintmodPOTSHIFT) ||
2270 (fr->vdw_modifier == eintmodPOTSWITCH) ||
2271 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2272 (fr->vdwtype == evdwSHIFT) ||
2273 (fr->vdwtype == evdwSWITCH))
2275 if (((fr->vdw_modifier == eintmodPOTSWITCH) ||
2276 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2277 (fr->vdwtype == evdwSWITCH)) && fr->rvdw_switch == 0)
2280 "With dispersion correction rvdw-switch can not be zero "
2281 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2284 scale = fr->nblists[0].table_vdw.scale;
2285 vdwtab = fr->nblists[0].table_vdw.data;
2287 /* Round the cut-offs to exact table values for precision */
2288 ri0 = floor(fr->rvdw_switch*scale);
2289 ri1 = ceil(fr->rvdw*scale);
2291 /* The code below has some support for handling force-switching, i.e.
2292 * when the force (instead of potential) is switched over a limited
2293 * region. This leads to a constant shift in the potential inside the
2294 * switching region, which we can handle by adding a constant energy
2295 * term in the force-switch case just like when we do potential-shift.
2297 * For now this is not enabled, but to keep the functionality in the
2298 * code we check separately for switch and shift. When we do force-switch
2299 * the shifting point is rvdw_switch, while it is the cutoff when we
2300 * have a classical potential-shift.
2302 * For a pure potential-shift the potential has a constant shift
2303 * all the way out to the cutoff, and that is it. For other forms
2304 * we need to calculate the constant shift up to the point where we
2305 * start modifying the potential.
2307 ri0 = (fr->vdw_modifier == eintmodPOTSHIFT) ? ri1 : ri0;
2314 if ((fr->vdw_modifier == eintmodFORCESWITCH) ||
2315 (fr->vdwtype == evdwSHIFT))
2317 /* Determine the constant energy shift below rvdw_switch.
2318 * Table has a scale factor since we have scaled it down to compensate
2319 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2321 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2322 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2324 else if (fr->vdw_modifier == eintmodPOTSHIFT)
2326 fr->enershiftsix = (real)(-1.0/(rc3*rc3));
2327 fr->enershifttwelve = (real)( 1.0/(rc9*rc3));
2330 /* Add the constant part from 0 to rvdw_switch.
2331 * This integration from 0 to rvdw_switch overcounts the number
2332 * of interactions by 1, as it also counts the self interaction.
2333 * We will correct for this later.
2335 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2336 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2338 /* Calculate the contribution in the range [r0,r1] where we
2339 * modify the potential. For a pure potential-shift modifier we will
2340 * have ri0==ri1, and there will not be any contribution here.
2342 for (i = 0; i < 2; i++)
2346 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2347 eners[i] -= enersum;
2351 /* Alright: Above we compensated by REMOVING the parts outside r0
2352 * corresponding to the ideal VdW 1/r6 and /r12 potentials.
2354 * Regardless of whether r0 is the point where we start switching,
2355 * or the cutoff where we calculated the constant shift, we include
2356 * all the parts we are missing out to infinity from r0 by
2357 * calculating the analytical dispersion correction.
2359 eners[0] += -4.0*M_PI/(3.0*rc3);
2360 eners[1] += 4.0*M_PI/(9.0*rc9);
2361 virs[0] += 8.0*M_PI/rc3;
2362 virs[1] += -16.0*M_PI/(3.0*rc9);
2364 else if (fr->vdwtype == evdwCUT ||
2365 EVDW_PME(fr->vdwtype) ||
2366 fr->vdwtype == evdwUSER)
2368 if (fr->vdwtype == evdwUSER && fplog)
2371 "WARNING: using dispersion correction with user tables\n");
2374 /* Note that with LJ-PME, the dispersion correction is multiplied
2375 * by the difference between the actual C6 and the value of C6
2376 * that would produce the combination rule.
2377 * This means the normal energy and virial difference formulas
2381 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2383 /* Contribution beyond the cut-off */
2384 eners[0] += -4.0*M_PI/(3.0*rc3);
2385 eners[1] += 4.0*M_PI/(9.0*rc9);
2386 if (fr->vdw_modifier == eintmodPOTSHIFT)
2388 /* Contribution within the cut-off */
2389 eners[0] += -4.0*M_PI/(3.0*rc3);
2390 eners[1] += 4.0*M_PI/(3.0*rc9);
2392 /* Contribution beyond the cut-off */
2393 virs[0] += 8.0*M_PI/rc3;
2394 virs[1] += -16.0*M_PI/(3.0*rc9);
2399 "Dispersion correction is not implemented for vdw-type = %s",
2400 evdw_names[fr->vdwtype]);
2403 /* When we deprecate the group kernels the code below can go too */
2404 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2406 /* Calculate self-interaction coefficient (assuming that
2407 * the reciprocal-space contribution is constant in the
2408 * region that contributes to the self-interaction).
2410 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2412 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2413 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2416 fr->enerdiffsix = eners[0];
2417 fr->enerdifftwelve = eners[1];
2418 /* The 0.5 is due to the Gromacs definition of the virial */
2419 fr->virdiffsix = 0.5*virs[0];
2420 fr->virdifftwelve = 0.5*virs[1];
2424 void calc_dispcorr(t_inputrec *ir, t_forcerec *fr,
2426 matrix box, real lambda, tensor pres, tensor virial,
2427 real *prescorr, real *enercorr, real *dvdlcorr)
2429 gmx_bool bCorrAll, bCorrPres;
2430 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2440 if (ir->eDispCorr != edispcNO)
2442 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2443 ir->eDispCorr == edispcAllEnerPres);
2444 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2445 ir->eDispCorr == edispcAllEnerPres);
2447 invvol = 1/det(box);
2450 /* Only correct for the interactions with the inserted molecule */
2451 dens = (natoms - fr->n_tpi)*invvol;
2456 dens = natoms*invvol;
2457 ninter = 0.5*natoms;
2460 if (ir->efep == efepNO)
2462 avcsix = fr->avcsix[0];
2463 avctwelve = fr->avctwelve[0];
2467 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2468 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2471 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2472 *enercorr += avcsix*enerdiff;
2474 if (ir->efep != efepNO)
2476 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2480 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2481 *enercorr += avctwelve*enerdiff;
2482 if (fr->efep != efepNO)
2484 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2490 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2491 if (ir->eDispCorr == edispcAllEnerPres)
2493 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2495 /* The factor 2 is because of the Gromacs virial definition */
2496 spres = -2.0*invvol*svir*PRESFAC;
2498 for (m = 0; m < DIM; m++)
2500 virial[m][m] += svir;
2501 pres[m][m] += spres;
2506 /* Can't currently control when it prints, for now, just print when degugging */
2511 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2517 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2518 *enercorr, spres, svir);
2522 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2526 if (fr->efep != efepNO)
2528 *dvdlcorr += dvdlambda;
2533 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2534 t_graph *graph, rvec x[])
2538 fprintf(fplog, "Removing pbc first time\n");
2540 calc_shifts(box, fr->shift_vec);
2543 mk_mshift(fplog, graph, fr->ePBC, box, x);
2546 p_graph(debug, "do_pbc_first 1", graph);
2548 shift_self(graph, box, x);
2549 /* By doing an extra mk_mshift the molecules that are broken
2550 * because they were e.g. imported from another software
2551 * will be made whole again. Such are the healing powers
2554 mk_mshift(fplog, graph, fr->ePBC, box, x);
2557 p_graph(debug, "do_pbc_first 2", graph);
2562 fprintf(fplog, "Done rmpbc\n");
2566 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2567 gmx_mtop_t *mtop, rvec x[],
2572 gmx_molblock_t *molb;
2574 if (bFirst && fplog)
2576 fprintf(fplog, "Removing pbc first time\n");
2581 for (mb = 0; mb < mtop->nmolblock; mb++)
2583 molb = &mtop->molblock[mb];
2584 if (molb->natoms_mol == 1 ||
2585 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2587 /* Just one atom or charge group in the molecule, no PBC required */
2588 as += molb->nmol*molb->natoms_mol;
2592 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2593 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2594 0, molb->natoms_mol, FALSE, FALSE, graph);
2596 for (mol = 0; mol < molb->nmol; mol++)
2598 mk_mshift(fplog, graph, ePBC, box, x+as);
2600 shift_self(graph, box, x+as);
2601 /* The molecule is whole now.
2602 * We don't need the second mk_mshift call as in do_pbc_first,
2603 * since we no longer need this graph.
2606 as += molb->natoms_mol;
2614 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2615 gmx_mtop_t *mtop, rvec x[])
2617 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2620 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2621 gmx_mtop_t *mtop, rvec x[])
2623 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2626 void finish_run(FILE *fplog, t_commrec *cr,
2627 t_inputrec *inputrec,
2628 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2629 gmx_walltime_accounting_t walltime_accounting,
2630 nonbonded_verlet_t *nbv,
2631 gmx_bool bWriteStat)
2634 t_nrnb *nrnb_tot = NULL;
2637 double elapsed_time,
2638 elapsed_time_over_all_ranks,
2639 elapsed_time_over_all_threads,
2640 elapsed_time_over_all_threads_over_all_ranks;
2641 wallcycle_sum(cr, wcycle);
2647 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2648 cr->mpi_comm_mysim);
2656 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2657 elapsed_time_over_all_ranks = elapsed_time;
2658 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2659 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2663 /* reduce elapsed_time over all MPI ranks in the current simulation */
2664 MPI_Allreduce(&elapsed_time,
2665 &elapsed_time_over_all_ranks,
2666 1, MPI_DOUBLE, MPI_SUM,
2667 cr->mpi_comm_mysim);
2668 elapsed_time_over_all_ranks /= cr->nnodes;
2669 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2670 * current simulation. */
2671 MPI_Allreduce(&elapsed_time_over_all_threads,
2672 &elapsed_time_over_all_threads_over_all_ranks,
2673 1, MPI_DOUBLE, MPI_SUM,
2674 cr->mpi_comm_mysim);
2680 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2687 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2689 print_dd_statistics(cr, inputrec, fplog);
2694 wallclock_gpu_t* gputimes = use_GPU(nbv) ?
2695 nbnxn_cuda_get_timings(nbv->cu_nbv) : NULL;
2696 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2697 elapsed_time_over_all_ranks,
2700 if (EI_DYNAMICS(inputrec->eI))
2702 delta_t = inputrec->delta_t;
2711 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2712 elapsed_time_over_all_ranks,
2713 walltime_accounting_get_nsteps_done(walltime_accounting),
2714 delta_t, nbfs, mflop);
2718 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2719 elapsed_time_over_all_ranks,
2720 walltime_accounting_get_nsteps_done(walltime_accounting),
2721 delta_t, nbfs, mflop);
2726 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2728 /* this function works, but could probably use a logic rewrite to keep all the different
2729 types of efep straight. */
2732 t_lambda *fep = ir->fepvals;
2734 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2736 for (i = 0; i < efptNR; i++)
2748 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2749 if checkpoint is set -- a kludge is in for now
2751 for (i = 0; i < efptNR; i++)
2753 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2754 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2756 lambda[i] = fep->init_lambda;
2759 lam0[i] = lambda[i];
2764 lambda[i] = fep->all_lambda[i][*fep_state];
2767 lam0[i] = lambda[i];
2773 /* need to rescale control temperatures to match current state */
2774 for (i = 0; i < ir->opts.ngtc; i++)
2776 if (ir->opts.ref_t[i] > 0)
2778 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2784 /* Send to the log the information on the current lambdas */
2787 fprintf(fplog, "Initial vector of lambda components:[ ");
2788 for (i = 0; i < efptNR; i++)
2790 fprintf(fplog, "%10.4f ", lambda[i]);
2792 fprintf(fplog, "]\n");
2798 void init_md(FILE *fplog,
2799 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2800 double *t, double *t0,
2801 real *lambda, int *fep_state, double *lam0,
2802 t_nrnb *nrnb, gmx_mtop_t *mtop,
2804 int nfile, const t_filenm fnm[],
2805 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2806 tensor force_vir, tensor shake_vir, rvec mu_tot,
2807 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags)
2812 /* Initial values */
2813 *t = *t0 = ir->init_t;
2816 for (i = 0; i < ir->opts.ngtc; i++)
2818 /* set bSimAnn if any group is being annealed */
2819 if (ir->opts.annealing[i] != eannNO)
2826 update_annealing_target_temp(&(ir->opts), ir->init_t);
2829 /* Initialize lambda variables */
2830 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2834 *upd = init_update(ir);
2840 *vcm = init_vcm(fplog, &mtop->groups, ir);
2843 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2845 if (ir->etc == etcBERENDSEN)
2847 please_cite(fplog, "Berendsen84a");
2849 if (ir->etc == etcVRESCALE)
2851 please_cite(fplog, "Bussi2007a");
2853 if (ir->eI == eiSD1)
2855 please_cite(fplog, "Goga2012");
2863 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv);
2865 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2866 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2871 please_cite(fplog, "Fritsch12");
2872 please_cite(fplog, "Junghans10");
2874 /* Initiate variables */
2875 clear_mat(force_vir);
2876 clear_mat(shake_vir);