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 "mpelogging.h"
66 #include "gmx_omp_nthreads.h"
80 gmx_grppairener_t *grppener,
82 gmx_bool bDoLongRangeNS)
87 GMX_MPE_LOG(ev_ns_start);
88 if (!fr->ns.nblist_initialized)
90 init_neighbor_list(fp, fr, md->homenr);
96 nsearch = search_neighbours(fp,fr,x,box,top,groups,cr,nrnb,md,
97 lambda,dvdlambda,grppener,
98 bFillGrid,bDoLongRangeNS,TRUE);
100 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)
108 dump_nblist(fp,cr,fr,fr->ns.dump_nl);
110 GMX_MPE_LOG(ev_ns_finish);
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)
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 do_force_lowlevel(FILE *fplog, gmx_large_int_t step,
139 t_forcerec *fr, t_inputrec *ir,
140 t_idef *idef, t_commrec *cr,
141 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
144 rvec x[], history_t *hist,
147 gmx_enerdata_t *enerd,
165 gmx_bool bDoEpot,bSepDVDL,bSB;
169 real Vsr,Vlr,Vcorr=0;
173 gmx_enerdata_t ed_lam;
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) fprintf(fplog,sepdvdlformat,s,v,dvdlambda);
184 GMX_MPE_LOG(ev_force_start);
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,md);
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. */
216 GMX_MPE_LOG(ev_do_fnbf_start);
219 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
220 #define TAKETIME FALSE
223 MPI_Barrier(cr->mpi_comm_mygroup);
230 /* foreign lambda component for walls */
231 dvdl = do_walls(ir,fr,box,md,x,f,lambda[efptVDW],
232 enerd->grpp.ener[egLJSR],nrnb);
233 PRINT_SEPDVDL("Walls",0.0,dvdl);
234 enerd->dvdl_lin[efptVDW] += dvdl;
237 /* If doing GB, reset dvda and calculate the Born radii */
238 if (ir->implicit_solvent)
240 wallcycle_sub_start(wcycle, ewcsNONBONDED);
242 for(i=0;i<born->nr;i++)
249 calc_gb_rad(cr,fr,ir,top,atype,x,&(fr->gblist),born,md,nrnb);
252 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
256 /* We only do non-bonded calculation with group scheme here, the verlet
257 * calls are done from do_force_cutsVERLET(). */
258 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
261 /* Add short-range interactions */
262 donb_flags |= GMX_NONBONDED_DO_SR;
264 if (flags & GMX_FORCE_FORCES)
266 donb_flags |= GMX_NONBONDED_DO_FORCE;
268 if (flags & GMX_FORCE_ENERGY)
270 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
272 if (flags & GMX_FORCE_DO_LR)
274 donb_flags |= GMX_NONBONDED_DO_LR;
277 wallcycle_sub_start(wcycle, ewcsNONBONDED);
278 do_nonbonded(cr,fr,x,f,f_longrange,md,excl,
279 &enerd->grpp,box_size,nrnb,
280 lambda,dvdl_nb,-1,-1,donb_flags);
281 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
284 /* If we do foreign lambda and we have soft-core interactions
285 * we have to recalculate the (non-linear) energies contributions.
287 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
289 wallcycle_sub_start(wcycle, ewcsNONBONDED);
290 init_enerdata(mtop->groups.grps[egcENER].nr,fepvals->n_lambda,&ed_lam);
292 for(i=0; i<enerd->n_lambda; i++)
294 for (j=0;j<efptNR;j++)
296 lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
298 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
299 do_nonbonded(cr,fr,x,f,f_longrange,md,excl,
300 &(ed_lam.grpp), box_size,nrnb,
301 lam_i,dvdl_dum,-1,-1,
302 GMX_NONBONDED_DO_FOREIGNLAMBDA | GMX_NONBONDED_DO_SR);
303 sum_epot(&ir->opts,&ed_lam);
304 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
306 destroy_enerdata(&ed_lam);
307 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
311 /* If we are doing GB, calculate bonded forces and apply corrections
312 * to the solvation forces */
313 /* MRS: Eventually, many need to include free energy contribution here! */
314 if (ir->implicit_solvent)
316 wallcycle_sub_start(wcycle, ewcsBONDED);
317 calc_gb_forces(cr,md,born,top,atype,x,f,fr,idef,
318 ir->gb_algorithm,ir->sa_algorithm,nrnb,bBornRadii,&pbc,graph,enerd);
319 wallcycle_sub_stop(wcycle, ewcsBONDED);
330 if (fepvals->sc_alpha!=0)
332 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
336 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
339 if (fepvals->sc_alpha!=0)
341 /* even though coulomb part is linear, we already added it, beacuse we
342 need to go through the vdw calculation anyway */
344 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
348 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
354 for(i=0; i<enerd->grpp.nener; i++)
358 enerd->grpp.ener[egBHAMSR][i] :
359 enerd->grpp.ener[egLJSR][i])
360 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
362 dvdlsum = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
363 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",Vsr,dvdlsum);
367 GMX_MPE_LOG(ev_do_fnbf_finish);
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 GMX_MPE_LOG(ev_calc_bonds_start);
410 wallcycle_sub_start(wcycle, ewcsBONDED);
411 calc_bonds(fplog,cr->ms,
412 idef,x,hist,f,fr,&pbc,graph,enerd,nrnb,lambda,md,fcd,
413 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
415 fr->bSepDVDL && do_per_step(step,ir->nstlog),step);
417 /* Check if we have to determine energy differences
418 * at foreign lambda's.
420 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
421 idef->ilsort != ilsortNO_FE)
423 if (idef->ilsort != ilsortFE_SORTED)
425 gmx_incons("The bonded interactions are not sorted for free energy");
427 init_enerdata(mtop->groups.grps[egcENER].nr,fepvals->n_lambda,&ed_lam);
429 for(i=0; i<enerd->n_lambda; i++)
431 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
432 for (j=0;j<efptNR;j++)
434 lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
436 calc_bonds_lambda(fplog,idef,x,fr,&pbc,graph,&ed_lam,nrnb,lam_i,md,
437 fcd,DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
438 sum_epot(&ir->opts,&ed_lam);
439 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
441 destroy_enerdata(&ed_lam);
444 GMX_MPE_LOG(ev_calc_bonds_finish);
445 wallcycle_sub_stop(wcycle, ewcsBONDED);
451 if (EEL_FULL(fr->eeltype))
453 bSB = (ir->nwall == 2);
457 svmul(ir->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
458 box_size[ZZ] *= ir->wall_ewald_zfac;
461 clear_mat(fr->vir_el_recip);
468 /* With the Verlet scheme exclusion forces are calculated
469 * in the non-bonded kernel.
471 /* The TPI molecule does not have exclusions with the rest
472 * of the system and no intra-molecular PME grid contributions
473 * will be calculated in gmx_pme_calc_energy.
475 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
476 ir->ewald_geometry != eewg3D ||
477 ir->epsilon_surface != 0)
481 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
485 gmx_fatal(FARGS,"TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
488 nthreads = gmx_omp_nthreads_get(emntBonded);
489 #pragma omp parallel for num_threads(nthreads) schedule(static)
490 for(t=0; t<nthreads; t++)
498 fnv = fr->f_novirsum;
499 vir = &fr->vir_el_recip;
506 vir = &fr->f_t[t].vir;
507 Vcorrt = &fr->f_t[t].Vcorr;
508 dvdlt = &fr->f_t[t].dvdl[efptCOUL];
509 for(i=0; i<fr->natoms_force; i++)
517 ewald_LRcorrection(fplog,
518 fr->excl_load[t],fr->excl_load[t+1],
521 md->nChargePerturbed ? md->chargeB : NULL,
522 ir->cutoff_scheme != ecutsVERLET,
523 excl,x,bSB ? boxs : box,mu_tot,
527 lambda[efptCOUL],dvdlt);
531 reduce_thread_forces(fr->natoms_force,fr->f_novirsum,
533 &Vcorr,efptCOUL,&dvdl,
537 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
542 Vcorr += ewald_charge_correction(cr,fr,lambda[efptCOUL],box,
543 &dvdl,fr->vir_el_recip);
546 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.",Vcorr,dvdl);
547 enerd->dvdl_lin[efptCOUL] += dvdl;
558 case eelPMEUSERSWITCH:
560 if (cr->duty & DUTY_PME)
562 assert(fr->n_tpi >= 0);
563 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
565 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
566 if (flags & GMX_FORCE_FORCES)
568 pme_flags |= GMX_PME_CALC_F;
570 if (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY))
572 pme_flags |= GMX_PME_CALC_ENER_VIR;
576 /* We don't calculate f, but we do want the potential */
577 pme_flags |= GMX_PME_CALC_POT;
579 wallcycle_start(wcycle,ewcPMEMESH);
580 status = gmx_pme_do(fr->pmedata,
581 md->start,md->homenr - fr->n_tpi,
583 md->chargeA,md->chargeB,
585 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
586 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
588 fr->vir_el_recip,fr->ewaldcoeff,
589 &Vlr,lambda[efptCOUL],&dvdl,
591 *cycles_pme = wallcycle_stop(wcycle,ewcPMEMESH);
593 /* We should try to do as little computation after
594 * this as possible, because parallel PME synchronizes
595 * the nodes, so we want all load imbalance of the rest
596 * of the force calculation to be before the PME call.
597 * DD load balancing is done on the whole time of
598 * the force call (without PME).
603 /* Determine the PME grid energy of the test molecule
604 * with the PME grid potential of the other charges.
606 gmx_pme_calc_energy(fr->pmedata,fr->n_tpi,
607 x + md->homenr - fr->n_tpi,
608 md->chargeA + md->homenr - fr->n_tpi,
611 PRINT_SEPDVDL("PME mesh",Vlr,dvdl);
615 Vlr = do_ewald(fplog,FALSE,ir,x,fr->f_novirsum,
616 md->chargeA,md->chargeB,
617 box_size,cr,md->homenr,
618 fr->vir_el_recip,fr->ewaldcoeff,
619 lambda[efptCOUL],&dvdl,fr->ewald_table);
620 PRINT_SEPDVDL("Ewald long-range",Vlr,dvdl);
623 gmx_fatal(FARGS,"No such electrostatics method implemented %s",
624 eel_names[fr->eeltype]);
628 gmx_fatal(FARGS,"Error %d in long range electrostatics routine %s",
629 status,EELTYPE(fr->eeltype));
631 /* Note that with separate PME nodes we get the real energies later */
632 enerd->dvdl_lin[efptCOUL] += dvdl;
633 enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
636 fprintf(debug,"Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
637 Vlr,Vcorr,enerd->term[F_COUL_RECIP]);
638 pr_rvecs(debug,0,"vir_el_recip after corr",fr->vir_el_recip,DIM);
639 pr_rvecs(debug,0,"fshift after LR Corrections",fr->fshift,SHIFTS);
644 if (EEL_RF(fr->eeltype))
646 /* With the Verlet scheme exclusion forces are calculated
647 * in the non-bonded kernel.
649 if (ir->cutoff_scheme != ecutsVERLET && fr->eeltype != eelRF_NEC)
652 enerd->term[F_RF_EXCL] =
653 RF_excl_correction(fplog,fr,graph,md,excl,x,f,
654 fr->fshift,&pbc,lambda[efptCOUL],&dvdl);
657 enerd->dvdl_lin[efptCOUL] += dvdl;
658 PRINT_SEPDVDL("RF exclusion correction",
659 enerd->term[F_RF_EXCL],dvdl);
667 print_nrnb(debug,nrnb);
675 MPI_Barrier(cr->mpi_comm_mygroup);
678 if (fr->timesteps == 11)
680 fprintf(stderr,"* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
681 cr->nodeid, gmx_step_str(fr->timesteps,buf),
682 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
683 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
691 pr_rvecs(debug,0,"fshift after bondeds",fr->fshift,SHIFTS);
694 GMX_MPE_LOG(ev_force_finish);
698 void init_enerdata(int ngener,int n_lambda,gmx_enerdata_t *enerd)
702 for(i=0; i<F_NRE; i++)
708 for(i=0; i<efptNR; i++) {
709 enerd->dvdl_lin[i] = 0;
710 enerd->dvdl_nonlin[i] = 0;
716 fprintf(debug,"Creating %d sized group matrix for energies\n",n2);
718 enerd->grpp.nener = n2;
719 for(i=0; (i<egNR); i++)
721 snew(enerd->grpp.ener[i],n2);
726 enerd->n_lambda = 1 + n_lambda;
727 snew(enerd->enerpart_lambda,enerd->n_lambda);
735 void destroy_enerdata(gmx_enerdata_t *enerd)
739 for(i=0; (i<egNR); i++)
741 sfree(enerd->grpp.ener[i]);
746 sfree(enerd->enerpart_lambda);
750 static real sum_v(int n,real v[])
762 void sum_epot(t_grpopts *opts,gmx_enerdata_t *enerd)
764 gmx_grppairener_t *grpp;
771 /* Accumulate energies */
772 epot[F_COUL_SR] = sum_v(grpp->nener,grpp->ener[egCOULSR]);
773 epot[F_LJ] = sum_v(grpp->nener,grpp->ener[egLJSR]);
774 epot[F_LJ14] = sum_v(grpp->nener,grpp->ener[egLJ14]);
775 epot[F_COUL14] = sum_v(grpp->nener,grpp->ener[egCOUL14]);
776 epot[F_COUL_LR] = sum_v(grpp->nener,grpp->ener[egCOULLR]);
777 epot[F_LJ_LR] = sum_v(grpp->nener,grpp->ener[egLJLR]);
778 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
779 epot[F_GBPOL] += sum_v(grpp->nener,grpp->ener[egGB]);
781 /* lattice part of LR doesnt belong to any group
782 * and has been added earlier
784 epot[F_BHAM] = sum_v(grpp->nener,grpp->ener[egBHAMSR]);
785 epot[F_BHAM_LR] = sum_v(grpp->nener,grpp->ener[egBHAMLR]);
788 for(i=0; (i<F_EPOT); i++)
790 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
792 epot[F_EPOT] += epot[i];
797 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
802 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
803 enerd->term[F_DVDL] = 0.0;
804 for (i=0;i<efptNR;i++)
806 if (fepvals->separate_dvdl[i])
808 /* could this be done more readably/compactly? */
817 index = F_DVDL_BONDED;
819 case (efptRESTRAINT):
820 index = F_DVDL_RESTRAINT;
829 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
832 fprintf(debug,"dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
833 efpt_names[i],i,enerd->term[index],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
838 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
841 fprintf(debug,"dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
842 efpt_names[0],i,enerd->term[F_DVDL],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
847 /* Notes on the foreign lambda free energy difference evaluation:
848 * Adding the potential and ekin terms that depend linearly on lambda
849 * as delta lam * dvdl to the energy differences is exact.
850 * For the constraints this is not exact, but we have no other option
851 * without literally changing the lengths and reevaluating the energies at each step.
852 * (try to remedy this post 4.6 - MRS)
853 * For the non-bonded LR term we assume that the soft-core (if present)
854 * no longer affects the energy beyond the short-range cut-off,
855 * which is a very good approximation (except for exotic settings).
856 * (investigate how to overcome this post 4.6 - MRS)
859 for(i=0; i<fepvals->n_lambda; i++)
860 { /* note we are iterating over fepvals here!
861 For the current lam, dlam = 0 automatically,
862 so we don't need to add anything to the
863 enerd->enerpart_lambda[0] */
865 /* we don't need to worry about dvdl contributions to the current lambda, because
866 it's automatically zero */
868 /* first kinetic energy term */
869 dlam = (fepvals->all_lambda[efptMASS][i] - lambda[efptMASS]);
871 enerd->enerpart_lambda[i+1] += enerd->term[F_DKDL]*dlam;
873 for (j=0;j<efptNR;j++)
875 if (j==efptMASS) {continue;} /* no other mass term to worry about */
877 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
878 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
881 fprintf(debug,"enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
882 fepvals->all_lambda[j][i],efpt_names[j],
883 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
884 dlam,enerd->dvdl_lin[j]);
890 void reset_enerdata(t_grpopts *opts,
891 t_forcerec *fr,gmx_bool bNS,
892 gmx_enerdata_t *enerd,
898 /* First reset all energy components, except for the long range terms
899 * on the master at non neighbor search steps, since the long range
900 * terms have already been summed at the last neighbor search step.
902 bKeepLR = (fr->bTwinRange && !bNS);
903 for(i=0; (i<egNR); i++) {
904 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR))) {
905 for(j=0; (j<enerd->grpp.nener); j++)
906 enerd->grpp.ener[i][j] = 0.0;
909 for (i=0;i<efptNR;i++)
911 enerd->dvdl_lin[i] = 0.0;
912 enerd->dvdl_nonlin[i] = 0.0;
915 /* Normal potential energy components */
916 for(i=0; (i<=F_EPOT); i++) {
917 enerd->term[i] = 0.0;
919 /* Initialize the dVdlambda term with the long range contribution */
920 /* Initialize the dvdl term with the long range contribution */
921 enerd->term[F_DVDL] = 0.0;
922 enerd->term[F_DVDL_COUL] = 0.0;
923 enerd->term[F_DVDL_VDW] = 0.0;
924 enerd->term[F_DVDL_BONDED] = 0.0;
925 enerd->term[F_DVDL_RESTRAINT] = 0.0;
926 enerd->term[F_DKDL] = 0.0;
927 if (enerd->n_lambda > 0)
929 for(i=0; i<enerd->n_lambda; i++)
931 enerd->enerpart_lambda[i] = 0.0;