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"
80 gmx_grppairener_t *grppener,
82 gmx_bool bDoLongRangeNS)
88 if (!fr->ns.nblist_initialized)
90 init_neighbor_list(fp, fr, md->homenr);
98 nsearch = search_neighbours(fp, fr, x, box, top, groups, cr, nrnb, md,
99 lambda, dvdlambda, grppener,
100 bFillGrid, bDoLongRangeNS, TRUE);
103 fprintf(debug, "nsearch = %d\n", nsearch);
106 /* Check whether we have to do dynamic load balancing */
107 /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
108 count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
109 &(top->idef),opts->ngener);
111 if (fr->ns.dump_nl > 0)
113 dump_nblist(fp, cr, fr, fr->ns.dump_nl);
117 static void reduce_thread_forces(int n, rvec *f,
120 int efpt_ind, real *dvdl,
121 int nthreads, f_thread_t *f_t)
125 /* This reduction can run over any number of threads */
126 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntBonded)) private(t) schedule(static)
127 for (i = 0; i < n; i++)
129 for (t = 1; t < nthreads; t++)
131 rvec_inc(f[i], f_t[t].f[i]);
134 for (t = 1; t < nthreads; t++)
136 *Vcorr += f_t[t].Vcorr;
137 *dvdl += f_t[t].dvdl[efpt_ind];
138 m_add(vir, f_t[t].vir, vir);
142 void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda)
144 fprintf(fplog, " %-30s V %12.5e dVdl %12.5e\n", s, v, dvdlambda);
147 void do_force_lowlevel(FILE *fplog, gmx_large_int_t step,
148 t_forcerec *fr, t_inputrec *ir,
149 t_idef *idef, t_commrec *cr,
150 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
153 rvec x[], history_t *hist,
156 gmx_enerdata_t *enerd,
174 gmx_bool bDoEpot, bSepDVDL, bSB;
178 real Vsr, Vlr, Vcorr = 0;
182 double clam_i, vlam_i;
183 real dvdl_dum[efptNR], dvdl, dvdl_nb[efptNR], lam_i[efptNR];
187 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
190 #define PRINT_SEPDVDL(s, v, dvdlambda) if (bSepDVDL) { gmx_print_sepdvdl(fplog, s, v, dvdlambda); }
192 set_pbc(&pbc, fr->ePBC, box);
194 /* reset free energy components */
195 for (i = 0; i < efptNR; i++)
202 for (i = 0; (i < DIM); i++)
204 box_size[i] = box[i][i];
207 bSepDVDL = (fr->bSepDVDL && do_per_step(step, ir->nstlog));
210 /* do QMMM first if requested */
213 enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr, md);
218 fprintf(fplog, "Step %s: non-bonded V and dVdl for node %d:\n",
219 gmx_step_str(step, buf), cr->nodeid);
222 /* Call the short range functions all in one go. */
225 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
226 #define TAKETIME FALSE
229 MPI_Barrier(cr->mpi_comm_mygroup);
236 /* foreign lambda component for walls */
237 dvdl = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
238 enerd->grpp.ener[egLJSR], nrnb);
239 PRINT_SEPDVDL("Walls", 0.0, dvdl);
240 enerd->dvdl_lin[efptVDW] += dvdl;
243 /* If doing GB, reset dvda and calculate the Born radii */
244 if (ir->implicit_solvent)
246 wallcycle_sub_start(wcycle, ewcsNONBONDED);
248 for (i = 0; i < born->nr; i++)
255 calc_gb_rad(cr, fr, ir, top, atype, x, &(fr->gblist), born, md, nrnb);
258 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
262 /* We only do non-bonded calculation with group scheme here, the verlet
263 * calls are done from do_force_cutsVERLET(). */
264 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
267 /* Add short-range interactions */
268 donb_flags |= GMX_NONBONDED_DO_SR;
270 if (flags & GMX_FORCE_FORCES)
272 donb_flags |= GMX_NONBONDED_DO_FORCE;
274 if (flags & GMX_FORCE_ENERGY)
276 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
278 if (flags & GMX_FORCE_DO_LR)
280 donb_flags |= GMX_NONBONDED_DO_LR;
283 wallcycle_sub_start(wcycle, ewcsNONBONDED);
284 do_nonbonded(cr, fr, x, f, f_longrange, md, excl,
285 &enerd->grpp, box_size, nrnb,
286 lambda, dvdl_nb, -1, -1, donb_flags);
288 /* If we do foreign lambda and we have soft-core interactions
289 * we have to recalculate the (non-linear) energies contributions.
291 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
293 for (i = 0; i < enerd->n_lambda; i++)
295 for (j = 0; j < efptNR; j++)
297 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
299 reset_foreign_enerdata(enerd);
300 do_nonbonded(cr, fr, x, f, f_longrange, md, excl,
301 &(enerd->foreign_grpp), box_size, nrnb,
302 lam_i, dvdl_dum, -1, -1,
303 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
304 sum_epot(&ir->opts, &(enerd->foreign_grpp), enerd->foreign_term);
305 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
308 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
312 /* If we are doing GB, calculate bonded forces and apply corrections
313 * to the solvation forces */
314 /* MRS: Eventually, many need to include free energy contribution here! */
315 if (ir->implicit_solvent)
317 wallcycle_sub_start(wcycle, ewcsBONDED);
318 calc_gb_forces(cr, md, born, top, atype, x, f, fr, idef,
319 ir->gb_algorithm, ir->sa_algorithm, nrnb, bBornRadii, &pbc, graph, enerd);
320 wallcycle_sub_stop(wcycle, ewcsBONDED);
331 if (fepvals->sc_alpha != 0)
333 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
337 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
340 if (fepvals->sc_alpha != 0)
342 /* even though coulomb part is linear, we already added it, beacuse we
343 need to go through the vdw calculation anyway */
345 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
349 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
355 for (i = 0; i < enerd->grpp.nener; i++)
359 enerd->grpp.ener[egBHAMSR][i] :
360 enerd->grpp.ener[egLJSR][i])
361 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
363 dvdlsum = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
364 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.", Vsr, dvdlsum);
371 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
374 /* Shift the coordinates. Must be done before bonded forces and PPPM,
375 * but is also necessary for SHAKE and update, therefore it can NOT
376 * go when no bonded forces have to be evaluated.
379 /* Here sometimes we would not need to shift with NBFonly,
380 * but we do so anyhow for consistency of the returned coordinates.
384 shift_self(graph, box, x);
387 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
391 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
394 /* Check whether we need to do bondeds or correct for exclusions */
396 ((flags & GMX_FORCE_BONDED)
397 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
399 /* Since all atoms are in the rectangular or triclinic unit-cell,
400 * only single box vector shifts (2 in x) are required.
402 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
406 if (flags & GMX_FORCE_BONDED)
408 wallcycle_sub_start(wcycle, ewcsBONDED);
409 calc_bonds(fplog, cr->ms,
410 idef, x, hist, f, fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
411 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
413 fr->bSepDVDL && do_per_step(step, ir->nstlog), step);
415 /* Check if we have to determine energy differences
416 * at foreign lambda's.
418 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
419 idef->ilsort != ilsortNO_FE)
421 if (idef->ilsort != ilsortFE_SORTED)
423 gmx_incons("The bonded interactions are not sorted for free energy");
425 for (i = 0; i < enerd->n_lambda; i++)
427 reset_foreign_enerdata(enerd);
428 for (j = 0; j < efptNR; j++)
430 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
432 calc_bonds_lambda(fplog, idef, x, fr, &pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
433 fcd, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
434 sum_epot(&ir->opts, &(enerd->foreign_grpp), enerd->foreign_term);
435 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
440 wallcycle_sub_stop(wcycle, ewcsBONDED);
446 if (EEL_FULL(fr->eeltype))
448 bSB = (ir->nwall == 2);
452 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
453 box_size[ZZ] *= ir->wall_ewald_zfac;
456 clear_mat(fr->vir_el_recip);
463 /* With the Verlet scheme exclusion forces are calculated
464 * in the non-bonded kernel.
466 /* The TPI molecule does not have exclusions with the rest
467 * of the system and no intra-molecular PME grid contributions
468 * will be calculated in gmx_pme_calc_energy.
470 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
471 ir->ewald_geometry != eewg3D ||
472 ir->epsilon_surface != 0)
476 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
480 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
483 nthreads = gmx_omp_nthreads_get(emntBonded);
484 #pragma omp parallel for num_threads(nthreads) schedule(static)
485 for (t = 0; t < nthreads; t++)
490 real *Vcorrt, *dvdlt;
493 fnv = fr->f_novirsum;
494 vir = &fr->vir_el_recip;
501 vir = &fr->f_t[t].vir;
502 Vcorrt = &fr->f_t[t].Vcorr;
503 dvdlt = &fr->f_t[t].dvdl[efptCOUL];
504 for (i = 0; i < fr->natoms_force; i++)
512 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
515 md->nChargePerturbed ? md->chargeB : NULL,
516 ir->cutoff_scheme != ecutsVERLET,
517 excl, x, bSB ? boxs : box, mu_tot,
521 lambda[efptCOUL], dvdlt);
525 reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
527 &Vcorr, efptCOUL, &dvdl,
531 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
536 Vcorr += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
537 &dvdl, fr->vir_el_recip);
540 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.", Vcorr, dvdl);
541 enerd->dvdl_lin[efptCOUL] += dvdl;
552 case eelPMEUSERSWITCH:
554 if (cr->duty & DUTY_PME)
556 assert(fr->n_tpi >= 0);
557 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
559 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
560 if (flags & GMX_FORCE_FORCES)
562 pme_flags |= GMX_PME_CALC_F;
564 if (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY))
566 pme_flags |= GMX_PME_CALC_ENER_VIR;
570 /* We don't calculate f, but we do want the potential */
571 pme_flags |= GMX_PME_CALC_POT;
573 wallcycle_start(wcycle, ewcPMEMESH);
574 status = gmx_pme_do(fr->pmedata,
575 md->start, md->homenr - fr->n_tpi,
577 md->chargeA, md->chargeB,
578 bSB ? boxs : box, cr,
579 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
580 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
582 fr->vir_el_recip, fr->ewaldcoeff,
583 &Vlr, lambda[efptCOUL], &dvdl,
585 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
587 /* We should try to do as little computation after
588 * this as possible, because parallel PME synchronizes
589 * the nodes, so we want all load imbalance of the rest
590 * of the force calculation to be before the PME call.
591 * DD load balancing is done on the whole time of
592 * the force call (without PME).
597 /* Determine the PME grid energy of the test molecule
598 * with the PME grid potential of the other charges.
600 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
601 x + md->homenr - fr->n_tpi,
602 md->chargeA + md->homenr - fr->n_tpi,
605 PRINT_SEPDVDL("PME mesh", Vlr, dvdl);
609 Vlr = do_ewald(fplog, FALSE, ir, x, fr->f_novirsum,
610 md->chargeA, md->chargeB,
611 box_size, cr, md->homenr,
612 fr->vir_el_recip, fr->ewaldcoeff,
613 lambda[efptCOUL], &dvdl, fr->ewald_table);
614 PRINT_SEPDVDL("Ewald long-range", Vlr, dvdl);
617 gmx_fatal(FARGS, "No such electrostatics method implemented %s",
618 eel_names[fr->eeltype]);
622 gmx_fatal(FARGS, "Error %d in long range electrostatics routine %s",
623 status, EELTYPE(fr->eeltype));
625 /* Note that with separate PME nodes we get the real energies later */
626 enerd->dvdl_lin[efptCOUL] += dvdl;
627 enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
630 fprintf(debug, "Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
631 Vlr, Vcorr, enerd->term[F_COUL_RECIP]);
632 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
633 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
638 if (EEL_RF(fr->eeltype))
640 /* With the Verlet scheme exclusion forces are calculated
641 * in the non-bonded kernel.
643 if (ir->cutoff_scheme != ecutsVERLET && fr->eeltype != eelRF_NEC)
646 enerd->term[F_RF_EXCL] =
647 RF_excl_correction(fplog, fr, graph, md, excl, x, f,
648 fr->fshift, &pbc, lambda[efptCOUL], &dvdl);
651 enerd->dvdl_lin[efptCOUL] += dvdl;
652 PRINT_SEPDVDL("RF exclusion correction",
653 enerd->term[F_RF_EXCL], dvdl);
661 print_nrnb(debug, nrnb);
669 MPI_Barrier(cr->mpi_comm_mygroup);
672 if (fr->timesteps == 11)
674 fprintf(stderr, "* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
675 cr->nodeid, gmx_step_str(fr->timesteps, buf),
676 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
677 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
685 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
690 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
694 for (i = 0; i < F_NRE; i++)
697 enerd->foreign_term[i] = 0;
701 for (i = 0; i < efptNR; i++)
703 enerd->dvdl_lin[i] = 0;
704 enerd->dvdl_nonlin[i] = 0;
710 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
712 enerd->grpp.nener = n2;
713 enerd->foreign_grpp.nener = n2;
714 for (i = 0; (i < egNR); i++)
716 snew(enerd->grpp.ener[i], n2);
717 snew(enerd->foreign_grpp.ener[i], n2);
722 enerd->n_lambda = 1 + n_lambda;
723 snew(enerd->enerpart_lambda, enerd->n_lambda);
731 void destroy_enerdata(gmx_enerdata_t *enerd)
735 for (i = 0; (i < egNR); i++)
737 sfree(enerd->grpp.ener[i]);
740 for (i = 0; (i < egNR); i++)
742 sfree(enerd->foreign_grpp.ener[i]);
747 sfree(enerd->enerpart_lambda);
751 static real sum_v(int n, real v[])
757 for (i = 0; (i < n); i++)
765 void sum_epot(t_grpopts *opts, gmx_grppairener_t *grpp, real *epot)
769 /* Accumulate energies */
770 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
771 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
772 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
773 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
774 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
775 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
776 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
777 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
779 /* lattice part of LR doesnt belong to any group
780 * and has been added earlier
782 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
783 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
786 for (i = 0; (i < F_EPOT); i++)
788 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
790 epot[F_EPOT] += epot[i];
795 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
800 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
801 enerd->term[F_DVDL] = 0.0;
802 for (i = 0; i < efptNR; i++)
804 if (fepvals->separate_dvdl[i])
806 /* could this be done more readably/compactly? */
819 index = F_DVDL_BONDED;
821 case (efptRESTRAINT):
822 index = F_DVDL_RESTRAINT;
828 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
831 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
832 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
837 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
840 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
841 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
846 /* Notes on the foreign lambda free energy difference evaluation:
847 * Adding the potential and ekin terms that depend linearly on lambda
848 * as delta lam * dvdl to the energy differences is exact.
849 * For the constraints this is not exact, but we have no other option
850 * without literally changing the lengths and reevaluating the energies at each step.
851 * (try to remedy this post 4.6 - MRS)
852 * For the non-bonded LR term we assume that the soft-core (if present)
853 * no longer affects the energy beyond the short-range cut-off,
854 * which is a very good approximation (except for exotic settings).
855 * (investigate how to overcome this post 4.6 - MRS)
857 if (fepvals->separate_dvdl[efptBONDED])
859 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
863 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
865 enerd->term[F_DVDL_CONSTR] = 0;
867 for (i = 0; i < fepvals->n_lambda; i++)
868 { /* note we are iterating over fepvals here!
869 For the current lam, dlam = 0 automatically,
870 so we don't need to add anything to the
871 enerd->enerpart_lambda[0] */
873 /* we don't need to worry about dvdl_lin contributions to dE at
874 current lambda, because the contributions to the current
875 lambda are automatically zeroed */
877 for (j = 0; j < efptNR; j++)
879 /* Note that this loop is over all dhdl components, not just the separated ones */
880 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
881 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
884 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
885 fepvals->all_lambda[j][i], efpt_names[j],
886 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
887 dlam, enerd->dvdl_lin[j]);
894 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
898 /* First reset all foreign energy components. Foreign energies always called on
899 neighbor search steps */
900 for (i = 0; (i < egNR); i++)
902 for (j = 0; (j < enerd->grpp.nener); j++)
904 enerd->foreign_grpp.ener[i][j] = 0.0;
908 /* potential energy components */
909 for (i = 0; (i <= F_EPOT); i++)
911 enerd->foreign_term[i] = 0.0;
915 void reset_enerdata(t_grpopts *opts,
916 t_forcerec *fr, gmx_bool bNS,
917 gmx_enerdata_t *enerd,
923 /* First reset all energy components, except for the long range terms
924 * on the master at non neighbor search steps, since the long range
925 * terms have already been summed at the last neighbor search step.
927 bKeepLR = (fr->bTwinRange && !bNS);
928 for (i = 0; (i < egNR); i++)
930 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
932 for (j = 0; (j < enerd->grpp.nener); j++)
934 enerd->grpp.ener[i][j] = 0.0;
938 for (i = 0; i < efptNR; i++)
940 enerd->dvdl_lin[i] = 0.0;
941 enerd->dvdl_nonlin[i] = 0.0;
944 /* Normal potential energy components */
945 for (i = 0; (i <= F_EPOT); i++)
947 enerd->term[i] = 0.0;
949 /* Initialize the dVdlambda term with the long range contribution */
950 /* Initialize the dvdl term with the long range contribution */
951 enerd->term[F_DVDL] = 0.0;
952 enerd->term[F_DVDL_COUL] = 0.0;
953 enerd->term[F_DVDL_VDW] = 0.0;
954 enerd->term[F_DVDL_BONDED] = 0.0;
955 enerd->term[F_DVDL_RESTRAINT] = 0.0;
956 enerd->term[F_DKDL] = 0.0;
957 if (enerd->n_lambda > 0)
959 for (i = 0; i < enerd->n_lambda; i++)
961 enerd->enerpart_lambda[i] = 0.0;
964 /* reset foreign energy data - separate function since we also call it elsewhere */
965 reset_foreign_enerdata(enerd);