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,
163 const DDBalanceRegionHandler &ddBalanceRegionHandler)
169 real dvdl_dum[efptNR], dvdl_nb[efptNR];
172 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
175 set_pbc(&pbc, fr->ePBC, box);
177 /* reset free energy components */
178 for (i = 0; i < efptNR; i++)
184 /* do QMMM first if requested */
187 enerd->term[F_EQM] = calculate_QMMM(cr, forceForUseWithShiftForces, fr);
190 /* Call the short range functions all in one go. */
193 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
194 #define TAKETIME FALSE
197 MPI_Barrier(cr->mpi_comm_mygroup);
204 /* foreign lambda component for walls */
205 real dvdl_walls = do_walls(*ir, *fr, box, *md, x,
206 forceWithVirial, lambda[efptVDW],
207 enerd->grpp.ener[egLJSR], nrnb);
208 enerd->dvdl_lin[efptVDW] += dvdl_walls;
211 /* We only do non-bonded calculation with group scheme here, the verlet
212 * calls are done from do_force_cutsVERLET(). */
213 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
216 /* Add short-range interactions */
217 donb_flags |= GMX_NONBONDED_DO_SR;
219 /* Currently all group scheme kernels always calculate (shift-)forces */
220 if (flags & GMX_FORCE_FORCES)
222 donb_flags |= GMX_NONBONDED_DO_FORCE;
224 if (flags & GMX_FORCE_VIRIAL)
226 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
228 if (flags & GMX_FORCE_ENERGY)
230 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
233 wallcycle_sub_start(wcycle, ewcsNONBONDED);
234 do_nonbonded(fr, x, forceForUseWithShiftForces, md, excl,
236 lambda, dvdl_nb, -1, -1, donb_flags);
238 /* If we do foreign lambda and we have soft-core interactions
239 * we have to recalculate the (non-linear) energies contributions.
241 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
243 for (i = 0; i < enerd->n_lambda; i++)
247 for (j = 0; j < efptNR; j++)
249 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
251 reset_foreign_enerdata(enerd);
252 do_nonbonded(fr, x, forceForUseWithShiftForces, md, excl,
253 &(enerd->foreign_grpp), nrnb,
254 lam_i, dvdl_dum, -1, -1,
255 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
256 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
257 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
260 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
271 if (fepvals->sc_alpha != 0)
273 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
277 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
280 if (fepvals->sc_alpha != 0)
282 /* even though coulomb part is linear, we already added it, beacuse we
283 need to go through the vdw calculation anyway */
285 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
289 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
294 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
297 /* Shift the coordinates. Must be done before listed forces and PPPM,
298 * but is also necessary for SHAKE and update, therefore it can NOT
299 * go when no listed forces have to be evaluated.
301 * The shifting and PBC code is deliberately not timed, since with
302 * the Verlet scheme it only takes non-zero time with triclinic
303 * boxes, and even then the time is around a factor of 100 less
304 * than the next smallest counter.
308 /* Here sometimes we would not need to shift with NBFonly,
309 * but we do so anyhow for consistency of the returned coordinates.
313 shift_self(graph, box, x);
316 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
320 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
323 /* Check whether we need to do listed interactions or correct for exclusions */
325 ((flags & GMX_FORCE_LISTED)
326 || EEL_RF(fr->ic->eeltype) || EEL_FULL(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype)))
328 /* TODO There are no electrostatics methods that require this
329 transformation, when using the Verlet scheme, so update the
330 above conditional. */
331 /* Since all atoms are in the rectangular or triclinic unit-cell,
332 * only single box vector shifts (2 in x) are required.
334 set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd->nc : nullptr,
338 do_force_listed(wcycle, box, ir->fepvals, cr, ms,
340 forceForUseWithShiftForces, forceWithVirial,
341 fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
342 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr,
348 /* Do long-range electrostatics and/or LJ-PME, including related short-range
351 if (EEL_FULL(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
354 real Vlr_q = 0, Vlr_lj = 0;
356 /* We reduce all virial, dV/dlambda and energy contributions, except
357 * for the reciprocal energies (Vlr_q, Vlr_lj) into the same struct.
359 ewald_corr_thread_t &ewaldOutput = fr->ewc_t[0];
360 clearEwaldThreadOutput(&ewaldOutput);
362 if (EEL_PME_EWALD(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
364 /* With the Verlet scheme exclusion forces are calculated
365 * in the non-bonded kernel.
367 /* The TPI molecule does not have exclusions with the rest
368 * of the system and no intra-molecular PME grid
369 * contributions will be calculated in
370 * gmx_pme_calc_energy.
372 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
373 ir->ewald_geometry != eewg3D ||
374 ir->epsilon_surface != 0)
378 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
382 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
385 nthreads = fr->nthread_ewc;
386 #pragma omp parallel for num_threads(nthreads) schedule(static)
387 for (t = 0; t < nthreads; t++)
391 ewald_corr_thread_t &ewc_t = fr->ewc_t[t];
394 clearEwaldThreadOutput(&ewc_t);
397 /* Threading is only supported with the Verlet cut-off
398 * scheme and then only single particle forces (no
399 * exclusion forces) are calculated, so we can store
400 * the forces in the normal, single forceWithVirial->force_ array.
402 ewald_LRcorrection(md->homenr, cr, nthreads, t, fr, ir,
403 md->chargeA, md->chargeB,
404 md->sqrt_c6A, md->sqrt_c6B,
405 md->sigmaA, md->sigmaB,
406 md->sigma3A, md->sigma3B,
407 (md->nChargePerturbed != 0) || (md->nTypePerturbed != 0),
408 ir->cutoff_scheme != ecutsVERLET,
409 excl, x, box, mu_tot,
412 as_rvec_array(forceWithVirial->force_.data()),
413 ewc_t.vir_q, ewc_t.vir_lj,
414 &ewc_t.Vcorr_q, &ewc_t.Vcorr_lj,
415 lambda[efptCOUL], lambda[efptVDW],
416 &ewc_t.dvdl[efptCOUL], &ewc_t.dvdl[efptVDW]);
418 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
422 reduceEwaldThreadOuput(nthreads, fr->ewc_t);
424 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
427 if (EEL_PME_EWALD(fr->ic->eeltype) && fr->n_tpi == 0)
429 /* This is not in a subcounter because it takes a
430 negligible and constant-sized amount of time */
431 ewaldOutput.Vcorr_q +=
432 ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
433 &ewaldOutput.dvdl[efptCOUL],
437 if ((EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype)) &&
438 thisRankHasDuty(cr, DUTY_PME) && (pme_run_mode(fr->pmedata) == PmeRunMode::CPU))
440 /* Do reciprocal PME for Coulomb and/or LJ. */
441 assert(fr->n_tpi >= 0);
442 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
444 pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
446 if (flags & GMX_FORCE_FORCES)
448 pme_flags |= GMX_PME_CALC_F;
450 if (flags & GMX_FORCE_VIRIAL)
452 pme_flags |= GMX_PME_CALC_ENER_VIR;
456 /* We don't calculate f, but we do want the potential */
457 pme_flags |= GMX_PME_CALC_POT;
460 /* With domain decomposition we close the CPU side load
461 * balancing region here, because PME does global
462 * communication that acts as a global barrier.
464 ddBalanceRegionHandler.closeAfterForceComputationCpu();
466 wallcycle_start(wcycle, ewcPMEMESH);
467 status = gmx_pme_do(fr->pmedata,
468 0, md->homenr - fr->n_tpi,
470 as_rvec_array(forceWithVirial->force_.data()),
471 md->chargeA, md->chargeB,
472 md->sqrt_c6A, md->sqrt_c6B,
473 md->sigmaA, md->sigmaB,
475 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
476 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
478 ewaldOutput.vir_q, ewaldOutput.vir_lj,
480 lambda[efptCOUL], lambda[efptVDW],
481 &ewaldOutput.dvdl[efptCOUL],
482 &ewaldOutput.dvdl[efptVDW],
484 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
487 gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
490 /* We should try to do as little computation after
491 * this as possible, because parallel PME synchronizes
492 * the nodes, so we want all load imbalance of the
493 * rest of the force calculation to be before the PME
494 * call. DD load balancing is done on the whole time
495 * of the force call (without PME).
500 if (EVDW_PME(ir->vdwtype))
503 gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
505 /* Determine the PME grid energy of the test molecule
506 * with the PME grid potential of the other charges.
508 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
509 x + md->homenr - fr->n_tpi,
510 md->chargeA + md->homenr - fr->n_tpi,
516 if (!EEL_PME(fr->ic->eeltype) && EEL_PME_EWALD(fr->ic->eeltype))
518 Vlr_q = do_ewald(ir, x, as_rvec_array(forceWithVirial->force_.data()),
519 md->chargeA, md->chargeB,
521 ewaldOutput.vir_q, fr->ic->ewaldcoeff_q,
522 lambda[efptCOUL], &ewaldOutput.dvdl[efptCOUL],
526 /* Note that with separate PME nodes we get the real energies later */
527 // TODO it would be simpler if we just accumulated a single
528 // long-range virial contribution.
529 forceWithVirial->addVirialContribution(ewaldOutput.vir_q);
530 forceWithVirial->addVirialContribution(ewaldOutput.vir_lj);
531 enerd->dvdl_lin[efptCOUL] += ewaldOutput.dvdl[efptCOUL];
532 enerd->dvdl_lin[efptVDW] += ewaldOutput.dvdl[efptVDW];
533 enerd->term[F_COUL_RECIP] = Vlr_q + ewaldOutput.Vcorr_q;
534 enerd->term[F_LJ_RECIP] = Vlr_lj + ewaldOutput.Vcorr_lj;
538 fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
539 Vlr_q, ewaldOutput.Vcorr_q, enerd->term[F_COUL_RECIP]);
540 pr_rvecs(debug, 0, "vir_el_recip after corr", ewaldOutput.vir_q, DIM);
541 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
542 fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
543 Vlr_lj, ewaldOutput.Vcorr_lj, enerd->term[F_LJ_RECIP]);
544 pr_rvecs(debug, 0, "vir_lj_recip after corr", ewaldOutput.vir_lj, DIM);
549 /* Is there a reaction-field exclusion correction needed?
550 * With the Verlet scheme, exclusion forces are calculated
551 * in the non-bonded kernel.
553 if (ir->cutoff_scheme != ecutsVERLET && EEL_RF(fr->ic->eeltype))
555 real dvdl_rf_excl = 0;
556 enerd->term[F_RF_EXCL] =
557 RF_excl_correction(fr, graph, md, excl, DOMAINDECOMP(cr),
558 x, forceForUseWithShiftForces,
559 fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
561 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
567 print_nrnb(debug, nrnb);
574 MPI_Barrier(cr->mpi_comm_mygroup);
577 if (fr->timesteps == 11)
580 fprintf(stderr, "* PP load balancing info: rank %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
581 cr->nodeid, gmx_step_str(fr->timesteps, buf),
582 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
583 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
591 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
596 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
600 for (i = 0; i < F_NRE; i++)
603 enerd->foreign_term[i] = 0;
607 for (i = 0; i < efptNR; i++)
609 enerd->dvdl_lin[i] = 0;
610 enerd->dvdl_nonlin[i] = 0;
616 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
618 enerd->grpp.nener = n2;
619 enerd->foreign_grpp.nener = n2;
620 for (i = 0; (i < egNR); i++)
622 snew(enerd->grpp.ener[i], n2);
623 snew(enerd->foreign_grpp.ener[i], n2);
628 enerd->n_lambda = 1 + n_lambda;
629 snew(enerd->enerpart_lambda, enerd->n_lambda);
637 void destroy_enerdata(gmx_enerdata_t *enerd)
641 for (i = 0; (i < egNR); i++)
643 sfree(enerd->grpp.ener[i]);
646 for (i = 0; (i < egNR); i++)
648 sfree(enerd->foreign_grpp.ener[i]);
653 sfree(enerd->enerpart_lambda);
657 static real sum_v(int n, const real v[])
663 for (i = 0; (i < n); i++)
671 void sum_epot(gmx_grppairener_t *grpp, real *epot)
675 /* Accumulate energies */
676 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
677 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
678 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
679 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
681 /* lattice part of LR doesnt belong to any group
682 * and has been added earlier
684 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
687 for (i = 0; (i < F_EPOT); i++)
689 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
691 epot[F_EPOT] += epot[i];
696 void sum_dhdl(gmx_enerdata_t *enerd, gmx::ArrayRef<const real> lambda, t_lambda *fepvals)
700 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
701 enerd->term[F_DVDL] = 0.0;
702 for (int i = 0; i < efptNR; i++)
704 if (fepvals->separate_dvdl[i])
706 /* could this be done more readably/compactly? */
719 index = F_DVDL_BONDED;
721 case (efptRESTRAINT):
722 index = F_DVDL_RESTRAINT;
728 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
731 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
732 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
737 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
740 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
741 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
746 if (fepvals->separate_dvdl[efptBONDED])
748 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
752 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
755 for (int i = 0; i < fepvals->n_lambda; i++)
757 /* note we are iterating over fepvals here!
758 For the current lam, dlam = 0 automatically,
759 so we don't need to add anything to the
760 enerd->enerpart_lambda[0] */
762 /* we don't need to worry about dvdl_lin contributions to dE at
763 current lambda, because the contributions to the current
764 lambda are automatically zeroed */
766 double &enerpart_lambda = enerd->enerpart_lambda[i + 1];
768 for (gmx::index j = 0; j < lambda.ssize(); j++)
770 /* Note that this loop is over all dhdl components, not just the separated ones */
771 const double dlam = fepvals->all_lambda[j][i] - lambda[j];
773 enerpart_lambda += dlam*enerd->dvdl_lin[j];
775 /* Constraints can not be evaluated at foreign lambdas, so we add
776 * a linear extrapolation. This is an approximation, but usually
777 * quite accurate since constraints change little between lambdas.
779 if ((j == efptBONDED && fepvals->separate_dvdl[efptBONDED]) ||
780 (j == efptFEP && !fepvals->separate_dvdl[efptBONDED]))
782 enerpart_lambda += dlam*enerd->term[F_DVDL_CONSTR];
785 if (j == efptMASS && !fepvals->separate_dvdl[j])
787 enerpart_lambda += dlam*enerd->term[F_DKDL];
792 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
793 fepvals->all_lambda[j][i], efpt_names[j],
794 enerpart_lambda - enerd->enerpart_lambda[0],
795 dlam, enerd->dvdl_lin[j]);
800 /* The constrain contribution is now included in other terms, so clear it */
801 enerd->term[F_DVDL_CONSTR] = 0;
805 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
809 /* First reset all foreign energy components. Foreign energies always called on
810 neighbor search steps */
811 for (i = 0; (i < egNR); i++)
813 for (j = 0; (j < enerd->grpp.nener); j++)
815 enerd->foreign_grpp.ener[i][j] = 0.0;
819 /* potential energy components */
820 for (i = 0; (i <= F_EPOT); i++)
822 enerd->foreign_term[i] = 0.0;
826 void reset_enerdata(gmx_enerdata_t *enerd)
830 /* First reset all energy components. */
831 for (i = 0; (i < egNR); i++)
833 for (j = 0; (j < enerd->grpp.nener); j++)
835 enerd->grpp.ener[i][j] = 0.0;
838 for (i = 0; i < efptNR; i++)
840 enerd->dvdl_lin[i] = 0.0;
841 enerd->dvdl_nonlin[i] = 0.0;
844 /* Normal potential energy components */
845 for (i = 0; (i <= F_EPOT); i++)
847 enerd->term[i] = 0.0;
849 enerd->term[F_DVDL] = 0.0;
850 enerd->term[F_DVDL_COUL] = 0.0;
851 enerd->term[F_DVDL_VDW] = 0.0;
852 enerd->term[F_DVDL_BONDED] = 0.0;
853 enerd->term[F_DVDL_RESTRAINT] = 0.0;
854 enerd->term[F_DKDL] = 0.0;
855 if (enerd->n_lambda > 0)
857 for (i = 0; i < enerd->n_lambda; i++)
859 enerd->enerpart_lambda[i] = 0.0;
862 /* reset foreign energy data - separate function since we also call it elsewhere */
863 reset_foreign_enerdata(enerd);