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.
42 #ifdef HAVE_SYS_TIME_H
54 #include "chargegroup.h"
76 #include "nbnxn_atomdata.h"
77 #include "nbnxn_search.h"
78 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
79 #include "nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
80 #include "nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
81 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
82 #include "nonbonded.h"
83 #include "../gmxlib/nonbonded/nb_kernel.h"
84 #include "../gmxlib/nonbonded/nb_free_energy.h"
86 #include "gromacs/timing/wallcycle.h"
87 #include "gromacs/timing/walltime_accounting.h"
88 #include "gromacs/utility/gmxmpi.h"
89 #include "gromacs/essentialdynamics/edsam.h"
90 #include "gromacs/pulling/pull.h"
91 #include "gromacs/pulling/pull_rotation.h"
96 #include "gmx_omp_nthreads.h"
98 #include "nbnxn_cuda_data_mgmt.h"
99 #include "nbnxn_cuda/nbnxn_cuda.h"
101 void print_time(FILE *out,
102 gmx_walltime_accounting_t walltime_accounting,
105 t_commrec gmx_unused *cr)
108 char timebuf[STRLEN];
109 double dt, elapsed_seconds, time_per_step;
112 #ifndef GMX_THREAD_MPI
118 fprintf(out, "step %s", gmx_step_str(step, buf));
119 if ((step >= ir->nstlist))
121 double seconds_since_epoch = gmx_gettime();
122 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
123 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
124 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
130 finish = (time_t) (seconds_since_epoch + dt);
131 gmx_ctime_r(&finish, timebuf, STRLEN);
132 sprintf(buf, "%s", timebuf);
133 buf[strlen(buf)-1] = '\0';
134 fprintf(out, ", will finish %s", buf);
138 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
143 fprintf(out, " performance: %.1f ns/day ",
144 ir->delta_t/1000*24*60*60/time_per_step);
147 #ifndef GMX_THREAD_MPI
157 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
158 const gmx_walltime_accounting_t walltime_accounting)
161 char timebuf[STRLEN];
162 char time_string[STRLEN];
167 if (walltime_accounting != NULL)
169 tmptime = (time_t) walltime_accounting_get_start_time_stamp(walltime_accounting);
170 gmx_ctime_r(&tmptime, timebuf, STRLEN);
174 tmptime = (time_t) gmx_gettime();
175 gmx_ctime_r(&tmptime, timebuf, STRLEN);
177 for (i = 0; timebuf[i] >= ' '; i++)
179 time_string[i] = timebuf[i];
181 time_string[i] = '\0';
183 fprintf(fplog, "%s on node %d %s\n", title, nodeid, time_string);
187 void print_start(FILE *fplog, t_commrec *cr,
188 gmx_walltime_accounting_t walltime_accounting,
193 sprintf(buf, "Started %s", name);
194 print_date_and_time(fplog, cr->nodeid, buf, 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(FILE *fplog,
324 matrix box, rvec x[],
325 gmx_enerdata_t *enerd,
333 /* Position restraints always require full pbc */
334 set_pbc(&pbc, ir->ePBC, box);
336 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
337 top->idef.iparams_posres,
338 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
339 ir->ePBC == epbcNONE ? NULL : &pbc,
340 lambda[efptRESTRAINT], &dvdl,
341 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
344 gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
346 enerd->term[F_POSRES] += v;
347 /* If just the force constant changes, the FEP term is linear,
348 * but if k changes, it is not.
350 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
351 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
353 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
355 for (i = 0; i < enerd->n_lambda; i++)
357 real dvdl_dum, lambda_dum;
359 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
360 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
361 top->idef.iparams_posres,
362 (const rvec*)x, NULL, NULL,
363 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
364 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
365 enerd->enerpart_lambda[i] += v;
370 static void fbposres_wrapper(t_inputrec *ir,
373 matrix box, rvec x[],
374 gmx_enerdata_t *enerd,
380 /* Flat-bottomed position restraints always require full pbc */
381 set_pbc(&pbc, ir->ePBC, box);
382 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
383 top->idef.iparams_fbposres,
384 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
385 ir->ePBC == epbcNONE ? NULL : &pbc,
386 fr->rc_scaling, fr->ePBC, fr->posres_com);
387 enerd->term[F_FBPOSRES] += v;
388 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
391 static void pull_potential_wrapper(FILE *fplog,
395 matrix box, rvec x[],
399 gmx_enerdata_t *enerd,
406 /* Calculate the center of mass forces, this requires communication,
407 * which is why pull_potential is called close to other communication.
408 * The virial contribution is calculated directly,
409 * which is why we call pull_potential after calc_virial.
411 set_pbc(&pbc, ir->ePBC, box);
413 enerd->term[F_COM_PULL] +=
414 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
415 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
418 gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
420 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
423 static void pme_receive_force_ener(FILE *fplog,
426 gmx_wallcycle_t wcycle,
427 gmx_enerdata_t *enerd,
430 real e_q, e_lj, v, dvdl_q, dvdl_lj;
431 float cycles_ppdpme, cycles_seppme;
433 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
434 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
436 /* In case of node-splitting, the PP nodes receive the long-range
437 * forces, virial and energy from the PME nodes here.
439 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
442 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
443 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
447 gmx_print_sepdvdl(fplog, "Electrostatic PME mesh", e_q, dvdl_q);
448 gmx_print_sepdvdl(fplog, "Lennard-Jones PME mesh", e_lj, dvdl_lj);
450 enerd->term[F_COUL_RECIP] += e_q;
451 enerd->term[F_LJ_RECIP] += e_lj;
452 enerd->dvdl_lin[efptCOUL] += dvdl_q;
453 enerd->dvdl_lin[efptVDW] += dvdl_lj;
457 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
459 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
462 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
463 gmx_int64_t step, real pforce, rvec *x, rvec *f)
467 char buf[STEPSTRSIZE];
470 for (i = 0; i < md->homenr; i++)
473 /* We also catch NAN, if the compiler does not optimize this away. */
474 if (fn2 >= pf2 || fn2 != fn2)
476 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
477 gmx_step_str(step, buf),
478 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
483 static void post_process_forces(t_commrec *cr,
485 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
487 matrix box, rvec x[],
492 t_forcerec *fr, gmx_vsite_t *vsite,
499 /* Spread the mesh force on virtual sites to the other particles...
500 * This is parallellized. MPI communication is performed
501 * if the constructing atoms aren't local.
503 wallcycle_start(wcycle, ewcVSITESPREAD);
504 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
505 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
507 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
508 wallcycle_stop(wcycle, ewcVSITESPREAD);
510 if (flags & GMX_FORCE_VIRIAL)
512 /* Now add the forces, this is local */
515 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
519 sum_forces(0, mdatoms->homenr,
522 if (EEL_FULL(fr->eeltype))
524 /* Add the mesh contribution to the virial */
525 m_add(vir_force, fr->vir_el_recip, vir_force);
527 if (EVDW_PME(fr->vdwtype))
529 /* Add the mesh contribution to the virial */
530 m_add(vir_force, fr->vir_lj_recip, vir_force);
534 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
539 if (fr->print_force >= 0)
541 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
545 static void do_nb_verlet(t_forcerec *fr,
546 interaction_const_t *ic,
547 gmx_enerdata_t *enerd,
548 int flags, int ilocality,
551 gmx_wallcycle_t wcycle)
553 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
555 nonbonded_verlet_group_t *nbvg;
558 if (!(flags & GMX_FORCE_NONBONDED))
560 /* skip non-bonded calculation */
564 nbvg = &fr->nbv->grp[ilocality];
566 /* CUDA kernel launch overhead is already timed separately */
567 if (fr->cutoff_scheme != ecutsVERLET)
569 gmx_incons("Invalid cut-off scheme passed!");
572 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
576 wallcycle_sub_start(wcycle, ewcsNONBONDED);
578 switch (nbvg->kernel_type)
580 case nbnxnk4x4_PlainC:
581 nbnxn_kernel_ref(&nbvg->nbl_lists,
587 enerd->grpp.ener[egCOULSR],
589 enerd->grpp.ener[egBHAMSR] :
590 enerd->grpp.ener[egLJSR]);
593 case nbnxnk4xN_SIMD_4xN:
594 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
601 enerd->grpp.ener[egCOULSR],
603 enerd->grpp.ener[egBHAMSR] :
604 enerd->grpp.ener[egLJSR]);
606 case nbnxnk4xN_SIMD_2xNN:
607 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
614 enerd->grpp.ener[egCOULSR],
616 enerd->grpp.ener[egBHAMSR] :
617 enerd->grpp.ener[egLJSR]);
620 case nbnxnk8x8x8_CUDA:
621 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
624 case nbnxnk8x8x8_PlainC:
625 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
630 nbvg->nbat->out[0].f,
632 enerd->grpp.ener[egCOULSR],
634 enerd->grpp.ener[egBHAMSR] :
635 enerd->grpp.ener[egLJSR]);
639 gmx_incons("Invalid nonbonded kernel type passed!");
644 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
647 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
649 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
651 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
652 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
654 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
658 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
660 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
661 if (flags & GMX_FORCE_ENERGY)
663 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
664 enr_nbnxn_kernel_ljc += 1;
665 enr_nbnxn_kernel_lj += 1;
668 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
669 nbvg->nbl_lists.natpair_ljq);
670 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
671 nbvg->nbl_lists.natpair_lj);
672 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
673 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
674 nbvg->nbl_lists.natpair_q);
676 if (ic->vdw_modifier == eintmodFORCESWITCH)
678 /* We add up the switch cost separately */
679 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
680 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
682 if (ic->vdw_modifier == eintmodPOTSWITCH)
684 /* We add up the switch cost separately */
685 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
686 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
688 if (ic->vdwtype == evdwPME)
690 /* We add up the LJ Ewald cost separately */
691 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
692 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
696 static void do_nb_verlet_fep(nbnxn_pairlist_set_t *nbl_lists,
703 gmx_enerdata_t *enerd,
706 gmx_wallcycle_t wcycle)
709 nb_kernel_data_t kernel_data;
711 real dvdl_nb[efptNR];
716 /* Add short-range interactions */
717 donb_flags |= GMX_NONBONDED_DO_SR;
719 /* Currently all group scheme kernels always calculate (shift-)forces */
720 if (flags & GMX_FORCE_FORCES)
722 donb_flags |= GMX_NONBONDED_DO_FORCE;
724 if (flags & GMX_FORCE_VIRIAL)
726 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
728 if (flags & GMX_FORCE_ENERGY)
730 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
732 if (flags & GMX_FORCE_DO_LR)
734 donb_flags |= GMX_NONBONDED_DO_LR;
737 kernel_data.flags = donb_flags;
738 kernel_data.lambda = lambda;
739 kernel_data.dvdl = dvdl_nb;
741 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
742 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
744 /* reset free energy components */
745 for (i = 0; i < efptNR; i++)
750 assert(gmx_omp_nthreads_get(emntNonbonded) == nbl_lists->nnbl);
752 wallcycle_sub_start(wcycle, ewcsNONBONDED);
753 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
754 for (th = 0; th < nbl_lists->nnbl; th++)
756 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
757 x, f, fr, mdatoms, &kernel_data, nrnb);
760 if (fepvals->sc_alpha != 0)
762 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
763 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
767 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
768 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
771 /* If we do foreign lambda and we have soft-core interactions
772 * we have to recalculate the (non-linear) energies contributions.
774 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
776 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
777 kernel_data.lambda = lam_i;
778 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
779 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
780 /* Note that we add to kernel_data.dvdl, but ignore the result */
782 for (i = 0; i < enerd->n_lambda; i++)
784 for (j = 0; j < efptNR; j++)
786 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
788 reset_foreign_enerdata(enerd);
789 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
790 for (th = 0; th < nbl_lists->nnbl; th++)
792 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
793 x, f, fr, mdatoms, &kernel_data, nrnb);
796 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
797 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
801 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
804 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
805 t_inputrec *inputrec,
806 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
808 gmx_groups_t gmx_unused *groups,
809 matrix box, rvec x[], history_t *hist,
813 gmx_enerdata_t *enerd, t_fcdata *fcd,
814 real *lambda, t_graph *graph,
815 t_forcerec *fr, interaction_const_t *ic,
816 gmx_vsite_t *vsite, rvec mu_tot,
817 double t, FILE *field, gmx_edsam_t ed,
825 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
826 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
827 gmx_bool bDiffKernels = FALSE;
829 rvec vzero, box_diag;
831 float cycles_pme, cycles_force, cycles_wait_gpu;
832 nonbonded_verlet_t *nbv;
837 nb_kernel_type = fr->nbv->grp[0].kernel_type;
840 homenr = mdatoms->homenr;
842 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
844 clear_mat(vir_force);
847 if (DOMAINDECOMP(cr))
849 cg1 = cr->dd->ncg_tot;
860 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
861 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
862 bFillGrid = (bNS && bStateChanged);
863 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
864 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
865 bDoForces = (flags & GMX_FORCE_FORCES);
866 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
867 bUseGPU = fr->nbv->bUseGPU;
868 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
872 update_forcerec(fr, box);
874 if (NEED_MUTOT(*inputrec))
876 /* Calculate total (local) dipole moment in a temporary common array.
877 * This makes it possible to sum them over nodes faster.
879 calc_mu(start, homenr,
880 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
885 if (fr->ePBC != epbcNONE)
887 /* Compute shift vectors every step,
888 * because of pressure coupling or box deformation!
890 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
892 calc_shifts(box, fr->shift_vec);
897 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
898 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
900 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
902 unshift_self(graph, box, x);
906 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
907 fr->shift_vec, nbv->grp[0].nbat);
910 if (!(cr->duty & DUTY_PME))
912 /* Send particle coordinates to the pme nodes.
913 * Since this is only implemented for domain decomposition
914 * and domain decomposition does not use the graph,
915 * we do not need to worry about shifting.
920 wallcycle_start(wcycle, ewcPP_PMESENDX);
922 bBS = (inputrec->nwall == 2);
926 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
929 if (EEL_PME(fr->eeltype))
931 pme_flags |= GMX_PME_DO_COULOMB;
934 if (EVDW_PME(fr->vdwtype))
936 pme_flags |= GMX_PME_DO_LJ;
939 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
940 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
941 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
944 wallcycle_stop(wcycle, ewcPP_PMESENDX);
948 /* do gridding for pair search */
951 if (graph && bStateChanged)
953 /* Calculate intramolecular shift vectors to make molecules whole */
954 mk_mshift(fplog, graph, fr->ePBC, box, x);
958 box_diag[XX] = box[XX][XX];
959 box_diag[YY] = box[YY][YY];
960 box_diag[ZZ] = box[ZZ][ZZ];
962 wallcycle_start(wcycle, ewcNS);
965 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
966 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
968 0, mdatoms->homenr, -1, fr->cginfo, x,
970 nbv->grp[eintLocal].kernel_type,
971 nbv->grp[eintLocal].nbat);
972 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
976 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
977 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
979 nbv->grp[eintNonlocal].kernel_type,
980 nbv->grp[eintNonlocal].nbat);
981 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
984 if (nbv->ngrp == 1 ||
985 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
987 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
988 nbv->nbs, mdatoms, fr->cginfo);
992 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
993 nbv->nbs, mdatoms, fr->cginfo);
994 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
995 nbv->nbs, mdatoms, fr->cginfo);
997 wallcycle_stop(wcycle, ewcNS);
1000 /* initialize the GPU atom data and copy shift vector */
1005 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1006 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
1007 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1010 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1011 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
1012 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1015 /* do local pair search */
1018 wallcycle_start_nocount(wcycle, ewcNS);
1019 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
1020 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
1023 nbv->min_ci_balanced,
1024 &nbv->grp[eintLocal].nbl_lists,
1026 nbv->grp[eintLocal].kernel_type,
1028 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
1032 /* initialize local pair-list on the GPU */
1033 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1034 nbv->grp[eintLocal].nbl_lists.nbl[0],
1037 wallcycle_stop(wcycle, ewcNS);
1041 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1042 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1043 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
1044 nbv->grp[eintLocal].nbat);
1045 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1046 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1051 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1052 /* launch local nonbonded F on GPU */
1053 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
1055 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1058 /* Communicate coordinates and sum dipole if necessary +
1059 do non-local pair search */
1060 if (DOMAINDECOMP(cr))
1062 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
1063 nbv->grp[eintLocal].kernel_type);
1067 /* With GPU+CPU non-bonded calculations we need to copy
1068 * the local coordinates to the non-local nbat struct
1069 * (in CPU format) as the non-local kernel call also
1070 * calculates the local - non-local interactions.
1072 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1073 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1074 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
1075 nbv->grp[eintNonlocal].nbat);
1076 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1077 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1082 wallcycle_start_nocount(wcycle, ewcNS);
1083 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1087 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
1090 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
1093 nbv->min_ci_balanced,
1094 &nbv->grp[eintNonlocal].nbl_lists,
1096 nbv->grp[eintNonlocal].kernel_type,
1099 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1101 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1103 /* initialize non-local pair-list on the GPU */
1104 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1105 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1108 wallcycle_stop(wcycle, ewcNS);
1112 wallcycle_start(wcycle, ewcMOVEX);
1113 dd_move_x(cr->dd, box, x);
1115 /* When we don't need the total dipole we sum it in global_stat */
1116 if (bStateChanged && NEED_MUTOT(*inputrec))
1118 gmx_sumd(2*DIM, mu, cr);
1120 wallcycle_stop(wcycle, ewcMOVEX);
1122 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1123 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1124 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1125 nbv->grp[eintNonlocal].nbat);
1126 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1127 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1130 if (bUseGPU && !bDiffKernels)
1132 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1133 /* launch non-local nonbonded F on GPU */
1134 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1136 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1142 /* launch D2H copy-back F */
1143 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1144 if (DOMAINDECOMP(cr) && !bDiffKernels)
1146 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1147 flags, eatNonlocal);
1149 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1151 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1154 if (bStateChanged && NEED_MUTOT(*inputrec))
1158 gmx_sumd(2*DIM, mu, cr);
1161 for (i = 0; i < 2; i++)
1163 for (j = 0; j < DIM; j++)
1165 fr->mu_tot[i][j] = mu[i*DIM + j];
1169 if (fr->efep == efepNO)
1171 copy_rvec(fr->mu_tot[0], mu_tot);
1175 for (j = 0; j < DIM; j++)
1178 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1179 lambda[efptCOUL]*fr->mu_tot[1][j];
1183 /* Reset energies */
1184 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1185 clear_rvecs(SHIFTS, fr->fshift);
1187 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1189 wallcycle_start(wcycle, ewcPPDURINGPME);
1190 dd_force_flop_start(cr->dd, nrnb);
1195 /* Enforced rotation has its own cycle counter that starts after the collective
1196 * coordinates have been communicated. It is added to ddCyclF to allow
1197 * for proper load-balancing */
1198 wallcycle_start(wcycle, ewcROT);
1199 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1200 wallcycle_stop(wcycle, ewcROT);
1203 /* Start the force cycle counter.
1204 * This counter is stopped in do_forcelow_level.
1205 * No parallel communication should occur while this counter is running,
1206 * since that will interfere with the dynamic load balancing.
1208 wallcycle_start(wcycle, ewcFORCE);
1211 /* Reset forces for which the virial is calculated separately:
1212 * PME/Ewald forces if necessary */
1213 if (fr->bF_NoVirSum)
1215 if (flags & GMX_FORCE_VIRIAL)
1217 fr->f_novirsum = fr->f_novirsum_alloc;
1220 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1224 clear_rvecs(homenr, fr->f_novirsum+start);
1229 /* We are not calculating the pressure so we do not need
1230 * a separate array for forces that do not contribute
1237 /* Clear the short- and long-range forces */
1238 clear_rvecs(fr->natoms_force_constr, f);
1239 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1241 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1244 clear_rvec(fr->vir_diag_posres);
1247 if (inputrec->ePull == epullCONSTRAINT)
1249 clear_pull_forces(inputrec->pull);
1252 /* We calculate the non-bonded forces, when done on the CPU, here.
1253 * We do this before calling do_force_lowlevel, as in there bondeds
1254 * forces are calculated before PME, which does communication.
1255 * With this order, non-bonded and bonded force calculation imbalance
1256 * can be balanced out by the domain decomposition load balancing.
1261 /* Maybe we should move this into do_force_lowlevel */
1262 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1266 if (fr->efep != efepNO)
1268 /* Calculate the local and non-local free energy interactions here.
1269 * Happens here on the CPU both with and without GPU.
1271 if (fr->nbv->grp[eintLocal].nbl_lists.nbl_fep[0]->nrj > 0)
1273 do_nb_verlet_fep(&fr->nbv->grp[eintLocal].nbl_lists,
1275 inputrec->fepvals, lambda,
1276 enerd, flags, nrnb, wcycle);
1279 if (DOMAINDECOMP(cr) &&
1280 fr->nbv->grp[eintNonlocal].nbl_lists.nbl_fep[0]->nrj > 0)
1282 do_nb_verlet_fep(&fr->nbv->grp[eintNonlocal].nbl_lists,
1284 inputrec->fepvals, lambda,
1285 enerd, flags, nrnb, wcycle);
1289 if (!bUseOrEmulGPU || bDiffKernels)
1293 if (DOMAINDECOMP(cr))
1295 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1296 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1306 aloc = eintNonlocal;
1309 /* Add all the non-bonded force to the normal force array.
1310 * This can be split into a local a non-local part when overlapping
1311 * communication with calculation with domain decomposition.
1313 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1314 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1315 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1316 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1317 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1318 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1319 wallcycle_start_nocount(wcycle, ewcFORCE);
1321 /* if there are multiple fshift output buffers reduce them */
1322 if ((flags & GMX_FORCE_VIRIAL) &&
1323 nbv->grp[aloc].nbl_lists.nnbl > 1)
1325 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1330 /* update QMMMrec, if necessary */
1333 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1336 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1338 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1342 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1344 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1347 /* Compute the bonded and non-bonded energies and optionally forces */
1348 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1349 cr, nrnb, wcycle, mdatoms,
1350 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1351 &(top->atomtypes), bBornRadii, box,
1352 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1353 flags, &cycles_pme);
1357 if (do_per_step(step, inputrec->nstcalclr))
1359 /* Add the long range forces to the short range forces */
1360 for (i = 0; i < fr->natoms_force_constr; i++)
1362 rvec_add(fr->f_twin[i], f[i], f[i]);
1367 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1371 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1374 if (bUseOrEmulGPU && !bDiffKernels)
1376 /* wait for non-local forces (or calculate in emulation mode) */
1377 if (DOMAINDECOMP(cr))
1383 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1384 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1385 nbv->grp[eintNonlocal].nbat,
1387 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1389 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1390 cycles_wait_gpu += cycles_tmp;
1391 cycles_force += cycles_tmp;
1395 wallcycle_start_nocount(wcycle, ewcFORCE);
1396 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1398 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1400 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1401 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1402 /* skip the reduction if there was no non-local work to do */
1403 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1405 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1406 nbv->grp[eintNonlocal].nbat, f);
1408 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1409 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1413 if (bDoForces && DOMAINDECOMP(cr))
1415 /* Communicate the forces */
1416 wallcycle_start(wcycle, ewcMOVEF);
1417 dd_move_f(cr->dd, f, fr->fshift);
1418 /* Do we need to communicate the separate force array
1419 * for terms that do not contribute to the single sum virial?
1420 * Position restraints and electric fields do not introduce
1421 * inter-cg forces, only full electrostatics methods do.
1422 * When we do not calculate the virial, fr->f_novirsum = f,
1423 * so we have already communicated these forces.
1425 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1426 (flags & GMX_FORCE_VIRIAL))
1428 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1432 /* We should not update the shift forces here,
1433 * since f_twin is already included in f.
1435 dd_move_f(cr->dd, fr->f_twin, NULL);
1437 wallcycle_stop(wcycle, ewcMOVEF);
1442 /* wait for local forces (or calculate in emulation mode) */
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_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1453 /* now clear the GPU outputs while we finish the step on the CPU */
1455 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1456 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1457 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1461 wallcycle_start_nocount(wcycle, ewcFORCE);
1462 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1463 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1465 wallcycle_stop(wcycle, ewcFORCE);
1467 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1468 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1469 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1471 /* skip the reduction if there was no non-local work to do */
1472 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1473 nbv->grp[eintLocal].nbat, f);
1475 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1476 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1479 if (DOMAINDECOMP(cr))
1481 dd_force_flop_stop(cr->dd, nrnb);
1484 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1487 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1494 if (IR_ELEC_FIELD(*inputrec))
1496 /* Compute forces due to electric field */
1497 calc_f_el(MASTER(cr) ? field : NULL,
1498 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1499 inputrec->ex, inputrec->et, t);
1502 /* If we have NoVirSum forces, but we do not calculate the virial,
1503 * we sum fr->f_novirum=f later.
1505 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1507 wallcycle_start(wcycle, ewcVSITESPREAD);
1508 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1509 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1510 wallcycle_stop(wcycle, ewcVSITESPREAD);
1514 wallcycle_start(wcycle, ewcVSITESPREAD);
1515 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1517 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1518 wallcycle_stop(wcycle, ewcVSITESPREAD);
1522 if (flags & GMX_FORCE_VIRIAL)
1524 /* Calculation of the virial must be done after vsites! */
1525 calc_virial(0, mdatoms->homenr, x, f,
1526 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1530 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1532 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1533 f, vir_force, mdatoms, enerd, lambda, t);
1536 /* Add the forces from enforced rotation potentials (if any) */
1539 wallcycle_start(wcycle, ewcROTadd);
1540 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1541 wallcycle_stop(wcycle, ewcROTadd);
1544 if (PAR(cr) && !(cr->duty & DUTY_PME))
1546 /* In case of node-splitting, the PP nodes receive the long-range
1547 * forces, virial and energy from the PME nodes here.
1549 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1554 post_process_forces(cr, step, nrnb, wcycle,
1555 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1559 /* Sum the potential energy terms from group contributions */
1560 sum_epot(&(enerd->grpp), enerd->term);
1563 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1564 t_inputrec *inputrec,
1565 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1566 gmx_localtop_t *top,
1567 gmx_groups_t *groups,
1568 matrix box, rvec x[], history_t *hist,
1572 gmx_enerdata_t *enerd, t_fcdata *fcd,
1573 real *lambda, t_graph *graph,
1574 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1575 double t, FILE *field, gmx_edsam_t ed,
1576 gmx_bool bBornRadii,
1582 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1583 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1584 gmx_bool bDoAdressWF;
1586 rvec vzero, box_diag;
1587 real e, v, dvdlambda[efptNR];
1589 float cycles_pme, cycles_force;
1592 homenr = mdatoms->homenr;
1594 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1596 clear_mat(vir_force);
1599 if (DOMAINDECOMP(cr))
1601 cg1 = cr->dd->ncg_tot;
1612 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1613 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1614 /* Should we update the long-range neighborlists at this step? */
1615 bDoLongRangeNS = fr->bTwinRange && bNS;
1616 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1617 bFillGrid = (bNS && bStateChanged);
1618 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1619 bDoForces = (flags & GMX_FORCE_FORCES);
1620 bDoPotential = (flags & GMX_FORCE_ENERGY);
1621 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1622 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1624 /* should probably move this to the forcerec since it doesn't change */
1625 bDoAdressWF = ((fr->adress_type != eAdressOff));
1629 update_forcerec(fr, box);
1631 if (NEED_MUTOT(*inputrec))
1633 /* Calculate total (local) dipole moment in a temporary common array.
1634 * This makes it possible to sum them over nodes faster.
1636 calc_mu(start, homenr,
1637 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1642 if (fr->ePBC != epbcNONE)
1644 /* Compute shift vectors every step,
1645 * because of pressure coupling or box deformation!
1647 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1649 calc_shifts(box, fr->shift_vec);
1654 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1655 &(top->cgs), x, fr->cg_cm);
1656 inc_nrnb(nrnb, eNR_CGCM, homenr);
1657 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1659 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1661 unshift_self(graph, box, x);
1666 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1667 inc_nrnb(nrnb, eNR_CGCM, homenr);
1670 if (bCalcCGCM && gmx_debug_at)
1672 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1676 if (!(cr->duty & DUTY_PME))
1678 /* Send particle coordinates to the pme nodes.
1679 * Since this is only implemented for domain decomposition
1680 * and domain decomposition does not use the graph,
1681 * we do not need to worry about shifting.
1686 wallcycle_start(wcycle, ewcPP_PMESENDX);
1688 bBS = (inputrec->nwall == 2);
1691 copy_mat(box, boxs);
1692 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1695 if (EEL_PME(fr->eeltype))
1697 pme_flags |= GMX_PME_DO_COULOMB;
1700 if (EVDW_PME(fr->vdwtype))
1702 pme_flags |= GMX_PME_DO_LJ;
1705 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1706 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1707 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1710 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1712 #endif /* GMX_MPI */
1714 /* Communicate coordinates and sum dipole if necessary */
1715 if (DOMAINDECOMP(cr))
1717 wallcycle_start(wcycle, ewcMOVEX);
1718 dd_move_x(cr->dd, box, x);
1719 wallcycle_stop(wcycle, ewcMOVEX);
1722 /* update adress weight beforehand */
1723 if (bStateChanged && bDoAdressWF)
1725 /* need pbc for adress weight calculation with pbc_dx */
1726 set_pbc(&pbc, inputrec->ePBC, box);
1727 if (fr->adress_site == eAdressSITEcog)
1729 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1730 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1732 else if (fr->adress_site == eAdressSITEcom)
1734 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1735 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1737 else if (fr->adress_site == eAdressSITEatomatom)
1739 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1740 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1744 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1745 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1749 if (NEED_MUTOT(*inputrec))
1756 gmx_sumd(2*DIM, mu, cr);
1758 for (i = 0; i < 2; i++)
1760 for (j = 0; j < DIM; j++)
1762 fr->mu_tot[i][j] = mu[i*DIM + j];
1766 if (fr->efep == efepNO)
1768 copy_rvec(fr->mu_tot[0], mu_tot);
1772 for (j = 0; j < DIM; j++)
1775 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1780 /* Reset energies */
1781 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1782 clear_rvecs(SHIFTS, fr->fshift);
1786 wallcycle_start(wcycle, ewcNS);
1788 if (graph && bStateChanged)
1790 /* Calculate intramolecular shift vectors to make molecules whole */
1791 mk_mshift(fplog, graph, fr->ePBC, box, x);
1794 /* Do the actual neighbour searching */
1796 groups, top, mdatoms,
1797 cr, nrnb, bFillGrid,
1800 wallcycle_stop(wcycle, ewcNS);
1803 if (inputrec->implicit_solvent && bNS)
1805 make_gb_nblist(cr, inputrec->gb_algorithm,
1806 x, box, fr, &top->idef, graph, fr->born);
1809 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1811 wallcycle_start(wcycle, ewcPPDURINGPME);
1812 dd_force_flop_start(cr->dd, nrnb);
1817 /* Enforced rotation has its own cycle counter that starts after the collective
1818 * coordinates have been communicated. It is added to ddCyclF to allow
1819 * for proper load-balancing */
1820 wallcycle_start(wcycle, ewcROT);
1821 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1822 wallcycle_stop(wcycle, ewcROT);
1825 /* Start the force cycle counter.
1826 * This counter is stopped in do_forcelow_level.
1827 * No parallel communication should occur while this counter is running,
1828 * since that will interfere with the dynamic load balancing.
1830 wallcycle_start(wcycle, ewcFORCE);
1834 /* Reset forces for which the virial is calculated separately:
1835 * PME/Ewald forces if necessary */
1836 if (fr->bF_NoVirSum)
1838 if (flags & GMX_FORCE_VIRIAL)
1840 fr->f_novirsum = fr->f_novirsum_alloc;
1843 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1847 clear_rvecs(homenr, fr->f_novirsum+start);
1852 /* We are not calculating the pressure so we do not need
1853 * a separate array for forces that do not contribute
1860 /* Clear the short- and long-range forces */
1861 clear_rvecs(fr->natoms_force_constr, f);
1862 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1864 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1867 clear_rvec(fr->vir_diag_posres);
1869 if (inputrec->ePull == epullCONSTRAINT)
1871 clear_pull_forces(inputrec->pull);
1874 /* update QMMMrec, if necessary */
1877 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1880 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1882 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1886 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1888 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1891 /* Compute the bonded and non-bonded energies and optionally forces */
1892 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1893 cr, nrnb, wcycle, mdatoms,
1894 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1895 &(top->atomtypes), bBornRadii, box,
1896 inputrec->fepvals, lambda,
1897 graph, &(top->excls), fr->mu_tot,
1903 if (do_per_step(step, inputrec->nstcalclr))
1905 /* Add the long range forces to the short range forces */
1906 for (i = 0; i < fr->natoms_force_constr; i++)
1908 rvec_add(fr->f_twin[i], f[i], f[i]);
1913 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1917 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1920 if (DOMAINDECOMP(cr))
1922 dd_force_flop_stop(cr->dd, nrnb);
1925 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1931 if (IR_ELEC_FIELD(*inputrec))
1933 /* Compute forces due to electric field */
1934 calc_f_el(MASTER(cr) ? field : NULL,
1935 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1936 inputrec->ex, inputrec->et, t);
1939 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1941 /* Compute thermodynamic force in hybrid AdResS region */
1942 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1943 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1946 /* Communicate the forces */
1947 if (DOMAINDECOMP(cr))
1949 wallcycle_start(wcycle, ewcMOVEF);
1950 dd_move_f(cr->dd, f, fr->fshift);
1951 /* Do we need to communicate the separate force array
1952 * for terms that do not contribute to the single sum virial?
1953 * Position restraints and electric fields do not introduce
1954 * inter-cg forces, only full electrostatics methods do.
1955 * When we do not calculate the virial, fr->f_novirsum = f,
1956 * so we have already communicated these forces.
1958 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1959 (flags & GMX_FORCE_VIRIAL))
1961 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1965 /* We should not update the shift forces here,
1966 * since f_twin is already included in f.
1968 dd_move_f(cr->dd, fr->f_twin, NULL);
1970 wallcycle_stop(wcycle, ewcMOVEF);
1973 /* If we have NoVirSum forces, but we do not calculate the virial,
1974 * we sum fr->f_novirum=f later.
1976 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1978 wallcycle_start(wcycle, ewcVSITESPREAD);
1979 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1980 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1981 wallcycle_stop(wcycle, ewcVSITESPREAD);
1985 wallcycle_start(wcycle, ewcVSITESPREAD);
1986 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1988 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1989 wallcycle_stop(wcycle, ewcVSITESPREAD);
1993 if (flags & GMX_FORCE_VIRIAL)
1995 /* Calculation of the virial must be done after vsites! */
1996 calc_virial(0, mdatoms->homenr, x, f,
1997 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
2001 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
2003 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
2004 f, vir_force, mdatoms, enerd, lambda, t);
2007 /* Add the forces from enforced rotation potentials (if any) */
2010 wallcycle_start(wcycle, ewcROTadd);
2011 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
2012 wallcycle_stop(wcycle, ewcROTadd);
2015 if (PAR(cr) && !(cr->duty & DUTY_PME))
2017 /* In case of node-splitting, the PP nodes receive the long-range
2018 * forces, virial and energy from the PME nodes here.
2020 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
2025 post_process_forces(cr, step, nrnb, wcycle,
2026 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
2030 /* Sum the potential energy terms from group contributions */
2031 sum_epot(&(enerd->grpp), enerd->term);
2034 void do_force(FILE *fplog, t_commrec *cr,
2035 t_inputrec *inputrec,
2036 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2037 gmx_localtop_t *top,
2038 gmx_groups_t *groups,
2039 matrix box, rvec x[], history_t *hist,
2043 gmx_enerdata_t *enerd, t_fcdata *fcd,
2044 real *lambda, t_graph *graph,
2046 gmx_vsite_t *vsite, rvec mu_tot,
2047 double t, FILE *field, gmx_edsam_t ed,
2048 gmx_bool bBornRadii,
2051 /* modify force flag if not doing nonbonded */
2052 if (!fr->bNonbonded)
2054 flags &= ~GMX_FORCE_NONBONDED;
2057 switch (inputrec->cutoff_scheme)
2060 do_force_cutsVERLET(fplog, cr, inputrec,
2076 do_force_cutsGROUP(fplog, cr, inputrec,
2091 gmx_incons("Invalid cut-off scheme passed!");
2096 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2097 t_inputrec *ir, t_mdatoms *md,
2098 t_state *state, t_commrec *cr, t_nrnb *nrnb,
2099 t_forcerec *fr, gmx_localtop_t *top)
2101 int i, m, start, end;
2103 real dt = ir->delta_t;
2107 snew(savex, state->natoms);
2114 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2115 start, md->homenr, end);
2117 /* Do a first constrain to reset particles... */
2118 step = ir->init_step;
2121 char buf[STEPSTRSIZE];
2122 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2123 gmx_step_str(step, buf));
2127 /* constrain the current position */
2128 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2129 ir, NULL, cr, step, 0, md,
2130 state->x, state->x, NULL,
2131 fr->bMolPBC, state->box,
2132 state->lambda[efptBONDED], &dvdl_dum,
2133 NULL, NULL, nrnb, econqCoord,
2134 ir->epc == epcMTTK, state->veta, state->veta);
2137 /* constrain the inital velocity, and save it */
2138 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2139 /* might not yet treat veta correctly */
2140 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2141 ir, NULL, cr, step, 0, md,
2142 state->x, state->v, state->v,
2143 fr->bMolPBC, state->box,
2144 state->lambda[efptBONDED], &dvdl_dum,
2145 NULL, NULL, nrnb, econqVeloc,
2146 ir->epc == epcMTTK, state->veta, state->veta);
2148 /* constrain the inital velocities at t-dt/2 */
2149 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2151 for (i = start; (i < end); i++)
2153 for (m = 0; (m < DIM); m++)
2155 /* Reverse the velocity */
2156 state->v[i][m] = -state->v[i][m];
2157 /* Store the position at t-dt in buf */
2158 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2161 /* Shake the positions at t=-dt with the positions at t=0
2162 * as reference coordinates.
2166 char buf[STEPSTRSIZE];
2167 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2168 gmx_step_str(step, buf));
2171 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2172 ir, NULL, cr, step, -1, md,
2173 state->x, savex, NULL,
2174 fr->bMolPBC, state->box,
2175 state->lambda[efptBONDED], &dvdl_dum,
2176 state->v, NULL, nrnb, econqCoord,
2177 ir->epc == epcMTTK, state->veta, state->veta);
2179 for (i = start; i < end; i++)
2181 for (m = 0; m < DIM; m++)
2183 /* Re-reverse the velocities */
2184 state->v[i][m] = -state->v[i][m];
2193 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2194 double *enerout, double *virout)
2196 double enersum, virsum;
2197 double invscale, invscale2, invscale3;
2198 double r, ea, eb, ec, pa, pb, pc, pd;
2200 int ri, offset, tabfactor;
2202 invscale = 1.0/scale;
2203 invscale2 = invscale*invscale;
2204 invscale3 = invscale*invscale2;
2206 /* Following summation derived from cubic spline definition,
2207 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2208 * the cubic spline. We first calculate the negative of the
2209 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2210 * add the more standard, abrupt cutoff correction to that result,
2211 * yielding the long-range correction for a switched function. We
2212 * perform both the pressure and energy loops at the same time for
2213 * simplicity, as the computational cost is low. */
2217 /* Since the dispersion table has been scaled down a factor
2218 * 6.0 and the repulsion a factor 12.0 to compensate for the
2219 * c6/c12 parameters inside nbfp[] being scaled up (to save
2220 * flops in kernels), we need to correct for this.
2231 for (ri = rstart; ri < rend; ++ri)
2235 eb = 2.0*invscale2*r;
2239 pb = 3.0*invscale2*r;
2240 pc = 3.0*invscale*r*r;
2243 /* this "8" is from the packing in the vdwtab array - perhaps
2244 should be defined? */
2246 offset = 8*ri + offstart;
2247 y0 = vdwtab[offset];
2248 f = vdwtab[offset+1];
2249 g = vdwtab[offset+2];
2250 h = vdwtab[offset+3];
2252 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);
2253 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);
2255 *enerout = 4.0*M_PI*enersum*tabfactor;
2256 *virout = 4.0*M_PI*virsum*tabfactor;
2259 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2261 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2262 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2263 double invscale, invscale2, invscale3;
2264 int ri0, ri1, ri, i, offstart, offset;
2265 real scale, *vdwtab, tabfactor, tmp;
2267 fr->enershiftsix = 0;
2268 fr->enershifttwelve = 0;
2269 fr->enerdiffsix = 0;
2270 fr->enerdifftwelve = 0;
2272 fr->virdifftwelve = 0;
2274 if (eDispCorr != edispcNO)
2276 for (i = 0; i < 2; i++)
2281 if (fr->vdwtype == evdwSWITCH || fr->vdwtype == evdwSHIFT ||
2282 fr->vdw_modifier == eintmodPOTSWITCH ||
2283 fr->vdw_modifier == eintmodFORCESWITCH)
2285 if (fr->rvdw_switch == 0)
2288 "With dispersion correction rvdw-switch can not be zero "
2289 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2292 scale = fr->nblists[0].table_elec_vdw.scale;
2293 vdwtab = fr->nblists[0].table_vdw.data;
2295 /* Round the cut-offs to exact table values for precision */
2296 ri0 = floor(fr->rvdw_switch*scale);
2297 ri1 = ceil(fr->rvdw*scale);
2303 if (fr->vdwtype == evdwSHIFT ||
2304 fr->vdw_modifier == eintmodFORCESWITCH)
2306 /* Determine the constant energy shift below rvdw_switch.
2307 * Table has a scale factor since we have scaled it down to compensate
2308 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2310 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2311 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2313 /* Add the constant part from 0 to rvdw_switch.
2314 * This integration from 0 to rvdw_switch overcounts the number
2315 * of interactions by 1, as it also counts the self interaction.
2316 * We will correct for this later.
2318 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2319 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2320 for (i = 0; i < 2; i++)
2324 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2325 eners[i] -= enersum;
2329 /* now add the correction for rvdw_switch to infinity */
2330 eners[0] += -4.0*M_PI/(3.0*rc3);
2331 eners[1] += 4.0*M_PI/(9.0*rc9);
2332 virs[0] += 8.0*M_PI/rc3;
2333 virs[1] += -16.0*M_PI/(3.0*rc9);
2335 else if (fr->vdwtype == evdwCUT ||
2336 EVDW_PME(fr->vdwtype) ||
2337 fr->vdwtype == evdwUSER)
2339 if (fr->vdwtype == evdwUSER && fplog)
2342 "WARNING: using dispersion correction with user tables\n");
2345 /* Note that with LJ-PME, the dispersion correction is multiplied
2346 * by the difference between the actual C6 and the value of C6
2347 * that would produce the combination rule.
2348 * This means the normal energy and virial difference formulas
2352 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2354 /* Contribution beyond the cut-off */
2355 eners[0] += -4.0*M_PI/(3.0*rc3);
2356 eners[1] += 4.0*M_PI/(9.0*rc9);
2357 if (fr->vdw_modifier == eintmodPOTSHIFT)
2359 /* Contribution within the cut-off */
2360 eners[0] += -4.0*M_PI/(3.0*rc3);
2361 eners[1] += 4.0*M_PI/(3.0*rc9);
2363 /* Contribution beyond the cut-off */
2364 virs[0] += 8.0*M_PI/rc3;
2365 virs[1] += -16.0*M_PI/(3.0*rc9);
2370 "Dispersion correction is not implemented for vdw-type = %s",
2371 evdw_names[fr->vdwtype]);
2374 /* TODO: remove this code once we have group LJ-PME kernels
2375 * that calculate the exact, full LJ param C6/r^6 within the cut-off,
2376 * as the current nbnxn kernels do.
2378 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2380 /* Calculate self-interaction coefficient (assuming that
2381 * the reciprocal-space contribution is constant in the
2382 * region that contributes to the self-interaction).
2384 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2386 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2387 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2390 fr->enerdiffsix = eners[0];
2391 fr->enerdifftwelve = eners[1];
2392 /* The 0.5 is due to the Gromacs definition of the virial */
2393 fr->virdiffsix = 0.5*virs[0];
2394 fr->virdifftwelve = 0.5*virs[1];
2398 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2399 gmx_int64_t step, int natoms,
2400 matrix box, real lambda, tensor pres, tensor virial,
2401 real *prescorr, real *enercorr, real *dvdlcorr)
2403 gmx_bool bCorrAll, bCorrPres;
2404 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2414 if (ir->eDispCorr != edispcNO)
2416 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2417 ir->eDispCorr == edispcAllEnerPres);
2418 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2419 ir->eDispCorr == edispcAllEnerPres);
2421 invvol = 1/det(box);
2424 /* Only correct for the interactions with the inserted molecule */
2425 dens = (natoms - fr->n_tpi)*invvol;
2430 dens = natoms*invvol;
2431 ninter = 0.5*natoms;
2434 if (ir->efep == efepNO)
2436 avcsix = fr->avcsix[0];
2437 avctwelve = fr->avctwelve[0];
2441 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2442 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2445 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2446 *enercorr += avcsix*enerdiff;
2448 if (ir->efep != efepNO)
2450 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2454 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2455 *enercorr += avctwelve*enerdiff;
2456 if (fr->efep != efepNO)
2458 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2464 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2465 if (ir->eDispCorr == edispcAllEnerPres)
2467 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2469 /* The factor 2 is because of the Gromacs virial definition */
2470 spres = -2.0*invvol*svir*PRESFAC;
2472 for (m = 0; m < DIM; m++)
2474 virial[m][m] += svir;
2475 pres[m][m] += spres;
2480 /* Can't currently control when it prints, for now, just print when degugging */
2485 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2491 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2492 *enercorr, spres, svir);
2496 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2500 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2502 gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
2504 if (fr->efep != efepNO)
2506 *dvdlcorr += dvdlambda;
2511 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2512 t_graph *graph, rvec x[])
2516 fprintf(fplog, "Removing pbc first time\n");
2518 calc_shifts(box, fr->shift_vec);
2521 mk_mshift(fplog, graph, fr->ePBC, box, x);
2524 p_graph(debug, "do_pbc_first 1", graph);
2526 shift_self(graph, box, x);
2527 /* By doing an extra mk_mshift the molecules that are broken
2528 * because they were e.g. imported from another software
2529 * will be made whole again. Such are the healing powers
2532 mk_mshift(fplog, graph, fr->ePBC, box, x);
2535 p_graph(debug, "do_pbc_first 2", graph);
2540 fprintf(fplog, "Done rmpbc\n");
2544 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2545 gmx_mtop_t *mtop, rvec x[],
2550 gmx_molblock_t *molb;
2552 if (bFirst && fplog)
2554 fprintf(fplog, "Removing pbc first time\n");
2559 for (mb = 0; mb < mtop->nmolblock; mb++)
2561 molb = &mtop->molblock[mb];
2562 if (molb->natoms_mol == 1 ||
2563 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2565 /* Just one atom or charge group in the molecule, no PBC required */
2566 as += molb->nmol*molb->natoms_mol;
2570 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2571 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2572 0, molb->natoms_mol, FALSE, FALSE, graph);
2574 for (mol = 0; mol < molb->nmol; mol++)
2576 mk_mshift(fplog, graph, ePBC, box, x+as);
2578 shift_self(graph, box, x+as);
2579 /* The molecule is whole now.
2580 * We don't need the second mk_mshift call as in do_pbc_first,
2581 * since we no longer need this graph.
2584 as += molb->natoms_mol;
2592 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2593 gmx_mtop_t *mtop, rvec x[])
2595 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2598 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2599 gmx_mtop_t *mtop, rvec x[])
2601 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2604 void finish_run(FILE *fplog, t_commrec *cr,
2605 t_inputrec *inputrec,
2606 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2607 gmx_walltime_accounting_t walltime_accounting,
2608 wallclock_gpu_t *gputimes,
2609 gmx_bool bWriteStat)
2612 t_nrnb *nrnb_tot = NULL;
2615 double elapsed_time,
2616 elapsed_time_over_all_ranks,
2617 elapsed_time_over_all_threads,
2618 elapsed_time_over_all_threads_over_all_ranks;
2619 wallcycle_sum(cr, wcycle);
2625 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2626 cr->mpi_comm_mysim);
2634 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2635 elapsed_time_over_all_ranks = elapsed_time;
2636 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2637 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2641 /* reduce elapsed_time over all MPI ranks in the current simulation */
2642 MPI_Allreduce(&elapsed_time,
2643 &elapsed_time_over_all_ranks,
2644 1, MPI_DOUBLE, MPI_SUM,
2645 cr->mpi_comm_mysim);
2646 elapsed_time_over_all_ranks /= cr->nnodes;
2647 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2648 * current simulation. */
2649 MPI_Allreduce(&elapsed_time_over_all_threads,
2650 &elapsed_time_over_all_threads_over_all_ranks,
2651 1, MPI_DOUBLE, MPI_SUM,
2652 cr->mpi_comm_mysim);
2658 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2665 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2667 print_dd_statistics(cr, inputrec, fplog);
2672 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2673 elapsed_time_over_all_ranks,
2676 if (EI_DYNAMICS(inputrec->eI))
2678 delta_t = inputrec->delta_t;
2687 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2688 elapsed_time_over_all_ranks,
2689 walltime_accounting_get_nsteps_done(walltime_accounting),
2690 delta_t, nbfs, mflop);
2694 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2695 elapsed_time_over_all_ranks,
2696 walltime_accounting_get_nsteps_done(walltime_accounting),
2697 delta_t, nbfs, mflop);
2702 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2704 /* this function works, but could probably use a logic rewrite to keep all the different
2705 types of efep straight. */
2708 t_lambda *fep = ir->fepvals;
2710 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2712 for (i = 0; i < efptNR; i++)
2724 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2725 if checkpoint is set -- a kludge is in for now
2727 for (i = 0; i < efptNR; i++)
2729 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2730 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2732 lambda[i] = fep->init_lambda;
2735 lam0[i] = lambda[i];
2740 lambda[i] = fep->all_lambda[i][*fep_state];
2743 lam0[i] = lambda[i];
2749 /* need to rescale control temperatures to match current state */
2750 for (i = 0; i < ir->opts.ngtc; i++)
2752 if (ir->opts.ref_t[i] > 0)
2754 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2760 /* Send to the log the information on the current lambdas */
2763 fprintf(fplog, "Initial vector of lambda components:[ ");
2764 for (i = 0; i < efptNR; i++)
2766 fprintf(fplog, "%10.4f ", lambda[i]);
2768 fprintf(fplog, "]\n");
2774 void init_md(FILE *fplog,
2775 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2776 double *t, double *t0,
2777 real *lambda, int *fep_state, double *lam0,
2778 t_nrnb *nrnb, gmx_mtop_t *mtop,
2780 int nfile, const t_filenm fnm[],
2781 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2782 tensor force_vir, tensor shake_vir, rvec mu_tot,
2783 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags)
2788 /* Initial values */
2789 *t = *t0 = ir->init_t;
2792 for (i = 0; i < ir->opts.ngtc; i++)
2794 /* set bSimAnn if any group is being annealed */
2795 if (ir->opts.annealing[i] != eannNO)
2802 update_annealing_target_temp(&(ir->opts), ir->init_t);
2805 /* Initialize lambda variables */
2806 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2810 *upd = init_update(ir);
2816 *vcm = init_vcm(fplog, &mtop->groups, ir);
2819 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2821 if (ir->etc == etcBERENDSEN)
2823 please_cite(fplog, "Berendsen84a");
2825 if (ir->etc == etcVRESCALE)
2827 please_cite(fplog, "Bussi2007a");
2835 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv);
2837 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2838 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2843 please_cite(fplog, "Fritsch12");
2844 please_cite(fplog, "Junghans10");
2846 /* Initiate variables */
2847 clear_mat(force_vir);
2848 clear_mat(shake_vir);