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"
77 gmx_bool bDoLongRangeNS)
83 if (!fr->ns.nblist_initialized)
85 init_neighbor_list(fp, fr, md->homenr);
93 nsearch = search_neighbours(fp, fr, box, top, groups, cr, nrnb, md,
94 bFillGrid, bDoLongRangeNS);
97 fprintf(debug, "nsearch = %d\n", nsearch);
100 /* Check whether we have to do dynamic load balancing */
101 /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
102 count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
103 &(top->idef),opts->ngener);
105 if (fr->ns.dump_nl > 0)
107 dump_nblist(fp, cr, fr, fr->ns.dump_nl);
111 static void reduce_thread_forces(int n, rvec *f,
114 int efpt_ind, real *dvdl,
115 int nthreads, f_thread_t *f_t)
119 /* This reduction can run over any number of threads */
120 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntBonded)) private(t) schedule(static)
121 for (i = 0; i < n; i++)
123 for (t = 1; t < nthreads; t++)
125 rvec_inc(f[i], f_t[t].f[i]);
128 for (t = 1; t < nthreads; t++)
130 *Vcorr += f_t[t].Vcorr;
131 *dvdl += f_t[t].dvdl[efpt_ind];
132 m_add(vir, f_t[t].vir, vir);
136 void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda)
138 fprintf(fplog, " %-30s V %12.5e dVdl %12.5e\n", s, v, dvdlambda);
141 void do_force_lowlevel(FILE *fplog, gmx_large_int_t step,
142 t_forcerec *fr, t_inputrec *ir,
143 t_idef *idef, t_commrec *cr,
144 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
146 rvec x[], history_t *hist,
149 gmx_enerdata_t *enerd,
166 gmx_bool bDoEpot, bSepDVDL, bSB;
170 real Vsr, Vlr, Vcorr = 0;
174 double clam_i, vlam_i;
175 real dvdl_dum[efptNR], dvdl, dvdl_nb[efptNR], lam_i[efptNR];
179 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
182 #define PRINT_SEPDVDL(s, v, dvdlambda) if (bSepDVDL) { gmx_print_sepdvdl(fplog, s, v, dvdlambda); }
184 set_pbc(&pbc, fr->ePBC, box);
186 /* reset free energy components */
187 for (i = 0; i < efptNR; i++)
194 for (i = 0; (i < DIM); i++)
196 box_size[i] = box[i][i];
199 bSepDVDL = (fr->bSepDVDL && do_per_step(step, ir->nstlog));
202 /* do QMMM first if requested */
205 enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr);
210 fprintf(fplog, "Step %s: non-bonded V and dVdl for node %d:\n",
211 gmx_step_str(step, buf), cr->nodeid);
214 /* Call the short range functions all in one go. */
217 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
218 #define TAKETIME FALSE
221 MPI_Barrier(cr->mpi_comm_mygroup);
228 /* foreign lambda component for walls */
229 dvdl = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
230 enerd->grpp.ener[egLJSR], nrnb);
231 PRINT_SEPDVDL("Walls", 0.0, dvdl);
232 enerd->dvdl_lin[efptVDW] += dvdl;
235 /* If doing GB, reset dvda and calculate the Born radii */
236 if (ir->implicit_solvent)
238 wallcycle_sub_start(wcycle, ewcsNONBONDED);
240 for (i = 0; i < born->nr; i++)
247 calc_gb_rad(cr, fr, ir, top, x, &(fr->gblist), born, md, nrnb);
250 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
254 /* We only do non-bonded calculation with group scheme here, the verlet
255 * calls are done from do_force_cutsVERLET(). */
256 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
259 /* Add short-range interactions */
260 donb_flags |= GMX_NONBONDED_DO_SR;
262 if (flags & GMX_FORCE_FORCES)
264 donb_flags |= GMX_NONBONDED_DO_FORCE;
266 if (flags & GMX_FORCE_ENERGY)
268 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
270 if (flags & GMX_FORCE_DO_LR)
272 donb_flags |= GMX_NONBONDED_DO_LR;
275 wallcycle_sub_start(wcycle, ewcsNONBONDED);
276 do_nonbonded(fr, x, f, f_longrange, md, excl,
278 lambda, dvdl_nb, -1, -1, donb_flags);
280 /* If we do foreign lambda and we have soft-core interactions
281 * we have to recalculate the (non-linear) energies contributions.
283 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
285 for (i = 0; i < enerd->n_lambda; i++)
287 for (j = 0; j < efptNR; j++)
289 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
291 reset_foreign_enerdata(enerd);
292 do_nonbonded(fr, x, f, f_longrange, md, excl,
293 &(enerd->foreign_grpp), nrnb,
294 lam_i, dvdl_dum, -1, -1,
295 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
296 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
297 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
300 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
304 /* If we are doing GB, calculate bonded forces and apply corrections
305 * to the solvation forces */
306 /* MRS: Eventually, many need to include free energy contribution here! */
307 if (ir->implicit_solvent)
309 wallcycle_sub_start(wcycle, ewcsBONDED);
310 calc_gb_forces(cr, md, born, top, x, f, fr, idef,
311 ir->gb_algorithm, ir->sa_algorithm, nrnb, &pbc, graph, enerd);
312 wallcycle_sub_stop(wcycle, ewcsBONDED);
323 if (fepvals->sc_alpha != 0)
325 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
329 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
332 if (fepvals->sc_alpha != 0)
334 /* even though coulomb part is linear, we already added it, beacuse we
335 need to go through the vdw calculation anyway */
337 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
341 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
347 for (i = 0; i < enerd->grpp.nener; i++)
351 enerd->grpp.ener[egBHAMSR][i] :
352 enerd->grpp.ener[egLJSR][i])
353 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
355 dvdlsum = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
356 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.", Vsr, dvdlsum);
363 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
366 /* Shift the coordinates. Must be done before bonded forces and PPPM,
367 * but is also necessary for SHAKE and update, therefore it can NOT
368 * go when no bonded forces have to be evaluated.
371 /* Here sometimes we would not need to shift with NBFonly,
372 * but we do so anyhow for consistency of the returned coordinates.
376 shift_self(graph, box, x);
379 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
383 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
386 /* Check whether we need to do bondeds or correct for exclusions */
388 ((flags & GMX_FORCE_BONDED)
389 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
391 /* Since all atoms are in the rectangular or triclinic unit-cell,
392 * only single box vector shifts (2 in x) are required.
394 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
398 if (flags & GMX_FORCE_BONDED)
400 wallcycle_sub_start(wcycle, ewcsBONDED);
401 calc_bonds(fplog, cr->ms,
402 idef, x, hist, f, fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
403 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
405 fr->bSepDVDL && do_per_step(step, ir->nstlog), step);
407 /* Check if we have to determine energy differences
408 * at foreign lambda's.
410 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
411 idef->ilsort != ilsortNO_FE)
413 if (idef->ilsort != ilsortFE_SORTED)
415 gmx_incons("The bonded interactions are not sorted for free energy");
417 for (i = 0; i < enerd->n_lambda; i++)
419 reset_foreign_enerdata(enerd);
420 for (j = 0; j < efptNR; j++)
422 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
424 calc_bonds_lambda(fplog, idef, x, fr, &pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
425 fcd, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
426 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
427 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
432 wallcycle_sub_stop(wcycle, ewcsBONDED);
438 if (EEL_FULL(fr->eeltype))
440 bSB = (ir->nwall == 2);
444 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
445 box_size[ZZ] *= ir->wall_ewald_zfac;
448 clear_mat(fr->vir_el_recip);
455 /* With the Verlet scheme exclusion forces are calculated
456 * in the non-bonded kernel.
458 /* The TPI molecule does not have exclusions with the rest
459 * of the system and no intra-molecular PME grid contributions
460 * will be calculated in gmx_pme_calc_energy.
462 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
463 ir->ewald_geometry != eewg3D ||
464 ir->epsilon_surface != 0)
468 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
472 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
475 nthreads = gmx_omp_nthreads_get(emntBonded);
476 #pragma omp parallel for num_threads(nthreads) schedule(static)
477 for (t = 0; t < nthreads; t++)
482 real *Vcorrt, *dvdlt;
485 fnv = fr->f_novirsum;
486 vir = &fr->vir_el_recip;
493 vir = &fr->f_t[t].vir;
494 Vcorrt = &fr->f_t[t].Vcorr;
495 dvdlt = &fr->f_t[t].dvdl[efptCOUL];
496 for (i = 0; i < fr->natoms_force; i++)
504 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
507 md->nChargePerturbed ? md->chargeB : NULL,
508 ir->cutoff_scheme != ecutsVERLET,
509 excl, x, bSB ? boxs : box, mu_tot,
513 lambda[efptCOUL], dvdlt);
517 reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
519 &Vcorr, efptCOUL, &dvdl,
523 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
528 Vcorr += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
529 &dvdl, fr->vir_el_recip);
532 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.", Vcorr, dvdl);
533 enerd->dvdl_lin[efptCOUL] += dvdl;
544 case eelPMEUSERSWITCH:
546 if (cr->duty & DUTY_PME)
548 assert(fr->n_tpi >= 0);
549 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
551 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
552 if (flags & GMX_FORCE_FORCES)
554 pme_flags |= GMX_PME_CALC_F;
556 if (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY))
558 pme_flags |= GMX_PME_CALC_ENER_VIR;
562 /* We don't calculate f, but we do want the potential */
563 pme_flags |= GMX_PME_CALC_POT;
565 wallcycle_start(wcycle, ewcPMEMESH);
566 status = gmx_pme_do(fr->pmedata,
567 md->start, md->homenr - fr->n_tpi,
569 md->chargeA, md->chargeB,
570 bSB ? boxs : box, cr,
571 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
572 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
574 fr->vir_el_recip, fr->ewaldcoeff,
575 &Vlr, lambda[efptCOUL], &dvdl,
577 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
579 /* We should try to do as little computation after
580 * this as possible, because parallel PME synchronizes
581 * the nodes, so we want all load imbalance of the rest
582 * of the force calculation to be before the PME call.
583 * DD load balancing is done on the whole time of
584 * the force call (without PME).
589 /* Determine the PME grid energy of the test molecule
590 * with the PME grid potential of the other charges.
592 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
593 x + md->homenr - fr->n_tpi,
594 md->chargeA + md->homenr - fr->n_tpi,
597 PRINT_SEPDVDL("PME mesh", Vlr, dvdl);
601 Vlr = do_ewald(ir, x, fr->f_novirsum,
602 md->chargeA, md->chargeB,
603 box_size, cr, md->homenr,
604 fr->vir_el_recip, fr->ewaldcoeff,
605 lambda[efptCOUL], &dvdl, fr->ewald_table);
606 PRINT_SEPDVDL("Ewald long-range", Vlr, dvdl);
609 gmx_fatal(FARGS, "No such electrostatics method implemented %s",
610 eel_names[fr->eeltype]);
614 gmx_fatal(FARGS, "Error %d in long range electrostatics routine %s",
615 status, EELTYPE(fr->eeltype));
617 /* Note that with separate PME nodes we get the real energies later */
618 enerd->dvdl_lin[efptCOUL] += dvdl;
619 enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
622 fprintf(debug, "Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
623 Vlr, Vcorr, enerd->term[F_COUL_RECIP]);
624 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
625 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
630 if (EEL_RF(fr->eeltype))
632 /* With the Verlet scheme exclusion forces are calculated
633 * in the non-bonded kernel.
635 if (ir->cutoff_scheme != ecutsVERLET && fr->eeltype != eelRF_NEC)
638 enerd->term[F_RF_EXCL] =
639 RF_excl_correction(fr, graph, md, excl, x, f,
640 fr->fshift, &pbc, lambda[efptCOUL], &dvdl);
643 enerd->dvdl_lin[efptCOUL] += dvdl;
644 PRINT_SEPDVDL("RF exclusion correction",
645 enerd->term[F_RF_EXCL], dvdl);
653 print_nrnb(debug, nrnb);
661 MPI_Barrier(cr->mpi_comm_mygroup);
664 if (fr->timesteps == 11)
666 fprintf(stderr, "* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
667 cr->nodeid, gmx_step_str(fr->timesteps, buf),
668 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
669 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
677 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
682 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
686 for (i = 0; i < F_NRE; i++)
689 enerd->foreign_term[i] = 0;
693 for (i = 0; i < efptNR; i++)
695 enerd->dvdl_lin[i] = 0;
696 enerd->dvdl_nonlin[i] = 0;
702 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
704 enerd->grpp.nener = n2;
705 enerd->foreign_grpp.nener = n2;
706 for (i = 0; (i < egNR); i++)
708 snew(enerd->grpp.ener[i], n2);
709 snew(enerd->foreign_grpp.ener[i], n2);
714 enerd->n_lambda = 1 + n_lambda;
715 snew(enerd->enerpart_lambda, enerd->n_lambda);
723 void destroy_enerdata(gmx_enerdata_t *enerd)
727 for (i = 0; (i < egNR); i++)
729 sfree(enerd->grpp.ener[i]);
732 for (i = 0; (i < egNR); i++)
734 sfree(enerd->foreign_grpp.ener[i]);
739 sfree(enerd->enerpart_lambda);
743 static real sum_v(int n, real v[])
749 for (i = 0; (i < n); i++)
757 void sum_epot(gmx_grppairener_t *grpp, real *epot)
761 /* Accumulate energies */
762 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
763 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
764 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
765 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
766 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
767 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
768 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
769 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
771 /* lattice part of LR doesnt belong to any group
772 * and has been added earlier
774 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
775 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
778 for (i = 0; (i < F_EPOT); i++)
780 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
782 epot[F_EPOT] += epot[i];
787 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
792 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
793 enerd->term[F_DVDL] = 0.0;
794 for (i = 0; i < efptNR; i++)
796 if (fepvals->separate_dvdl[i])
798 /* could this be done more readably/compactly? */
811 index = F_DVDL_BONDED;
813 case (efptRESTRAINT):
814 index = F_DVDL_RESTRAINT;
820 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
823 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
824 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
829 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
832 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
833 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
838 /* Notes on the foreign lambda free energy difference evaluation:
839 * Adding the potential and ekin terms that depend linearly on lambda
840 * as delta lam * dvdl to the energy differences is exact.
841 * For the constraints this is not exact, but we have no other option
842 * without literally changing the lengths and reevaluating the energies at each step.
843 * (try to remedy this post 4.6 - MRS)
844 * For the non-bonded LR term we assume that the soft-core (if present)
845 * no longer affects the energy beyond the short-range cut-off,
846 * which is a very good approximation (except for exotic settings).
847 * (investigate how to overcome this post 4.6 - MRS)
849 if (fepvals->separate_dvdl[efptBONDED])
851 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
855 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
857 enerd->term[F_DVDL_CONSTR] = 0;
859 for (i = 0; i < fepvals->n_lambda; i++)
861 /* note we are iterating over fepvals here!
862 For the current lam, dlam = 0 automatically,
863 so we don't need to add anything to the
864 enerd->enerpart_lambda[0] */
866 /* we don't need to worry about dvdl_lin contributions to dE at
867 current lambda, because the contributions to the current
868 lambda are automatically zeroed */
870 for (j = 0; j < efptNR; j++)
872 /* Note that this loop is over all dhdl components, not just the separated ones */
873 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
874 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
877 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
878 fepvals->all_lambda[j][i], efpt_names[j],
879 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
880 dlam, enerd->dvdl_lin[j]);
887 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
891 /* First reset all foreign energy components. Foreign energies always called on
892 neighbor search steps */
893 for (i = 0; (i < egNR); i++)
895 for (j = 0; (j < enerd->grpp.nener); j++)
897 enerd->foreign_grpp.ener[i][j] = 0.0;
901 /* potential energy components */
902 for (i = 0; (i <= F_EPOT); i++)
904 enerd->foreign_term[i] = 0.0;
908 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
909 gmx_enerdata_t *enerd,
915 /* First reset all energy components, except for the long range terms
916 * on the master at non neighbor search steps, since the long range
917 * terms have already been summed at the last neighbor search step.
919 bKeepLR = (fr->bTwinRange && !bNS);
920 for (i = 0; (i < egNR); i++)
922 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
924 for (j = 0; (j < enerd->grpp.nener); j++)
926 enerd->grpp.ener[i][j] = 0.0;
930 for (i = 0; i < efptNR; i++)
932 enerd->dvdl_lin[i] = 0.0;
933 enerd->dvdl_nonlin[i] = 0.0;
936 /* Normal potential energy components */
937 for (i = 0; (i <= F_EPOT); i++)
939 enerd->term[i] = 0.0;
941 /* Initialize the dVdlambda term with the long range contribution */
942 /* Initialize the dvdl term with the long range contribution */
943 enerd->term[F_DVDL] = 0.0;
944 enerd->term[F_DVDL_COUL] = 0.0;
945 enerd->term[F_DVDL_VDW] = 0.0;
946 enerd->term[F_DVDL_BONDED] = 0.0;
947 enerd->term[F_DVDL_RESTRAINT] = 0.0;
948 enerd->term[F_DKDL] = 0.0;
949 if (enerd->n_lambda > 0)
951 for (i = 0; i < enerd->n_lambda; i++)
953 enerd->enerpart_lambda[i] = 0.0;
956 /* reset foreign energy data - separate function since we also call it elsewhere */
957 reset_foreign_enerdata(enerd);