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, 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.
51 #include "nonbonded.h"
66 #include "gmx_omp_nthreads.h"
68 #include "gromacs/timing/wallcycle.h"
79 gmx_bool bDoLongRangeNS)
85 if (!fr->ns.nblist_initialized)
87 init_neighbor_list(fp, fr, md->homenr);
95 nsearch = search_neighbours(fp, fr, box, top, groups, cr, nrnb, md,
96 bFillGrid, bDoLongRangeNS);
99 fprintf(debug, "nsearch = %d\n", nsearch);
102 /* Check whether we have to do dynamic load balancing */
103 /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
104 count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
105 &(top->idef),opts->ngener);
107 if (fr->ns.dump_nl > 0)
109 dump_nblist(fp, cr, fr, fr->ns.dump_nl);
113 static void reduce_thread_forces(int n, rvec *f,
116 int efpt_ind, real *dvdl,
117 int nthreads, f_thread_t *f_t)
121 /* This reduction can run over any number of threads */
122 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntBonded)) private(t) schedule(static)
123 for (i = 0; i < n; i++)
125 for (t = 1; t < nthreads; t++)
127 rvec_inc(f[i], f_t[t].f[i]);
130 for (t = 1; t < nthreads; t++)
132 *Vcorr += f_t[t].Vcorr;
133 *dvdl += f_t[t].dvdl[efpt_ind];
134 m_add(vir, f_t[t].vir, vir);
138 void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda)
140 fprintf(fplog, " %-30s V %12.5e dVdl %12.5e\n", s, v, dvdlambda);
143 void do_force_lowlevel(FILE *fplog, gmx_large_int_t step,
144 t_forcerec *fr, t_inputrec *ir,
145 t_idef *idef, t_commrec *cr,
146 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
148 rvec x[], history_t *hist,
151 gmx_enerdata_t *enerd,
168 gmx_bool bDoEpot, bSepDVDL, bSB;
174 double clam_i, vlam_i;
175 real dvdl_dum[efptNR], dvdl_nb[efptNR], lam_i[efptNR];
178 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
181 #define PRINT_SEPDVDL(s, v, dvdlambda) if (bSepDVDL) { gmx_print_sepdvdl(fplog, s, v, dvdlambda); }
183 set_pbc(&pbc, fr->ePBC, box);
185 /* reset free energy components */
186 for (i = 0; i < efptNR; i++)
193 for (i = 0; (i < DIM); i++)
195 box_size[i] = box[i][i];
198 bSepDVDL = (fr->bSepDVDL && do_per_step(step, ir->nstlog));
201 /* do QMMM first if requested */
204 enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr);
209 fprintf(fplog, "Step %s: non-bonded V and dVdl for node %d:\n",
210 gmx_step_str(step, buf), cr->nodeid);
213 /* Call the short range functions all in one go. */
216 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
217 #define TAKETIME FALSE
220 MPI_Barrier(cr->mpi_comm_mygroup);
227 /* foreign lambda component for walls */
228 real dvdl_walls = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
229 enerd->grpp.ener[egLJSR], nrnb);
230 PRINT_SEPDVDL("Walls", 0.0, dvdl_walls);
231 enerd->dvdl_lin[efptVDW] += dvdl_walls;
234 /* If doing GB, reset dvda and calculate the Born radii */
235 if (ir->implicit_solvent)
237 wallcycle_sub_start(wcycle, ewcsNONBONDED);
239 for (i = 0; i < born->nr; i++)
246 calc_gb_rad(cr, fr, ir, top, x, &(fr->gblist), born, md, nrnb);
249 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
253 /* We only do non-bonded calculation with group scheme here, the verlet
254 * calls are done from do_force_cutsVERLET(). */
255 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
258 /* Add short-range interactions */
259 donb_flags |= GMX_NONBONDED_DO_SR;
261 if (flags & GMX_FORCE_FORCES)
263 donb_flags |= GMX_NONBONDED_DO_FORCE;
265 if (flags & GMX_FORCE_ENERGY)
267 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
269 if (flags & GMX_FORCE_DO_LR)
271 donb_flags |= GMX_NONBONDED_DO_LR;
274 wallcycle_sub_start(wcycle, ewcsNONBONDED);
275 do_nonbonded(fr, x, f, f_longrange, md, excl,
277 lambda, dvdl_nb, -1, -1, donb_flags);
279 /* If we do foreign lambda and we have soft-core interactions
280 * we have to recalculate the (non-linear) energies contributions.
282 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
284 for (i = 0; i < enerd->n_lambda; i++)
286 for (j = 0; j < efptNR; j++)
288 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
290 reset_foreign_enerdata(enerd);
291 do_nonbonded(fr, x, f, f_longrange, md, excl,
292 &(enerd->foreign_grpp), nrnb,
293 lam_i, dvdl_dum, -1, -1,
294 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
295 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
296 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
299 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
303 /* If we are doing GB, calculate bonded forces and apply corrections
304 * to the solvation forces */
305 /* MRS: Eventually, many need to include free energy contribution here! */
306 if (ir->implicit_solvent)
308 wallcycle_sub_start(wcycle, ewcsBONDED);
309 calc_gb_forces(cr, md, born, top, x, f, fr, idef,
310 ir->gb_algorithm, ir->sa_algorithm, nrnb, &pbc, graph, enerd);
311 wallcycle_sub_stop(wcycle, ewcsBONDED);
322 if (fepvals->sc_alpha != 0)
324 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
328 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
331 if (fepvals->sc_alpha != 0)
333 /* even though coulomb part is linear, we already added it, beacuse we
334 need to go through the vdw calculation anyway */
336 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
340 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
345 real V_short_range = 0;
346 real dvdl_short_range;
348 for (i = 0; i < enerd->grpp.nener; i++)
352 enerd->grpp.ener[egBHAMSR][i] :
353 enerd->grpp.ener[egLJSR][i])
354 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
356 dvdl_short_range = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
357 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",
366 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
369 /* Shift the coordinates. Must be done before bonded forces and PPPM,
370 * but is also necessary for SHAKE and update, therefore it can NOT
371 * go when no bonded forces have to be evaluated.
374 /* Here sometimes we would not need to shift with NBFonly,
375 * but we do so anyhow for consistency of the returned coordinates.
379 shift_self(graph, box, x);
382 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
386 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
389 /* Check whether we need to do bondeds or correct for exclusions */
391 ((flags & GMX_FORCE_BONDED)
392 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
394 /* Since all atoms are in the rectangular or triclinic unit-cell,
395 * only single box vector shifts (2 in x) are required.
397 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
401 if (flags & GMX_FORCE_BONDED)
403 wallcycle_sub_start(wcycle, ewcsBONDED);
404 calc_bonds(fplog, cr->ms,
405 idef, x, hist, f, fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
406 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
408 fr->bSepDVDL && do_per_step(step, ir->nstlog), step);
410 /* Check if we have to determine energy differences
411 * at foreign lambda's.
413 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
414 idef->ilsort != ilsortNO_FE)
416 if (idef->ilsort != ilsortFE_SORTED)
418 gmx_incons("The bonded interactions are not sorted for free energy");
420 for (i = 0; i < enerd->n_lambda; i++)
422 reset_foreign_enerdata(enerd);
423 for (j = 0; j < efptNR; j++)
425 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
427 calc_bonds_lambda(fplog, idef, x, fr, &pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
428 fcd, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
429 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
430 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
435 wallcycle_sub_stop(wcycle, ewcsBONDED);
441 if (EEL_FULL(fr->eeltype))
443 real Vlr = 0, Vcorr = 0;
444 real dvdl_long_range = 0;
447 bSB = (ir->nwall == 2);
451 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
452 box_size[ZZ] *= ir->wall_ewald_zfac;
455 clear_mat(fr->vir_el_recip);
459 real dvdl_long_range_correction = 0;
461 /* With the Verlet scheme exclusion forces are calculated
462 * in the non-bonded kernel.
464 /* The TPI molecule does not have exclusions with the rest
465 * of the system and no intra-molecular PME grid contributions
466 * will be calculated in gmx_pme_calc_energy.
468 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
469 ir->ewald_geometry != eewg3D ||
470 ir->epsilon_surface != 0)
474 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
478 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
481 nthreads = gmx_omp_nthreads_get(emntBonded);
482 #pragma omp parallel for num_threads(nthreads) schedule(static)
483 for (t = 0; t < nthreads; t++)
488 real *Vcorrt, *dvdlt;
491 fnv = fr->f_novirsum;
492 vir = &fr->vir_el_recip;
494 dvdlt = &dvdl_long_range_correction;
499 vir = &fr->f_t[t].vir;
500 Vcorrt = &fr->f_t[t].Vcorr;
501 dvdlt = &fr->f_t[t].dvdl[efptCOUL];
502 for (i = 0; i < fr->natoms_force; i++)
510 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
513 md->nChargePerturbed ? md->chargeB : NULL,
514 ir->cutoff_scheme != ecutsVERLET,
515 excl, x, bSB ? boxs : box, mu_tot,
519 lambda[efptCOUL], dvdlt);
523 reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
525 &Vcorr, efptCOUL, &dvdl_long_range_correction,
529 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
534 Vcorr += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
535 &dvdl_long_range_correction,
539 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.", Vcorr,
540 dvdl_long_range_correction);
541 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_correction;
549 case eelPMEUSERSWITCH:
551 if (cr->duty & DUTY_PME)
553 assert(fr->n_tpi >= 0);
554 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
556 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
557 if (flags & GMX_FORCE_FORCES)
559 pme_flags |= GMX_PME_CALC_F;
561 if (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY))
563 pme_flags |= GMX_PME_CALC_ENER_VIR;
567 /* We don't calculate f, but we do want the potential */
568 pme_flags |= GMX_PME_CALC_POT;
570 wallcycle_start(wcycle, ewcPMEMESH);
571 status = gmx_pme_do(fr->pmedata,
572 md->start, md->homenr - fr->n_tpi,
574 md->chargeA, md->chargeB,
575 bSB ? boxs : box, cr,
576 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
577 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
579 fr->vir_el_recip, fr->ewaldcoeff,
580 &Vlr, lambda[efptCOUL], &dvdl_long_range,
582 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
584 /* We should try to do as little computation after
585 * this as possible, because parallel PME synchronizes
586 * the nodes, so we want all load imbalance of the rest
587 * of the force calculation to be before the PME call.
588 * DD load balancing is done on the whole time of
589 * the force call (without PME).
594 /* Determine the PME grid energy of the test molecule
595 * with the PME grid potential of the other charges.
597 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
598 x + md->homenr - fr->n_tpi,
599 md->chargeA + md->homenr - fr->n_tpi,
602 PRINT_SEPDVDL("PME mesh", Vlr, dvdl_long_range);
606 Vlr = do_ewald(ir, x, fr->f_novirsum,
607 md->chargeA, md->chargeB,
608 box_size, cr, md->homenr,
609 fr->vir_el_recip, fr->ewaldcoeff,
610 lambda[efptCOUL], &dvdl_long_range, fr->ewald_table);
611 PRINT_SEPDVDL("Ewald long-range", Vlr, dvdl_long_range);
614 gmx_fatal(FARGS, "No such electrostatics method implemented %s",
615 eel_names[fr->eeltype]);
619 gmx_fatal(FARGS, "Error %d in long range electrostatics routine %s",
620 status, EELTYPE(fr->eeltype));
622 /* Note that with separate PME nodes we get the real energies later */
623 enerd->dvdl_lin[efptCOUL] += dvdl_long_range;
624 enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
627 fprintf(debug, "Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
628 Vlr, Vcorr, enerd->term[F_COUL_RECIP]);
629 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
630 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
635 /* Is there a reaction-field exclusion correction needed? */
636 if (EEL_RF(fr->eeltype) && eelRF_NEC != fr->eeltype)
638 /* With the Verlet scheme, exclusion forces are calculated
639 * in the non-bonded kernel.
641 if (ir->cutoff_scheme != ecutsVERLET)
643 real dvdl_rf_excl = 0;
644 enerd->term[F_RF_EXCL] =
645 RF_excl_correction(fr, graph, md, excl, x, f,
646 fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
648 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
649 PRINT_SEPDVDL("RF exclusion correction",
650 enerd->term[F_RF_EXCL], dvdl_rf_excl);
659 print_nrnb(debug, nrnb);
667 MPI_Barrier(cr->mpi_comm_mygroup);
670 if (fr->timesteps == 11)
672 fprintf(stderr, "* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
673 cr->nodeid, gmx_step_str(fr->timesteps, buf),
674 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
675 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
683 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
688 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
692 for (i = 0; i < F_NRE; i++)
695 enerd->foreign_term[i] = 0;
699 for (i = 0; i < efptNR; i++)
701 enerd->dvdl_lin[i] = 0;
702 enerd->dvdl_nonlin[i] = 0;
708 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
710 enerd->grpp.nener = n2;
711 enerd->foreign_grpp.nener = n2;
712 for (i = 0; (i < egNR); i++)
714 snew(enerd->grpp.ener[i], n2);
715 snew(enerd->foreign_grpp.ener[i], n2);
720 enerd->n_lambda = 1 + n_lambda;
721 snew(enerd->enerpart_lambda, enerd->n_lambda);
729 void destroy_enerdata(gmx_enerdata_t *enerd)
733 for (i = 0; (i < egNR); i++)
735 sfree(enerd->grpp.ener[i]);
738 for (i = 0; (i < egNR); i++)
740 sfree(enerd->foreign_grpp.ener[i]);
745 sfree(enerd->enerpart_lambda);
749 static real sum_v(int n, real v[])
755 for (i = 0; (i < n); i++)
763 void sum_epot(gmx_grppairener_t *grpp, real *epot)
767 /* Accumulate energies */
768 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
769 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
770 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
771 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
772 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
773 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
774 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
775 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
777 /* lattice part of LR doesnt belong to any group
778 * and has been added earlier
780 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
781 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
784 for (i = 0; (i < F_EPOT); i++)
786 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
788 epot[F_EPOT] += epot[i];
793 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
798 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
799 enerd->term[F_DVDL] = 0.0;
800 for (i = 0; i < efptNR; i++)
802 if (fepvals->separate_dvdl[i])
804 /* could this be done more readably/compactly? */
817 index = F_DVDL_BONDED;
819 case (efptRESTRAINT):
820 index = F_DVDL_RESTRAINT;
826 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
829 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
830 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
835 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
838 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
839 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
844 /* Notes on the foreign lambda free energy difference evaluation:
845 * Adding the potential and ekin terms that depend linearly on lambda
846 * as delta lam * dvdl to the energy differences is exact.
847 * For the constraints this is not exact, but we have no other option
848 * without literally changing the lengths and reevaluating the energies at each step.
849 * (try to remedy this post 4.6 - MRS)
850 * For the non-bonded LR term we assume that the soft-core (if present)
851 * no longer affects the energy beyond the short-range cut-off,
852 * which is a very good approximation (except for exotic settings).
853 * (investigate how to overcome this post 4.6 - MRS)
855 if (fepvals->separate_dvdl[efptBONDED])
857 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
861 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
863 enerd->term[F_DVDL_CONSTR] = 0;
865 for (i = 0; i < fepvals->n_lambda; i++)
867 /* note we are iterating over fepvals here!
868 For the current lam, dlam = 0 automatically,
869 so we don't need to add anything to the
870 enerd->enerpart_lambda[0] */
872 /* we don't need to worry about dvdl_lin contributions to dE at
873 current lambda, because the contributions to the current
874 lambda are automatically zeroed */
876 for (j = 0; j < efptNR; j++)
878 /* Note that this loop is over all dhdl components, not just the separated ones */
879 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
880 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
883 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
884 fepvals->all_lambda[j][i], efpt_names[j],
885 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
886 dlam, enerd->dvdl_lin[j]);
893 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
897 /* First reset all foreign energy components. Foreign energies always called on
898 neighbor search steps */
899 for (i = 0; (i < egNR); i++)
901 for (j = 0; (j < enerd->grpp.nener); j++)
903 enerd->foreign_grpp.ener[i][j] = 0.0;
907 /* potential energy components */
908 for (i = 0; (i <= F_EPOT); i++)
910 enerd->foreign_term[i] = 0.0;
914 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
915 gmx_enerdata_t *enerd,
921 /* First reset all energy components, except for the long range terms
922 * on the master at non neighbor search steps, since the long range
923 * terms have already been summed at the last neighbor search step.
925 bKeepLR = (fr->bTwinRange && !bNS);
926 for (i = 0; (i < egNR); i++)
928 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
930 for (j = 0; (j < enerd->grpp.nener); j++)
932 enerd->grpp.ener[i][j] = 0.0;
936 for (i = 0; i < efptNR; i++)
938 enerd->dvdl_lin[i] = 0.0;
939 enerd->dvdl_nonlin[i] = 0.0;
942 /* Normal potential energy components */
943 for (i = 0; (i <= F_EPOT); i++)
945 enerd->term[i] = 0.0;
947 /* Initialize the dVdlambda term with the long range contribution */
948 /* Initialize the dvdl term with the long range contribution */
949 enerd->term[F_DVDL] = 0.0;
950 enerd->term[F_DVDL_COUL] = 0.0;
951 enerd->term[F_DVDL_VDW] = 0.0;
952 enerd->term[F_DVDL_BONDED] = 0.0;
953 enerd->term[F_DVDL_RESTRAINT] = 0.0;
954 enerd->term[F_DKDL] = 0.0;
955 if (enerd->n_lambda > 0)
957 for (i = 0; i < enerd->n_lambda; i++)
959 enerd->enerpart_lambda[i] = 0.0;
962 /* reset foreign energy data - separate function since we also call it elsewhere */
963 reset_foreign_enerdata(enerd);