2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
51 #include "nonbonded.h"
66 #include "gmx_omp_nthreads.h"
68 #include "gromacs/timing/wallcycle.h"
79 gmx_bool bDoLongRangeNS)
85 if (!fr->ns.nblist_initialized)
87 init_neighbor_list(fp, fr, md->homenr);
95 nsearch = search_neighbours(fp, fr, box, top, groups, cr, nrnb, md,
96 bFillGrid, bDoLongRangeNS);
99 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)
109 dump_nblist(fp, cr, fr, fr->ns.dump_nl);
113 static void reduce_thread_forces(int n, rvec *f,
114 tensor vir_q, tensor vir_lj,
115 real *Vcorr_q, real *Vcorr_lj,
116 real *dvdl_q, real *dvdl_lj,
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)
123 for (i = 0; i < n; i++)
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_q += f_t[t].Vcorr_q;
133 *Vcorr_lj += f_t[t].Vcorr_lj;
134 *dvdl_q += f_t[t].dvdl[efptCOUL];
135 *dvdl_lj += f_t[t].dvdl[efptVDW];
136 m_add(vir_q, f_t[t].vir_q, vir_q);
137 m_add(vir_lj, f_t[t].vir_lj, vir_lj);
141 void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda)
143 fprintf(fplog, " %-30s V %12.5e dVdl %12.5e\n", s, v, dvdlambda);
146 void do_force_lowlevel(FILE *fplog, gmx_int64_t step,
147 t_forcerec *fr, t_inputrec *ir,
148 t_idef *idef, t_commrec *cr,
149 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
151 rvec x[], history_t *hist,
154 gmx_enerdata_t *enerd,
171 gmx_bool bDoEpot, bSepDVDL, bSB;
177 double clam_i, vlam_i;
178 real dvdl_dum[efptNR], dvdl_nb[efptNR], lam_i[efptNR];
179 real dvdl_q, dvdl_lj;
182 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
185 #define PRINT_SEPDVDL(s, v, dvdlambda) if (bSepDVDL) { gmx_print_sepdvdl(fplog, s, v, dvdlambda); }
187 set_pbc(&pbc, fr->ePBC, box);
189 /* reset free energy components */
190 for (i = 0; i < efptNR; i++)
197 for (i = 0; (i < DIM); i++)
199 box_size[i] = box[i][i];
202 bSepDVDL = (fr->bSepDVDL && do_per_step(step, ir->nstlog));
205 /* do QMMM first if requested */
208 enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr);
213 fprintf(fplog, "Step %s: non-bonded V and dVdl for node %d:\n",
214 gmx_step_str(step, buf), cr->nodeid);
217 /* Call the short range functions all in one go. */
220 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
221 #define TAKETIME FALSE
224 MPI_Barrier(cr->mpi_comm_mygroup);
231 /* foreign lambda component for walls */
232 real dvdl_walls = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
233 enerd->grpp.ener[egLJSR], nrnb);
234 PRINT_SEPDVDL("Walls", 0.0, dvdl_walls);
235 enerd->dvdl_lin[efptVDW] += dvdl_walls;
238 /* If doing GB, reset dvda and calculate the Born radii */
239 if (ir->implicit_solvent)
241 wallcycle_sub_start(wcycle, ewcsNONBONDED);
243 for (i = 0; i < born->nr; i++)
250 calc_gb_rad(cr, fr, ir, top, x, &(fr->gblist), born, md, nrnb);
253 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
257 /* We only do non-bonded calculation with group scheme here, the verlet
258 * calls are done from do_force_cutsVERLET(). */
259 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
262 /* Add short-range interactions */
263 donb_flags |= GMX_NONBONDED_DO_SR;
265 if (flags & GMX_FORCE_FORCES)
267 donb_flags |= GMX_NONBONDED_DO_FORCE;
269 if (flags & GMX_FORCE_ENERGY)
271 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
273 if (flags & GMX_FORCE_DO_LR)
275 donb_flags |= GMX_NONBONDED_DO_LR;
278 wallcycle_sub_start(wcycle, ewcsNONBONDED);
279 do_nonbonded(fr, x, f, f_longrange, md, excl,
281 lambda, dvdl_nb, -1, -1, donb_flags);
283 /* If we do foreign lambda and we have soft-core interactions
284 * we have to recalculate the (non-linear) energies contributions.
286 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
288 for (i = 0; i < enerd->n_lambda; i++)
290 for (j = 0; j < efptNR; j++)
292 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
294 reset_foreign_enerdata(enerd);
295 do_nonbonded(fr, x, f, f_longrange, md, excl,
296 &(enerd->foreign_grpp), nrnb,
297 lam_i, dvdl_dum, -1, -1,
298 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
299 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
300 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
303 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
307 /* If we are doing GB, calculate bonded forces and apply corrections
308 * to the solvation forces */
309 /* MRS: Eventually, many need to include free energy contribution here! */
310 if (ir->implicit_solvent)
312 wallcycle_sub_start(wcycle, ewcsBONDED);
313 calc_gb_forces(cr, md, born, top, x, f, fr, idef,
314 ir->gb_algorithm, ir->sa_algorithm, nrnb, &pbc, graph, enerd);
315 wallcycle_sub_stop(wcycle, ewcsBONDED);
326 if (fepvals->sc_alpha != 0)
328 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
332 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
335 if (fepvals->sc_alpha != 0)
337 /* even though coulomb part is linear, we already added it, beacuse we
338 need to go through the vdw calculation anyway */
340 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
344 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
349 real V_short_range = 0;
350 real dvdl_short_range = 0;
352 for (i = 0; i < enerd->grpp.nener; i++)
356 enerd->grpp.ener[egBHAMSR][i] :
357 enerd->grpp.ener[egLJSR][i])
358 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
360 dvdl_short_range = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
361 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",
370 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
373 /* Shift the coordinates. Must be done before bonded forces and PPPM,
374 * but is also necessary for SHAKE and update, therefore it can NOT
375 * go when no bonded forces have to be evaluated.
378 /* Here sometimes we would not need to shift with NBFonly,
379 * but we do so anyhow for consistency of the returned coordinates.
383 shift_self(graph, box, x);
386 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
390 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
393 /* Check whether we need to do bondeds or correct for exclusions */
395 ((flags & GMX_FORCE_BONDED)
396 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype)))
398 /* Since all atoms are in the rectangular or triclinic unit-cell,
399 * only single box vector shifts (2 in x) are required.
401 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
405 if (flags & GMX_FORCE_BONDED)
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 for (i = 0; i < enerd->n_lambda; i++)
426 reset_foreign_enerdata(enerd);
427 for (j = 0; j < efptNR; j++)
429 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
431 calc_bonds_lambda(fplog, idef, x, fr, &pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
432 fcd, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
433 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
434 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
439 wallcycle_sub_stop(wcycle, ewcsBONDED);
445 if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
447 real Vlr = 0, Vcorr = 0;
448 real dvdl_long_range = 0;
451 bSB = (ir->nwall == 2);
455 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
456 box_size[ZZ] *= ir->wall_ewald_zfac;
460 /* Do long-range electrostatics and/or LJ-PME, including related short-range
464 clear_mat(fr->vir_el_recip);
465 clear_mat(fr->vir_lj_recip);
467 if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
469 real Vlr_q = 0, Vlr_lj = 0, Vcorr_q = 0, Vcorr_lj = 0;
470 real dvdl_long_range_q = 0, dvdl_long_range_lj = 0;
473 if (EEL_EWALD(fr->eeltype) || EVDW_PME(fr->vdwtype))
475 real dvdl_long_range_correction_q = 0;
476 real dvdl_long_range_correction_lj = 0;
477 /* With the Verlet scheme exclusion forces are calculated
478 * in the non-bonded kernel.
480 /* The TPI molecule does not have exclusions with the rest
481 * of the system and no intra-molecular PME grid
482 * contributions will be calculated in
483 * gmx_pme_calc_energy.
485 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
486 ir->ewald_geometry != eewg3D ||
487 ir->epsilon_surface != 0)
491 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
495 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
498 nthreads = gmx_omp_nthreads_get(emntBonded);
499 #pragma omp parallel for num_threads(nthreads) schedule(static)
500 for (t = 0; t < nthreads; t++)
504 tensor *vir_q, *vir_lj;
505 real *Vcorrt_q, *Vcorrt_lj, *dvdlt_q, *dvdlt_lj;
508 fnv = fr->f_novirsum;
509 vir_q = &fr->vir_el_recip;
510 vir_lj = &fr->vir_lj_recip;
512 Vcorrt_lj = &Vcorr_lj;
513 dvdlt_q = &dvdl_long_range_correction_q;
514 dvdlt_lj = &dvdl_long_range_correction_lj;
519 vir_q = &fr->f_t[t].vir_q;
520 vir_lj = &fr->f_t[t].vir_lj;
521 Vcorrt_q = &fr->f_t[t].Vcorr_q;
522 Vcorrt_lj = &fr->f_t[t].Vcorr_lj;
523 dvdlt_q = &fr->f_t[t].dvdl[efptCOUL];
524 dvdlt_lj = &fr->f_t[t].dvdl[efptVDW];
525 for (i = 0; i < fr->natoms_force; i++)
534 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
537 md->nChargePerturbed ? md->chargeB : NULL,
539 md->nChargePerturbed ? md->c6B : NULL,
541 md->nChargePerturbed ? md->sigmaB : NULL,
543 md->nChargePerturbed ? md->sigma3B : NULL,
544 ir->cutoff_scheme != ecutsVERLET,
545 excl, x, bSB ? boxs : box, mu_tot,
548 fnv, *vir_q, *vir_lj,
550 lambda[efptCOUL], lambda[efptVDW],
555 reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
556 fr->vir_el_recip, fr->vir_lj_recip,
558 &dvdl_long_range_correction_q,
559 &dvdl_long_range_correction_lj,
562 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
567 Vcorr_q += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
568 &dvdl_long_range_correction_q,
572 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.", Vcorr_q, dvdl_long_range_correction_q);
573 PRINT_SEPDVDL("Ewald excl. corr. LJ", Vcorr_lj, dvdl_long_range_correction_lj);
574 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_correction_q;
575 enerd->dvdl_lin[efptVDW] += dvdl_long_range_correction_lj;
578 if ((EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype)))
580 if (cr->duty & DUTY_PME)
582 /* Do reciprocal PME for Coulomb and/or LJ. */
583 assert(fr->n_tpi >= 0);
584 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
586 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
587 if (EEL_PME(fr->eeltype))
589 pme_flags |= GMX_PME_DO_COULOMB;
591 if (EVDW_PME(fr->vdwtype))
593 pme_flags |= GMX_PME_DO_LJ;
594 if (fr->ljpme_combination_rule == eljpmeLB)
596 /*Lorentz-Berthelot Comb. Rules in LJ-PME*/
597 pme_flags |= GMX_PME_LJ_LB;
600 if (flags & GMX_FORCE_FORCES)
602 pme_flags |= GMX_PME_CALC_F;
604 if (flags & GMX_FORCE_VIRIAL)
606 pme_flags |= GMX_PME_CALC_ENER_VIR;
610 /* We don't calculate f, but we do want the potential */
611 pme_flags |= GMX_PME_CALC_POT;
613 wallcycle_start(wcycle, ewcPMEMESH);
614 status = gmx_pme_do(fr->pmedata,
615 md->start, md->homenr - fr->n_tpi,
617 md->chargeA, md->chargeB,
619 md->sigmaA, md->sigmaB,
620 bSB ? boxs : box, cr,
621 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
622 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
624 fr->vir_el_recip, fr->ewaldcoeff_q,
625 fr->vir_lj_recip, fr->ewaldcoeff_lj,
627 lambda[efptCOUL], lambda[efptVDW],
628 &dvdl_long_range_q, &dvdl_long_range_lj, pme_flags);
629 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
632 gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
634 /* We should try to do as little computation after
635 * this as possible, because parallel PME synchronizes
636 * the nodes, so we want all load imbalance of the
637 * rest of the force calculation to be before the PME
638 * call. DD load balancing is done on the whole time
639 * of the force call (without PME).
644 if (EVDW_PME(ir->vdwtype))
647 gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
649 /* Determine the PME grid energy of the test molecule
650 * with the PME grid potential of the other charges.
652 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
653 x + md->homenr - fr->n_tpi,
654 md->chargeA + md->homenr - fr->n_tpi,
657 PRINT_SEPDVDL("PME mesh", Vlr_q + Vlr_lj, dvdl_long_range_q+dvdl_long_range_lj);
661 if (!EEL_PME(fr->eeltype) && EEL_EWALD(fr->eeltype))
663 Vlr_q = do_ewald(ir, x, fr->f_novirsum,
664 md->chargeA, md->chargeB,
665 box_size, cr, md->homenr,
666 fr->vir_el_recip, fr->ewaldcoeff_q,
667 lambda[efptCOUL], &dvdl_long_range_q, fr->ewald_table);
668 PRINT_SEPDVDL("Ewald long-range", Vlr_q, dvdl_long_range_q);
670 else if (!EEL_EWALD(fr->eeltype))
672 gmx_fatal(FARGS, "No such electrostatics method implemented %s",
673 eel_names[fr->eeltype]);
675 /* Note that with separate PME nodes we get the real energies later */
676 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_q;
677 enerd->dvdl_lin[efptVDW] += dvdl_long_range_lj;
678 enerd->term[F_COUL_RECIP] = Vlr_q + Vcorr_q;
679 enerd->term[F_LJ_RECIP] = Vlr_lj + Vcorr_lj;
682 fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
683 Vlr_q, Vcorr_q, enerd->term[F_COUL_RECIP]);
684 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
685 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
686 fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
687 Vlr_lj, Vcorr_lj, enerd->term[F_LJ_RECIP]);
688 pr_rvecs(debug, 0, "vir_lj_recip after corr", fr->vir_lj_recip, DIM);
693 /* Is there a reaction-field exclusion correction needed? */
694 if (EEL_RF(fr->eeltype) && eelRF_NEC != fr->eeltype)
696 /* With the Verlet scheme, exclusion forces are calculated
697 * in the non-bonded kernel.
699 if (ir->cutoff_scheme != ecutsVERLET)
701 real dvdl_rf_excl = 0;
702 enerd->term[F_RF_EXCL] =
703 RF_excl_correction(fr, graph, md, excl, x, f,
704 fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
706 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
707 PRINT_SEPDVDL("RF exclusion correction",
708 enerd->term[F_RF_EXCL], dvdl_rf_excl);
717 print_nrnb(debug, nrnb);
725 MPI_Barrier(cr->mpi_comm_mygroup);
728 if (fr->timesteps == 11)
730 fprintf(stderr, "* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
731 cr->nodeid, gmx_step_str(fr->timesteps, buf),
732 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
733 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
741 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
746 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
750 for (i = 0; i < F_NRE; i++)
753 enerd->foreign_term[i] = 0;
757 for (i = 0; i < efptNR; i++)
759 enerd->dvdl_lin[i] = 0;
760 enerd->dvdl_nonlin[i] = 0;
766 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
768 enerd->grpp.nener = n2;
769 enerd->foreign_grpp.nener = n2;
770 for (i = 0; (i < egNR); i++)
772 snew(enerd->grpp.ener[i], n2);
773 snew(enerd->foreign_grpp.ener[i], n2);
778 enerd->n_lambda = 1 + n_lambda;
779 snew(enerd->enerpart_lambda, enerd->n_lambda);
787 void destroy_enerdata(gmx_enerdata_t *enerd)
791 for (i = 0; (i < egNR); i++)
793 sfree(enerd->grpp.ener[i]);
796 for (i = 0; (i < egNR); i++)
798 sfree(enerd->foreign_grpp.ener[i]);
803 sfree(enerd->enerpart_lambda);
807 static real sum_v(int n, real v[])
813 for (i = 0; (i < n); i++)
821 void sum_epot(gmx_grppairener_t *grpp, real *epot)
825 /* Accumulate energies */
826 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
827 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
828 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
829 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
830 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
831 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
832 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
833 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
835 /* lattice part of LR doesnt belong to any group
836 * and has been added earlier
838 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
839 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
842 for (i = 0; (i < F_EPOT); i++)
844 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
846 epot[F_EPOT] += epot[i];
851 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
856 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
857 enerd->term[F_DVDL] = 0.0;
858 for (i = 0; i < efptNR; i++)
860 if (fepvals->separate_dvdl[i])
862 /* could this be done more readably/compactly? */
875 index = F_DVDL_BONDED;
877 case (efptRESTRAINT):
878 index = F_DVDL_RESTRAINT;
884 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
887 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
888 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
893 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
896 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
897 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
902 /* Notes on the foreign lambda free energy difference evaluation:
903 * Adding the potential and ekin terms that depend linearly on lambda
904 * as delta lam * dvdl to the energy differences is exact.
905 * For the constraints this is not exact, but we have no other option
906 * without literally changing the lengths and reevaluating the energies at each step.
907 * (try to remedy this post 4.6 - MRS)
908 * For the non-bonded LR term we assume that the soft-core (if present)
909 * no longer affects the energy beyond the short-range cut-off,
910 * which is a very good approximation (except for exotic settings).
911 * (investigate how to overcome this post 4.6 - MRS)
913 if (fepvals->separate_dvdl[efptBONDED])
915 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
919 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
921 enerd->term[F_DVDL_CONSTR] = 0;
923 for (i = 0; i < fepvals->n_lambda; i++)
925 /* note we are iterating over fepvals here!
926 For the current lam, dlam = 0 automatically,
927 so we don't need to add anything to the
928 enerd->enerpart_lambda[0] */
930 /* we don't need to worry about dvdl_lin contributions to dE at
931 current lambda, because the contributions to the current
932 lambda are automatically zeroed */
934 for (j = 0; j < efptNR; j++)
936 /* Note that this loop is over all dhdl components, not just the separated ones */
937 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
938 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
941 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
942 fepvals->all_lambda[j][i], efpt_names[j],
943 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
944 dlam, enerd->dvdl_lin[j]);
951 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
955 /* First reset all foreign energy components. Foreign energies always called on
956 neighbor search steps */
957 for (i = 0; (i < egNR); i++)
959 for (j = 0; (j < enerd->grpp.nener); j++)
961 enerd->foreign_grpp.ener[i][j] = 0.0;
965 /* potential energy components */
966 for (i = 0; (i <= F_EPOT); i++)
968 enerd->foreign_term[i] = 0.0;
972 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
973 gmx_enerdata_t *enerd,
979 /* First reset all energy components, except for the long range terms
980 * on the master at non neighbor search steps, since the long range
981 * terms have already been summed at the last neighbor search step.
983 bKeepLR = (fr->bTwinRange && !bNS);
984 for (i = 0; (i < egNR); i++)
986 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
988 for (j = 0; (j < enerd->grpp.nener); j++)
990 enerd->grpp.ener[i][j] = 0.0;
994 for (i = 0; i < efptNR; i++)
996 enerd->dvdl_lin[i] = 0.0;
997 enerd->dvdl_nonlin[i] = 0.0;
1000 /* Normal potential energy components */
1001 for (i = 0; (i <= F_EPOT); i++)
1003 enerd->term[i] = 0.0;
1005 /* Initialize the dVdlambda term with the long range contribution */
1006 /* Initialize the dvdl term with the long range contribution */
1007 enerd->term[F_DVDL] = 0.0;
1008 enerd->term[F_DVDL_COUL] = 0.0;
1009 enerd->term[F_DVDL_VDW] = 0.0;
1010 enerd->term[F_DVDL_BONDED] = 0.0;
1011 enerd->term[F_DVDL_RESTRAINT] = 0.0;
1012 enerd->term[F_DKDL] = 0.0;
1013 if (enerd->n_lambda > 0)
1015 for (i = 0; i < enerd->n_lambda; i++)
1017 enerd->enerpart_lambda[i] = 0.0;
1020 /* reset foreign energy data - separate function since we also call it elsewhere */
1021 reset_foreign_enerdata(enerd);