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 if (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(cr,fr,x,f,f_longrange,md,excl,
277 &enerd->grpp,box_size,nrnb,
278 lambda,dvdl_nb,-1,-1,donb_flags);
279 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
282 /* If we do foreign lambda and we have soft-core interactions
283 * we have to recalculate the (non-linear) energies contributions.
285 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
287 wallcycle_sub_start(wcycle, ewcsNONBONDED);
288 init_enerdata(mtop->groups.grps[egcENER].nr,fepvals->n_lambda,&ed_lam);
290 for(i=0; i<enerd->n_lambda; i++)
292 for (j=0;j<efptNR;j++)
294 lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
296 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
297 do_nonbonded(cr,fr,x,f,f_longrange,md,excl,
298 &(ed_lam.grpp), box_size,nrnb,
299 lam_i,dvdl_dum,-1,-1,
300 GMX_NONBONDED_DO_FOREIGNLAMBDA | GMX_NONBONDED_DO_SR);
301 sum_epot(&ir->opts,&ed_lam);
302 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
304 destroy_enerdata(&ed_lam);
305 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
309 /* If we are doing GB, calculate bonded forces and apply corrections
310 * to the solvation forces */
311 /* MRS: Eventually, many need to include free energy contribution here! */
312 if (ir->implicit_solvent)
314 calc_gb_forces(cr,md,born,top,atype,x,f,fr,idef,
315 ir->gb_algorithm,ir->sa_algorithm,nrnb,bBornRadii,&pbc,graph,enerd);
316 wallcycle_sub_stop(wcycle, ewcsBONDED);
327 if (fepvals->sc_alpha!=0)
329 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
333 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
336 if (fepvals->sc_alpha!=0)
338 /* even though coulomb part is linear, we already added it, beacuse we
339 need to go through the vdw calculation anyway */
341 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
345 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
351 for(i=0; i<enerd->grpp.nener; i++)
355 enerd->grpp.ener[egBHAMSR][i] :
356 enerd->grpp.ener[egLJSR][i])
357 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
359 dvdlsum = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
360 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",Vsr,dvdlsum);
364 GMX_MPE_LOG(ev_do_fnbf_finish);
368 pr_rvecs(debug,0,"fshift after SR",fr->fshift,SHIFTS);
371 /* Shift the coordinates. Must be done before bonded forces and PPPM,
372 * but is also necessary for SHAKE and update, therefore it can NOT
373 * go when no bonded forces have to be evaluated.
376 /* Here sometimes we would not need to shift with NBFonly,
377 * but we do so anyhow for consistency of the returned coordinates.
381 shift_self(graph,box,x);
384 inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
388 inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
391 /* Check whether we need to do bondeds or correct for exclusions */
393 ((flags & GMX_FORCE_BONDED)
394 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
396 /* Since all atoms are in the rectangular or triclinic unit-cell,
397 * only single box vector shifts (2 in x) are required.
399 set_pbc_dd(&pbc,fr->ePBC,cr->dd,TRUE,box);
403 if (flags & GMX_FORCE_BONDED)
405 GMX_MPE_LOG(ev_calc_bonds_start);
407 wallcycle_sub_start(wcycle, ewcsBONDED);
408 calc_bonds(fplog,cr->ms,
409 idef,x,hist,f,fr,&pbc,graph,enerd,nrnb,lambda,md,fcd,
410 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
412 fr->bSepDVDL && do_per_step(step,ir->nstlog),step);
414 /* Check if we have to determine energy differences
415 * at foreign lambda's.
417 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
418 idef->ilsort != ilsortNO_FE)
420 if (idef->ilsort != ilsortFE_SORTED)
422 gmx_incons("The bonded interactions are not sorted for free energy");
424 init_enerdata(mtop->groups.grps[egcENER].nr,fepvals->n_lambda,&ed_lam);
426 for(i=0; i<enerd->n_lambda; i++)
428 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
429 for (j=0;j<efptNR;j++)
431 lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
433 calc_bonds_lambda(fplog,idef,x,fr,&pbc,graph,&ed_lam,nrnb,lam_i,md,
434 fcd,DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
435 sum_epot(&ir->opts,&ed_lam);
436 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
438 destroy_enerdata(&ed_lam);
441 GMX_MPE_LOG(ev_calc_bonds_finish);
442 wallcycle_sub_stop(wcycle, ewcsBONDED);
448 if (EEL_FULL(fr->eeltype))
450 bSB = (ir->nwall == 2);
454 svmul(ir->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
455 box_size[ZZ] *= ir->wall_ewald_zfac;
458 clear_mat(fr->vir_el_recip);
465 /* With the Verlet scheme exclusion forces are calculated
466 * in the non-bonded kernel.
468 /* The TPI molecule does not have exclusions with the rest
469 * of the system and no intra-molecular PME grid contributions
470 * will be calculated in gmx_pme_calc_energy.
472 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
473 ir->ewald_geometry != eewg3D ||
474 ir->epsilon_surface != 0)
478 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
482 gmx_fatal(FARGS,"TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
485 nthreads = gmx_omp_nthreads_get(emntBonded);
486 #pragma omp parallel for num_threads(nthreads) schedule(static)
487 for(t=0; t<nthreads; t++)
495 fnv = fr->f_novirsum;
496 vir = &fr->vir_el_recip;
503 vir = &fr->f_t[t].vir;
504 Vcorrt = &fr->f_t[t].Vcorr;
505 dvdlt = &fr->f_t[t].dvdl[efptCOUL];
506 for(i=0; i<fr->natoms_force; i++)
514 ewald_LRcorrection(fplog,
515 fr->excl_load[t],fr->excl_load[t+1],
518 md->nChargePerturbed ? md->chargeB : NULL,
519 ir->cutoff_scheme != ecutsVERLET,
520 excl,x,bSB ? boxs : box,mu_tot,
524 lambda[efptCOUL],dvdlt);
528 reduce_thread_forces(fr->natoms_force,fr->f_novirsum,
530 &Vcorr,efptCOUL,&dvdl,
534 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
539 Vcorr += ewald_charge_correction(cr,fr,lambda[efptCOUL],box,
540 &dvdl,fr->vir_el_recip);
543 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.",Vcorr,dvdl);
544 enerd->dvdl_lin[efptCOUL] += dvdl;
555 case eelPMEUSERSWITCH:
557 if (cr->duty & DUTY_PME)
559 assert(fr->n_tpi >= 0);
560 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
562 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
563 if (flags & GMX_FORCE_FORCES)
565 pme_flags |= GMX_PME_CALC_F;
567 if (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY))
569 pme_flags |= GMX_PME_CALC_ENER_VIR;
573 /* We don't calculate f, but we do want the potential */
574 pme_flags |= GMX_PME_CALC_POT;
576 wallcycle_start(wcycle,ewcPMEMESH);
577 status = gmx_pme_do(fr->pmedata,
578 md->start,md->homenr - fr->n_tpi,
580 md->chargeA,md->chargeB,
582 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
583 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
585 fr->vir_el_recip,fr->ewaldcoeff,
586 &Vlr,lambda[efptCOUL],&dvdl,
588 *cycles_pme = wallcycle_stop(wcycle,ewcPMEMESH);
590 /* We should try to do as little computation after
591 * this as possible, because parallel PME synchronizes
592 * the nodes, so we want all load imbalance of the rest
593 * of the force calculation to be before the PME call.
594 * DD load balancing is done on the whole time of
595 * the force call (without PME).
600 /* Determine the PME grid energy of the test molecule
601 * with the PME grid potential of the other charges.
603 gmx_pme_calc_energy(fr->pmedata,fr->n_tpi,
604 x + md->homenr - fr->n_tpi,
605 md->chargeA + md->homenr - fr->n_tpi,
608 PRINT_SEPDVDL("PME mesh",Vlr,dvdl);
612 Vlr = do_ewald(fplog,FALSE,ir,x,fr->f_novirsum,
613 md->chargeA,md->chargeB,
614 box_size,cr,md->homenr,
615 fr->vir_el_recip,fr->ewaldcoeff,
616 lambda[efptCOUL],&dvdl,fr->ewald_table);
617 PRINT_SEPDVDL("Ewald long-range",Vlr,dvdl);
620 gmx_fatal(FARGS,"No such electrostatics method implemented %s",
621 eel_names[fr->eeltype]);
625 gmx_fatal(FARGS,"Error %d in long range electrostatics routine %s",
626 status,EELTYPE(fr->eeltype));
628 /* Note that with separate PME nodes we get the real energies later */
629 enerd->dvdl_lin[efptCOUL] += dvdl;
630 enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
633 fprintf(debug,"Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
634 Vlr,Vcorr,enerd->term[F_COUL_RECIP]);
635 pr_rvecs(debug,0,"vir_el_recip after corr",fr->vir_el_recip,DIM);
636 pr_rvecs(debug,0,"fshift after LR Corrections",fr->fshift,SHIFTS);
641 if (EEL_RF(fr->eeltype))
643 /* With the Verlet scheme exclusion forces are calculated
644 * in the non-bonded kernel.
646 if (ir->cutoff_scheme != ecutsVERLET && fr->eeltype != eelRF_NEC)
649 enerd->term[F_RF_EXCL] =
650 RF_excl_correction(fplog,fr,graph,md,excl,x,f,
651 fr->fshift,&pbc,lambda[efptCOUL],&dvdl);
654 enerd->dvdl_lin[efptCOUL] += dvdl;
655 PRINT_SEPDVDL("RF exclusion correction",
656 enerd->term[F_RF_EXCL],dvdl);
664 print_nrnb(debug,nrnb);
672 MPI_Barrier(cr->mpi_comm_mygroup);
675 if (fr->timesteps == 11)
677 fprintf(stderr,"* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
678 cr->nodeid, gmx_step_str(fr->timesteps,buf),
679 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
680 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
688 pr_rvecs(debug,0,"fshift after bondeds",fr->fshift,SHIFTS);
691 GMX_MPE_LOG(ev_force_finish);
695 void init_enerdata(int ngener,int n_lambda,gmx_enerdata_t *enerd)
699 for(i=0; i<F_NRE; i++)
705 for(i=0; i<efptNR; i++) {
706 enerd->dvdl_lin[i] = 0;
707 enerd->dvdl_nonlin[i] = 0;
713 fprintf(debug,"Creating %d sized group matrix for energies\n",n2);
715 enerd->grpp.nener = n2;
716 for(i=0; (i<egNR); i++)
718 snew(enerd->grpp.ener[i],n2);
723 enerd->n_lambda = 1 + n_lambda;
724 snew(enerd->enerpart_lambda,enerd->n_lambda);
732 void destroy_enerdata(gmx_enerdata_t *enerd)
736 for(i=0; (i<egNR); i++)
738 sfree(enerd->grpp.ener[i]);
743 sfree(enerd->enerpart_lambda);
747 static real sum_v(int n,real v[])
759 void sum_epot(t_grpopts *opts,gmx_enerdata_t *enerd)
761 gmx_grppairener_t *grpp;
768 /* Accumulate energies */
769 epot[F_COUL_SR] = sum_v(grpp->nener,grpp->ener[egCOULSR]);
770 epot[F_LJ] = sum_v(grpp->nener,grpp->ener[egLJSR]);
771 epot[F_LJ14] = sum_v(grpp->nener,grpp->ener[egLJ14]);
772 epot[F_COUL14] = sum_v(grpp->nener,grpp->ener[egCOUL14]);
773 epot[F_COUL_LR] = sum_v(grpp->nener,grpp->ener[egCOULLR]);
774 epot[F_LJ_LR] = sum_v(grpp->nener,grpp->ener[egLJLR]);
775 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
776 epot[F_GBPOL] += sum_v(grpp->nener,grpp->ener[egGB]);
778 /* lattice part of LR doesnt belong to any group
779 * and has been added earlier
781 epot[F_BHAM] = sum_v(grpp->nener,grpp->ener[egBHAMSR]);
782 epot[F_BHAM_LR] = sum_v(grpp->nener,grpp->ener[egBHAMLR]);
785 for(i=0; (i<F_EPOT); i++)
787 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
789 epot[F_EPOT] += epot[i];
794 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
799 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
800 enerd->term[F_DVDL] = 0.0;
801 for (i=0;i<efptNR;i++)
803 if (fepvals->separate_dvdl[i])
805 /* could this be done more readably/compactly? */
814 index = F_DVDL_BONDED;
816 case (efptRESTRAINT):
817 index = F_DVDL_RESTRAINT;
826 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
829 fprintf(debug,"dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
830 efpt_names[i],i,enerd->term[index],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
835 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
838 fprintf(debug,"dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
839 efpt_names[0],i,enerd->term[F_DVDL],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
844 /* Notes on the foreign lambda free energy difference evaluation:
845 * Adding the potential and ekin terms that depend linearly on lambda
846 * as delta lam * dvdl to the energy differences is exact.
847 * For the constraints this is not exact, but we have no other option
848 * without literally changing the lengths and reevaluating the energies at each step.
849 * (try to remedy this post 4.6 - MRS)
850 * For the non-bonded LR term we assume that the soft-core (if present)
851 * no longer affects the energy beyond the short-range cut-off,
852 * which is a very good approximation (except for exotic settings).
853 * (investigate how to overcome this post 4.6 - MRS)
856 for(i=0; i<fepvals->n_lambda; i++)
857 { /* note we are iterating over fepvals here!
858 For the current lam, dlam = 0 automatically,
859 so we don't need to add anything to the
860 enerd->enerpart_lambda[0] */
862 /* we don't need to worry about dvdl contributions to the current lambda, because
863 it's automatically zero */
865 /* first kinetic energy term */
866 dlam = (fepvals->all_lambda[efptMASS][i] - lambda[efptMASS]);
868 enerd->enerpart_lambda[i+1] += enerd->term[F_DKDL]*dlam;
870 for (j=0;j<efptNR;j++)
872 if (j==efptMASS) {continue;} /* no other mass term to worry about */
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]);
887 void reset_enerdata(t_grpopts *opts,
888 t_forcerec *fr,gmx_bool bNS,
889 gmx_enerdata_t *enerd,
895 /* First reset all energy components, except for the long range terms
896 * on the master at non neighbor search steps, since the long range
897 * terms have already been summed at the last neighbor search step.
899 bKeepLR = (fr->bTwinRange && !bNS);
900 for(i=0; (i<egNR); i++) {
901 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR))) {
902 for(j=0; (j<enerd->grpp.nener); j++)
903 enerd->grpp.ener[i][j] = 0.0;
906 for (i=0;i<efptNR;i++)
908 enerd->dvdl_lin[i] = 0.0;
909 enerd->dvdl_nonlin[i] = 0.0;
912 /* Normal potential energy components */
913 for(i=0; (i<=F_EPOT); i++) {
914 enerd->term[i] = 0.0;
916 /* Initialize the dVdlambda term with the long range contribution */
917 /* Initialize the dvdl term with the long range contribution */
918 enerd->term[F_DVDL] = 0.0;
919 enerd->term[F_DVDL_COUL] = 0.0;
920 enerd->term[F_DVDL_VDW] = 0.0;
921 enerd->term[F_DVDL_BONDED] = 0.0;
922 enerd->term[F_DVDL_RESTRAINT] = 0.0;
923 enerd->term[F_DKDL] = 0.0;
924 if (enerd->n_lambda > 0)
926 for(i=0; i<enerd->n_lambda; i++)
928 enerd->enerpart_lambda[i] = 0.0;