2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 #include "gromacs/legacyheaders/sim_util.h"
48 #ifdef HAVE_SYS_TIME_H
52 #include "gromacs/essentialdynamics/edsam.h"
53 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
54 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
55 #include "gromacs/imd/imd.h"
56 #include "gromacs/legacyheaders/calcmu.h"
57 #include "gromacs/legacyheaders/chargegroup.h"
58 #include "gromacs/legacyheaders/constr.h"
59 #include "gromacs/legacyheaders/copyrite.h"
60 #include "gromacs/legacyheaders/disre.h"
61 #include "gromacs/legacyheaders/domdec.h"
62 #include "gromacs/legacyheaders/force.h"
63 #include "gromacs/legacyheaders/genborn.h"
64 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
65 #include "gromacs/legacyheaders/mdatoms.h"
66 #include "gromacs/legacyheaders/mdrun.h"
67 #include "gromacs/legacyheaders/names.h"
68 #include "gromacs/legacyheaders/network.h"
69 #include "gromacs/legacyheaders/nonbonded.h"
70 #include "gromacs/legacyheaders/nrnb.h"
71 #include "gromacs/legacyheaders/orires.h"
72 #include "gromacs/legacyheaders/pme.h"
73 #include "gromacs/legacyheaders/qmmm.h"
74 #include "gromacs/legacyheaders/txtdump.h"
75 #include "gromacs/legacyheaders/typedefs.h"
76 #include "gromacs/legacyheaders/update.h"
77 #include "gromacs/legacyheaders/types/commrec.h"
78 #include "gromacs/listed-forces/bonded.h"
79 #include "gromacs/math/units.h"
80 #include "gromacs/math/vec.h"
81 #include "gromacs/mdlib/nb_verlet.h"
82 #include "gromacs/mdlib/nbnxn_atomdata.h"
83 #include "gromacs/mdlib/nbnxn_search.h"
84 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda.h"
85 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_data_mgmt.h"
86 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
87 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref.h"
88 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
89 #include "gromacs/mdlib/nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
90 #include "gromacs/pbcutil/ishift.h"
91 #include "gromacs/pbcutil/mshift.h"
92 #include "gromacs/pbcutil/pbc.h"
93 #include "gromacs/pulling/pull.h"
94 #include "gromacs/pulling/pull_rotation.h"
95 #include "gromacs/timing/wallcycle.h"
96 #include "gromacs/timing/walltime_accounting.h"
97 #include "gromacs/utility/cstringutil.h"
98 #include "gromacs/utility/gmxmpi.h"
99 #include "gromacs/utility/smalloc.h"
103 void print_time(FILE *out,
104 gmx_walltime_accounting_t walltime_accounting,
107 t_commrec gmx_unused *cr)
110 char timebuf[STRLEN];
111 double dt, elapsed_seconds, time_per_step;
114 #ifndef GMX_THREAD_MPI
120 fprintf(out, "step %s", gmx_step_str(step, buf));
121 if ((step >= ir->nstlist))
123 double seconds_since_epoch = gmx_gettime();
124 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
125 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
126 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
132 finish = (time_t) (seconds_since_epoch + dt);
133 gmx_ctime_r(&finish, timebuf, STRLEN);
134 sprintf(buf, "%s", timebuf);
135 buf[strlen(buf)-1] = '\0';
136 fprintf(out, ", will finish %s", buf);
140 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
145 fprintf(out, " performance: %.1f ns/day ",
146 ir->delta_t/1000*24*60*60/time_per_step);
149 #ifndef GMX_THREAD_MPI
159 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
162 char time_string[STRLEN];
171 char timebuf[STRLEN];
172 time_t temp_time = (time_t) the_time;
174 gmx_ctime_r(&temp_time, timebuf, STRLEN);
175 for (i = 0; timebuf[i] >= ' '; i++)
177 time_string[i] = timebuf[i];
179 time_string[i] = '\0';
182 fprintf(fplog, "%s on rank %d %s\n", title, nodeid, time_string);
185 void print_start(FILE *fplog, t_commrec *cr,
186 gmx_walltime_accounting_t walltime_accounting,
191 sprintf(buf, "Started %s", name);
192 print_date_and_time(fplog, cr->nodeid, buf,
193 walltime_accounting_get_start_time_stamp(walltime_accounting));
196 static void sum_forces(int start, int end, rvec f[], rvec flr[])
202 pr_rvecs(debug, 0, "fsr", f+start, end-start);
203 pr_rvecs(debug, 0, "flr", flr+start, end-start);
205 for (i = start; (i < end); i++)
207 rvec_inc(f[i], flr[i]);
212 * calc_f_el calculates forces due to an electric field.
214 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
216 * Et[] contains the parameters for the time dependent
217 * part of the field (not yet used).
218 * Ex[] contains the parameters for
219 * the spatial dependent part of the field. You can have cool periodic
220 * fields in principle, but only a constant field is supported
222 * The function should return the energy due to the electric field
223 * (if any) but for now returns 0.
226 * There can be problems with the virial.
227 * Since the field is not self-consistent this is unavoidable.
228 * For neutral molecules the virial is correct within this approximation.
229 * For neutral systems with many charged molecules the error is small.
230 * But for systems with a net charge or a few charged molecules
231 * the error can be significant when the field is high.
232 * Solution: implement a self-consitent electric field into PME.
234 static void calc_f_el(FILE *fp, int start, int homenr,
235 real charge[], rvec f[],
236 t_cosines Ex[], t_cosines Et[], double t)
242 for (m = 0; (m < DIM); m++)
249 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
253 Ext[m] = cos(Et[m].a[0]*t);
262 /* Convert the field strength from V/nm to MD-units */
263 Ext[m] *= Ex[m].a[0]*FIELDFAC;
264 for (i = start; (i < start+homenr); i++)
266 f[i][m] += charge[i]*Ext[m];
276 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
277 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
281 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
282 tensor vir_part, t_graph *graph, matrix box,
283 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
288 /* The short-range virial from surrounding boxes */
290 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
291 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
293 /* Calculate partial virial, for local atoms only, based on short range.
294 * Total virial is computed in global_stat, called from do_md
296 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
297 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
299 /* Add position restraint contribution */
300 for (i = 0; i < DIM; i++)
302 vir_part[i][i] += fr->vir_diag_posres[i];
305 /* Add wall contribution */
306 for (i = 0; i < DIM; i++)
308 vir_part[i][ZZ] += fr->vir_wall_z[i];
313 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
317 static void posres_wrapper(int flags,
321 matrix box, rvec x[],
322 gmx_enerdata_t *enerd,
330 /* Position restraints always require full pbc */
331 set_pbc(&pbc, ir->ePBC, box);
333 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
334 top->idef.iparams_posres,
335 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
336 ir->ePBC == epbcNONE ? NULL : &pbc,
337 lambda[efptRESTRAINT], &dvdl,
338 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
339 enerd->term[F_POSRES] += v;
340 /* If just the force constant changes, the FEP term is linear,
341 * but if k changes, it is not.
343 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
344 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
346 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
348 for (i = 0; i < enerd->n_lambda; i++)
350 real dvdl_dum, lambda_dum;
352 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
353 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
354 top->idef.iparams_posres,
355 (const rvec*)x, NULL, NULL,
356 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
357 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
358 enerd->enerpart_lambda[i] += v;
363 static void fbposres_wrapper(t_inputrec *ir,
366 matrix box, rvec x[],
367 gmx_enerdata_t *enerd,
373 /* Flat-bottomed position restraints always require full pbc */
374 set_pbc(&pbc, ir->ePBC, box);
375 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
376 top->idef.iparams_fbposres,
377 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
378 ir->ePBC == epbcNONE ? NULL : &pbc,
379 fr->rc_scaling, fr->ePBC, fr->posres_com);
380 enerd->term[F_FBPOSRES] += v;
381 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
384 static void pull_potential_wrapper(t_commrec *cr,
386 matrix box, rvec x[],
390 gmx_enerdata_t *enerd,
393 gmx_wallcycle_t wcycle)
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 wallcycle_start(wcycle, ewcPULLPOT);
404 set_pbc(&pbc, ir->ePBC, box);
406 enerd->term[F_COM_PULL] +=
407 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
408 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
409 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
410 wallcycle_stop(wcycle, ewcPULLPOT);
413 static void pme_receive_force_ener(t_commrec *cr,
414 gmx_wallcycle_t wcycle,
415 gmx_enerdata_t *enerd,
418 real e_q, e_lj, v, dvdl_q, dvdl_lj;
419 float cycles_ppdpme, cycles_seppme;
421 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
422 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
424 /* In case of node-splitting, the PP nodes receive the long-range
425 * forces, virial and energy from the PME nodes here.
427 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
430 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
431 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
433 enerd->term[F_COUL_RECIP] += e_q;
434 enerd->term[F_LJ_RECIP] += e_lj;
435 enerd->dvdl_lin[efptCOUL] += dvdl_q;
436 enerd->dvdl_lin[efptVDW] += dvdl_lj;
440 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
442 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
445 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
446 gmx_int64_t step, real pforce, rvec *x, rvec *f)
450 char buf[STEPSTRSIZE];
453 for (i = 0; i < md->homenr; i++)
456 /* We also catch NAN, if the compiler does not optimize this away. */
457 if (fn2 >= pf2 || fn2 != fn2)
459 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
460 gmx_step_str(step, buf),
461 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
466 static void post_process_forces(t_commrec *cr,
468 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
470 matrix box, rvec x[],
475 t_forcerec *fr, gmx_vsite_t *vsite,
482 /* Spread the mesh force on virtual sites to the other particles...
483 * This is parallellized. MPI communication is performed
484 * if the constructing atoms aren't local.
486 wallcycle_start(wcycle, ewcVSITESPREAD);
487 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
488 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
490 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
491 wallcycle_stop(wcycle, ewcVSITESPREAD);
493 if (flags & GMX_FORCE_VIRIAL)
495 /* Now add the forces, this is local */
498 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
502 sum_forces(0, mdatoms->homenr,
505 if (EEL_FULL(fr->eeltype))
507 /* Add the mesh contribution to the virial */
508 m_add(vir_force, fr->vir_el_recip, vir_force);
510 if (EVDW_PME(fr->vdwtype))
512 /* Add the mesh contribution to the virial */
513 m_add(vir_force, fr->vir_lj_recip, vir_force);
517 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
522 if (fr->print_force >= 0)
524 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
528 static void do_nb_verlet(t_forcerec *fr,
529 interaction_const_t *ic,
530 gmx_enerdata_t *enerd,
531 int flags, int ilocality,
534 gmx_wallcycle_t wcycle)
536 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
538 nonbonded_verlet_group_t *nbvg;
541 if (!(flags & GMX_FORCE_NONBONDED))
543 /* skip non-bonded calculation */
547 nbvg = &fr->nbv->grp[ilocality];
549 /* CUDA kernel launch overhead is already timed separately */
550 if (fr->cutoff_scheme != ecutsVERLET)
552 gmx_incons("Invalid cut-off scheme passed!");
555 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
559 wallcycle_sub_start(wcycle, ewcsNONBONDED);
561 switch (nbvg->kernel_type)
563 case nbnxnk4x4_PlainC:
564 nbnxn_kernel_ref(&nbvg->nbl_lists,
570 enerd->grpp.ener[egCOULSR],
572 enerd->grpp.ener[egBHAMSR] :
573 enerd->grpp.ener[egLJSR]);
576 case nbnxnk4xN_SIMD_4xN:
577 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
584 enerd->grpp.ener[egCOULSR],
586 enerd->grpp.ener[egBHAMSR] :
587 enerd->grpp.ener[egLJSR]);
589 case nbnxnk4xN_SIMD_2xNN:
590 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
597 enerd->grpp.ener[egCOULSR],
599 enerd->grpp.ener[egBHAMSR] :
600 enerd->grpp.ener[egLJSR]);
603 case nbnxnk8x8x8_CUDA:
604 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
607 case nbnxnk8x8x8_PlainC:
608 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
613 nbvg->nbat->out[0].f,
615 enerd->grpp.ener[egCOULSR],
617 enerd->grpp.ener[egBHAMSR] :
618 enerd->grpp.ener[egLJSR]);
622 gmx_incons("Invalid nonbonded kernel type passed!");
627 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
630 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
632 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
634 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
635 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
637 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
641 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
643 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
644 if (flags & GMX_FORCE_ENERGY)
646 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
647 enr_nbnxn_kernel_ljc += 1;
648 enr_nbnxn_kernel_lj += 1;
651 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
652 nbvg->nbl_lists.natpair_ljq);
653 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
654 nbvg->nbl_lists.natpair_lj);
655 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
656 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
657 nbvg->nbl_lists.natpair_q);
659 if (ic->vdw_modifier == eintmodFORCESWITCH)
661 /* We add up the switch cost separately */
662 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
663 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
665 if (ic->vdw_modifier == eintmodPOTSWITCH)
667 /* We add up the switch cost separately */
668 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
669 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
671 if (ic->vdwtype == evdwPME)
673 /* We add up the LJ Ewald cost separately */
674 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
675 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
679 static void do_nb_verlet_fep(nbnxn_pairlist_set_t *nbl_lists,
686 gmx_enerdata_t *enerd,
689 gmx_wallcycle_t wcycle)
692 nb_kernel_data_t kernel_data;
694 real dvdl_nb[efptNR];
699 /* Add short-range interactions */
700 donb_flags |= GMX_NONBONDED_DO_SR;
702 /* Currently all group scheme kernels always calculate (shift-)forces */
703 if (flags & GMX_FORCE_FORCES)
705 donb_flags |= GMX_NONBONDED_DO_FORCE;
707 if (flags & GMX_FORCE_VIRIAL)
709 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
711 if (flags & GMX_FORCE_ENERGY)
713 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
715 if (flags & GMX_FORCE_DO_LR)
717 donb_flags |= GMX_NONBONDED_DO_LR;
720 kernel_data.flags = donb_flags;
721 kernel_data.lambda = lambda;
722 kernel_data.dvdl = dvdl_nb;
724 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
725 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
727 /* reset free energy components */
728 for (i = 0; i < efptNR; i++)
733 assert(gmx_omp_nthreads_get(emntNonbonded) == nbl_lists->nnbl);
735 wallcycle_sub_start(wcycle, ewcsNONBONDED);
736 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
737 for (th = 0; th < nbl_lists->nnbl; th++)
739 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
740 x, f, fr, mdatoms, &kernel_data, nrnb);
743 if (fepvals->sc_alpha != 0)
745 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
746 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
750 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
751 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
754 /* If we do foreign lambda and we have soft-core interactions
755 * we have to recalculate the (non-linear) energies contributions.
757 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
759 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
760 kernel_data.lambda = lam_i;
761 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
762 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
763 /* Note that we add to kernel_data.dvdl, but ignore the result */
765 for (i = 0; i < enerd->n_lambda; i++)
767 for (j = 0; j < efptNR; j++)
769 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
771 reset_foreign_enerdata(enerd);
772 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
773 for (th = 0; th < nbl_lists->nnbl; th++)
775 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
776 x, f, fr, mdatoms, &kernel_data, nrnb);
779 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
780 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
784 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
787 gmx_bool use_GPU(const nonbonded_verlet_t *nbv)
789 return nbv != NULL && nbv->bUseGPU;
792 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
793 t_inputrec *inputrec,
794 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
796 gmx_groups_t gmx_unused *groups,
797 matrix box, rvec x[], history_t *hist,
801 gmx_enerdata_t *enerd, t_fcdata *fcd,
802 real *lambda, t_graph *graph,
803 t_forcerec *fr, interaction_const_t *ic,
804 gmx_vsite_t *vsite, rvec mu_tot,
805 double t, FILE *field, gmx_edsam_t ed,
813 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
814 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
815 gmx_bool bDiffKernels = FALSE;
817 rvec vzero, box_diag;
819 float cycles_pme, cycles_force, cycles_wait_gpu;
820 nonbonded_verlet_t *nbv;
825 nb_kernel_type = fr->nbv->grp[0].kernel_type;
828 homenr = mdatoms->homenr;
830 clear_mat(vir_force);
833 if (DOMAINDECOMP(cr))
835 cg1 = cr->dd->ncg_tot;
846 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
847 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
848 bFillGrid = (bNS && bStateChanged);
849 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
850 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
851 bDoForces = (flags & GMX_FORCE_FORCES);
852 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
853 bUseGPU = fr->nbv->bUseGPU;
854 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
858 update_forcerec(fr, box);
860 if (NEED_MUTOT(*inputrec))
862 /* Calculate total (local) dipole moment in a temporary common array.
863 * This makes it possible to sum them over nodes faster.
865 calc_mu(start, homenr,
866 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
871 if (fr->ePBC != epbcNONE)
873 /* Compute shift vectors every step,
874 * because of pressure coupling or box deformation!
876 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
878 calc_shifts(box, fr->shift_vec);
883 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
884 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
886 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
888 unshift_self(graph, box, x);
892 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
893 fr->shift_vec, nbv->grp[0].nbat);
896 if (!(cr->duty & DUTY_PME))
898 /* Send particle coordinates to the pme nodes.
899 * Since this is only implemented for domain decomposition
900 * and domain decomposition does not use the graph,
901 * we do not need to worry about shifting.
906 wallcycle_start(wcycle, ewcPP_PMESENDX);
908 bBS = (inputrec->nwall == 2);
912 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
915 if (EEL_PME(fr->eeltype))
917 pme_flags |= GMX_PME_DO_COULOMB;
920 if (EVDW_PME(fr->vdwtype))
922 pme_flags |= GMX_PME_DO_LJ;
925 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
926 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
927 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
930 wallcycle_stop(wcycle, ewcPP_PMESENDX);
934 /* do gridding for pair search */
937 if (graph && bStateChanged)
939 /* Calculate intramolecular shift vectors to make molecules whole */
940 mk_mshift(fplog, graph, fr->ePBC, box, x);
944 box_diag[XX] = box[XX][XX];
945 box_diag[YY] = box[YY][YY];
946 box_diag[ZZ] = box[ZZ][ZZ];
948 wallcycle_start(wcycle, ewcNS);
951 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
952 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
954 0, mdatoms->homenr, -1, fr->cginfo, x,
956 nbv->grp[eintLocal].kernel_type,
957 nbv->grp[eintLocal].nbat);
958 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
962 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
963 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
965 nbv->grp[eintNonlocal].kernel_type,
966 nbv->grp[eintNonlocal].nbat);
967 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
970 if (nbv->ngrp == 1 ||
971 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
973 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
974 nbv->nbs, mdatoms, fr->cginfo);
978 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
979 nbv->nbs, mdatoms, fr->cginfo);
980 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
981 nbv->nbs, mdatoms, fr->cginfo);
983 wallcycle_stop(wcycle, ewcNS);
986 /* initialize the GPU atom data and copy shift vector */
991 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
992 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
993 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
996 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
997 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
998 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1001 /* do local pair search */
1004 wallcycle_start_nocount(wcycle, ewcNS);
1005 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1006 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
1009 nbv->min_ci_balanced,
1010 &nbv->grp[eintLocal].nbl_lists,
1012 nbv->grp[eintLocal].kernel_type,
1014 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1018 /* initialize local pair-list on the GPU */
1019 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1020 nbv->grp[eintLocal].nbl_lists.nbl[0],
1023 wallcycle_stop(wcycle, ewcNS);
1027 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1028 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1029 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
1030 nbv->grp[eintLocal].nbat);
1031 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1032 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1037 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1038 /* launch local nonbonded F on GPU */
1039 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
1041 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1044 /* Communicate coordinates and sum dipole if necessary +
1045 do non-local pair search */
1046 if (DOMAINDECOMP(cr))
1048 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
1049 nbv->grp[eintLocal].kernel_type);
1053 /* With GPU+CPU non-bonded calculations we need to copy
1054 * the local coordinates to the non-local nbat struct
1055 * (in CPU format) as the non-local kernel call also
1056 * calculates the local - non-local interactions.
1058 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1059 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1060 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
1061 nbv->grp[eintNonlocal].nbat);
1062 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1063 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1068 wallcycle_start_nocount(wcycle, ewcNS);
1069 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1073 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
1076 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
1079 nbv->min_ci_balanced,
1080 &nbv->grp[eintNonlocal].nbl_lists,
1082 nbv->grp[eintNonlocal].kernel_type,
1085 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1087 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1089 /* initialize non-local pair-list on the GPU */
1090 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1091 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1094 wallcycle_stop(wcycle, ewcNS);
1098 wallcycle_start(wcycle, ewcMOVEX);
1099 dd_move_x(cr->dd, box, x);
1101 /* When we don't need the total dipole we sum it in global_stat */
1102 if (bStateChanged && NEED_MUTOT(*inputrec))
1104 gmx_sumd(2*DIM, mu, cr);
1106 wallcycle_stop(wcycle, ewcMOVEX);
1108 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1109 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1110 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1111 nbv->grp[eintNonlocal].nbat);
1112 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1113 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1116 if (bUseGPU && !bDiffKernels)
1118 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1119 /* launch non-local nonbonded F on GPU */
1120 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1122 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1128 /* launch D2H copy-back F */
1129 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1130 if (DOMAINDECOMP(cr) && !bDiffKernels)
1132 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1133 flags, eatNonlocal);
1135 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1137 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1140 if (bStateChanged && NEED_MUTOT(*inputrec))
1144 gmx_sumd(2*DIM, mu, cr);
1147 for (i = 0; i < 2; i++)
1149 for (j = 0; j < DIM; j++)
1151 fr->mu_tot[i][j] = mu[i*DIM + j];
1155 if (fr->efep == efepNO)
1157 copy_rvec(fr->mu_tot[0], mu_tot);
1161 for (j = 0; j < DIM; j++)
1164 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1165 lambda[efptCOUL]*fr->mu_tot[1][j];
1169 /* Reset energies */
1170 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1171 clear_rvecs(SHIFTS, fr->fshift);
1173 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1175 wallcycle_start(wcycle, ewcPPDURINGPME);
1176 dd_force_flop_start(cr->dd, nrnb);
1181 /* Enforced rotation has its own cycle counter that starts after the collective
1182 * coordinates have been communicated. It is added to ddCyclF to allow
1183 * for proper load-balancing */
1184 wallcycle_start(wcycle, ewcROT);
1185 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1186 wallcycle_stop(wcycle, ewcROT);
1189 /* Start the force cycle counter.
1190 * This counter is stopped in do_forcelow_level.
1191 * No parallel communication should occur while this counter is running,
1192 * since that will interfere with the dynamic load balancing.
1194 wallcycle_start(wcycle, ewcFORCE);
1197 /* Reset forces for which the virial is calculated separately:
1198 * PME/Ewald forces if necessary */
1199 if (fr->bF_NoVirSum)
1201 if (flags & GMX_FORCE_VIRIAL)
1203 fr->f_novirsum = fr->f_novirsum_alloc;
1206 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1210 clear_rvecs(homenr, fr->f_novirsum+start);
1215 /* We are not calculating the pressure so we do not need
1216 * a separate array for forces that do not contribute
1223 /* Clear the short- and long-range forces */
1224 clear_rvecs(fr->natoms_force_constr, f);
1225 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1227 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1230 clear_rvec(fr->vir_diag_posres);
1233 if (inputrec->ePull == epullCONSTRAINT)
1235 clear_pull_forces(inputrec->pull);
1238 /* We calculate the non-bonded forces, when done on the CPU, here.
1239 * We do this before calling do_force_lowlevel, because in that
1240 * function, the listed forces are calculated before PME, which
1241 * does communication. With this order, non-bonded and listed
1242 * force calculation imbalance can be balanced out by the domain
1243 * decomposition load balancing.
1248 /* Maybe we should move this into do_force_lowlevel */
1249 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1253 if (fr->efep != efepNO)
1255 /* Calculate the local and non-local free energy interactions here.
1256 * Happens here on the CPU both with and without GPU.
1258 if (fr->nbv->grp[eintLocal].nbl_lists.nbl_fep[0]->nrj > 0)
1260 do_nb_verlet_fep(&fr->nbv->grp[eintLocal].nbl_lists,
1262 inputrec->fepvals, lambda,
1263 enerd, flags, nrnb, wcycle);
1266 if (DOMAINDECOMP(cr) &&
1267 fr->nbv->grp[eintNonlocal].nbl_lists.nbl_fep[0]->nrj > 0)
1269 do_nb_verlet_fep(&fr->nbv->grp[eintNonlocal].nbl_lists,
1271 inputrec->fepvals, lambda,
1272 enerd, flags, nrnb, wcycle);
1276 if (!bUseOrEmulGPU || bDiffKernels)
1280 if (DOMAINDECOMP(cr))
1282 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1283 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1293 aloc = eintNonlocal;
1296 /* Add all the non-bonded force to the normal force array.
1297 * This can be split into a local a non-local part when overlapping
1298 * communication with calculation with domain decomposition.
1300 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1301 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1302 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1303 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1304 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1305 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1306 wallcycle_start_nocount(wcycle, ewcFORCE);
1308 /* if there are multiple fshift output buffers reduce them */
1309 if ((flags & GMX_FORCE_VIRIAL) &&
1310 nbv->grp[aloc].nbl_lists.nnbl > 1)
1312 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1317 /* update QMMMrec, if necessary */
1320 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1323 if ((flags & GMX_FORCE_LISTED) && top->idef.il[F_POSRES].nr > 0)
1325 posres_wrapper(flags, inputrec, nrnb, top, box, x,
1329 if ((flags & GMX_FORCE_LISTED) && top->idef.il[F_FBPOSRES].nr > 0)
1331 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1334 /* Compute the bonded and non-bonded energies and optionally forces */
1335 do_force_lowlevel(fr, inputrec, &(top->idef),
1336 cr, nrnb, wcycle, mdatoms,
1337 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1339 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1340 flags, &cycles_pme);
1344 if (do_per_step(step, inputrec->nstcalclr))
1346 /* Add the long range forces to the short range forces */
1347 for (i = 0; i < fr->natoms_force_constr; i++)
1349 rvec_add(fr->f_twin[i], f[i], f[i]);
1354 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1358 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1361 if (bUseOrEmulGPU && !bDiffKernels)
1363 /* wait for non-local forces (or calculate in emulation mode) */
1364 if (DOMAINDECOMP(cr))
1370 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1371 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1372 nbv->grp[eintNonlocal].nbat,
1374 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1376 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1377 cycles_wait_gpu += cycles_tmp;
1378 cycles_force += cycles_tmp;
1382 wallcycle_start_nocount(wcycle, ewcFORCE);
1383 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1385 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1387 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1388 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1389 /* skip the reduction if there was no non-local work to do */
1390 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1392 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1393 nbv->grp[eintNonlocal].nbat, f);
1395 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1396 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1400 if (bDoForces && DOMAINDECOMP(cr))
1404 /* We are done with the CPU compute, but the GPU local non-bonded
1405 * kernel can still be running while we communicate the forces.
1406 * We start a counter here, so we can, hopefully, time the rest
1407 * of the GPU kernel execution and data transfer.
1409 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L_EST);
1412 /* Communicate the forces */
1413 wallcycle_start(wcycle, ewcMOVEF);
1414 dd_move_f(cr->dd, f, fr->fshift);
1415 /* Do we need to communicate the separate force array
1416 * for terms that do not contribute to the single sum virial?
1417 * Position restraints and electric fields do not introduce
1418 * inter-cg forces, only full electrostatics methods do.
1419 * When we do not calculate the virial, fr->f_novirsum = f,
1420 * so we have already communicated these forces.
1422 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1423 (flags & GMX_FORCE_VIRIAL))
1425 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1429 /* We should not update the shift forces here,
1430 * since f_twin is already included in f.
1432 dd_move_f(cr->dd, fr->f_twin, NULL);
1434 wallcycle_stop(wcycle, ewcMOVEF);
1439 /* wait for local forces (or calculate in emulation mode) */
1442 float cycles_tmp, cycles_wait_est;
1443 const float cuda_api_overhead_margin = 50000.0f; /* cycles */
1445 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1446 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1447 nbv->grp[eintLocal].nbat,
1449 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1451 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1453 if (bDoForces && DOMAINDECOMP(cr))
1455 cycles_wait_est = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L_EST);
1457 if (cycles_tmp < cuda_api_overhead_margin)
1459 /* We measured few cycles, it could be that the kernel
1460 * and transfer finished earlier and there was no actual
1461 * wait time, only API call overhead.
1462 * Then the actual time could be anywhere between 0 and
1463 * cycles_wait_est. As a compromise, we use half the time.
1465 cycles_wait_est *= 0.5f;
1470 /* No force communication so we actually timed the wait */
1471 cycles_wait_est = cycles_tmp;
1473 /* Even though this is after dd_move_f, the actual task we are
1474 * waiting for runs asynchronously with dd_move_f and we usually
1475 * have nothing to balance it with, so we can and should add
1476 * the time to the force time for load balancing.
1478 cycles_force += cycles_wait_est;
1479 cycles_wait_gpu += cycles_wait_est;
1481 /* now clear the GPU outputs while we finish the step on the CPU */
1483 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1484 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1485 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1489 wallcycle_start_nocount(wcycle, ewcFORCE);
1490 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1491 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1493 wallcycle_stop(wcycle, ewcFORCE);
1495 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1496 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1497 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1499 /* skip the reduction if there was no non-local work to do */
1500 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1501 nbv->grp[eintLocal].nbat, f);
1503 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1504 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1507 if (DOMAINDECOMP(cr))
1509 dd_force_flop_stop(cr->dd, nrnb);
1512 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1515 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1522 if (IR_ELEC_FIELD(*inputrec))
1524 /* Compute forces due to electric field */
1525 calc_f_el(MASTER(cr) ? field : NULL,
1526 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1527 inputrec->ex, inputrec->et, t);
1530 /* If we have NoVirSum forces, but we do not calculate the virial,
1531 * we sum fr->f_novirum=f later.
1533 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1535 wallcycle_start(wcycle, ewcVSITESPREAD);
1536 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1537 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1538 wallcycle_stop(wcycle, ewcVSITESPREAD);
1542 wallcycle_start(wcycle, ewcVSITESPREAD);
1543 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1545 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1546 wallcycle_stop(wcycle, ewcVSITESPREAD);
1550 if (flags & GMX_FORCE_VIRIAL)
1552 /* Calculation of the virial must be done after vsites! */
1553 calc_virial(0, mdatoms->homenr, x, f,
1554 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1558 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1560 /* Since the COM pulling is always done mass-weighted, no forces are
1561 * applied to vsites and this call can be done after vsite spreading.
1563 pull_potential_wrapper(cr, inputrec, box, x,
1564 f, vir_force, mdatoms, enerd, lambda, t,
1568 /* Add the forces from enforced rotation potentials (if any) */
1571 wallcycle_start(wcycle, ewcROTadd);
1572 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1573 wallcycle_stop(wcycle, ewcROTadd);
1576 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1577 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1579 if (PAR(cr) && !(cr->duty & DUTY_PME))
1581 /* In case of node-splitting, the PP nodes receive the long-range
1582 * forces, virial and energy from the PME nodes here.
1584 pme_receive_force_ener(cr, wcycle, enerd, fr);
1589 post_process_forces(cr, step, nrnb, wcycle,
1590 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1594 /* Sum the potential energy terms from group contributions */
1595 sum_epot(&(enerd->grpp), enerd->term);
1598 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1599 t_inputrec *inputrec,
1600 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1601 gmx_localtop_t *top,
1602 gmx_groups_t *groups,
1603 matrix box, rvec x[], history_t *hist,
1607 gmx_enerdata_t *enerd, t_fcdata *fcd,
1608 real *lambda, t_graph *graph,
1609 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1610 double t, FILE *field, gmx_edsam_t ed,
1611 gmx_bool bBornRadii,
1617 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1618 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1619 gmx_bool bDoAdressWF;
1621 rvec vzero, box_diag;
1622 real e, v, dvdlambda[efptNR];
1624 float cycles_pme, cycles_force;
1627 homenr = mdatoms->homenr;
1629 clear_mat(vir_force);
1632 if (DOMAINDECOMP(cr))
1634 cg1 = cr->dd->ncg_tot;
1645 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1646 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1647 /* Should we update the long-range neighborlists at this step? */
1648 bDoLongRangeNS = fr->bTwinRange && bNS;
1649 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1650 bFillGrid = (bNS && bStateChanged);
1651 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1652 bDoForces = (flags & GMX_FORCE_FORCES);
1653 bDoPotential = (flags & GMX_FORCE_ENERGY);
1654 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1655 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1657 /* should probably move this to the forcerec since it doesn't change */
1658 bDoAdressWF = ((fr->adress_type != eAdressOff));
1662 update_forcerec(fr, box);
1664 if (NEED_MUTOT(*inputrec))
1666 /* Calculate total (local) dipole moment in a temporary common array.
1667 * This makes it possible to sum them over nodes faster.
1669 calc_mu(start, homenr,
1670 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1675 if (fr->ePBC != epbcNONE)
1677 /* Compute shift vectors every step,
1678 * because of pressure coupling or box deformation!
1680 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1682 calc_shifts(box, fr->shift_vec);
1687 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1688 &(top->cgs), x, fr->cg_cm);
1689 inc_nrnb(nrnb, eNR_CGCM, homenr);
1690 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1692 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1694 unshift_self(graph, box, x);
1699 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1700 inc_nrnb(nrnb, eNR_CGCM, homenr);
1703 if (bCalcCGCM && gmx_debug_at)
1705 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1709 if (!(cr->duty & DUTY_PME))
1711 /* Send particle coordinates to the pme nodes.
1712 * Since this is only implemented for domain decomposition
1713 * and domain decomposition does not use the graph,
1714 * we do not need to worry about shifting.
1719 wallcycle_start(wcycle, ewcPP_PMESENDX);
1721 bBS = (inputrec->nwall == 2);
1724 copy_mat(box, boxs);
1725 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1728 if (EEL_PME(fr->eeltype))
1730 pme_flags |= GMX_PME_DO_COULOMB;
1733 if (EVDW_PME(fr->vdwtype))
1735 pme_flags |= GMX_PME_DO_LJ;
1738 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1739 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1740 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1743 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1745 #endif /* GMX_MPI */
1747 /* Communicate coordinates and sum dipole if necessary */
1748 if (DOMAINDECOMP(cr))
1750 wallcycle_start(wcycle, ewcMOVEX);
1751 dd_move_x(cr->dd, box, x);
1752 wallcycle_stop(wcycle, ewcMOVEX);
1755 /* update adress weight beforehand */
1756 if (bStateChanged && bDoAdressWF)
1758 /* need pbc for adress weight calculation with pbc_dx */
1759 set_pbc(&pbc, inputrec->ePBC, box);
1760 if (fr->adress_site == eAdressSITEcog)
1762 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1763 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1765 else if (fr->adress_site == eAdressSITEcom)
1767 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1768 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1770 else if (fr->adress_site == eAdressSITEatomatom)
1772 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1773 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1777 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1778 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1782 if (NEED_MUTOT(*inputrec))
1789 gmx_sumd(2*DIM, mu, cr);
1791 for (i = 0; i < 2; i++)
1793 for (j = 0; j < DIM; j++)
1795 fr->mu_tot[i][j] = mu[i*DIM + j];
1799 if (fr->efep == efepNO)
1801 copy_rvec(fr->mu_tot[0], mu_tot);
1805 for (j = 0; j < DIM; j++)
1808 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1813 /* Reset energies */
1814 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1815 clear_rvecs(SHIFTS, fr->fshift);
1819 wallcycle_start(wcycle, ewcNS);
1821 if (graph && bStateChanged)
1823 /* Calculate intramolecular shift vectors to make molecules whole */
1824 mk_mshift(fplog, graph, fr->ePBC, box, x);
1827 /* Do the actual neighbour searching */
1829 groups, top, mdatoms,
1830 cr, nrnb, bFillGrid,
1833 wallcycle_stop(wcycle, ewcNS);
1836 if (inputrec->implicit_solvent && bNS)
1838 make_gb_nblist(cr, inputrec->gb_algorithm,
1839 x, box, fr, &top->idef, graph, fr->born);
1842 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1844 wallcycle_start(wcycle, ewcPPDURINGPME);
1845 dd_force_flop_start(cr->dd, nrnb);
1850 /* Enforced rotation has its own cycle counter that starts after the collective
1851 * coordinates have been communicated. It is added to ddCyclF to allow
1852 * for proper load-balancing */
1853 wallcycle_start(wcycle, ewcROT);
1854 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1855 wallcycle_stop(wcycle, ewcROT);
1858 /* Start the force cycle counter.
1859 * This counter is stopped in do_forcelow_level.
1860 * No parallel communication should occur while this counter is running,
1861 * since that will interfere with the dynamic load balancing.
1863 wallcycle_start(wcycle, ewcFORCE);
1867 /* Reset forces for which the virial is calculated separately:
1868 * PME/Ewald forces if necessary */
1869 if (fr->bF_NoVirSum)
1871 if (flags & GMX_FORCE_VIRIAL)
1873 fr->f_novirsum = fr->f_novirsum_alloc;
1876 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1880 clear_rvecs(homenr, fr->f_novirsum+start);
1885 /* We are not calculating the pressure so we do not need
1886 * a separate array for forces that do not contribute
1893 /* Clear the short- and long-range forces */
1894 clear_rvecs(fr->natoms_force_constr, f);
1895 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1897 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1900 clear_rvec(fr->vir_diag_posres);
1902 if (inputrec->ePull == epullCONSTRAINT)
1904 clear_pull_forces(inputrec->pull);
1907 /* update QMMMrec, if necessary */
1910 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1913 if ((flags & GMX_FORCE_LISTED) && top->idef.il[F_POSRES].nr > 0)
1915 posres_wrapper(flags, inputrec, nrnb, top, box, x,
1919 if ((flags & GMX_FORCE_LISTED) && top->idef.il[F_FBPOSRES].nr > 0)
1921 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1924 /* Compute the bonded and non-bonded energies and optionally forces */
1925 do_force_lowlevel(fr, inputrec, &(top->idef),
1926 cr, nrnb, wcycle, mdatoms,
1927 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1929 inputrec->fepvals, lambda,
1930 graph, &(top->excls), fr->mu_tot,
1936 if (do_per_step(step, inputrec->nstcalclr))
1938 /* Add the long range forces to the short range forces */
1939 for (i = 0; i < fr->natoms_force_constr; i++)
1941 rvec_add(fr->f_twin[i], f[i], f[i]);
1946 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1950 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1953 if (DOMAINDECOMP(cr))
1955 dd_force_flop_stop(cr->dd, nrnb);
1958 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1964 if (IR_ELEC_FIELD(*inputrec))
1966 /* Compute forces due to electric field */
1967 calc_f_el(MASTER(cr) ? field : NULL,
1968 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1969 inputrec->ex, inputrec->et, t);
1972 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1974 /* Compute thermodynamic force in hybrid AdResS region */
1975 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1976 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1979 /* Communicate the forces */
1980 if (DOMAINDECOMP(cr))
1982 wallcycle_start(wcycle, ewcMOVEF);
1983 dd_move_f(cr->dd, f, fr->fshift);
1984 /* Do we need to communicate the separate force array
1985 * for terms that do not contribute to the single sum virial?
1986 * Position restraints and electric fields do not introduce
1987 * inter-cg forces, only full electrostatics methods do.
1988 * When we do not calculate the virial, fr->f_novirsum = f,
1989 * so we have already communicated these forces.
1991 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1992 (flags & GMX_FORCE_VIRIAL))
1994 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1998 /* We should not update the shift forces here,
1999 * since f_twin is already included in f.
2001 dd_move_f(cr->dd, fr->f_twin, NULL);
2003 wallcycle_stop(wcycle, ewcMOVEF);
2006 /* If we have NoVirSum forces, but we do not calculate the virial,
2007 * we sum fr->f_novirum=f later.
2009 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
2011 wallcycle_start(wcycle, ewcVSITESPREAD);
2012 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
2013 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
2014 wallcycle_stop(wcycle, ewcVSITESPREAD);
2018 wallcycle_start(wcycle, ewcVSITESPREAD);
2019 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
2021 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
2022 wallcycle_stop(wcycle, ewcVSITESPREAD);
2026 if (flags & GMX_FORCE_VIRIAL)
2028 /* Calculation of the virial must be done after vsites! */
2029 calc_virial(0, mdatoms->homenr, x, f,
2030 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
2034 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
2036 pull_potential_wrapper(cr, inputrec, box, x,
2037 f, vir_force, mdatoms, enerd, lambda, t,
2041 /* Add the forces from enforced rotation potentials (if any) */
2044 wallcycle_start(wcycle, ewcROTadd);
2045 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
2046 wallcycle_stop(wcycle, ewcROTadd);
2049 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
2050 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
2052 if (PAR(cr) && !(cr->duty & DUTY_PME))
2054 /* In case of node-splitting, the PP nodes receive the long-range
2055 * forces, virial and energy from the PME nodes here.
2057 pme_receive_force_ener(cr, wcycle, enerd, fr);
2062 post_process_forces(cr, step, nrnb, wcycle,
2063 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
2067 /* Sum the potential energy terms from group contributions */
2068 sum_epot(&(enerd->grpp), enerd->term);
2071 void do_force(FILE *fplog, t_commrec *cr,
2072 t_inputrec *inputrec,
2073 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2074 gmx_localtop_t *top,
2075 gmx_groups_t *groups,
2076 matrix box, rvec x[], history_t *hist,
2080 gmx_enerdata_t *enerd, t_fcdata *fcd,
2081 real *lambda, t_graph *graph,
2083 gmx_vsite_t *vsite, rvec mu_tot,
2084 double t, FILE *field, gmx_edsam_t ed,
2085 gmx_bool bBornRadii,
2088 /* modify force flag if not doing nonbonded */
2089 if (!fr->bNonbonded)
2091 flags &= ~GMX_FORCE_NONBONDED;
2094 switch (inputrec->cutoff_scheme)
2097 do_force_cutsVERLET(fplog, cr, inputrec,
2113 do_force_cutsGROUP(fplog, cr, inputrec,
2128 gmx_incons("Invalid cut-off scheme passed!");
2133 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2134 t_inputrec *ir, t_mdatoms *md,
2135 t_state *state, t_commrec *cr, t_nrnb *nrnb,
2136 t_forcerec *fr, gmx_localtop_t *top)
2138 int i, m, start, end;
2140 real dt = ir->delta_t;
2144 snew(savex, state->natoms);
2151 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2152 start, md->homenr, end);
2154 /* Do a first constrain to reset particles... */
2155 step = ir->init_step;
2158 char buf[STEPSTRSIZE];
2159 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2160 gmx_step_str(step, buf));
2164 /* constrain the current position */
2165 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2166 ir, NULL, cr, step, 0, 1.0, md,
2167 state->x, state->x, NULL,
2168 fr->bMolPBC, state->box,
2169 state->lambda[efptBONDED], &dvdl_dum,
2170 NULL, NULL, nrnb, econqCoord,
2171 ir->epc == epcMTTK, state->veta, state->veta);
2174 /* constrain the inital velocity, and save it */
2175 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2176 /* might not yet treat veta correctly */
2177 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2178 ir, NULL, cr, step, 0, 1.0, md,
2179 state->x, state->v, state->v,
2180 fr->bMolPBC, state->box,
2181 state->lambda[efptBONDED], &dvdl_dum,
2182 NULL, NULL, nrnb, econqVeloc,
2183 ir->epc == epcMTTK, state->veta, state->veta);
2185 /* constrain the inital velocities at t-dt/2 */
2186 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2188 for (i = start; (i < end); i++)
2190 for (m = 0; (m < DIM); m++)
2192 /* Reverse the velocity */
2193 state->v[i][m] = -state->v[i][m];
2194 /* Store the position at t-dt in buf */
2195 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2198 /* Shake the positions at t=-dt with the positions at t=0
2199 * as reference coordinates.
2203 char buf[STEPSTRSIZE];
2204 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2205 gmx_step_str(step, buf));
2208 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2209 ir, NULL, cr, step, -1, 1.0, md,
2210 state->x, savex, NULL,
2211 fr->bMolPBC, state->box,
2212 state->lambda[efptBONDED], &dvdl_dum,
2213 state->v, NULL, nrnb, econqCoord,
2214 ir->epc == epcMTTK, state->veta, state->veta);
2216 for (i = start; i < end; i++)
2218 for (m = 0; m < DIM; m++)
2220 /* Re-reverse the velocities */
2221 state->v[i][m] = -state->v[i][m];
2230 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2231 double *enerout, double *virout)
2233 double enersum, virsum;
2234 double invscale, invscale2, invscale3;
2235 double r, ea, eb, ec, pa, pb, pc, pd;
2237 int ri, offset, tabfactor;
2239 invscale = 1.0/scale;
2240 invscale2 = invscale*invscale;
2241 invscale3 = invscale*invscale2;
2243 /* Following summation derived from cubic spline definition,
2244 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2245 * the cubic spline. We first calculate the negative of the
2246 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2247 * add the more standard, abrupt cutoff correction to that result,
2248 * yielding the long-range correction for a switched function. We
2249 * perform both the pressure and energy loops at the same time for
2250 * simplicity, as the computational cost is low. */
2254 /* Since the dispersion table has been scaled down a factor
2255 * 6.0 and the repulsion a factor 12.0 to compensate for the
2256 * c6/c12 parameters inside nbfp[] being scaled up (to save
2257 * flops in kernels), we need to correct for this.
2268 for (ri = rstart; ri < rend; ++ri)
2272 eb = 2.0*invscale2*r;
2276 pb = 3.0*invscale2*r;
2277 pc = 3.0*invscale*r*r;
2280 /* this "8" is from the packing in the vdwtab array - perhaps
2281 should be defined? */
2283 offset = 8*ri + offstart;
2284 y0 = vdwtab[offset];
2285 f = vdwtab[offset+1];
2286 g = vdwtab[offset+2];
2287 h = vdwtab[offset+3];
2289 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);
2290 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);
2292 *enerout = 4.0*M_PI*enersum*tabfactor;
2293 *virout = 4.0*M_PI*virsum*tabfactor;
2296 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2298 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2299 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2300 double invscale, invscale2, invscale3;
2301 int ri0, ri1, ri, i, offstart, offset;
2302 real scale, *vdwtab, tabfactor, tmp;
2304 fr->enershiftsix = 0;
2305 fr->enershifttwelve = 0;
2306 fr->enerdiffsix = 0;
2307 fr->enerdifftwelve = 0;
2309 fr->virdifftwelve = 0;
2311 if (eDispCorr != edispcNO)
2313 for (i = 0; i < 2; i++)
2318 if ((fr->vdw_modifier == eintmodPOTSHIFT) ||
2319 (fr->vdw_modifier == eintmodPOTSWITCH) ||
2320 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2321 (fr->vdwtype == evdwSHIFT) ||
2322 (fr->vdwtype == evdwSWITCH))
2324 if (((fr->vdw_modifier == eintmodPOTSWITCH) ||
2325 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2326 (fr->vdwtype == evdwSWITCH)) && fr->rvdw_switch == 0)
2329 "With dispersion correction rvdw-switch can not be zero "
2330 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2333 scale = fr->nblists[0].table_vdw.scale;
2334 vdwtab = fr->nblists[0].table_vdw.data;
2336 /* Round the cut-offs to exact table values for precision */
2337 ri0 = floor(fr->rvdw_switch*scale);
2338 ri1 = ceil(fr->rvdw*scale);
2340 /* The code below has some support for handling force-switching, i.e.
2341 * when the force (instead of potential) is switched over a limited
2342 * region. This leads to a constant shift in the potential inside the
2343 * switching region, which we can handle by adding a constant energy
2344 * term in the force-switch case just like when we do potential-shift.
2346 * For now this is not enabled, but to keep the functionality in the
2347 * code we check separately for switch and shift. When we do force-switch
2348 * the shifting point is rvdw_switch, while it is the cutoff when we
2349 * have a classical potential-shift.
2351 * For a pure potential-shift the potential has a constant shift
2352 * all the way out to the cutoff, and that is it. For other forms
2353 * we need to calculate the constant shift up to the point where we
2354 * start modifying the potential.
2356 ri0 = (fr->vdw_modifier == eintmodPOTSHIFT) ? ri1 : ri0;
2363 if ((fr->vdw_modifier == eintmodFORCESWITCH) ||
2364 (fr->vdwtype == evdwSHIFT))
2366 /* Determine the constant energy shift below rvdw_switch.
2367 * Table has a scale factor since we have scaled it down to compensate
2368 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2370 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2371 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2373 else if (fr->vdw_modifier == eintmodPOTSHIFT)
2375 fr->enershiftsix = (real)(-1.0/(rc3*rc3));
2376 fr->enershifttwelve = (real)( 1.0/(rc9*rc3));
2379 /* Add the constant part from 0 to rvdw_switch.
2380 * This integration from 0 to rvdw_switch overcounts the number
2381 * of interactions by 1, as it also counts the self interaction.
2382 * We will correct for this later.
2384 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2385 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2387 /* Calculate the contribution in the range [r0,r1] where we
2388 * modify the potential. For a pure potential-shift modifier we will
2389 * have ri0==ri1, and there will not be any contribution here.
2391 for (i = 0; i < 2; i++)
2395 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2396 eners[i] -= enersum;
2400 /* Alright: Above we compensated by REMOVING the parts outside r0
2401 * corresponding to the ideal VdW 1/r6 and /r12 potentials.
2403 * Regardless of whether r0 is the point where we start switching,
2404 * or the cutoff where we calculated the constant shift, we include
2405 * all the parts we are missing out to infinity from r0 by
2406 * calculating the analytical dispersion correction.
2408 eners[0] += -4.0*M_PI/(3.0*rc3);
2409 eners[1] += 4.0*M_PI/(9.0*rc9);
2410 virs[0] += 8.0*M_PI/rc3;
2411 virs[1] += -16.0*M_PI/(3.0*rc9);
2413 else if (fr->vdwtype == evdwCUT ||
2414 EVDW_PME(fr->vdwtype) ||
2415 fr->vdwtype == evdwUSER)
2417 if (fr->vdwtype == evdwUSER && fplog)
2420 "WARNING: using dispersion correction with user tables\n");
2423 /* Note that with LJ-PME, the dispersion correction is multiplied
2424 * by the difference between the actual C6 and the value of C6
2425 * that would produce the combination rule.
2426 * This means the normal energy and virial difference formulas
2430 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2432 /* Contribution beyond the cut-off */
2433 eners[0] += -4.0*M_PI/(3.0*rc3);
2434 eners[1] += 4.0*M_PI/(9.0*rc9);
2435 if (fr->vdw_modifier == eintmodPOTSHIFT)
2437 /* Contribution within the cut-off */
2438 eners[0] += -4.0*M_PI/(3.0*rc3);
2439 eners[1] += 4.0*M_PI/(3.0*rc9);
2441 /* Contribution beyond the cut-off */
2442 virs[0] += 8.0*M_PI/rc3;
2443 virs[1] += -16.0*M_PI/(3.0*rc9);
2448 "Dispersion correction is not implemented for vdw-type = %s",
2449 evdw_names[fr->vdwtype]);
2452 /* When we deprecate the group kernels the code below can go too */
2453 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2455 /* Calculate self-interaction coefficient (assuming that
2456 * the reciprocal-space contribution is constant in the
2457 * region that contributes to the self-interaction).
2459 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2461 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2462 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2465 fr->enerdiffsix = eners[0];
2466 fr->enerdifftwelve = eners[1];
2467 /* The 0.5 is due to the Gromacs definition of the virial */
2468 fr->virdiffsix = 0.5*virs[0];
2469 fr->virdifftwelve = 0.5*virs[1];
2473 void calc_dispcorr(t_inputrec *ir, t_forcerec *fr,
2475 matrix box, real lambda, tensor pres, tensor virial,
2476 real *prescorr, real *enercorr, real *dvdlcorr)
2478 gmx_bool bCorrAll, bCorrPres;
2479 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2489 if (ir->eDispCorr != edispcNO)
2491 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2492 ir->eDispCorr == edispcAllEnerPres);
2493 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2494 ir->eDispCorr == edispcAllEnerPres);
2496 invvol = 1/det(box);
2499 /* Only correct for the interactions with the inserted molecule */
2500 dens = (natoms - fr->n_tpi)*invvol;
2505 dens = natoms*invvol;
2506 ninter = 0.5*natoms;
2509 if (ir->efep == efepNO)
2511 avcsix = fr->avcsix[0];
2512 avctwelve = fr->avctwelve[0];
2516 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2517 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2520 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2521 *enercorr += avcsix*enerdiff;
2523 if (ir->efep != efepNO)
2525 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2529 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2530 *enercorr += avctwelve*enerdiff;
2531 if (fr->efep != efepNO)
2533 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2539 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2540 if (ir->eDispCorr == edispcAllEnerPres)
2542 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2544 /* The factor 2 is because of the Gromacs virial definition */
2545 spres = -2.0*invvol*svir*PRESFAC;
2547 for (m = 0; m < DIM; m++)
2549 virial[m][m] += svir;
2550 pres[m][m] += spres;
2555 /* Can't currently control when it prints, for now, just print when degugging */
2560 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2566 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2567 *enercorr, spres, svir);
2571 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2575 if (fr->efep != efepNO)
2577 *dvdlcorr += dvdlambda;
2582 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2583 t_graph *graph, rvec x[])
2587 fprintf(fplog, "Removing pbc first time\n");
2589 calc_shifts(box, fr->shift_vec);
2592 mk_mshift(fplog, graph, fr->ePBC, box, x);
2595 p_graph(debug, "do_pbc_first 1", graph);
2597 shift_self(graph, box, x);
2598 /* By doing an extra mk_mshift the molecules that are broken
2599 * because they were e.g. imported from another software
2600 * will be made whole again. Such are the healing powers
2603 mk_mshift(fplog, graph, fr->ePBC, box, x);
2606 p_graph(debug, "do_pbc_first 2", graph);
2611 fprintf(fplog, "Done rmpbc\n");
2615 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2616 gmx_mtop_t *mtop, rvec x[],
2621 gmx_molblock_t *molb;
2623 if (bFirst && fplog)
2625 fprintf(fplog, "Removing pbc first time\n");
2630 for (mb = 0; mb < mtop->nmolblock; mb++)
2632 molb = &mtop->molblock[mb];
2633 if (molb->natoms_mol == 1 ||
2634 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2636 /* Just one atom or charge group in the molecule, no PBC required */
2637 as += molb->nmol*molb->natoms_mol;
2641 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2642 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2643 0, molb->natoms_mol, FALSE, FALSE, graph);
2645 for (mol = 0; mol < molb->nmol; mol++)
2647 mk_mshift(fplog, graph, ePBC, box, x+as);
2649 shift_self(graph, box, x+as);
2650 /* The molecule is whole now.
2651 * We don't need the second mk_mshift call as in do_pbc_first,
2652 * since we no longer need this graph.
2655 as += molb->natoms_mol;
2663 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2664 gmx_mtop_t *mtop, rvec x[])
2666 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2669 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2670 gmx_mtop_t *mtop, rvec x[])
2672 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2675 void finish_run(FILE *fplog, t_commrec *cr,
2676 t_inputrec *inputrec,
2677 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2678 gmx_walltime_accounting_t walltime_accounting,
2679 nonbonded_verlet_t *nbv,
2680 gmx_bool bWriteStat)
2683 t_nrnb *nrnb_tot = NULL;
2686 double elapsed_time,
2687 elapsed_time_over_all_ranks,
2688 elapsed_time_over_all_threads,
2689 elapsed_time_over_all_threads_over_all_ranks;
2690 wallcycle_sum(cr, wcycle);
2696 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2697 cr->mpi_comm_mysim);
2705 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2706 elapsed_time_over_all_ranks = elapsed_time;
2707 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2708 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2712 /* reduce elapsed_time over all MPI ranks in the current simulation */
2713 MPI_Allreduce(&elapsed_time,
2714 &elapsed_time_over_all_ranks,
2715 1, MPI_DOUBLE, MPI_SUM,
2716 cr->mpi_comm_mysim);
2717 elapsed_time_over_all_ranks /= cr->nnodes;
2718 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2719 * current simulation. */
2720 MPI_Allreduce(&elapsed_time_over_all_threads,
2721 &elapsed_time_over_all_threads_over_all_ranks,
2722 1, MPI_DOUBLE, MPI_SUM,
2723 cr->mpi_comm_mysim);
2729 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2736 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2738 print_dd_statistics(cr, inputrec, fplog);
2743 wallclock_gpu_t* gputimes = use_GPU(nbv) ?
2744 nbnxn_cuda_get_timings(nbv->cu_nbv) : NULL;
2745 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2746 elapsed_time_over_all_ranks,
2749 if (EI_DYNAMICS(inputrec->eI))
2751 delta_t = inputrec->delta_t;
2760 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2761 elapsed_time_over_all_ranks,
2762 walltime_accounting_get_nsteps_done(walltime_accounting),
2763 delta_t, nbfs, mflop);
2767 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2768 elapsed_time_over_all_ranks,
2769 walltime_accounting_get_nsteps_done(walltime_accounting),
2770 delta_t, nbfs, mflop);
2775 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2777 /* this function works, but could probably use a logic rewrite to keep all the different
2778 types of efep straight. */
2781 t_lambda *fep = ir->fepvals;
2783 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2785 for (i = 0; i < efptNR; i++)
2797 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2798 if checkpoint is set -- a kludge is in for now
2800 for (i = 0; i < efptNR; i++)
2802 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2803 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2805 lambda[i] = fep->init_lambda;
2808 lam0[i] = lambda[i];
2813 lambda[i] = fep->all_lambda[i][*fep_state];
2816 lam0[i] = lambda[i];
2822 /* need to rescale control temperatures to match current state */
2823 for (i = 0; i < ir->opts.ngtc; i++)
2825 if (ir->opts.ref_t[i] > 0)
2827 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2833 /* Send to the log the information on the current lambdas */
2836 fprintf(fplog, "Initial vector of lambda components:[ ");
2837 for (i = 0; i < efptNR; i++)
2839 fprintf(fplog, "%10.4f ", lambda[i]);
2841 fprintf(fplog, "]\n");
2847 void init_md(FILE *fplog,
2848 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2849 double *t, double *t0,
2850 real *lambda, int *fep_state, double *lam0,
2851 t_nrnb *nrnb, gmx_mtop_t *mtop,
2853 int nfile, const t_filenm fnm[],
2854 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2855 tensor force_vir, tensor shake_vir, rvec mu_tot,
2856 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags,
2857 gmx_wallcycle_t wcycle)
2862 /* Initial values */
2863 *t = *t0 = ir->init_t;
2866 for (i = 0; i < ir->opts.ngtc; i++)
2868 /* set bSimAnn if any group is being annealed */
2869 if (ir->opts.annealing[i] != eannNO)
2876 update_annealing_target_temp(&(ir->opts), ir->init_t);
2879 /* Initialize lambda variables */
2880 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2884 *upd = init_update(ir);
2890 *vcm = init_vcm(fplog, &mtop->groups, ir);
2893 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2895 if (ir->etc == etcBERENDSEN)
2897 please_cite(fplog, "Berendsen84a");
2899 if (ir->etc == etcVRESCALE)
2901 please_cite(fplog, "Bussi2007a");
2903 if (ir->eI == eiSD1)
2905 please_cite(fplog, "Goga2012");
2913 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv, wcycle);
2915 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2916 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2921 please_cite(fplog, "Fritsch12");
2922 please_cite(fplog, "Junghans10");
2924 /* Initiate variables */
2925 clear_mat(force_vir);
2926 clear_mat(shake_vir);