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,2014, 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.
46 #include "nonbonded.h"
58 #include "gmx_omp_nthreads.h"
60 #include "gromacs/legacyheaders/types/commrec.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/pbcutil/ishift.h"
63 #include "gromacs/pbcutil/mshift.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/timing/wallcycle.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/smalloc.h"
78 gmx_bool bDoLongRangeNS)
84 if (!fr->ns.nblist_initialized)
86 init_neighbor_list(fp, fr, md->homenr);
94 nsearch = search_neighbours(fp, fr, box, top, groups, cr, nrnb, md,
95 bFillGrid, bDoLongRangeNS);
98 fprintf(debug, "nsearch = %d\n", nsearch);
101 /* Check whether we have to do dynamic load balancing */
102 /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
103 count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
104 &(top->idef),opts->ngener);
106 if (fr->ns.dump_nl > 0)
108 dump_nblist(fp, cr, fr, fr->ns.dump_nl);
112 static void reduce_thread_forces(int n, rvec *f,
113 tensor vir_q, tensor vir_lj,
114 real *Vcorr_q, real *Vcorr_lj,
115 real *dvdl_q, real *dvdl_lj,
116 int nthreads, f_thread_t *f_t)
120 /* This reduction can run over any number of threads */
121 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntBonded)) private(t) schedule(static)
122 for (i = 0; i < n; i++)
124 for (t = 1; t < nthreads; t++)
126 rvec_inc(f[i], f_t[t].f[i]);
129 for (t = 1; t < nthreads; t++)
131 *Vcorr_q += f_t[t].Vcorr_q;
132 *Vcorr_lj += f_t[t].Vcorr_lj;
133 *dvdl_q += f_t[t].dvdl[efptCOUL];
134 *dvdl_lj += f_t[t].dvdl[efptVDW];
135 m_add(vir_q, f_t[t].vir_q, vir_q);
136 m_add(vir_lj, f_t[t].vir_lj, vir_lj);
140 void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda)
142 fprintf(fplog, " %-30s V %12.5e dVdl %12.5e\n", s, v, dvdlambda);
145 void do_force_lowlevel(FILE *fplog, gmx_int64_t step,
146 t_forcerec *fr, t_inputrec *ir,
147 t_idef *idef, t_commrec *cr,
148 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
150 rvec x[], history_t *hist,
153 gmx_enerdata_t *enerd,
170 gmx_bool bDoEpot, bSepDVDL, bSB;
176 double clam_i, vlam_i;
177 real dvdl_dum[efptNR], dvdl_nb[efptNR], lam_i[efptNR];
178 real dvdl_q, dvdl_lj;
181 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
184 #define PRINT_SEPDVDL(s, v, dvdlambda) if (bSepDVDL) { gmx_print_sepdvdl(fplog, s, v, dvdlambda); }
186 set_pbc(&pbc, fr->ePBC, box);
188 /* reset free energy components */
189 for (i = 0; i < efptNR; i++)
196 for (i = 0; (i < DIM); i++)
198 box_size[i] = box[i][i];
201 bSepDVDL = (fr->bSepDVDL && do_per_step(step, ir->nstlog));
204 /* do QMMM first if requested */
207 enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr);
212 fprintf(fplog, "Step %s: non-bonded V and dVdl for rank %d:\n",
213 gmx_step_str(step, buf), cr->nodeid);
216 /* Call the short range functions all in one go. */
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 real dvdl_walls = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
232 enerd->grpp.ener[egLJSR], nrnb);
233 PRINT_SEPDVDL("Walls", 0.0, dvdl_walls);
234 enerd->dvdl_lin[efptVDW] += dvdl_walls;
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, 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 /* Currently all group scheme kernels always calculate (shift-)forces */
265 if (flags & GMX_FORCE_FORCES)
267 donb_flags |= GMX_NONBONDED_DO_FORCE;
269 if (flags & GMX_FORCE_VIRIAL)
271 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
273 if (flags & GMX_FORCE_ENERGY)
275 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
277 if (flags & GMX_FORCE_DO_LR)
279 donb_flags |= GMX_NONBONDED_DO_LR;
282 wallcycle_sub_start(wcycle, ewcsNONBONDED);
283 do_nonbonded(fr, x, f, f_longrange, md, excl,
285 lambda, dvdl_nb, -1, -1, donb_flags);
287 /* If we do foreign lambda and we have soft-core interactions
288 * we have to recalculate the (non-linear) energies contributions.
290 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
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_foreign_enerdata(enerd);
299 do_nonbonded(fr, x, f, f_longrange, md, excl,
300 &(enerd->foreign_grpp), nrnb,
301 lam_i, dvdl_dum, -1, -1,
302 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
303 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
304 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
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, x, f, fr, idef,
318 ir->gb_algorithm, ir->sa_algorithm, nrnb, &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];
353 real V_short_range = 0;
354 real dvdl_short_range = 0;
356 for (i = 0; i < enerd->grpp.nener; i++)
360 enerd->grpp.ener[egBHAMSR][i] :
361 enerd->grpp.ener[egLJSR][i])
362 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
364 dvdl_short_range = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
365 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",
374 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
377 /* Shift the coordinates. Must be done before bonded forces and PPPM,
378 * but is also necessary for SHAKE and update, therefore it can NOT
379 * go when no bonded forces have to be evaluated.
382 /* Here sometimes we would not need to shift with NBFonly,
383 * but we do so anyhow for consistency of the returned coordinates.
387 shift_self(graph, box, x);
390 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
394 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
397 /* Check whether we need to do bondeds or correct for exclusions */
399 ((flags & GMX_FORCE_BONDED)
400 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype)))
402 /* Since all atoms are in the rectangular or triclinic unit-cell,
403 * only single box vector shifts (2 in x) are required.
405 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
409 if (flags & GMX_FORCE_BONDED)
411 wallcycle_sub_start(wcycle, ewcsBONDED);
412 calc_bonds(fplog, cr->ms,
413 idef, x, hist, f, fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
414 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
416 fr->bSepDVDL && do_per_step(step, ir->nstlog), step);
418 /* Check if we have to determine energy differences
419 * at foreign lambda's.
421 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
422 idef->ilsort != ilsortNO_FE)
424 if (idef->ilsort != ilsortFE_SORTED)
426 gmx_incons("The bonded interactions are not sorted for free energy");
428 for (i = 0; i < enerd->n_lambda; i++)
430 reset_foreign_enerdata(enerd);
431 for (j = 0; j < efptNR; j++)
433 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
435 calc_bonds_lambda(fplog, idef, x, fr, &pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
436 fcd, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
437 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
438 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
443 wallcycle_sub_stop(wcycle, ewcsBONDED);
449 clear_mat(fr->vir_el_recip);
450 clear_mat(fr->vir_lj_recip);
452 /* Do long-range electrostatics and/or LJ-PME, including related short-range
455 if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
457 real Vlr = 0, Vcorr = 0;
458 real dvdl_long_range = 0;
460 real Vlr_q = 0, Vlr_lj = 0, Vcorr_q = 0, Vcorr_lj = 0;
461 real dvdl_long_range_q = 0, dvdl_long_range_lj = 0;
463 bSB = (ir->nwall == 2);
467 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
468 box_size[ZZ] *= ir->wall_ewald_zfac;
471 if (EEL_PME_EWALD(fr->eeltype) || EVDW_PME(fr->vdwtype))
473 real dvdl_long_range_correction_q = 0;
474 real dvdl_long_range_correction_lj = 0;
475 /* With the Verlet scheme exclusion forces are calculated
476 * in the non-bonded kernel.
478 /* The TPI molecule does not have exclusions with the rest
479 * of the system and no intra-molecular PME grid
480 * contributions will be calculated in
481 * gmx_pme_calc_energy.
483 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
484 ir->ewald_geometry != eewg3D ||
485 ir->epsilon_surface != 0)
489 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
493 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
496 nthreads = gmx_omp_nthreads_get(emntBonded);
497 #pragma omp parallel for num_threads(nthreads) schedule(static)
498 for (t = 0; t < nthreads; t++)
502 tensor *vir_q, *vir_lj;
503 real *Vcorrt_q, *Vcorrt_lj, *dvdlt_q, *dvdlt_lj;
506 fnv = fr->f_novirsum;
507 vir_q = &fr->vir_el_recip;
508 vir_lj = &fr->vir_lj_recip;
510 Vcorrt_lj = &Vcorr_lj;
511 dvdlt_q = &dvdl_long_range_correction_q;
512 dvdlt_lj = &dvdl_long_range_correction_lj;
517 vir_q = &fr->f_t[t].vir_q;
518 vir_lj = &fr->f_t[t].vir_lj;
519 Vcorrt_q = &fr->f_t[t].Vcorr_q;
520 Vcorrt_lj = &fr->f_t[t].Vcorr_lj;
521 dvdlt_q = &fr->f_t[t].dvdl[efptCOUL];
522 dvdlt_lj = &fr->f_t[t].dvdl[efptVDW];
523 for (i = 0; i < fr->natoms_force; i++)
533 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
536 md->nChargePerturbed ? md->chargeB : NULL,
538 md->nTypePerturbed ? md->sqrt_c6B : NULL,
540 md->nTypePerturbed ? md->sigmaB : NULL,
542 md->nTypePerturbed ? md->sigma3B : NULL,
543 ir->cutoff_scheme != ecutsVERLET,
544 excl, x, bSB ? boxs : box, mu_tot,
547 fnv, *vir_q, *vir_lj,
549 lambda[efptCOUL], lambda[efptVDW],
554 reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
555 fr->vir_el_recip, fr->vir_lj_recip,
557 &dvdl_long_range_correction_q,
558 &dvdl_long_range_correction_lj,
561 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
564 if (EEL_PME_EWALD(fr->eeltype) && fr->n_tpi == 0)
566 Vcorr_q += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
567 &dvdl_long_range_correction_q,
571 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.", Vcorr_q, dvdl_long_range_correction_q);
572 PRINT_SEPDVDL("Ewald excl. corr. LJ", Vcorr_lj, dvdl_long_range_correction_lj);
573 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_correction_q;
574 enerd->dvdl_lin[efptVDW] += dvdl_long_range_correction_lj;
576 if ((EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype)) && (cr->duty & DUTY_PME))
578 /* Do reciprocal PME for Coulomb and/or LJ. */
579 assert(fr->n_tpi >= 0);
580 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
582 pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
583 if (EEL_PME(fr->eeltype))
585 pme_flags |= GMX_PME_DO_COULOMB;
587 if (EVDW_PME(fr->vdwtype))
589 pme_flags |= GMX_PME_DO_LJ;
591 if (flags & GMX_FORCE_FORCES)
593 pme_flags |= GMX_PME_CALC_F;
595 if (flags & GMX_FORCE_VIRIAL)
597 pme_flags |= GMX_PME_CALC_ENER_VIR;
601 /* We don't calculate f, but we do want the potential */
602 pme_flags |= GMX_PME_CALC_POT;
604 wallcycle_start(wcycle, ewcPMEMESH);
605 status = gmx_pme_do(fr->pmedata,
606 0, md->homenr - fr->n_tpi,
608 md->chargeA, md->chargeB,
609 md->sqrt_c6A, md->sqrt_c6B,
610 md->sigmaA, md->sigmaB,
611 bSB ? boxs : box, cr,
612 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
613 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
615 fr->vir_el_recip, fr->ewaldcoeff_q,
616 fr->vir_lj_recip, fr->ewaldcoeff_lj,
618 lambda[efptCOUL], lambda[efptVDW],
619 &dvdl_long_range_q, &dvdl_long_range_lj, pme_flags);
620 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
623 gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
625 /* We should try to do as little computation after
626 * this as possible, because parallel PME synchronizes
627 * the nodes, so we want all load imbalance of the
628 * rest of the force calculation to be before the PME
629 * call. DD load balancing is done on the whole time
630 * of the force call (without PME).
635 if (EVDW_PME(ir->vdwtype))
638 gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
640 /* Determine the PME grid energy of the test molecule
641 * with the PME grid potential of the other charges.
643 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
644 x + md->homenr - fr->n_tpi,
645 md->chargeA + md->homenr - fr->n_tpi,
648 PRINT_SEPDVDL("PME mesh", Vlr_q + Vlr_lj, dvdl_long_range_q+dvdl_long_range_lj);
652 if (!EEL_PME(fr->eeltype) && EEL_PME_EWALD(fr->eeltype))
654 Vlr_q = do_ewald(ir, x, fr->f_novirsum,
655 md->chargeA, md->chargeB,
656 box_size, cr, md->homenr,
657 fr->vir_el_recip, fr->ewaldcoeff_q,
658 lambda[efptCOUL], &dvdl_long_range_q, fr->ewald_table);
659 PRINT_SEPDVDL("Ewald long-range", Vlr_q, dvdl_long_range_q);
662 /* Note that with separate PME nodes we get the real energies later */
663 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_q;
664 enerd->dvdl_lin[efptVDW] += dvdl_long_range_lj;
665 enerd->term[F_COUL_RECIP] = Vlr_q + Vcorr_q;
666 enerd->term[F_LJ_RECIP] = Vlr_lj + Vcorr_lj;
669 fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
670 Vlr_q, Vcorr_q, enerd->term[F_COUL_RECIP]);
671 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
672 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
673 fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
674 Vlr_lj, Vcorr_lj, enerd->term[F_LJ_RECIP]);
675 pr_rvecs(debug, 0, "vir_lj_recip after corr", fr->vir_lj_recip, DIM);
680 /* Is there a reaction-field exclusion correction needed? */
681 if (EEL_RF(fr->eeltype) && eelRF_NEC != fr->eeltype)
683 /* With the Verlet scheme, exclusion forces are calculated
684 * in the non-bonded kernel.
686 if (ir->cutoff_scheme != ecutsVERLET)
688 real dvdl_rf_excl = 0;
689 enerd->term[F_RF_EXCL] =
690 RF_excl_correction(fr, graph, md, excl, x, f,
691 fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
693 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
694 PRINT_SEPDVDL("RF exclusion correction",
695 enerd->term[F_RF_EXCL], dvdl_rf_excl);
704 print_nrnb(debug, nrnb);
712 MPI_Barrier(cr->mpi_comm_mygroup);
715 if (fr->timesteps == 11)
717 fprintf(stderr, "* PP load balancing info: rank %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
718 cr->nodeid, gmx_step_str(fr->timesteps, buf),
719 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
720 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
728 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
733 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
737 for (i = 0; i < F_NRE; i++)
740 enerd->foreign_term[i] = 0;
744 for (i = 0; i < efptNR; i++)
746 enerd->dvdl_lin[i] = 0;
747 enerd->dvdl_nonlin[i] = 0;
753 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
755 enerd->grpp.nener = n2;
756 enerd->foreign_grpp.nener = n2;
757 for (i = 0; (i < egNR); i++)
759 snew(enerd->grpp.ener[i], n2);
760 snew(enerd->foreign_grpp.ener[i], n2);
765 enerd->n_lambda = 1 + n_lambda;
766 snew(enerd->enerpart_lambda, enerd->n_lambda);
774 void destroy_enerdata(gmx_enerdata_t *enerd)
778 for (i = 0; (i < egNR); i++)
780 sfree(enerd->grpp.ener[i]);
783 for (i = 0; (i < egNR); i++)
785 sfree(enerd->foreign_grpp.ener[i]);
790 sfree(enerd->enerpart_lambda);
794 static real sum_v(int n, real v[])
800 for (i = 0; (i < n); i++)
808 void sum_epot(gmx_grppairener_t *grpp, real *epot)
812 /* Accumulate energies */
813 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
814 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
815 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
816 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
817 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
818 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
819 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
820 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
822 /* lattice part of LR doesnt belong to any group
823 * and has been added earlier
825 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
826 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
829 for (i = 0; (i < F_EPOT); i++)
831 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
833 epot[F_EPOT] += epot[i];
838 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
843 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
844 enerd->term[F_DVDL] = 0.0;
845 for (i = 0; i < efptNR; i++)
847 if (fepvals->separate_dvdl[i])
849 /* could this be done more readably/compactly? */
862 index = F_DVDL_BONDED;
864 case (efptRESTRAINT):
865 index = F_DVDL_RESTRAINT;
871 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
874 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
875 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
880 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
883 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
884 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
889 /* Notes on the foreign lambda free energy difference evaluation:
890 * Adding the potential and ekin terms that depend linearly on lambda
891 * as delta lam * dvdl to the energy differences is exact.
892 * For the constraints this is not exact, but we have no other option
893 * without literally changing the lengths and reevaluating the energies at each step.
894 * (try to remedy this post 4.6 - MRS)
895 * For the non-bonded LR term we assume that the soft-core (if present)
896 * no longer affects the energy beyond the short-range cut-off,
897 * which is a very good approximation (except for exotic settings).
898 * (investigate how to overcome this post 4.6 - MRS)
900 if (fepvals->separate_dvdl[efptBONDED])
902 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
906 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
908 enerd->term[F_DVDL_CONSTR] = 0;
910 for (i = 0; i < fepvals->n_lambda; i++)
912 /* note we are iterating over fepvals here!
913 For the current lam, dlam = 0 automatically,
914 so we don't need to add anything to the
915 enerd->enerpart_lambda[0] */
917 /* we don't need to worry about dvdl_lin contributions to dE at
918 current lambda, because the contributions to the current
919 lambda are automatically zeroed */
921 for (j = 0; j < efptNR; j++)
923 /* Note that this loop is over all dhdl components, not just the separated ones */
924 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
925 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
928 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
929 fepvals->all_lambda[j][i], efpt_names[j],
930 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
931 dlam, enerd->dvdl_lin[j]);
938 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
942 /* First reset all foreign energy components. Foreign energies always called on
943 neighbor search steps */
944 for (i = 0; (i < egNR); i++)
946 for (j = 0; (j < enerd->grpp.nener); j++)
948 enerd->foreign_grpp.ener[i][j] = 0.0;
952 /* potential energy components */
953 for (i = 0; (i <= F_EPOT); i++)
955 enerd->foreign_term[i] = 0.0;
959 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
960 gmx_enerdata_t *enerd,
966 /* First reset all energy components, except for the long range terms
967 * on the master at non neighbor search steps, since the long range
968 * terms have already been summed at the last neighbor search step.
970 bKeepLR = (fr->bTwinRange && !bNS);
971 for (i = 0; (i < egNR); i++)
973 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
975 for (j = 0; (j < enerd->grpp.nener); j++)
977 enerd->grpp.ener[i][j] = 0.0;
981 for (i = 0; i < efptNR; i++)
983 enerd->dvdl_lin[i] = 0.0;
984 enerd->dvdl_nonlin[i] = 0.0;
987 /* Normal potential energy components */
988 for (i = 0; (i <= F_EPOT); i++)
990 enerd->term[i] = 0.0;
992 /* Initialize the dVdlambda term with the long range contribution */
993 /* Initialize the dvdl term with the long range contribution */
994 enerd->term[F_DVDL] = 0.0;
995 enerd->term[F_DVDL_COUL] = 0.0;
996 enerd->term[F_DVDL_VDW] = 0.0;
997 enerd->term[F_DVDL_BONDED] = 0.0;
998 enerd->term[F_DVDL_RESTRAINT] = 0.0;
999 enerd->term[F_DKDL] = 0.0;
1000 if (enerd->n_lambda > 0)
1002 for (i = 0; i < enerd->n_lambda; i++)
1004 enerd->enerpart_lambda[i] = 0.0;
1007 /* reset foreign energy data - separate function since we also call it elsewhere */
1008 reset_foreign_enerdata(enerd);