1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
50 #include "nonbonded.h"
65 #include "gmx_omp_nthreads.h"
67 #include "gromacs/timing/wallcycle.h"
78 gmx_bool bDoLongRangeNS)
84 if (!fr->ns.nblist_initialized)
86 init_neighbor_list(fp, fr, md->homenr);
94 nsearch = search_neighbours(fp, fr, box, top, groups, cr, nrnb, md,
95 bFillGrid, bDoLongRangeNS);
98 fprintf(debug, "nsearch = %d\n", nsearch);
101 /* Check whether we have to do dynamic load balancing */
102 /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
103 count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
104 &(top->idef),opts->ngener);
106 if (fr->ns.dump_nl > 0)
108 dump_nblist(fp, cr, fr, fr->ns.dump_nl);
112 static void reduce_thread_forces(int n, rvec *f,
115 int efpt_ind, real *dvdl,
116 int nthreads, f_thread_t *f_t)
120 /* This reduction can run over any number of threads */
121 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntBonded)) private(t) schedule(static)
122 for (i = 0; i < n; i++)
124 for (t = 1; t < nthreads; t++)
126 rvec_inc(f[i], f_t[t].f[i]);
129 for (t = 1; t < nthreads; t++)
131 *Vcorr += f_t[t].Vcorr;
132 *dvdl += f_t[t].dvdl[efpt_ind];
133 m_add(vir, f_t[t].vir, vir);
137 void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda)
139 fprintf(fplog, " %-30s V %12.5e dVdl %12.5e\n", s, v, dvdlambda);
142 void do_force_lowlevel(FILE *fplog, gmx_large_int_t step,
143 t_forcerec *fr, t_inputrec *ir,
144 t_idef *idef, t_commrec *cr,
145 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
147 rvec x[], history_t *hist,
150 gmx_enerdata_t *enerd,
167 gmx_bool bDoEpot, bSepDVDL, bSB;
171 real Vsr, Vlr, Vcorr = 0;
175 double clam_i, vlam_i;
176 real dvdl_dum[efptNR], dvdl, dvdl_nb[efptNR], lam_i[efptNR];
180 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
183 #define PRINT_SEPDVDL(s, v, dvdlambda) if (bSepDVDL) { gmx_print_sepdvdl(fplog, s, v, dvdlambda); }
185 set_pbc(&pbc, fr->ePBC, box);
187 /* reset free energy components */
188 for (i = 0; i < efptNR; i++)
195 for (i = 0; (i < DIM); i++)
197 box_size[i] = box[i][i];
200 bSepDVDL = (fr->bSepDVDL && do_per_step(step, ir->nstlog));
203 /* do QMMM first if requested */
206 enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr);
211 fprintf(fplog, "Step %s: non-bonded V and dVdl for node %d:\n",
212 gmx_step_str(step, buf), cr->nodeid);
215 /* Call the short range functions all in one go. */
218 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
219 #define TAKETIME FALSE
222 MPI_Barrier(cr->mpi_comm_mygroup);
229 /* foreign lambda component for walls */
230 dvdl = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
231 enerd->grpp.ener[egLJSR], nrnb);
232 PRINT_SEPDVDL("Walls", 0.0, dvdl);
233 enerd->dvdl_lin[efptVDW] += dvdl;
236 /* If doing GB, reset dvda and calculate the Born radii */
237 if (ir->implicit_solvent)
239 wallcycle_sub_start(wcycle, ewcsNONBONDED);
241 for (i = 0; i < born->nr; i++)
248 calc_gb_rad(cr, fr, ir, top, x, &(fr->gblist), born, md, nrnb);
251 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
255 /* We only do non-bonded calculation with group scheme here, the verlet
256 * calls are done from do_force_cutsVERLET(). */
257 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
260 /* Add short-range interactions */
261 donb_flags |= GMX_NONBONDED_DO_SR;
263 if (flags & GMX_FORCE_FORCES)
265 donb_flags |= GMX_NONBONDED_DO_FORCE;
267 if (flags & GMX_FORCE_ENERGY)
269 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
271 if (flags & GMX_FORCE_DO_LR)
273 donb_flags |= GMX_NONBONDED_DO_LR;
276 wallcycle_sub_start(wcycle, ewcsNONBONDED);
277 do_nonbonded(fr, x, f, f_longrange, md, excl,
279 lambda, dvdl_nb, -1, -1, donb_flags);
281 /* If we do foreign lambda and we have soft-core interactions
282 * we have to recalculate the (non-linear) energies contributions.
284 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
286 for (i = 0; i < enerd->n_lambda; i++)
288 for (j = 0; j < efptNR; j++)
290 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
292 reset_foreign_enerdata(enerd);
293 do_nonbonded(fr, x, f, f_longrange, md, excl,
294 &(enerd->foreign_grpp), nrnb,
295 lam_i, dvdl_dum, -1, -1,
296 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
297 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
298 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
301 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
305 /* If we are doing GB, calculate bonded forces and apply corrections
306 * to the solvation forces */
307 /* MRS: Eventually, many need to include free energy contribution here! */
308 if (ir->implicit_solvent)
310 wallcycle_sub_start(wcycle, ewcsBONDED);
311 calc_gb_forces(cr, md, born, top, x, f, fr, idef,
312 ir->gb_algorithm, ir->sa_algorithm, nrnb, &pbc, graph, enerd);
313 wallcycle_sub_stop(wcycle, ewcsBONDED);
324 if (fepvals->sc_alpha != 0)
326 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
330 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
333 if (fepvals->sc_alpha != 0)
335 /* even though coulomb part is linear, we already added it, beacuse we
336 need to go through the vdw calculation anyway */
338 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
342 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
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 dvdlsum = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
357 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.", Vsr, dvdlsum);
364 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
367 /* Shift the coordinates. Must be done before bonded forces and PPPM,
368 * but is also necessary for SHAKE and update, therefore it can NOT
369 * go when no bonded forces have to be evaluated.
372 /* Here sometimes we would not need to shift with NBFonly,
373 * but we do so anyhow for consistency of the returned coordinates.
377 shift_self(graph, box, x);
380 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
384 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
387 /* Check whether we need to do bondeds or correct for exclusions */
389 ((flags & GMX_FORCE_BONDED)
390 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
392 /* Since all atoms are in the rectangular or triclinic unit-cell,
393 * only single box vector shifts (2 in x) are required.
395 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
399 if (flags & GMX_FORCE_BONDED)
401 wallcycle_sub_start(wcycle, ewcsBONDED);
402 calc_bonds(fplog, cr->ms,
403 idef, x, hist, f, fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
404 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
406 fr->bSepDVDL && do_per_step(step, ir->nstlog), step);
408 /* Check if we have to determine energy differences
409 * at foreign lambda's.
411 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
412 idef->ilsort != ilsortNO_FE)
414 if (idef->ilsort != ilsortFE_SORTED)
416 gmx_incons("The bonded interactions are not sorted for free energy");
418 for (i = 0; i < enerd->n_lambda; i++)
420 reset_foreign_enerdata(enerd);
421 for (j = 0; j < efptNR; j++)
423 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
425 calc_bonds_lambda(fplog, idef, x, fr, &pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
426 fcd, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
427 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
428 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
433 wallcycle_sub_stop(wcycle, ewcsBONDED);
439 if (EEL_FULL(fr->eeltype))
441 bSB = (ir->nwall == 2);
445 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
446 box_size[ZZ] *= ir->wall_ewald_zfac;
449 clear_mat(fr->vir_el_recip);
456 /* With the Verlet scheme exclusion forces are calculated
457 * in the non-bonded kernel.
459 /* The TPI molecule does not have exclusions with the rest
460 * of the system and no intra-molecular PME grid contributions
461 * will be calculated in gmx_pme_calc_energy.
463 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
464 ir->ewald_geometry != eewg3D ||
465 ir->epsilon_surface != 0)
469 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
473 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
476 nthreads = gmx_omp_nthreads_get(emntBonded);
477 #pragma omp parallel for num_threads(nthreads) schedule(static)
478 for (t = 0; t < nthreads; t++)
483 real *Vcorrt, *dvdlt;
486 fnv = fr->f_novirsum;
487 vir = &fr->vir_el_recip;
494 vir = &fr->f_t[t].vir;
495 Vcorrt = &fr->f_t[t].Vcorr;
496 dvdlt = &fr->f_t[t].dvdl[efptCOUL];
497 for (i = 0; i < fr->natoms_force; i++)
505 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
508 md->nChargePerturbed ? md->chargeB : NULL,
509 ir->cutoff_scheme != ecutsVERLET,
510 excl, x, bSB ? boxs : box, mu_tot,
514 lambda[efptCOUL], dvdlt);
518 reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
520 &Vcorr, efptCOUL, &dvdl,
524 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
529 Vcorr += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
530 &dvdl, fr->vir_el_recip);
533 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.", Vcorr, dvdl);
534 enerd->dvdl_lin[efptCOUL] += dvdl;
545 case eelPMEUSERSWITCH:
547 if (cr->duty & DUTY_PME)
549 assert(fr->n_tpi >= 0);
550 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
552 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
553 if (flags & GMX_FORCE_FORCES)
555 pme_flags |= GMX_PME_CALC_F;
557 if (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY))
559 pme_flags |= GMX_PME_CALC_ENER_VIR;
563 /* We don't calculate f, but we do want the potential */
564 pme_flags |= GMX_PME_CALC_POT;
566 wallcycle_start(wcycle, ewcPMEMESH);
567 status = gmx_pme_do(fr->pmedata,
568 md->start, md->homenr - fr->n_tpi,
570 md->chargeA, md->chargeB,
571 bSB ? boxs : box, cr,
572 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
573 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
575 fr->vir_el_recip, fr->ewaldcoeff,
576 &Vlr, lambda[efptCOUL], &dvdl,
578 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
580 /* We should try to do as little computation after
581 * this as possible, because parallel PME synchronizes
582 * the nodes, so we want all load imbalance of the rest
583 * of the force calculation to be before the PME call.
584 * DD load balancing is done on the whole time of
585 * the force call (without PME).
590 /* Determine the PME grid energy of the test molecule
591 * with the PME grid potential of the other charges.
593 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
594 x + md->homenr - fr->n_tpi,
595 md->chargeA + md->homenr - fr->n_tpi,
598 PRINT_SEPDVDL("PME mesh", Vlr, dvdl);
602 Vlr = do_ewald(ir, x, fr->f_novirsum,
603 md->chargeA, md->chargeB,
604 box_size, cr, md->homenr,
605 fr->vir_el_recip, fr->ewaldcoeff,
606 lambda[efptCOUL], &dvdl, fr->ewald_table);
607 PRINT_SEPDVDL("Ewald long-range", Vlr, dvdl);
610 gmx_fatal(FARGS, "No such electrostatics method implemented %s",
611 eel_names[fr->eeltype]);
615 gmx_fatal(FARGS, "Error %d in long range electrostatics routine %s",
616 status, EELTYPE(fr->eeltype));
618 /* Note that with separate PME nodes we get the real energies later */
619 enerd->dvdl_lin[efptCOUL] += dvdl;
620 enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
623 fprintf(debug, "Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
624 Vlr, Vcorr, enerd->term[F_COUL_RECIP]);
625 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
626 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
631 if (EEL_RF(fr->eeltype))
633 /* With the Verlet scheme exclusion forces are calculated
634 * in the non-bonded kernel.
636 if (ir->cutoff_scheme != ecutsVERLET && fr->eeltype != eelRF_NEC)
639 enerd->term[F_RF_EXCL] =
640 RF_excl_correction(fr, graph, md, excl, x, f,
641 fr->fshift, &pbc, lambda[efptCOUL], &dvdl);
644 enerd->dvdl_lin[efptCOUL] += dvdl;
645 PRINT_SEPDVDL("RF exclusion correction",
646 enerd->term[F_RF_EXCL], dvdl);
654 print_nrnb(debug, nrnb);
662 MPI_Barrier(cr->mpi_comm_mygroup);
665 if (fr->timesteps == 11)
667 fprintf(stderr, "* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
668 cr->nodeid, gmx_step_str(fr->timesteps, buf),
669 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
670 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
678 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
683 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
687 for (i = 0; i < F_NRE; i++)
690 enerd->foreign_term[i] = 0;
694 for (i = 0; i < efptNR; i++)
696 enerd->dvdl_lin[i] = 0;
697 enerd->dvdl_nonlin[i] = 0;
703 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
705 enerd->grpp.nener = n2;
706 enerd->foreign_grpp.nener = n2;
707 for (i = 0; (i < egNR); i++)
709 snew(enerd->grpp.ener[i], n2);
710 snew(enerd->foreign_grpp.ener[i], n2);
715 enerd->n_lambda = 1 + n_lambda;
716 snew(enerd->enerpart_lambda, enerd->n_lambda);
724 void destroy_enerdata(gmx_enerdata_t *enerd)
728 for (i = 0; (i < egNR); i++)
730 sfree(enerd->grpp.ener[i]);
733 for (i = 0; (i < egNR); i++)
735 sfree(enerd->foreign_grpp.ener[i]);
740 sfree(enerd->enerpart_lambda);
744 static real sum_v(int n, real v[])
750 for (i = 0; (i < n); i++)
758 void sum_epot(gmx_grppairener_t *grpp, real *epot)
762 /* Accumulate energies */
763 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
764 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
765 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
766 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
767 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
768 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
769 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
770 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
772 /* lattice part of LR doesnt belong to any group
773 * and has been added earlier
775 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
776 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
779 for (i = 0; (i < F_EPOT); i++)
781 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
783 epot[F_EPOT] += epot[i];
788 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
793 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
794 enerd->term[F_DVDL] = 0.0;
795 for (i = 0; i < efptNR; i++)
797 if (fepvals->separate_dvdl[i])
799 /* could this be done more readably/compactly? */
812 index = F_DVDL_BONDED;
814 case (efptRESTRAINT):
815 index = F_DVDL_RESTRAINT;
821 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
824 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
825 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
830 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
833 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
834 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
839 /* Notes on the foreign lambda free energy difference evaluation:
840 * Adding the potential and ekin terms that depend linearly on lambda
841 * as delta lam * dvdl to the energy differences is exact.
842 * For the constraints this is not exact, but we have no other option
843 * without literally changing the lengths and reevaluating the energies at each step.
844 * (try to remedy this post 4.6 - MRS)
845 * For the non-bonded LR term we assume that the soft-core (if present)
846 * no longer affects the energy beyond the short-range cut-off,
847 * which is a very good approximation (except for exotic settings).
848 * (investigate how to overcome this post 4.6 - MRS)
850 if (fepvals->separate_dvdl[efptBONDED])
852 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
856 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
858 enerd->term[F_DVDL_CONSTR] = 0;
860 for (i = 0; i < fepvals->n_lambda; i++)
862 /* note we are iterating over fepvals here!
863 For the current lam, dlam = 0 automatically,
864 so we don't need to add anything to the
865 enerd->enerpart_lambda[0] */
867 /* we don't need to worry about dvdl_lin contributions to dE at
868 current lambda, because the contributions to the current
869 lambda are automatically zeroed */
871 for (j = 0; j < efptNR; j++)
873 /* Note that this loop is over all dhdl components, not just the separated ones */
874 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
875 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
878 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
879 fepvals->all_lambda[j][i], efpt_names[j],
880 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
881 dlam, enerd->dvdl_lin[j]);
888 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
892 /* First reset all foreign energy components. Foreign energies always called on
893 neighbor search steps */
894 for (i = 0; (i < egNR); i++)
896 for (j = 0; (j < enerd->grpp.nener); j++)
898 enerd->foreign_grpp.ener[i][j] = 0.0;
902 /* potential energy components */
903 for (i = 0; (i <= F_EPOT); i++)
905 enerd->foreign_term[i] = 0.0;
909 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
910 gmx_enerdata_t *enerd,
916 /* First reset all energy components, except for the long range terms
917 * on the master at non neighbor search steps, since the long range
918 * terms have already been summed at the last neighbor search step.
920 bKeepLR = (fr->bTwinRange && !bNS);
921 for (i = 0; (i < egNR); i++)
923 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
925 for (j = 0; (j < enerd->grpp.nener); j++)
927 enerd->grpp.ener[i][j] = 0.0;
931 for (i = 0; i < efptNR; i++)
933 enerd->dvdl_lin[i] = 0.0;
934 enerd->dvdl_nonlin[i] = 0.0;
937 /* Normal potential energy components */
938 for (i = 0; (i <= F_EPOT); i++)
940 enerd->term[i] = 0.0;
942 /* Initialize the dVdlambda term with the long range contribution */
943 /* Initialize the dvdl term with the long range contribution */
944 enerd->term[F_DVDL] = 0.0;
945 enerd->term[F_DVDL_COUL] = 0.0;
946 enerd->term[F_DVDL_VDW] = 0.0;
947 enerd->term[F_DVDL_BONDED] = 0.0;
948 enerd->term[F_DVDL_RESTRAINT] = 0.0;
949 enerd->term[F_DKDL] = 0.0;
950 if (enerd->n_lambda > 0)
952 for (i = 0; i < enerd->n_lambda; i++)
954 enerd->enerpart_lambda[i] = 0.0;
957 /* reset foreign energy data - separate function since we also call it elsewhere */
958 reset_foreign_enerdata(enerd);