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,2015,2016,2017,2018,2019, 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.
47 #include "gromacs/domdec/dlbtiming.h"
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/domdec/domdec_struct.h"
50 #include "gromacs/ewald/ewald.h"
51 #include "gromacs/ewald/long_range_correction.h"
52 #include "gromacs/ewald/pme.h"
53 #include "gromacs/gmxlib/network.h"
54 #include "gromacs/gmxlib/nrnb.h"
55 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
56 #include "gromacs/listed_forces/listed_forces.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/math/vecdump.h"
59 #include "gromacs/mdlib/force_flags.h"
60 #include "gromacs/mdlib/forcerec_threading.h"
61 #include "gromacs/mdlib/ns.h"
62 #include "gromacs/mdlib/qmmm.h"
63 #include "gromacs/mdlib/rf_util.h"
64 #include "gromacs/mdlib/wall.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/enerdata.h"
67 #include "gromacs/mdtypes/forceoutput.h"
68 #include "gromacs/mdtypes/forcerec.h"
69 #include "gromacs/mdtypes/inputrec.h"
70 #include "gromacs/mdtypes/md_enums.h"
71 #include "gromacs/pbcutil/ishift.h"
72 #include "gromacs/pbcutil/mshift.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/timing/wallcycle.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/exceptions.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/smalloc.h"
83 const gmx_groups_t *groups,
93 if (!fr->ns->nblist_initialized)
95 init_neighbor_list(fp, fr, md->homenr);
98 nsearch = search_neighbours(fp, fr, box, top, groups, cr, nrnb, md,
102 fprintf(debug, "nsearch = %d\n", nsearch);
105 /* Check whether we have to do dynamic load balancing */
106 /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
107 count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
108 &(top->idef),opts->ngener);
110 if (fr->ns->dump_nl > 0)
112 dump_nblist(fp, cr, fr, fr->ns->dump_nl);
116 static void clearEwaldThreadOutput(ewald_corr_thread_t *ewc_t)
120 ewc_t->dvdl[efptCOUL] = 0;
121 ewc_t->dvdl[efptVDW] = 0;
122 clear_mat(ewc_t->vir_q);
123 clear_mat(ewc_t->vir_lj);
126 static void reduceEwaldThreadOuput(int nthreads, ewald_corr_thread_t *ewc_t)
128 ewald_corr_thread_t &dest = ewc_t[0];
130 for (int t = 1; t < nthreads; t++)
132 dest.Vcorr_q += ewc_t[t].Vcorr_q;
133 dest.Vcorr_lj += ewc_t[t].Vcorr_lj;
134 dest.dvdl[efptCOUL] += ewc_t[t].dvdl[efptCOUL];
135 dest.dvdl[efptVDW] += ewc_t[t].dvdl[efptVDW];
136 m_add(dest.vir_q, ewc_t[t].vir_q, dest.vir_q);
137 m_add(dest.vir_lj, ewc_t[t].vir_lj, dest.vir_lj);
141 void do_force_lowlevel(t_forcerec *fr,
142 const t_inputrec *ir,
145 const gmx_multisim_t *ms,
147 gmx_wallcycle_t wcycle,
151 rvec *forceForUseWithShiftForces,
152 gmx::ForceWithVirial *forceWithVirial,
153 gmx_enerdata_t *enerd,
158 const t_graph *graph,
159 const t_blocka *excl,
162 const DDBalanceRegionHandler &ddBalanceRegionHandler)
167 real dvdl_dum[efptNR], dvdl_nb[efptNR];
170 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
173 /* reset free energy components */
174 for (i = 0; i < efptNR; i++)
180 /* do QMMM first if requested */
183 enerd->term[F_EQM] = calculate_QMMM(cr, forceForUseWithShiftForces, fr);
186 /* Call the short range functions all in one go. */
189 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
190 #define TAKETIME FALSE
193 MPI_Barrier(cr->mpi_comm_mygroup);
200 /* foreign lambda component for walls */
201 real dvdl_walls = do_walls(*ir, *fr, box, *md, x,
202 forceWithVirial, lambda[efptVDW],
203 enerd->grpp.ener[egLJSR], nrnb);
204 enerd->dvdl_lin[efptVDW] += dvdl_walls;
207 /* We only do non-bonded calculation with group scheme here, the verlet
208 * calls are done from do_force_cutsVERLET(). */
209 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
212 /* Add short-range interactions */
213 donb_flags |= GMX_NONBONDED_DO_SR;
215 /* Currently all group scheme kernels always calculate (shift-)forces */
216 if (flags & GMX_FORCE_FORCES)
218 donb_flags |= GMX_NONBONDED_DO_FORCE;
220 if (flags & GMX_FORCE_VIRIAL)
222 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
224 if (flags & GMX_FORCE_ENERGY)
226 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
229 wallcycle_sub_start(wcycle, ewcsNONBONDED);
230 do_nonbonded(fr, x, forceForUseWithShiftForces, md, excl,
232 lambda, dvdl_nb, -1, -1, donb_flags);
234 /* If we do foreign lambda and we have soft-core interactions
235 * we have to recalculate the (non-linear) energies contributions.
237 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
239 for (i = 0; i < enerd->n_lambda; i++)
243 for (j = 0; j < efptNR; j++)
245 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
247 reset_foreign_enerdata(enerd);
248 do_nonbonded(fr, x, forceForUseWithShiftForces, md, excl,
249 &(enerd->foreign_grpp), nrnb,
250 lam_i, dvdl_dum, -1, -1,
251 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
252 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
253 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
256 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
267 if (fepvals->sc_alpha != 0)
269 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
273 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
276 if (fepvals->sc_alpha != 0)
278 /* even though coulomb part is linear, we already added it, beacuse we
279 need to go through the vdw calculation anyway */
281 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
285 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
290 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
293 /* Shift the coordinates. Must be done before listed forces and PPPM,
294 * but is also necessary for SHAKE and update, therefore it can NOT
295 * go when no listed forces have to be evaluated.
297 * The shifting and PBC code is deliberately not timed, since with
298 * the Verlet scheme it only takes non-zero time with triclinic
299 * boxes, and even then the time is around a factor of 100 less
300 * than the next smallest counter.
304 /* Here sometimes we would not need to shift with NBFonly,
305 * but we do so anyhow for consistency of the returned coordinates.
309 shift_self(graph, box, x);
312 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
316 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
321 const bool useRfWithGroupScheme = (fr->cutoff_scheme == ecutsGROUP) && EEL_RF(fr->ic->eeltype);
322 /* Check whether we need to take into account PBC in the following force tasks:
323 * listed interactions or when correcting for exclusions (for the Group scheme with RF) */
326 const auto needPbcForListedForces = bool(flags & GMX_FORCE_LISTED) && haveCpuListedForces(*fr, *idef, *fcd);
327 if (needPbcForListedForces || useRfWithGroupScheme)
329 /* Since all atoms are in the rectangular or triclinic unit-cell,
330 * only single box vector shifts (2 in x) are required.
332 set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd->nc : nullptr,
337 do_force_listed(wcycle, box, ir->fepvals, cr, ms,
339 forceForUseWithShiftForces, forceWithVirial,
340 fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
341 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr,
345 /* Do long-range electrostatics and/or LJ-PME, including related short-range
348 if (EEL_FULL(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
351 real Vlr_q = 0, Vlr_lj = 0;
353 /* We reduce all virial, dV/dlambda and energy contributions, except
354 * for the reciprocal energies (Vlr_q, Vlr_lj) into the same struct.
356 ewald_corr_thread_t &ewaldOutput = fr->ewc_t[0];
357 clearEwaldThreadOutput(&ewaldOutput);
359 if (EEL_PME_EWALD(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
361 /* With the Verlet scheme exclusion forces are calculated
362 * in the non-bonded kernel.
364 /* The TPI molecule does not have exclusions with the rest
365 * of the system and no intra-molecular PME grid
366 * contributions will be calculated in
367 * gmx_pme_calc_energy.
369 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
370 ir->ewald_geometry != eewg3D ||
371 ir->epsilon_surface != 0)
375 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
379 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
382 nthreads = fr->nthread_ewc;
383 #pragma omp parallel for num_threads(nthreads) schedule(static)
384 for (t = 0; t < nthreads; t++)
388 ewald_corr_thread_t &ewc_t = fr->ewc_t[t];
391 clearEwaldThreadOutput(&ewc_t);
394 /* Threading is only supported with the Verlet cut-off
395 * scheme and then only single particle forces (no
396 * exclusion forces) are calculated, so we can store
397 * the forces in the normal, single forceWithVirial->force_ array.
399 ewald_LRcorrection(md->homenr, cr, nthreads, t, fr, ir,
400 md->chargeA, md->chargeB,
401 md->sqrt_c6A, md->sqrt_c6B,
402 md->sigmaA, md->sigmaB,
403 md->sigma3A, md->sigma3B,
404 (md->nChargePerturbed != 0) || (md->nTypePerturbed != 0),
405 ir->cutoff_scheme != ecutsVERLET,
406 excl, x, box, mu_tot,
409 as_rvec_array(forceWithVirial->force_.data()),
410 ewc_t.vir_q, ewc_t.vir_lj,
411 &ewc_t.Vcorr_q, &ewc_t.Vcorr_lj,
412 lambda[efptCOUL], lambda[efptVDW],
413 &ewc_t.dvdl[efptCOUL], &ewc_t.dvdl[efptVDW]);
415 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
419 reduceEwaldThreadOuput(nthreads, fr->ewc_t);
421 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
424 if (EEL_PME_EWALD(fr->ic->eeltype) && fr->n_tpi == 0)
426 /* This is not in a subcounter because it takes a
427 negligible and constant-sized amount of time */
428 ewaldOutput.Vcorr_q +=
429 ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
430 &ewaldOutput.dvdl[efptCOUL],
434 if ((EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype)) &&
435 thisRankHasDuty(cr, DUTY_PME) && (pme_run_mode(fr->pmedata) == PmeRunMode::CPU))
437 /* Do reciprocal PME for Coulomb and/or LJ. */
438 assert(fr->n_tpi >= 0);
439 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
441 pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
443 if (flags & GMX_FORCE_FORCES)
445 pme_flags |= GMX_PME_CALC_F;
447 if (flags & GMX_FORCE_VIRIAL)
449 pme_flags |= GMX_PME_CALC_ENER_VIR;
453 /* We don't calculate f, but we do want the potential */
454 pme_flags |= GMX_PME_CALC_POT;
457 /* With domain decomposition we close the CPU side load
458 * balancing region here, because PME does global
459 * communication that acts as a global barrier.
461 ddBalanceRegionHandler.closeAfterForceComputationCpu();
463 wallcycle_start(wcycle, ewcPMEMESH);
464 status = gmx_pme_do(fr->pmedata,
465 0, md->homenr - fr->n_tpi,
467 as_rvec_array(forceWithVirial->force_.data()),
468 md->chargeA, md->chargeB,
469 md->sqrt_c6A, md->sqrt_c6B,
470 md->sigmaA, md->sigmaB,
472 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
473 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
475 ewaldOutput.vir_q, ewaldOutput.vir_lj,
477 lambda[efptCOUL], lambda[efptVDW],
478 &ewaldOutput.dvdl[efptCOUL],
479 &ewaldOutput.dvdl[efptVDW],
481 wallcycle_stop(wcycle, ewcPMEMESH);
484 gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
487 /* We should try to do as little computation after
488 * this as possible, because parallel PME synchronizes
489 * the nodes, so we want all load imbalance of the
490 * rest of the force calculation to be before the PME
491 * call. DD load balancing is done on the whole time
492 * of the force call (without PME).
497 if (EVDW_PME(ir->vdwtype))
500 gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
502 /* Determine the PME grid energy of the test molecule
503 * with the PME grid potential of the other charges.
505 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
506 x + md->homenr - fr->n_tpi,
507 md->chargeA + md->homenr - fr->n_tpi,
513 if (!EEL_PME(fr->ic->eeltype) && EEL_PME_EWALD(fr->ic->eeltype))
515 Vlr_q = do_ewald(ir, x, as_rvec_array(forceWithVirial->force_.data()),
516 md->chargeA, md->chargeB,
518 ewaldOutput.vir_q, fr->ic->ewaldcoeff_q,
519 lambda[efptCOUL], &ewaldOutput.dvdl[efptCOUL],
523 /* Note that with separate PME nodes we get the real energies later */
524 // TODO it would be simpler if we just accumulated a single
525 // long-range virial contribution.
526 forceWithVirial->addVirialContribution(ewaldOutput.vir_q);
527 forceWithVirial->addVirialContribution(ewaldOutput.vir_lj);
528 enerd->dvdl_lin[efptCOUL] += ewaldOutput.dvdl[efptCOUL];
529 enerd->dvdl_lin[efptVDW] += ewaldOutput.dvdl[efptVDW];
530 enerd->term[F_COUL_RECIP] = Vlr_q + ewaldOutput.Vcorr_q;
531 enerd->term[F_LJ_RECIP] = Vlr_lj + ewaldOutput.Vcorr_lj;
535 fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
536 Vlr_q, ewaldOutput.Vcorr_q, enerd->term[F_COUL_RECIP]);
537 pr_rvecs(debug, 0, "vir_el_recip after corr", ewaldOutput.vir_q, DIM);
538 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
539 fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
540 Vlr_lj, ewaldOutput.Vcorr_lj, enerd->term[F_LJ_RECIP]);
541 pr_rvecs(debug, 0, "vir_lj_recip after corr", ewaldOutput.vir_lj, DIM);
546 /* Is there a reaction-field exclusion correction needed?
547 * With the Verlet scheme, exclusion forces are calculated
548 * in the non-bonded kernel.
550 if (useRfWithGroupScheme)
552 real dvdl_rf_excl = 0;
553 enerd->term[F_RF_EXCL] =
554 RF_excl_correction(fr, graph, md, excl, DOMAINDECOMP(cr),
555 x, forceForUseWithShiftForces,
556 fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
558 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
564 print_nrnb(debug, nrnb);
571 MPI_Barrier(cr->mpi_comm_mygroup);
574 if (fr->timesteps == 11)
577 fprintf(stderr, "* PP load balancing info: rank %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
578 cr->nodeid, gmx_step_str(fr->timesteps, buf),
579 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
580 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
588 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
593 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
597 for (i = 0; i < F_NRE; i++)
600 enerd->foreign_term[i] = 0;
604 for (i = 0; i < efptNR; i++)
606 enerd->dvdl_lin[i] = 0;
607 enerd->dvdl_nonlin[i] = 0;
613 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
615 enerd->grpp.nener = n2;
616 enerd->foreign_grpp.nener = n2;
617 for (i = 0; (i < egNR); i++)
619 snew(enerd->grpp.ener[i], n2);
620 snew(enerd->foreign_grpp.ener[i], n2);
625 enerd->n_lambda = 1 + n_lambda;
626 snew(enerd->enerpart_lambda, enerd->n_lambda);
634 void destroy_enerdata(gmx_enerdata_t *enerd)
638 for (i = 0; (i < egNR); i++)
640 sfree(enerd->grpp.ener[i]);
643 for (i = 0; (i < egNR); i++)
645 sfree(enerd->foreign_grpp.ener[i]);
650 sfree(enerd->enerpart_lambda);
654 static real sum_v(int n, const real v[])
660 for (i = 0; (i < n); i++)
668 void sum_epot(gmx_grppairener_t *grpp, real *epot)
672 /* Accumulate energies */
673 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
674 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
675 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
676 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
678 /* lattice part of LR doesnt belong to any group
679 * and has been added earlier
681 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
684 for (i = 0; (i < F_EPOT); i++)
686 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
688 epot[F_EPOT] += epot[i];
693 void sum_dhdl(gmx_enerdata_t *enerd, gmx::ArrayRef<const real> lambda, t_lambda *fepvals)
697 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
698 enerd->term[F_DVDL] = 0.0;
699 for (int i = 0; i < efptNR; i++)
701 if (fepvals->separate_dvdl[i])
703 /* could this be done more readably/compactly? */
716 index = F_DVDL_BONDED;
718 case (efptRESTRAINT):
719 index = F_DVDL_RESTRAINT;
725 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
728 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
729 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
734 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
737 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
738 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
743 if (fepvals->separate_dvdl[efptBONDED])
745 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
749 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
752 for (int i = 0; i < fepvals->n_lambda; i++)
754 /* note we are iterating over fepvals here!
755 For the current lam, dlam = 0 automatically,
756 so we don't need to add anything to the
757 enerd->enerpart_lambda[0] */
759 /* we don't need to worry about dvdl_lin contributions to dE at
760 current lambda, because the contributions to the current
761 lambda are automatically zeroed */
763 double &enerpart_lambda = enerd->enerpart_lambda[i + 1];
765 for (gmx::index j = 0; j < lambda.ssize(); j++)
767 /* Note that this loop is over all dhdl components, not just the separated ones */
768 const double dlam = fepvals->all_lambda[j][i] - lambda[j];
770 enerpart_lambda += dlam*enerd->dvdl_lin[j];
772 /* Constraints can not be evaluated at foreign lambdas, so we add
773 * a linear extrapolation. This is an approximation, but usually
774 * quite accurate since constraints change little between lambdas.
776 if ((j == efptBONDED && fepvals->separate_dvdl[efptBONDED]) ||
777 (j == efptFEP && !fepvals->separate_dvdl[efptBONDED]))
779 enerpart_lambda += dlam*enerd->term[F_DVDL_CONSTR];
782 if (j == efptMASS && !fepvals->separate_dvdl[j])
784 enerpart_lambda += dlam*enerd->term[F_DKDL];
789 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
790 fepvals->all_lambda[j][i], efpt_names[j],
791 enerpart_lambda - enerd->enerpart_lambda[0],
792 dlam, enerd->dvdl_lin[j]);
797 /* The constrain contribution is now included in other terms, so clear it */
798 enerd->term[F_DVDL_CONSTR] = 0;
802 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
806 /* First reset all foreign energy components. Foreign energies always called on
807 neighbor search steps */
808 for (i = 0; (i < egNR); i++)
810 for (j = 0; (j < enerd->grpp.nener); j++)
812 enerd->foreign_grpp.ener[i][j] = 0.0;
816 /* potential energy components */
817 for (i = 0; (i <= F_EPOT); i++)
819 enerd->foreign_term[i] = 0.0;
823 void reset_enerdata(gmx_enerdata_t *enerd)
827 /* First reset all energy components. */
828 for (i = 0; (i < egNR); i++)
830 for (j = 0; (j < enerd->grpp.nener); j++)
832 enerd->grpp.ener[i][j] = 0.0;
835 for (i = 0; i < efptNR; i++)
837 enerd->dvdl_lin[i] = 0.0;
838 enerd->dvdl_nonlin[i] = 0.0;
841 /* Normal potential energy components */
842 for (i = 0; (i <= F_EPOT); i++)
844 enerd->term[i] = 0.0;
846 enerd->term[F_DVDL] = 0.0;
847 enerd->term[F_DVDL_COUL] = 0.0;
848 enerd->term[F_DVDL_VDW] = 0.0;
849 enerd->term[F_DVDL_BONDED] = 0.0;
850 enerd->term[F_DVDL_RESTRAINT] = 0.0;
851 enerd->term[F_DKDL] = 0.0;
852 if (enerd->n_lambda > 0)
854 for (i = 0; i < enerd->n_lambda; i++)
856 enerd->enerpart_lambda[i] = 0.0;
859 /* reset foreign energy data - separate function since we also call it elsewhere */
860 reset_foreign_enerdata(enerd);