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);
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);
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)
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 do_force_lowlevel(FILE *fplog, gmx_large_int_t step,
137 t_forcerec *fr, t_inputrec *ir,
138 t_idef *idef, t_commrec *cr,
139 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
142 rvec x[], history_t *hist,
145 gmx_enerdata_t *enerd,
163 gmx_bool bDoEpot,bSepDVDL,bSB;
167 real Vsr,Vlr,Vcorr=0;
171 double clam_i,vlam_i;
172 real dvdl_dum[efptNR], dvdl, dvdl_nb[efptNR], lam_i[efptNR];
176 double t0=0.0,t1,t2,t3; /* time measurement for coarse load balancing */
179 #define PRINT_SEPDVDL(s,v,dvdlambda) if (bSepDVDL) fprintf(fplog,sepdvdlformat,s,v,dvdlambda);
182 set_pbc(&pbc,fr->ePBC,box);
184 /* reset free energy components */
185 for (i=0;i<efptNR;i++)
192 for(i=0; (i<DIM); i++)
194 box_size[i]=box[i][i];
197 bSepDVDL=(fr->bSepDVDL && do_per_step(step,ir->nstlog));
200 /* do QMMM first if requested */
203 enerd->term[F_EQM] = calculate_QMMM(cr,x,f,fr,md);
208 fprintf(fplog,"Step %s: non-bonded V and dVdl for node %d:\n",
209 gmx_step_str(step,buf),cr->nodeid);
212 /* Call the short range functions all in one go. */
215 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
216 #define TAKETIME FALSE
219 MPI_Barrier(cr->mpi_comm_mygroup);
226 /* foreign lambda component for walls */
227 dvdl = do_walls(ir,fr,box,md,x,f,lambda[efptVDW],
228 enerd->grpp.ener[egLJSR],nrnb);
229 PRINT_SEPDVDL("Walls",0.0,dvdl);
230 enerd->dvdl_lin[efptVDW] += dvdl;
233 /* If doing GB, reset dvda and calculate the Born radii */
234 if (ir->implicit_solvent)
236 wallcycle_sub_start(wcycle, ewcsNONBONDED);
238 for(i=0;i<born->nr;i++)
245 calc_gb_rad(cr,fr,ir,top,atype,x,&(fr->gblist),born,md,nrnb);
248 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
252 /* We only do non-bonded calculation with group scheme here, the verlet
253 * calls are done from do_force_cutsVERLET(). */
254 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
257 /* Add short-range interactions */
258 donb_flags |= GMX_NONBONDED_DO_SR;
260 if (flags & GMX_FORCE_FORCES)
262 donb_flags |= GMX_NONBONDED_DO_FORCE;
264 if (flags & GMX_FORCE_ENERGY)
266 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
268 if (flags & GMX_FORCE_DO_LR)
270 donb_flags |= GMX_NONBONDED_DO_LR;
273 wallcycle_sub_start(wcycle, ewcsNONBONDED);
274 do_nonbonded(cr,fr,x,f,f_longrange,md,excl,
275 &enerd->grpp,box_size,nrnb,
276 lambda,dvdl_nb,-1,-1,donb_flags);
278 /* If we do foreign lambda and we have soft-core interactions
279 * we have to recalculate the (non-linear) energies contributions.
281 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
283 for(i=0; i<enerd->n_lambda; i++)
285 for (j=0;j<efptNR;j++)
287 lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
289 reset_foreign_enerdata(enerd);
290 do_nonbonded(cr,fr,x,f,f_longrange,md,excl,
291 &(enerd->foreign_grpp),box_size,nrnb,
292 lam_i,dvdl_dum,-1,-1,
293 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
294 sum_epot(&ir->opts,&(enerd->foreign_grpp),enerd->foreign_term);
295 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
298 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
302 /* If we are doing GB, calculate bonded forces and apply corrections
303 * to the solvation forces */
304 /* MRS: Eventually, many need to include free energy contribution here! */
305 if (ir->implicit_solvent)
307 wallcycle_sub_start(wcycle, ewcsBONDED);
308 calc_gb_forces(cr,md,born,top,atype,x,f,fr,idef,
309 ir->gb_algorithm,ir->sa_algorithm,nrnb,bBornRadii,&pbc,graph,enerd);
310 wallcycle_sub_stop(wcycle, ewcsBONDED);
321 if (fepvals->sc_alpha!=0)
323 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
327 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
330 if (fepvals->sc_alpha!=0)
332 /* even though coulomb part is linear, we already added it, beacuse we
333 need to go through the vdw calculation anyway */
335 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
339 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
345 for(i=0; i<enerd->grpp.nener; i++)
349 enerd->grpp.ener[egBHAMSR][i] :
350 enerd->grpp.ener[egLJSR][i])
351 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
353 dvdlsum = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
354 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",Vsr,dvdlsum);
361 pr_rvecs(debug,0,"fshift after SR",fr->fshift,SHIFTS);
364 /* Shift the coordinates. Must be done before bonded forces and PPPM,
365 * but is also necessary for SHAKE and update, therefore it can NOT
366 * go when no bonded forces have to be evaluated.
369 /* Here sometimes we would not need to shift with NBFonly,
370 * but we do so anyhow for consistency of the returned coordinates.
374 shift_self(graph,box,x);
377 inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
381 inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
384 /* Check whether we need to do bondeds or correct for exclusions */
386 ((flags & GMX_FORCE_BONDED)
387 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
389 /* Since all atoms are in the rectangular or triclinic unit-cell,
390 * only single box vector shifts (2 in x) are required.
392 set_pbc_dd(&pbc,fr->ePBC,cr->dd,TRUE,box);
396 if (flags & GMX_FORCE_BONDED)
398 wallcycle_sub_start(wcycle, ewcsBONDED);
399 calc_bonds(fplog,cr->ms,
400 idef,x,hist,f,fr,&pbc,graph,enerd,nrnb,lambda,md,fcd,
401 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
403 fr->bSepDVDL && do_per_step(step,ir->nstlog),step);
405 /* Check if we have to determine energy differences
406 * at foreign lambda's.
408 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
409 idef->ilsort != ilsortNO_FE)
411 if (idef->ilsort != ilsortFE_SORTED)
413 gmx_incons("The bonded interactions are not sorted for free energy");
415 for(i=0; i<enerd->n_lambda; i++)
417 reset_foreign_enerdata(enerd);
418 for (j=0;j<efptNR;j++)
420 lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
422 calc_bonds_lambda(fplog,idef,x,fr,&pbc,graph,&(enerd->foreign_grpp),enerd->foreign_term,nrnb,lam_i,md,
423 fcd,DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
424 sum_epot(&ir->opts,&(enerd->foreign_grpp),enerd->foreign_term);
425 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
430 wallcycle_sub_stop(wcycle, ewcsBONDED);
436 if (EEL_FULL(fr->eeltype))
438 bSB = (ir->nwall == 2);
442 svmul(ir->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
443 box_size[ZZ] *= ir->wall_ewald_zfac;
446 clear_mat(fr->vir_el_recip);
453 /* With the Verlet scheme exclusion forces are calculated
454 * in the non-bonded kernel.
456 /* The TPI molecule does not have exclusions with the rest
457 * of the system and no intra-molecular PME grid contributions
458 * will be calculated in gmx_pme_calc_energy.
460 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
461 ir->ewald_geometry != eewg3D ||
462 ir->epsilon_surface != 0)
466 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
470 gmx_fatal(FARGS,"TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
473 nthreads = gmx_omp_nthreads_get(emntBonded);
474 #pragma omp parallel for num_threads(nthreads) schedule(static)
475 for(t=0; t<nthreads; t++)
483 fnv = fr->f_novirsum;
484 vir = &fr->vir_el_recip;
491 vir = &fr->f_t[t].vir;
492 Vcorrt = &fr->f_t[t].Vcorr;
493 dvdlt = &fr->f_t[t].dvdl[efptCOUL];
494 for(i=0; i<fr->natoms_force; i++)
502 ewald_LRcorrection(fplog,
503 fr->excl_load[t],fr->excl_load[t+1],
506 md->nChargePerturbed ? md->chargeB : NULL,
507 ir->cutoff_scheme != ecutsVERLET,
508 excl,x,bSB ? boxs : box,mu_tot,
512 lambda[efptCOUL],dvdlt);
516 reduce_thread_forces(fr->natoms_force,fr->f_novirsum,
518 &Vcorr,efptCOUL,&dvdl,
522 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
527 Vcorr += ewald_charge_correction(cr,fr,lambda[efptCOUL],box,
528 &dvdl,fr->vir_el_recip);
531 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.",Vcorr,dvdl);
532 enerd->dvdl_lin[efptCOUL] += dvdl;
543 case eelPMEUSERSWITCH:
545 if (cr->duty & DUTY_PME)
547 assert(fr->n_tpi >= 0);
548 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
550 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
551 if (flags & GMX_FORCE_FORCES)
553 pme_flags |= GMX_PME_CALC_F;
555 if (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY))
557 pme_flags |= GMX_PME_CALC_ENER_VIR;
561 /* We don't calculate f, but we do want the potential */
562 pme_flags |= GMX_PME_CALC_POT;
564 wallcycle_start(wcycle,ewcPMEMESH);
565 status = gmx_pme_do(fr->pmedata,
566 md->start,md->homenr - fr->n_tpi,
568 md->chargeA,md->chargeB,
570 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
571 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
573 fr->vir_el_recip,fr->ewaldcoeff,
574 &Vlr,lambda[efptCOUL],&dvdl,
576 *cycles_pme = wallcycle_stop(wcycle,ewcPMEMESH);
578 /* We should try to do as little computation after
579 * this as possible, because parallel PME synchronizes
580 * the nodes, so we want all load imbalance of the rest
581 * of the force calculation to be before the PME call.
582 * DD load balancing is done on the whole time of
583 * the force call (without PME).
588 /* Determine the PME grid energy of the test molecule
589 * with the PME grid potential of the other charges.
591 gmx_pme_calc_energy(fr->pmedata,fr->n_tpi,
592 x + md->homenr - fr->n_tpi,
593 md->chargeA + md->homenr - fr->n_tpi,
596 PRINT_SEPDVDL("PME mesh",Vlr,dvdl);
600 Vlr = do_ewald(fplog,FALSE,ir,x,fr->f_novirsum,
601 md->chargeA,md->chargeB,
602 box_size,cr,md->homenr,
603 fr->vir_el_recip,fr->ewaldcoeff,
604 lambda[efptCOUL],&dvdl,fr->ewald_table);
605 PRINT_SEPDVDL("Ewald long-range",Vlr,dvdl);
608 gmx_fatal(FARGS,"No such electrostatics method implemented %s",
609 eel_names[fr->eeltype]);
613 gmx_fatal(FARGS,"Error %d in long range electrostatics routine %s",
614 status,EELTYPE(fr->eeltype));
616 /* Note that with separate PME nodes we get the real energies later */
617 enerd->dvdl_lin[efptCOUL] += dvdl;
618 enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
621 fprintf(debug,"Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
622 Vlr,Vcorr,enerd->term[F_COUL_RECIP]);
623 pr_rvecs(debug,0,"vir_el_recip after corr",fr->vir_el_recip,DIM);
624 pr_rvecs(debug,0,"fshift after LR Corrections",fr->fshift,SHIFTS);
629 if (EEL_RF(fr->eeltype))
631 /* With the Verlet scheme exclusion forces are calculated
632 * in the non-bonded kernel.
634 if (ir->cutoff_scheme != ecutsVERLET && fr->eeltype != eelRF_NEC)
637 enerd->term[F_RF_EXCL] =
638 RF_excl_correction(fplog,fr,graph,md,excl,x,f,
639 fr->fshift,&pbc,lambda[efptCOUL],&dvdl);
642 enerd->dvdl_lin[efptCOUL] += dvdl;
643 PRINT_SEPDVDL("RF exclusion correction",
644 enerd->term[F_RF_EXCL],dvdl);
652 print_nrnb(debug,nrnb);
660 MPI_Barrier(cr->mpi_comm_mygroup);
663 if (fr->timesteps == 11)
665 fprintf(stderr,"* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
666 cr->nodeid, gmx_step_str(fr->timesteps,buf),
667 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
668 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
676 pr_rvecs(debug,0,"fshift after bondeds",fr->fshift,SHIFTS);
681 void init_enerdata(int ngener,int n_lambda,gmx_enerdata_t *enerd)
685 for(i=0; i<F_NRE; i++)
688 enerd->foreign_term[i] = 0;
692 for(i=0; i<efptNR; i++) {
693 enerd->dvdl_lin[i] = 0;
694 enerd->dvdl_nonlin[i] = 0;
700 fprintf(debug,"Creating %d sized group matrix for energies\n",n2);
702 enerd->grpp.nener = n2;
703 enerd->foreign_grpp.nener = n2;
704 for(i=0; (i<egNR); i++)
706 snew(enerd->grpp.ener[i],n2);
707 snew(enerd->foreign_grpp.ener[i],n2);
712 enerd->n_lambda = 1 + n_lambda;
713 snew(enerd->enerpart_lambda,enerd->n_lambda);
721 void destroy_enerdata(gmx_enerdata_t *enerd)
725 for(i=0; (i<egNR); i++)
727 sfree(enerd->grpp.ener[i]);
730 for(i=0; (i<egNR); i++)
732 sfree(enerd->foreign_grpp.ener[i]);
737 sfree(enerd->enerpart_lambda);
741 static real sum_v(int n,real v[])
753 void sum_epot(t_grpopts *opts, gmx_grppairener_t *grpp, real *epot)
757 /* Accumulate energies */
758 epot[F_COUL_SR] = sum_v(grpp->nener,grpp->ener[egCOULSR]);
759 epot[F_LJ] = sum_v(grpp->nener,grpp->ener[egLJSR]);
760 epot[F_LJ14] = sum_v(grpp->nener,grpp->ener[egLJ14]);
761 epot[F_COUL14] = sum_v(grpp->nener,grpp->ener[egCOUL14]);
762 epot[F_COUL_LR] = sum_v(grpp->nener,grpp->ener[egCOULLR]);
763 epot[F_LJ_LR] = sum_v(grpp->nener,grpp->ener[egLJLR]);
764 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
765 epot[F_GBPOL] += sum_v(grpp->nener,grpp->ener[egGB]);
767 /* lattice part of LR doesnt belong to any group
768 * and has been added earlier
770 epot[F_BHAM] = sum_v(grpp->nener,grpp->ener[egBHAMSR]);
771 epot[F_BHAM_LR] = sum_v(grpp->nener,grpp->ener[egBHAMLR]);
774 for(i=0; (i<F_EPOT); i++)
776 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
778 epot[F_EPOT] += epot[i];
783 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
788 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
789 enerd->term[F_DVDL] = 0.0;
790 for (i=0;i<efptNR;i++)
792 if (fepvals->separate_dvdl[i])
794 /* could this be done more readably/compactly? */
803 index = F_DVDL_BONDED;
805 case (efptRESTRAINT):
806 index = F_DVDL_RESTRAINT;
815 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
818 fprintf(debug,"dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
819 efpt_names[i],i,enerd->term[index],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
824 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
827 fprintf(debug,"dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
828 efpt_names[0],i,enerd->term[F_DVDL],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
833 /* Notes on the foreign lambda free energy difference evaluation:
834 * Adding the potential and ekin terms that depend linearly on lambda
835 * as delta lam * dvdl to the energy differences is exact.
836 * For the constraints this is not exact, but we have no other option
837 * without literally changing the lengths and reevaluating the energies at each step.
838 * (try to remedy this post 4.6 - MRS)
839 * For the non-bonded LR term we assume that the soft-core (if present)
840 * no longer affects the energy beyond the short-range cut-off,
841 * which is a very good approximation (except for exotic settings).
842 * (investigate how to overcome this post 4.6 - MRS)
845 for(i=0; i<fepvals->n_lambda; i++)
846 { /* note we are iterating over fepvals here!
847 For the current lam, dlam = 0 automatically,
848 so we don't need to add anything to the
849 enerd->enerpart_lambda[0] */
851 /* we don't need to worry about dvdl contributions to the current lambda, because
852 it's automatically zero */
854 /* first kinetic energy term */
855 dlam = (fepvals->all_lambda[efptMASS][i] - lambda[efptMASS]);
857 enerd->enerpart_lambda[i+1] += enerd->term[F_DKDL]*dlam;
859 for (j=0;j<efptNR;j++)
861 if (j==efptMASS) {continue;} /* no other mass term to worry about */
863 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
864 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
867 fprintf(debug,"enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
868 fepvals->all_lambda[j][i],efpt_names[j],
869 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
870 dlam,enerd->dvdl_lin[j]);
877 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
881 /* First reset all foreign energy components. Foreign energies always called on
882 neighbor search steps */
883 for(i=0; (i<egNR); i++)
885 for(j=0; (j<enerd->grpp.nener); j++)
887 enerd->foreign_grpp.ener[i][j] = 0.0;
891 /* potential energy components */
892 for(i=0; (i<=F_EPOT); i++)
894 enerd->foreign_term[i] = 0.0;
898 void reset_enerdata(t_grpopts *opts,
899 t_forcerec *fr,gmx_bool bNS,
900 gmx_enerdata_t *enerd,
906 /* First reset all energy components, except for the long range terms
907 * on the master at non neighbor search steps, since the long range
908 * terms have already been summed at the last neighbor search step.
910 bKeepLR = (fr->bTwinRange && !bNS);
911 for(i=0; (i<egNR); i++) {
912 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR))) {
913 for(j=0; (j<enerd->grpp.nener); j++)
914 enerd->grpp.ener[i][j] = 0.0;
917 for (i=0;i<efptNR;i++)
919 enerd->dvdl_lin[i] = 0.0;
920 enerd->dvdl_nonlin[i] = 0.0;
923 /* Normal potential energy components */
924 for(i=0; (i<=F_EPOT); i++) {
925 enerd->term[i] = 0.0;
927 /* Initialize the dVdlambda term with the long range contribution */
928 /* Initialize the dvdl term with the long range contribution */
929 enerd->term[F_DVDL] = 0.0;
930 enerd->term[F_DVDL_COUL] = 0.0;
931 enerd->term[F_DVDL_VDW] = 0.0;
932 enerd->term[F_DVDL_BONDED] = 0.0;
933 enerd->term[F_DVDL_RESTRAINT] = 0.0;
934 enerd->term[F_DKDL] = 0.0;
935 if (enerd->n_lambda > 0)
937 for(i=0; i<enerd->n_lambda; i++)
939 enerd->enerpart_lambda[i] = 0.0;
942 /* reset foreign energy data - separate function since we also call it elsewhere */
943 reset_foreign_enerdata(enerd);