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 double clam_i,vlam_i;
174 real dvdl_dum[efptNR], dvdl, 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) fprintf(fplog,sepdvdlformat,s,v,dvdlambda);
183 GMX_MPE_LOG(ev_force_start);
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,md);
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. */
215 GMX_MPE_LOG(ev_do_fnbf_start);
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,atype,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(cr,fr,x,f,f_longrange,md,excl,
278 &enerd->grpp,box_size,nrnb,
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(cr,fr,x,f,f_longrange,md,excl,
294 &(enerd->foreign_grpp),box_size,nrnb,
295 lam_i,dvdl_dum,-1,-1,
296 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
297 sum_epot(&ir->opts,&(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,atype,x,f,fr,idef,
312 ir->gb_algorithm,ir->sa_algorithm,nrnb,bBornRadii,&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);
361 GMX_MPE_LOG(ev_do_fnbf_finish);
365 pr_rvecs(debug,0,"fshift after SR",fr->fshift,SHIFTS);
368 /* Shift the coordinates. Must be done before bonded forces and PPPM,
369 * but is also necessary for SHAKE and update, therefore it can NOT
370 * go when no bonded forces have to be evaluated.
373 /* Here sometimes we would not need to shift with NBFonly,
374 * but we do so anyhow for consistency of the returned coordinates.
378 shift_self(graph,box,x);
381 inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
385 inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
388 /* Check whether we need to do bondeds or correct for exclusions */
390 ((flags & GMX_FORCE_BONDED)
391 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
393 /* Since all atoms are in the rectangular or triclinic unit-cell,
394 * only single box vector shifts (2 in x) are required.
396 set_pbc_dd(&pbc,fr->ePBC,cr->dd,TRUE,box);
400 if (flags & GMX_FORCE_BONDED)
402 GMX_MPE_LOG(ev_calc_bonds_start);
404 wallcycle_sub_start(wcycle, ewcsBONDED);
405 calc_bonds(fplog,cr->ms,
406 idef,x,hist,f,fr,&pbc,graph,enerd,nrnb,lambda,md,fcd,
407 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
409 fr->bSepDVDL && do_per_step(step,ir->nstlog),step);
411 /* Check if we have to determine energy differences
412 * at foreign lambda's.
414 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
415 idef->ilsort != ilsortNO_FE)
417 if (idef->ilsort != ilsortFE_SORTED)
419 gmx_incons("The bonded interactions are not sorted for free energy");
421 for(i=0; i<enerd->n_lambda; i++)
423 reset_foreign_enerdata(enerd);
424 for (j=0;j<efptNR;j++)
426 lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
428 calc_bonds_lambda(fplog,idef,x,fr,&pbc,graph,&(enerd->foreign_grpp),enerd->foreign_term,nrnb,lam_i,md,
429 fcd,DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
430 sum_epot(&ir->opts,&(enerd->foreign_grpp),enerd->foreign_term);
431 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
435 GMX_MPE_LOG(ev_calc_bonds_finish);
436 wallcycle_sub_stop(wcycle, ewcsBONDED);
442 if (EEL_FULL(fr->eeltype))
444 bSB = (ir->nwall == 2);
448 svmul(ir->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
449 box_size[ZZ] *= ir->wall_ewald_zfac;
452 clear_mat(fr->vir_el_recip);
459 /* With the Verlet scheme exclusion forces are calculated
460 * in the non-bonded kernel.
462 /* The TPI molecule does not have exclusions with the rest
463 * of the system and no intra-molecular PME grid contributions
464 * will be calculated in gmx_pme_calc_energy.
466 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
467 ir->ewald_geometry != eewg3D ||
468 ir->epsilon_surface != 0)
472 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
476 gmx_fatal(FARGS,"TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
479 nthreads = gmx_omp_nthreads_get(emntBonded);
480 #pragma omp parallel for num_threads(nthreads) schedule(static)
481 for(t=0; t<nthreads; t++)
489 fnv = fr->f_novirsum;
490 vir = &fr->vir_el_recip;
497 vir = &fr->f_t[t].vir;
498 Vcorrt = &fr->f_t[t].Vcorr;
499 dvdlt = &fr->f_t[t].dvdl[efptCOUL];
500 for(i=0; i<fr->natoms_force; i++)
508 ewald_LRcorrection(fplog,
509 fr->excl_load[t],fr->excl_load[t+1],
512 md->nChargePerturbed ? md->chargeB : NULL,
513 ir->cutoff_scheme != ecutsVERLET,
514 excl,x,bSB ? boxs : box,mu_tot,
518 lambda[efptCOUL],dvdlt);
522 reduce_thread_forces(fr->natoms_force,fr->f_novirsum,
524 &Vcorr,efptCOUL,&dvdl,
528 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
533 Vcorr += ewald_charge_correction(cr,fr,lambda[efptCOUL],box,
534 &dvdl,fr->vir_el_recip);
537 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.",Vcorr,dvdl);
538 enerd->dvdl_lin[efptCOUL] += dvdl;
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,
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,
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);
606 Vlr = do_ewald(fplog,FALSE,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,fr->ewald_table);
611 PRINT_SEPDVDL("Ewald long-range",Vlr,dvdl);
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;
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 if (EEL_RF(fr->eeltype))
637 /* With the Verlet scheme exclusion forces are calculated
638 * in the non-bonded kernel.
640 if (ir->cutoff_scheme != ecutsVERLET && fr->eeltype != eelRF_NEC)
643 enerd->term[F_RF_EXCL] =
644 RF_excl_correction(fplog,fr,graph,md,excl,x,f,
645 fr->fshift,&pbc,lambda[efptCOUL],&dvdl);
648 enerd->dvdl_lin[efptCOUL] += dvdl;
649 PRINT_SEPDVDL("RF exclusion correction",
650 enerd->term[F_RF_EXCL],dvdl);
658 print_nrnb(debug,nrnb);
666 MPI_Barrier(cr->mpi_comm_mygroup);
669 if (fr->timesteps == 11)
671 fprintf(stderr,"* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
672 cr->nodeid, gmx_step_str(fr->timesteps,buf),
673 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
674 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
682 pr_rvecs(debug,0,"fshift after bondeds",fr->fshift,SHIFTS);
685 GMX_MPE_LOG(ev_force_finish);
689 void init_enerdata(int ngener,int n_lambda,gmx_enerdata_t *enerd)
693 for(i=0; i<F_NRE; i++)
696 enerd->foreign_term[i] = 0;
700 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[])
761 void sum_epot(t_grpopts *opts, gmx_grppairener_t *grpp, real *epot)
765 /* Accumulate energies */
766 epot[F_COUL_SR] = sum_v(grpp->nener,grpp->ener[egCOULSR]);
767 epot[F_LJ] = sum_v(grpp->nener,grpp->ener[egLJSR]);
768 epot[F_LJ14] = sum_v(grpp->nener,grpp->ener[egLJ14]);
769 epot[F_COUL14] = sum_v(grpp->nener,grpp->ener[egCOUL14]);
770 epot[F_COUL_LR] = sum_v(grpp->nener,grpp->ener[egCOULLR]);
771 epot[F_LJ_LR] = sum_v(grpp->nener,grpp->ener[egLJLR]);
772 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
773 epot[F_GBPOL] += sum_v(grpp->nener,grpp->ener[egGB]);
775 /* lattice part of LR doesnt belong to any group
776 * and has been added earlier
778 epot[F_BHAM] = sum_v(grpp->nener,grpp->ener[egBHAMSR]);
779 epot[F_BHAM_LR] = sum_v(grpp->nener,grpp->ener[egBHAMLR]);
782 for(i=0; (i<F_EPOT); i++)
784 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
786 epot[F_EPOT] += epot[i];
791 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
796 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
797 enerd->term[F_DVDL] = 0.0;
798 for (i=0;i<efptNR;i++)
800 if (fepvals->separate_dvdl[i])
802 /* could this be done more readably/compactly? */
811 index = F_DVDL_BONDED;
813 case (efptRESTRAINT):
814 index = F_DVDL_RESTRAINT;
823 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
826 fprintf(debug,"dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
827 efpt_names[i],i,enerd->term[index],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
832 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
835 fprintf(debug,"dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
836 efpt_names[0],i,enerd->term[F_DVDL],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
841 /* Notes on the foreign lambda free energy difference evaluation:
842 * Adding the potential and ekin terms that depend linearly on lambda
843 * as delta lam * dvdl to the energy differences is exact.
844 * For the constraints this is not exact, but we have no other option
845 * without literally changing the lengths and reevaluating the energies at each step.
846 * (try to remedy this post 4.6 - MRS)
847 * For the non-bonded LR term we assume that the soft-core (if present)
848 * no longer affects the energy beyond the short-range cut-off,
849 * which is a very good approximation (except for exotic settings).
850 * (investigate how to overcome this post 4.6 - MRS)
853 for(i=0; i<fepvals->n_lambda; i++)
854 { /* note we are iterating over fepvals here!
855 For the current lam, dlam = 0 automatically,
856 so we don't need to add anything to the
857 enerd->enerpart_lambda[0] */
859 /* we don't need to worry about dvdl contributions to the current lambda, because
860 it's automatically zero */
862 /* first kinetic energy term */
863 dlam = (fepvals->all_lambda[efptMASS][i] - lambda[efptMASS]);
865 enerd->enerpart_lambda[i+1] += enerd->term[F_DKDL]*dlam;
867 for (j=0;j<efptNR;j++)
869 if (j==efptMASS) {continue;} /* no other mass term to worry about */
871 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
872 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
875 fprintf(debug,"enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
876 fepvals->all_lambda[j][i],efpt_names[j],
877 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
878 dlam,enerd->dvdl_lin[j]);
885 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
889 /* First reset all foreign energy components. Foreign energies always called on
890 neighbor search steps */
891 for(i=0; (i<egNR); i++)
893 for(j=0; (j<enerd->grpp.nener); j++)
895 enerd->foreign_grpp.ener[i][j] = 0.0;
899 /* potential energy components */
900 for(i=0; (i<=F_EPOT); i++)
902 enerd->foreign_term[i] = 0.0;
906 void reset_enerdata(t_grpopts *opts,
907 t_forcerec *fr,gmx_bool bNS,
908 gmx_enerdata_t *enerd,
914 /* First reset all energy components, except for the long range terms
915 * on the master at non neighbor search steps, since the long range
916 * terms have already been summed at the last neighbor search step.
918 bKeepLR = (fr->bTwinRange && !bNS);
919 for(i=0; (i<egNR); i++) {
920 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR))) {
921 for(j=0; (j<enerd->grpp.nener); j++)
922 enerd->grpp.ener[i][j] = 0.0;
925 for (i=0;i<efptNR;i++)
927 enerd->dvdl_lin[i] = 0.0;
928 enerd->dvdl_nonlin[i] = 0.0;
931 /* Normal potential energy components */
932 for(i=0; (i<=F_EPOT); i++) {
933 enerd->term[i] = 0.0;
935 /* Initialize the dVdlambda term with the long range contribution */
936 /* Initialize the dvdl term with the long range contribution */
937 enerd->term[F_DVDL] = 0.0;
938 enerd->term[F_DVDL_COUL] = 0.0;
939 enerd->term[F_DVDL_VDW] = 0.0;
940 enerd->term[F_DVDL_BONDED] = 0.0;
941 enerd->term[F_DVDL_RESTRAINT] = 0.0;
942 enerd->term[F_DKDL] = 0.0;
943 if (enerd->n_lambda > 0)
945 for(i=0; i<enerd->n_lambda; i++)
947 enerd->enerpart_lambda[i] = 0.0;
950 /* reset foreign energy data - separate function since we also call it elsewhere */
951 reset_foreign_enerdata(enerd);